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