Merge release-5-0 into master
[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 /*! \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 "gromacs/legacyheaders/typedefs.h"
64 #include "../fileio/filenm.h"
65 #include "../timing/wallcycle.h"
66
67 #ifdef GMX_NATIVE_WINDOWS
68 #include <Windows.h>
69 #define NOFLAGS 0
70 #endif
71
72
73 #ifdef __cplusplus
74 extern "C" {
75 #endif
76
77 static const char IMDstr[] = "IMD:";  /**< Tag output from the IMD module with this string. */
78
79
80 /*! \brief Writes out the group of atoms selected for interactive manipulation.
81  *
82  * Called by grompp.
83  * The resulting file has to be read in by VMD if one wants it to connect to mdrun.
84  *
85  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
86  * \param ir      Structure containing MD input parameters, among those
87  *                the IMD data structure.
88  * \param state   The current state of the MD system.
89  * \param sys     The global, complete system topology.
90  * \param nfile   Number of files.
91  * \param fnm     Filename struct.
92  */
93 extern void write_IMDgroup_to_file(gmx_bool bIMD, t_inputrec *ir, t_state *state,
94                                    gmx_mtop_t *sys, int nfile, const t_filenm fnm[]);
95
96
97 /*! \brief Make a selection of the home atoms for the IMD group.
98  *
99  * Should be called at every domain decomposition.
100  * Each node checks which of the atoms from "ind" are local and puts its local
101  * atom numbers into the "ind_local" array. Furthermore, in "xa_ind" it is
102  * stored at which position each local atom belongs in the assembled/collective
103  * array, so that on the master node all positions can be merged into the
104  * assembled array correctly.
105  *
106  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
107  * \param dd      Structure containing domain decomposition data.
108  * \param imd     The IMD group of atoms.
109  */
110 extern void dd_make_local_IMD_atoms(gmx_bool bIMD, gmx_domdec_t *dd, t_IMD *imd);
111
112
113 /*! \brief Initializes (or disables) IMD.
114  *
115  * This function is called before the main MD loop over time steps,
116  * and it must be called prior to any call to dd_partition_system if in parallel.
117  *
118  * \param ir           The inputrec structure containing the MD input parameters
119  *                     including a pointer to the IMD data structure.
120  * \param cr           Information structure for MPI communication.
121  * \param top_global   The topology of the whole system.
122  * \param fplog        General output file, normally md.log.
123  * \param defnstimd    Default IMD update (=communication) frequency.
124  * \param x            The starting positions of the atoms.
125  * \param nfile        Number of files.
126  * \param fnm          Struct containing file names etc.
127  * \param oenv         Output options.
128  * \param imdport      Port to use for IMD connections.
129  * \param Flags        Flags passed over from main, used to determine
130  *                     whether or not we are appending.
131  */
132 extern void init_IMD(t_inputrec *ir, t_commrec *cr, gmx_mtop_t *top_global,
133                      FILE *fplog, int defnstimd, rvec x[],
134                      int nfile, const t_filenm fnm[], output_env_t oenv,
135                      int imdport, unsigned long  Flags);
136
137
138 /*! \brief IMD required in this time step?
139  * Also checks for new IMD connection and syncs the nodes.
140  *
141  * \param bIMD         Only springs into action if bIMD is TRUE. Otherwise returns directly.
142  * \param step         The time step.
143  * \param cr           Information structure for MPI communication.
144  * \param bNS          Is this a neighbor searching step?
145  * \param box          The simulation box.
146  * \param x            The local atomic positions on this node.
147  * \param ir           The inputrec structure containing the MD input parameters
148  *                     including a pointer to the IMD data structure.
149  * \param t            The time.
150  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
151  *
152  * \returns            Whether or not we have to do IMD communication at this step.
153  */
154 extern gmx_bool do_IMD(gmx_bool bIMD, gmx_int64_t step, t_commrec *cr, gmx_bool bNS,
155                        matrix box, rvec x[], t_inputrec *ir, double t,
156                        gmx_wallcycle_t wcycle);
157
158
159 /*! \brief Get the IMD update frequency.
160  *
161  * \param IMDsetup     Opaque pointer to IMD private data.
162  *
163  * \returns            The current IMD update/communication frequency
164  */
165 extern int IMD_get_step(t_gmx_IMD IMDsetup);
166
167
168 /*! \brief Add external forces from a running interactive molecular dynamics session.
169  *
170  * \param bIMD         Returns directly if bIMD is FALSE.
171  * \param imd          The IMD data structure.
172  * \param cr           Information structure for MPI communication.
173  * \param f            The forces.
174  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
175  */
176 extern void IMD_apply_forces(gmx_bool bIMD, t_IMD *imd, t_commrec *cr, rvec *f,
177                              gmx_wallcycle_t wcycle);
178
179
180 /*! \brief Copy energies and convert to float from enerdata to the IMD energy record.
181  *
182  * We do no conversion, so units in client are SI!
183  *
184  * \param bIMD             Only springs into action if bIMD is TRUE. Otherwise returns directly.
185  * \param imd              The IMD data structure.
186  * \param enerd            Contains the GROMACS energies for the different interaction types.
187  * \param step             The time step.
188  * \param bHaveNewEnergies Only copy energies if we have done global summing of them before.
189  *
190  */
191 extern void IMD_fill_energy_record(gmx_bool bIMD, t_IMD *imd, gmx_enerdata_t *enerd,
192                                    gmx_int64_t step, gmx_bool bHaveNewEnergies);
193
194
195 /*! \brief Send positions and energies to the client.
196  *
197  * \param imd              The IMD data structure.
198  */
199 extern void IMD_send_positions(t_IMD *imd);
200
201
202 /*! \brief Calls IMD_prepare_energies() and then IMD_send_positions().
203  *
204  * \param bIMD             Returns directly if bIMD is FALSE.
205  * \param bIMDstep         If true, transfer the positions. Otherwise just update the time step and potentially the energy record.
206  * \param imd              The IMD data structure.
207  * \param enerd            Contains the GROMACS energies for the different interaction types.
208  * \param step             The time step.
209  * \param bHaveNewEnergies Only update the energy record if we have done global summing of the energies.
210  * \param wcycle           Count wallcycles of IMD routines for diagnostic output.
211  *
212  */
213 extern void IMD_prep_energies_send_positions(gmx_bool bIMD, gmx_bool bIMDstep,
214                                              t_IMD *imd, gmx_enerdata_t *enerd,
215                                              gmx_int64_t step, gmx_bool bHaveNewEnergies,
216                                              gmx_wallcycle_t wcycle);
217
218 /*! \brief Finalize IMD and do some cleaning up.
219  *
220  * Currently, IMD finalize closes the force output file.
221  *
222  * \param bIMD         Returns directly if bIMD is FALSE.
223  * \param imd          The IMD data structure.
224  */
225 extern void IMD_finalize(gmx_bool bIMD, t_IMD *imd);
226
227
228 #ifdef __cplusplus
229 }
230 #endif
231
232 #endif