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