17bd7df8e357e29b58c49ff6521802a11b5fa93c
[alexxy/gromacs.git] / src / gromacs / mdlib / resethandler.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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
37  * Defines the reset handler class.
38  *
39  * \author Pascal Merz <pascal.merz@colorado.edu>
40  * \ingroup module_mdlib
41  */
42
43 #include "gmxpre.h"
44
45 #include "resethandler.h"
46
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/ewald/pme.h"
49 #include "gromacs/ewald/pme-load-balancing.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/gpu_utils/gpu_utils.h"
52 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
53 #include "gromacs/mdlib/sim_util.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/timing/walltime_accounting.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58
59 namespace gmx
60 {
61
62 /*! \brief Convert signed char (as used by SimulationSignal) to ResetSignal enum
63  *
64  * Expected values are
65  *   \p sig == 0 -- no signal
66  *   \p sig >= 1 -- signal received
67  */
68 static inline ResetSignal convertToResetSignal(signed char sig)
69 {
70     GMX_ASSERT(sig >= 0, "Unexpected reset signal < 0 received");
71     return sig >= 1 ? ResetSignal::doResetCounters : ResetSignal::noSignal;
72 }
73
74 ResetHandler::ResetHandler(
75         compat::not_null<SimulationSignal*> signal,
76         bool                                simulationsShareState,
77         int64_t                             nsteps,
78         bool                                isMaster,
79         bool                                resetHalfway,
80         real                                maximumHoursToRun,
81         const MDLogger                     &mdlog,
82         gmx_wallcycle_t                     wcycle,
83         gmx_walltime_accounting_t           walltime_accounting) :
84     signal_(*signal),
85     rankCanSetSignal_(false),
86     simulationNeedsReset_(false),
87     maximumHoursToRun_(maximumHoursToRun)
88 {
89     if (simulationsShareState)
90     {
91         signal_.isLocal = false;
92     }
93     if (resetHalfway)
94     {
95         GMX_LOG(mdlog.info).asParagraph().
96             appendText(
97                 "The -resethway functionality is deprecated, and may be removed in a future version.");
98         if (nsteps > 0)
99         {
100             /* Signal to reset the counters half the simulation steps. */
101             wcycle_set_reset_counters(wcycle, nsteps / 2);
102         }
103         simulationNeedsReset_ = true;
104
105         if (isMaster && (maximumHoursToRun > 0))
106         {
107             rankCanSetSignal_ = true;
108         }
109     }
110     else if (wcycle_get_reset_counters(wcycle) > 0)
111     {
112         simulationNeedsReset_ = true;
113     }
114     else
115     {
116         // if no reset is happening, this will always be valid
117         walltime_accounting_set_valid_finish(walltime_accounting);
118     }
119 }
120
121 bool ResetHandler::setSignalImpl(gmx_walltime_accounting_t walltime_accounting)
122 {
123     const double secondsSinceStart = walltime_accounting_get_time_since_start(walltime_accounting);
124     if (secondsSinceStart > maximumHoursToRun_ * 60.0 * 60.0 * 0.495)
125     {
126         /* Set flag that will communicate the signal to all ranks in the simulation */
127         signal_.sig = static_cast<signed char>(ResetSignal::doResetCounters);
128         /* Let handler know that we did signal a reset */
129         return true;
130     }
131     /* Let handler know that we did not signal a reset */
132     return false;
133 }
134
135 bool ResetHandler::resetCountersImpl(
136         int64_t                     step,
137         int64_t                     step_rel,
138         const MDLogger             &mdlog,
139         FILE                       *fplog,
140         const t_commrec            *cr,
141         nonbonded_verlet_t         *nbv,
142         t_nrnb                     *nrnb,
143         const gmx_pme_t            *pme,
144         const pme_load_balancing_t *pme_loadbal,
145         gmx_wallcycle_t             wcycle,
146         gmx_walltime_accounting_t   walltime_accounting)
147 {
148     /* Reset either if signal has been passed, or if reset step has been reached */
149     if (convertToResetSignal(signal_.set) == ResetSignal::doResetCounters ||
150         step_rel == wcycle_get_reset_counters(wcycle))
151     {
152         if (pme_loadbal_is_active(pme_loadbal))
153         {
154             /* Do not permit counter reset while PME load
155              * balancing is active. The only purpose for resetting
156              * counters is to measure reliable performance data,
157              * and that can't be done before balancing
158              * completes.
159              *
160              * TODO consider fixing this by delaying the reset
161              * until after load balancing completes,
162              * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
163             gmx_fatal(FARGS, "PME tuning was still active when attempting to "
164                       "reset mdrun counters at step %" PRId64 ". Try "
165                       "resetting counters later in the run, e.g. with gmx "
166                       "mdrun -resetstep.", step);
167         }
168
169         char sbuf[STEPSTRSIZE];
170
171         /* Reset all the counters related to performance over the run */
172         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
173                 "step %s: resetting all time and cycle counters",
174                 gmx_step_str(step, sbuf));
175
176         if (use_GPU(nbv))
177         {
178             nbnxn_gpu_reset_timings(nbv);
179         }
180
181         if (pme_gpu_task_enabled(pme))
182         {
183             pme_gpu_reset_timings(pme);
184         }
185
186         if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
187         {
188             resetGpuProfiler();
189         }
190
191         wallcycle_stop(wcycle, ewcRUN);
192         wallcycle_reset_all(wcycle);
193         if (DOMAINDECOMP(cr))
194         {
195             reset_dd_statistics_counters(cr->dd);
196         }
197         init_nrnb(nrnb);
198         wallcycle_start(wcycle, ewcRUN);
199         walltime_accounting_reset_time(walltime_accounting, step);
200         print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
201
202         wcycle_set_reset_counters(wcycle, -1);
203         if (!thisRankHasDuty(cr, DUTY_PME))
204         {
205             /* Tell our PME node to reset its counters */
206             gmx_pme_send_resetcounters(cr, step);
207         }
208         /* Reset can only happen once, so clear the triggering flag. */
209         signal_.set = static_cast<signed char>(ResetSignal::noSignal);
210         /* We have done a reset, so the finish will be valid. */
211         walltime_accounting_set_valid_finish(walltime_accounting);
212         /* Let handler know that we handled a reset */
213         return true;
214     }
215
216     /* Let handler know that we did not handle a reset */
217     return false;
218 }
219
220 }  // namespace gmx