703e632b34cf4bf18e1f35d9f6d8f85e0709945d
[alexxy/gromacs.git] / src / gromacs / mdrun / runner.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2017,2018,2019,2020, 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 /*! \libinternal \file
36  *
37  * \brief Declares the routine running the inetgrators.
38  *
39  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40  * \ingroup module_mdrun
41  */
42 #ifndef GMX_MDRUN_RUNNER_H
43 #define GMX_MDRUN_RUNNER_H
44
45 #include <cstdio>
46
47 #include <array>
48 #include <memory>
49
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/compat/pointers.h"
52 #include "gromacs/domdec/options.h"
53 #include "gromacs/hardware/hw_info.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdrun/mdmodules.h"
56 #include "gromacs/mdrun/simulationinputhandle.h"
57 #include "gromacs/mdrunutility/handlerestart.h"
58 #include "gromacs/mdtypes/mdrunoptions.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/gmxmpi.h"
62 #include "gromacs/utility/real.h"
63
64 #include "replicaexchange.h"
65
66 struct gmx_multisim_t;
67 struct gmx_output_env_t;
68 struct ReplicaExchangeParameters;
69 struct t_fileio;
70
71 namespace gmx
72 {
73
74 // Todo: move to forward declaration headers...
75 class MDModules;
76 class IRestraintPotential; // defined in restraint/restraintpotential.h
77 class RestraintManager;
78 class SimulationContext;
79 class StopHandlerBuilder;
80
81 /*! \libinternal \brief Runner object for supporting setup and execution of mdrun.
82  *
83  * This class has responsibility for the lifetime of data structures
84  * that exist for the life of the simulation, e.g. for logging and
85  * communication.
86  *
87  * It is also responsible for initializing data members that
88  * e.g. correspond to values potentially set by commmand-line
89  * options. Later these will be obtained directly from modules, and
90  * the results of command-line option handling returned directly to
91  * the modules, rather than propagated to them by data members of this
92  * class.
93  *
94  * \todo Most of the attributes should be declared by specific modules
95  * as command-line options. Accordingly, they do not conform to the
96  * naming scheme, because that would make for a lot of noise in the
97  * diff, only to have it change again when the options move to their
98  * modules.
99  *
100  * \todo Preparing logging and MPI contexts could probably be a
101  * higher-level responsibility, so that an Mdrunner would get made
102  * without needing to re-initialize these components (as currently
103  * happens always for the master rank, and differently for the spawned
104  * ranks with thread-MPI).
105  *
106  * \ingroup module_mdrun
107  */
108 class Mdrunner
109 {
110 public:
111     /*! \brief Builder class to manage object creation.
112      *
113      * This class is a member of gmx::Mdrunner so that it can initialize
114      * private members of gmx::Mdrunner.
115      *
116      * It is non-trivial to establish an initialized gmx::Mdrunner invariant,
117      * so objects can be obtained by clients using a Builder or a move.
118      * Clients cannot default initialize or copy gmx::Mdrunner.
119      */
120     class BuilderImplementation;
121
122     ~Mdrunner();
123
124     /*!
125      * \brief Copy not allowed.
126      *
127      * An Mdrunner has unique resources and it is not clear whether any of
128      * one of those resources should be duplicated or shared unless the
129      * specific use case is known. Either build a fresh runner or use a
130      * helper function for clearly indicated behavior. API clarification may
131      * allow unambiguous initialization by copy in future versions.
132      *
133      * \{
134      */
135     Mdrunner(const Mdrunner&) = delete;
136     Mdrunner& operator=(const Mdrunner&) = delete;
137     /* \} */
138
139     /*!
140      * \brief Mdrunner objects can be passed by value via move semantics.
141      *
142      * \param handle runner instance to be moved from.
143      * \{
144      */
145     Mdrunner(Mdrunner&& handle) noexcept;
146     Mdrunner& operator=(Mdrunner&& handle) noexcept;
147     /* \} */
148
149     /*! \brief Driver routine, that calls the different simulation methods. */
150     /*!
151      * Currently, thread-MPI does not spawn threads until during mdrunner() and parallelism
152      * is not initialized until some time during this call...
153      */
154     int mdrunner();
155
156     /*!
157      * \brief Add a potential to be evaluated during MD integration.
158      *
159      * \param restraint MD restraint potential to apply
160      * \param name User-friendly plain-text name to uniquely identify the puller
161      *
162      * This implementation attaches an object providing the gmx::IRestraintPotential
163      * interface.
164      * \todo Mdrunner should fetch such resources from the SimulationContext
165      * rather than offering this public interface.
166      */
167     void addPotential(std::shared_ptr<IRestraintPotential> restraint, const std::string& name);
168
169     /*! \brief Prepare the thread-MPI communicator to have \c
170      * numThreadsToLaunch ranks, by spawning new thread-MPI
171      * threads.
172      *
173      * Called by mdrunner() to start a specific number of threads
174      * (including the main thread) for thread-parallel runs. This
175      * in turn calls mdrunner() for each thread. */
176     void spawnThreads(int numThreadsToLaunch);
177
178     /*! \brief Initializes a new Mdrunner from the master.
179      *
180      * Run this method in a new thread from a master runner to get additional
181      * workers on spawned threads.
182      *
183      * \returns New Mdrunner instance suitable for thread-MPI work on new ranks.
184      *
185      * \internal
186      * \todo clarify (multiple) invariants during MD runner start-up.
187      * The runner state before and after launching threads is distinct enough that
188      * it should be codified in the invariants of different classes. That would
189      * mean that the object returned by this method would be of a different type
190      * than the object held by the client up to the point of call, and its name
191      * would be changed to "launchOnSpawnedThread" or something not including the
192      * word "clone".
193      */
194     Mdrunner cloneOnSpawnedThread() const;
195
196 private:
197     /*! \brief Constructor. */
198     explicit Mdrunner(std::unique_ptr<MDModules> mdModules);
199
200     //! Parallelism-related user options.
201     gmx_hw_opt_t hw_opt;
202
203     //! Filenames and properties from command-line argument values.
204     ArrayRef<const t_filenm> filenames;
205
206     /*! \brief Output context for writing text files
207      *
208      * \internal
209      * \todo push this data member down when the information can be queried from an encapsulated resource.
210      */
211     gmx_output_env_t* oenv = nullptr;
212     //! Ongoing collection of mdrun options
213     MdrunOptions mdrunOptions;
214     //! Options for the domain decomposition.
215     DomdecOptions domdecOptions;
216
217     /*! \brief Target short-range interations for "cpu", "gpu", or "auto". Default is "auto".
218      *
219      * \internal
220      * \todo replace with string or enum class and initialize with sensible value.
221      */
222     const char* nbpu_opt = nullptr;
223
224     /*! \brief Target long-range interactions for "cpu", "gpu", or "auto". Default is "auto".
225      *
226      * \internal
227      * \todo replace with string or enum class and initialize with sensible value.
228      */
229     const char* pme_opt = nullptr;
230
231     /*! \brief Target long-range interactions FFT/solve stages for "cpu", "gpu", or "auto". Default is "auto".
232      *
233      * \internal
234      * \todo replace with string or enum class and initialize with sensible value.
235      */
236     const char* pme_fft_opt = nullptr;
237
238     /*! \brief Target bonded interations for "cpu", "gpu", or "auto". Default is "auto".
239      *
240      * \internal
241      * \todo replace with string or enum class and initialize with sensible value.
242      */
243     const char* bonded_opt = nullptr;
244
245     /*! \brief Target update calculation for "cpu", "gpu", or "auto". Default is "auto".
246      *
247      * \internal
248      * \todo replace with string or enum class and initialize with sensible value.
249      */
250     const char* update_opt = nullptr;
251
252     //! Command-line override for the duration of a neighbor list with the Verlet scheme.
253     int nstlist_cmdline = 0;
254     //! Parameters for replica-exchange simulations.
255     ReplicaExchangeParameters replExParams;
256     //! Print a warning if any force is larger than this (in kJ/mol nm).
257     real pforce = -1;
258
259     //! Handle to file used for logging.
260     LogFilePtr logFileGuard = nullptr;
261     //! \brief Non-owning handle to file used for logging.
262     t_fileio* logFileHandle = nullptr;
263
264     /*! \brief Non-owning handle to world communication data structure for task assigment.
265      *
266      * With real MPI, gets a value from the SimulationContext
267      * supplied to the MdrunnerBuilder. With thread-MPI gets a
268      * value after threads have been spawned. */
269     MPI_Comm libraryWorldCommunicator = MPI_COMM_NULL;
270
271     /*! \brief Non-owning handle to communication data structure for the current simulation.
272      *
273      * With real MPI, gets a value from the SimulationContext
274      * supplied to the MdrunnerBuilder. With thread-MPI gets a
275      * value after threads have been spawned. */
276     MPI_Comm simulationCommunicator = MPI_COMM_NULL;
277
278     //! \brief Non-owning handle to multi-simulation handler.
279     gmx_multisim_t* ms = nullptr;
280
281     //! Whether the simulation will start afresh, or restart with/without appending.
282     StartingBehavior startingBehavior = StartingBehavior::NewSimulation;
283
284     /*!
285      * \brief Handle to restraints manager for the current process.
286      *
287      * \internal
288      * Use opaque pointer for this implementation detail.
289      */
290     std::unique_ptr<RestraintManager> restraintManager_;
291
292     /*!
293      * \brief Builder for stop signal handler
294      *
295      * Optionally provided through MdrunnerBuilder. Client may create a
296      * StopHandlerBuilder and register any number of signal providers before
297      * launching the Mdrunner.
298      *
299      * Default is an empty signal handler that will have local signal issuers
300      * added after being passed into the integrator.
301      *
302      * \internal
303      * We do not need a full type specification here, so we use an opaque pointer.
304      */
305     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_;
306     //! The modules that comprise mdrun.
307     std::unique_ptr<MDModules> mdModules_;
308
309     //! Non-owning handle to the results of the hardware detection.
310     const gmx_hw_info_t* hwinfo_ = nullptr;
311
312     /*!
313      * \brief Holds simulation input specification provided by client, if any.
314      *
315      * If present on any instance (rank) of a simulation runner, an identical
316      * (or compatible) SimulationInput must be held on all cooperating instances.
317      */
318     SimulationInputHandle inputHolder_;
319 };
320
321 /*! \libinternal
322  * \brief Build a gmx::Mdrunner.
323  *
324  * Client code (such as `gmx mdrun`) uses this builder to get an initialized Mdrunner.
325  *
326  * A builder allows the library to ensure that client code cannot obtain an
327  * uninitialized or partially initialized runner by refusing to build() if the
328  * client has not provided sufficient or self-consistent direction. Director
329  * code can be implemented for different user interfaces, encapsulating any
330  * run-time functionality that does not belong in the library MD code, such
331  * as command-line option processing or interfacing to external libraries.
332  *
333  * \ingroup module_mdrun
334  *
335  * \internal
336  *
337  * The initial Builder implementation is neither extensible at run time nor
338  * at compile time. Future implementations should evolve to compose the runner,
339  * rather than just consolidating the parameters for initialization, but there
340  * is not yet a firm design for how flexibly module code will be coupled to
341  * the builder and how much of the client interface will be in this Builder
342  * versus Builders provided by the various modules.
343  *
344  * The named components for the initial builder implementation are descriptive
345  * of the state of mdrun at the time, and are not intended to be prescriptive of
346  * future design.
347  * The probable course of GROMACS development is for the modular components that
348  * support MD simulation to independently express their input parameters (required
349  * and optional) and to provide some sort of help to the UI for input preparation.
350  * If each module provides or aids the instantiation of a Director
351  * for the client code, the Directors could be constructed with a handle to this
352  * Builder and it would not need a public interface.
353  *
354  * As the modules are more clearly encapsulated, each module can provide its own
355  * builder, user interface helpers, and/or composable Director code.
356  * The runner and client code will also have to be updated as appropriate
357  * default behavior is clarified for
358  * (a) default behavior of client when user does not provide input,
359  * (b) default behavior of builder when client does not provide input, and
360  * (c) default behavior of runner when builder does not provide input.
361  */
362 class MdrunnerBuilder final
363 {
364 public:
365     /*!
366      * \brief Constructor requires a handle to a SimulationContext to share.
367      *
368      * \param mdModules  The handle to the set of modules active in mdrun
369      * \param context    Required handle to simulation context
370      *
371      * The calling code must guarantee that the
372      * pointer remains valid for the lifetime of the builder, and that the
373      * resources retrieved from the context remain valid for the lifetime of
374      * the runner produced.
375      */
376     explicit MdrunnerBuilder(std::unique_ptr<MDModules>           mdModules,
377                              compat::not_null<SimulationContext*> context);
378
379     //! \cond
380     MdrunnerBuilder()                       = delete;
381     MdrunnerBuilder(const MdrunnerBuilder&) = delete;
382     MdrunnerBuilder& operator=(const MdrunnerBuilder&) = delete;
383     //! \endcond
384
385     /*! \brief Allow transfer of ownership with move semantics.
386      *
387      * \param builder source object to transfer.
388      *
389      * \{
390      */
391     MdrunnerBuilder(MdrunnerBuilder&& builder) noexcept;
392     MdrunnerBuilder& operator=(MdrunnerBuilder&& builder) noexcept;
393     //! \}
394
395     /*!
396      * \brief Get ownership of an initialized gmx::Mdrunner.
397      *
398      * After build() is called, the Builder object should not be used
399      * again. It is an error to call build without first calling all builder
400      * methods described as "required."
401      *
402      * \return A new Mdrunner.
403      *
404      * \throws APIError if a required component has not been added before calling build().
405      */
406     Mdrunner build();
407
408     /*!
409      * \brief Supply the result of hardware detection to the gmx::Mdrunner
410      *
411      * \param hwinfo  Non-owning not-null handle to result of hardware detection.
412      *
413      * \todo It would be better to express this as either a not-null const pointer or
414      * a const reference, but neither of those is consistent with incremental
415      * building of an object. This motivates future work to be able to make a deep copy
416      * of the detection result. See https://gitlab.com/gromacs/gromacs/-/issues/3650 */
417     MdrunnerBuilder& addHardwareDetectionResult(const gmx_hw_info_t* hwinfo);
418
419     /*!
420      * \brief Set up non-bonded short-range force calculations.
421      *
422      * Required. Director code must provide valid options for the non-bonded
423      * interaction code. The builder does not apply any defaults.
424      *
425      * \param nbpu_opt Target short-range interactions for "cpu", "gpu", or "auto".
426      *
427      * Calling must guarantee that the pointed-to C string is valid through
428      * simulation launch.
429      *
430      * \internal
431      * \todo Replace with string or enum that we can have sensible defaults for.
432      * \todo Either the Builder or modular Director code should provide sensible defaults.
433      */
434     MdrunnerBuilder& addNonBonded(const char* nbpu_opt);
435
436     /*!
437      * \brief Set up long-range electrostatics calculations.
438      *
439      * Required. Director code should provide valid options for PME electrostatics,
440      * whether or not PME electrostatics are used. The builder does not apply
441      * any defaults, so client code should be prepared to provide (e.g.) "auto"
442      * in the event no user input or logic provides an alternative argument.
443      *
444      * \param pme_opt Target long-range interactions for "cpu", "gpu", or "auto".
445      * \param pme_fft_opt Target long-range interactions FFT/solve stages for "cpu", "gpu", or "auto".
446      *
447      * Calling must guarantee that the pointed-to C strings are valid through
448      * simulation launch.
449      *
450      * \internal
451      * The arguments are passed as references to elements of arrays of C strings.
452      * \todo Replace with modern strings or (better) enum classes.
453      * \todo Make optional and/or encapsulate into electrostatics module.
454      */
455     MdrunnerBuilder& addElectrostatics(const char* pme_opt, const char* pme_fft_opt);
456
457     /*!
458      * \brief Assign responsibility for tasks for bonded interactions.
459      *
460      * Required. Director code should provide valid options for
461      * bonded interaction task assignment, whether or not such
462      * interactions are present. The builder does not apply any
463      * defaults, so client code should be prepared to provide
464      * (e.g.) "auto" in the event no user input or logic provides
465      * an alternative argument.
466      *
467      * \param bonded_opt Target bonded interactions for "cpu", "gpu", or "auto".
468      *
469      * Calling must guarantee that the pointed-to C strings are valid through
470      * simulation launch.
471      *
472      * \internal
473      * The arguments are passed as references to elements of arrays of C strings.
474      * \todo Replace with modern strings or (better) enum classes.
475      * \todo Make optional and/or encapsulate into task assignment module.
476      */
477     MdrunnerBuilder& addBondedTaskAssignment(const char* bonded_opt);
478
479     /*! \brief
480      * Assign responsibility for tasks for update and constrain calculation.
481      *
482      * Required. Director code should provide valid options for
483      * update and constraint task assignment. The builder does not apply any
484      * defaults, so client code should be prepared to provide
485      * (e.g.) "auto" in the event no user input or logic provides
486      * an alternative argument.
487      *
488      * \param[in] update_opt Target update calculation for "cpu", "gpu", or "auto".
489      *
490      * Calling must guarantee that the pointed-to C strings are valid through
491      * simulation launch.
492      *
493      * \internal
494      * The arguments are passed as references to elements of arrays of C strings.
495      * \todo Replace with modern strings or (better) enum classes.
496      * \todo Make optional and/or encapsulate into task assignment module.
497      */
498     MdrunnerBuilder& addUpdateTaskAssignment(const char* update_opt);
499
500     /*!
501      * \brief Set MD options not owned by some other module.
502      *
503      * Optional. Override simulation parameters
504      *
505      * \param options structure to copy
506      * \param forceWarningThreshold Print a warning if any force is larger than this (in kJ/mol nm) (default -1)
507      * \param startingBehavior Whether the simulation will start afresh, or restart with/without appending.
508      *
509      * \internal
510      * \todo Map these parameters to more appropriate encapsulating types.
511      * Find a better way to indicate "unspecified" than a magic value of the parameter type.
512      */
513     MdrunnerBuilder& addSimulationMethod(const MdrunOptions& options,
514                                          real                forceWarningThreshold,
515                                          StartingBehavior    startingBehavior);
516
517     /*!
518      * \brief Set the domain decomposition module.
519      *
520      * Optional. Overrides default constructed DomdecOptions if provided.
521      *
522      * \param options options with which to construct domain decomposition.
523      *
524      * \internal
525      * \todo revisit whether we should be passing this parameter struct or a higher-level handle of some sort.
526      */
527     MdrunnerBuilder& addDomainDecomposition(const DomdecOptions& options);
528
529     /*!
530      * \brief Set Verlet list manager.
531      *
532      * Optional. Neighbor list existence, type, and parameters are mostly determined
533      * by the simulation parameters loaded elsewhere. This is just an override.
534      *
535      * \param rebuildInterval override for the duration of a neighbor list with the Verlet scheme.
536      */
537     MdrunnerBuilder& addNeighborList(int rebuildInterval);
538
539     /*!
540      * \brief Set replica exchange manager.
541      *
542      * Optional. For guidance on preparing a valid ReplicaExchangeParameters
543      * value, refer to the details in mdrun.cpp, the `t_pargs pa[]` defined there,
544      * and the action of parse_common_args() with regards to that structure.
545      * If not provided by client, a default constructed ReplicaExchangeParameters
546      * is used.
547      *
548      * \param params parameters with which to set up replica exchange.
549      *
550      * \internal
551      * \todo revisit whether we should be passing this parameter struct or a higher-level handle of some sort.
552      */
553     MdrunnerBuilder& addReplicaExchange(const ReplicaExchangeParameters& params);
554
555     /*!
556      * \brief Specify parameters determining hardware resource allocation.
557      *
558      * Optional. If not provided, default-constructed gmx_hw_opt_t will be used.
559      *
560      * \param hardwareOptions Parallelism-related user options.
561      */
562     MdrunnerBuilder& addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
563
564     /*!
565      * \brief Provide the filenames options structure with option values chosen
566      *
567      * Required. The object is assumed to have been updated by
568      * parse_common_args or equivalent.
569      *
570      * \param filenames Filenames and properties from command-line argument values or defaults.
571      *
572      * \internal
573      * \todo Modules should manage their own filename options and defaults.
574      */
575     MdrunnerBuilder& addFilenames(ArrayRef<const t_filenm> filenames);
576
577     /*!
578      * \brief Provide parameters for setting up output environment.
579      *
580      * Required. Handle is assumed to have been produced by output_env_init
581      * as in parse_common_args.
582      *
583      * \param outputEnvironment Output context for writing text files.
584      *
585      * \internal
586      * \todo Allow client code to set up output environment and provide as a resource.
587      * This parameter is used to set up resources that are dependent on the execution
588      * environment and API context. Such resources should be retrieved by the simulator
589      * from a client-provided resource, but currently the resources are only fully
590      * initialized in Mdrunner.
591      */
592     MdrunnerBuilder& addOutputEnvironment(gmx_output_env_t* outputEnvironment);
593
594     /*!
595      * \brief Provide the filehandle pointer to be used for the MD log.
596      *
597      * Required. Either nullptr if no log should be written, or
598      * valid and open reading for writing.
599      *
600      * \param logFileHandle Non-owning handle to file used for logging.
601      * \internal
602      */
603     MdrunnerBuilder& addLogFile(t_fileio* logFileHandle);
604
605     /*!
606      * \brief Provide a StopHandlerBuilder for the MD stop signal handling.
607      *
608      * Optional. Defaults to empty.
609      *
610      * Client may provide additional (non-default) issuers of simulation stop
611      * signals by preconfiguring the StopHandlerBuilder used later when the
612      * simulation runs.
613      *
614      * \param builder
615      */
616     MdrunnerBuilder& addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
617
618     /*!
619      * \brief Acquire a handle to the SimulationInput.
620      *
621      * Required. SimulationInput will be taking responsibility for some of the
622      * input provided through other methods, such as addFilenames.
623      *
624      * See also issue https://gitlab.com/gromacs/gromacs/-/issues/3374
625      *
626      * \param input Shared ownership of a SimulationInput.
627      */
628     MdrunnerBuilder& addInput(SimulationInputHandle input);
629
630     ~MdrunnerBuilder();
631
632 private:
633     std::unique_ptr<Mdrunner::BuilderImplementation> impl_;
634 };
635
636 } // namespace gmx
637
638 #endif // GMX_MDRUN_RUNNER_H