1a038a180119397a5da6c1a56d16e941440a6728
[alexxy/gromacs.git] / src / gromacs / mdlib / simulationsignal.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011,2014,2015,2016,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \file
38  *
39  * \brief This file declares functions for inter-rank signalling by mdrun
40  *
41  * This handles details of responding to termination conditions,
42  * coordinating checkpoints, and coordinating multi-simulations.
43  *
44  * \todo Move this to mdrunutility module alongside gathering
45  * multi-simulation communication infrastructure there.
46  *
47  * \author Berk Hess <hess@kth.se>
48  * \author Mark Abraham <mark.j.abraham@gmail.com>
49  * \inlibraryapi
50  * \ingroup module_mdlib
51  */
52 #ifndef GMX_MDLIB_SIMULATIONSIGNAL_H
53 #define GMX_MDLIB_SIMULATIONSIGNAL_H
54
55 #include <array>
56
57 #include "gromacs/utility/real.h"
58
59 struct gmx_multisim_t;
60 struct t_commrec;
61
62 //! Kinds of simulation conditions to signal about.
63 enum
64 {
65     eglsCHKPT,
66     eglsSTOPCOND,
67     eglsRESETCOUNTERS,
68     eglsNR
69 };
70
71 namespace gmx
72 {
73
74 template<typename T>
75 class ArrayRef;
76
77 /*!
78  * \brief POD-style object used by mdrun ranks to set and
79  * receive signals within and between simulations.
80  *
81  * Keep in mind that the values of signals are transmitted to other
82  * ranks through an MPI_Reduce after casting them to a real (so the
83  * signals can be sent together with other data). This means that the
84  * only meaningful values are positive, negative or zero.
85  *
86  * isLocal permits (for example) replica-exchange to require that any
87  * checkpointing is synchronized across all simulations, by setting
88  * isLocal to false, so that the trigger for action is set only when
89  * inter-simulation signalling happens. Replica-exchange can
90  * coordinate this at run time when a SimulationSignaller is made. */
91 class SimulationSignal
92 {
93 public:
94     //! Constructor
95     SimulationSignal(bool isSignalLocal = true) : sig(0), set(0), isLocal(isSignalLocal) {}
96     //! The signal set by this rank in do_md().
97     signed char sig;
98     //! The communicated signal that triggers action, which will be equal for all ranks, once communication has occured.
99     signed char set;
100     //! Is the signal in one simulation independent of other simulations?
101     bool isLocal;
102 };
103
104 //! Convenience typedef for the group of signals used.
105 typedef std::array<SimulationSignal, eglsNR> SimulationSignals;
106
107 /*!
108  * \brief Object used by mdrun ranks to signal to each other at this step.
109  *
110  * This object has responsibility to read signal values from \c gs,
111  * coordinate communication within and perhaps between simulations,
112  * and set result signal values in \c gs as appropriate.
113  *
114  * It is intended to have a very short lifetime, so should remain easy
115  * to construct and destruct on the stack just when the global
116  * communication occurs. */
117 class SimulationSignaller
118 {
119 public:
120     //! Constructor
121     SimulationSignaller(SimulationSignals*    signals,
122                         const t_commrec*      cr,
123                         const gmx_multisim_t* ms,
124                         bool                  doInterSim,
125                         bool                  doIntraSim);
126     /*! \brief Return a reference to an array of signal values to communicate.
127      *
128      * \return If intra-sim signalling will take place, fill and
129      * return a reference to the array of reals in which signals
130      * will be communicated with the signal values to be
131      * sent. Otherwise return a EmptyArrayRef. */
132     gmx::ArrayRef<real> getCommunicationBuffer();
133     /*! \brief Handle inter-simulation signal communication.
134      *
135      * If an inter-simulation signal should be handled, communicate between
136      * simulation-master ranks, then propagate from the masters to the
137      * rest of the ranks for each simulation. It is the responsibility of
138      * the calling code to ensure that any necessary intra-simulation
139      * signalling has already occurred, e.g. in global_stat(). */
140     void signalInterSim();
141     /*! \brief Propagate signals when appropriate.
142      *
143      * Always propagate an mdrun signal value when doing
144      * inter-simulation signalling; otherwise, propagate it only
145      * if should be propagated within this simulation,
146      * ie. locally. See documentation of SimulationSignal for
147      * details. */
148     void setSignals();
149     //! Convenience wrapper that calls signalInterSim() then setSignals().
150     void finalizeSignals();
151
152 private:
153     //! Source and sink for mdrun signals
154     SimulationSignals* signals_;
155     //! Communication object.
156     const t_commrec* cr_;
157     //! Multi-sim handler.
158     const gmx_multisim_t* ms_;
159     //! Do inter-sim communication at this step.
160     bool doInterSim_;
161     //! Do intra-sim communication at this step.
162     bool doIntraSim_;
163     //! Buffer for MPI communication.
164     std::array<real, eglsNR> mpiBuffer_;
165 };
166
167 } // namespace gmx
168
169 #endif