Replace compat::make_unique with std::make_unique
[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, 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(
56         compat::not_null<SimulationSignal*>         signal,
57         bool                                        simulationShareState,
58         std::vector < std::function<StopSignal()> > stopConditions,
59         bool                                        neverUpdateNeighborList) :
60     signal_(*signal),
61     stopConditions_(std::move(stopConditions)),
62     neverUpdateNeighborlist_(neverUpdateNeighborList)
63 {
64     if (simulationShareState)
65     {
66         signal_.isLocal = false;
67     }
68 }
69
70 StopConditionSignal::StopConditionSignal(
71         int  nstList,
72         bool makeBinaryReproducibleSimulation,
73         int  nstSignalComm) :
74     handledStopCondition_(gmx_stop_cond_none),
75     makeBinaryReproducibleSimulation_(makeBinaryReproducibleSimulation),
76     nstSignalComm_(nstSignalComm),
77     nstList_(nstList)
78 {}
79
80 StopSignal StopConditionSignal::getSignal(FILE *fplog)
81 {
82     StopSignal signal = StopSignal::noSignal;
83
84     /* Check whether everything is still alright */
85     if (static_cast<int>(gmx_get_stop_condition()) > handledStopCondition_)
86     {
87         int nsteps_stop   = -1;
88
89         /* this just makes signals[].sig compatible with the hack
90            of sending signals around by MPI_Reduce together with
91            other floats */
92         if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
93             (makeBinaryReproducibleSimulation_ && gmx_get_stop_condition() == gmx_stop_cond_next))
94         {
95             /* We need at least two global communication steps to pass
96              * around the signal. We stop at a pair-list creation step
97              * to allow for exact continuation, when possible.
98              */
99             signal      = StopSignal::stopAtNextNSStep;
100             nsteps_stop = std::max(nstList_, 2 * nstSignalComm_);
101         }
102         else if (gmx_get_stop_condition() == gmx_stop_cond_next)
103         {
104             /* Stop directly after the next global communication step.
105              * This breaks exact continuation.
106              */
107             signal      = StopSignal::stopImmediately;
108             nsteps_stop = nstSignalComm_ + 1;
109         }
110         if (fplog)
111         {
112             fprintf(fplog,
113                     "\n\nReceived the %s signal, stopping within %d steps\n\n",
114                     gmx_get_signal_name(), nsteps_stop);
115             fflush(fplog);
116         }
117         fprintf(stderr,
118                 "\n\nReceived the %s signal, stopping within %d steps\n\n",
119                 gmx_get_signal_name(), nsteps_stop);
120         fflush(stderr);
121         handledStopCondition_ = static_cast<int>(gmx_get_stop_condition());
122     }
123
124     return signal;
125 }
126
127 StopConditionTime::StopConditionTime(
128         int  nstList,
129         real maximumHoursToRun,
130         int  nstSignalComm) :
131     signalSent_(false),
132     maximumHoursToRun_(maximumHoursToRun),
133     nstList_(nstList),
134     nstSignalComm_(nstSignalComm),
135     neverUpdateNeighborlist_(nstList <= 0)
136 {}
137
138 StopSignal StopConditionTime::getSignal(
139         bool                      bNS,
140         int64_t                   step,
141         FILE                     *fplog,
142         gmx_walltime_accounting_t walltime_accounting)
143 {
144     if (signalSent_)
145     {
146         // We only want to send it once, but might be called again before run is terminated
147         return StopSignal::noSignal;
148     }
149     if ((bNS || neverUpdateNeighborlist_) &&
150         walltime_accounting_get_time_since_start(walltime_accounting) > maximumHoursToRun_ * 60.0 * 60.0 * 0.99)
151     {
152         /* Signal to terminate the run */
153         char sbuf[STEPSTRSIZE];
154         int  nsteps_stop = std::max(nstList_, 2 * nstSignalComm_);
155         if (fplog)
156         {
157             fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, "
158                     "will terminate the run within %d steps\n",
159                     gmx_step_str(step, sbuf), maximumHoursToRun_ * 0.99, nsteps_stop);
160         }
161         fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, "
162                 "will terminate the run within %d steps\n",
163                 gmx_step_str(step, sbuf), maximumHoursToRun_ * 0.99, 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 (
176         compat::not_null<SimulationSignal*> signal,
177         bool                                simulationShareState,
178         bool                                isMaster,
179         int                                 nstList,
180         bool                                makeBinaryReproducibleSimulation,
181         int                                 nstSignalComm,
182         real                                maximumHoursToRun,
183         bool                                neverUpdateNeighborList,
184         FILE                               *fplog,
185         const int64_t                      &step,
186         const gmx_bool                     &bNS,
187         gmx_walltime_accounting_t           walltime_accounting)
188 {
189     if (!GMX_THREAD_MPI || isMaster)
190     {
191         // Using shared ptr because move-only callable not supported by std::function.
192         // Would require replacement such as fu2::function or cxx_function.
193         auto stopConditionSignal = std::make_shared<StopConditionSignal>(
194                     nstList, makeBinaryReproducibleSimulation, nstSignalComm);
195         registerStopCondition(
196                 [stopConditionSignal, fplog]()
197                 {return stopConditionSignal->getSignal(fplog); });
198     }
199
200     if (isMaster && maximumHoursToRun > 0)
201     {
202         auto stopConditionTime = std::make_shared<StopConditionTime>(
203                     nstList, maximumHoursToRun, nstSignalComm);
204         registerStopCondition(
205                 [stopConditionTime, &bNS, &step, fplog, walltime_accounting]()
206                 {return stopConditionTime->getSignal(bNS, step, fplog, walltime_accounting); });
207     }
208
209     return std::make_unique<StopHandler>(
210             signal, simulationShareState, stopConditions_, neverUpdateNeighborList);
211 }
212
213 } // namespace gmx