Introduce PpForceWorkload
[gromacs.git] / src / gromacs / mdrun / rerun.cpp
blob4a29fc9b97cb8d2eecbed5109e3ee6c5fb1dfc07
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
37 * \brief Implements the loop for simulation reruns
39 * \author Pascal Merz <pascal.merz@colorado.edu>
40 * \ingroup module_mdrun
42 #include "gmxpre.h"
44 #include "config.h"
46 #include <cinttypes>
47 #include <cmath>
48 #include <cstdio>
49 #include <cstdlib>
51 #include <algorithm>
52 #include <memory>
54 #include "gromacs/awh/awh.h"
55 #include "gromacs/commandline/filenm.h"
56 #include "gromacs/compat/make_unique.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_network.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/partition.h"
62 #include "gromacs/essentialdynamics/edsam.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/ewald/pme-load-balancing.h"
65 #include "gromacs/fileio/trxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/gpu_utils/gpu_utils.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/listed-forces/manage-threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/utilities.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/math/vectypes.h"
75 #include "gromacs/mdlib/checkpointhandler.h"
76 #include "gromacs/mdlib/compute_io.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/ebin.h"
79 #include "gromacs/mdlib/expanded.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdebin.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/mdrun.h"
88 #include "gromacs/mdlib/mdsetup.h"
89 #include "gromacs/mdlib/membed.h"
90 #include "gromacs/mdlib/nb_verlet.h"
91 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
92 #include "gromacs/mdlib/ns.h"
93 #include "gromacs/mdlib/resethandler.h"
94 #include "gromacs/mdlib/shellfc.h"
95 #include "gromacs/mdlib/sighandler.h"
96 #include "gromacs/mdlib/sim_util.h"
97 #include "gromacs/mdlib/simulationsignal.h"
98 #include "gromacs/mdlib/stophandler.h"
99 #include "gromacs/mdlib/tgroup.h"
100 #include "gromacs/mdlib/trajectory_writing.h"
101 #include "gromacs/mdlib/update.h"
102 #include "gromacs/mdlib/vcm.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdtypes/awh-history.h"
105 #include "gromacs/mdtypes/awh-params.h"
106 #include "gromacs/mdtypes/commrec.h"
107 #include "gromacs/mdtypes/df_history.h"
108 #include "gromacs/mdtypes/energyhistory.h"
109 #include "gromacs/mdtypes/fcdata.h"
110 #include "gromacs/mdtypes/forcerec.h"
111 #include "gromacs/mdtypes/group.h"
112 #include "gromacs/mdtypes/inputrec.h"
113 #include "gromacs/mdtypes/interaction_const.h"
114 #include "gromacs/mdtypes/md_enums.h"
115 #include "gromacs/mdtypes/mdatom.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/state.h"
118 #include "gromacs/mimic/MimicUtils.h"
119 #include "gromacs/pbcutil/mshift.h"
120 #include "gromacs/pbcutil/pbc.h"
121 #include "gromacs/pulling/pull.h"
122 #include "gromacs/swap/swapcoords.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/timing/walltime_accounting.h"
125 #include "gromacs/topology/atoms.h"
126 #include "gromacs/topology/idef.h"
127 #include "gromacs/topology/mtop_util.h"
128 #include "gromacs/topology/topology.h"
129 #include "gromacs/trajectory/trajectoryframe.h"
130 #include "gromacs/utility/basedefinitions.h"
131 #include "gromacs/utility/cstringutil.h"
132 #include "gromacs/utility/fatalerror.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/real.h"
135 #include "gromacs/utility/smalloc.h"
137 #include "integrator.h"
138 #include "replicaexchange.h"
140 using gmx::SimulationSignaller;
142 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
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
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)
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);
167 if (constructVsites)
169 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
171 if (graph)
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
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()));
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)
184 unshift_self(graph, globalState->box, globalState->x.rvec_array());
189 void gmx::Integrator::do_rerun()
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 gmx_mdoutf *outf = nullptr;
199 int64_t step, step_rel;
200 double t, lam0[efptNR];
201 bool isLastStep = false;
202 bool doFreeEnergyPerturbation = false;
203 int force_flags;
204 tensor force_vir, shake_vir, total_vir, pres;
205 t_trxstatus *status;
206 rvec mu_tot;
207 t_trxframe rerun_fr;
208 gmx_localtop_t *top;
209 t_mdebin *mdebin = nullptr;
210 gmx_enerdata_t *enerd;
211 PaddedVector<gmx::RVec> f {};
212 gmx_global_stat_t gstat;
213 t_graph *graph = nullptr;
214 gmx_groups_t *groups;
215 gmx_ekindata_t *ekind;
216 gmx_shellfc_t *shellfc;
218 double cycles;
220 /* Domain decomposition could incorrectly miss a bonded
221 interaction, but checking for that requires a global
222 communication stage, which does not otherwise happen in DD
223 code. So we do that alongside the first global energy reduction
224 after a new DD is made. These variables handle whether the
225 check happens, and the result it returns. */
226 bool shouldCheckNumberOfBondedInteractions = false;
227 int totalNumberOfBondedInteractions = -1;
229 SimulationSignals signals;
230 // Most global communnication stages don't propagate mdrun
231 // signals, and will use this object to achieve that.
232 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
234 GMX_LOG(mdlog.info).asParagraph().
235 appendText("Note that it is planned that the command gmx mdrun -rerun will "
236 "be available in a different form in a future version of GROMACS, "
237 "e.g. gmx rerun -f.");
239 if (ir->bExpanded)
241 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
243 if (ir->bSimTemp)
245 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
247 if (ir->bDoAwh)
249 gmx_fatal(FARGS, "AWH not supported by rerun.");
251 if (replExParams.exchangeInterval > 0)
253 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
255 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
257 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
259 if (ir->bIMD)
261 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
263 if (isMultiSim(ms))
265 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
267 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
268 [](int i){return i != eannNO; }))
270 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
273 /* Rerun can't work if an output file name is the same as the input file name.
274 * If this is the case, the user will get an error telling them what the issue is.
276 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0 ||
277 strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
279 gmx_fatal(FARGS, "When using mdrun -rerun, the name of the input trajectory file "
280 "%s cannot be identical to the name of an output file (whether "
281 "given explicitly with -o or -x, or by default)",
282 opt2fn("-rerun", nfile, fnm));
285 /* Settings for rerun */
286 ir->nstlist = 1;
287 ir->nstcalcenergy = 1;
288 int nstglobalcomm = 1;
289 const bool bNS = true;
291 ir->nstxout_compressed = 0;
292 groups = &top_global->groups;
293 if (ir->eI == eiMimic)
295 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
298 /* Initial values */
299 init_rerun(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
300 state_global, lam0, nrnb, top_global,
301 nfile, fnm, &outf, &mdebin, wcycle);
303 /* Energy terms and groups */
304 snew(enerd, 1);
305 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
306 enerd);
308 /* Kinetic energy data */
309 snew(ekind, 1);
310 init_ekindata(fplog, top_global, &(ir->opts), ekind);
311 /* Copy the cos acceleration to the groups struct */
312 ekind->cosacc.cos_accel = ir->cos_accel;
314 gstat = global_stat_init(ir);
316 /* Check for polarizable models and flexible constraints */
317 shellfc = init_shell_flexcon(fplog,
318 top_global, constr ? constr->numFlexibleConstraints() : 0,
319 ir->nstcalcenergy, DOMAINDECOMP(cr));
322 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
323 if ((io > 2000) && MASTER(cr))
325 fprintf(stderr,
326 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
327 io);
331 // Local state only becomes valid now.
332 std::unique_ptr<t_state> stateInstance;
333 t_state * state;
335 if (DOMAINDECOMP(cr))
337 top = dd_init_local_top(top_global);
339 stateInstance = compat::make_unique<t_state>();
340 state = stateInstance.get();
341 dd_init_local_state(cr->dd, state_global, state);
343 /* Distribute the charge groups over the nodes from the master node */
344 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
345 state_global, top_global, ir,
346 state, &f, mdAtoms, top, fr,
347 vsite, constr,
348 nrnb, nullptr, FALSE);
349 shouldCheckNumberOfBondedInteractions = true;
351 else
353 state_change_natoms(state_global, state_global->natoms);
354 /* We need to allocate one element extra, since we might use
355 * (unaligned) 4-wide SIMD loads to access rvec entries.
357 f.resizeWithPadding(state_global->natoms);
358 /* Copy the pointer to the global state */
359 state = state_global;
361 snew(top, 1);
362 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
363 &graph, mdAtoms, constr, vsite, shellfc);
366 auto mdatoms = mdAtoms->mdatoms();
368 // NOTE: The global state is no longer used at this point.
369 // But state_global is still used as temporary storage space for writing
370 // the global state to file and potentially for replica exchange.
371 // (Global topology should persist.)
373 update_mdatoms(mdatoms, state->lambda[efptMASS]);
375 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
377 doFreeEnergyPerturbation = true;
381 int cglo_flags = (CGLO_INITIALIZATION | CGLO_GSTAT |
382 (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
383 bool bSumEkinhOld = false;
384 t_vcm *vcm = nullptr;
385 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, 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);
390 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
391 top_global, top, state,
392 &shouldCheckNumberOfBondedInteractions);
394 if (MASTER(cr))
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)
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");
405 fprintf(fplog, "\n");
408 walltime_accounting_start_time(walltime_accounting);
409 wallcycle_start(wcycle, ewcRUN);
410 print_start(fplog, cr, walltime_accounting, "mdrun");
412 /***********************************************************
414 * Loop over MD steps
416 ************************************************************/
418 if (constr)
420 GMX_LOG(mdlog.info).asParagraph().
421 appendText("Simulations has constraints. Rerun does not recalculate constraints.");
424 rerun_fr.natoms = 0;
425 if (MASTER(cr))
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)
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);
438 if (ir->ePBC != epbcNONE)
440 if (!rerun_fr.bBox)
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);
446 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
448 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
449 "has too small box dimensions", rerun_fr.step, rerun_fr.time);
454 GMX_LOG(mdlog.info).asParagraph().
455 appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
457 if (PAR(cr))
459 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
462 if (ir->ePBC != epbcNONE)
464 /* Set the shift vectors.
465 * Necessary here when have a static box different from the tpr box.
467 calc_shifts(rerun_fr.box, fr->shift_vec);
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);
475 // we don't do counter resetting in rerun - finish will always be valid
476 walltime_accounting_set_valid_finish(walltime_accounting);
478 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
479 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
481 step = ir->init_step;
482 step_rel = 0;
484 /* and stop now if we should */
485 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
486 while (!isLastStep)
488 wallcycle_start(wcycle, ewcSTEP);
490 if (rerun_fr.bStep)
492 step = rerun_fr.step;
493 step_rel = step - ir->init_step;
495 if (rerun_fr.bTime)
497 t = rerun_fr.time;
499 else
501 t = step;
504 if (ir->efep != efepNO && MASTER(cr))
506 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
509 if (MASTER(cr))
511 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
512 if (constructVsites && DOMAINDECOMP(cr))
514 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
515 "use a single rank");
517 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph);
520 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
522 if (DOMAINDECOMP(cr))
524 /* Repartition the domain decomposition */
525 const bool bMasterState = true;
526 dd_partition_system(fplog, mdlog, step, cr,
527 bMasterState, nstglobalcomm,
528 state_global, top_global, ir,
529 state, &f, mdAtoms, top, fr,
530 vsite, constr,
531 nrnb, wcycle,
532 mdrunOptions.verbose);
533 shouldCheckNumberOfBondedInteractions = true;
536 if (MASTER(cr))
538 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
541 if (ir->efep != efepNO)
543 update_mdatoms(mdatoms, state->lambda[efptMASS]);
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));
553 if (shellfc)
555 /* Now is the time to relax the shells */
556 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
557 enforcedRotation, step,
558 ir, bNS, force_flags, top,
559 constr, enerd, fcd,
560 state, f.arrayRefWithPadding(), force_vir, mdatoms,
561 nrnb, wcycle, graph, groups,
562 shellfc, fr, ppForceWorkload, t, mu_tot,
563 vsite,
564 ddOpenBalanceRegion, ddCloseBalanceRegion);
566 else
568 /* The coordinates (x) are shifted (to get whole molecules)
569 * in do_force.
570 * This is parallellized as well, and does communication too.
571 * Check comments in sim_util.c
573 Awh *awh = nullptr;
574 gmx_edsam *ed = nullptr;
575 do_force(fplog, cr, ms, ir, awh, enforcedRotation,
576 step, nrnb, wcycle, top, groups,
577 state->box, state->x.arrayRefWithPadding(), &state->hist,
578 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
579 state->lambda, graph,
580 fr, ppForceWorkload, vsite, mu_tot, t, ed,
581 GMX_FORCE_NS | force_flags,
582 ddOpenBalanceRegion, ddCloseBalanceRegion);
585 /* Now we have the energies and forces corresponding to the
586 * coordinates at time t.
589 const bool isCheckpointingStep = false;
590 const bool doRerun = true;
591 const bool bSumEkinhOld = false;
592 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
593 ir, state, state_global, observablesHistory,
594 top_global, fr,
595 outf, mdebin, ekind, f,
596 isCheckpointingStep, doRerun, isLastStep,
597 mdrunOptions.writeConfout,
598 bSumEkinhOld);
601 stopHandler->setSignal();
603 if (graph)
605 /* Need to unshift here */
606 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
609 if (vsite != nullptr)
611 wallcycle_start(wcycle, ewcVSITECONSTR);
612 if (graph != nullptr)
614 shift_self(graph, state->box, as_rvec_array(state->x.data()));
616 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
617 top->idef.iparams, top->idef.il,
618 fr->ePBC, fr->bMolPBC, cr, state->box);
620 if (graph != nullptr)
622 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
624 wallcycle_stop(wcycle, ewcVSITECONSTR);
628 const bool doInterSimSignal = false;
629 const bool doIntraSimSignal = true;
630 bool bSumEkinhOld = false;
631 t_vcm *vcm = nullptr;
632 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
634 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
635 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
636 constr, &signaller,
637 state->box,
638 &totalNumberOfBondedInteractions, &bSumEkinhOld,
639 CGLO_GSTAT | CGLO_ENERGY
640 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
642 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
643 top_global, top, state,
644 &shouldCheckNumberOfBondedInteractions);
647 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
648 the virial that should probably be addressed eventually. state->veta has better properies,
649 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
650 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
652 if (ir->efep != efepNO)
654 /* Sum up the foreign energy and dhdl terms for md and sd.
655 Currently done every step so that dhdl is correct in the .edr */
656 sum_dhdl(enerd, state->lambda, ir->fepvals);
659 /* Output stuff */
660 if (MASTER(cr))
662 const bool bCalcEnerStep = true;
663 upd_mdebin(mdebin, doFreeEnergyPerturbation, bCalcEnerStep,
664 t, mdatoms->tmass, enerd, state,
665 ir->fepvals, ir->expandedvals, rerun_fr.box,
666 shake_vir, force_vir, total_vir, pres,
667 ekind, mu_tot, constr);
669 const bool do_ene = true;
670 const bool do_log = true;
671 Awh *awh = nullptr;
672 const bool do_dr = ir->nstdisreout != 0;
673 const bool do_or = ir->nstorireout != 0;
675 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
676 step, t,
677 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh);
679 if (do_per_step(step, ir->nstlog))
681 if (fflush(fplog) != 0)
683 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
688 /* Print the remaining wall clock time for the run */
689 if (isMasterSimMasterRank(ms, cr) &&
690 (mdrunOptions.verbose || gmx_got_usr_signal()))
692 if (shellfc)
694 fprintf(stderr, "\n");
696 print_time(stderr, walltime_accounting, step, ir, cr);
699 /* Ion/water position swapping.
700 * Not done in last step since trajectory writing happens before this call
701 * in the MD loop and exchanges would be lost anyway. */
702 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep &&
703 do_per_step(step, ir->swap->nstswap))
705 const bool doRerun = true;
706 do_swapcoords(cr, step, t, ir, wcycle,
707 rerun_fr.x,
708 rerun_fr.box,
709 MASTER(cr) && mdrunOptions.verbose,
710 doRerun);
713 if (MASTER(cr))
715 /* read next frame from input trajectory */
716 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
719 if (PAR(cr))
721 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
724 cycles = wallcycle_stop(wcycle, ewcSTEP);
725 if (DOMAINDECOMP(cr) && wcycle)
727 dd_cycles_add(cr->dd, cycles, ddCyclStep);
730 if (!rerun_fr.bStep)
732 /* increase the MD step number */
733 step++;
734 step_rel++;
737 /* End of main MD loop */
739 /* Closing TNG files can include compressing data. Therefore it is good to do that
740 * before stopping the time measurements. */
741 mdoutf_tng_close(outf);
743 /* Stop measuring walltime */
744 walltime_accounting_end_time(walltime_accounting);
746 if (MASTER(cr))
748 close_trx(status);
751 if (!thisRankHasDuty(cr, DUTY_PME))
753 /* Tell the PME only node to finish */
754 gmx_pme_send_finish(cr);
757 done_mdebin(mdebin);
758 done_mdoutf(outf);
760 done_shellfc(fplog, shellfc, step_rel);
762 // Clean up swapcoords
763 if (ir->eSwapCoords != eswapNO)
765 finish_swapcoords(ir->swap);
768 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
770 destroy_enerdata(enerd);
771 sfree(enerd);
772 sfree(top);