Gromacs  2016.6
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Typedefs | Functions
edsam.h File Reference
#include "gromacs/math/vectypes.h"
#include "gromacs/utility/basedefinitions.h"
+ Include dependency graph for edsam.h:
+ This graph shows which files directly or indirectly include this file:

Description

Declares functions to calculate both essential dynamics constraints as well as flooding potentials and forces.

Authors
Bert de Groot bgroo.nosp@m.t@gw.nosp@m.dg.de, Oliver Lange olive.nosp@m.r.la.nosp@m.nge@t.nosp@m.um.d.nosp@m.e, Carsten Kutzner ckutz.nosp@m.ne@g.nosp@m.wdg.d.nosp@m.e

Typedefs

typedef struct gmx_edsam * gmx_edsam_t
 Abstract type for essential dynamics. More...
 

Functions

void do_edsam (const t_inputrec *ir, gmx_int64_t step, t_commrec *cr, rvec xs[], rvec v[], matrix box, gmx_edsam_t ed)
 Applies essential dynamics constrains as defined in the .edi input file. More...
 
gmx_edsam_t ed_open (int natoms, edsamstate_t *EDstate, int nfile, const t_filenm fnm[], unsigned long Flags, const gmx_output_env_t *oenv, t_commrec *cr)
 Reads in the .edi file containing the essential dynamics and flooding data. More...
 
void init_edsam (const gmx_mtop_t *mtop, const t_inputrec *ir, t_commrec *cr, gmx_edsam_t ed, rvec x[], matrix box, edsamstate_t *EDstate)
 Initializes the essential dynamics and flooding module. More...
 
void dd_make_local_ed_indices (gmx_domdec_t *dd, gmx_edsam_t ed)
 Make a selection of the home atoms for the ED groups. More...
 
void do_flood (t_commrec *cr, const t_inputrec *ir, rvec x[], rvec force[], gmx_edsam_t ed, matrix box, gmx_int64_t step, gmx_bool bNS)
 Evaluate the flooding potential(s) and forces as requested in the .edi input file. More...
 
void done_ed (gmx_edsam_t *ed)
 Clean up. More...
 

Typedef Documentation

typedef struct gmx_edsam* gmx_edsam_t

Abstract type for essential dynamics.

The main type is defined only in edsam.c

Function Documentation

void dd_make_local_ed_indices ( gmx_domdec_t *  dd,
gmx_edsam_t  ed 
)

Make a selection of the home atoms for the ED groups.

Should be called at every domain decomposition.

Parameters
ddDomain decomposition data.
edEssential dynamics and flooding data.
void do_edsam ( const t_inputrec *  ir,
gmx_int64_t  step,
t_commrec *  cr,
rvec  xs[],
rvec  v[],
matrix  box,
gmx_edsam_t  ed 
)

Applies essential dynamics constrains as defined in the .edi input file.

Parameters
irMD input parameter record.
stepNumber of the time step.
crData needed for MPI communication.
xsThe local positions on this processor.
vThe local velocities.
boxThe simulation box.
edThe essential dynamics data.
void do_flood ( t_commrec *  cr,
const t_inputrec *  ir,
rvec  x[],
rvec  force[],
gmx_edsam_t  ed,
matrix  box,
gmx_int64_t  step,
gmx_bool  bNS 
)

Evaluate the flooding potential(s) and forces as requested in the .edi input file.

Parameters
crData needed for MPI communication.
irMD input parameter record.
xPositions on the local processor.
forceForcefield forces to which the flooding forces are added.
edThe essential dynamics data.
boxThe simulation box.
stepNumber of the time step.
bNSAre we in a neighbor searching step?
void done_ed ( gmx_edsam_t ed)

Clean up.

Parameters
edThe essential dynamics data
gmx_edsam_t ed_open ( int  natoms,
edsamstate_t *  EDstate,
int  nfile,
const t_filenm  fnm[],
unsigned long  Flags,
const gmx_output_env_t *  oenv,
t_commrec *  cr 
)

Reads in the .edi file containing the essential dynamics and flooding data.

This function opens the ED input and output files, reads in all datasets it finds in the input file, and cross-checks whether the .edi file information is consistent with the essential dynamics data found in the checkpoint file (if present). gmx make_edi can be used to create an .edi input file.

Parameters
natomsNumber of atoms of the whole MD system.
EDstateEssential dynamics and flooding data stored in the checkpoint file.
nfileNumber of entries (files) in the fnm structure.
fnmThe filenames struct; it contains also the names of the essential dynamics and flooding in + output files.
FlagsFlags passed over from main, used to determine whether we are appending.
oenvNeeded to open the output xvgr file.
crData needed for MPI communication.
Returns
Pointer to the initialized essential dynamics / flooding data.
void init_edsam ( const gmx_mtop_t *  mtop,
const t_inputrec *  ir,
t_commrec *  cr,
gmx_edsam_t  ed,
rvec  x[],
matrix  box,
edsamstate_t *  EDstate 
)

Initializes the essential dynamics and flooding module.

Parameters
mtopMolecular topology.
irMD input parameter record.
crData needed for MPI communication.
edThe essential dynamics data.
xPositions of the whole MD system.
boxThe simulation box.
EDstateED data stored in the checkpoint file.