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