Gromacs  2018.8
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Macros | Typedefs | Enumerations | Functions | Variables
#include "gmxpre.h"
#include "imd.h"
#include "config.h"
#include <errno.h>
#include <string.h>
#include <unistd.h>
#include "gromacs/commandline/filenm.h"
#include "gromacs/domdec/domdec_struct.h"
#include "gromacs/domdec/ga2la.h"
#include "gromacs/fileio/confio.h"
#include "gromacs/fileio/gmxfio.h"
#include "gromacs/fileio/xvgr.h"
#include "gromacs/gmxlib/network.h"
#include "gromacs/imd/imdsocket.h"
#include "gromacs/math/units.h"
#include "gromacs/math/vec.h"
#include "gromacs/mdlib/broadcaststructs.h"
#include "gromacs/mdlib/groupcoord.h"
#include "gromacs/mdlib/mdrun.h"
#include "gromacs/mdlib/sighandler.h"
#include "gromacs/mdlib/sim_util.h"
#include "gromacs/mdtypes/inputrec.h"
#include "gromacs/mdtypes/md_enums.h"
#include "gromacs/mdtypes/state.h"
#include "gromacs/pbcutil/pbc.h"
#include "gromacs/timing/wallcycle.h"
#include "gromacs/topology/mtop_util.h"
#include "gromacs/topology/topology.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/smalloc.h"
+ Include dependency graph for imd.cpp:

Description

Implements functions of imd.h.

Re-implementation of basic IMD functions from NAMD/VMD from scratch, see imdsocket.h for references to the IMD API.

Author
Martin Hoefling, Carsten Kutzner ckutz.nosp@m.ne@g.nosp@m.wdg.d.nosp@m.e

Classes

struct  IMDEnergyBlock
 IMD (interactive molecular dynamics) energy record. More...
 
struct  IMDHeader
 IMD (interactive molecular dynamics) communication structure. More...
 
struct  t_gmx_IMD
 IMD (interactive molecular dynamics) main data structure. More...
 

Macros

#define IMDLOOPWAIT   1
 How long shall we wait in seconds until we check for a connection again?
 
#define IMDCONNECTWAIT   2
 How long shall we check for the IMD_GO?
 
#define HEADERSIZE   8
 IMD Header Size.
 
#define IMDVERSION   2
 IMD Protocol Version.
 

Typedefs

typedef struct t_gmx_IMD t_gmx_IMD_setup
 IMD (interactive molecular dynamics) main data structure. More...
 
typedef enum IMDType_t IMDMessageType
 Enum for types of IMD messages. More...
 

Enumerations

enum  IMDType_t {
  IMD_DISCONNECT, IMD_ENERGIES, IMD_FCOORDS, IMD_GO,
  IMD_HANDSHAKE, IMD_KILL, IMD_MDCOMM, IMD_PAUSE,
  IMD_TRATE, IMD_IOERROR, IMD_NR
}
 Enum for types of IMD messages. More...
 

Functions

static void fill_header (IMDHeader *header, IMDMessageType type, gmx_int32_t length)
 Fills the header with message and the length argument.
 
static void swap_header (IMDHeader *header)
 Swaps the endianess of the header.
 
static gmx_int32_t imd_read_multiple (IMDSocket *socket, char *datptr, gmx_int32_t toread)
 Reads multiple bytes from socket.
 
static gmx_int32_t imd_write_multiple (IMDSocket *socket, const char *datptr, gmx_int32_t towrite)
 Writes multiple bytes to socket in analogy to imd_read_multiple.
 
static int imd_handshake (IMDSocket *socket)
 Handshake with IMD client.
 
static int imd_send_energies (IMDSocket *socket, const IMDEnergyBlock *energies, char *buffer)
 Send energies using the energy block and the send buffer.
 
static IMDMessageType imd_recv_header (IMDSocket *socket, gmx_int32_t *length)
 Receive IMD header from socket, sets the length and returns the IMD message.
 
static int imd_recv_mdcomm (IMDSocket *socket, gmx_int32_t nforces, gmx_int32_t *forcendx, float *forces)
 Receive force indices and forces. More...
 
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...
 
static int imd_send_rvecs (IMDSocket *socket, int nat, rvec *x, char *buffer)
 Send positions from rvec. More...
 
static t_gmx_IMD_setupimd_create (int imdatoms, int nstimddef, int imdport)
 Initializes the IMD private data.
 
