Gromacs  2022-beta1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Functions | Variables
pull.h File Reference
#include <cstdio>
#include <optional>
#include "gromacs/math/vectypes.h"
#include "gromacs/mdtypes/pull_params.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/basedefinitions.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for pull.h:
+ This graph shows which files directly or indirectly include this file:

Description

This file contains datatypes and function declarations necessary for mdrun to interface with the pull code.

Author
Berk Hess

Classes

class  gmx::ArrayRef< typename >
 STL-like interface to a C array of T (or part of a std container of T). More...
 

Functions

const char * pull_coordinate_units (const t_pull_coord &pcrd)
 Returns the units of the pull coordinate. More...
 
double pull_conversion_factor_userinput2internal (const t_pull_coord &pcrd)
 Returns the conversion factor from the pull coord init/rate unit to internal value unit. More...
 
double pull_conversion_factor_internal2userinput (const t_pull_coord &pcrd)
 Returns the conversion factor from the pull coord internal value unit to the init/rate unit. More...
 
double get_pull_coord_value (pull_t *pull, int coordIndex, const t_pbc &pbc)
 Get the value for pull coord coord_ind. More...
 
void register_external_pull_potential (struct pull_t *pull, int coord_index, const char *provider)
 Registers the provider of an external potential for a coordinate. More...
 
void apply_external_pull_coord_force (struct pull_t *pull, int coord_index, double coord_force, gmx::ArrayRef< const real > masses, gmx::ForceWithVirial *forceWithVirial)
 Apply forces of an external potential to a pull coordinate. More...
 
void clear_pull_forces (pull_t *pull)
 Set the all the pull forces to zero. More...
 
real pull_potential (pull_t *pull, gmx::ArrayRef< const real > masses, const t_pbc &pbc, const t_commrec *cr, double t, real lambda, gmx::ArrayRef< const gmx::RVec > x, gmx::ForceWithVirial *force, real *dvdlambda)
 Determine the COM pull forces and add them to f, return the potential. More...
 
void pull_constraint (struct pull_t *pull, gmx::ArrayRef< const real > masses, const t_pbc &pbc, const t_commrec *cr, double dt, double t, gmx::ArrayRef< gmx::RVec > x, gmx::ArrayRef< gmx::RVec > xp, gmx::ArrayRef< gmx::RVec > v, tensor vir)
 Constrain the coordinates xp in the directions in x and also constrain v when v != NULL. More...
 
void dd_make_local_pull_groups (const t_commrec *cr, pull_t *pull)
 Make a selection of the home atoms for all pull groups. Should be called at every domain decomposition. More...
 
struct pull_t * init_pull (FILE *fplog, const pull_params_t *pull_params, const t_inputrec *ir, const gmx_mtop_t &mtop, const t_commrec *cr, gmx::LocalAtomSetManager *atomSets, real lambda)
 Allocate, initialize and return a pull work struct. More...
 
void finish_pull (struct pull_t *pull)
 Close the pull output files and delete pull. More...
 
void pull_calc_coms (const t_commrec *cr, pull_t *pull, gmx::ArrayRef< const real > masses, const t_pbc &pbc, double t, gmx::ArrayRef< const gmx::RVec > x, gmx::ArrayRef< gmx::RVec > xp)
 Calculates centers of mass all pull groups. More...
 
int pullCheckPbcWithinGroups (const pull_t &pull, gmx::ArrayRef< const gmx::RVec > x, const t_pbc &pbc, real pbcMargin)
 Checks whether all groups that use a reference atom are within PBC restrictions. More...
 
bool pullCheckPbcWithinGroup (const pull_t &pull, gmx::ArrayRef< const gmx::RVec > x, const t_pbc &pbc, int groupNr, real pbcMargin)
 Checks whether a specific group that uses a reference atom is within PBC restrictions. More...
 
bool pull_have_potential (const pull_t &pull)
 Returns if we have pull coordinates with potential pulling. More...
 
bool pull_have_constraint (const pull_t &pull)
 Returns if we have pull coordinates with constraint pulling. More...
 
bool pull_have_constraint (const pull_params_t &pullParameters)
 Returns if inputrec has pull coordinates with constraint pulling. More...
 
real max_pull_distance2 (const pull_coord_work_t &pcrd, const t_pbc &pbc)
 Returns the maxing distance for pulling. More...
 
void updatePrevStepPullCom (pull_t *pull, std::optional< gmx::ArrayRef< double >> comPreviousStep)
 Sets the previous step COM in pull to the current COM, and optionally updates it in the provided ArrayRef. More...
 
