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