static void imd_prepare_master_socket (t_gmx_IMD_setup *IMDsetup)
 Prepare the socket on the MASTER.
 
static void imd_disconnect (t_gmx_IMD_setup *IMDsetup)
 Disconnect the client.
 
static void imd_fatal (t_gmx_IMD_setup *IMDsetup, const char *msg)
 Prints an error message and disconnects the client. More...
 
static gmx_bool imd_tryconnect (t_gmx_IMD_setup *IMDsetup)
 Check whether we got an incoming connection.
 
static void imd_blockconnect (t_gmx_IMD_setup *IMDsetup)
 Wrap imd_tryconnect in order to make it blocking. More...
 
static void imd_prepare_vmd_Forces (t_gmx_IMD_setup *IMDsetup)
 Make sure that our array holding the forces received via IMD is large enough.
 
static void imd_read_vmd_Forces (t_gmx_IMD_setup *IMDsetup)
 Reads forces received via IMD.
 
static void imd_prepare_MD_Forces (t_gmx_IMD_setup *IMDsetup)
 Prepares the MD force arrays.
 
static void imd_copyto_MD_Forces (t_gmx_IMD_setup *IMDsetup)
 Copy IMD forces to MD forces. More...
 
static gmx_bool bForcesChanged (t_gmx_IMD_setup *IMDsetup)
 Return TRUE if any of the forces or indices changed.
 
static void keep_old_values (t_gmx_IMD_setup *IMDsetup)
 Fill the old_f_ind and old_forces arrays with the new, old values.
 
static gmx_bool rvecs_differ (const rvec v1, const rvec v2)
 Returns TRUE if any component of the two rvecs differs.
 
static void output_imd_forces (t_inputrec *ir, double time)
 Write the applied pull forces to logfile. More...
 
static void imd_sync_nodes (t_inputrec *ir, t_commrec *cr, double t)
 Synchronize the nodes.
 
static void imd_readcommand (t_gmx_IMD_setup *IMDsetup)
 Reads header from the client and decides what to do.
 
static FILE * open_imd_out (const char *fn, t_gmx_IMD_setup *IMDsetup, int nat_total, const gmx_output_env_t *oenv, const ContinuationOptions &continuationOptions)
 Open IMD output file and write header information. More...
 
void IMD_finalize (gmx_bool bIMD, t_IMD *imd)
 Finalize IMD and do some cleaning up. More...
 
static void init_imd_prepare_mols_in_imdgroup (t_gmx_IMD_setup *IMDsetup, gmx_mtop_t *top_global)
 Creates the molecule start-end position array of molecules in the IMD group.
 
static void shift_positions (matrix box, rvec x[], ivec is, int nr)
 Copied and modified from groupcoord.c shift_positions_group().
 
static void imd_remove_molshifts (t_gmx_IMD_setup *IMDsetup, matrix box)
 Removes shifts of molecules diffused outside of the box.
 
static void init_imd_prepare_for_x_assembly (t_commrec *cr, rvec x[], t_gmx_IMD_setup *IMDsetup)
 Initialize arrays used to assemble the positions from the other nodes.
 
static void imd_check_integrator_parallel (t_inputrec *ir, t_commrec *cr)
 Check for non-working integrator / parallel options.
 
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...
 
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...
 
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...
 

Variables

const char * eIMDType_names [IMD_NR+1]
 Names of the IMDType for error messages. More...
 

Typedef Documentation

typedef enum IMDType_t IMDMessageType

Enum for types of IMD messages.

We use the same records as the NAMD/VMD IMD implementation.

typedef struct t_gmx_IMD t_gmx_IMD_setup

IMD (interactive molecular dynamics) main data structure.

Contains private IMD data

Enumeration Type Documentation

enum IMDType_t

Enum for types of IMD messages.

We use the same records as the NAMD/VMD IMD implementation.

Enumerator
IMD_DISCONNECT 

client disconnect

IMD_ENERGIES 

energy data

IMD_FCOORDS 

atomic coordinates

IMD_GO 

start command for the simulation

IMD_HANDSHAKE 

handshake to determine little/big endianness

IMD_KILL 

terminates the simulation

IMD_MDCOMM 

force data

IMD_PAUSE 

pauses the simulation

IMD_TRATE 

sets the IMD transmission and processing rate

IMD_IOERROR 

I/O error.

IMD_NR 

number of entries

Function Documentation

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.

Parameters
bIMDOnly springs into action if bIMD is TRUE. Otherwise returns directly.
ddStructure containing domain decomposition data.
imdThe 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.