std::vector< double > prevStepPullCom (const pull_t *pull)
 Returns a copy of the previous step pull COM as flat vector. More...
 
void setPrevStepPullCom (pull_t *pull, gmx::ArrayRef< const double > prevStepPullCom)
 Set the previous step pull COM from a flat vector. More...
 
void preparePrevStepPullCom (const t_inputrec *ir, pull_t *pull_work, gmx::ArrayRef< const real > masses, t_state *state, const t_state *state_global, const t_commrec *cr, bool startingFromCheckpoint)
 Allocates, initializes and communicates the previous step pull COM (if that option is set to true). More...
 
void initPullComFromPrevStep (const t_commrec *cr, pull_t *pull, gmx::ArrayRef< const real > masses, const t_pbc &pbc, gmx::ArrayRef< const gmx::RVec > x)
 Initializes the COM of the previous step (set to initial COM) More...
 
void preparePrevStepPullComNewSimulation (const t_commrec *cr, pull_t *pull_work, gmx::ArrayRef< const real > masses, gmx::ArrayRef< const gmx::RVec > x, const matrix box, PbcType pbcType, std::optional< gmx::ArrayRef< double >> &&comPreviousStep)
 Initializes the previous step pull COM for new simulations (no reading from checkpoint). More...
 

Variables

static constexpr real c_pullGroupPbcMargin = 0.9
 Margin for checking pull group PBC distances compared to half the box size.
 
static constexpr real c_pullGroupSmallGroupThreshold = 0.5
 Threshold (as a factor of half the box size) for accepting pull groups without explicitly set refatom.
 

Function Documentation

void apply_external_pull_coord_force ( struct pull_t *  pull,
int  coord_index,
double  coord_force,
gmx::ArrayRef< const real masses,
gmx::ForceWithVirial forceWithVirial 
)

Apply forces of an external potential to a pull coordinate.

This function applies the external scalar force coord_force to the pull coordinate, distributing it over the atoms in the groups involved in the pull coordinate. The corresponding potential energy value should be added to the pull or the module's potential energy term separately by the module itself. This function should be called after pull_potential has been called and, obviously, before the coordinates are updated uses the forces.

Parameters
[in,out]pullThe pull struct.
[in]coord_indexThe pull coordinate index to set the force for.
[in]coord_forceThe scalar force for the pull coordinate.
[in]massesAtoms masses.
[in,out]forceWithVirialForce and virial buffers.
void clear_pull_forces ( pull_t *  pull)

Set the all the pull forces to zero.

Parameters
pullThe pull group.
void dd_make_local_pull_groups ( const t_commrec *  cr,
pull_t *  pull 
)

Make a selection of the home atoms for all pull groups. Should be called at every domain decomposition.

Parameters
crStructure for communication info.
pullThe pull group.
void finish_pull ( struct pull_t *  pull)

Close the pull output files and delete pull.

Parameters
pullThe pull data structure.
double get_pull_coord_value ( pull_t *  pull,
int  coordIndex,
const t_pbc pbc 
)

Get the value for pull coord coord_ind.

Parameters
[in,out]pullThe pull struct.
[in]coordIndexIndex of the pull coordinate in the list of coordinates
[in]pbcInformation structure about periodicity.
Returns
the value of the pull coordinate.
struct pull_t* init_pull ( FILE *  fplog,
const pull_params_t *  pull_params,
const t_inputrec *  ir,
const gmx_mtop_t &  mtop,
const t_commrec *  cr,
gmx::LocalAtomSetManager atomSets,
real  lambda 
)

Allocate, initialize and return a pull work struct.

Parameters
fplogGeneral output file, normally md.log.
pull_paramsThe pull input parameters containing all pull settings.
irThe inputrec.
mtopThe topology of the whole system.
crStruct for communication info.
atomSetsThe manager that handles the pull atom sets
lambdaFEP lambda.
void initPullComFromPrevStep ( const t_commrec *  cr,
pull_t *  pull,
gmx::ArrayRef< const real masses,
const t_pbc pbc,
gmx::ArrayRef< const gmx::RVec x 
)

Initializes the COM of the previous step (set to initial COM)

Parameters
[in]crStruct for communication info.
[in]pullThe pull data structure.
[in]massesAtoms masses.
[in]pbcInformation struct about periodicity.
[in]xThe local positions.
real max_pull_distance2 ( const pull_coord_work_t pcrd,
const t_pbc pbc 
)

Returns the maxing distance for pulling.

For distance geometries, only dimensions with pcrd->params[dim]=1 are included in the distance calculation. For directional geometries, only dimensions with pcrd->vec[dim]!=0 are included in the distance calculation.

