Update clang-tidy to clang version 8
[gromacs.git] / src / gromacs / mdrun / runner.h
blob30b44147b7553b0f1fdb1780e89097a736c9860b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \libinternal \file
37 * \brief Declares the routine running the inetgrators.
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_mdrun
42 #ifndef GMX_MDRUN_RUNNER_H
43 #define GMX_MDRUN_RUNNER_H
45 #include <cstdio>
47 #include <array>
48 #include <memory>
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/real.h"
61 #include "replicaexchange.h"
63 struct gmx_output_env_t;
64 struct ReplicaExchangeParameters;
65 struct t_commrec;
66 struct t_fileio;
68 namespace gmx
71 // Todo: move to forward declaration headers...
72 class MDModules;
73 class IRestraintPotential; // defined in restraint/restraintpotential.h
74 class RestraintManager;
75 class SimulationContext;
76 class StopHandlerBuilder;
78 //! Work-around for GCC bug 58265
79 constexpr bool BUGFREE_NOEXCEPT_STRING = std::is_nothrow_move_assignable<std::string>::value;
81 /*! \libinternal \brief Runner object for supporting setup and execution of mdrun.
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.
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.
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.
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).
106 * \ingroup module_mdrun
108 class Mdrunner
110 public:
111 /*! \brief Builder class to manage object creation.
113 * This class is a member of gmx::Mdrunner so that it can initialize
114 * private members of gmx::Mdrunner.
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.
120 class BuilderImplementation;
122 ~Mdrunner();
125 * \brief Copy not allowed.
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.
133 * \{
135 Mdrunner(const Mdrunner &) = delete;
136 Mdrunner &operator=(const Mdrunner &) = delete;
137 /* \} */
140 * \brief Mdrunner objects can be passed by value via move semantics.
142 * \param handle runner instance to be moved from.
143 * \{
145 Mdrunner(Mdrunner &&handle) noexcept;
146 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
147 Mdrunner &operator=(Mdrunner &&handle) noexcept(BUGFREE_NOEXCEPT_STRING);
148 /* \} */
150 /*! \brief Driver routine, that calls the different simulation methods. */
152 * Currently, thread-MPI does not spawn threads until during mdrunner() and parallelism
153 * is not initialized until some time during this call...
155 int mdrunner();
158 * \brief Add a potential to be evaluated during MD integration.
160 * \param restraint MD restraint potential to apply
161 * \param name User-friendly plain-text name to uniquely identify the puller
163 * This implementation attaches an object providing the gmx::IRestraintPotential
164 * interface.
165 * \todo Mdrunner should fetch such resources from the SimulationContext
166 * rather than offering this public interface.
168 void addPotential(std::shared_ptr<IRestraintPotential> restraint,
169 const std::string &name);
171 //! Called when thread-MPI spawns threads.
172 t_commrec *spawnThreads(int numThreadsToLaunch) const;
174 /*! \brief Initializes a new Mdrunner from the master.
176 * Run this method in a new thread from a master runner to get additional
177 * workers on spawned threads.
179 * \returns New Mdrunner instance suitable for thread-MPI work on new ranks.
181 * \internal
182 * \todo clarify (multiple) invariants during MD runner start-up.
183 * The runner state before and after launching threads is distinct enough that
184 * it should be codified in the invariants of different classes. That would
185 * mean that the object returned by this method would be of a different type
186 * than the object held by the client up to the point of call, and its name
187 * would be changed to "launchOnSpawnedThread" or something not including the
188 * word "clone".
190 Mdrunner cloneOnSpawnedThread() const;
192 private:
193 /*! \brief Constructor. */
194 explicit Mdrunner(std::unique_ptr<MDModules> mdModules);
196 //! Parallelism-related user options.
197 gmx_hw_opt_t hw_opt;
199 //! Filenames and properties from command-line argument values.
200 ArrayRef<const t_filenm> filenames;
202 /*! \brief Output context for writing text files
204 * \internal
205 * \todo push this data member down when the information can be queried from an encapsulated resource.
207 gmx_output_env_t *oenv = nullptr;
208 //! Ongoing collection of mdrun options
209 MdrunOptions mdrunOptions;
210 //! Options for the domain decomposition.
211 DomdecOptions domdecOptions;
213 /*! \brief Target short-range interations for "cpu", "gpu", or "auto". Default is "auto".
215 * \internal
216 * \todo replace with string or enum class and initialize with sensible value.
218 const char *nbpu_opt = nullptr;
220 /*! \brief Target long-range interactions for "cpu", "gpu", or "auto". Default is "auto".
222 * \internal
223 * \todo replace with string or enum class and initialize with sensible value.
225 const char *pme_opt = nullptr;
227 /*! \brief Target long-range interactions FFT/solve stages for "cpu", "gpu", or "auto". Default is "auto".
229 * \internal
230 * \todo replace with string or enum class and initialize with sensible value.
232 const char *pme_fft_opt = nullptr;
234 /*! \brief Target bonded interations for "cpu", "gpu", or "auto". Default is "auto".
236 * \internal
237 * \todo replace with string or enum class and initialize with sensible value.
239 const char *bonded_opt = nullptr;
241 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
242 int nstlist_cmdline = 0;
243 //! Parameters for replica-exchange simulations.
244 ReplicaExchangeParameters replExParams;
245 //! Print a warning if any force is larger than this (in kJ/mol nm).
246 real pforce = -1;
248 //! Handle to file used for logging.
249 LogFilePtr logFileGuard = nullptr;
250 //! \brief Non-owning handle to file used for logging.
251 t_fileio *logFileHandle = nullptr;
253 //! \brief Non-owning handle to communication data structure.
254 t_commrec *cr = nullptr;
256 //! \brief Non-owning handle to multi-simulation handler.
257 gmx_multisim_t *ms = nullptr;
259 //! Whether the simulation will start afresh, or restart with/without appending.
260 StartingBehavior startingBehavior = StartingBehavior::NewSimulation;
263 * \brief Handle to restraints manager for the current process.
265 * \internal
266 * Use opaque pointer for this implementation detail.
268 std::unique_ptr<RestraintManager> restraintManager_;
271 * \brief Builder for stop signal handler
273 * Optionally provided through MdrunnerBuilder. Client may create a
274 * StopHandlerBuilder and register any number of signal providers before
275 * launching the Mdrunner.
277 * Default is an empty signal handler that will have local signal issuers
278 * added after being passed into the integrator.
280 * \internal
281 * We do not need a full type specification here, so we use an opaque pointer.
283 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_;
284 //! The modules that comprise mdrun.
285 std::unique_ptr<MDModules> mdModules_;
288 /*! \libinternal
289 * \brief Build a gmx::Mdrunner.
291 * Client code (such as `gmx mdrun`) uses this builder to get an initialized Mdrunner.
293 * A builder allows the library to ensure that client code cannot obtain an
294 * uninitialized or partially initialized runner by refusing to build() if the
295 * client has not provided sufficient or self-consistent direction. Director
296 * code can be implemented for different user interfaces, encapsulating any
297 * run-time functionality that does not belong in the library MD code, such
298 * as command-line option processing or interfacing to external libraries.
300 * \ingroup module_mdrun
302 * \internal
304 * The initial Builder implementation is neither extensible at run time nor
305 * at compile time. Future implementations should evolve to compose the runner,
306 * rather than just consolidating the parameters for initialization, but there
307 * is not yet a firm design for how flexibly module code will be coupled to
308 * the builder and how much of the client interface will be in this Builder
309 * versus Builders provided by the various modules.
311 * The named components for the initial builder implementation are descriptive
312 * of the state of mdrun at the time, and are not intended to be prescriptive of
313 * future design.
314 * The probable course of GROMACS development is for the modular components that
315 * support MD simulation to independently express their input parameters (required
316 * and optional) and to provide some sort of help to the UI for input preparation.
317 * If each module provides or aids the instantiation of a Director
318 * for the client code, the Directors could be constructed with a handle to this
319 * Builder and it would not need a public interface.
321 * As the modules are more clearly encapsulated, each module can provide its own
322 * builder, user interface helpers, and/or composable Director code.
323 * The runner and client code will also have to be updated as appropriate
324 * default behavior is clarified for
325 * (a) default behavior of client when user does not provide input,
326 * (b) default behavior of builder when client does not provide input, and
327 * (c) default behavior of runner when builder does not provide input.
329 class MdrunnerBuilder final
331 public:
333 * \brief Constructor requires a handle to a SimulationContext to share.
335 * \param mdModules The handle to the set of modules active in mdrun
336 * \param context The handle to run-time resources and execution environment details.
338 * The calling code must guarantee that the
339 * pointer remains valid for the lifetime of the builder, and that the
340 * resources retrieved from the context remain valid for the lifetime of
341 * the runner produced.
343 * \internal
344 * \todo Find and implement appropriate abstraction layers for SimulationContext.
345 * At some point this parameter should become a constant reference or value
346 * instead of a pointer.
347 * Ref e.g. https://redmine.gromacs.org/issues/2587
349 explicit MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
350 compat::not_null<SimulationContext*> context);
352 //! \cond
353 MdrunnerBuilder() = delete;
354 MdrunnerBuilder(const MdrunnerBuilder&) = delete;
355 MdrunnerBuilder &operator=(const MdrunnerBuilder &) = delete;
356 //! \endcond
358 /*! \brief Allow transfer of ownership with move semantics.
360 * \param builder source object to transfer.
362 * \{
364 MdrunnerBuilder(MdrunnerBuilder && builder) noexcept;
365 MdrunnerBuilder &operator=(MdrunnerBuilder &&builder) noexcept;
366 //! \}
369 * \brief Get ownership of an initialized gmx::Mdrunner.
371 * After build() is called, the Builder object should not be used
372 * again. It is an error to call build without first calling all builder
373 * methods described as "required."
375 * \return A new Mdrunner.
377 * \throws APIError if a required component has not been added before calling build().
379 Mdrunner build();
382 * \brief Set up non-bonded short-range force calculations.
384 * Required. Director code must provide valid options for the non-bonded
385 * interaction code. The builder does not apply any defaults.
387 * \param nbpu_opt Target short-range interactions for "cpu", "gpu", or "auto".
389 * Calling must guarantee that the pointed-to C string is valid through
390 * simulation launch.
392 * \internal
393 * \todo Replace with string or enum that we can have sensible defaults for.
394 * \todo Either the Builder or modular Director code should provide sensible defaults.
396 MdrunnerBuilder &addNonBonded(const char* nbpu_opt);
399 * \brief Set up long-range electrostatics calculations.
401 * Required. Director code should provide valid options for PME electrostatics,
402 * whether or not PME electrostatics are used. The builder does not apply
403 * any defaults, so client code should be prepared to provide (e.g.) "auto"
404 * in the event no user input or logic provides an alternative argument.
406 * \param pme_opt Target long-range interactions for "cpu", "gpu", or "auto".
407 * \param pme_fft_opt Target long-range interactions FFT/solve stages for "cpu", "gpu", or "auto".
409 * Calling must guarantee that the pointed-to C strings are valid through
410 * simulation launch.
412 * \internal
413 * The arguments are passed as references to elements of arrays of C strings.
414 * \todo Replace with modern strings or (better) enum classes.
415 * \todo Make optional and/or encapsulate into electrostatics module.
417 MdrunnerBuilder &addElectrostatics(const char* pme_opt,
418 const char* pme_fft_opt);
421 * \brief Assign responsibility for tasks for bonded interactions.
423 * Required. Director code should provide valid options for
424 * bonded interaction task assignment, whether or not such
425 * interactions are present. The builder does not apply any
426 * defaults, so client code should be prepared to provide
427 * (e.g.) "auto" in the event no user input or logic provides
428 * an alternative argument.
430 * \param bonded_opt Target bonded interactions for "cpu", "gpu", or "auto".
432 * Calling must guarantee that the pointed-to C strings are valid through
433 * simulation launch.
435 * \internal
436 * The arguments are passed as references to elements of arrays of C strings.
437 * \todo Replace with modern strings or (better) enum classes.
438 * \todo Make optional and/or encapsulate into task assignment module.
440 MdrunnerBuilder &addBondedTaskAssignment(const char *bonded_opt);
443 * \brief Provide access to the multisim communicator to use.
445 * \param multisim borrowed handle to multisim record.
447 * Required. Client should get a `gmx_multisim_t*` value from init_multisystem().
448 * Client is responsible for calling done_multisim() on the handle after
449 * simulation.
451 * \internal Pointer is copied, but calling code retains ownership and
452 * responsibility for multisim. Mdrunner must not do anything that would
453 * invalidate the original copied-from pointer.
455 * \todo Clarify semantics of specifying optional multisim work
456 * \todo Clarify ownership and management of multisim resources.
458 MdrunnerBuilder &addMultiSim(gmx_multisim_t* multisim);
461 * \brief Set MD options not owned by some other module.
463 * Optional. Override simulation parameters
465 * \param options structure to copy
466 * \param forceWarningThreshold Print a warning if any force is larger than this (in kJ/mol nm) (default -1)
467 * \param startingBehavior Whether the simulation will start afresh, or restart with/without appending.
469 * \internal
470 * \todo Map these parameters to more appropriate encapsulating types.
471 * Find a better way to indicate "unspecified" than a magic value of the parameter type.
473 MdrunnerBuilder &addSimulationMethod(const MdrunOptions &options,
474 real forceWarningThreshold,
475 StartingBehavior startingBehavior);
478 * \brief Set the domain decomposition module.
480 * Optional. Overrides default constructed DomdecOptions if provided.
482 * \param options options with which to construct domain decomposition.
484 * \internal
485 * \todo revisit whether we should be passing this parameter struct or a higher-level handle of some sort.
487 MdrunnerBuilder &addDomainDecomposition(const DomdecOptions &options);
490 * \brief Set Verlet list manager.
492 * Optional. Neighbor list existence, type, and parameters are mostly determined
493 * by the simulation parameters loaded elsewhere. This is just an override.
495 * \param rebuildInterval override for the duration of a neighbor list with the Verlet scheme.
497 MdrunnerBuilder &addNeighborList(int rebuildInterval);
500 * \brief Set replica exchange manager.
502 * Optional. For guidance on preparing a valid ReplicaExchangeParameters
503 * value, refer to the details in mdrun.cpp, the `t_pargs pa[]` defined there,
504 * and the action of parse_common_args() with regards to that structure.
505 * If not provided by client, a default constructed ReplicaExchangeParameters
506 * is used.
508 * \param params parameters with which to set up replica exchange.
510 * \internal
511 * \todo revisit whether we should be passing this parameter struct or a higher-level handle of some sort.
513 MdrunnerBuilder &addReplicaExchange(const ReplicaExchangeParameters &params);
516 * \brief Specify parameters determining hardware resource allocation.
518 * Optional. If not provided, default-constructed gmx_hw_opt_t will be used.
520 * \param hardwareOptions Parallelism-related user options.
522 MdrunnerBuilder &addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
525 * \brief Provide the filenames options structure with option values chosen
527 * Required. The object is assumed to have been updated by
528 * parse_common_args or equivalent.
530 * \param filenames Filenames and properties from command-line argument values or defaults.
532 * \internal
533 * \todo Modules should manage their own filename options and defaults.
535 MdrunnerBuilder &addFilenames(ArrayRef<const t_filenm> filenames);
538 * \brief Provide parameters for setting up output environment.
540 * Required. Handle is assumed to have been produced by output_env_init
541 * as in parse_common_args.
543 * \param outputEnvironment Output context for writing text files.
545 * \internal
546 * \todo Allow client code to set up output environment and provide as a resource.
547 * This parameter is used to set up resources that are dependent on the execution
548 * environment and API context. Such resources should be retrieved by the simulator
549 * from a client-provided resource, but currently the resources are only fully
550 * initialized in Mdrunner.
552 MdrunnerBuilder &addOutputEnvironment(gmx_output_env_t* outputEnvironment);
555 * \brief Provide the filehandle pointer to be used for the MD log.
557 * Required. Either nullptr if no log should be written, or
558 * valid and open reading for writing.
560 * \param logFileHandle Non-owning handle to file used for logging.
561 * \internal
563 MdrunnerBuilder &addLogFile(t_fileio *logFileHandle);
566 * \brief Provide a StopHandlerBuilder for the MD stop signal handling.
568 * Optional. Defaults to empty.
570 * Client may provide additional (non-default) issuers of simulation stop
571 * signals by preconfiguring the StopHandlerBuilder used later when the
572 * simulation runs.
574 * \param builder
576 MdrunnerBuilder &addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
578 ~MdrunnerBuilder();
580 private:
581 std::unique_ptr<Mdrunner::BuilderImplementation> impl_;
584 } // namespace gmx
586 #endif // GMX_MDRUN_RUNNER_H