Gromacs  2024.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Classes | Public Member Functions
gmx::Update Class Reference

#include <gromacs/mdlib/update.h>

Description

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...
 

Constructor & Destructor Documentation

Update::Update ( const t_inputrec &  inputRecord,
const gmx_ekindata_t &  ekind,
BoxDeformation *  boxDeformation 
)

Constructor.

Parameters
[in]inputRecordInput record, used to construct SD object.
[in]ekindKinetic energy data
[in]boxDeformationPeriodic box deformation object.

Member Function Documentation

BoxDeformation * Update::deform ( ) const

Getter to local copy of box deformation class.

Returns
handle to 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.

Parameters
[in]inputRecordInput record.
[in]havePartiallyFrozenAtomsWhether atoms are frozen along 1 or 2 (not 3) dimensions?
[in]homenrThe number of atoms on this processor.
[in]stateSystem state object.
[in]wcycleWall-clock cycle counter.
[in]haveConstraintsIf the system has constraints.
const std::vector< bool > & Update::getAndersenRandomizeGroup ( ) const

Getter for the list of the randomize groups.

Needed for Andersen temperature control.

Returns
Reference to the groups from the SD data object.
const std::vector< real > & Update::getBoltzmanFactor ( ) const

Getter for the list of the Boltzmann factors.

Needed for Andersen temperature control.

Returns
Reference to the Boltzmann factors from the SD data object.
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.

Parameters
[in]inputRecordInput record.
[in]stepCurrent timestep.
[in]homenrThe number of atoms on this processor.
[in]havePartiallyFrozenAtomsWhether atoms are frozen along 1 or 2 (not 3) dimensions?
[in]ptypeThe list of particle types.
[in]invMassInverse atomic mass per atom, 0 for vsites and shells.
[in]invMassPerDimInverse atomic mass per atom and dimension, 0 for vsites, shells and frozen dimensions
[in]stateSystem state object.
[in]fBuffer with atomic forces for home particles.
[in]fcdataForce calculation data to update distance and orientation restraints.
[in]ekindKinetic energy data (for temperature coupling, energy groups, etc.).
[in]parrinelloRahmanMParrinello-Rahman velocity scaling matrix.
[in]updatePartWhat should be updated, coordinates or velocities. This enum only used in VV integrator.
[in]crComunication record (Old comment: these shouldn't be here – need to think about it).
[in]haveConstraintsIf 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.

Parameters
[in]inputRecordInput record.
[in]stepCurrent timestep.
[in]dvdlambdaFree energy derivative. Contribution to be added to the bonded interactions.
[in]homenrThe number of atoms on this processor.
[in]ptypeThe list of particle types.
[in]invMassInverse atomic mass per atom, 0 for vsites and shells.
[in]stateSystem state object.
[in]crComunication record.
[in]nrnbCycle counters.
[in]wcycleWall-clock cycle counter.
[in]constrConstraints object. The constraints are applied on coordinates after update.
[in]do_logIf this is logging step.
[in]do_eneIf 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.

Parameters
[in]inputRecordThe input record
[in]ekindKinetic 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.

Parameters
[in]numAtomsUpdated number of atoms.
[in]cFREEZEGroup index for freezing
[in]cTCGroup index for center of mass motion removal
[in]cAccelerationGroup 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.

Returns
The pointer to the intermediate coordinates buffer.

The documentation for this class was generated from the following files: