Use more modern clang for linters and sanitizers
[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,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  *
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/applied_forces/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/listed_forces.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/freeenergyparameters.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/resethandler.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/simulationsignal.h"
91 #include "gromacs/mdlib/stat.h"
92 #include "gromacs/mdlib/stophandler.h"
93 #include "gromacs/mdlib/tgroup.h"
94 #include "gromacs/mdlib/trajectory_writing.h"
95 #include "gromacs/mdlib/update.h"
96 #include "gromacs/mdlib/vcm.h"
97 #include "gromacs/mdlib/vsite.h"
98 #include "gromacs/mdrunutility/handlerestart.h"
99 #include "gromacs/mdrunutility/multisim.h"
100 #include "gromacs/mdrunutility/printtime.h"
101 #include "gromacs/mdtypes/awh_history.h"
102 #include "gromacs/mdtypes/awh_params.h"
103 #include "gromacs/mdtypes/commrec.h"
104 #include "gromacs/mdtypes/df_history.h"
105 #include "gromacs/mdtypes/energyhistory.h"
106 #include "gromacs/mdtypes/forcebuffers.h"
107 #include "gromacs/mdtypes/forcerec.h"
108 #include "gromacs/mdtypes/group.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/interaction_const.h"
111 #include "gromacs/mdtypes/md_enums.h"
112 #include "gromacs/mdtypes/mdatom.h"
113 #include "gromacs/mdtypes/mdrunoptions.h"
114 #include "gromacs/mdtypes/observableshistory.h"
115 #include "gromacs/mdtypes/simulation_workload.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/mimic/utilities.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/timing/wallcycle.h"
122 #include "gromacs/timing/walltime_accounting.h"
123 #include "gromacs/topology/atoms.h"
124 #include "gromacs/topology/idef.h"
125 #include "gromacs/topology/mtop_util.h"
126 #include "gromacs/topology/topology.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basedefinitions.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/fatalerror.h"
131 #include "gromacs/utility/logger.h"
132 #include "gromacs/utility/real.h"
133
134 #include "legacysimulator.h"
135 #include "replicaexchange.h"
136 #include "shellfc.h"
137
138 using gmx::SimulationSignaller;
139 using gmx::VirtualSitesHandler;
140
141 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
142  *
143  * \param[in]     rerunFrame      The trajectory frame to compute energy/forces for
144  * \param[in,out] globalState     The global state container
145  * \param[in]     constructVsites When true, vsite coordinates are constructed
146  * \param[in]     vsite           Vsite setup, can be nullptr when \p constructVsites = false
147  */
148 static void prepareRerunState(const t_trxframe&          rerunFrame,
149                               t_state*                   globalState,
150                               bool                       constructVsites,
151                               const VirtualSitesHandler* vsite)
152 {
153     auto x      = makeArrayRef(globalState->x);
154     auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec*>(rerunFrame.x), globalState->natoms);
155     std::copy(rerunX.begin(), rerunX.end(), x.begin());
156     copy_mat(rerunFrame.box, globalState->box);
157
158     if (constructVsites)
159     {
160         GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
161
162         vsite->construct(globalState->x, globalState->v, globalState->box, gmx::VSiteOperation::PositionsAndVelocities);
163     }
164 }
165
166 void gmx::LegacySimulator::do_rerun()
167 {
168     // TODO Historically, the EM and MD "integrators" used different
169     // names for the t_inputrec *parameter, but these must have the
170     // same name, now that it's a member of a struct. We use this ir
171     // alias to avoid a large ripple of nearly useless changes.
172     // t_inputrec is being replaced by IMdpOptionsProvider, so this
173     // will go away eventually.
174     const t_inputrec* ir = inputrec;
175     int64_t           step, step_rel;
176     double            t;
177     bool              isLastStep               = false;
178     bool              doFreeEnergyPerturbation = false;
179     unsigned int      force_flags;
180     tensor            force_vir, shake_vir, total_vir, pres;
181     t_trxstatus*      status = nullptr;
182     rvec              mu_tot;
183     t_trxframe        rerun_fr;
184     gmx_localtop_t    top(top_global.ffparams);
185     ForceBuffers      f;
186     gmx_global_stat_t gstat;
187     gmx_shellfc_t*    shellfc;
188
189     double cycles;
190
191     SimulationSignals signals;
192     // Most global communnication stages don't propagate mdrun
193     // signals, and will use this object to achieve that.
194     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
195
196     GMX_LOG(mdlog.info)
197             .asParagraph()
198             .appendText(
199                     "Note that it is planned that the command gmx mdrun -rerun will "
200                     "be available in a different form in a future version of GROMACS, "
201                     "e.g. gmx rerun -f.");
202
203     if (ir->efep != FreeEnergyPerturbationType::No
204         && (mdAtoms->mdatoms()->nMassPerturbed > 0 || (constr && constr->havePerturbedConstraints())))
205     {
206         gmx_fatal(FARGS,
207                   "Perturbed masses or constraints are not supported by rerun. "
208                   "Either make a .tpr without mass and constraint perturbation, "
209                   "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
210     }
211     if (ir->bExpanded)
212     {
213         gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
214     }
215     if (ir->bSimTemp)
216     {
217         gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
218     }
219     if (ir->bDoAwh)
220     {
221         gmx_fatal(FARGS, "AWH not supported by rerun.");
222     }
223     if (replExParams.exchangeInterval > 0)
224     {
225         gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
226     }
227     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
228     {
229         gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
230     }
231     if (ir->bIMD)
232     {
233         gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
234     }
235     if (isMultiSim(ms))
236     {
237         gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
238     }
239     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc, [](SimulatedAnnealing i) {
240             return i != SimulatedAnnealing::No;
241         }))
242     {
243         gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
244     }
245
246     /* Rerun can't work if an output file name is the same as the input file name.
247      * If this is the case, the user will get an error telling them what the issue is.
248      */
249     if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0
250         || strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
251     {
252         gmx_fatal(FARGS,
253                   "When using mdrun -rerun, the name of the input trajectory file "
254                   "%s cannot be identical to the name of an output file (whether "
255                   "given explicitly with -o or -x, or by default)",
256                   opt2fn("-rerun", nfile, fnm));
257     }
258
259     /* Settings for rerun */
260     {
261         // TODO: Avoid changing inputrec (#3854)
262         auto* nonConstInputrec               = const_cast<t_inputrec*>(inputrec);
263         nonConstInputrec->nstlist            = 1;
264         nonConstInputrec->nstcalcenergy      = 1;
265         nonConstInputrec->nstxout_compressed = 0;
266     }
267     int        nstglobalcomm = 1;
268     const bool bNS           = true;
269
270     const SimulationGroups* groups = &top_global.groups;
271     if (ir->eI == IntegrationAlgorithm::Mimic)
272     {
273         auto* nonConstGlobalTopology                         = const_cast<gmx_mtop_t*>(&top_global);
274         nonConstGlobalTopology->intermolecularExclusionGroup = genQmmmIndices(top_global);
275     }
276     int*                fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
277     gmx::ArrayRef<real> lambda    = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
278     initialize_lambdas(fplog,
279                        ir->efep,
280                        ir->bSimTemp,
281                        *ir->fepvals,
282                        ir->simtempvals->temperatures,
283                        gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
284                        MASTER(cr),
285                        fep_state,
286                        lambda);
287     const bool        simulationsShareState = false;
288     gmx_mdoutf*       outf                  = init_mdoutf(fplog,
289                                    nfile,
290                                    fnm,
291                                    mdrunOptions,
292                                    cr,
293                                    outputProvider,
294                                    mdModulesNotifiers,
295                                    ir,
296                                    top_global,
297                                    oenv,
298                                    wcycle,
299                                    StartingBehavior::NewSimulation,
300                                    simulationsShareState,
301                                    ms);
302     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
303                                    top_global,
304                                    *ir,
305                                    pull_work,
306                                    mdoutf_get_fp_dhdl(outf),
307                                    true,
308                                    StartingBehavior::NewSimulation,
309                                    simulationsShareState,
310                                    mdModulesNotifiers);
311
312     gstat = global_stat_init(ir);
313
314     /* Check for polarizable models and flexible constraints */
315     shellfc = init_shell_flexcon(fplog,
316                                  top_global,
317                                  constr ? constr->numFlexibleConstraints() : 0,
318                                  ir->nstcalcenergy,
319                                  DOMAINDECOMP(cr),
320                                  runScheduleWork->simulationWork.useGpuPme);
321
322     {
323         double io = compute_io(ir, top_global.natoms, *groups, energyOutput.numEnergyTerms(), 1);
324         if ((io > 2000) && MASTER(cr))
325         {
326             fprintf(stderr, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io);
327         }
328     }
329
330     // Local state only becomes valid now.
331     std::unique_ptr<t_state> stateInstance;
332     t_state*                 state;
333
334     if (DOMAINDECOMP(cr))
335     {
336         stateInstance = std::make_unique<t_state>();
337         state         = stateInstance.get();
338         dd_init_local_state(*cr->dd, state_global, state);
339
340         /* Distribute the charge groups over the nodes from the master node */
341         dd_partition_system(fplog,
342                             mdlog,
343                             ir->init_step,
344                             cr,
345                             TRUE,
346                             1,
347                             state_global,
348                             top_global,
349                             *ir,
350                             imdSession,
351                             pull_work,
352                             state,
353                             &f,
354                             mdAtoms,
355                             &top,
356                             fr,
357                             vsite,
358                             constr,
359                             nrnb,
360                             nullptr,
361                             FALSE);
362     }
363     else
364     {
365         state_change_natoms(state_global, state_global->natoms);
366         /* Copy the pointer to the global state */
367         state = state_global;
368
369         mdAlgorithmsSetupAtomData(cr, *ir, top_global, &top, fr, &f, mdAtoms, constr, vsite, shellfc);
370     }
371
372     auto* mdatoms = mdAtoms->mdatoms();
373
374     // NOTE: The global state is no longer used at this point.
375     // But state_global is still used as temporary storage space for writing
376     // the global state to file and potentially for replica exchange.
377     // (Global topology should persist.)
378
379     update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
380
381     if (ir->efep != FreeEnergyPerturbationType::No && ir->fepvals->nstdhdl != 0)
382     {
383         doFreeEnergyPerturbation = true;
384     }
385
386     {
387         int cglo_flags = CGLO_GSTAT;
388         if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
389         {
390             cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
391         }
392         bool   bSumEkinhOld = false;
393         t_vcm* vcm          = nullptr;
394         compute_globals(gstat,
395                         cr,
396                         ir,
397                         fr,
398                         ekind,
399                         makeConstArrayRef(state->x),
400                         makeConstArrayRef(state->v),
401                         state->box,
402                         mdatoms,
403                         nrnb,
404                         vcm,
405                         nullptr,
406                         enerd,
407                         force_vir,
408                         shake_vir,
409                         total_vir,
410                         pres,
411                         gmx::ArrayRef<real>{},
412                         &nullSignaller,
413                         state->box,
414                         &bSumEkinhOld,
415                         cglo_flags);
416         if (DOMAINDECOMP(cr))
417         {
418             checkNumberOfBondedInteractions(
419                     mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
420         }
421     }
422
423     if (MASTER(cr))
424     {
425         fprintf(stderr,
426                 "starting md rerun '%s', reading coordinates from"
427                 " input trajectory '%s'\n\n",
428                 *(top_global.name),
429                 opt2fn("-rerun", nfile, fnm));
430         if (mdrunOptions.verbose)
431         {
432             fprintf(stderr,
433                     "Calculated time to finish depends on nsteps from "
434                     "run input file,\nwhich may not correspond to the time "
435                     "needed to process input trajectory.\n\n");
436         }
437         fprintf(fplog, "\n");
438     }
439
440     walltime_accounting_start_time(walltime_accounting);
441     wallcycle_start(wcycle, WallCycleCounter::Run);
442     print_start(fplog, cr, walltime_accounting, "mdrun");
443
444     /***********************************************************
445      *
446      *             Loop over MD steps
447      *
448      ************************************************************/
449
450     if (constr)
451     {
452         GMX_LOG(mdlog.info)
453                 .asParagraph()
454                 .appendText("Simulations has constraints. Rerun does not recalculate constraints.");
455     }
456
457     rerun_fr.natoms = 0;
458     if (MASTER(cr))
459     {
460         isLastStep = !read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm), &rerun_fr, TRX_NEED_X);
461         if (rerun_fr.natoms != top_global.natoms)
462         {
463             gmx_fatal(FARGS,
464                       "Number of atoms in trajectory (%d) does not match the "
465                       "run input file (%d)\n",
466                       rerun_fr.natoms,
467                       top_global.natoms);
468         }
469
470         if (ir->pbcType != PbcType::No)
471         {
472             if (!rerun_fr.bBox)
473             {
474                 gmx_fatal(FARGS,
475                           "Rerun trajectory frame step %" PRId64
476                           " time %f "
477                           "does not contain a box, while pbc is used",
478                           rerun_fr.step,
479                           rerun_fr.time);
480             }
481             if (max_cutoff2(ir->pbcType, rerun_fr.box) < gmx::square(fr->rlist))
482             {
483                 gmx_fatal(FARGS,
484                           "Rerun trajectory frame step %" PRId64
485                           " time %f "
486                           "has too small box dimensions",
487                           rerun_fr.step,
488                           rerun_fr.time);
489             }
490         }
491     }
492
493     GMX_LOG(mdlog.info)
494             .asParagraph()
495             .appendText(
496                     "Rerun does not report kinetic energy, total energy, temperature, virial and "
497                     "pressure.");
498
499     if (PAR(cr))
500     {
501         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
502     }
503
504     if (ir->pbcType != PbcType::No)
505     {
506         /* Set the shift vectors.
507          * Necessary here when have a static box different from the tpr box.
508          */
509         calc_shifts(rerun_fr.box, fr->shift_vec);
510     }
511
512     step     = ir->init_step;
513     step_rel = 0;
514
515     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
516             compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]),
517             false,
518             MASTER(cr),
519             ir->nstlist,
520             mdrunOptions.reproducible,
521             nstglobalcomm,
522             mdrunOptions.maximumHoursToRun,
523             ir->nstlist == 0,
524             fplog,
525             step,
526             bNS,
527             walltime_accounting);
528
529     // we don't do counter resetting in rerun - finish will always be valid
530     walltime_accounting_set_valid_finish(walltime_accounting);
531
532     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
533
534     /* and stop now if we should */
535     isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
536     while (!isLastStep)
537     {
538         wallcycle_start(wcycle, WallCycleCounter::Step);
539
540         if (rerun_fr.bStep)
541         {
542             step     = rerun_fr.step;
543             step_rel = step - ir->init_step;
544         }
545         if (rerun_fr.bTime)
546         {
547             t = rerun_fr.time;
548         }
549         else
550         {
551             t = step;
552         }
553
554         if (ir->efep != FreeEnergyPerturbationType::No && MASTER(cr))
555         {
556             if (rerun_fr.bLambda)
557             {
558                 ir->fepvals->init_lambda = rerun_fr.lambda;
559             }
560             else
561             {
562                 if (rerun_fr.bFepState)
563                 {
564                     state->fep_state = rerun_fr.fep_state;
565                 }
566             }
567
568             state_global->lambda = currentLambdas(step, *(ir->fepvals), state->fep_state);
569         }
570
571         if (MASTER(cr))
572         {
573             const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
574             if (constructVsites && DOMAINDECOMP(cr))
575             {
576                 gmx_fatal(FARGS,
577                           "Vsite recalculation with -rerun is not implemented with domain "
578                           "decomposition, "
579                           "use a single rank");
580             }
581             prepareRerunState(rerun_fr, state_global, constructVsites, vsite);
582         }
583
584         isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
585
586         if (DOMAINDECOMP(cr))
587         {
588             /* Repartition the domain decomposition */
589             const bool bMasterState = true;
590             dd_partition_system(fplog,
591                                 mdlog,
592                                 step,
593                                 cr,
594                                 bMasterState,
595                                 nstglobalcomm,
596                                 state_global,
597                                 top_global,
598                                 *ir,
599                                 imdSession,
600                                 pull_work,
601                                 state,
602                                 &f,
603                                 mdAtoms,
604                                 &top,
605                                 fr,
606                                 vsite,
607                                 constr,
608                                 nrnb,
609                                 wcycle,
610                                 mdrunOptions.verbose);
611         }
612
613         if (MASTER(cr))
614         {
615             EnergyOutput::printHeader(fplog, step, t); /* can we improve the information printed here? */
616         }
617
618         if (ir->efep != FreeEnergyPerturbationType::No)
619         {
620             update_mdatoms(mdatoms, state->lambda[FreeEnergyPerturbationCouplingType::Mass]);
621         }
622
623         force_flags = (GMX_FORCE_STATECHANGED | GMX_FORCE_DYNAMICBOX | GMX_FORCE_ALLFORCES
624                        | GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 and #3400 are solved
625                        GMX_FORCE_ENERGY | (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
626
627         if (shellfc)
628         {
629             /* Now is the time to relax the shells */
630             relax_shell_flexcon(fplog,
631                                 cr,
632                                 ms,
633                                 mdrunOptions.verbose,
634                                 enforcedRotation,
635                                 step,
636                                 ir,
637                                 imdSession,
638                                 pull_work,
639                                 bNS,
640                                 force_flags,
641                                 &top,
642                                 constr,
643                                 enerd,
644                                 state->natoms,
645                                 state->x.arrayRefWithPadding(),
646                                 state->v.arrayRefWithPadding(),
647                                 state->box,
648                                 state->lambda,
649                                 &state->hist,
650                                 &f.view(),
651                                 force_vir,
652                                 *mdatoms,
653                                 nrnb,
654                                 wcycle,
655                                 shellfc,
656                                 fr,
657                                 runScheduleWork,
658                                 t,
659                                 mu_tot,
660                                 vsite,
661                                 ddBalanceRegionHandler);
662         }
663         else
664         {
665             /* The coordinates (x) are shifted (to get whole molecules)
666              * in do_force.
667              * This is parallellized as well, and does communication too.
668              * Check comments in sim_util.c
669              */
670             Awh*       awh = nullptr;
671             gmx_edsam* ed  = nullptr;
672             do_force(fplog,
673                      cr,
674                      ms,
675                      *ir,
676                      awh,
677                      enforcedRotation,
678                      imdSession,
679                      pull_work,
680                      step,
681                      nrnb,
682                      wcycle,
683                      &top,
684                      state->box,
685                      state->x.arrayRefWithPadding(),
686                      &state->hist,
687                      &f.view(),
688                      force_vir,
689                      mdatoms,
690                      enerd,
691                      state->lambda,
692                      fr,
693                      runScheduleWork,
694                      vsite,
695                      mu_tot,
696                      t,
697                      ed,
698                      GMX_FORCE_NS | force_flags,
699                      ddBalanceRegionHandler);
700         }
701
702         /* Now we have the energies and forces corresponding to the
703          * coordinates at time t.
704          */
705         {
706             const bool isCheckpointingStep = false;
707             const bool doRerun             = true;
708             const bool bSumEkinhOld        = false;
709             do_md_trajectory_writing(fplog,
710                                      cr,
711                                      nfile,
712                                      fnm,
713                                      step,
714                                      step_rel,
715                                      t,
716                                      ir,
717                                      state,
718                                      state_global,
719                                      observablesHistory,
720                                      top_global,
721                                      fr,
722                                      outf,
723                                      energyOutput,
724                                      ekind,
725                                      f.view().force(),
726                                      isCheckpointingStep,
727                                      doRerun,
728                                      isLastStep,
729                                      mdrunOptions.writeConfout,
730                                      bSumEkinhOld);
731         }
732
733         stopHandler->setSignal();
734
735         {
736             const bool          doInterSimSignal = false;
737             const bool          doIntraSimSignal = true;
738             bool                bSumEkinhOld     = false;
739             t_vcm*              vcm              = nullptr;
740             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
741
742             int cglo_flags = CGLO_GSTAT | CGLO_ENERGY;
743             if (DOMAINDECOMP(cr) && shouldCheckNumberOfBondedInteractions(*cr->dd))
744             {
745                 cglo_flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
746             }
747             compute_globals(gstat,
748                             cr,
749                             ir,
750                             fr,
751                             ekind,
752                             makeConstArrayRef(state->x),
753                             makeConstArrayRef(state->v),
754                             state->box,
755                             mdatoms,
756                             nrnb,
757                             vcm,
758                             wcycle,
759                             enerd,
760                             force_vir,
761                             shake_vir,
762                             total_vir,
763                             pres,
764                             constr != nullptr ? constr->rmsdData() : gmx::ArrayRef<real>{},
765                             &signaller,
766                             state->box,
767                             &bSumEkinhOld,
768                             cglo_flags);
769             if (DOMAINDECOMP(cr))
770             {
771                 checkNumberOfBondedInteractions(
772                         mdlog, cr, top_global, &top, makeConstArrayRef(state->x), state->box);
773             }
774         }
775
776         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
777            the virial that should probably be addressed eventually. state->veta has better properies,
778            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
779            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
780
781         /* Output stuff */
782         if (MASTER(cr))
783         {
784             const bool bCalcEnerStep = true;
785             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation,
786                                              bCalcEnerStep,
787                                              t,
788                                              mdatoms->tmass,
789                                              enerd,
790                                              ir->fepvals.get(),
791                                              ir->expandedvals.get(),
792                                              state->box,
793                                              PTCouplingArrays({ state->boxv,
794                                                                 state->nosehoover_xi,
795                                                                 state->nosehoover_vxi,
796                                                                 state->nhpres_xi,
797                                                                 state->nhpres_vxi }),
798                                              state->fep_state,
799                                              shake_vir,
800                                              force_vir,
801                                              total_vir,
802                                              pres,
803                                              ekind,
804                                              mu_tot,
805                                              constr);
806
807             const bool do_ene = true;
808             const bool do_log = true;
809             Awh*       awh    = nullptr;
810             const bool do_dr  = ir->nstdisreout != 0;
811             const bool do_or  = ir->nstorireout != 0;
812
813             EnergyOutput::printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
814             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
815                                                do_ene,
816                                                do_dr,
817                                                do_or,
818                                                do_log ? fplog : nullptr,
819                                                step,
820                                                t,
821                                                fr->fcdata.get(),
822                                                awh);
823
824             if (do_per_step(step, ir->nstlog))
825             {
826                 if (fflush(fplog) != 0)
827                 {
828                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
829                 }
830             }
831         }
832
833         /* Print the remaining wall clock time for the run */
834         if (isMasterSimMasterRank(ms, MASTER(cr)) && (mdrunOptions.verbose || gmx_got_usr_signal()))
835         {
836             if (shellfc)
837             {
838                 fprintf(stderr, "\n");
839             }
840             print_time(stderr, walltime_accounting, step, ir, cr);
841         }
842
843         /* Ion/water position swapping.
844          * Not done in last step since trajectory writing happens before this call
845          * in the MD loop and exchanges would be lost anyway. */
846         if ((ir->eSwapCoords != SwapType::No) && (step > 0) && !isLastStep
847             && do_per_step(step, ir->swap->nstswap))
848         {
849             const bool doRerun = true;
850             do_swapcoords(cr,
851                           step,
852                           t,
853                           ir,
854                           swap,
855                           wcycle,
856                           rerun_fr.x,
857                           rerun_fr.box,
858                           MASTER(cr) && mdrunOptions.verbose,
859                           doRerun);
860         }
861
862         if (MASTER(cr))
863         {
864             /* read next frame from input trajectory */
865             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
866         }
867
868         if (PAR(cr))
869         {
870             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
871         }
872
873         cycles = wallcycle_stop(wcycle, WallCycleCounter::Step);
874         if (DOMAINDECOMP(cr) && wcycle)
875         {
876             dd_cycles_add(cr->dd, cycles, ddCyclStep);
877         }
878
879         if (!rerun_fr.bStep)
880         {
881             /* increase the MD step number */
882             step++;
883             step_rel++;
884         }
885     }
886     /* End of main MD loop */
887
888     /* Closing TNG files can include compressing data. Therefore it is good to do that
889      * before stopping the time measurements. */
890     mdoutf_tng_close(outf);
891
892     /* Stop measuring walltime */
893     walltime_accounting_end_time(walltime_accounting);
894
895     if (MASTER(cr))
896     {
897         close_trx(status);
898     }
899
900     if (!thisRankHasDuty(cr, DUTY_PME))
901     {
902         /* Tell the PME only node to finish */
903         gmx_pme_send_finish(cr);
904     }
905
906     done_mdoutf(outf);
907
908     done_shellfc(fplog, shellfc, step_rel);
909
910     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
911 }