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