Parameters
[in]pcrdPulling data structure
[in]pbcInformation on periodic boundary conditions
Returns
The maximume distance
void preparePrevStepPullCom ( const t_inputrec *  ir,
pull_t *  pull_work,
gmx::ArrayRef< const real masses,
t_state state,
const t_state state_global,
const t_commrec *  cr,
bool  startingFromCheckpoint 
)

Allocates, initializes and communicates the previous step pull COM (if that option is set to true).

If ir->pull->bSetPbcRefToPrevStepCOM is not true nothing is done.

Parameters
[in]irThe input options/settings of the simulation.
[in]pull_workThe COM pull force calculation data structure
[in]massesAtoms masses.
[in]stateThe local (to this rank) state.
[in]state_globalThe global state.
[in]crStruct for communication info.
[in]startingFromCheckpointIs the simulation starting from a checkpoint?
void preparePrevStepPullComNewSimulation ( const t_commrec *  cr,
pull_t *  pull_work,
gmx::ArrayRef< const real masses,
gmx::ArrayRef< const gmx::RVec x,
const matrix  box,
PbcType  pbcType,
std::optional< gmx::ArrayRef< double >> &&  comPreviousStep 
)

Initializes the previous step pull COM for new simulations (no reading from checkpoint).

Parameters
[in]crStruct for communication info.
[in]pull_workThe COM pull force calculation data structure.
[in]massesAtoms masses.
[in]xThe local positions.
[in]boxThe current box matrix.
[in]pbcTypeThe type of periodic boundary conditions.
[in]comPreviousStepThe COM of the previous step of each pull group.
std::vector<double> prevStepPullCom ( const pull_t *  pull)

Returns a copy of the previous step pull COM as flat vector.

Used for modular simulator checkpointing. Allows to keep the implementation details of pull_t hidden from its users.

Parameters
[in]pullThe COM pull force calculation data structure
Returns
A copy of the previous step COM
void pull_calc_coms ( const t_commrec *  cr,
pull_t *  pull,
gmx::ArrayRef< const real masses,
const t_pbc pbc,
double  t,
gmx::ArrayRef< const gmx::RVec x,
gmx::ArrayRef< gmx::RVec xp 
)

Calculates centers of mass all pull groups.

Parameters
[in]crStruct for communication info.
[in]pullThe pull data structure.
[in]massesAtoms masses.
[in]pbcInformation struct about periodicity.
[in]tTime, only used for cylinder ref.
[in]xThe local positions.
[in,out]xpUpdated x, can be NULL.
void pull_constraint ( struct pull_t *  pull,
gmx::ArrayRef< const real masses,
const t_pbc pbc,
const t_commrec *  cr,
double  dt,
double  t,
gmx::ArrayRef< gmx::RVec x,
gmx::ArrayRef< gmx::RVec xp,
gmx::ArrayRef< gmx::RVec v,
tensor  vir 
)

Constrain the coordinates xp in the directions in x and also constrain v when v != NULL.

Parameters
[in,out]pullThe pull data.
[in]massesAtoms masses.
[in]pbcInformation struct about periodicity.
[in]crStruct for communication info.
[in]dtThe time step length.
[in]tThe time.
[in]xPositions.
[in,out]xpUpdated x, can be NULL.
[in,out]vVelocities, which may get a pull correction.
[in,out]virThe virial, which, if != NULL, gets a pull correction.
double pull_conversion_factor_internal2userinput ( const t_pull_coord &  pcrd)

Returns the conversion factor from the pull coord internal value unit to the init/rate unit.

Parameters
[in]pcrdThe pull coordinate to get the conversion factor for.
Returns
the conversion factor.
double pull_conversion_factor_userinput2internal ( const t_pull_coord &  pcrd)

Returns the conversion factor from the pull coord init/rate unit to internal value unit.

Parameters
[in]pcrdThe pull coordinate to get the conversion factor for.
Returns
the conversion factor.
const char* pull_coordinate_units ( const t_pull_coord &  pcrd)

Returns the units of the pull coordinate.

Parameters
[in]pcrdThe pull coordinate to query the units for.
Returns
a string with the units of the coordinate.
bool pull_have_constraint ( const pull_t &  pull)

Returns if we have pull coordinates with constraint pulling.

Parameters
[in]pullThe pull data structure.
bool pull_have_constraint ( const pull_params_t &  pullParameters)

Returns if inputrec has pull coordinates with constraint pulling.

Parameters
[in]pullParametersPulling input parameters from input record.
bool pull_have_potential ( const pull_t &  pull)

