Gromacs
2018.8
|
#include <cstdio>
#include "gromacs/math/vectypes.h"
#include "gromacs/mdtypes/pull-params.h"
#include "gromacs/pulling/pull_internal.h"
#include "gromacs/utility/basedefinitions.h"
#include "gromacs/utility/real.h"
This file contains datatypes and function declarations necessary for mdrun to interface with the pull code.
Functions | |
bool | pull_coordinate_is_angletype (const t_pull_coord *pcrd) |
Returns if the pull coordinate is an angle. More... | |
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 (struct pull_t *pull, int coord_ind, const struct 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, const t_mdatoms *mdatoms, gmx::ForceWithVirial *forceWithVirial) |
Apply forces of an external potential to a pull coordinate. More... | |
void | clear_pull_forces (struct pull_t *pull) |
Set the all the pull forces to zero. More... | |
real | pull_potential (struct pull_t *pull, t_mdatoms *md, struct t_pbc *pbc, t_commrec *cr, double t, real lambda, 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, t_mdatoms *md, struct t_pbc *pbc, t_commrec *cr, double dt, double t, rvec *x, rvec *xp, 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 (t_commrec *cr, struct pull_t *pull, t_mdatoms *md) |
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, int nfile, const t_filenm fnm[], const gmx_mtop_t *mtop, t_commrec *cr, const gmx_output_env_t *oenv, real lambda, gmx_bool bOutFile, const ContinuationOptions &continuationOptions) |
Allocate, initialize and return a pull work struct. More... | |
void | finish_pull (struct pull_t *pull) |
Close the pull output files. More... | |
void | pull_print_output (struct pull_t *pull, gmx_int64_t step, double time) |
Print the pull output (x and/or f) More... | |
void | pull_calc_coms (t_commrec *cr, struct pull_t *pull, t_mdatoms *md, struct t_pbc *pbc, double t, rvec x[], rvec *xp) |
Calculates centers of mass all pull groups. More... | |
int | pullCheckPbcWithinGroups (const pull_t &pull, const rvec *x, const t_pbc &pbc) |
Checks whether all groups that use a reference atom are within PBC restrictions. More... | |
gmx_bool | pull_have_potential (const struct pull_t *pull) |
Returns if we have pull coordinates with potential pulling. More... | |
gmx_bool | pull_have_constraint (const struct pull_t *pull) |
Returns if we have 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... | |
Variables | |
static constexpr real | c_pullGroupPbcMargin = 0.9 |
Margin for checking PBC distances compared to half the box size in pullCheckPbcWithinGroups() | |
void apply_external_pull_coord_force | ( | struct pull_t * | pull, |
int | coord_index, | ||
double | coord_force, | ||
const t_mdatoms * | mdatoms, | ||
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.
[in,out] | pull | The pull struct. |
[in] | coord_index | The pull coordinate index to set the force for. |
[in] | coord_force | The scalar force for the pull coordinate. |
[in] | mdatoms | Atom properties, only masses are used. |
[in,out] | forceWithVirial | Force and virial buffers. |
void clear_pull_forces | ( | struct pull_t * | pull | ) |
Set the all the pull forces to zero.
pull | The pull group. |
void dd_make_local_pull_groups | ( | t_commrec * | cr, |
struct pull_t * | pull, | ||
t_mdatoms * | md | ||
) |
Make a selection of the home atoms for all pull groups. Should be called at every domain decomposition.
cr | Structure for communication info. |
pull | The pull group. |
md | All atoms. |
void finish_pull | ( | struct pull_t * | pull | ) |
Close the pull output files.
pull | The pull group. |
double get_pull_coord_value | ( | struct pull_t * | pull, |
int | coord_ind, | ||
const struct t_pbc * | pbc | ||
) |
Get the value for pull coord coord_ind.
[in,out] | pull | The pull struct. |
[in] | coord_ind | Number of the pull coordinate. |
[in] | pbc | Information structure about periodicity. |
struct pull_t* init_pull | ( | FILE * | fplog, |
const pull_params_t * | pull_params, | ||
const t_inputrec * | ir, | ||
int | nfile, | ||
const t_filenm | fnm[], | ||
const gmx_mtop_t * | mtop, | ||
t_commrec * | cr, | ||
const gmx_output_env_t * | oenv, | ||
real | lambda, | ||
gmx_bool | bOutFile, | ||
const ContinuationOptions & | continuationOptions | ||
) |
Allocate, initialize and return a pull work struct.
fplog | General output file, normally md.log. |
pull_params | The pull input parameters containing all pull settings. |
ir | The inputrec. |
nfile | Number of files. |
fnm | Standard filename struct. |
mtop | The topology of the whole system. |
cr | Struct for communication info. |
oenv | Output options. |
lambda | FEP lambda. |
bOutFile | Open output files? |
continuationOptions | Options for continuing from checkpoint file |
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.
[in] | pcrd | Pulling data structure |
[in] | pbc | Information on periodic boundary conditions |
void pull_calc_coms | ( | t_commrec * | cr, |
struct pull_t * | pull, | ||
t_mdatoms * | md, | ||
struct t_pbc * | pbc, | ||
double | t, | ||
rvec | x[], | ||
rvec * | xp | ||
) |
Calculates centers of mass all pull groups.
[in] | cr | Struct for communication info. |
[in] | pull | The pull data structure. |
[in] | md | All atoms. |
[in] | pbc | Information struct about periodicity. |
[in] | t | Time, only used for cylinder ref. |
[in] | x | The local positions. |
[in,out] | xp | Updated x, can be NULL. |
void pull_constraint | ( | struct pull_t * | pull, |
t_mdatoms * | md, | ||
struct t_pbc * | pbc, | ||
t_commrec * | cr, | ||
double | dt, | ||
double | t, | ||
rvec * | x, | ||
rvec * | xp, | ||
rvec * | v, | ||
tensor | vir | ||
) |
Constrain the coordinates xp in the directions in x and also constrain v when v != NULL.
[in,out] | pull | The pull data. |
[in] | md | All atoms. |
[in] | pbc | Information struct about periodicity. |
[in] | cr | Struct for communication info. |
[in] | dt | The time step length. |
[in] | t | The time. |
[in] | x | Positions. |
[in,out] | xp | Updated x, can be NULL. |
[in,out] | v | Velocities, which may get a pull correction. |
[in,out] | vir | The 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.
[in] | pcrd | The pull coordinate to get the conversion factor for. |
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.
[in] | pcrd | The pull coordinate to get the conversion factor for. |
bool pull_coordinate_is_angletype | ( | const t_pull_coord * | pcrd | ) |
Returns if the pull coordinate is an angle.
[in] | pcrd | The pull coordinate to query the type for. |
const char* pull_coordinate_units | ( | const t_pull_coord * | pcrd | ) |
Returns the units of the pull coordinate.
[in] | pcrd | The pull coordinate to query the units for. |
gmx_bool pull_have_constraint | ( | const struct pull_t * | pull | ) |
Returns if we have pull coordinates with constraint pulling.
[in] | pull | The pull data structure. |
gmx_bool pull_have_potential | ( | const struct pull_t * | pull | ) |
Returns if we have pull coordinates with potential pulling.
[in] | pull | The pull data structure. |
real pull_potential | ( | struct pull_t * | pull, |
t_mdatoms * | md, | ||
struct t_pbc * | pbc, | ||
t_commrec * | cr, | ||
double | t, | ||
real | lambda, | ||
rvec * | x, | ||
gmx::ForceWithVirial * | force, | ||
real * | dvdlambda | ||
) |
Determine the COM pull forces and add them to f, return the potential.
[in,out] | pull | The pull struct. |
[in] | md | All atoms. |
[in] | pbc | Information struct about periodicity. |
[in] | cr | Struct for communication info. |
[in] | t | Time. |
[in] | lambda | The value of lambda in FEP calculations. |
[in] | x | Positions. |
[in,out] | force | Forces and virial. |
[out] | dvdlambda | Pull contribution to dV/d(lambda). |
void pull_print_output | ( | struct pull_t * | pull, |
gmx_int64_t | step, | ||
double | time | ||
) |
Print the pull output (x and/or f)
pull | The pull data structure. |
step | Time step number. |
time | Time. |
int pullCheckPbcWithinGroups | ( | const pull_t & | pull, |
const rvec * | x, | ||
const t_pbc & | pbc | ||
) |
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.
Should be called without MPI parallelization and after pull_calc_coms() has been called at least once.
[in] | pull | The pull data structure |
[in] | x | The coordinates |
[in] | pbc | Information struct about periodicity |
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.
[in,out] | pull | The pull struct. |
[in] | coord_index | The pull coordinate index to register the external potential for. |
[in] | provider | Provider string, should match the potential-provider pull coordinate mdp option. |