Start making IMD and swap model IMDModule
[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,2019, 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  * \todo Rename the directory, source and test files to
47  * interactive_md, and prefer type names like
48  * InteractiveMDSession. Avoid ambiguity with IMDModule.
49  */
50
51 /*! \libinternal \file
52  *
53  * \brief This file declares the class that coordinates with VMD via
54  * the Interactive Molecular Dynamics protocol, along with supporting
55  * free functions.
56  *
57  * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
58  *
59  * \inlibraryapi
60  * \ingroup module_imd
61  */
62
63 #ifndef GMX_IMD_IMD_H
64 #define GMX_IMD_IMD_H
65
66 #include <memory>
67
68 #include "gromacs/math/vectypes.h"
69 #include "gromacs/utility/basedefinitions.h"
70 #include "gromacs/utility/classhelpers.h"
71
72 struct gmx_domdec_t;
73 struct gmx_enerdata_t;
74 struct gmx_mtop_t;
75 struct gmx_multisim_t;
76 struct gmx_output_env_t;
77 struct gmx_wallcycle;
78 struct t_commrec;
79 struct t_filenm;
80 struct t_IMD;
81 struct t_inputrec;
82 class t_state;
83
84 namespace gmx
85 {
86 class IMDModule;
87 class ImdSession;
88 class InteractiveMolecularDynamics;
89 class MDLogger;
90 struct MdrunOptions;
91
92 /*! \brief
93  * Creates a module for interactive molecular dynamics.
94  */
95 std::unique_ptr<IMDModule> createInteractiveMolecularDynamicsModule();
96
97 static const char IMDstr[] = "IMD:"; /**< Tag output from the IMD module with this string. */
98
99 /*! \brief Writes out the group of atoms selected for interactive manipulation.
100  *
101  * Called by grompp.
102  * The resulting file has to be read in by VMD if one wants it to connect to mdrun.
103  *
104  * \param bIMD    Only springs into action if bIMD is TRUE. Otherwise returns directly.
105  * \param ir      Structure containing MD input parameters, among those
106  *                the IMD data structure.
107  * \param state   The current state of the MD system.
108  * \param sys     The global, complete system topology.
109  * \param nfile   Number of files.
110  * \param fnm     Filename struct.
111  */
112 void write_IMDgroup_to_file(bool bIMD, t_inputrec *ir, const t_state *state,
113                             const gmx_mtop_t *sys, int nfile, const t_filenm fnm[]);
114
115
116 /*! \brief Makes and returns an initialized IMD session, which may be inactive.
117  *
118  * This function is called before the main MD loop over time steps.
119  *
120  * \param ir           The inputrec structure containing the MD input parameters
121  * \param cr           Information structure for MPI communication.
122  * \param wcycle       Count wallcycles of IMD routines for diagnostic output.
123  * \param enerd        Contains the GROMACS energies for the different interaction types.
124  * \param ms           Handler for multi-simulations.
125  * \param top_global   The topology of the whole system.
126  * \param mdlog        Logger
127  * \param x            The starting positions of the atoms.
128  * \param nfile        Number of files.
129  * \param fnm          Struct containing file names etc.
130  * \param oenv         Output options.
131  * \param mdrunOptions Options for mdrun.
132  */
133 std::unique_ptr<ImdSession>
134 makeImdSession(const t_inputrec *ir,
135                const t_commrec *cr,
136                gmx_wallcycle *wcycle,
137                gmx_enerdata_t *enerd,
138                const gmx_multisim_t *ms,
139                const gmx_mtop_t *top_global,
140                const MDLogger &mdlog,
141                const rvec x[],
142                int nfile, const t_filenm fnm[], const gmx_output_env_t *oenv,
143                const MdrunOptions &mdrunOptions);
144
145 class ImdSession
146 {
147     private:
148         //! Private constructor, to force the use of makeImdSession()
149         ImdSession(const MDLogger &mdlog);
150     public:
151         ~ImdSession();
152
153         /*! \brief Make a selection of the home atoms for the IMD group.
154          *
155          * Should be called at every domain decomposition.
156          * Each node checks which of the atoms from "ind" are local and puts its local
157          * atom numbers into the "ind_local" array. Furthermore, in "xa_ind" it is
158          * stored at which position each local atom belongs in the assembled/collective
159          * array, so that on the master node all positions can be merged into the
160          * assembled array correctly.
161          *
162          * \param dd          Structure containing domain decomposition data.
163          */
164         void dd_make_local_IMD_atoms(const gmx_domdec_t *dd);
165
166         /*! \brief Prepare to do IMD if required in this time step
167          * Also checks for new IMD connection and syncs the nodes.
168          *
169          * \param step         The time step.
170          * \param bNS          Is this a neighbor searching step?
171          * \param box          The simulation box.
172          * \param x            The local atomic positions on this node.
173          * \param t            The time.
174          *
175          * \returns            Whether or not we have to do IMD communication at this step.
176          */
177         bool run(int64_t      step,
178                  bool         bNS,
179                  const matrix box,
180                  const rvec   x[],
181                  double       t);
182
183         /*! \brief Add external forces from a running interactive molecular dynamics session.
184          *
185          * \param f            The forces.
186          */
187         void applyForces(rvec *f);
188
189         /*! \brief Copy energies and convert to float from enerdata to the IMD energy record.
190          *
191          * We do no conversion, so units in client are SI!
192          *
193          * \param step             The time step.
194          * \param bHaveNewEnergies Only copy energies if we have done global summing of them before.
195          */
196         void fillEnergyRecord(int64_t step,
197                               bool    bHaveNewEnergies);
198
199         /*! \brief Send positions and energies to the client. */
200         void sendPositionsAndEnergies();
201
202         /*! \brief Updates the energy record sent to VMD if needed, and sends positions and energies.
203          *
204          * \param bIMDstep         If true, transfer the positions. Otherwise just update the time step and potentially the energy record.
205          * \param step             The time step.
206          * \param bHaveNewEnergies Update the energy record if we have done global summing of the energies.
207          */
208         void updateEnergyRecordAndSendPositionsAndEnergies(bool    bIMDstep,
209                                                            int64_t step,
210                                                            bool    bHaveNewEnergies);
211
212     private:
213         //! Implementation type.
214         class Impl;
215         //! Implementation object.
216         PrivateImplPointer<Impl> impl_;
217
218     public:
219         // Befriend the factory function.
220         friend std::unique_ptr<ImdSession>
221         makeImdSession(const t_inputrec       *ir,
222                        const t_commrec        *cr,
223                        gmx_wallcycle          *wcycle,
224                        gmx_enerdata_t         *enerd,
225                        const gmx_multisim_t   *ms,
226                        const gmx_mtop_t       *top_global,
227                        const MDLogger         &mdlog,
228                        const rvec              x[],
229                        int                     nfile,
230                        const t_filenm          fnm[],
231                        const gmx_output_env_t *oenv,
232                        const MdrunOptions     &mdrunOptions);
233 };
234
235 } // namespace gmx
236
237 #endif