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