Fixes for clang-tidy-9
[alexxy/gromacs.git] / src / gromacs / mdrun / rerun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020, 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  *
37  * \brief Implements the loop for simulation reruns
38  *
39  * \author Pascal Merz <pascal.merz@colorado.edu>
40  * \ingroup module_mdrun
41  */
42 #include "gmxpre.h"
43
44 #include <cinttypes>
45 #include <cmath>
46 #include <cstdio>
47 #include <cstdlib>
48
49 #include <algorithm>
50 #include <memory>
51
52 #include "gromacs/awh/awh.h"
53 #include "gromacs/commandline/filenm.h"
54 #include "gromacs/domdec/collect.h"
55 #include "gromacs/domdec/dlbtiming.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/mdsetup.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/ewald/pme_load_balancing.h"
63 #include "gromacs/ewald/pme_pp.h"
64 #include "gromacs/fileio/trxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/listed_forces/manage_threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/mdlib/checkpointhandler.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/enerdata_utils.h"
78 #include "gromacs/mdlib/energyoutput.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdoutf.h"
86 #include "gromacs/mdlib/membed.h"
87 #include "gromacs/mdlib/resethandler.h"
88 #include "gromacs/mdlib/sighandler.h"
89 #include "gromacs/mdlib/simulationsignal.h"
90 #include "gromacs/mdlib/stat.h"
91 #include "gromacs/mdlib/stophandler.h"
92 #include "gromacs/mdlib/tgroup.h"
93 #include "gromacs/mdlib/trajectory_writing.h"
94 #include "gromacs/mdlib/update.h"
95 #include "gromacs/mdlib/vcm.h"
96 #include "gromacs/mdlib/vsite.h"
97 #include "gromacs/mdrunutility/handlerestart.h"
98 #include "gromacs/mdrunutility/multisim.h"
99 #include "gromacs/mdrunutility/printtime.h"
100 #include "gromacs/mdtypes/awh_history.h"
101 #include "gromacs/mdtypes/awh_params.h"
102 #include "gromacs/mdtypes/commrec.h"
103 #include "gromacs/mdtypes/df_history.h"
104 #include "gromacs/mdtypes/energyhistory.h"
105 #include "gromacs/mdtypes/fcdata.h"
106 #include "gromacs/mdtypes/forcerec.h"
107 #include "gromacs/mdtypes/group.h"
108 #include "gromacs/mdtypes/inputrec.h"
109 #include "gromacs/mdtypes/interaction_const.h"
110 #include "gromacs/mdtypes/md_enums.h"
111 #include "gromacs/mdtypes/mdatom.h"
112 #include "gromacs/mdtypes/mdrunoptions.h"
113 #include "gromacs/mdtypes/observableshistory.h"
114 #include "gromacs/mdtypes/state.h"
115 #include "gromacs/mimic/utilities.h"
116 #include "gromacs/pbcutil/pbc.h"
117 #include "gromacs/pulling/pull.h"
118 #include "gromacs/swap/swapcoords.h"
119 #include "gromacs/timing/wallcycle.h"
120 #include "gromacs/timing/walltime_accounting.h"
121 #include "gromacs/topology/atoms.h"
122 #include "gromacs/topology/idef.h"
123 #include "gromacs/topology/mtop_util.h"
124 #include "gromacs/topology/topology.h"
125 #include "gromacs/trajectory/trajectoryframe.h"
126 #include "gromacs/utility/basedefinitions.h"
127 #include "gromacs/utility/cstringutil.h"
128 #include "gromacs/utility/fatalerror.h"
129 #include "gromacs/utility/logger.h"
130 #include "gromacs/utility/real.h"
131 #include "gromacs/utility/smalloc.h"
132
133 #include "legacysimulator.h"
134 #include "replicaexchange.h"
135 #include "shellfc.h"
136
137 using gmx::SimulationSignaller;
138
139 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
140  *
141  * \param[in]     rerunFrame      The trajectory frame to compute energy/forces for
142  * \param[in,out] globalState     The global state container
143  * \param[in]     constructVsites When true, vsite coordinates are constructed
144  * \param[in]     vsite           Vsite setup, can be nullptr when \p constructVsites = false
145  * \param[in]     idef            Topology parameters, used for constructing vsites
146  * \param[in]     timeStep        Time step, used for constructing vsites
147  * \param[in]     forceRec        Force record, used for constructing vsites
148  */
149 static void prepareRerunState(const t_trxframe&             rerunFrame,
150                               t_state*                      globalState,
151                               bool                          constructVsites,
152                               const gmx_vsite_t*            vsite,
153                               const InteractionDefinitions& idef,
154                               double                        timeStep,
155                               const t_forcerec&             forceRec)
156 {
157     auto x      = makeArrayRef(globalState->x);
158     auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
159     std::copy(rerunX.begin(), rerunX.end(), x.begin());
160     copy_mat(rerunFrame.box, globalState->box);
161
162     if (constructVsites)
163     {
164         GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
165
166         construct_vsites(vsite, globalState->x.rvec_array(), timeStep, globalState->v.rvec_array(),
167                          idef.iparams, idef.il, forceRec.pbcType, forceRec.bMolPBC, nullptr,
168                          globalState->box);
169     }
170 }
171
172 void gmx::LegacySimulator::do_rerun()
173 {
174     // TODO Historically, the EM and MD "integrators" used different
175     // names for the t_inputrec *parameter, but these must have the
176     // same name, now that it's a member of a struct. We use this ir
177     // alias to avoid a large ripple of nearly useless changes.
178     // t_inputrec is being replaced by IMdpOptionsProvider, so this
179     // will go away eventually.
180     t_inputrec*                 ir = inputrec;
181     int64_t                     step, step_rel;
182     double                      t, lam0[efptNR];
183     bool                        isLastStep               = false;
184     bool                        doFreeEnergyPerturbation = false;
185     unsigned int                force_flags;
186     tensor                      force_vir, shake_vir, total_vir, pres;
187     t_trxstatus*                status = nullptr;
188     rvec                        mu_tot;
189     t_trxframe                  rerun_fr;
190     gmx_localtop_t              top(top_global->ffparams);
191     PaddedHostVector<gmx::RVec> f{};
192     gmx_global_stat_t           gstat;
193     gmx_shellfc_t*              shellfc;
194
195     double cycles;
196
197     /* Domain decomposition could incorrectly miss a bonded
198        interaction, but checking for that requires a global
199        communication stage, which does not otherwise happen in DD
200        code. So we do that alongside the first global energy reduction
201        after a new DD is made. These variables handle whether the
202        check happens, and the result it returns. */
203     bool shouldCheckNumberOfBondedInteractions = false;
204     int  totalNumberOfBondedInteractions       = -1;
205
206     SimulationSignals signals;
207     // Most global communnication stages don't propagate mdrun
208     // signals, and will use this object to achieve that.
209     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
210
211     GMX_LOG(mdlog.info)
212             .asParagraph()
213             .appendText(
214                     "Note that it is planned that the command gmx mdrun -rerun will "
215                     "be available in a different form in a future version of GROMACS, "
216                     "e.g. gmx rerun -f.");
217
218     if (ir->efep != efepNO
219         && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
220     {
221         gmx_fatal(FARGS,
222                   "Perturbed masses or constraints are not supported by rerun. "
223                   "Either make a .tpr without mass and constraint perturbation, "
224                   "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
225     }
226     if (ir->bExpanded)
227     {
228         gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
229     }
230     if (ir->bSimTemp)
231     {
232         gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
233     }
234     if (ir->bDoAwh)
235     {
236         gmx_fatal(FARGS, "AWH not supported by rerun.");
237     }
238     if (replExParams.exchangeInterval > 0)
239     {
240         gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
241     }
242     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
243     {
244         gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
245     }
246     if (ir->bIMD)
247     {
248         gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
249     }
250     if (isMultiSim(ms))
251     {
252         gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
253     }
254     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
255                     [](int i) { return i != eannNO; }))
256     {
257         gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
258     }
259
260     /* Rerun can't work if an output file name is the same as the input file name.
261      * If this is the case, the user will get an error telling them what the issue is.
262      */
263     if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
264         || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
265     {
266         gmx_fatal(FARGS,
267                   "When using mdrun -rerun, the name of the input trajectory file "
268                   "%s cannot be identical to the name of an output file (whether "
269                   "given explicitly with -o or -x, or by default)",
270                   opt2fn("-rerun", nfile, fnm));
271     }
272
273     /* Settings for rerun */
274     ir->nstlist              = 1;
275     ir->nstcalcenergy        = 1;
276     int        nstglobalcomm = 1;
277     const bool bNS           = true;
278
279     ir->nstxout_compressed         = 0;
280     const SimulationGroups* groups = &top_global->groups;
281     if (ir->eI == eiMimic)
282     {
283         top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
284     }
285
286     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
287     const bool        simulationsShareState = false;
288     gmx_mdoutf*       outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
289                                    mdModulesNotifier, ir, top_global, oenv, wcycle,
290                                    StartingBehavior::NewSimulation, simulationsShareState, ms);
291     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work,
292                                    mdoutf_get_fp_dhdl(outf), true, StartingBehavior::NewSimulation,
293                                    mdModulesNotifier);
294
295     gstat = global_stat_init(ir);
296
297     /* Check for polarizable models and flexible constraints */
298     shellfc = init_shell_flexcon(fplog, top_global, constr ? constr->numFlexibleConstraints() : 0,
299                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
300
301     {
302         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
303         if ((io > 2000) && MASTER(cr))
304         {
305             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
306         }
307     }
308
309     // Local state only becomes valid now.
310     std::unique_ptr<t_state> stateInstance;
311     t_state*                 state;
312
313     if (DOMAINDECOMP(cr))
314     {
315         stateInstance = std::make_unique<t_state>();
316         state         = stateInstance.get();
317         dd_init_local_state(cr->dd, state_global, state);
318
319         /* Distribute the charge groups over the nodes from the master node */
320         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
321                             imdSession, pull_work, state, &f, mdAtoms, &top, fr, vsite, constr,
322                             nrnb, nullptr, FALSE);
323         shouldCheckNumberOfBondedInteractions = true;
324     }
325     else
326     {
327         state_change_natoms(state_global, state_global->natoms);
328         /* Copy the pointer to the global state */
329         state = state_global;
330
331         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
332     }
333
334     auto mdatoms = mdAtoms->mdatoms();
335
336     // NOTE: The global state is no longer used at this point.
337     // But state_global is still used as temporary storage space for writing
338     // the global state to file and potentially for replica exchange.
339     // (Global topology should persist.)
340
341     update_mdatoms(mdatoms, state->lambda[efptMASS]);
342
343     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
344     {
345         doFreeEnergyPerturbation = true;
346     }
347
348     {
349         int cglo_flags =
350                 (CGLO_GSTAT
351                  | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
352         bool   bSumEkinhOld = false;
353         t_vcm* vcm          = nullptr;
354         compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
355                         makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms, nrnb,
356                         vcm, nullptr, enerd, force_vir, shake_vir, total_vir, pres, constr, &nullSignaller,
357                         state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
358     }
359     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global, &top,
360                                     makeConstArrayRef(state->x), state->box,
361                                     &shouldCheckNumberOfBondedInteractions);
362
363     if (MASTER(cr))
364     {
365         fprintf(stderr,
366                 "starting md rerun '%s', reading coordinates from"
367                 " input trajectory '%s'\n\n",
368                 *(top_global->name), opt2fn("-rerun", nfile, fnm));
369         if (mdrunOptions.verbose)
370         {
371             fprintf(stderr,
372                     "Calculated time to finish depends on nsteps from "
373                     "run input file,\nwhich may not correspond to the time "
374                     "needed to process input trajectory.\n\n");
375         }
376         fprintf(fplog, "\n");
377     }
378
379     walltime_accounting_start_time(walltime_accounting);
380     wallcycle_start(wcycle, ewcRUN);
381     print_start(fplog, cr, walltime_accounting, "mdrun");
382
383     /***********************************************************
384      *
385      *             Loop over MD steps
386      *
387      ************************************************************/
388
389     if (constr)
390     {
391         GMX_LOG(mdlog.info)
392                 .asParagraph()
393                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
394     }
395
396     rerun_fr.natoms = 0;
397     if (MASTER(cr))
398     {
399         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
400         if (rerun_fr.natoms != top_global->natoms)
401         {
402             gmx_fatal(FARGS,
403                       "Number of atoms in trajectory (%d) does not match the "
404                       "run input file (%d)\n",
405                       rerun_fr.natoms, top_global->natoms);
406         }
407
408         if (ir->pbcType != PbcType::No)
409         {
410             if (!rerun_fr.bBox)
411             {
412                 gmx_fatal(FARGS,
413                           "Rerun trajectory frame step %" PRId64
414                           " time %f "
415                           "does not contain a box, while pbc is used",
416                           rerun_fr.step, rerun_fr.time);
417             }
418             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
419             {
420                 gmx_fatal(FARGS,
421                           "Rerun trajectory frame step %" PRId64
422                           " time %f "
423                           "has too small box dimensions",
424                           rerun_fr.step, rerun_fr.time);
425             }
426         }
427     }
428
429     GMX_LOG(mdlog.info)
430             .asParagraph()
431             .appendText(
432                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
433                     "pressure.");
434
435     if (PAR(cr))
436     {
437         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
438     }
439
440     if (ir->pbcType != PbcType::No)
441     {
442         /* Set the shift vectors.
443          * Necessary here when have a static box different from the tpr box.
444          */
445         calc_shifts(rerun_fr.box, fr->shift_vec);
446     }
447
448     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
449             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false, MASTER(cr),
450             ir->nstlist, mdrunOptions.reproducible, nstglobalcomm, mdrunOptions.maximumHoursToRun,
451             ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
452
453     // we don't do counter resetting in rerun - finish will always be valid
454     walltime_accounting_set_valid_finish(walltime_accounting);
455
456     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
457
458     step     = ir->init_step;
459     step_rel = 0;
460
461     /* and stop now if we should */
462     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
463     while (!isLastStep)
464     {
465         wallcycle_start(wcycle, ewcSTEP);
466
467         if (rerun_fr.bStep)
468         {
469             step     = rerun_fr.step;
470             step_rel = step - ir->init_step;
471         }
472         if (rerun_fr.bTime)
473         {
474             t = rerun_fr.time;
475         }
476         else
477         {
478             t = step;
479         }
480
481         if (ir->efep != efepNO && MASTER(cr))
482         {
483             setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
484         }
485
486         if (MASTER(cr))
487         {
488             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
489             if (constructVsites && DOMAINDECOMP(cr))
490             {
491                 gmx_fatal(FARGS,
492                           "Vsite recalculation with -rerun is not implemented with domain "
493                           "decomposition, "
494                           "use a single rank");
495             }
496             prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top.idef, ir->delta_t, *fr);
497         }
498
499         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
500
501         if (DOMAINDECOMP(cr))
502         {
503             /* Repartition the domain decomposition */
504             const bool bMasterState = true;
505             dd_partition_system(fplog, mdlog, step, cr, bMasterState, nstglobalcomm, state_global,
506                                 *top_global, ir, imdSession, pull_work, state, &f, mdAtoms, &top,
507                                 fr, vsite, constr, nrnb, wcycle, mdrunOptions.verbose);
508             shouldCheckNumberOfBondedInteractions = true;
509         }
510
511         if (MASTER(cr))
512         {
513             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
514         }
515
516         if (ir->efep != efepNO)
517         {
518             update_mdatoms(mdatoms, state->lambda[efptMASS]);
519         }
520
521         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
522                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
523                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
524
525         if (shellfc)
526         {
527             /* Now is the time to relax the shells */
528             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, enforcedRotation, step, ir,
529                                 imdSession, pull_work, bNS, force_flags, &top, constr, enerd, fcd,
530                                 state->natoms, state->x.arrayRefWithPadding(),
531                                 state->v.arrayRefWithPadding(), state->box, state->lambda, &state->hist,
532                                 f.arrayRefWithPadding(), force_vir, mdatoms, nrnb, wcycle, shellfc,
533                                 fr, runScheduleWork, t, mu_tot, vsite, ddBalanceRegionHandler);
534         }
535         else
536         {
537             /* The coordinates (x) are shifted (to get whole molecules)
538              * in do_force.
539              * This is parallellized as well, and does communication too.
540              * Check comments in sim_util.c
541              */
542             Awh*       awh = nullptr;
543             gmx_edsam* ed  = nullptr;
544             do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession, pull_work, step, nrnb,
545                      wcycle, &top, state->box, state->x.arrayRefWithPadding(), &state->hist,
546                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd, state->lambda, fr,
547                      runScheduleWork, vsite, mu_tot, t, ed, GMX_FORCE_NS | force_flags,
548                      ddBalanceRegionHandler);
549         }
550
551         /* Now we have the energies and forces corresponding to the
552          * coordinates at time t.
553          */
554         {
555             const bool isCheckpointingStep = false;
556             const bool doRerun             = true;
557             const bool bSumEkinhOld        = false;
558             do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t, ir, state,
559                                      state_global, observablesHistory, top_global, fr, outf,
560                                      energyOutput, ekind, f, isCheckpointingStep, doRerun,
561                                      isLastStep, mdrunOptions.writeConfout, bSumEkinhOld);
562         }
563
564         stopHandler->setSignal();
565
566         if (vsite != nullptr)
567         {
568             wallcycle_start(wcycle, ewcVSITECONSTR);
569             construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t,
570                              as_rvec_array(state->v.data()), top.idef.iparams, top.idef.il,
571                              fr->pbcType, fr->bMolPBC, cr, state->box);
572             wallcycle_stop(wcycle, ewcVSITECONSTR);
573         }
574
575         {
576             const bool          doInterSimSignal = false;
577             const bool          doIntraSimSignal = true;
578             bool                bSumEkinhOld     = false;
579             t_vcm*              vcm              = nullptr;
580             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
581
582             compute_globals(gstat, cr, ir, fr, ekind, makeConstArrayRef(state->x),
583                             makeConstArrayRef(state->v), state->box, state->lambda[efptVDW], mdatoms,
584                             nrnb, vcm, wcycle, enerd, force_vir, shake_vir, total_vir, pres, constr,
585                             &signaller, state->box, &totalNumberOfBondedInteractions, &bSumEkinhOld,
586                             CGLO_GSTAT | CGLO_ENERGY
587                                     | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
588                                                                              : 0));
589             checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions, top_global,
590                                             &top, makeConstArrayRef(state->x), state->box,
591                                             &shouldCheckNumberOfBondedInteractions);
592         }
593
594         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
595            the virial that should probably be addressed eventually. state->veta has better properies,
596            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
597            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
598
599         if (ir->efep != efepNO)
600         {
601             /* Sum up the foreign energy and dhdl terms for md and sd.
602                Currently done every step so that dhdl is correct in the .edr */
603             sum_dhdl(enerd, state->lambda, *ir->fepvals);
604         }
605
606         /* Output stuff */
607         if (MASTER(cr))
608         {
609             const bool bCalcEnerStep = true;
610             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep, t,
611                                              mdatoms->tmass, enerd, state, ir->fepvals,
612                                              ir->expandedvals, state->box, shake_vir, force_vir,
613                                              total_vir, pres, ekind, mu_tot, constr);
614
615             const bool do_ene = true;
616             const bool do_log = true;
617             Awh*       awh    = nullptr;
618             const bool do_dr  = ir->nstdisreout != 0;
619             const bool do_or  = ir->nstorireout != 0;
620
621             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
622             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
623                                                do_log ? fplog : nullptr, step, t, fcd, awh);
624
625             if (do_per_step(step, ir->nstlog))
626             {
627                 if (fflush(fplog) != 0)
628                 {
629                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
630                 }
631             }
632         }
633
634         /* Print the remaining wall clock time for the run */
635         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
636         {
637             if (shellfc)
638             {
639                 fprintf(stderr, "\n");
640             }
641             print_time(stderr, walltime_accounting, step, ir, cr);
642         }
643
644         /* Ion/water position swapping.
645          * Not done in last step since trajectory writing happens before this call
646          * in the MD loop and exchanges would be lost anyway. */
647         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep && do_per_step(step, ir->swap->nstswap))
648         {
649             const bool doRerun = true;
650             do_swapcoords(cr, step, t, ir, swap, wcycle, rerun_fr.x, rerun_fr.box,
651                           MASTER(cr) && mdrunOptions.verbose, doRerun);
652         }
653
654         if (MASTER(cr))
655         {
656             /* read next frame from input trajectory */
657             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
658         }
659
660         if (PAR(cr))
661         {
662             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
663         }
664
665         cycles = wallcycle_stop(wcycle, ewcSTEP);
666         if (DOMAINDECOMP(cr) && wcycle)
667         {
668             dd_cycles_add(cr->dd, cycles, ddCyclStep);
669         }
670
671         if (!rerun_fr.bStep)
672         {
673             /* increase the MD step number */
674             step++;
675             step_rel++;
676         }
677     }
678     /* End of main MD loop */
679
680     /* Closing TNG files can include compressing data. Therefore it is good to do that
681      * before stopping the time measurements. */
682     mdoutf_tng_close(outf);
683
684     /* Stop measuring walltime */
685     walltime_accounting_end_time(walltime_accounting);
686
687     if (MASTER(cr))
688     {
689         close_trx(status);
690     }
691
692     if (!thisRankHasDuty(cr, DUTY_PME))
693     {
694         /* Tell the PME only node to finish */
695         gmx_pme_send_finish(cr);
696     }
697
698     done_mdoutf(outf);
699
700     done_shellfc(fplog, shellfc, step_rel);
701
702     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
703 }