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.
37 * \brief Implements the loop for simulation reruns
39 * \author Pascal Merz <pascal.merz@colorado.edu>
40 * \ingroup module_mdrun
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/simulationsignal.h"
94 #include "gromacs/mdlib/stat.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
,
158 const t_forcerec
&forceRec
,
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
);
168 GMX_ASSERT(vsite
, "Need valid vsite for constructing vsites");
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
);
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;
202 tensor force_vir
, shake_vir
, total_vir
, pres
;
207 PaddedVector
<gmx::RVec
> f
{};
208 gmx_global_stat_t gstat
;
209 t_graph
*graph
= nullptr;
210 gmx_shellfc_t
*shellfc
;
214 /* Domain decomposition could incorrectly miss a bonded
215 interaction, but checking for that requires a global
216 communication stage, which does not otherwise happen in DD
217 code. So we do that alongside the first global energy reduction
218 after a new DD is made. These variables handle whether the
219 check happens, and the result it returns. */
220 bool shouldCheckNumberOfBondedInteractions
= false;
221 int totalNumberOfBondedInteractions
= -1;
223 SimulationSignals signals
;
224 // Most global communnication stages don't propagate mdrun
225 // signals, and will use this object to achieve that.
226 SimulationSignaller
nullSignaller(nullptr, nullptr, nullptr, false, false);
228 GMX_LOG(mdlog
.info
).asParagraph().
229 appendText("Note that it is planned that the command gmx mdrun -rerun will "
230 "be available in a different form in a future version of GROMACS, "
231 "e.g. gmx rerun -f.");
233 if (ir
->efep
!= efepNO
&& (mdAtoms
->mdatoms()->nMassPerturbed
> 0 ||
234 (constr
&& constr
->havePerturbedConstraints())))
236 gmx_fatal(FARGS
, "Perturbed masses or constraints are not supported by rerun. "
237 "Either make a .tpr without mass and constraint perturbation, "
238 "or use GROMACS 2018.4, 2018.5 or later 2018 version.");
242 gmx_fatal(FARGS
, "Expanded ensemble not supported by rerun.");
246 gmx_fatal(FARGS
, "Simulated tempering not supported by rerun.");
250 gmx_fatal(FARGS
, "AWH not supported by rerun.");
252 if (replExParams
.exchangeInterval
> 0)
254 gmx_fatal(FARGS
, "Replica exchange not supported by rerun.");
256 if (opt2bSet("-ei", nfile
, fnm
) || observablesHistory
->edsamHistory
!= nullptr)
258 gmx_fatal(FARGS
, "Essential dynamics not supported by rerun.");
262 gmx_fatal(FARGS
, "Interactive MD not supported by rerun.");
266 gmx_fatal(FARGS
, "Multiple simulations not supported by rerun.");
268 if (std::any_of(ir
->opts
.annealing
, ir
->opts
.annealing
+ ir
->opts
.ngtc
,
269 [](int i
){return i
!= eannNO
; }))
271 gmx_fatal(FARGS
, "Simulated annealing not supported by rerun.");
274 /* Rerun can't work if an output file name is the same as the input file name.
275 * If this is the case, the user will get an error telling them what the issue is.
277 if (strcmp(opt2fn("-rerun", nfile
, fnm
), opt2fn("-o", nfile
, fnm
)) == 0 ||
278 strcmp(opt2fn("-rerun", nfile
, fnm
), opt2fn("-x", nfile
, fnm
)) == 0)
280 gmx_fatal(FARGS
, "When using mdrun -rerun, the name of the input trajectory file "
281 "%s cannot be identical to the name of an output file (whether "
282 "given explicitly with -o or -x, or by default)",
283 opt2fn("-rerun", nfile
, fnm
));
286 /* Settings for rerun */
288 ir
->nstcalcenergy
= 1;
289 int nstglobalcomm
= 1;
290 const bool bNS
= true;
292 ir
->nstxout_compressed
= 0;
293 SimulationGroups
*groups
= &top_global
->groups
;
294 if (ir
->eI
== eiMimic
)
296 top_global
->intermolecularExclusionGroup
= genQmmmIndices(*top_global
);
299 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 gmx::EnergyOutput energyOutput
;
303 energyOutput
.prepare(mdoutf_get_fp_ene(outf
), top_global
, ir
, mdoutf_get_fp_dhdl(outf
), true);
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
;
312 gstat
= global_stat_init(ir
);
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
));
320 double io
= compute_io(ir
, top_global
->natoms
, *groups
, energyOutput
.numEnergyTerms(), 1);
321 if ((io
> 2000) && MASTER(cr
))
324 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
329 // Local state only becomes valid now.
330 std::unique_ptr
<t_state
> stateInstance
;
333 if (DOMAINDECOMP(cr
))
335 dd_init_local_top(*top_global
, &top
);
337 stateInstance
= std::make_unique
<t_state
>();
338 state
= stateInstance
.get();
339 dd_init_local_state(cr
->dd
, state_global
, state
);
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
,
344 state
, &f
, mdAtoms
, &top
, fr
,
346 nrnb
, nullptr, FALSE
);
347 shouldCheckNumberOfBondedInteractions
= true;
351 state_change_natoms(state_global
, state_global
->natoms
);
352 /* We need to allocate one element extra, since we might use
353 * (unaligned) 4-wide SIMD loads to access rvec entries.
355 f
.resizeWithPadding(state_global
->natoms
);
356 /* Copy the pointer to the global state */
357 state
= state_global
;
359 mdAlgorithmsSetupAtomData(cr
, ir
, *top_global
, &top
, fr
,
360 &graph
, mdAtoms
, constr
, vsite
, shellfc
);
363 auto mdatoms
= mdAtoms
->mdatoms();
365 // NOTE: The global state is no longer used at this point.
366 // But state_global is still used as temporary storage space for writing
367 // the global state to file and potentially for replica exchange.
368 // (Global topology should persist.)
370 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
372 if (ir
->efep
!= efepNO
&& ir
->fepvals
->nstdhdl
!= 0)
374 doFreeEnergyPerturbation
= true;
378 int cglo_flags
= (CGLO_INITIALIZATION
| CGLO_GSTAT
|
379 (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
380 bool bSumEkinhOld
= false;
381 t_vcm
*vcm
= nullptr;
382 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
383 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
384 constr
, &nullSignaller
, state
->box
,
385 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags
);
387 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
388 top_global
, &top
, state
,
389 &shouldCheckNumberOfBondedInteractions
);
393 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
394 " input trajectory '%s'\n\n",
395 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
396 if (mdrunOptions
.verbose
)
398 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
399 "run input file,\nwhich may not correspond to the time "
400 "needed to process input trajectory.\n\n");
402 fprintf(fplog
, "\n");
405 walltime_accounting_start_time(walltime_accounting
);
406 wallcycle_start(wcycle
, ewcRUN
);
407 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
409 /***********************************************************
413 ************************************************************/
417 GMX_LOG(mdlog
.info
).asParagraph().
418 appendText("Simulations has constraints. Rerun does not recalculate constraints.");
424 isLastStep
= !read_first_frame(oenv
, &status
,
425 opt2fn("-rerun", nfile
, fnm
),
426 &rerun_fr
, TRX_NEED_X
);
427 if (rerun_fr
.natoms
!= top_global
->natoms
)
430 "Number of atoms in trajectory (%d) does not match the "
431 "run input file (%d)\n",
432 rerun_fr
.natoms
, top_global
->natoms
);
435 if (ir
->ePBC
!= epbcNONE
)
439 gmx_fatal(FARGS
, "Rerun trajectory frame step %" PRId64
" time %f "
440 "does not contain a box, while pbc is used",
441 rerun_fr
.step
, rerun_fr
.time
);
443 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < gmx::square(fr
->rlist
))
445 gmx_fatal(FARGS
, "Rerun trajectory frame step %" PRId64
" time %f "
446 "has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
451 GMX_LOG(mdlog
.info
).asParagraph().
452 appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
456 rerun_parallel_comm(cr
, &rerun_fr
, &isLastStep
);
459 if (ir
->ePBC
!= epbcNONE
)
461 /* Set the shift vectors.
462 * Necessary here when have a static box different from the tpr box.
464 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
467 auto stopHandler
= stopHandlerBuilder
->getStopHandlerMD(
468 compat::not_null
<SimulationSignal
*>(&signals
[eglsSTOPCOND
]), false,
469 MASTER(cr
), ir
->nstlist
, mdrunOptions
.reproducible
, nstglobalcomm
,
470 mdrunOptions
.maximumHoursToRun
, ir
->nstlist
== 0, fplog
, step
, bNS
, walltime_accounting
);
472 // we don't do counter resetting in rerun - finish will always be valid
473 walltime_accounting_set_valid_finish(walltime_accounting
);
475 const DDBalanceRegionHandler
ddBalanceRegionHandler(cr
);
477 step
= ir
->init_step
;
480 /* and stop now if we should */
481 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
484 wallcycle_start(wcycle
, ewcSTEP
);
488 step
= rerun_fr
.step
;
489 step_rel
= step
- ir
->init_step
;
500 if (ir
->efep
!= efepNO
&& MASTER(cr
))
502 setCurrentLambdasRerun(step
, ir
->fepvals
, &rerun_fr
, lam0
, state_global
);
507 const bool constructVsites
= ((vsite
!= nullptr) && mdrunOptions
.rerunConstructVsites
);
508 if (constructVsites
&& DOMAINDECOMP(cr
))
510 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
511 "use a single rank");
513 prepareRerunState(rerun_fr
, state_global
, constructVsites
, vsite
, top
.idef
, ir
->delta_t
, *fr
, graph
);
516 isLastStep
= isLastStep
|| stopHandler
->stoppingAfterCurrentStep(bNS
);
518 if (DOMAINDECOMP(cr
))
520 /* Repartition the domain decomposition */
521 const bool bMasterState
= true;
522 dd_partition_system(fplog
, mdlog
, step
, cr
,
523 bMasterState
, nstglobalcomm
,
524 state_global
, *top_global
, ir
,
525 state
, &f
, mdAtoms
, &top
, fr
,
528 mdrunOptions
.verbose
);
529 shouldCheckNumberOfBondedInteractions
= true;
534 print_ebin_header(fplog
, step
, t
); /* can we improve the information printed here? */
537 if (ir
->efep
!= efepNO
)
539 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
542 force_flags
= (GMX_FORCE_STATECHANGED
|
543 GMX_FORCE_DYNAMICBOX
|
544 GMX_FORCE_ALLFORCES
|
545 (GMX_GPU
? GMX_FORCE_VIRIAL
: 0) | // TODO: Get rid of this once #2649 is solved
547 (doFreeEnergyPerturbation
? GMX_FORCE_DHDL
: 0));
551 /* Now is the time to relax the shells */
552 relax_shell_flexcon(fplog
, cr
, ms
, mdrunOptions
.verbose
,
553 enforcedRotation
, step
,
554 ir
, bNS
, force_flags
, &top
,
556 state
, f
.arrayRefWithPadding(), force_vir
, mdatoms
,
558 shellfc
, fr
, ppForceWorkload
, t
, mu_tot
,
560 ddBalanceRegionHandler
);
564 /* The coordinates (x) are shifted (to get whole molecules)
566 * This is parallellized as well, and does communication too.
567 * Check comments in sim_util.c
570 gmx_edsam
*ed
= nullptr;
571 do_force(fplog
, cr
, ms
, ir
, awh
, enforcedRotation
,
572 step
, nrnb
, wcycle
, &top
,
573 state
->box
, state
->x
.arrayRefWithPadding(), &state
->hist
,
574 f
.arrayRefWithPadding(), force_vir
, mdatoms
, enerd
, fcd
,
575 state
->lambda
, graph
,
576 fr
, ppForceWorkload
, vsite
, mu_tot
, t
, ed
,
577 GMX_FORCE_NS
| force_flags
,
578 ddBalanceRegionHandler
);
581 /* Now we have the energies and forces corresponding to the
582 * coordinates at time t.
585 const bool isCheckpointingStep
= false;
586 const bool doRerun
= true;
587 const bool bSumEkinhOld
= false;
588 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
589 ir
, state
, state_global
, observablesHistory
,
591 outf
, energyOutput
, ekind
, f
,
592 isCheckpointingStep
, doRerun
, isLastStep
,
593 mdrunOptions
.writeConfout
,
597 stopHandler
->setSignal();
601 /* Need to unshift here */
602 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
605 if (vsite
!= nullptr)
607 wallcycle_start(wcycle
, ewcVSITECONSTR
);
608 if (graph
!= nullptr)
610 shift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
612 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, as_rvec_array(state
->v
.data()),
613 top
.idef
.iparams
, top
.idef
.il
,
614 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
616 if (graph
!= nullptr)
618 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
620 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
624 const bool doInterSimSignal
= false;
625 const bool doIntraSimSignal
= true;
626 bool bSumEkinhOld
= false;
627 t_vcm
*vcm
= nullptr;
628 SimulationSignaller
signaller(&signals
, cr
, ms
, doInterSimSignal
, doIntraSimSignal
);
630 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
631 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
634 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
635 CGLO_GSTAT
| CGLO_ENERGY
636 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
638 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
639 top_global
, &top
, state
,
640 &shouldCheckNumberOfBondedInteractions
);
643 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
644 the virial that should probably be addressed eventually. state->veta has better properies,
645 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
646 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
648 if (ir
->efep
!= efepNO
)
650 /* Sum up the foreign energy and dhdl terms for md and sd.
651 Currently done every step so that dhdl is correct in the .edr */
652 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
658 const bool bCalcEnerStep
= true;
659 energyOutput
.addDataAtEnergyStep(doFreeEnergyPerturbation
, bCalcEnerStep
,
660 t
, mdatoms
->tmass
, enerd
, state
,
661 ir
->fepvals
, ir
->expandedvals
, state
->box
,
662 shake_vir
, force_vir
, total_vir
, pres
,
663 ekind
, mu_tot
, constr
);
665 const bool do_ene
= true;
666 const bool do_log
= true;
668 const bool do_dr
= ir
->nstdisreout
!= 0;
669 const bool do_or
= ir
->nstorireout
!= 0;
671 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
,
672 do_log
? fplog
: nullptr,
674 eprNORMAL
, fcd
, groups
, &(ir
->opts
), awh
);
676 if (do_per_step(step
, ir
->nstlog
))
678 if (fflush(fplog
) != 0)
680 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
685 /* Print the remaining wall clock time for the run */
686 if (isMasterSimMasterRank(ms
, cr
) &&
687 (mdrunOptions
.verbose
|| gmx_got_usr_signal()))
691 fprintf(stderr
, "\n");
693 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
696 /* Ion/water position swapping.
697 * Not done in last step since trajectory writing happens before this call
698 * in the MD loop and exchanges would be lost anyway. */
699 if ((ir
->eSwapCoords
!= eswapNO
) && (step
> 0) && !isLastStep
&&
700 do_per_step(step
, ir
->swap
->nstswap
))
702 const bool doRerun
= true;
703 do_swapcoords(cr
, step
, t
, ir
, wcycle
,
706 MASTER(cr
) && mdrunOptions
.verbose
,
712 /* read next frame from input trajectory */
713 isLastStep
= !read_next_frame(oenv
, status
, &rerun_fr
);
718 rerun_parallel_comm(cr
, &rerun_fr
, &isLastStep
);
721 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
722 if (DOMAINDECOMP(cr
) && wcycle
)
724 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
729 /* increase the MD step number */
734 /* End of main MD loop */
736 /* Closing TNG files can include compressing data. Therefore it is good to do that
737 * before stopping the time measurements. */
738 mdoutf_tng_close(outf
);
740 /* Stop measuring walltime */
741 walltime_accounting_end_time(walltime_accounting
);
748 if (!thisRankHasDuty(cr
, DUTY_PME
))
750 /* Tell the PME only node to finish */
751 gmx_pme_send_finish(cr
);
756 done_shellfc(fplog
, shellfc
, step_rel
);
758 // Clean up swapcoords
759 if (ir
->eSwapCoords
!= eswapNO
)
761 finish_swapcoords(ir
->swap
);
764 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);