OCL: Make variables const
[gromacs.git] / src / gromacs / mdlib / stophandler.cpp
blob39babd5e072b452e3c3b26473ea308bf04b349c5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief Defines StopHandler, a helper class and two stop conditions.
38 * \author Pascal Merz <pascal.merz@colorado.edu>
39 * \ingroup module_mdlib
41 #include "gmxpre.h"
43 #include "stophandler.h"
45 #include "config.h"
47 #include "gromacs/compat/make_unique.h"
48 #include "gromacs/timing/walltime_accounting.h"
49 #include "gromacs/utility/cstringutil.h"
51 namespace gmx
54 StopHandler::StopHandler(
55 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)
63 if (simulationShareState)
65 signal_.isLocal = false;
69 StopConditionSignal::StopConditionSignal(
70 int nstList,
71 bool makeBinaryReproducibleSimulation,
72 int nstSignalComm) :
73 handledStopCondition_(gmx_stop_cond_none),
74 makeBinaryReproducibleSimulation_(makeBinaryReproducibleSimulation),
75 nstSignalComm_(nstSignalComm),
76 nstList_(nstList)
79 StopSignal StopConditionSignal::getSignal(FILE *fplog)
81 StopSignal signal = StopSignal::noSignal;
83 /* Check whether everything is still alright */
84 if (static_cast<int>(gmx_get_stop_condition()) > handledStopCondition_)
86 int nsteps_stop = -1;
88 /* this just makes signals[].sig compatible with the hack
89 of sending signals around by MPI_Reduce together with
90 other floats */
91 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
92 (makeBinaryReproducibleSimulation_ && gmx_get_stop_condition() == gmx_stop_cond_next))
94 /* We need at least two global communication steps to pass
95 * around the signal. We stop at a pair-list creation step
96 * to allow for exact continuation, when possible.
98 signal = StopSignal::stopAtNextNSStep;
99 nsteps_stop = std::max(nstList_, 2 * nstSignalComm_);
101 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
103 /* Stop directly after the next global communication step.
104 * This breaks exact continuation.
106 signal = StopSignal::stopImmediately;
107 nsteps_stop = nstSignalComm_ + 1;
109 if (fplog)
111 fprintf(fplog,
112 "\n\nReceived the %s signal, stopping within %d steps\n\n",
113 gmx_get_signal_name(), nsteps_stop);
114 fflush(fplog);
116 fprintf(stderr,
117 "\n\nReceived the %s signal, stopping within %d steps\n\n",
118 gmx_get_signal_name(), nsteps_stop);
119 fflush(stderr);
120 handledStopCondition_ = static_cast<int>(gmx_get_stop_condition());
123 return signal;
126 StopConditionTime::StopConditionTime(
127 int nstList,
128 real maximumHoursToRun,
129 int nstSignalComm) :
130 signalSent_(false),
131 maximumHoursToRun_(maximumHoursToRun),
132 nstList_(nstList),
133 nstSignalComm_(nstSignalComm),
134 neverUpdateNeighborlist_(nstList <= 0)
137 StopSignal StopConditionTime::getSignal(
138 bool bNS,
139 int64_t step,
140 FILE *fplog,
141 gmx_walltime_accounting_t walltime_accounting)
143 if (signalSent_)
145 // We only want to send it once, but might be called again before run is terminated
146 return StopSignal::noSignal;
148 if ((bNS || neverUpdateNeighborlist_) &&
149 walltime_accounting_get_time_since_start(walltime_accounting) > maximumHoursToRun_ * 60.0 * 60.0 * 0.99)
151 /* Signal to terminate the run */
152 char sbuf[STEPSTRSIZE];
153 int nsteps_stop = std::max(nstList_, 2 * nstSignalComm_);
154 if (fplog)
156 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, "
157 "will terminate the run within %d steps\n",
158 gmx_step_str(step, sbuf), maximumHoursToRun_ * 0.99, nsteps_stop);
160 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, "
161 "will terminate the run within %d steps\n",
162 gmx_step_str(step, sbuf), maximumHoursToRun_ * 0.99, nsteps_stop);
163 signalSent_ = true;
164 return StopSignal::stopAtNextNSStep;
166 return StopSignal::noSignal;
169 void StopHandlerBuilder::registerStopCondition(std::function<StopSignal()> stopCondition)
171 stopConditions_.emplace_back(std::move(stopCondition));
174 std::unique_ptr<StopHandler> StopHandlerBuilder::getStopHandlerMD (
175 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)
188 if (!GMX_THREAD_MPI || isMaster)
190 // TODO: Use unique_ptr once we switch to C++14 (unique_ptr can not easily be
191 // captured in lambda functions in C++11)
192 auto stopConditionSignal = std::make_shared<StopConditionSignal>(
193 nstList, makeBinaryReproducibleSimulation, nstSignalComm);
194 registerStopCondition(
195 [stopConditionSignal, fplog]()
196 {return stopConditionSignal->getSignal(fplog); });
199 if (isMaster && maximumHoursToRun > 0)
201 // TODO: Use unique_ptr once we switch to C++14 (unique_ptr can not easily be
202 // captured in lambda functions in C++11)
203 auto stopConditionTime = std::make_shared<StopConditionTime>(
204 nstList, maximumHoursToRun, nstSignalComm);
205 registerStopCondition(
206 [stopConditionTime, &bNS, &step, fplog, walltime_accounting]()
207 {return stopConditionTime->getSignal(bNS, step, fplog, walltime_accounting); });
210 return compat::make_unique<StopHandler>(
211 signal, simulationShareState, stopConditions_, neverUpdateNeighborList);
214 } // namespace gmx