Parameters
bIMDOnly springs into action if bIMD is TRUE. Otherwise returns directly.
stepThe time step.
crInformation structure for MPI communication.
bNSIs this a neighbor searching step?
boxThe simulation box.
xThe local atomic positions on this node.
irThe inputrec structure containing the MD input parameters including a pointer to the IMD data structure.
tThe time.
wcycleCount wallcycles of IMD routines for diagnostic output.
Returns
Whether or not we have to do IMD communication at this step.
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.

Parameters
bIMDReturns directly if bIMD is FALSE.
imdThe IMD data structure.
crInformation structure for MPI communication.
fThe forces.
wcycleCount wallcycles of IMD routines for diagnostic output.
static void imd_blockconnect ( t_gmx_IMD_setup IMDsetup)
static

Wrap imd_tryconnect in order to make it blocking.

Used when the simulation should wait for an incoming connection.

static void imd_copyto_MD_Forces ( t_gmx_IMD_setup IMDsetup)
static

Copy IMD forces to MD forces.

Do conversion from Cal->Joule and from Angstrom -> nm and from a pointer array to arrays to 3*N array.

static void imd_fatal ( t_gmx_IMD_setup IMDsetup,
const char *  msg 
)
static

Prints an error message and disconnects the client.

Does not terminate mdrun!

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!

Parameters
bIMDOnly springs into action if bIMD is TRUE. Otherwise returns directly.
imdThe IMD data structure.
enerdContains the GROMACS energies for the different interaction types.
stepThe time step.
bHaveNewEnergiesOnly 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.

Parameters
bIMDReturns directly if bIMD is FALSE.
imdThe IMD data structure.
int IMD_get_step ( t_gmx_IMD IMDsetup)

Get the IMD update frequency.

Parameters
IMDsetupOpaque pointer to IMD private data.
Returns
The current IMD update/communication frequency
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().

Parameters
bIMDReturns directly if bIMD is FALSE.
bIMDstepIf true, transfer the positions. Otherwise just update the time step and potentially the energy record.
imdThe IMD data structure.
enerdContains the GROMACS energies for the different interaction types.
stepThe time step.
bHaveNewEnergiesOnly update the energy record if we have done global summing of the energies.
wcycleCount wallcycles of IMD routines for diagnostic output.
static int imd_recv_mdcomm ( IMDSocket socket,
gmx_int32_t  nforces,
gmx_int32_t *  forcendx,
float *  forces 
)
static

Receive force indices and forces.

The number of forces was previously communicated via the header.

void IMD_send_positions ( t_IMD *  imd)

Send positions and energies to the client.

Parameters
imdThe IMD data structure.
static int imd_send_rvecs ( IMDSocket socket,
int  nat,
rvec *  x,
char *  buffer 
)
static

Send positions from rvec.

We need a separate send buffer and conversion to Angstrom.

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.

Parameters
irThe inputrec structure containing the MD input parameters including a pointer to the IMD data structure.
crInformation structure for MPI communication.
top_globalThe topology of the whole system.
fplogGeneral output file, normally md.log.
defnstimdDefault IMD update (=communication) frequency.
xThe starting positions of the atoms.
nfileNumber of files.
fnmStruct containing file names etc.
oenvOutput options.
mdrunOptionsOptions for mdrun.
static FILE* open_imd_out ( const char *  fn,
t_gmx_IMD_setup IMDsetup,
int  nat_total,
const gmx_output_env_t *  oenv,
const ContinuationOptions continuationOptions 
)
static

Open IMD output file and write header information.

Call on master only.

static void output_imd_forces ( t_inputrec *  ir,
double  time 
)
static

Write the applied pull forces to logfile.

Call on master only!

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.

Parameters
bIMDOnly springs into action if bIMD is TRUE. Otherwise returns directly.
irStructure containing MD input parameters, among those the IMD data structure.
stateThe current state of the MD system.
sysThe global, complete system topology.
nfileNumber of files.
fnmFilename struct.

Variable Documentation

const char* eIMDType_names[IMD_NR+1]
Initial value:
= {
"IMD_DISCONNECT",
"IMD_ENERGIES",
"IMD_FCOORDS",
"IMD_GO",
"IMD_HANDSHAKE",
"IMD_KILL",
"IMD_MDCOMM",
"IMD_PAUSE",
"IMD_TRATE",
"IMD_IOERROR",
nullptr
}

Names of the IMDType for error messages.