Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / simulationsignal.cpp
blob15f44ccbf4cfadfec6671f09190c109e2c05e9d9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /*! \internal \file
40 * \brief This file defines functions for inter- and intra-simulation
41 * signalling by mdrun.
43 * This handles details of responding to termination conditions,
44 * coordinating checkpoints, and coordinating multi-simulations.
46 * \todo Move this to mdrunutility module alongside gathering
47 * multi-simulation communication infrastructure there.
49 * \author Berk Hess <hess@kth.se>
50 * \author Mark Abraham <mark.j.abraham@gmail.com>
51 * \ingroup module_mdlib
54 #include "gmxpre.h"
56 #include "simulationsignal.h"
58 #include <algorithm>
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/mdrunutility/multisim.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/utility/arrayref.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/real.h"
67 namespace gmx
70 SimulationSignaller::SimulationSignaller(SimulationSignals* signals,
71 const t_commrec* cr,
72 const gmx_multisim_t* ms,
73 bool doInterSim,
74 bool doIntraSim) :
75 signals_(signals),
76 cr_(cr),
77 ms_(ms),
78 doInterSim_(doInterSim),
79 doIntraSim_(doInterSim || doIntraSim),
80 mpiBuffer_{}
84 gmx::ArrayRef<real> SimulationSignaller::getCommunicationBuffer()
86 if (doIntraSim_)
88 std::transform(std::begin(*signals_), std::end(*signals_), std::begin(mpiBuffer_),
89 [](const SimulationSignals::value_type& s) { return s.sig; });
91 return mpiBuffer_;
93 else
95 return {};
99 void SimulationSignaller::signalInterSim()
101 if (!doInterSim_)
103 return;
105 // The situations that lead to doInterSim_ == true without a
106 // multi-simulation begin active should already have issued an
107 // error at mdrun time in release mode, so there's no need for a
108 // release-mode assertion.
109 GMX_ASSERT(isMultiSim(ms_), "Cannot do inter-simulation signalling without a multi-simulation");
110 if (MASTER(cr_))
112 // Communicate the signals between the simulations.
113 gmx_sum_sim(eglsNR, mpiBuffer_.data(), ms_);
115 if (DOMAINDECOMP(cr_))
117 // Communicate the signals from the master to the others.
118 gmx_bcast(eglsNR * sizeof(mpiBuffer_[0]), mpiBuffer_.data(), cr_->mpi_comm_mygroup);
122 void SimulationSignaller::setSignals()
124 if (!doIntraSim_)
126 return;
129 SimulationSignals& s = *signals_;
130 for (size_t i = 0; i < s.size(); i++)
132 if (doInterSim_ || s[i].isLocal)
134 // Floating-point reduction of integer values is always
135 // exact, so we can use a simple cast.
136 signed char gsi = static_cast<signed char>(mpiBuffer_[i]);
138 /* Set the communicated signal only when it is non-zero,
139 * since a previously set signal might not have been
140 * handled immediately. */
141 if (gsi != 0)
143 s[i].set = gsi;
145 // Turn off any local signal now that it has been processed.
146 s[i].sig = 0;
151 void SimulationSignaller::finalizeSignals()
153 signalInterSim();
154 setSignals();
157 } // namespace gmx