Returns if we have pull coordinates with potential pulling.

Parameters
[in]pullThe pull data structure.
real pull_potential ( pull_t *  pull,
gmx::ArrayRef< const real masses,
const t_pbc pbc,
const t_commrec *  cr,
double  t,
real  lambda,
gmx::ArrayRef< const gmx::RVec x,
gmx::ForceWithVirial force,
real dvdlambda 
)

Determine the COM pull forces and add them to f, return the potential.

Parameters
[in,out]pullThe pull struct.
[in]massesAtoms masses.
[in]pbcInformation struct about periodicity.
[in]crStruct for communication info.
[in]tTime.
[in]lambdaThe value of lambda in FEP calculations.
[in]xPositions.
[in,out]forceForces and virial.
[out]dvdlambdaPull contribution to dV/d(lambda).
Returns
The pull potential energy.
bool pullCheckPbcWithinGroup ( const pull_t &  pull,
gmx::ArrayRef< const gmx::RVec x,
const t_pbc pbc,
int  groupNr,
real  pbcMargin 
)

Checks whether a specific group that uses a reference atom is within PBC restrictions.

Groups that use a reference atom for determining PBC should have all their atoms within half the box size from the PBC atom. The box size is used per dimension for rectangular boxes, but can be a combination of dimensions for triclinic boxes, depending on which dimensions are involved in the pull coordinates a group is involved in. A margin is specified to ensure that atoms are not too close to the maximum distance. Only one group is checked.

Should be called without MPI parallelization and after pull_calc_coms() has been called at least once.

Parameters
[in]pullThe pull data structure
[in]xThe coordinates
[in]pbcInformation struct about periodicity
[in]groupNrThe index of the group (in pull.group[]) to check
[in]pbcMarginThe minimum margin (as a fraction) to half the box size
Returns
true if the group obeys PBC otherwise false
int pullCheckPbcWithinGroups ( const pull_t &  pull,
gmx::ArrayRef< const gmx::RVec x,
const t_pbc pbc,
real  pbcMargin 
)

Checks whether all groups that use a reference atom are within PBC restrictions.

Groups that use a reference atom for determining PBC should have all their atoms within half the box size from the PBC atom. The box size is used per dimension for rectangular boxes, but can be a combination of dimensions for triclinic boxes, depending on which dimensions are involved in the pull coordinates a group is involved in. A margin is specified to ensure that atoms are not too close to the maximum distance.

Should be called without MPI parallelization and after pull_calc_coms() has been called at least once.

Parameters
[in]pullThe pull data structure
[in]xThe coordinates
[in]pbcInformation struct about periodicity
[in]pbcMarginThe minimum margin (as a fraction) to half the box size
Returns
-1 when all groups obey PBC or the first group index that fails PBC
void register_external_pull_potential ( struct pull_t *  pull,
int  coord_index,
const char *  provider 
)

Registers the provider of an external potential for a coordinate.

This function is only used for checking the consistency of the pull setup. For each pull coordinate of type external-potential, selected by the user in the mdp file, there has to be a module that provides this potential. The module registers itself as the provider by calling this function. The passed provider string has to match the string that the user passed with the potential-provider pull coordinate mdp option. This function should be called after init_pull has been called and before pull_potential is called for the first time. This function does many consistency checks and when it returns and the first call to do_potential passes, the pull setup is guaranteed to be correct (unless the module doesn't call apply_external_pull_coord_force every step or calls it with incorrect forces). This registering function will exit with a (release) assertion failure when used incorrely or with a fatal error when the user (mdp) input in inconsistent.

Thread-safe for simultaneous registration from multiple threads.

Parameters
[in,out]pullThe pull struct.
[in]coord_indexThe pull coordinate index to register the external potential for.
[in]providerProvider string, should match the potential-provider pull coordinate mdp option.
void setPrevStepPullCom ( pull_t *  pull,
gmx::ArrayRef< const double >  prevStepPullCom 
)

Set the previous step pull COM from a flat vector.

Used to restore modular simulator checkpoints. Allows to keep the implementation details of pull_t hidden from its users.

Parameters
[in]pullThe COM pull force calculation data structure
[in]prevStepPullComThe previous step COM to set
void updatePrevStepPullCom ( pull_t *  pull,
std::optional< gmx::ArrayRef< double >>  comPreviousStep 
)

Sets the previous step COM in pull to the current COM, and optionally updates it in the provided ArrayRef.

Parameters
[in]pullThe COM pull force calculation data structure
[in]comPreviousStepThe COM of the previous step of each pull group