Interactive Molecular Dynamics (IMD)
[alexxy/gromacs.git] / src / gromacs / imd / imd.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 /*! \defgroup module_imd Interactive molecular dynamics (IMD)
37  * \ingroup group_mdrun
38  *
39  * \brief
40  * Allows mdrun to interface with VMD via the interactive molecular dynamics
41  * (IMD) protocol.
42  *
43  * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
44  *
45  */
46
47 /*! \libinternal \file
48  *
49  * \brief
50  * This file contains datatypes and function declarations necessary for mdrun
51  * to interface with VMD via the interactive molecular dynamics protocol.
52  *
53  * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
54  *
55  * \inlibraryapi
56  * \ingroup module_imd
57  */
58
59 #ifndef GMX_IMD_IMD_H
60 #define GMX_IMD_IMD_H
61
62 #include "typedefs.h"
63 #include "../fileio/filenm.h"
64
65 #ifdef GMX_NATIVE_WINDOWS
66 #include <Windows.h>
67 #define NOFLAGS 0
68 #endif
69
70
71 #ifdef __cplusplus
72 extern "C" {
73 #endif
74
75 static const char IMDstr[] = "IMD:";  /**< Tag output from the IMD module with this string. */
76
77
78 /*! \brief Writes out the group of atoms selected for interactive manipulation.
79  *
80  * Called by grompp.
81  * The resulting file has to be read in by VMD if one wants it to connect to mdrun.
82  *
83  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
84  * \param ir      Structure containing MD input parameters, among those
85  *                the IMD data structure.
86  * \param state   The current state of the MD system.
87  * \param sys     The global, complete system topology.
88  * \param nfile   Number of files.
89  * \param fnm     Filename struct.
90  */
91 extern void write_IMDgroup_to_file(gmx_bool bIMD, t_inputrec *ir, t_state *state,
92                                    gmx_mtop_t *sys, int nfile, const t_filenm fnm[]);
93
94
95 /*! \brief Make a selection of the home atoms for the IMD group.
96  *
97  * Should be called at every domain decomposition.
98  * Each node checks which of the atoms from "ind" are local and puts its local
99  * atom numbers into the "ind_local" array. Furthermore, in "xa_ind" it is
100  * stored at which position each local atom belongs in the assembled/collective
101  * array, so that on the master node all positions can be merged into the
102  * assembled array correctly.
103  *
104  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
105  * \param dd      Structure containing domain decomposition data.
106  * \param imd     The IMD group of atoms.
107  */
108 extern void dd_make_local_IMD_atoms(gmx_bool bIMD, gmx_domdec_t *dd, t_IMD *imd);
109
110
111 /*! \brief Initializes (or disables) IMD.
112  *
113  * This function is called before the main MD loop over time steps,
114  * and it must be called prior to any call to dd_partition_system if in parallel.
115  *
116  * \param ir           The inputrec structure containing the MD input parameters
117  *                     including a pointer to the IMD data structure.
118  * \param cr           Information structure for MPI communication.
119  * \param top_global   The topology of the whole system.
120  * \param fplog        General output file, normally md.log.
121  * \param defnstimd    Default IMD update (=communication) frequency.
122  * \param x            The starting positions of the atoms.
123  * \param nfile        Number of files.
124  * \param fnm          Struct containing file names etc.
125  * \param oenv         Output options.
126  * \param imdport      Port to use for IMD connections.
127  * \param Flags        Flags passed over from main, used to determine
128  *                     whether or not we are appending.
129  */
130 extern void init_IMD(t_inputrec *ir, t_commrec *cr, gmx_mtop_t *top_global,
131                      FILE *fplog, int defnstimd, rvec x[],
132                      int nfile, const t_filenm fnm[], output_env_t oenv,
133                      int imdport, unsigned long  Flags);
134
135
136 /*! \brief IMD required in this time step?
137  * Also checks for new IMD connection and syncs the nodes.
138  *
139  * \param bIMD         Only springs into action if bIMD is TRUE. Otherwise returns directly.
140  * \param step         The time step.
141  * \param cr           Information structure for MPI communication.
142  * \param bNS          Is this a neighbor searching step?
143  * \param box          The simulation box.
144  * \param x            The local atomic positions on this node.
145  * \param ir           The inputrec structure containing the MD input parameters
146  *                     including a pointer to the IMD data structure.
147  * \param t            The time.
148  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
149  *
150  * \returns            Whether or not we have to do IMD communication at this step.
151  */
152 extern gmx_bool do_IMD(gmx_bool bIMD, gmx_int64_t step, t_commrec *cr, gmx_bool bNS,
153                        matrix box, rvec x[], t_inputrec *ir, double t,
154                        gmx_wallcycle_t wcycle);
155
156
157 /*! \brief Get the IMD update frequency.
158  *
159  * \param IMDsetup     Opaque pointer to IMD private data.
160  *
161  * \returns            The current IMD update/communication frequency
162  */
163 extern int IMD_get_step(t_gmx_IMD IMDsetup);
164
165
166 /*! \brief Add external forces from a running interactive molecular dynamics session.
167  *
168  * \param bIMD         Returns directly if bIMD is FALSE.
169  * \param imd          The IMD data structure.
170  * \param cr           Information structure for MPI communication.
171  * \param f            The forces.
172  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
173  */
174 extern void IMD_apply_forces(gmx_bool bIMD, t_IMD *imd, t_commrec *cr, rvec *f,
175                              gmx_wallcycle_t wcycle);
176
177
178 /*! \brief Copy energies and convert to float from enerdata to the IMD energy record.
179  *
180  * We do no conversion, so units in client are SI!
181  *
182  * \param bIMD             Only springs into action if bIMD is TRUE. Otherwise returns directly.
183  * \param imd              The IMD data structure.
184  * \param enerd            Contains the GROMACS energies for the different interaction types.
185  * \param step             The time step.
186  * \param bHaveNewEnergies Only copy energies if we have done global summing of them before.
187  *
188  */
189 extern void IMD_fill_energy_record(gmx_bool bIMD, t_IMD *imd, gmx_enerdata_t *enerd,
190                                    gmx_int64_t step, gmx_bool bHaveNewEnergies);
191
192
193 /*! \brief Send positions and energies to the client.
194  *
195  * \param imd              The IMD data structure.
196  */
197 extern void IMD_send_positions(t_IMD *imd);
198
199
200 /*! \brief Calls IMD_prepare_energies() and then IMD_send_positions().
201  *
202  * \param bIMD             Returns directly if bIMD is FALSE.
203  * \param bIMDstep         If true, transfer the positions. Otherwise just update the time step and potentially the energy record.
204  * \param imd              The IMD data structure.
205  * \param enerd            Contains the GROMACS energies for the different interaction types.
206  * \param step             The time step.
207  * \param bHaveNewEnergies Only update the energy record if we have done global summing of the energies.
208  * \param wcycle           Count wallcycles of IMD routines for diagnostic output.
209  *
210  */
211 extern void IMD_prep_energies_send_positions(gmx_bool bIMD, gmx_bool bIMDstep,
212                                              t_IMD *imd, gmx_enerdata_t *enerd,
213                                              gmx_int64_t step, gmx_bool bHaveNewEnergies,
214                                              gmx_wallcycle_t wcycle);
215
216 /*! \brief Finalize IMD and do some cleaning up.
217  *
218  * Currently, IMD finalize closes the force output file.
219  *
220  * \param bIMD         Returns directly if bIMD is FALSE.
221  * \param imd          The IMD data structure.
222  */
223 extern void IMD_finalize(gmx_bool bIMD, t_IMD *imd);
224
225
226 #ifdef __cplusplus
227 }
228 #endif
229
230 #endif