Split simulationWork.useGpuBufferOps into separate x and f flags
[alexxy/gromacs.git] / src / gromacs / mdrun / simulationcontext.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 #ifndef GROMACS_SIMULATIONCONTEXT_H
36 #define GROMACS_SIMULATIONCONTEXT_H
37
38 /*! \libinternal
39  * \file
40  * \brief Provide ways for client code to own simulation resources.
41  *
42  * For `gmx mdrun` to be implemented as a client program, public API needs to
43  * provide a way to create and manipulate the SimulationContext.
44  *
45  * \author M. Eric Irrgang <ericirrgang@gmail.com>
46  * \ingroup module_mdrun
47  * \inlibraryapi
48  */
49
50 #include <memory>
51 #include <string>
52
53 #include "gromacs/mdrunutility/multisim.h"
54 #include "gromacs/utility/gmxmpi.h"
55
56 struct gmx_multisim_t;
57
58 namespace gmx
59 {
60
61 template<typename T>
62 class ArrayRef;
63
64 /*! \cond libapi
65  * \libinternal
66  * \brief Simulation environment and configuration.
67  *
68  * SimulationContext allows a simulation module (\relates gmx::mdrun) to retrieve
69  * runtime parameters and resources from client code. The client retains ownership
70  * of the context and its resources, with exceptions as noted.
71  *
72  * A client must share ownership of some resources and transfer ownership of
73  * other resources to create or configure the context. The simulation may in
74  * turn consume or borrow some resources from the context. This is a new
75  * framework that will evolve in the contexts of
76  * https://gitlab.com/gromacs/gromacs/-/issues/2375
77  * https://gitlab.com/gromacs/gromacs/-/issues/2587
78  *
79  * The public interface of SimulationContext is not yet well-specified.
80  * Client code can create an instance with gmx::createSimulationContext()
81  *
82  * \warning The SimulationContext does not yet encapsulate the resources allocated
83  * to the simulator. See https://gitlab.com/gromacs/gromacs/-/issues/3650
84  *
85  * \warning On thread-MPI code paths, the SimulationContext *does not* represent
86  * constraints on the number of simulation ranks, nor does it represent initialized
87  * communication resources.
88  *
89  * \todo Clarify and strengthen the invariant represented by SimulationContext.
90  *
91  * \todo This class should also handle aspects of simulation
92  * environment such as working directory and environment variables.
93  *
94  * \ingroup module_mdrun
95  * \inlibraryapi
96  *
97  * \internal
98  * This is a minimal placeholder for a more complete implementation.
99  * Interfaces for different API levels are not yet final.
100  * \todo Impose sensible access restrictions.
101  * Either the SimulationContext should be passed to the Mdrunner as logically constant or
102  * a separate handle class can provide access to resources that have been
103  * allocated by (negotiated with) the client for the current simulation
104  * (or simulation segment).
105  * Non-const accessors to shared resources need to be associated with update
106  * signals that the simulation components (modules and runner) can subscribe to.
107  *
108  * Also reference https://gitlab.com/gromacs/gromacs/-/issues/2587
109  */
110 class SimulationContext final
111 {
112 public:
113     /*!
114      * \brief Object must be initialized with non-default constructor.
115      */
116     SimulationContext() = delete;
117     /*!
118      * \brief Initialize from borrowed communicator.
119      *
120      * \param communicator Communicator for all collaborating simulation processes.
121      * \param multiSimDirectoryNames  Names of any directories used with -multidir
122      *
123      * Caller is responsible for keeping *communicator* valid for the life of SimulationContext,
124      * and then freeing the communicator, if appropriate.
125      *
126      * With an external MPI library (non-thread-MPI chosen when configuring with CMake),
127      * the client promises that MPI has been initialized (such as by calling gmx::init()).
128      * This communicator is "borrowed" (not duplicated) from the caller.
129      * Additional communicators may be split from the provided communicator
130      * during the life of the SimulationContext or its consumers.
131      *
132      * With thread-MPI, *communicator* must be MPI_COMM_NULL.
133      * The communicator is set up later
134      * during the process of spawning the threads that will be the MPI
135      * ranks. (Multi-simulation is not supported with thread-MPI.)
136      *
137      * \todo Refine distribution of environment management.
138      *       There should be a context level at which only the current simulator directory matters,
139      *       and a level above which encapsulates multisim details in a specialized type.
140      */
141     explicit SimulationContext(MPI_Comm communicator, ArrayRef<const std::string> multiSimDirectoryNames);
142
143     /*!
144      * \brief MPI resources for the entire simulation context.
145      *
146      * With an external MPI library (non-thread-MPI chosen when configuring with CMake),
147      * gmx::init() has called MPI_Init and the provided communicator is valid to use.
148      * The communicator is "borrowed" (not duplicated) from the caller.
149      *
150      * With thread-MPI, the communicator is set up later
151      * during the process of spawning the threads that will be the MPI
152      * ranks.
153      */
154     MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
155
156     /*!
157      * \brief MPI communicator object for this simulation.
158      *
159      * With an external MPI library (non-thread-MPI chosen when configuring with CMake),
160      * gmx::init() has called MPI_Init and the provided communicator is valid to use.
161      * The communicator is "borrowed" (not duplicated) from the caller.
162      *
163      * With thread-MPI, the communicator is set up later
164      * during the process of spawning the threads that will be the MPI
165      * ranks.
166      */
167     MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
168
169     /*!
170      * \brief Multi-sim handler (if required by e.g. gmx mdrun
171      * -multidir; only supported with real MPI)
172      *
173      * If a multi-simulation is in use, then its
174      * communicator(s) are found in multiSimulation_. This
175      * communicator is that of all ranks from all simulations, and
176      * will later be split into one for each simulation.
177      * TODO: Perhaps (for simplicity) communicator splitting
178      *       can be undertaken during multi-sim setup (acquisition of the multisim resource).
179      *
180      * Multi-simulation is not supported with thread-MPI.
181      */
182     std::unique_ptr<gmx_multisim_t> multiSimulation_;
183 };
184 //! \endcond
185
186 } // end namespace gmx
187
188 #endif // GROMACS_SIMULATIONCONTEXT_H