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