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