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

#include <gromacs/mdlib/constr.h>

Inherited by gmx::Constraints::CreationHelper.

Description

Handles constraints.

Public Member Functions

int numFlexibleConstraints () const
 Returns the total number of flexible constraints in the system.
 
bool havePerturbedConstraints () const
 Returns whether the system contains perturbed constraints.
 
void setConstraints (gmx_localtop_t *top, int numAtoms, int numHomeAtoms, gmx::ArrayRef< const real > masses, gmx::ArrayRef< const real > inverseMasses, bool hasMassPerturbedAtoms, real lambda, gmx::ArrayRef< const unsigned short > cFREEZE)
 Set up all the local constraints for the domain. More...
 
bool apply (bool bLog, bool bEner, int64_t step, int delta_step, real step_scaling, ArrayRefWithPadding< RVec > x, ArrayRefWithPadding< RVec > xprime, ArrayRef< RVec > min_proj, const matrix box, real lambda, real *dvdlambda, ArrayRefWithPadding< RVec > v, bool computeVirial, tensor constraintsVirial, ConstraintVariable econq)
 Applies constraints to coordinates. More...
 
void saveEdsamPointer (gmx_edsam *ed)
 Links the essentialdynamics and constraint code.
 
ArrayRef< const ListOfLists
< int > > 
atom2constraints_moltype () const
 Getter for use by domain decomposition.
 
ArrayRef< const std::vector
< int > > 
atom2settle_moltype () const
 Getter for use by domain decomposition.
 
real rmsd () const
 Return the RMSD of the constraints when available.
 
int numConstraintsTotal ()
 Get the total number of constraints. More...
 

Member Function Documentation

bool gmx::Constraints::apply ( bool  bLog,
bool  bEner,
int64_t  step,
int  delta_step,
real  step_scaling,
ArrayRefWithPadding< RVec x,
ArrayRefWithPadding< RVec xprime,
ArrayRef< RVec min_proj,
const matrix  box,
real  lambda,
real dvdlambda,
ArrayRefWithPadding< RVec v,
bool  computeVirial,
tensor  constraintsVirial,
ConstraintVariable  econq 
)

Applies constraints to coordinates.

When econq=ConstraintVariable::Positions constrains coordinates xprime using the directions in x, min_proj is not used.

When econq=ConstraintVariable::Derivative, calculates the components xprime in the constraint directions and subtracts these components from min_proj. So when min_proj=xprime, the constraint components are projected out.

When econq=ConstraintVariable::Deriv_FlexCon, the same is done as with ConstraintVariable::Derivative, but only the components of the flexible constraints are stored.

delta_step is used for determining the constraint reference lengths when lenA != lenB or will the pull code with a pulling rate. step + delta_step is the step at which the final configuration is meant to be; for update delta_step = 1.

step_scaling can be used to update coordinates based on the time step multiplied by this factor. Thus, normally 1.0 is passed. The SD1 integrator uses 0.5 in one of its calls, to correct positions for half a step of changed velocities.

If v!=NULL also constrain v by adding the constraint corrections / dt.

If computeVirial is true, calculate the constraint virial.

Return whether the application of constraints succeeded without error.

/note x is non-const, because non-local atoms need to be communicated.

int gmx::Constraints::numConstraintsTotal ( )

Get the total number of constraints.

Returns
Total number of constraints, including SETTLE and LINCS/SHAKE constraints.
void gmx::Constraints::setConstraints ( gmx_localtop_t top,
int  numAtoms,
int  numHomeAtoms,
gmx::ArrayRef< const real masses,
gmx::ArrayRef< const real inverseMasses,
bool  hasMassPerturbedAtoms,
real  lambda,
gmx::ArrayRef< const unsigned short >  cFREEZE 
)

Set up all the local constraints for the domain.

Todo:
Make this a callback that is called automatically once a new domain has been made.

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