Merge branch 'origin/release-5-1'
[alexxy/gromacs.git] / src / gromacs / mdlib / mdrun_signalling.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, 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 /*! \internal \file
38  *
39  * \brief This file defines 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  * \author Berk Hess <hess@kth.se>
45  * \author Mark Abraham <mark.j.abraham@gmail.com>
46  * \ingroup module_mdlib
47  */
48
49 #include "gmxpre.h"
50
51 #include "mdrun_signalling.h"
52
53 #include <algorithm>
54
55 #include "gromacs/gmxlib/network.h"
56 #include "gromacs/mdlib/md_support.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/real.h"
62
63 void init_global_signals(struct gmx_signalling_t *gs, const t_commrec *cr,
64                          const t_inputrec *ir, int repl_ex_nst)
65 {
66     int i;
67
68     if (MULTISIM(cr))
69     {
70         gs->nstms = multisim_nstsimsync(cr, ir, repl_ex_nst);
71         if (debug)
72         {
73             fprintf(debug, "Syncing simulations for checkpointing and termination every %d steps\n", gs->nstms);
74         }
75     }
76     else
77     {
78         gs->nstms = 1;
79     }
80
81     for (i = 0; i < eglsNR; i++)
82     {
83         gs->sig[i] = 0;
84         gs->set[i] = 0;
85     }
86 }
87
88 gmx::ArrayRef<real>
89 prepareSignalBuffer(struct gmx_signalling_t *gs)
90 {
91     if (gs)
92     {
93         gmx::ArrayRef<signed char> sig(gs->sig);
94         gmx::ArrayRef<real>        temp(gs->mpiBuffer);
95
96         std::copy(sig.begin(), sig.end(), temp.begin());
97
98         return temp;
99     }
100     else
101     {
102         return gmx::EmptyArrayRef();
103     }
104 }
105
106 void
107 handleSignals(struct gmx_signalling_t  *gs,
108               const t_commrec          *cr,
109               bool                      bInterSimGS)
110 {
111     /* Is the signal in one simulation independent of other simulations? */
112     bool bIsSignalLocal[eglsNR] = { false, false, true };
113
114     if (!gs)
115     {
116         return;
117     }
118
119     if (MULTISIM(cr) && bInterSimGS)
120     {
121         if (MASTER(cr))
122         {
123             /* Communicate the signals between the simulations */
124             gmx_sum_sim(eglsNR, gs->mpiBuffer, cr->ms);
125         }
126         /* Communicate the signals from the master to the others */
127         gmx_bcast(eglsNR*sizeof(gs->mpiBuffer[0]), gs->mpiBuffer, cr);
128     }
129     for (int i = 0; i < eglsNR; i++)
130     {
131         if (bInterSimGS || bIsSignalLocal[i])
132         {
133             /* Set the communicated signal only when it is non-zero,
134              * since signals might not be processed at each MD step.
135              */
136             signed char gsi = (gs->mpiBuffer[i] >= 0.0 ?
137                                static_cast<signed char>(gs->mpiBuffer[i] + 0.5) :
138                                static_cast<signed char>(gs->mpiBuffer[i] - 0.5));
139             if (gsi != 0)
140             {
141                 gs->set[i] = gsi;
142             }
143             /* Turn off the local signal */
144             gs->sig[i] = 0;
145         }
146     }
147 }