Make stepWorkload.useGpuXBufferOps flag consistent
[alexxy/gromacs.git] / src / gromacs / mdlib / stophandler.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief Defines StopHandler, a helper class and two stop conditions.
37  *
38  * \author Pascal Merz <pascal.merz@colorado.edu>
39  * \ingroup module_mdlib
40  */
41 #include "gmxpre.h"
42
43 #include "stophandler.h"
44
45 #include "config.h"
46
47 #include <memory>
48
49 #include "gromacs/timing/walltime_accounting.h"
50 #include "gromacs/utility/cstringutil.h"
51
52 namespace gmx
53 {
54
55 StopHandler::StopHandler(compat::not_null<SimulationSignal*>      signal,
56                          bool                                     simulationShareState,
57                          std::vector<std::function<StopSignal()>> stopConditions,
58                          bool                                     neverUpdateNeighborList) :
59     signal_(*signal),
60     stopConditions_(std::move(stopConditions)),
61     neverUpdateNeighborlist_(neverUpdateNeighborList)
62 {
63     if (simulationShareState)
64     {
65         signal_.isLocal = false;
66     }
67 }
68
69 StopConditionSignal::StopConditionSignal(int nstList, bool makeBinaryReproducibleSimulation, int nstSignalComm) :
70     handledStopCondition_(StopCondition::None),
71     makeBinaryReproducibleSimulation_(makeBinaryReproducibleSimulation),
72     nstSignalComm_(nstSignalComm),
73     nstList_(nstList)
74 {
75 }
76
77 StopSignal StopConditionSignal::getSignal(FILE* fplog)
78 {
79     StopSignal signal = StopSignal::noSignal;
80
81     /* Check whether everything is still alright */
82     if (gmx_get_stop_condition() > handledStopCondition_)
83     {
84         int nsteps_stop = -1;
85
86         /* this just makes signals[].sig compatible with the hack
87            of sending signals around by MPI_Reduce together with
88            other floats */
89         if ((gmx_get_stop_condition() == StopCondition::NextNS)
90             || (makeBinaryReproducibleSimulation_ && gmx_get_stop_condition() == StopCondition::Next))
91         {
92             /* We need at least two global communication steps to pass
93              * around the signal. We stop at a pair-list creation step
94              * to allow for exact continuation, when possible.
95              */
96             signal      = StopSignal::stopAtNextNSStep;
97             nsteps_stop = std::max(nstList_, 2 * nstSignalComm_);
98         }
99         else if (gmx_get_stop_condition() == StopCondition::Next)
100         {
101             /* Stop directly after the next global communication step.
102              * This breaks exact continuation.
103              */
104             signal      = StopSignal::stopImmediately;
105             nsteps_stop = nstSignalComm_ + 1;
106         }
107         if (fplog)
108         {
109             fprintf(fplog,
110                     "\n\nReceived the %s signal, stopping within %d steps\n\n",
111                     gmx_get_signal_name(),
112                     nsteps_stop);
113             fflush(fplog);
114         }
115         fprintf(stderr,
116                 "\n\nReceived the %s signal, stopping within %d steps\n\n",
117                 gmx_get_signal_name(),
118                 nsteps_stop);
119         fflush(stderr);
120         handledStopCondition_ = gmx_get_stop_condition();
121     }
122
123     return signal;
124 }
125
126 StopConditionTime::StopConditionTime(int nstList, real maximumHoursToRun, int nstSignalComm) :
127     signalSent_(false),
128     maximumHoursToRun_(maximumHoursToRun),
129     nstList_(nstList),
130     nstSignalComm_(nstSignalComm),
131     neverUpdateNeighborlist_(nstList <= 0)
132 {
133 }
134
135 StopSignal StopConditionTime::getSignal(bool bNS, int64_t step, FILE* fplog, gmx_walltime_accounting_t walltime_accounting)
136 {
137     if (signalSent_)
138     {
139         // We only want to send it once, but might be called again before run is terminated
140         return StopSignal::noSignal;
141     }
142     if ((bNS || neverUpdateNeighborlist_)
143         && walltime_accounting_get_time_since_start(walltime_accounting)
144                    > maximumHoursToRun_ * 60.0 * 60.0 * 0.99)
145     {
146         /* Signal to terminate the run */
147         char sbuf[STEPSTRSIZE];
148         int  nsteps_stop = std::max(nstList_, 2 * nstSignalComm_);
149         if (fplog)
150         {
151             fprintf(fplog,
152                     "\nStep %s: Run time exceeded %.3f hours, "
153                     "will terminate the run within %d steps\n",
154                     gmx_step_str(step, sbuf),
155                     maximumHoursToRun_ * 0.99,
156                     nsteps_stop);
157         }
158         fprintf(stderr,
159                 "\nStep %s: Run time exceeded %.3f hours, "
160                 "will terminate the run within %d steps\n",
161                 gmx_step_str(step, sbuf),
162                 maximumHoursToRun_ * 0.99,
163                 nsteps_stop);
164         signalSent_ = true;
165         return StopSignal::stopAtNextNSStep;
166     }
167     return StopSignal::noSignal;
168 }
169
170 void StopHandlerBuilder::registerStopCondition(std::function<StopSignal()> stopCondition)
171 {
172     stopConditions_.emplace_back(std::move(stopCondition));
173 };
174
175 std::unique_ptr<StopHandler> StopHandlerBuilder::getStopHandlerMD(compat::not_null<SimulationSignal*> signal,
176                                                                   bool simulationShareState,
177                                                                   bool isMaster,
178                                                                   int  nstList,
179                                                                   bool makeBinaryReproducibleSimulation,
180                                                                   int   nstSignalComm,
181                                                                   real  maximumHoursToRun,
182                                                                   bool  neverUpdateNeighborList,
183                                                                   FILE* fplog,
184                                                                   const int64_t&  step,
185                                                                   const gmx_bool& bNS,
186                                                                   gmx_walltime_accounting_t walltime_accounting)
187 {
188     if (!GMX_THREAD_MPI || isMaster)
189     {
190         // Using shared ptr because move-only callable not supported by std::function.
191         // Would require replacement such as fu2::function or cxx_function.
192         auto stopConditionSignal = std::make_shared<StopConditionSignal>(
193                 nstList, makeBinaryReproducibleSimulation, nstSignalComm);
194         registerStopCondition(
195                 [stopConditionSignal, fplog]() { return stopConditionSignal->getSignal(fplog); });
196     }
197
198     if (isMaster && maximumHoursToRun > 0)
199     {
200         auto stopConditionTime =
201                 std::make_shared<StopConditionTime>(nstList, maximumHoursToRun, nstSignalComm);
202         registerStopCondition([stopConditionTime, &bNS, &step, fplog, walltime_accounting]() {
203             return stopConditionTime->getSignal(bNS, step, fplog, walltime_accounting);
204         });
205     }
206
207     return std::make_unique<StopHandler>(
208             signal, simulationShareState, stopConditions_, neverUpdateNeighborList);
209 }
210
211 } // namespace gmx