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