Gromacs
2024.3
|
#include <gromacs/mdlib/update.h>
Contains data for update phase.
Classes | |
class | Impl |
pImpled implementation for Update More... | |
Public Member Functions | |
Update (const t_inputrec &inputRecord, const gmx_ekindata_t &ekind, BoxDeformation *boxDeformation) | |
Constructor. More... | |
~Update () | |
Destructor. | |
PaddedVector< gmx::RVec > * | xp () |
Get the pointer to updated coordinates. More... | |
BoxDeformation * | deform () const |
Getter to local copy of box deformation class. More... | |
void | updateAfterPartition (int numAtoms, gmx::ArrayRef< const unsigned short > cFREEZE, gmx::ArrayRef< const unsigned short > cTC, gmx::ArrayRef< const unsigned short > cAcceleration) |
Sets data that changes only at domain decomposition time. More... | |
void | update_coords (const t_inputrec &inputRecord, int64_t step, int homenr, bool havePartiallyFrozenAtoms, gmx::ArrayRef< const ParticleType > ptype, gmx::ArrayRef< const real > invMass, gmx::ArrayRef< const gmx::RVec > invMassPerDim, t_state *state, const gmx::ArrayRefWithPadding< const gmx::RVec > &f, t_fcdata *fcdata, const gmx_ekindata_t *ekind, const Matrix3x3 &parrinelloRahmanM, int updatePart, const t_commrec *cr, bool haveConstraints) |
Perform numerical integration step. More... | |
void | finish_update (const t_inputrec &inputRecord, bool havePartiallyFrozenAtoms, int homenr, t_state *state, gmx_wallcycle *wcycle, bool haveConstraints) |
Finalize the coordinate update. More... | |
void | update_sd_second_half (const t_inputrec &inputRecord, int64_t step, real *dvdlambda, int homenr, gmx::ArrayRef< const ParticleType > ptype, gmx::ArrayRef< const real > invMass, t_state *state, const t_commrec *cr, t_nrnb *nrnb, gmx_wallcycle *wcycle, gmx::Constraints *constr, bool do_log, bool do_ene) |
Secong part of the SD integrator. More... | |
void | update_for_constraint_virial (const t_inputrec &inputRecord, int homenr, bool havePartiallyFrozenAtoms, gmx::ArrayRef< const real > invmass, gmx::ArrayRef< const gmx::RVec > invMassPerDim, const t_state &state, const gmx::ArrayRefWithPadding< const gmx::RVec > &f, const gmx_ekindata_t &ekind) |
Performs a leap-frog update without updating state so the constrain virial can be computed. | |
void | update_temperature_constants (const t_inputrec &inputRecord, const gmx_ekindata_t &ekind) |
Update pre-computed constants that depend on the reference temperature for coupling. More... | |
const std::vector< bool > & | getAndersenRandomizeGroup () const |
Getter for the list of the randomize groups. More... | |
const std::vector< real > & | getBoltzmanFactor () const |
Getter for the list of the Boltzmann factors. More... | |
Update::Update | ( | const t_inputrec & | inputRecord, |
const gmx_ekindata_t & | ekind, | ||
BoxDeformation * | boxDeformation | ||
) |
Constructor.
[in] | inputRecord | Input record, used to construct SD object. |
[in] | ekind | Kinetic energy data |
[in] | boxDeformation | Periodic box deformation object. |
BoxDeformation * Update::deform | ( | ) | const |
Getter to local copy of box deformation class.
void Update::finish_update | ( | const t_inputrec & | inputRecord, |
bool | havePartiallyFrozenAtoms, | ||
int | homenr, | ||
t_state * | state, | ||
gmx_wallcycle * | wcycle, | ||
bool | haveConstraints | ||
) |
Finalize the coordinate update.
Copy the updated coordinates to the main coordinates buffer for the atoms that are not frozen.
[in] | inputRecord | Input record. |
[in] | havePartiallyFrozenAtoms | Whether atoms are frozen along 1 or 2 (not 3) dimensions? |
[in] | homenr | The number of atoms on this processor. |
[in] | state | System state object. |
[in] | wcycle | Wall-clock cycle counter. |
[in] | haveConstraints | If the system has constraints. |
const std::vector< bool > & Update::getAndersenRandomizeGroup | ( | ) | const |
Getter for the list of the randomize groups.
Needed for Andersen temperature control.
const std::vector< real > & Update::getBoltzmanFactor | ( | ) | const |
Getter for the list of the Boltzmann factors.
Needed for Andersen temperature control.
void Update::update_coords | ( | const t_inputrec & | inputRecord, |
int64_t | step, | ||
int | homenr, | ||
bool | havePartiallyFrozenAtoms, | ||
gmx::ArrayRef< const ParticleType > | ptype, | ||
gmx::ArrayRef< const real > | invMass, | ||
gmx::ArrayRef< const gmx::RVec > | invMassPerDim, | ||
t_state * | state, | ||
const gmx::ArrayRefWithPadding< const gmx::RVec > & | f, | ||
t_fcdata * | fcdata, | ||
const gmx_ekindata_t * | ekind, | ||
const Matrix3x3 & | parrinelloRahmanM, | ||
int | updatePart, | ||
const t_commrec * | cr, | ||
bool | haveConstraints | ||
) |
Perform numerical integration step.
Selects the appropriate integrator, based on the input record and performs a numerical integration step.
[in] | inputRecord | Input record. |
[in] | step | Current timestep. |
[in] | homenr | The number of atoms on this processor. |
[in] | havePartiallyFrozenAtoms | Whether atoms are frozen along 1 or 2 (not 3) dimensions? |
[in] | ptype | The list of particle types. |
[in] | invMass | Inverse atomic mass per atom, 0 for vsites and shells. |
[in] | invMassPerDim | Inverse atomic mass per atom and dimension, 0 for vsites, shells and frozen dimensions |
[in] | state | System state object. |
[in] | f | Buffer with atomic forces for home particles. |
[in] | fcdata | Force calculation data to update distance and orientation restraints. |
[in] | ekind | Kinetic energy data (for temperature coupling, energy groups, etc.). |
[in] | parrinelloRahmanM | Parrinello-Rahman velocity scaling matrix. |
[in] | updatePart | What should be updated, coordinates or velocities. This enum only used in VV integrator. |
[in] | cr | Comunication record (Old comment: these shouldn't be here – need to think about it). |
[in] | haveConstraints | If the system has constraints. |
void Update::update_sd_second_half | ( | const t_inputrec & | inputRecord, |
int64_t | step, | ||
real * | dvdlambda, | ||
int | homenr, | ||
gmx::ArrayRef< const ParticleType > | ptype, | ||
gmx::ArrayRef< const real > | invMass, | ||
t_state * | state, | ||
const t_commrec * | cr, | ||
t_nrnb * | nrnb, | ||
gmx_wallcycle * | wcycle, | ||
gmx::Constraints * | constr, | ||
bool | do_log, | ||
bool | do_ene | ||
) |
Secong part of the SD integrator.
The first part of integration is performed in the update_coords(...) method.
[in] | inputRecord | Input record. |
[in] | step | Current timestep. |
[in] | dvdlambda | Free energy derivative. Contribution to be added to the bonded interactions. |
[in] | homenr | The number of atoms on this processor. |
[in] | ptype | The list of particle types. |
[in] | invMass | Inverse atomic mass per atom, 0 for vsites and shells. |
[in] | state | System state object. |
[in] | cr | Comunication record. |
[in] | nrnb | Cycle counters. |
[in] | wcycle | Wall-clock cycle counter. |
[in] | constr | Constraints object. The constraints are applied on coordinates after update. |
[in] | do_log | If this is logging step. |
[in] | do_ene | If this is an energy evaluation step. |
void Update::update_temperature_constants | ( | const t_inputrec & | inputRecord, |
const gmx_ekindata_t & | ekind | ||
) |
Update pre-computed constants that depend on the reference temperature for coupling.
This could change e.g. in simulated annealing.
[in] | inputRecord | The input record |
[in] | ekind | Kinetic energy data |
void Update::updateAfterPartition | ( | int | numAtoms, |
gmx::ArrayRef< const unsigned short > | cFREEZE, | ||
gmx::ArrayRef< const unsigned short > | cTC, | ||
gmx::ArrayRef< const unsigned short > | cAcceleration | ||
) |
Sets data that changes only at domain decomposition time.
[in] | numAtoms | Updated number of atoms. |
[in] | cFREEZE | Group index for freezing |
[in] | cTC | Group index for center of mass motion removal |
[in] | cAcceleration | Group index for constant acceleration groups |
PaddedVector< RVec > * Update::xp | ( | ) |
Get the pointer to updated coordinates.
Update saves the updated coordinates into separate buffer, so that constraints will have access to both updated and not update coordinates. For that, update owns a separate buffer. See finish_update(...) for details.