Gromacs  2018.8
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions | Variables
pull.h File Reference
#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"
+ Include dependency graph for pull.h:

Description

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

Author
Berk Hess

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()
 

Function Documentation

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.

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]mdatomsAtom properties, only masses are used.
[in,out]forceWithVirialForce and virial buffers.
void clear_pull_forces ( struct pull_t *  pull)

Set the all the pull forces to zero.

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

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

Close the pull output files.

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

Parameters
[in,out]pullThe pull struct.
[in]coord_indNumber of the pull coordinate.
[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,
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.

Parameters
fplogGeneral output file, normally md.log.
pull_paramsThe pull input parameters containing all pull settings.
irThe inputrec.
nfileNumber of files.
fnmStandard filename struct.
mtopThe topology of the whole system.
crStruct for communication info.
oenvOutput options.
lambdaFEP lambda.
bOutFileOpen output files?
continuationOptionsOptions for continuing from checkpoint file
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 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.

Parameters
[in]crStruct for communication info.
[in]pullThe pull data structure.
[in]mdAll atoms.
[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,
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.

Parameters
[in,out]pullThe pull data.
[in]mdAll atoms.
[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.
bool pull_coordinate_is_angletype ( const t_pull_coord *  pcrd)

Returns if the pull coordinate is an angle.

Parameters
[in]pcrdThe pull coordinate to query the type for.
Returns
a boolean telling if the coordinate is of angle type.
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.
gmx_bool pull_have_constraint ( const struct pull_t *  pull)

Returns if we have pull coordinates with constraint pulling.

Parameters
[in]pullThe pull data structure.
gmx_bool pull_have_potential ( const struct pull_t *  pull)

Returns if we have pull coordinates with potential pulling.

Parameters
[in]pullThe 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.

Parameters
[in,out]pullThe pull struct.
[in]mdAll atoms.
[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.
void pull_print_output ( struct pull_t *  pull,
gmx_int64_t  step,
double  time 
)

Print the pull output (x and/or f)

Parameters
pullThe pull data structure.
stepTime step number.
timeTime.
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.

Parameters
[in]pullThe pull data structure
[in]xThe coordinates
[in]pbcInformation struct about periodicity
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.