Improve mimic module
[gromacs.git] / src / gromacs / mdrun / rerun.cpp
blobb025d6ab2805d8bf5c86db15349c85390a2b97ad
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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/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/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/energyoutput.h"
80 #include "gromacs/mdlib/expanded.h"
81 #include "gromacs/mdlib/force.h"
82 #include "gromacs/mdlib/force_flags.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/mdsetup.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/ns.h"
90 #include "gromacs/mdlib/resethandler.h"
91 #include "gromacs/mdlib/shellfc.h"
92 #include "gromacs/mdlib/sighandler.h"
93 #include "gromacs/mdlib/sim_util.h"
94 #include "gromacs/mdlib/simulationsignal.h"
95 #include "gromacs/mdlib/stophandler.h"
96 #include "gromacs/mdlib/tgroup.h"
97 #include "gromacs/mdlib/trajectory_writing.h"
98 #include "gromacs/mdlib/update.h"
99 #include "gromacs/mdlib/vcm.h"
100 #include "gromacs/mdlib/vsite.h"
101 #include "gromacs/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"
136 #include "integrator.h"
137 #include "replicaexchange.h"
139 using gmx::SimulationSignaller;
141 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
143 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
144 * \param[in,out] globalState The global state container
145 * \param[in] constructVsites When true, vsite coordinates are constructed
146 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
147 * \param[in] idef Topology parameters, used for constructing vsites
148 * \param[in] timeStep Time step, used for constructing vsites
149 * \param[in] forceRec Force record, used for constructing vsites
150 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
152 static void prepareRerunState(const t_trxframe &rerunFrame,
153 t_state *globalState,
154 bool constructVsites,
155 const gmx_vsite_t *vsite,
156 const t_idef &idef,
157 double timeStep,
158 const t_forcerec &forceRec,
159 t_graph *graph)
161 auto x = makeArrayRef(globalState->x);
162 auto rerunX = arrayRefFromArray(reinterpret_cast<gmx::RVec *>(rerunFrame.x), globalState->natoms);
163 std::copy(rerunX.begin(), rerunX.end(), x.begin());
164 copy_mat(rerunFrame.box, globalState->box);
166 if (constructVsites)
168 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
170 if (graph)
172 /* Following is necessary because the graph may get out of sync
173 * with the coordinates if we only have every N'th coordinate set
175 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, globalState->x.rvec_array());
176 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
178 construct_vsites(vsite, globalState->x.rvec_array(), timeStep, globalState->v.rvec_array(),
179 idef.iparams, idef.il,
180 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
181 if (graph)
183 unshift_self(graph, globalState->box, globalState->x.rvec_array());
188 void gmx::Integrator::do_rerun()
190 // TODO Historically, the EM and MD "integrators" used different
191 // names for the t_inputrec *parameter, but these must have the
192 // same name, now that it's a member of a struct. We use this ir
193 // alias to avoid a large ripple of nearly useless changes.
194 // t_inputrec is being replaced by IMdpOptionsProvider, so this
195 // will go away eventually.
196 t_inputrec *ir = inputrec;
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 gmx_enerdata_t *enerd;
208 PaddedVector<gmx::RVec> f {};
209 gmx_global_stat_t gstat;
210 t_graph *graph = nullptr;
211 gmx_groups_t *groups;
212 gmx_shellfc_t *shellfc;
214 double cycles;
216 /* Domain decomposition could incorrectly miss a bonded
217 interaction, but checking for that requires a global
218 communication stage, which does not otherwise happen in DD
219 code. So we do that alongside the first global energy reduction
220 after a new DD is made. These variables handle whether the
221 check happens, and the result it returns. */
222 bool shouldCheckNumberOfBondedInteractions = false;
223 int totalNumberOfBondedInteractions = -1;
225 SimulationSignals signals;
226 // Most global communnication stages don't propagate mdrun
227 // signals, and will use this object to achieve that.
228 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
230 GMX_LOG(mdlog.info).asParagraph().
231 appendText("Note that it is planned that the command gmx mdrun -rerun will "
232 "be available in a different form in a future version of GROMACS, "
233 "e.g. gmx rerun -f.");
235 if (ir->efep != efepNO && (mdAtoms->mdatoms()->nMassPerturbed > 0 ||
236 (constr && constr->havePerturbedConstraints())))
238 gmx_fatal(FARGS, "Perturbed masses or constraints are not supported by rerun. "
239 "Either make a .tpr without mass and constraint perturbation, "
240 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
242 if (ir->bExpanded)
244 gmx_fatal(FARGS, "Expanded ensemble not supported by rerun.");
246 if (ir->bSimTemp)
248 gmx_fatal(FARGS, "Simulated tempering not supported by rerun.");
250 if (ir->bDoAwh)
252 gmx_fatal(FARGS, "AWH not supported by rerun.");
254 if (replExParams.exchangeInterval > 0)
256 gmx_fatal(FARGS, "Replica exchange not supported by rerun.");
258 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
260 gmx_fatal(FARGS, "Essential dynamics not supported by rerun.");
262 if (ir->bIMD)
264 gmx_fatal(FARGS, "Interactive MD not supported by rerun.");
266 if (isMultiSim(ms))
268 gmx_fatal(FARGS, "Multiple simulations not supported by rerun.");
270 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
271 [](int i){return i != eannNO; }))
273 gmx_fatal(FARGS, "Simulated annealing not supported by rerun.");
276 /* Rerun can't work if an output file name is the same as the input file name.
277 * If this is the case, the user will get an error telling them what the issue is.
279 if (strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-o", nfile, fnm)) == 0 ||
280 strcmp(opt2fn("-rerun", nfile, fnm), opt2fn("-x", nfile, fnm)) == 0)
282 gmx_fatal(FARGS, "When using mdrun -rerun, the name of the input trajectory file "
283 "%s cannot be identical to the name of an output file (whether "
284 "given explicitly with -o or -x, or by default)",
285 opt2fn("-rerun", nfile, fnm));
288 /* Settings for rerun */
289 ir->nstlist = 1;
290 ir->nstcalcenergy = 1;
291 int nstglobalcomm = 1;
292 const bool bNS = true;
294 ir->nstxout_compressed = 0;
295 groups = &top_global->groups;
296 if (ir->eI == eiMimic)
298 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
301 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
302 init_nrnb(nrnb);
303 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle);
304 gmx::EnergyOutput energyOutput;
305 energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, ir, mdoutf_get_fp_dhdl(outf), true);
307 /* Energy terms and groups */
308 snew(enerd, 1);
309 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
310 enerd);
312 /* Kinetic energy data */
313 std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
314 gmx_ekindata_t *ekind = eKinData.get();
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;
319 gstat = global_stat_init(ir);
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));
327 double io = compute_io(ir, top_global->natoms, groups, energyOutput.numEnergyTerms(), 1);
328 if ((io > 2000) && MASTER(cr))
330 fprintf(stderr,
331 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
332 io);
336 // Local state only becomes valid now.
337 std::unique_ptr<t_state> stateInstance;
338 t_state * state;
340 if (DOMAINDECOMP(cr))
342 dd_init_local_top(*top_global, &top);
344 stateInstance = std::make_unique<t_state>();
345 state = stateInstance.get();
346 dd_init_local_state(cr->dd, state_global, state);
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;
356 else
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.
362 f.resizeWithPadding(state_global->natoms);
363 /* Copy the pointer to the global state */
364 state = state_global;
366 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
367 &graph, mdAtoms, constr, vsite, shellfc);
370 auto mdatoms = mdAtoms->mdatoms();
372 // NOTE: The global state is no longer used at this point.
373 // But state_global is still used as temporary storage space for writing
374 // the global state to file and potentially for replica exchange.
375 // (Global topology should persist.)
377 update_mdatoms(mdatoms, state->lambda[efptMASS]);
379 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
381 doFreeEnergyPerturbation = true;
385 int cglo_flags = (CGLO_INITIALIZATION | CGLO_GSTAT |
386 (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
387 bool bSumEkinhOld = false;
388 t_vcm *vcm = nullptr;
389 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
390 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
391 constr, &nullSignaller, state->box,
392 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
394 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
395 top_global, &top, state,
396 &shouldCheckNumberOfBondedInteractions);
398 if (MASTER(cr))
400 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
401 " input trajectory '%s'\n\n",
402 *(top_global->name), opt2fn("-rerun", nfile, fnm));
403 if (mdrunOptions.verbose)
405 fprintf(stderr, "Calculated time to finish depends on nsteps from "
406 "run input file,\nwhich may not correspond to the time "
407 "needed to process input trajectory.\n\n");
409 fprintf(fplog, "\n");
412 walltime_accounting_start_time(walltime_accounting);
413 wallcycle_start(wcycle, ewcRUN);
414 print_start(fplog, cr, walltime_accounting, "mdrun");
416 /***********************************************************
418 * Loop over MD steps
420 ************************************************************/
422 if (constr)
424 GMX_LOG(mdlog.info).asParagraph().
425 appendText("Simulations has constraints. Rerun does not recalculate constraints.");
428 rerun_fr.natoms = 0;
429 if (MASTER(cr))
431 isLastStep = !read_first_frame(oenv, &status,
432 opt2fn("-rerun", nfile, fnm),
433 &rerun_fr, TRX_NEED_X);
434 if (rerun_fr.natoms != top_global->natoms)
436 gmx_fatal(FARGS,
437 "Number of atoms in trajectory (%d) does not match the "
438 "run input file (%d)\n",
439 rerun_fr.natoms, top_global->natoms);
442 if (ir->ePBC != epbcNONE)
444 if (!rerun_fr.bBox)
446 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
447 "does not contain a box, while pbc is used",
448 rerun_fr.step, rerun_fr.time);
450 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
452 gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f "
453 "has too small box dimensions", rerun_fr.step, rerun_fr.time);
458 GMX_LOG(mdlog.info).asParagraph().
459 appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
461 if (PAR(cr))
463 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
466 if (ir->ePBC != epbcNONE)
468 /* Set the shift vectors.
469 * Necessary here when have a static box different from the tpr box.
471 calc_shifts(rerun_fr.box, fr->shift_vec);
474 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
475 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
476 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
477 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
479 // we don't do counter resetting in rerun - finish will always be valid
480 walltime_accounting_set_valid_finish(walltime_accounting);
482 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
484 step = ir->init_step;
485 step_rel = 0;
487 /* and stop now if we should */
488 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
489 while (!isLastStep)
491 wallcycle_start(wcycle, ewcSTEP);
493 if (rerun_fr.bStep)
495 step = rerun_fr.step;
496 step_rel = step - ir->init_step;
498 if (rerun_fr.bTime)
500 t = rerun_fr.time;
502 else
504 t = step;
507 if (ir->efep != efepNO && MASTER(cr))
509 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
512 if (MASTER(cr))
514 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
515 if (constructVsites && DOMAINDECOMP(cr))
517 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
518 "use a single rank");
520 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top.idef, ir->delta_t, *fr, graph);
523 isLastStep = isLastStep || stopHandler->stoppingAfterCurrentStep(bNS);
525 if (DOMAINDECOMP(cr))
527 /* Repartition the domain decomposition */
528 const bool bMasterState = true;
529 dd_partition_system(fplog, mdlog, step, cr,
530 bMasterState, nstglobalcomm,
531 state_global, *top_global, ir,
532 state, &f, mdAtoms, &top, fr,
533 vsite, constr,
534 nrnb, wcycle,
535 mdrunOptions.verbose);
536 shouldCheckNumberOfBondedInteractions = true;
539 if (MASTER(cr))
541 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
544 if (ir->efep != efepNO)
546 update_mdatoms(mdatoms, state->lambda[efptMASS]);
549 force_flags = (GMX_FORCE_STATECHANGED |
550 GMX_FORCE_DYNAMICBOX |
551 GMX_FORCE_ALLFORCES |
552 (GMX_GPU ? GMX_FORCE_VIRIAL : 0) | // TODO: Get rid of this once #2649 is solved
553 GMX_FORCE_ENERGY |
554 (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
556 if (shellfc)
558 /* Now is the time to relax the shells */
559 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
560 enforcedRotation, step,
561 ir, bNS, force_flags, &top,
562 constr, enerd, fcd,
563 state, f.arrayRefWithPadding(), force_vir, mdatoms,
564 nrnb, wcycle, graph, groups,
565 shellfc, fr, ppForceWorkload, t, mu_tot,
566 vsite,
567 ddBalanceRegionHandler);
569 else
571 /* The coordinates (x) are shifted (to get whole molecules)
572 * in do_force.
573 * This is parallellized as well, and does communication too.
574 * Check comments in sim_util.c
576 Awh *awh = nullptr;
577 gmx_edsam *ed = nullptr;
578 do_force(fplog, cr, ms, ir, awh, enforcedRotation,
579 step, nrnb, wcycle, &top, groups,
580 state->box, state->x.arrayRefWithPadding(), &state->hist,
581 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
582 state->lambda, graph,
583 fr, ppForceWorkload, vsite, mu_tot, t, ed,
584 GMX_FORCE_NS | force_flags,
585 ddBalanceRegionHandler);
588 /* Now we have the energies and forces corresponding to the
589 * coordinates at time t.
592 const bool isCheckpointingStep = false;
593 const bool doRerun = true;
594 const bool bSumEkinhOld = false;
595 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
596 ir, state, state_global, observablesHistory,
597 top_global, fr,
598 outf, energyOutput, ekind, f,
599 isCheckpointingStep, doRerun, isLastStep,
600 mdrunOptions.writeConfout,
601 bSumEkinhOld);
604 stopHandler->setSignal();
606 if (graph)
608 /* Need to unshift here */
609 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
612 if (vsite != nullptr)
614 wallcycle_start(wcycle, ewcVSITECONSTR);
615 if (graph != nullptr)
617 shift_self(graph, state->box, as_rvec_array(state->x.data()));
619 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
620 top.idef.iparams, top.idef.il,
621 fr->ePBC, fr->bMolPBC, cr, state->box);
623 if (graph != nullptr)
625 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
627 wallcycle_stop(wcycle, ewcVSITECONSTR);
631 const bool doInterSimSignal = false;
632 const bool doIntraSimSignal = true;
633 bool bSumEkinhOld = false;
634 t_vcm *vcm = nullptr;
635 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
637 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
638 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
639 constr, &signaller,
640 state->box,
641 &totalNumberOfBondedInteractions, &bSumEkinhOld,
642 CGLO_GSTAT | CGLO_ENERGY
643 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
645 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
646 top_global, &top, state,
647 &shouldCheckNumberOfBondedInteractions);
650 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
651 the virial that should probably be addressed eventually. state->veta has better properies,
652 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
653 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
655 if (ir->efep != efepNO)
657 /* Sum up the foreign energy and dhdl terms for md and sd.
658 Currently done every step so that dhdl is correct in the .edr */
659 sum_dhdl(enerd, state->lambda, ir->fepvals);
662 /* Output stuff */
663 if (MASTER(cr))
665 const bool bCalcEnerStep = true;
666 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep,
667 t, mdatoms->tmass, enerd, state,
668 ir->fepvals, ir->expandedvals, state->box,
669 shake_vir, force_vir, total_vir, pres,
670 ekind, mu_tot, constr);
672 const bool do_ene = true;
673 const bool do_log = true;
674 Awh *awh = nullptr;
675 const bool do_dr = ir->nstdisreout != 0;
676 const bool do_or = ir->nstorireout != 0;
678 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
679 do_log ? fplog : nullptr,
680 step, t,
681 eprNORMAL, fcd, groups, &(ir->opts), awh);
683 if (do_per_step(step, ir->nstlog))
685 if (fflush(fplog) != 0)
687 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
692 /* Print the remaining wall clock time for the run */
693 if (isMasterSimMasterRank(ms, cr) &&
694 (mdrunOptions.verbose || gmx_got_usr_signal()))
696 if (shellfc)
698 fprintf(stderr, "\n");
700 print_time(stderr, walltime_accounting, step, ir, cr);
703 /* Ion/water position swapping.
704 * Not done in last step since trajectory writing happens before this call
705 * in the MD loop and exchanges would be lost anyway. */
706 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !isLastStep &&
707 do_per_step(step, ir->swap->nstswap))
709 const bool doRerun = true;
710 do_swapcoords(cr, step, t, ir, wcycle,
711 rerun_fr.x,
712 rerun_fr.box,
713 MASTER(cr) && mdrunOptions.verbose,
714 doRerun);
717 if (MASTER(cr))
719 /* read next frame from input trajectory */
720 isLastStep = !read_next_frame(oenv, status, &rerun_fr);
723 if (PAR(cr))
725 rerun_parallel_comm(cr, &rerun_fr, &isLastStep);
728 cycles = wallcycle_stop(wcycle, ewcSTEP);
729 if (DOMAINDECOMP(cr) && wcycle)
731 dd_cycles_add(cr->dd, cycles, ddCyclStep);
734 if (!rerun_fr.bStep)
736 /* increase the MD step number */
737 step++;
738 step_rel++;
741 /* End of main MD loop */
743 /* Closing TNG files can include compressing data. Therefore it is good to do that
744 * before stopping the time measurements. */
745 mdoutf_tng_close(outf);
747 /* Stop measuring walltime */
748 walltime_accounting_end_time(walltime_accounting);
750 if (MASTER(cr))
752 close_trx(status);
755 if (!thisRankHasDuty(cr, DUTY_PME))
757 /* Tell the PME only node to finish */
758 gmx_pme_send_finish(cr);
761 done_mdoutf(outf);
763 done_shellfc(fplog, shellfc, step_rel);
765 // Clean up swapcoords
766 if (ir->eSwapCoords != eswapNO)
768 finish_swapcoords(ir->swap);
771 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
773 destroy_enerdata(enerd);
774 sfree(enerd);