Remove unused CGLO_INITIALIZATION
[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/dlbtiming.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/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme_load_balancing.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/gpu_utils/gpu_utils.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/enerdata_utils.h"
80 #include "gromacs/mdlib/energyoutput.h"
81 #include "gromacs/mdlib/expanded.h"
82 #include "gromacs/mdlib/force.h"
83 #include "gromacs/mdlib/force_flags.h"
84 #include "gromacs/mdlib/forcerec.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/fcdata.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/state.h"
117 #include "gromacs/mimic/utilities.h"
118 #include "gromacs/pbcutil/mshift.h"
119 #include "gromacs/pbcutil/pbc.h"
120 #include "gromacs/pulling/pull.h"
121 #include "gromacs/swap/swapcoords.h"
122 #include "gromacs/timing/wallcycle.h"
123 #include "gromacs/timing/walltime_accounting.h"
124 #include "gromacs/topology/atoms.h"
125 #include "gromacs/topology/idef.h"
126 #include "gromacs/topology/mtop_util.h"
127 #include "gromacs/topology/topology.h"
128 #include "gromacs/trajectory/trajectoryframe.h"
129 #include "gromacs/utility/basedefinitions.h"
130 #include "gromacs/utility/cstringutil.h"
131 #include "gromacs/utility/fatalerror.h"
132 #include "gromacs/utility/logger.h"
133 #include "gromacs/utility/real.h"
134 #include "gromacs/utility/smalloc.h"
135
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
138 #include "shellfc.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::LegacySimulator::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     unsigned 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     PaddedVector<gmx::RVec> f {};
209     gmx_global_stat_t       gstat;
210     t_graph                *graph = nullptr;
211     gmx_shellfc_t          *shellfc;
212
213     double                  cycles;
214
215     /* Domain decomposition could incorrectly miss a bonded
216        interaction, but checking for that requires a global
217        communication stage, which does not otherwise happen in DD
218        code. So we do that alongside the first global energy reduction
219        after a new DD is made. These variables handle whether the
220        check happens, and the result it returns. */
221     bool              shouldCheckNumberOfBondedInteractions = false;
222     int               totalNumberOfBondedInteractions       = -1;
223
224     SimulationSignals signals;
225     // Most global communnication stages don't propagate mdrun
226     // signals, and will use this object to achieve that.
227     SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
228
229     GMX_LOG(mdlog.info).asParagraph().
230         appendText("Note that it is planned that the command gmx mdrun -rerun will "
231                    "be available in a different form in a future version of GROMACS, "
232                    "e.g. gmx rerun -f.");
233
234     if (ir->efep != efepNO && (mdAtoms->mdatoms()->nMassPerturbed > 0 ||
235                                (constr && constr->havePerturbedConstraints())))
236     {
237         gmx_fatal(FARGS, "Perturbed masses or constraints are not supported by rerun. "
238                   "Either make a .tpr without mass and constraint perturbation, "
239                   "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
240     }
241     if (ir->bExpanded)
242     {
243         gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
244     }
245     if (ir->bSimTemp)
246     {
247         gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
248     }
249     if (ir->bDoAwh)
250     {
251         gmx_fatal(FARGS, "AWH not supported by rerun.");
252     }
253     if (replExParams.exchangeInterval > 0)
254     {
255         gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
256     }
257     if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
258     {
259         gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
260     }
261     if (ir->bIMD)
262     {
263         gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
264     }
265     if (isMultiSim(ms))
266     {
267         gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
268     }
269     if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
270                     [](int i){return i != eannNO; }))
271     {
272         gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
273     }
274
275     /* Rerun can't work if an output file name is the same as the input file name.
276      * If this is the case, the user will get an error telling them what the issue is.
277      */
278     if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0 ||
279         strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
280     {
281         gmx_fatal(FARGS, "When using mdrun -rerun, the name of the input trajectory file "
282                   "%s cannot be identical to the name of an output file (whether "
283                   "given explicitly with -o or -x, or by default)",
284                   opt2fn("-rerun", nfile, fnm));
285     }
286
287     /* Settings for rerun */
288     ir->nstlist       = 1;
289     ir->nstcalcenergy = 1;
290     int        nstglobalcomm = 1;
291     const bool bNS           = true;
292
293     ir->nstxout_compressed = 0;
294     SimulationGroups *groups                 = &top_global->groups;
295     if (ir->eI == eiMimic)
296     {
297         top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
298     }
299
300     initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
301     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle,
302                                          StartingBehavior::NewSimulation);
303     gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, ir, pull_work, mdoutf_get_fp_dhdl(outf), true);
304
305     /* Kinetic energy data */
306     std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
307     gmx_ekindata_t                 *ekind    = eKinData.get();
308     init_ekindata(fplog, top_global, &(ir->opts), ekind);
309     /* Copy the cos acceleration to the groups struct */
310     ekind->cosacc.cos_accel = ir->cos_accel;
311
312     gstat = global_stat_init(ir);
313
314     /* Check for polarizable models and flexible constraints */
315     shellfc = init_shell_flexcon(fplog,
316                                  top_global, constr ? constr->numFlexibleConstraints() : 0,
317                                  ir->nstcalcenergy, DOMAINDECOMP(cr));
318
319     {
320         double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
321         if ((io > 2000) && MASTER(cr))
322         {
323             fprintf(stderr,
324                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
325                     io);
326         }
327     }
328
329     // Local state only becomes valid now.
330     std::unique_ptr<t_state> stateInstance;
331     t_state *                state;
332
333     if (DOMAINDECOMP(cr))
334     {
335         dd_init_local_top(*top_global, &top);
336
337         stateInstance = std::make_unique<t_state>();
338         state         = stateInstance.get();
339         dd_init_local_state(cr->dd, state_global, state);
340
341         /* Distribute the charge groups over the nodes from the master node */
342         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
343                             state_global, *top_global, ir, imdSession,
344                             pull_work,
345                             state, &f, mdAtoms, &top, fr,
346                             vsite, constr,
347                             nrnb, nullptr, FALSE);
348         shouldCheckNumberOfBondedInteractions = true;
349     }
350     else
351     {
352         state_change_natoms(state_global, state_global->natoms);
353         /* We need to allocate one element extra, since we might use
354          * (unaligned) 4-wide SIMD loads to access rvec entries.
355          */
356         f.resizeWithPadding(state_global->natoms);
357         /* Copy the pointer to the global state */
358         state = state_global;
359
360         mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
361                                   &graph, mdAtoms, constr, vsite, shellfc);
362     }
363
364     auto mdatoms = mdAtoms->mdatoms();
365
366     // NOTE: The global state is no longer used at this point.
367     // But state_global is still used as temporary storage space for writing
368     // the global state to file and potentially for replica exchange.
369     // (Global topology should persist.)
370
371     update_mdatoms(mdatoms, state->lambda[efptMASS]);
372
373     if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
374     {
375         doFreeEnergyPerturbation = true;
376     }
377
378     {
379         int    cglo_flags = (CGLO_GSTAT |
380                              (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
381         bool   bSumEkinhOld = false;
382         t_vcm *vcm          = nullptr;
383         compute_globals(gstat, cr, ir, fr, ekind,
384                         state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
385                         mdatoms, nrnb, vcm,
386                         nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
387                         constr, &nullSignaller, state->box,
388                         &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
389     }
390     checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
391                                     top_global, &top, state->x.rvec_array(), state->box,
392                                     &shouldCheckNumberOfBondedInteractions);
393
394     if (MASTER(cr))
395     {
396         fprintf(stderr, "starting md rerun '%s', reading coordinates from"
397                 " input trajectory '%s'\n\n",
398                 *(top_global->name), opt2fn("-rerun", nfile, fnm));
399         if (mdrunOptions.verbose)
400         {
401             fprintf(stderr, "Calculated time to finish depends on nsteps from "
402                     "run input file,\nwhich may not correspond to the time "
403                     "needed to process input trajectory.\n\n");
404         }
405         fprintf(fplog, "\n");
406     }
407
408     walltime_accounting_start_time(walltime_accounting);
409     wallcycle_start(wcycle, ewcRUN);
410     print_start(fplog, cr, walltime_accounting, "mdrun");
411
412     /***********************************************************
413      *
414      *             Loop over MD steps
415      *
416      ************************************************************/
417
418     if (constr)
419     {
420         GMX_LOG(mdlog.info).asParagraph().
421             appendText("Simulations has constraints. Rerun does not recalculate constraints.");
422     }
423
424     rerun_fr.natoms = 0;
425     if (MASTER(cr))
426     {
427         isLastStep = !read_first_frame(oenv, &status,
428                                        opt2fn("-rerun", nfile, fnm),
429                                        &rerun_fr, TRX_NEED_X);
430         if (rerun_fr.natoms != top_global->natoms)
431         {
432             gmx_fatal(FARGS,
433                       "Number of atoms in trajectory (%d) does not match the "
434                       "run input file (%d)\n",
435                       rerun_fr.natoms, top_global->natoms);
436         }
437
438         if (ir->ePBC != epbcNONE)
439         {
440             if (!rerun_fr.bBox)
441             {
442                 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
443                           "does not contain a box, while pbc is used",
444                           rerun_fr.step, rerun_fr.time);
445             }
446             if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
447             {
448                 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
449                           "has too small box dimensions", rerun_fr.step, rerun_fr.time);
450             }
451         }
452     }
453
454     GMX_LOG(mdlog.info).asParagraph().
455         appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
456
457     if (PAR(cr))
458     {
459         rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
460     }
461
462     if (ir->ePBC != epbcNONE)
463     {
464         /* Set the shift vectors.
465          * Necessary here when have a static box different from the tpr box.
466          */
467         calc_shifts(rerun_fr.box, fr->shift_vec);
468     }
469
470     auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
471                 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
472                 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
473                 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
474
475     // we don't do counter resetting in rerun - finish will always be valid
476     walltime_accounting_set_valid_finish(walltime_accounting);
477
478     const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
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, imdSession,
528                                 pull_work,
529                                 state, &f, mdAtoms, &top, fr,
530                                 vsite, constr,
531                                 nrnb, wcycle,
532                                 mdrunOptions.verbose);
533             shouldCheckNumberOfBondedInteractions = true;
534         }
535
536         if (MASTER(cr))
537         {
538             energyOutput.printHeader(fplog, step, t); /* can we improve the information printed here? */
539         }
540
541         if (ir->efep != efepNO)
542         {
543             update_mdatoms(mdatoms, state->lambda[efptMASS]);
544         }
545
546         force_flags = (GMX_FORCE_STATECHANGED |
547                        GMX_FORCE_DYNAMICBOX |
548                        GMX_FORCE_ALLFORCES |
549                        (GMX_GPU ? GMX_FORCE_VIRIAL : 0) |  // TODO: Get rid of this once #2649 is solved
550                        GMX_FORCE_ENERGY |
551                        (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
552
553         if (shellfc)
554         {
555             /* Now is the time to relax the shells */
556             relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
557                                 enforcedRotation, step,
558                                 ir, imdSession, pull_work, bNS, force_flags, &top,
559                                 constr, enerd, fcd,
560                                 state->natoms,
561                                 state->x.arrayRefWithPadding(),
562                                 state->v.arrayRefWithPadding(),
563                                 state->box,
564                                 state->lambda,
565                                 &state->hist,
566                                 f.arrayRefWithPadding(), force_vir, mdatoms,
567                                 nrnb, wcycle, graph,
568                                 shellfc, fr, ppForceWorkload, t, mu_tot,
569                                 vsite,
570                                 ddBalanceRegionHandler);
571         }
572         else
573         {
574             /* The coordinates (x) are shifted (to get whole molecules)
575              * in do_force.
576              * This is parallellized as well, and does communication too.
577              * Check comments in sim_util.c
578              */
579             Awh       *awh = nullptr;
580             gmx_edsam *ed  = nullptr;
581             do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession,
582                      pull_work,
583                      step, nrnb, wcycle, &top,
584                      state->box, state->x.arrayRefWithPadding(), &state->hist,
585                      f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
586                      state->lambda, graph,
587                      fr, ppForceWorkload, vsite, mu_tot, t, ed,
588                      GMX_FORCE_NS | force_flags,
589                      ddBalanceRegionHandler);
590         }
591
592         /* Now we have the energies and forces corresponding to the
593          * coordinates at time t.
594          */
595         {
596             const bool isCheckpointingStep = false;
597             const bool doRerun             = true;
598             const bool bSumEkinhOld        = false;
599             do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
600                                      ir, state, state_global, observablesHistory,
601                                      top_global, fr,
602                                      outf, energyOutput, ekind, f,
603                                      isCheckpointingStep, doRerun, isLastStep,
604                                      mdrunOptions.writeConfout,
605                                      bSumEkinhOld);
606         }
607
608         stopHandler->setSignal();
609
610         if (graph)
611         {
612             /* Need to unshift here */
613             unshift_self(graph, state->box, as_rvec_array(state->x.data()));
614         }
615
616         if (vsite != nullptr)
617         {
618             wallcycle_start(wcycle, ewcVSITECONSTR);
619             if (graph != nullptr)
620             {
621                 shift_self(graph, state->box, as_rvec_array(state->x.data()));
622             }
623             construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
624                              top.idef.iparams, top.idef.il,
625                              fr->ePBC, fr->bMolPBC, cr, state->box);
626
627             if (graph != nullptr)
628             {
629                 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
630             }
631             wallcycle_stop(wcycle, ewcVSITECONSTR);
632         }
633
634         {
635             const bool          doInterSimSignal = false;
636             const bool          doIntraSimSignal = true;
637             bool                bSumEkinhOld     = false;
638             t_vcm              *vcm              = nullptr;
639             SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
640
641             compute_globals(gstat, cr, ir, fr, ekind,
642                             state->x.rvec_array(), state->v.rvec_array(), state->box, state->lambda[efptVDW],
643                             mdatoms, nrnb, vcm,
644                             wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
645                             constr, &signaller,
646                             state->box,
647                             &totalNumberOfBondedInteractions, &bSumEkinhOld,
648                             CGLO_GSTAT | CGLO_ENERGY
649                             | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
650                             );
651             checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
652                                             top_global, &top, state->x.rvec_array(), state->box,
653                                             &shouldCheckNumberOfBondedInteractions);
654         }
655
656         /* Note: this is OK, but there are some numerical precision issues with using the convergence of
657            the virial that should probably be addressed eventually. state->veta has better properies,
658            but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
659            generate the new shake_vir, but test the veta value for convergence.  This will take some thought. */
660
661         if (ir->efep != efepNO)
662         {
663             /* Sum up the foreign energy and dhdl terms for md and sd.
664                Currently done every step so that dhdl is correct in the .edr */
665             sum_dhdl(enerd, state->lambda, ir->fepvals);
666         }
667
668         /* Output stuff */
669         if (MASTER(cr))
670         {
671             const bool bCalcEnerStep = true;
672             energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep,
673                                              t, mdatoms->tmass, enerd, state,
674                                              ir->fepvals, ir->expandedvals, state->box,
675                                              shake_vir, force_vir, total_vir, pres,
676                                              ekind, mu_tot,
677                                              constr);
678
679             const bool do_ene = true;
680             const bool do_log = true;
681             Awh       *awh    = nullptr;
682             const bool do_dr  = ir->nstdisreout != 0;
683             const bool do_or  = ir->nstorireout != 0;
684
685             energyOutput.printAnnealingTemperatures(do_log ? fplog : nullptr, groups, &(ir->opts));
686             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
687                                                do_log ? fplog : nullptr,
688                                                step, t,
689                                                fcd, awh);
690
691             if (do_per_step(step, ir->nstlog))
692             {
693                 if (fflush(fplog) != 0)
694                 {
695                     gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
696                 }
697             }
698         }
699
700         /* Print the remaining wall clock time for the run */
701         if (isMasterSimMasterRank(ms, MASTER(cr)) &&
702             (mdrunOptions.verbose || gmx_got_usr_signal()))
703         {
704             if (shellfc)
705             {
706                 fprintf(stderr, "\n");
707             }
708             print_time(stderr, walltime_accounting, step, ir, cr);
709         }
710
711         /* Ion/water position swapping.
712          * Not done in last step since trajectory writing happens before this call
713          * in the MD loop and exchanges would be lost anyway. */
714         if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep &&
715             do_per_step(step, ir->swap->nstswap))
716         {
717             const bool doRerun = true;
718             do_swapcoords(cr, step, t, ir, swap, wcycle,
719                           rerun_fr.x,
720                           rerun_fr.box,
721                           MASTER(cr) && mdrunOptions.verbose,
722                           doRerun);
723         }
724
725         if (MASTER(cr))
726         {
727             /* read next frame from input trajectory */
728             isLastStep = !read_next_frame(oenv, status, &rerun_fr);
729         }
730
731         if (PAR(cr))
732         {
733             rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
734         }
735
736         cycles = wallcycle_stop(wcycle, ewcSTEP);
737         if (DOMAINDECOMP(cr) && wcycle)
738         {
739             dd_cycles_add(cr->dd, cycles, ddCyclStep);
740         }
741
742         if (!rerun_fr.bStep)
743         {
744             /* increase the MD step number */
745             step++;
746             step_rel++;
747         }
748     }
749     /* End of main MD loop */
750
751     /* Closing TNG files can include compressing data. Therefore it is good to do that
752      * before stopping the time measurements. */
753     mdoutf_tng_close(outf);
754
755     /* Stop measuring walltime */
756     walltime_accounting_end_time(walltime_accounting);
757
758     if (MASTER(cr))
759     {
760         close_trx(status);
761     }
762
763     if (!thisRankHasDuty(cr, DUTY_PME))
764     {
765         /* Tell the PME only node to finish */
766         gmx_pme_send_finish(cr);
767     }
768
769     done_mdoutf(outf);
770
771     done_shellfc(fplog, shellfc, step_rel);
772
773     walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
774 }