Revert "Eliminate mdlib/mdrun.h"
[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,2015,2016,2017,2018, 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 /*! \libinternal
37  * \defgroup module_imd Interactive molecular dynamics (IMD)
38  * \ingroup group_mdrun
39  *
40  * \brief
41  * Allows mdrun to interface with VMD via the interactive molecular dynamics
42  * (IMD) protocol.
43  *
44  * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
45  *
46  */
47
48 /*! \libinternal \file
49  *
50  * \brief
51  * This file contains datatypes and function declarations necessary for mdrun
52  * to interface with VMD via the interactive molecular dynamics protocol.
53  *
54  * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
55  *
56  * \inlibraryapi
57  * \ingroup module_imd
58  */
59
60 #ifndef GMX_IMD_IMD_H
61 #define GMX_IMD_IMD_H
62
63 #include "config.h"
64
65 #include <cstdio>
66
67 #include "gromacs/math/vectypes.h"
68 #include "gromacs/utility/basedefinitions.h"
69
70
71 #if GMX_NATIVE_WINDOWS
72 #include <Windows.h>
73 #define NOFLAGS 0
74 #endif
75
76 struct gmx_domdec_t;
77 struct gmx_enerdata_t;
78 struct gmx_mtop_t;
79 struct gmx_multisim_t;
80 struct gmx_output_env_t;
81 struct gmx_wallcycle;
82 struct MdrunOptions;
83 struct t_commrec;
84 struct t_filenm;
85 struct t_gmx_IMD;
86 struct t_IMD;
87 struct t_inputrec;
88 class t_state;
89
90 static const char IMDstr[] = "IMD:";  /**< Tag output from the IMD module with this string. */
91
92 /*! \brief Writes out the group of atoms selected for interactive manipulation.
93  *
94  * Called by grompp.
95  * The resulting file has to be read in by VMD if one wants it to connect to mdrun.
96  *
97  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
98  * \param ir      Structure containing MD input parameters, among those
99  *                the IMD data structure.
100  * \param state   The current state of the MD system.
101  * \param sys     The global, complete system topology.
102  * \param nfile   Number of files.
103  * \param fnm     Filename struct.
104  */
105 void write_IMDgroup_to_file(gmx_bool bIMD, t_inputrec *ir, const t_state *state,
106                             const gmx_mtop_t *sys, int nfile, const t_filenm fnm[]);
107
108
109 /*! \brief Make a selection of the home atoms for the IMD group.
110  *
111  * Should be called at every domain decomposition.
112  * Each node checks which of the atoms from "ind" are local and puts its local
113  * atom numbers into the "ind_local" array. Furthermore, in "xa_ind" it is
114  * stored at which position each local atom belongs in the assembled/collective
115  * array, so that on the master node all positions can be merged into the
116  * assembled array correctly.
117  *
118  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
119  * \param dd      Structure containing domain decomposition data.
120  * \param imd     The IMD group of atoms.
121  */
122 void dd_make_local_IMD_atoms(gmx_bool bIMD, const gmx_domdec_t *dd, t_IMD *imd);
123
124
125 /*! \brief Initializes (or disables) IMD.
126  *
127  * This function is called before the main MD loop over time steps.
128  *
129  * \param ir           The inputrec structure containing the MD input parameters
130  *                     including a pointer to the IMD data structure.
131  * \param cr           Information structure for MPI communication.
132  * \param ms           Handler for multi-simulations.
133  * \param top_global   The topology of the whole system.
134  * \param fplog        General output file, normally md.log.
135  * \param defnstimd    Default IMD update (=communication) frequency.
136  * \param x            The starting positions of the atoms.
137  * \param nfile        Number of files.
138  * \param fnm          Struct containing file names etc.
139  * \param oenv         Output options.
140  * \param mdrunOptions Options for mdrun.
141  */
142 void init_IMD(t_inputrec *ir, const t_commrec *cr,
143               const gmx_multisim_t *ms,
144               const gmx_mtop_t *top_global,
145               FILE *fplog, int defnstimd, const rvec x[],
146               int nfile, const t_filenm fnm[], const gmx_output_env_t *oenv,
147               const MdrunOptions &mdrunOptions);
148
149
150 /*! \brief IMD required in this time step?
151  * Also checks for new IMD connection and syncs the nodes.
152  *
153  * \param bIMD         Only springs into action if bIMD is TRUE. Otherwise returns directly.
154  * \param step         The time step.
155  * \param cr           Information structure for MPI communication.
156  * \param bNS          Is this a neighbor searching step?
157  * \param box          The simulation box.
158  * \param x            The local atomic positions on this node.
159  * \param ir           The inputrec structure containing the MD input parameters
160  *                     including a pointer to the IMD data structure.
161  * \param t            The time.
162  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
163  *
164  * \returns            Whether or not we have to do IMD communication at this step.
165  */
166 gmx_bool do_IMD(gmx_bool bIMD, int64_t step, const t_commrec *cr,
167                 gmx_bool bNS,
168                 const matrix box, const rvec x[], t_inputrec *ir, double t,
169                 gmx_wallcycle *wcycle);
170
171
172 /*! \brief Get the IMD update frequency.
173  *
174  * \param IMDsetup     Opaque pointer to IMD private data.
175  *
176  * \returns            The current IMD update/communication frequency
177  */
178 int IMD_get_step(t_gmx_IMD *IMDsetup);
179
180
181 /*! \brief Add external forces from a running interactive molecular dynamics session.
182  *
183  * \param bIMD         Returns directly if bIMD is FALSE.
184  * \param imd          The IMD data structure.
185  * \param cr           Information structure for MPI communication.
186  * \param f            The forces.
187  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
188  */
189 void IMD_apply_forces(gmx_bool bIMD, t_IMD *imd,
190                       const t_commrec *cr, rvec *f,
191                       gmx_wallcycle *wcycle);
192
193
194 /*! \brief Copy energies and convert to float from enerdata to the IMD energy record.
195  *
196  * We do no conversion, so units in client are SI!
197  *
198  * \param bIMD             Only springs into action if bIMD is TRUE. Otherwise returns directly.
199  * \param imd              The IMD data structure.
200  * \param enerd            Contains the GROMACS energies for the different interaction types.
201  * \param step             The time step.
202  * \param bHaveNewEnergies Only copy energies if we have done global summing of them before.
203  *
204  */
205 void IMD_fill_energy_record(gmx_bool bIMD, t_IMD *imd, const gmx_enerdata_t *enerd,
206                             int64_t step, gmx_bool bHaveNewEnergies);
207
208
209 /*! \brief Send positions and energies to the client.
210  *
211  * \param imd              The IMD data structure.
212  */
213 void IMD_send_positions(t_IMD *imd);
214
215
216 /*! \brief Calls IMD_prepare_energies() and then IMD_send_positions().
217  *
218  * \param bIMD             Returns directly if bIMD is FALSE.
219  * \param bIMDstep         If true, transfer the positions. Otherwise just update the time step and potentially the energy record.
220  * \param imd              The IMD data structure.
221  * \param enerd            Contains the GROMACS energies for the different interaction types.
222  * \param step             The time step.
223  * \param bHaveNewEnergies Only update the energy record if we have done global summing of the energies.
224  * \param wcycle           Count wallcycles of IMD routines for diagnostic output.
225  *
226  */
227 void IMD_prep_energies_send_positions(gmx_bool bIMD, gmx_bool bIMDstep,
228                                       t_IMD *imd, const gmx_enerdata_t *enerd,
229                                       int64_t step, gmx_bool bHaveNewEnergies,
230                                       gmx_wallcycle *wcycle);
231
232 /*! \brief Finalize IMD and do some cleaning up.
233  *
234  * Currently, IMD finalize closes the force output file.
235  *
236  * \param bIMD         Returns directly if bIMD is FALSE.
237  * \param imd          The IMD data structure.
238  */
239 void IMD_finalize(gmx_bool bIMD, t_IMD *imd);
240
241 #endif