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