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