Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / stophandler.h
blobd156b1c6d862a63a4e8de9207b2b2ce0c9750778
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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
36 * \brief Declares StopHandler, a helper class and two stop conditions.
38 * These classes encapsulate the setting and handling of stop signals.
40 * StopHandler lives during the lifetime of do_md. It checks via registered stop
41 * conditions whether the simulation should be stopped at the next possible step or
42 * at the next possible neighbor-searching step. It communicates this via signal to
43 * all ranks and communicates this to do_md via stoppingAfterCurrentStep().
45 * StopHandlerBuilder is owned by the runner, and allows to register stop conditions
46 * at a higher level, outside of do_md. Within do_md, it is creating a StopHandler
47 * object by binding local data and passing a reference to the stop conditions.
49 * Here, we are implementing two stop conditions: StopConditionTime sets a stop condition
50 * based on the elapsed time (only relevant if the -maxh flag was set), while
51 * StopConditionSignal sets stop conditions via signals received from the operating
52 * systems (SIGINT / SIGTERM).
54 * The stop conditions are stored as function pointers created by a lambda expression.
55 * They bind to required local data, in the case of StopConditionTime and StopConditionSignal
56 * these are partially owned by do_md. This requires these function pointers to be deleted
57 * at the end of do_md(). This is achieved by having the do_md() specific function pointers
58 * owned by StopHandler, which in turn is owned (via unique_ptr) by do_md().
60 * \author Pascal Merz <pascal.merz@colorado.edu>
61 * \inlibraryapi
62 * \ingroup module_mdlib
65 #ifndef GMX_MDLIB_STOPHANDLER_H
66 #define GMX_MDLIB_STOPHANDLER_H
68 #include <functional>
69 #include <memory>
70 #include <vector>
72 #include "gromacs/compat/pointers.h"
73 #include "gromacs/mdlib/sighandler.h"
74 #include "gromacs/mdlib/simulationsignal.h"
76 struct gmx_walltime_accounting;
78 namespace gmx
80 /*! \brief Stop signals
82 * Signals that stop conditions can send to all ranks. Possible signals include
83 * * nothing to signal
84 * * stop at the next neighbor-searching step
85 * * stop as soon as signal is received
87 enum class StopSignal
89 noSignal = 0,
90 stopAtNextNSStep = 1,
91 stopImmediately = -1
94 /*! \brief Convert signed char (as used by SimulationSignal) to StopSignal enum
96 * * Expected values are
97 * \p sig == 0 -- no signal
98 * \p sig >= 1 -- stop at next NS
99 * \p sig <= -1 -- stop asap
101 static inline StopSignal convertToStopSignal(signed char sig)
103 if (sig <= -1)
105 return StopSignal::stopImmediately;
107 else if (sig >= 1)
109 return StopSignal::stopAtNextNSStep;
111 else // sig == 0
113 return StopSignal::noSignal;
117 /*! \libinternal
118 * \brief Class handling the stop signal
120 * Loops over the registered stop conditions and sets a signal if
121 * requested (currently only done by master rank).
122 * All ranks receive the stop signal and set the respective flag.
123 * The functions are implemented within this header file to avoid leaving
124 * the translation unit unnecessarily.
126 class StopHandler final
128 public:
129 /*! \brief StopHandler constructor (will be called by StopHandlerBuilder)
131 * @param signal Non-null pointer to a signal used for reading and writing of signals
132 * @param simulationShareState Whether this signal needs to be shared across multiple simulations
133 * @param stopConditions Vector of callback functions setting the signal
134 * @param neverUpdateNeighborList Whether simulation keeps same neighbor list forever
136 * Note: As the StopHandler does not work without this signal, it keeps a non-const reference
137 * to it as a member variable.
139 StopHandler(compat::not_null<SimulationSignal*> signal,
140 bool simulationShareState,
141 std::vector<std::function<StopSignal()>> stopConditions,
142 bool neverUpdateNeighborList);
144 /*! \brief Decides whether a stop signal shall be sent
146 * Loops over the stopCondition vector passed at build time (consisting of conditions
147 * registered with StopHandlerBuilder, and conditions built by StopHandlerBuilder by
148 * default), and sets any signal obtained.
149 * Returns as soon as a StopSignal::stopImmediately signal was obtained, or after
150 * checking all registered stop conditions.
152 void setSignal() const
154 for (auto& condition : stopConditions_)
156 const StopSignal sig = condition();
157 if (sig != StopSignal::noSignal)
159 signal_.sig = static_cast<signed char>(sig);
160 if (sig == StopSignal::stopImmediately)
162 // We don't want this to be overwritten by a less urgent stop
163 break;
169 /*! \brief Decides whether the simulation shall be stopped after the current step
171 * The simulation is stopped after the current step if
172 * * the signal for immediate stop was received, or
173 * * the signal for stop at the next neighbor-searching step was received, and
174 * the current step is a neighbor-searching step.
176 bool stoppingAfterCurrentStep(bool bNS) const
178 return convertToStopSignal(signal_.set) == StopSignal::stopImmediately
179 || (convertToStopSignal(signal_.set) == StopSignal::stopAtNextNSStep
180 && (bNS || neverUpdateNeighborlist_));
183 private:
184 SimulationSignal& signal_;
185 const std::vector<std::function<StopSignal()>> stopConditions_;
186 const bool neverUpdateNeighborlist_;
189 /*! \libinternal
190 * \brief Class setting the stop signal based on gmx_get_stop_condition()
192 * Master rank sets the stop signal if required (generally due to SIGINT).
194 class StopConditionSignal final
196 public:
197 /*! \brief StopConditionSignal constructor
199 StopConditionSignal(int nstList, bool makeBinaryReproducibleSimulation, int nstSignalComm);
201 /*! \brief Decides whether a stopping signal needs to be set
203 * Stop signal is set based on the value of gmx_get_stop_condition(): Set signal for
204 * stop at the next neighbor-searching step at first SIGINT / SIGTERM, set signal
205 * for stop at the next step at second SIGINT / SIGTERM.
207 StopSignal getSignal(FILE* fplog);
209 private:
210 int handledStopCondition_;
211 const bool makeBinaryReproducibleSimulation_;
212 const int nstSignalComm_;
213 const int nstList_;
216 /*! \libinternal
217 * \brief Class setting the stop signal based on maximal run time
219 * Master rank sets the stop signal if run time exceeds maximal run time.
221 class StopConditionTime final
223 public:
224 /*! \brief StopConditionTime constructor
226 StopConditionTime(int nstList, real maximumHoursToRun, int nstSignalComm);
228 /*! \brief Decides whether a stopping signal needs to be set
230 * Stop signal is set if run time is greater than 99% of maximal run time. Signal will
231 * trigger stopping of the simulation at the next neighbor-searching step.
233 StopSignal getSignal(bool bNS, int64_t step, FILE* fplog, gmx_walltime_accounting* walltime_accounting);
235 private:
236 bool signalSent_;
238 const real maximumHoursToRun_;
239 const int nstList_;
240 const int nstSignalComm_;
241 const bool neverUpdateNeighborlist_;
244 /*! \libinternal
245 * \brief Class preparing the creation of a StopHandler
247 * An object of this helper class (owned by the runner) allows to register stop conditions
248 * outside of the actual simulation run via registerStopCondition(), accepting a std::function
249 * object. It then builds a StopHandler object inside do_md, once it can bind to required
250 * local data.
252 * The registered stop conditions plus the standard MD stop conditions (stop based on
253 * received signal from OS [SIGINT / SIGTERM] or based on maximal run time) are then called
254 * by the Stophandler every step to determine whether the simulation should be stopped.
255 * The registered functions need to be of type `std::function<StopSignal()>`, i.e. not
256 * taking any input arguments and returning a `StopSignal` signal which will get propagated
257 * to all ranks. If the function needs input arguments, these need to be bound (e.g. via
258 * lambda capturing) before being registered with the StopHandlerBuilder.
260 class StopHandlerBuilder final
262 public:
263 /*! \brief Register stop condition
265 * This allows code in the scope of the StopHandlerBuilder (runner level) to inject
266 * stop conditions in simulations. Stop conditions are defined as argument-less functions
267 * which return a StopSignal. The return value of this function is then propagated to all
268 * ranks, and allows to stop the simulation at the next global communication step (returned
269 * signal StopSignal::stopImmediately), or at the next NS step (returned signal
270 * StopSignal::stopAtNextNSStep, allows for exact continuation).
272 * Arguments needed by the stop condition function need to be bound / captured. If these
273 * arguments are captured by reference or using a pointer, it is the registrant's
274 * responsibility to ensure that these arguments do not go out of scope during the lifetime
275 * of the StopHandlerBuilder.
277 void registerStopCondition(std::function<StopSignal()> stopCondition);
279 /*! \brief Create StopHandler
281 * Gets called in the scope of the integrator (aka do_md()) to get a pointer to the
282 * StopHandler for the current simulations. Adds the standard MD stop conditions
283 * (e.g. gmx::StopConditionTime, gmx::StopConditionSignal) to the currently registered
284 * stop conditions. Initializes a new StopHandler with this extended vector of
285 * stop conditions. It is the caller's responsibility to make sure arguments passed by
286 * pointer or reference remain valid for the lifetime of the returned StopHandler.
288 std::unique_ptr<StopHandler> getStopHandlerMD(compat::not_null<SimulationSignal*> signal,
289 bool simulationShareState,
290 bool isMaster,
291 int nstList,
292 bool makeBinaryReproducibleSimulation,
293 int nstSignalComm,
294 real maximumHoursToRun,
295 bool neverUpdateNeighborList,
296 FILE* fplog,
297 const int64_t& step,
298 const gmx_bool& bNS,
299 gmx_walltime_accounting* walltime_accounting);
301 private:
302 /*! \brief Initial stopConditions
304 * StopConditions registered via registerStopCondition(). getStopHandlerMD will
305 * copy this vector and add additional conditions before passing the new vector
306 * to the built StopHandler object.
308 std::vector<std::function<StopSignal()>> stopConditions_;
311 } // namespace gmx
313 #endif // GMX_MDLIB_STOPHANDLER_H