Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / mdlib / simulationsignal.cpp
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) 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /*! \internal \file
39  *
40  * \brief This file defines functions for inter- and intra-simulation
41  * signalling by mdrun.
42  *
43  * This handles details of responding to termination conditions,
44  * coordinating checkpoints, and coordinating multi-simulations.
45  *
46  * \todo Move this to mdrunutility module alongside gathering
47  * multi-simulation communication infrastructure there.
48  *
49  * \author Berk Hess <hess@kth.se>
50  * \author Mark Abraham <mark.j.abraham@gmail.com>
51  * \ingroup module_mdlib
52  */
53
54 #include "gmxpre.h"
55
56 #include "simulationsignal.h"
57
58 #include <algorithm>
59
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"
66
67 namespace gmx
68 {
69
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_{}
81 {
82 }
83
84 gmx::ArrayRef<real> SimulationSignaller::getCommunicationBuffer()
85 {
86     if (doIntraSim_)
87     {
88         std::transform(std::begin(*signals_),
89                        std::end(*signals_),
90                        std::begin(mpiBuffer_),
91                        [](const SimulationSignals::value_type& s) { return s.sig; });
92
93         return mpiBuffer_;
94     }
95     else
96     {
97         return {};
98     }
99 }
100
101 void SimulationSignaller::signalInterSim()
102 {
103     if (!doInterSim_)
104     {
105         return;
106     }
107     // The situations that lead to doInterSim_ == true without a
108     // multi-simulation begin active should already have issued an
109     // error at mdrun time in release mode, so there's no need for a
110     // release-mode assertion.
111     GMX_ASSERT(isMultiSim(ms_), "Cannot do inter-simulation signalling without a multi-simulation");
112     if (MASTER(cr_))
113     {
114         // Communicate the signals between the simulations.
115         gmx_sum_sim(eglsNR, mpiBuffer_.data(), ms_);
116     }
117     if (DOMAINDECOMP(cr_))
118     {
119         // Communicate the signals from the master to the others.
120         gmx_bcast(eglsNR * sizeof(mpiBuffer_[0]), mpiBuffer_.data(), cr_->mpi_comm_mygroup);
121     }
122 }
123
124 void SimulationSignaller::setSignals()
125 {
126     if (!doIntraSim_)
127     {
128         return;
129     }
130
131     SimulationSignals& s = *signals_;
132     for (size_t i = 0; i < s.size(); i++)
133     {
134         if (doInterSim_ || s[i].isLocal)
135         {
136             // Floating-point reduction of integer values is always
137             // exact, so we can use a simple cast.
138             signed char gsi = static_cast<signed char>(mpiBuffer_[i]);
139
140             /* Set the communicated signal only when it is non-zero,
141              * since a previously set signal might not have been
142              * handled immediately. */
143             if (gsi != 0)
144             {
145                 s[i].set = gsi;
146             }
147             // Turn off any local signal now that it has been processed.
148             s[i].sig = 0;
149         }
150     }
151 }
152
153 void SimulationSignaller::finalizeSignals()
154 {
155     signalInterSim();
156     setSignals();
157 }
158
159 } // namespace gmx