Gromacs
2018.8
|
#include "config.h"
#include <cstdio>
#include "gromacs/math/vectypes.h"
#include "gromacs/utility/basedefinitions.h"
This file contains datatypes and function declarations necessary for mdrun to interface with VMD via the interactive molecular dynamics protocol.
Functions | |
void | write_IMDgroup_to_file (gmx_bool bIMD, t_inputrec *ir, t_state *state, gmx_mtop_t *sys, int nfile, const t_filenm fnm[]) |
Writes out the group of atoms selected for interactive manipulation. More... | |
void | dd_make_local_IMD_atoms (gmx_bool bIMD, gmx_domdec_t *dd, t_IMD *imd) |
Make a selection of the home atoms for the IMD group. More... | |
void | init_IMD (t_inputrec *ir, t_commrec *cr, gmx_mtop_t *top_global, FILE *fplog, int defnstimd, rvec x[], int nfile, const t_filenm fnm[], const gmx_output_env_t *oenv, const MdrunOptions &mdrunOptions) |
Initializes (or disables) IMD. More... | |
gmx_bool | do_IMD (gmx_bool bIMD, gmx_int64_t step, t_commrec *cr, gmx_bool bNS, matrix box, rvec x[], t_inputrec *ir, double t, gmx_wallcycle *wcycle) |
IMD required in this time step? Also checks for new IMD connection and syncs the nodes. More... | |
int | IMD_get_step (t_gmx_IMD *IMDsetup) |
Get the IMD update frequency. More... | |
void | IMD_apply_forces (gmx_bool bIMD, t_IMD *imd, t_commrec *cr, rvec *f, gmx_wallcycle *wcycle) |
Add external forces from a running interactive molecular dynamics session. More... | |
void | IMD_fill_energy_record (gmx_bool bIMD, t_IMD *imd, gmx_enerdata_t *enerd, gmx_int64_t step, gmx_bool bHaveNewEnergies) |
Copy energies and convert to float from enerdata to the IMD energy record. More... | |
void | IMD_send_positions (t_IMD *imd) |
Send positions and energies to the client. More... | |
void | IMD_prep_energies_send_positions (gmx_bool bIMD, gmx_bool bIMDstep, t_IMD *imd, gmx_enerdata_t *enerd, gmx_int64_t step, gmx_bool bHaveNewEnergies, gmx_wallcycle *wcycle) |
Calls IMD_prepare_energies() and then IMD_send_positions(). More... | |
void | IMD_finalize (gmx_bool bIMD, t_IMD *imd) |
Finalize IMD and do some cleaning up. More... | |
Variables | |
static const char | IMDstr [] = "IMD:" |
Tag output from the IMD module with this string. More... | |
void dd_make_local_IMD_atoms | ( | gmx_bool | bIMD, |
gmx_domdec_t * | dd, | ||
t_IMD * | imd | ||
) |
Make a selection of the home atoms for the IMD group.
Should be called at every domain decomposition. Each node checks which of the atoms from "ind" are local and puts its local atom numbers into the "ind_local" array. Furthermore, in "xa_ind" it is stored at which position each local atom belongs in the assembled/collective array, so that on the master node all positions can be merged into the assembled array correctly.
bIMD | Only springs into action if bIMD is TRUE. Otherwise returns directly. |
dd | Structure containing domain decomposition data. |
imd | The IMD group of atoms. |
gmx_bool do_IMD | ( | gmx_bool | bIMD, |
gmx_int64_t | step, | ||
t_commrec * | cr, | ||
gmx_bool | bNS, | ||
matrix | box, | ||
rvec | x[], | ||
t_inputrec * | ir, | ||
double | t, | ||
gmx_wallcycle * | wcycle | ||
) |
IMD required in this time step? Also checks for new IMD connection and syncs the nodes.
bIMD | Only springs into action if bIMD is TRUE. Otherwise returns directly. |
step | The time step. |
cr | Information structure for MPI communication. |
bNS | Is this a neighbor searching step? |
box | The simulation box. |
x | The local atomic positions on this node. |
ir | The inputrec structure containing the MD input parameters including a pointer to the IMD data structure. |
t | The time. |
wcycle | Count wallcycles of IMD routines for diagnostic output. |
void IMD_apply_forces | ( | gmx_bool | bIMD, |
t_IMD * | imd, | ||
t_commrec * | cr, | ||
rvec * | f, | ||
gmx_wallcycle * | wcycle | ||
) |
Add external forces from a running interactive molecular dynamics session.
bIMD | Returns directly if bIMD is FALSE. |
imd | The IMD data structure. |
cr | Information structure for MPI communication. |
f | The forces. |
wcycle | Count wallcycles of IMD routines for diagnostic output. |
void IMD_fill_energy_record | ( | gmx_bool | bIMD, |
t_IMD * | imd, | ||
gmx_enerdata_t * | enerd, | ||
gmx_int64_t | step, | ||
gmx_bool | bHaveNewEnergies | ||
) |
Copy energies and convert to float from enerdata to the IMD energy record.
We do no conversion, so units in client are SI!
bIMD | Only springs into action if bIMD is TRUE. Otherwise returns directly. |
imd | The IMD data structure. |
enerd | Contains the GROMACS energies for the different interaction types. |
step | The time step. |
bHaveNewEnergies | Only copy energies if we have done global summing of them before. |
void IMD_finalize | ( | gmx_bool | bIMD, |
t_IMD * | imd | ||
) |
Finalize IMD and do some cleaning up.
Currently, IMD finalize closes the force output file.
bIMD | Returns directly if bIMD is FALSE. |
imd | The IMD data structure. |
int IMD_get_step | ( | t_gmx_IMD * | IMDsetup | ) |
Get the IMD update frequency.
IMDsetup | Opaque pointer to IMD private data. |
void IMD_prep_energies_send_positions | ( | gmx_bool | bIMD, |
gmx_bool | bIMDstep, | ||
t_IMD * | imd, | ||
gmx_enerdata_t * | enerd, | ||
gmx_int64_t | step, | ||
gmx_bool | bHaveNewEnergies, | ||
gmx_wallcycle * | wcycle | ||
) |
Calls IMD_prepare_energies() and then IMD_send_positions().
bIMD | Returns directly if bIMD is FALSE. |
bIMDstep | If true, transfer the positions. Otherwise just update the time step and potentially the energy record. |
imd | The IMD data structure. |
enerd | Contains the GROMACS energies for the different interaction types. |
step | The time step. |
bHaveNewEnergies | Only update the energy record if we have done global summing of the energies. |
wcycle | Count wallcycles of IMD routines for diagnostic output. |
void IMD_send_positions | ( | t_IMD * | imd | ) |
Send positions and energies to the client.
imd | The IMD data structure. |
void init_IMD | ( | t_inputrec * | ir, |
t_commrec * | cr, | ||
gmx_mtop_t * | top_global, | ||
FILE * | fplog, | ||
int | defnstimd, | ||
rvec | x[], | ||
int | nfile, | ||
const t_filenm | fnm[], | ||
const gmx_output_env_t * | oenv, | ||
const MdrunOptions & | mdrunOptions | ||
) |
Initializes (or disables) IMD.
This function is called before the main MD loop over time steps, and it must be called prior to any call to dd_partition_system if in parallel.
ir | The inputrec structure containing the MD input parameters including a pointer to the IMD data structure. |
cr | Information structure for MPI communication. |
top_global | The topology of the whole system. |
fplog | General output file, normally md.log. |
defnstimd | Default IMD update (=communication) frequency. |
x | The starting positions of the atoms. |
nfile | Number of files. |
fnm | Struct containing file names etc. |
oenv | Output options. |
mdrunOptions | Options for mdrun. |
void write_IMDgroup_to_file | ( | gmx_bool | bIMD, |
t_inputrec * | ir, | ||
t_state * | state, | ||
gmx_mtop_t * | sys, | ||
int | nfile, | ||
const t_filenm | fnm[] | ||
) |
Writes out the group of atoms selected for interactive manipulation.
Called by grompp. The resulting file has to be read in by VMD if one wants it to connect to mdrun.
bIMD | Only springs into action if bIMD is TRUE. Otherwise returns directly. |
ir | Structure containing MD input parameters, among those the IMD data structure. |
state | The current state of the MD system. |
sys | The global, complete system topology. |
nfile | Number of files. |
fnm | Filename struct. |
|
static |
Tag output from the IMD module with this string.