Unify handling of GMX_ENABLE_GPU_TIMING and GMX_DISABLE_GPU_TIMING
[alexxy/gromacs.git] / src / gromacs / mdrun / simulatorbuilder.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 simulator builder for mdrun
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_mdrun
40  */
41 #ifndef GMX_MDRUN_SIMULATORBUILDER_H
42 #define GMX_MDRUN_SIMULATORBUILDER_H
43
44 #include <memory>
45
46
47 class energyhistory_t;
48 class gmx_ekindata_t;
49 struct gmx_enerdata_t;
50 struct gmx_enfrot;
51 struct gmx_localtop_t;
52 struct gmx_mtop_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_filenm;
62 struct t_forcerec;
63 struct t_inputrec;
64 struct t_nrnb;
65 class t_state;
66 struct t_swap;
67
68 namespace gmx
69 {
70 class BoxDeformation;
71 class Constraints;
72 class IMDOutputProvider;
73 class ImdSession;
74 class ISimulator;
75 class MdrunScheduleWorkload;
76 class MembedHolder;
77 class MDAtoms;
78 class MDLogger;
79 struct MDModulesNotifiers;
80 struct MdrunOptions;
81 class ObservablesReducerBuilder;
82 class ReadCheckpointDataHolder;
83 enum class StartingBehavior;
84 class StopHandlerBuilder;
85 class VirtualSitesHandler;
86
87 /*! \brief
88  * Simulation configuation settings.
89  */
90 struct SimulatorConfig
91 {
92 public:
93     //! Build from settings for this simulation.
94     SimulatorConfig(const MdrunOptions&    mdrunOptions,
95                     StartingBehavior       startingBehavior,
96                     MdrunScheduleWorkload* runScheduleWork) :
97         mdrunOptions_(mdrunOptions), startingBehavior_(startingBehavior), runScheduleWork_(runScheduleWork)
98     {
99     }
100     // TODO: Specify copy and move semantics.
101
102     //! Handle to user options.
103     const MdrunOptions& mdrunOptions_;
104     //! How are we starting the simulation.
105     StartingBehavior startingBehavior_;
106     //! How are we scheduling the tasks for this simulation.
107     MdrunScheduleWorkload* runScheduleWork_;
108 };
109
110
111 /*! \brief
112  * Data for a specific simulation state.
113  *
114  * \todo Think of a better name and annoy people that forget
115  *       to add documentation for their code.
116  */
117 struct SimulatorStateData
118 {
119     //! Build collection of current state data.
120     SimulatorStateData(t_state*            globalState,
121                        t_state*            localState,
122                        ObservablesHistory* observablesHistory,
123                        gmx_enerdata_t*     enerdata,
124                        gmx_ekindata_t*     ekindata) :
125         globalState_p(globalState),
126         localState_p(localState),
127         observablesHistory_p(observablesHistory),
128         enerdata_p(enerdata),
129         ekindata_p(ekindata)
130     {
131     }
132
133     //! Can perform copy of current state.
134     SimulatorStateData(const SimulatorStateData& simulatorStateData) = default;
135
136     //! Handle to global state of the simulation.
137     t_state* globalState_p;
138     //! Handle to local state of the simulation.
139     t_state* localState_p;
140     //! Handle to current simulation history.
141     ObservablesHistory* observablesHistory_p;
142     //! Handle to collected data for energy groups.
143     gmx_enerdata_t* enerdata_p;
144     //! Handle to collected data for kinectic energy.
145     gmx_ekindata_t* ekindata_p;
146 };
147
148 /*! \brief
149  * Collection of environmental information for a simulation.
150  *
151  * \todo Fix doxygen checks.
152  */
153 class SimulatorEnv
154 {
155 public:
156     //! Build from current simulation environment.
157     SimulatorEnv(FILE*                      fplog,
158                  t_commrec*                 commRec,
159                  gmx_multisim_t*            multisimCommRec,
160                  const MDLogger&            logger,
161                  gmx_output_env_t*          outputEnv,
162                  ObservablesReducerBuilder* observablesReducerBuilder) :
163         fplog_{ fplog }, commRec_{ commRec }, multisimCommRec_{ multisimCommRec }, logger_{ logger }, outputEnv_{ outputEnv }, observablesReducerBuilder_{ observablesReducerBuilder }
164     {
165     }
166
167     //! Handle to log file.
168     FILE* fplog_;
169     //! Handle to communication record.
170     t_commrec* commRec_;
171     //! Handle to multisim communication record.
172     const gmx_multisim_t* multisimCommRec_;
173     //! Handle to propper logging framework.
174     const MDLogger& logger_;
175     //! Handle to file output handling.
176     const gmx_output_env_t* outputEnv_;
177     //! Builder for coordinator of reduction for observables
178     ObservablesReducerBuilder* observablesReducerBuilder_;
179 };
180
181 /*! \brief
182  * Collection of profiling information.
183  */
184 class Profiling
185 {
186 public:
187     //! Build profiling information collection.
188     Profiling(t_nrnb* nrnb, gmx_walltime_accounting* walltimeAccounting, gmx_wallcycle* wallCycle) :
189         nrnb(nrnb), wallCycle(wallCycle), walltimeAccounting(walltimeAccounting)
190     {
191     }
192
193     //! Handle to datastructure.
194     t_nrnb* nrnb;
195     //! Handle to wallcycle stuff.
196     gmx_wallcycle* wallCycle;
197     //! Handle to wallcycle time accounting stuff.
198     gmx_walltime_accounting* walltimeAccounting;
199 };
200
201 /*! \brief
202  * Collection of constraint parameters.
203  */
204 class ConstraintsParam
205 {
206 public:
207     //! Build collection with handle to actual objects.
208     ConstraintsParam(Constraints* constraints, gmx_enfrot* enforcedRotation, VirtualSitesHandler* vSite) :
209         constr(constraints), enforcedRotation(enforcedRotation), vsite(vSite)
210     {
211     }
212
213     //! Handle to constraint object.
214     Constraints* constr;
215     //! Handle to information about using enforced rotation.
216     gmx_enfrot* enforcedRotation;
217     //! Handle to vsite stuff.
218     VirtualSitesHandler* vsite;
219 };
220
221 /*! \brief
222  * Collection of legacy input information.
223  */
224 class LegacyInput
225 {
226 public:
227     //! Build collection from legacy input data.
228     LegacyInput(int filenamesSize, const t_filenm* filenamesData, t_inputrec* inputRec, t_forcerec* forceRec) :
229         numFile(filenamesSize), filenames(filenamesData), inputrec(inputRec), forceRec(forceRec)
230     {
231     }
232
233     //! Number of input files.
234     int numFile;
235     //! File names.
236     const t_filenm* filenames;
237     //! Handle to simulation input record.
238     t_inputrec* inputrec;
239     //! Handle to simulation force record.
240     t_forcerec* forceRec;
241 };
242
243 /*! \brief SimulatorBuilder parameter type for InteractiveMD.
244  *
245  * Conveys a non-owning pointer to implementation details.
246  *
247  * \todo If adding doxygen stubs actual add the full stub.
248  */
249 class InteractiveMD
250 {
251 public:
252     //! Create handle to IMD information.
253     explicit InteractiveMD(ImdSession* imdSession) : imdSession(imdSession) {}
254
255     //! Internal handle to IMD info.
256     ImdSession* imdSession;
257 };
258
259 class SimulatorModules
260 {
261 public:
262     SimulatorModules(IMDOutputProvider* mdOutputProvider, const MDModulesNotifiers& notifiers) :
263         outputProvider(mdOutputProvider), mdModulesNotifiers(notifiers)
264     {
265     }
266
267     IMDOutputProvider*        outputProvider;
268     const MDModulesNotifiers& mdModulesNotifiers;
269 };
270
271 class CenterOfMassPulling
272 {
273 public:
274     explicit CenterOfMassPulling(pull_t* pullWork) : pull_work(pullWork) {}
275
276     pull_t* pull_work;
277 };
278
279 /*! \brief
280  * Parameter type for IonSwapping SimulatorBuilder component.
281  *
282  * Conveys a non-owning pointer to implementation details.
283  *
284  * \todo Add full information.
285  */
286 class IonSwapping
287 {
288 public:
289     //! Create handle.
290     IonSwapping(t_swap* ionSwap) : ionSwap(ionSwap) {}
291
292     //! Internal storage for handle.
293     t_swap* ionSwap;
294 };
295
296 /*! \brief
297  * Collection of handles to topology information.
298  */
299 class TopologyData
300 {
301 public:
302     //! Build collection from simulation data.
303     TopologyData(const gmx_mtop_t& globalTopology, gmx_localtop_t* localTopology, MDAtoms* mdAtoms) :
304         top_global(globalTopology), localTopology(localTopology), mdAtoms(mdAtoms)
305     {
306     }
307
308     //! Handle to global simulation topology.
309     const gmx_mtop_t& top_global;
310     //! Handle to local simulation topology.
311     gmx_localtop_t* localTopology;
312     //! Handle to information about MDAtoms.
313     MDAtoms* mdAtoms;
314 };
315
316 /*! \brief
317  * Handle to information about the box.
318  *
319  * Design note: The client may own the BoxDeformation via std::unique_ptr, but we are not
320  * transferring ownership at this time. (May be the subject of future changes.)
321  */
322 class BoxDeformationHandle
323 {
324 public:
325     //! Build handle to box stuff.
326     BoxDeformationHandle(BoxDeformation* boxDeformation) : deform(boxDeformation) {}
327
328     //! Internal storage for handle.
329     BoxDeformation* deform;
330 };
331
332 /*! \libinternal
333  * \brief Class preparing the creation of Simulator objects
334  *
335  * Objects of this class build Simulator objects, which in turn are used to
336  * run molecular simulations.
337  */
338 class SimulatorBuilder
339 {
340 public:
341     void add(MembedHolder&& membedHolder);
342
343     void add(std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder)
344     {
345         stopHandlerBuilder_ = std::move(stopHandlerBuilder);
346     }
347
348     void add(SimulatorStateData&& simulatorStateData)
349     {
350         simulatorStateData_ = std::make_unique<SimulatorStateData>(simulatorStateData);
351     }
352
353     void add(SimulatorConfig&& simulatorConfig)
354     {
355         // Note: SimulatorConfig appears to the compiler to be trivially copyable,
356         // but this may not be safe and may change in the future.
357         simulatorConfig_ = std::make_unique<SimulatorConfig>(simulatorConfig);
358     }
359
360     void add(SimulatorEnv&& simulatorEnv)
361     {
362         simulatorEnv_ = std::make_unique<SimulatorEnv>(simulatorEnv);
363     }
364
365     void add(Profiling&& profiling) { profiling_ = std::make_unique<Profiling>(profiling); }
366
367     void add(ConstraintsParam&& constraintsParam)
368     {
369         constraintsParam_ = std::make_unique<ConstraintsParam>(constraintsParam);
370     }
371
372     void add(LegacyInput&& legacyInput)
373     {
374         legacyInput_ = std::make_unique<LegacyInput>(legacyInput);
375     }
376
377     void add(ReplicaExchangeParameters&& replicaExchangeParameters);
378
379     void add(InteractiveMD&& interactiveMd)
380     {
381         interactiveMD_ = std::make_unique<InteractiveMD>(interactiveMd);
382     }
383
384     void add(SimulatorModules&& simulatorModules)
385     {
386         simulatorModules_ = std::make_unique<SimulatorModules>(simulatorModules);
387     }
388
389     void add(CenterOfMassPulling&& centerOfMassPulling)
390     {
391         centerOfMassPulling_ = std::make_unique<CenterOfMassPulling>(centerOfMassPulling);
392     }
393
394     void add(IonSwapping&& ionSwapping)
395     {
396         ionSwapping_ = std::make_unique<IonSwapping>(ionSwapping);
397     }
398
399     void add(TopologyData&& topologyData)
400     {
401         topologyData_ = std::make_unique<TopologyData>(topologyData);
402     }
403
404     void add(BoxDeformationHandle&& boxDeformation)
405     {
406         boxDeformation_ = std::make_unique<BoxDeformationHandle>(boxDeformation);
407     }
408
409     /*!
410      * \brief Pass the read checkpoint data for modular simulator
411      *
412      * Note that this is currently the point at which the ReadCheckpointDataHolder
413      * is fully filled. Consequently it stops being an object at which read
414      * operations from file are targeted, and becomes a read-only object from
415      * which elements read their data to recreate an earlier internal state.
416      *
417      * Currently, this behavior change is not enforced. Once input reading and
418      * simulator builder have matured, these restrictions could be imposed.
419      *
420      * See #3656
421      */
422     void add(std::unique_ptr<ReadCheckpointDataHolder> modularSimulatorCheckpointData);
423
424     /*! \brief Build a Simulator object based on input data
425      *
426      * Return a pointer to a simulation object. The use of a parameter
427      * pack insulates the builder from changes to the arguments of the
428      * Simulator objects.
429      *
430      * \throws gmx::APIError if expected set-up methods have not been called before build()
431      *
432      * \return  Unique pointer to a Simulator object
433      */
434     std::unique_ptr<ISimulator> build(bool useModularSimulator);
435
436 private:
437     // Note: we use std::unique_ptr instead of std::optional because we want to
438     // allow for opaque types at the discretion of the module developer.
439     /*! \brief Collection of handles to  individual information. */
440     /*! \{ */
441     std::unique_ptr<SimulatorConfig>           simulatorConfig_;
442     std::unique_ptr<MembedHolder>              membedHolder_;
443     std::unique_ptr<StopHandlerBuilder>        stopHandlerBuilder_;
444     std::unique_ptr<SimulatorStateData>        simulatorStateData_;
445     std::unique_ptr<SimulatorEnv>              simulatorEnv_;
446     std::unique_ptr<Profiling>                 profiling_;
447     std::unique_ptr<ConstraintsParam>          constraintsParam_;
448     std::unique_ptr<LegacyInput>               legacyInput_;
449     std::unique_ptr<ReplicaExchangeParameters> replicaExchangeParameters_;
450     std::unique_ptr<InteractiveMD>             interactiveMD_;
451     std::unique_ptr<SimulatorModules>          simulatorModules_;
452     std::unique_ptr<CenterOfMassPulling>       centerOfMassPulling_;
453     std::unique_ptr<IonSwapping>               ionSwapping_;
454     std::unique_ptr<TopologyData>              topologyData_;
455     std::unique_ptr<BoxDeformationHandle>      boxDeformation_;
456     //! Contains checkpointing data for the modular simulator
457     std::unique_ptr<ReadCheckpointDataHolder> modularSimulatorCheckpointData_;
458     /*! \} */
459 };
460
461 } // namespace gmx
462
463 #endif // GMX_MDRUN_SIMULATORBUILDER_SIMULATORBUILDER_H