86ea3279a8ae7ebc7e7aed3aae894f93034e5fd4
[alexxy/gromacs.git] / src / gromacs / mdrun / isimulator.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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 /*! \internal
36  * \brief Declares the general simulator interface
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_mdrun
40  */
41 #ifndef GMX_MDRUN_ISIMULATOR_H
42 #define GMX_MDRUN_ISIMULATOR_H
43
44 #include "gromacs/mdlib/stophandler.h"
45
46 class energyhistory_t;
47 struct gmx_ekindata_t;
48 struct gmx_enerdata_t;
49 struct gmx_enfrot;
50 struct gmx_mtop_t;
51 struct gmx_membed_t;
52 struct gmx_multisim_t;
53 struct gmx_output_env_t;
54 struct gmx_wallcycle;
55 struct gmx_walltime_accounting;
56 struct ObservablesHistory;
57 struct pull_t;
58 struct ReplicaExchangeParameters;
59 struct t_commrec;
60 struct t_forcerec;
61 struct t_filenm;
62 struct t_inputrec;
63 struct t_nrnb;
64 struct t_swap;
65 class t_state;
66
67 namespace gmx
68 {
69 enum class StartingBehavior;
70 class BoxDeformation;
71 class Constraints;
72 class MdrunScheduleWorkload;
73 class IMDOutputProvider;
74 struct MdModulesNotifier;
75 class ImdSession;
76 class MDLogger;
77 class MDAtoms;
78 class StopHandlerBuilder;
79 struct MdrunOptions;
80 class VirtualSitesHandler;
81
82 /*! \internal
83  * \brief The Simulator interface
84  *
85  * This is the general interface for any type of simulation type
86  * run with GROMACS. This allows having a builder return different
87  * Simulator objects based on user input.
88  */
89 class ISimulator
90 {
91 public:
92     /*! \brief The simulation run
93      *
94      * This will be called by the owner of the simulator object. To be redefined
95      * by the child classes. This function is expected to run the simulation.
96      */
97     virtual void run() = 0;
98     //! Standard destructor
99     virtual ~ISimulator() = default;
100 };
101
102 /*! \internal
103  * \brief The legacy simulator data
104  *
105  * This contains the data passed into the GROMACS simulators from
106  * the Mdrunner object.
107  */
108 class LegacySimulatorData
109 {
110 public:
111     //! The constructor
112     LegacySimulatorData(FILE*                               fplog,
113                         t_commrec*                          cr,
114                         const gmx_multisim_t*               ms,
115                         const MDLogger&                     mdlog,
116                         int                                 nfile,
117                         const t_filenm*                     fnm,
118                         const gmx_output_env_t*             oenv,
119                         const MdrunOptions&                 mdrunOptions,
120                         StartingBehavior                    startingBehavior,
121                         VirtualSitesHandler*                vsite,
122                         Constraints*                        constr,
123                         gmx_enfrot*                         enforcedRotation,
124                         BoxDeformation*                     deform,
125                         IMDOutputProvider*                  outputProvider,
126                         const MdModulesNotifier&            mdModulesNotifier,
127                         t_inputrec*                         inputrec,
128                         ImdSession*                         imdSession,
129                         pull_t*                             pull_work,
130                         t_swap*                             swap,
131                         gmx_mtop_t*                         top_global,
132                         t_state*                            state_global,
133                         ObservablesHistory*                 observablesHistory,
134                         MDAtoms*                            mdAtoms,
135                         t_nrnb*                             nrnb,
136                         gmx_wallcycle*                      wcycle,
137                         t_forcerec*                         fr,
138                         gmx_enerdata_t*                     enerd,
139                         gmx_ekindata_t*                     ekind,
140                         MdrunScheduleWorkload*              runScheduleWork,
141                         const ReplicaExchangeParameters&    replExParams,
142                         gmx_membed_t*                       membed,
143                         gmx_walltime_accounting*            walltime_accounting,
144                         std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder,
145                         bool                                doRerun) :
146         fplog(fplog),
147         cr(cr),
148         ms(ms),
149         mdlog(mdlog),
150         nfile(nfile),
151         fnm(fnm),
152         oenv(oenv),
153         mdrunOptions(mdrunOptions),
154         startingBehavior(startingBehavior),
155         vsite(vsite),
156         constr(constr),
157         enforcedRotation(enforcedRotation),
158         deform(deform),
159         outputProvider(outputProvider),
160         mdModulesNotifier(mdModulesNotifier),
161         inputrec(inputrec),
162         imdSession(imdSession),
163         pull_work(pull_work),
164         swap(swap),
165         top_global(top_global),
166         state_global(state_global),
167         observablesHistory(observablesHistory),
168         mdAtoms(mdAtoms),
169         nrnb(nrnb),
170         wcycle(wcycle),
171         fr(fr),
172         enerd(enerd),
173         ekind(ekind),
174         runScheduleWork(runScheduleWork),
175         replExParams(replExParams),
176         membed(membed),
177         walltime_accounting(walltime_accounting),
178         stopHandlerBuilder(std::move(stopHandlerBuilder)),
179         doRerun(doRerun)
180     {
181     }
182
183     //! Handles logging.
184     FILE* fplog;
185     //! Handles communication.
186     t_commrec* cr;
187     //! Coordinates multi-simulations.
188     const gmx_multisim_t* ms;
189     //! Handles logging.
190     const MDLogger& mdlog;
191     //! Count of input file options.
192     int nfile;
193     //! Content of input file options.
194     const t_filenm* fnm;
195     //! Handles writing text output.
196     const gmx_output_env_t* oenv;
197     //! Contains command-line options to mdrun.
198     const MdrunOptions& mdrunOptions;
199     //! Whether the simulation will start afresh, or restart with/without appending.
200     const StartingBehavior startingBehavior;
201     //! Handles virtual sites.
202     VirtualSitesHandler* vsite;
203     //! Handles constraints.
204     Constraints* constr;
205     //! Handles enforced rotation.
206     gmx_enfrot* enforcedRotation;
207     //! Handles box deformation.
208     BoxDeformation* deform;
209     //! Handles writing output files.
210     IMDOutputProvider* outputProvider;
211     //! Handles notifications to MdModules for checkpoint writing
212     const MdModulesNotifier& mdModulesNotifier;
213     //! Contains user input mdp options.
214     t_inputrec* inputrec;
215     //! The Interactive Molecular Dynamics session.
216     ImdSession* imdSession;
217     //! The pull work object.
218     pull_t* pull_work;
219     //! The coordinate-swapping session.
220     t_swap* swap;
221     //! Full system topology.
222     const gmx_mtop_t* top_global;
223     //! Full simulation state (only non-nullptr on master rank).
224     t_state* state_global;
225     //! History of simulation observables.
226     ObservablesHistory* observablesHistory;
227     //! Atom parameters for this domain.
228     MDAtoms* mdAtoms;
229     //! Manages flop accounting.
230     t_nrnb* nrnb;
231     //! Manages wall cycle accounting.
232     gmx_wallcycle* wcycle;
233     //! Parameters for force calculations.
234     t_forcerec* fr;
235     //! Data for energy output.
236     gmx_enerdata_t* enerd;
237     //! Kinetic energy data.
238     gmx_ekindata_t* ekind;
239     //! Schedule of work for each MD step for this task.
240     MdrunScheduleWorkload* runScheduleWork;
241     //! Parameters for replica exchange algorihtms.
242     const ReplicaExchangeParameters& replExParams;
243     //! Parameters for membrane embedding.
244     gmx_membed_t* membed;
245     //! Manages wall time accounting.
246     gmx_walltime_accounting* walltime_accounting;
247     //! Registers stop conditions
248     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder;
249     //! Whether we're doing a rerun.
250     bool doRerun;
251 };
252
253 } // namespace gmx
254
255 #endif // GMX_MDRUN_ISIMULATOR_H