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