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/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
,
159 const t_forcerec
&forceRec
,
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
);
169 GMX_ASSERT(vsite
, "Need valid vsite for constructing vsites");
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
);
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;
204 tensor force_vir
, shake_vir
, total_vir
, pres
;
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_shellfc_t
*shellfc
;
219 /* Domain decomposition could incorrectly miss a bonded
220 interaction, but checking for that requires a global
221 communication stage, which does not otherwise happen in DD
222 code. So we do that alongside the first global energy reduction
223 after a new DD is made. These variables handle whether the
224 check happens, and the result it returns. */
225 bool shouldCheckNumberOfBondedInteractions
= false;
226 int totalNumberOfBondedInteractions
= -1;
228 SimulationSignals signals
;
229 // Most global communnication stages don't propagate mdrun
230 // signals, and will use this object to achieve that.
231 SimulationSignaller
nullSignaller(nullptr, nullptr, nullptr, false, false);
233 GMX_LOG(mdlog
.info
).asParagraph().
234 appendText("Note that it is planned that the command gmx mdrun -rerun will "
235 "be available in a different form in a future version of GROMACS, "
236 "e.g. gmx rerun -f.");
240 gmx_fatal(FARGS
, "Expanded ensemble not supported by rerun.");
244 gmx_fatal(FARGS
, "Simulated tempering not supported by rerun.");
248 gmx_fatal(FARGS
, "AWH not supported by rerun.");
250 if (replExParams
.exchangeInterval
> 0)
252 gmx_fatal(FARGS
, "Replica exchange not supported by rerun.");
254 if (opt2bSet("-ei", nfile
, fnm
) || observablesHistory
->edsamHistory
!= nullptr)
256 gmx_fatal(FARGS
, "Essential dynamics not supported by rerun.");
260 gmx_fatal(FARGS
, "Interactive MD not supported by rerun.");
264 gmx_fatal(FARGS
, "Multiple simulations not supported by rerun.");
266 if (std::any_of(ir
->opts
.annealing
, ir
->opts
.annealing
+ ir
->opts
.ngtc
,
267 [](int i
){return i
!= eannNO
; }))
269 gmx_fatal(FARGS
, "Simulated annealing not supported by rerun.");
272 /* Rerun can't work if an output file name is the same as the input file name.
273 * If this is the case, the user will get an error telling them what the issue is.
275 if (strcmp(opt2fn("-rerun", nfile
, fnm
), opt2fn("-o", nfile
, fnm
)) == 0 ||
276 strcmp(opt2fn("-rerun", nfile
, fnm
), opt2fn("-x", nfile
, fnm
)) == 0)
278 gmx_fatal(FARGS
, "When using mdrun -rerun, the name of the input trajectory file "
279 "%s cannot be identical to the name of an output file (whether "
280 "given explicitly with -o or -x, or by default)",
281 opt2fn("-rerun", nfile
, fnm
));
284 /* Settings for rerun */
286 ir
->nstcalcenergy
= 1;
287 int nstglobalcomm
= 1;
288 const bool bNS
= true;
290 ir
->nstxout_compressed
= 0;
291 groups
= &top_global
->groups
;
292 if (ir
->eI
== eiMimic
)
294 top_global
->intermolecularExclusionGroup
= genQmmmIndices(*top_global
);
298 init_rerun(fplog
, cr
, outputProvider
, ir
, oenv
, mdrunOptions
,
299 state_global
, lam0
, nrnb
, top_global
,
300 nfile
, fnm
, &outf
, &mdebin
, wcycle
);
302 /* Energy terms and groups */
304 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
307 /* Kinetic energy data */
308 std::unique_ptr
<gmx_ekindata_t
> eKinData
= compat::make_unique
<gmx_ekindata_t
>();
309 gmx_ekindata_t
*ekind
= eKinData
.get();
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
))
326 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
331 // Local state only becomes valid now.
332 std::unique_ptr
<t_state
> stateInstance
;
335 if (DOMAINDECOMP(cr
))
337 dd_init_local_top(*top_global
, &top
);
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
,
348 nrnb
, nullptr, FALSE
);
349 shouldCheckNumberOfBondedInteractions
= true;
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 mdAlgorithmsSetupAtomData(cr
, ir
, *top_global
, &top
, fr
,
362 &graph
, mdAtoms
, constr
, vsite
, shellfc
);
365 auto mdatoms
= mdAtoms
->mdatoms();
367 // NOTE: The global state is no longer used at this point.
368 // But state_global is still used as temporary storage space for writing
369 // the global state to file and potentially for replica exchange.
370 // (Global topology should persist.)
372 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
374 if (ir
->efep
!= efepNO
&& ir
->fepvals
->nstdhdl
!= 0)
376 doFreeEnergyPerturbation
= true;
380 int cglo_flags
= (CGLO_INITIALIZATION
| CGLO_GSTAT
|
381 (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
382 bool bSumEkinhOld
= false;
383 t_vcm
*vcm
= nullptr;
384 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
385 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
386 constr
, &nullSignaller
, state
->box
,
387 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags
);
389 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
390 top_global
, &top
, state
,
391 &shouldCheckNumberOfBondedInteractions
);
395 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
396 " input trajectory '%s'\n\n",
397 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
398 if (mdrunOptions
.verbose
)
400 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
401 "run input file,\nwhich may not correspond to the time "
402 "needed to process input trajectory.\n\n");
404 fprintf(fplog
, "\n");
407 walltime_accounting_start_time(walltime_accounting
);
408 wallcycle_start(wcycle
, ewcRUN
);
409 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
411 /***********************************************************
415 ************************************************************/
419 GMX_LOG(mdlog
.info
).asParagraph().
420 appendText("Simulations has constraints. Rerun does not recalculate constraints.");
426 isLastStep
= !read_first_frame(oenv
, &status
,
427 opt2fn("-rerun", nfile
, fnm
),
428 &rerun_fr
, TRX_NEED_X
);
429 if (rerun_fr
.natoms
!= top_global
->natoms
)
432 "Number of atoms in trajectory (%d) does not match the "
433 "run input file (%d)\n",
434 rerun_fr
.natoms
, top_global
->natoms
);
437 if (ir
->ePBC
!= epbcNONE
)
441 gmx_fatal(FARGS
, "Rerun trajectory frame step %" PRId64
" time %f "
442 "does not contain a box, while pbc is used",
443 rerun_fr
.step
, rerun_fr
.time
);
445 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < gmx::square(fr
->rlist
))
447 gmx_fatal(FARGS
, "Rerun trajectory frame step %" PRId64
" time %f "
448 "has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
453 GMX_LOG(mdlog
.info
).asParagraph().
454 appendText("Rerun does not report kinetic energy, total energy, temperature, virial and pressure.");
458 rerun_parallel_comm(cr
, &rerun_fr
, &isLastStep
);
461 if (ir
->ePBC
!= epbcNONE
)
463 /* Set the shift vectors.
464 * Necessary here when have a static box different from the tpr box.
466 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
469 auto stopHandler
= stopHandlerBuilder
->getStopHandlerMD(
470 compat::not_null
<SimulationSignal
*>(&signals
[eglsSTOPCOND
]), false,
471 MASTER(cr
), ir
->nstlist
, mdrunOptions
.reproducible
, nstglobalcomm
,
472 mdrunOptions
.maximumHoursToRun
, ir
->nstlist
== 0, fplog
, step
, bNS
, walltime_accounting
);
474 // we don't do counter resetting in rerun - finish will always be valid
475 walltime_accounting_set_valid_finish(walltime_accounting
);
477 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion
= (DOMAINDECOMP(cr
) ? DdOpenBalanceRegionBeforeForceComputation::yes
: DdOpenBalanceRegionBeforeForceComputation::no
);
478 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion
= (DOMAINDECOMP(cr
) ? DdCloseBalanceRegionAfterForceComputation::yes
: DdCloseBalanceRegionAfterForceComputation::no
);
480 step
= ir
->init_step
;
483 /* and stop now if we should */
484 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
487 wallcycle_start(wcycle
, ewcSTEP
);
491 step
= rerun_fr
.step
;
492 step_rel
= step
- ir
->init_step
;
503 if (ir
->efep
!= efepNO
&& MASTER(cr
))
505 setCurrentLambdasRerun(step
, ir
->fepvals
, &rerun_fr
, lam0
, state_global
);
510 const bool constructVsites
= ((vsite
!= nullptr) && mdrunOptions
.rerunConstructVsites
);
511 if (constructVsites
&& DOMAINDECOMP(cr
))
513 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
514 "use a single rank");
516 prepareRerunState(rerun_fr
, state_global
, constructVsites
, vsite
, top
.idef
, ir
->delta_t
, *fr
, graph
);
519 isLastStep
= isLastStep
|| stopHandler
->stoppingAfterCurrentStep(bNS
);
521 if (DOMAINDECOMP(cr
))
523 /* Repartition the domain decomposition */
524 const bool bMasterState
= true;
525 dd_partition_system(fplog
, mdlog
, step
, cr
,
526 bMasterState
, nstglobalcomm
,
527 state_global
, *top_global
, ir
,
528 state
, &f
, mdAtoms
, &top
, fr
,
531 mdrunOptions
.verbose
);
532 shouldCheckNumberOfBondedInteractions
= true;
537 print_ebin_header(fplog
, step
, t
); /* can we improve the information printed here? */
540 if (ir
->efep
!= efepNO
)
542 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
545 force_flags
= (GMX_FORCE_STATECHANGED
|
546 GMX_FORCE_DYNAMICBOX
|
547 GMX_FORCE_ALLFORCES
|
548 (GMX_GPU
? GMX_FORCE_VIRIAL
: 0) | // TODO: Get rid of this once #2649 is solved
550 (doFreeEnergyPerturbation
? GMX_FORCE_DHDL
: 0));
554 /* Now is the time to relax the shells */
555 relax_shell_flexcon(fplog
, cr
, ms
, mdrunOptions
.verbose
,
556 enforcedRotation
, step
,
557 ir
, bNS
, force_flags
, &top
,
559 state
, f
.arrayRefWithPadding(), force_vir
, mdatoms
,
560 nrnb
, wcycle
, graph
, groups
,
561 shellfc
, fr
, ppForceWorkload
, t
, mu_tot
,
563 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
567 /* The coordinates (x) are shifted (to get whole molecules)
569 * This is parallellized as well, and does communication too.
570 * Check comments in sim_util.c
573 gmx_edsam
*ed
= nullptr;
574 do_force(fplog
, cr
, ms
, ir
, awh
, enforcedRotation
,
575 step
, nrnb
, wcycle
, &top
, groups
,
576 state
->box
, state
->x
.arrayRefWithPadding(), &state
->hist
,
577 f
.arrayRefWithPadding(), force_vir
, mdatoms
, enerd
, fcd
,
578 state
->lambda
, graph
,
579 fr
, ppForceWorkload
, vsite
, mu_tot
, t
, ed
,
580 GMX_FORCE_NS
| force_flags
,
581 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
584 /* Now we have the energies and forces corresponding to the
585 * coordinates at time t.
588 const bool isCheckpointingStep
= false;
589 const bool doRerun
= true;
590 const bool bSumEkinhOld
= false;
591 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
592 ir
, state
, state_global
, observablesHistory
,
594 outf
, mdebin
, ekind
, f
,
595 isCheckpointingStep
, doRerun
, isLastStep
,
596 mdrunOptions
.writeConfout
,
600 stopHandler
->setSignal();
604 /* Need to unshift here */
605 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
608 if (vsite
!= nullptr)
610 wallcycle_start(wcycle
, ewcVSITECONSTR
);
611 if (graph
!= nullptr)
613 shift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
615 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, as_rvec_array(state
->v
.data()),
616 top
.idef
.iparams
, top
.idef
.il
,
617 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
619 if (graph
!= nullptr)
621 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
623 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
627 const bool doInterSimSignal
= false;
628 const bool doIntraSimSignal
= true;
629 bool bSumEkinhOld
= false;
630 t_vcm
*vcm
= nullptr;
631 SimulationSignaller
signaller(&signals
, cr
, ms
, doInterSimSignal
, doIntraSimSignal
);
633 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
634 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
637 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
638 CGLO_GSTAT
| CGLO_ENERGY
639 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
641 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
642 top_global
, &top
, state
,
643 &shouldCheckNumberOfBondedInteractions
);
646 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
647 the virial that should probably be addressed eventually. state->veta has better properies,
648 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
649 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
651 if (ir
->efep
!= efepNO
)
653 /* Sum up the foreign energy and dhdl terms for md and sd.
654 Currently done every step so that dhdl is correct in the .edr */
655 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
661 const bool bCalcEnerStep
= true;
662 upd_mdebin(mdebin
, doFreeEnergyPerturbation
, bCalcEnerStep
,
663 t
, mdatoms
->tmass
, enerd
, state
,
664 ir
->fepvals
, ir
->expandedvals
, rerun_fr
.box
,
665 shake_vir
, force_vir
, total_vir
, pres
,
666 ekind
, mu_tot
, constr
);
668 const bool do_ene
= true;
669 const bool do_log
= true;
671 const bool do_dr
= ir
->nstdisreout
!= 0;
672 const bool do_or
= ir
->nstorireout
!= 0;
674 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
, do_log
? fplog
: nullptr,
676 eprNORMAL
, mdebin
, fcd
, groups
, &(ir
->opts
), awh
);
678 if (do_per_step(step
, ir
->nstlog
))
680 if (fflush(fplog
) != 0)
682 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
687 /* Print the remaining wall clock time for the run */
688 if (isMasterSimMasterRank(ms
, cr
) &&
689 (mdrunOptions
.verbose
|| gmx_got_usr_signal()))
693 fprintf(stderr
, "\n");
695 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
698 /* Ion/water position swapping.
699 * Not done in last step since trajectory writing happens before this call
700 * in the MD loop and exchanges would be lost anyway. */
701 if ((ir
->eSwapCoords
!= eswapNO
) && (step
> 0) && !isLastStep
&&
702 do_per_step(step
, ir
->swap
->nstswap
))
704 const bool doRerun
= true;
705 do_swapcoords(cr
, step
, t
, ir
, wcycle
,
708 MASTER(cr
) && mdrunOptions
.verbose
,
714 /* read next frame from input trajectory */
715 isLastStep
= !read_next_frame(oenv
, status
, &rerun_fr
);
720 rerun_parallel_comm(cr
, &rerun_fr
, &isLastStep
);
723 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
724 if (DOMAINDECOMP(cr
) && wcycle
)
726 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
731 /* increase the MD step number */
736 /* End of main MD loop */
738 /* Closing TNG files can include compressing data. Therefore it is good to do that
739 * before stopping the time measurements. */
740 mdoutf_tng_close(outf
);
742 /* Stop measuring walltime */
743 walltime_accounting_end_time(walltime_accounting
);
750 if (!thisRankHasDuty(cr
, DUTY_PME
))
752 /* Tell the PME only node to finish */
753 gmx_pme_send_finish(cr
);
759 done_shellfc(fplog
, shellfc
, step_rel
);
761 // Clean up swapcoords
762 if (ir
->eSwapCoords
!= eswapNO
)
764 finish_swapcoords(ir
->swap
);
767 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);
769 destroy_enerdata(enerd
);