Moved essential dynamics initialization call into do_md()
[gromacs.git] / src / programs / mdrun / md.cpp
blobb73092eaf81aab9a123441b7f2e4b046205af1c7
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "md.h"
41 #include "config.h"
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
47 #include <algorithm>
48 #include <memory>
50 #include "thread_mpi/threads.h"
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_network.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/essentialdynamics/edsam.h"
57 #include "gromacs/ewald/pme.h"
58 #include "gromacs/ewald/pme-load-balancing.h"
59 #include "gromacs/fileio/trxio.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/gmxlib/nrnb.h"
62 #include "gromacs/gpu_utils/gpu_utils.h"
63 #include "gromacs/imd/imd.h"
64 #include "gromacs/listed-forces/manage-threading.h"
65 #include "gromacs/math/functions.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/math/vectypes.h"
69 #include "gromacs/mdlib/compute_io.h"
70 #include "gromacs/mdlib/constr.h"
71 #include "gromacs/mdlib/ebin.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/force_flags.h"
74 #include "gromacs/mdlib/forcerec.h"
75 #include "gromacs/mdlib/md_support.h"
76 #include "gromacs/mdlib/mdatoms.h"
77 #include "gromacs/mdlib/mdebin.h"
78 #include "gromacs/mdlib/mdoutf.h"
79 #include "gromacs/mdlib/mdrun.h"
80 #include "gromacs/mdlib/mdsetup.h"
81 #include "gromacs/mdlib/nb_verlet.h"
82 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
83 #include "gromacs/mdlib/ns.h"
84 #include "gromacs/mdlib/shellfc.h"
85 #include "gromacs/mdlib/sighandler.h"
86 #include "gromacs/mdlib/sim_util.h"
87 #include "gromacs/mdlib/simulationsignal.h"
88 #include "gromacs/mdlib/tgroup.h"
89 #include "gromacs/mdlib/trajectory_writing.h"
90 #include "gromacs/mdlib/update.h"
91 #include "gromacs/mdlib/vcm.h"
92 #include "gromacs/mdlib/vsite.h"
93 #include "gromacs/mdtypes/commrec.h"
94 #include "gromacs/mdtypes/df_history.h"
95 #include "gromacs/mdtypes/energyhistory.h"
96 #include "gromacs/mdtypes/fcdata.h"
97 #include "gromacs/mdtypes/forcerec.h"
98 #include "gromacs/mdtypes/group.h"
99 #include "gromacs/mdtypes/inputrec.h"
100 #include "gromacs/mdtypes/interaction_const.h"
101 #include "gromacs/mdtypes/md_enums.h"
102 #include "gromacs/mdtypes/mdatom.h"
103 #include "gromacs/mdtypes/observableshistory.h"
104 #include "gromacs/mdtypes/state.h"
105 #include "gromacs/pbcutil/mshift.h"
106 #include "gromacs/pbcutil/pbc.h"
107 #include "gromacs/pulling/pull.h"
108 #include "gromacs/swap/swapcoords.h"
109 #include "gromacs/timing/wallcycle.h"
110 #include "gromacs/timing/walltime_accounting.h"
111 #include "gromacs/topology/atoms.h"
112 #include "gromacs/topology/idef.h"
113 #include "gromacs/topology/mtop_util.h"
114 #include "gromacs/topology/topology.h"
115 #include "gromacs/trajectory/trajectoryframe.h"
116 #include "gromacs/utility/basedefinitions.h"
117 #include "gromacs/utility/cstringutil.h"
118 #include "gromacs/utility/fatalerror.h"
119 #include "gromacs/utility/logger.h"
120 #include "gromacs/utility/real.h"
121 #include "gromacs/utility/smalloc.h"
123 #include "deform.h"
124 #include "membed.h"
125 #include "repl_ex.h"
127 #ifdef GMX_FAHCORE
128 #include "corewrap.h"
129 #endif
131 using gmx::SimulationSignaller;
133 /*! \brief Check whether bonded interactions are missing, if appropriate
135 * \param[in] fplog Log file pointer
136 * \param[in] cr Communication object
137 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
138 * \param[in] top_global Global topology for the error message
139 * \param[in] top_local Local topology for the error message
140 * \param[in] state Global state for the error message
141 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
143 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
144 * is always set to false after exit.
146 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
147 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
148 bool *shouldCheckNumberOfBondedInteractions)
150 if (*shouldCheckNumberOfBondedInteractions)
152 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
154 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
156 *shouldCheckNumberOfBondedInteractions = false;
160 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
161 gmx_int64_t step,
162 gmx_int64_t *step_rel, t_inputrec *ir,
163 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
164 gmx_walltime_accounting_t walltime_accounting,
165 struct nonbonded_verlet_t *nbv)
167 char sbuf[STEPSTRSIZE];
169 /* Reset all the counters related to performance over the run */
170 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
171 "step %s: resetting all time and cycle counters",
172 gmx_step_str(step, sbuf));
174 if (use_GPU(nbv))
176 nbnxn_gpu_reset_timings(nbv);
177 resetGpuProfiler();
180 wallcycle_stop(wcycle, ewcRUN);
181 wallcycle_reset_all(wcycle);
182 if (DOMAINDECOMP(cr))
184 reset_dd_statistics_counters(cr->dd);
186 init_nrnb(nrnb);
187 ir->init_step += *step_rel;
188 ir->nsteps -= *step_rel;
189 *step_rel = 0;
190 wallcycle_start(wcycle, ewcRUN);
191 walltime_accounting_start(walltime_accounting);
192 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
195 /*! \libinternal
196 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
197 int nfile, const t_filenm fnm[],
198 const gmx_output_env_t *oenv, gmx_bool bVerbose,
199 int nstglobalcomm,
200 gmx_vsite_t *vsite, gmx_constr_t constr,
201 int stepout,
202 gmx::IMDOutputProvider *outputProvider,
203 t_inputrec *inputrec,
204 gmx_mtop_t *top_global, t_fcdata *fcd,
205 t_state *state_global,
206 t_mdatoms *mdatoms,
207 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
208 t_forcerec *fr,
209 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
210 real cpt_period, real max_hours,
211 int imdport,
212 unsigned long Flags,
213 gmx_walltime_accounting_t walltime_accounting)
215 double gmx::do_md(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
216 int nfile, const t_filenm fnm[],
217 const gmx_output_env_t *oenv, gmx_bool bVerbose,
218 int nstglobalcomm,
219 gmx_vsite_t *vsite, gmx_constr_t constr,
220 int stepout, gmx::IMDOutputProvider *outputProvider,
221 t_inputrec *ir,
222 gmx_mtop_t *top_global,
223 t_fcdata *fcd,
224 t_state *state_global,
225 ObservablesHistory *observablesHistory,
226 t_mdatoms *mdatoms,
227 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
228 t_forcerec *fr,
229 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
230 gmx_membed_t *membed,
231 real cpt_period, real max_hours,
232 int imdport,
233 unsigned long Flags,
234 gmx_walltime_accounting_t walltime_accounting)
236 gmx_mdoutf_t outf = nullptr;
237 gmx_int64_t step, step_rel;
238 double elapsed_time;
239 double t, t0, lam0[efptNR];
240 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
241 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD,
242 bFirstStep, startingFromCheckpoint, bInitStep, bLastStep = FALSE,
243 bBornRadii, bUsingEnsembleRestraints;
244 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
245 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
246 bForceUpdate = FALSE, bCPT;
247 gmx_bool bMasterState;
248 int force_flags, cglo_flags;
249 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
250 int i, m;
251 t_trxstatus *status;
252 rvec mu_tot;
253 t_vcm *vcm;
254 matrix parrinellorahmanMu, M;
255 t_trxframe rerun_fr;
256 gmx_repl_ex_t repl_ex = nullptr;
257 int nchkpt = 1;
258 gmx_localtop_t *top;
259 t_mdebin *mdebin = nullptr;
260 gmx_enerdata_t *enerd;
261 PaddedRVecVector f {};
262 gmx_global_stat_t gstat;
263 gmx_update_t *upd = nullptr;
264 t_graph *graph = nullptr;
265 gmx_groups_t *groups;
266 gmx_ekindata_t *ekind;
267 gmx_shellfc_t *shellfc;
268 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
269 gmx_bool bResetCountersHalfMaxH = FALSE;
270 gmx_bool bTemp, bPres, bTrotter;
271 real dvdl_constr;
272 rvec *cbuf = nullptr;
273 int cbuf_nalloc = 0;
274 matrix lastbox;
275 int lamnew = 0;
276 /* for FEP */
277 int nstfep = 0;
278 double cycles;
279 real saved_conserved_quantity = 0;
280 real last_ekin = 0;
281 t_extmass MassQ;
282 int **trotter_seq;
283 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
284 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
287 /* PME load balancing data for GPU kernels */
288 pme_load_balancing_t *pme_loadbal = nullptr;
289 gmx_bool bPMETune = FALSE;
290 gmx_bool bPMETunePrinting = FALSE;
292 /* Interactive MD */
293 gmx_bool bIMDstep = FALSE;
295 #ifdef GMX_FAHCORE
296 /* Temporary addition for FAHCORE checkpointing */
297 int chkpt_ret;
298 #endif
299 /* Domain decomposition could incorrectly miss a bonded
300 interaction, but checking for that requires a global
301 communication stage, which does not otherwise happen in DD
302 code. So we do that alongside the first global energy reduction
303 after a new DD is made. These variables handle whether the
304 check happens, and the result it returns. */
305 bool shouldCheckNumberOfBondedInteractions = false;
306 int totalNumberOfBondedInteractions = -1;
308 SimulationSignals signals;
309 // Most global communnication stages don't propagate mdrun
310 // signals, and will use this object to achieve that.
311 SimulationSignaller nullSignaller(nullptr, nullptr, false, false);
313 /* Check for special mdrun options */
314 bRerunMD = (Flags & MD_RERUN);
315 if (Flags & MD_RESETCOUNTERSHALFWAY)
317 if (ir->nsteps > 0)
319 /* Signal to reset the counters half the simulation steps. */
320 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
322 /* Signal to reset the counters halfway the simulation time. */
323 bResetCountersHalfMaxH = (max_hours > 0);
326 /* md-vv uses averaged full step velocities for T-control
327 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
328 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
329 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
331 if (bRerunMD)
333 /* Since we don't know if the frames read are related in any way,
334 * rebuild the neighborlist at every step.
336 ir->nstlist = 1;
337 ir->nstcalcenergy = 1;
338 nstglobalcomm = 1;
341 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
342 bGStatEveryStep = (nstglobalcomm == 1);
344 if (bRerunMD)
346 ir->nstxout_compressed = 0;
348 groups = &top_global->groups;
350 gmx_edsam *ed = nullptr;
351 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
353 /* Initialize essential dynamics sampling */
354 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
355 top_global, ir, cr, constr,
356 as_rvec_array(state_global->x.data()),
357 state_global->box, observablesHistory,
358 oenv, Flags & MD_APPENDFILES);
361 if (ir->eSwapCoords != eswapNO)
363 /* Initialize ion swapping code */
364 init_swapcoords(fplog, bVerbose, ir, opt2fn_master("-swap", nfile, fnm, cr),
365 top_global, as_rvec_array(state_global->x.data()), state_global->box, observablesHistory, cr, oenv, Flags);
368 /* Initial values */
369 init_md(fplog, cr, outputProvider, ir, oenv, &t, &t0, state_global->lambda,
370 &(state_global->fep_state), lam0,
371 nrnb, top_global, &upd,
372 nfile, fnm, &outf, &mdebin,
373 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
375 clear_mat(total_vir);
376 clear_mat(pres);
377 /* Energy terms and groups */
378 snew(enerd, 1);
379 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
380 enerd);
382 /* Kinetic energy data */
383 snew(ekind, 1);
384 init_ekindata(fplog, top_global, &(ir->opts), ekind);
385 /* Copy the cos acceleration to the groups struct */
386 ekind->cosacc.cos_accel = ir->cos_accel;
388 gstat = global_stat_init(ir);
390 /* Check for polarizable models and flexible constraints */
391 shellfc = init_shell_flexcon(fplog,
392 top_global, n_flexible_constraints(constr),
393 ir->nstcalcenergy, DOMAINDECOMP(cr));
395 if (shellfc && ir->nstcalcenergy != 1)
397 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir->nstcalcenergy);
399 if (shellfc && DOMAINDECOMP(cr))
401 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
404 if (inputrecDeform(ir))
406 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
407 set_deform_reference_box(upd,
408 deform_init_init_step_tpx,
409 deform_init_box_tpx);
410 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
414 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
415 if ((io > 2000) && MASTER(cr))
417 fprintf(stderr,
418 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
419 io);
423 std::unique_ptr<t_state> stateInstance;
424 t_state * state;
426 if (DOMAINDECOMP(cr))
428 top = dd_init_local_top(top_global);
430 stateInstance = std::unique_ptr<t_state>(new t_state);
431 state = stateInstance.get();
432 dd_init_local_state(cr->dd, state_global, state);
434 else
436 state_change_natoms(state_global, state_global->natoms);
437 /* We need to allocate one element extra, since we might use
438 * (unaligned) 4-wide SIMD loads to access rvec entries.
440 f.resize(state_global->natoms + 1);
441 /* Copy the pointer to the global state */
442 state = state_global;
444 snew(top, 1);
445 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
446 &graph, mdatoms, vsite, shellfc);
448 update_realloc(upd, state->natoms);
451 /* Set up interactive MD (IMD) */
452 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, as_rvec_array(state_global->x.data()),
453 nfile, fnm, oenv, imdport, Flags);
455 if (DOMAINDECOMP(cr))
457 /* Distribute the charge groups over the nodes from the master node */
458 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
459 state_global, top_global, ir,
460 state, &f, mdatoms, top, fr,
461 vsite, constr,
462 nrnb, nullptr, FALSE);
463 shouldCheckNumberOfBondedInteractions = true;
464 update_realloc(upd, state->natoms);
467 update_mdatoms(mdatoms, state->lambda[efptMASS]);
469 startingFromCheckpoint = Flags & MD_STARTFROMCPT;
471 if (ir->bExpanded)
473 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
476 if (MASTER(cr))
478 if (startingFromCheckpoint)
480 /* Update mdebin with energy history if appending to output files */
481 if (Flags & MD_APPENDFILES)
483 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
485 else if (observablesHistory->energyHistory.get() != nullptr)
487 /* We might have read an energy history from checkpoint.
488 * As we are not appending, we want to restart the statistics.
489 * Free the allocated memory and reset the counts.
491 observablesHistory->energyHistory = {};
494 if (observablesHistory->energyHistory.get() == nullptr)
496 observablesHistory->energyHistory = std::unique_ptr<energyhistory_t>(new energyhistory_t {});
498 /* Set the initial energy history in state by updating once */
499 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
502 /* Initialize constraints */
503 if (constr && !DOMAINDECOMP(cr))
505 set_constraints(constr, top, ir, mdatoms, cr);
508 if (repl_ex_nst > 0 && MASTER(cr))
510 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
511 repl_ex_nst, repl_ex_nex, repl_ex_seed);
514 /* PME tuning is only supported with PME for Coulomb. Is is not supported
515 * with only LJ PME, or for reruns.
517 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
518 !(Flags & MD_REPRODUCIBLE));
519 if (bPMETune)
521 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
522 fr->ic, fr->pmedata, use_GPU(fr->nbv),
523 &bPMETunePrinting);
526 if (!ir->bContinuation && !bRerunMD)
528 if (state->flags & (1 << estV))
530 /* Set the velocities of vsites, shells and frozen atoms to zero */
531 for (i = 0; i < mdatoms->homenr; i++)
533 if (mdatoms->ptype[i] == eptVSite ||
534 mdatoms->ptype[i] == eptShell)
536 clear_rvec(state->v[i]);
538 else if (mdatoms->cFREEZE)
540 for (m = 0; m < DIM; m++)
542 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
544 state->v[i][m] = 0;
551 if (constr)
553 /* Constrain the initial coordinates and velocities */
554 do_constrain_first(fplog, constr, ir, mdatoms, state,
555 cr, nrnb, fr, top);
557 if (vsite)
559 /* Construct the virtual sites for the initial configuration */
560 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
561 top->idef.iparams, top->idef.il,
562 fr->ePBC, fr->bMolPBC, cr, state->box);
566 if (ir->efep != efepNO)
568 /* Set free energy calculation frequency as the greatest common
569 * denominator of nstdhdl and repl_ex_nst. */
570 nstfep = ir->fepvals->nstdhdl;
571 if (ir->bExpanded)
573 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
575 if (repl_ex_nst > 0)
577 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
581 /* Be REALLY careful about what flags you set here. You CANNOT assume
582 * this is the first step, since we might be restarting from a checkpoint,
583 * and in that case we should not do any modifications to the state.
585 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
587 if (Flags & MD_READ_EKIN)
589 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
592 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
593 | (bStopCM ? CGLO_STOPCM : 0)
594 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
595 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
596 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
598 bSumEkinhOld = FALSE;
599 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
600 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
601 constr, &nullSignaller, state->box,
602 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags
603 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
604 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
605 top_global, top, state,
606 &shouldCheckNumberOfBondedInteractions);
607 if (ir->eI == eiVVAK)
609 /* a second call to get the half step temperature initialized as well */
610 /* we do the same call as above, but turn the pressure off -- internally to
611 compute_globals, this is recognized as a velocity verlet half-step
612 kinetic energy calculation. This minimized excess variables, but
613 perhaps loses some logic?*/
615 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
616 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
617 constr, &nullSignaller, state->box,
618 nullptr, &bSumEkinhOld,
619 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
622 /* Calculate the initial half step temperature, and save the ekinh_old */
623 if (!(Flags & MD_STARTFROMCPT))
625 for (i = 0; (i < ir->opts.ngtc); i++)
627 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
630 if (ir->eI != eiVV)
632 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
633 and there is no previous step */
636 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
637 temperature control */
638 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
640 if (MASTER(cr))
642 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
644 fprintf(fplog,
645 "RMS relative constraint deviation after constraining: %.2e\n",
646 constr_rmsd(constr));
648 if (EI_STATE_VELOCITY(ir->eI))
650 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
652 if (bRerunMD)
654 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
655 " input trajectory '%s'\n\n",
656 *(top_global->name), opt2fn("-rerun", nfile, fnm));
657 if (bVerbose)
659 fprintf(stderr, "Calculated time to finish depends on nsteps from "
660 "run input file,\nwhich may not correspond to the time "
661 "needed to process input trajectory.\n\n");
664 else
666 char tbuf[20];
667 fprintf(stderr, "starting mdrun '%s'\n",
668 *(top_global->name));
669 if (ir->nsteps >= 0)
671 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
673 else
675 sprintf(tbuf, "%s", "infinite");
677 if (ir->init_step > 0)
679 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
680 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
681 gmx_step_str(ir->init_step, sbuf2),
682 ir->init_step*ir->delta_t);
684 else
686 fprintf(stderr, "%s steps, %s ps.\n",
687 gmx_step_str(ir->nsteps, sbuf), tbuf);
690 fprintf(fplog, "\n");
693 walltime_accounting_start(walltime_accounting);
694 wallcycle_start(wcycle, ewcRUN);
695 print_start(fplog, cr, walltime_accounting, "mdrun");
697 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
698 #ifdef GMX_FAHCORE
699 chkpt_ret = fcCheckPointParallel( cr->nodeid,
700 NULL, 0);
701 if (chkpt_ret == 0)
703 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
705 #endif
707 /***********************************************************
709 * Loop over MD steps
711 ************************************************************/
713 /* if rerunMD then read coordinates and velocities from input trajectory */
714 if (bRerunMD)
716 if (getenv("GMX_FORCE_UPDATE"))
718 bForceUpdate = TRUE;
721 rerun_fr.natoms = 0;
722 if (MASTER(cr))
724 bLastStep = !read_first_frame(oenv, &status,
725 opt2fn("-rerun", nfile, fnm),
726 &rerun_fr, TRX_NEED_X | TRX_READ_V);
727 if (rerun_fr.natoms != top_global->natoms)
729 gmx_fatal(FARGS,
730 "Number of atoms in trajectory (%d) does not match the "
731 "run input file (%d)\n",
732 rerun_fr.natoms, top_global->natoms);
734 if (ir->ePBC != epbcNONE)
736 if (!rerun_fr.bBox)
738 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time);
740 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
742 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
747 if (PAR(cr))
749 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
752 if (ir->ePBC != epbcNONE)
754 /* Set the shift vectors.
755 * Necessary here when have a static box different from the tpr box.
757 calc_shifts(rerun_fr.box, fr->shift_vec);
761 /* loop over MD steps or if rerunMD to end of input trajectory */
762 bFirstStep = TRUE;
763 /* Skip the first Nose-Hoover integration when we get the state from tpx */
764 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
765 bSumEkinhOld = FALSE;
766 bExchanged = FALSE;
767 bNeedRepartition = FALSE;
768 // TODO This implementation of ensemble orientation restraints is nasty because
769 // a user can't just do multi-sim with single-sim orientation restraints.
770 bUsingEnsembleRestraints = (fcd->disres.nsystems > 1) || (cr->ms && fcd->orires.nr);
773 // Replica exchange and ensemble restraints need all
774 // simulations to remain synchronized, so they need
775 // checkpoints and stop conditions to act on the same step, so
776 // the propagation of such signals must take place between
777 // simulations, not just within simulations.
778 bool checkpointIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
779 bool stopConditionIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
780 bool resetCountersIsLocal = true;
781 signals[eglsCHKPT] = SimulationSignal(checkpointIsLocal);
782 signals[eglsSTOPCOND] = SimulationSignal(stopConditionIsLocal);
783 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
786 step = ir->init_step;
787 step_rel = 0;
789 // TODO extract this to new multi-simulation module
790 if (MASTER(cr) && MULTISIM(cr) && (repl_ex_nst <= 0 ))
792 if (!multisim_int_all_are_equal(cr->ms, ir->nsteps))
794 GMX_LOG(mdlog.warning).appendText(
795 "Note: The number of steps is not consistent across multi simulations,\n"
796 "but we are proceeding anyway!");
798 if (!multisim_int_all_are_equal(cr->ms, ir->init_step))
800 GMX_LOG(mdlog.warning).appendText(
801 "Note: The initial step is not consistent across multi simulations,\n"
802 "but we are proceeding anyway!");
806 /* and stop now if we should */
807 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
808 while (!bLastStep)
811 /* Determine if this is a neighbor search step */
812 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
814 if (bPMETune && bNStList)
816 /* PME grid + cut-off optimization with GPUs or PME nodes */
817 pme_loadbal_do(pme_loadbal, cr,
818 (bVerbose && MASTER(cr)) ? stderr : nullptr,
819 fplog, mdlog,
820 ir, fr, state,
821 wcycle,
822 step, step_rel,
823 &bPMETunePrinting);
826 wallcycle_start(wcycle, ewcSTEP);
828 if (bRerunMD)
830 if (rerun_fr.bStep)
832 step = rerun_fr.step;
833 step_rel = step - ir->init_step;
835 if (rerun_fr.bTime)
837 t = rerun_fr.time;
839 else
841 t = step;
844 else
846 bLastStep = (step_rel == ir->nsteps);
847 t = t0 + step*ir->delta_t;
850 // TODO Refactor this, so that nstfep does not need a default value of zero
851 if (ir->efep != efepNO || ir->bSimTemp)
853 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
854 requiring different logic. */
856 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
857 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
858 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
859 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
860 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
863 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
864 do_per_step(step, repl_ex_nst));
866 if (bSimAnn)
868 update_annealing_target_temp(ir, t, upd);
871 if (bRerunMD)
873 if (!DOMAINDECOMP(cr) || MASTER(cr))
875 for (i = 0; i < state_global->natoms; i++)
877 copy_rvec(rerun_fr.x[i], state_global->x[i]);
879 if (rerun_fr.bV)
881 for (i = 0; i < state_global->natoms; i++)
883 copy_rvec(rerun_fr.v[i], state_global->v[i]);
886 else
888 for (i = 0; i < state_global->natoms; i++)
890 clear_rvec(state_global->v[i]);
892 if (bRerunWarnNoV)
894 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
895 " Ekin, temperature and pressure are incorrect,\n"
896 " the virial will be incorrect when constraints are present.\n"
897 "\n");
898 bRerunWarnNoV = FALSE;
902 copy_mat(rerun_fr.box, state_global->box);
903 copy_mat(state_global->box, state->box);
905 if (vsite && (Flags & MD_RERUN_VSITE))
907 if (DOMAINDECOMP(cr))
909 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
911 if (graph)
913 /* Following is necessary because the graph may get out of sync
914 * with the coordinates if we only have every N'th coordinate set
916 mk_mshift(fplog, graph, fr->ePBC, state->box, as_rvec_array(state->x.data()));
917 shift_self(graph, state->box, as_rvec_array(state->x.data()));
919 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
920 top->idef.iparams, top->idef.il,
921 fr->ePBC, fr->bMolPBC, cr, state->box);
922 if (graph)
924 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
929 /* Stop Center of Mass motion */
930 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
932 if (bRerunMD)
934 /* for rerun MD always do Neighbour Searching */
935 bNS = (bFirstStep || ir->nstlist != 0);
937 else
939 /* Determine whether or not to do Neighbour Searching */
940 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
943 /* < 0 means stop at next step, > 0 means stop at next NS step */
944 if ( (signals[eglsSTOPCOND].set < 0) ||
945 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
947 bLastStep = TRUE;
950 /* Determine whether or not to update the Born radii if doing GB */
951 bBornRadii = bFirstStep;
952 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
954 bBornRadii = TRUE;
957 /* do_log triggers energy and virial calculation. Because this leads
958 * to different code paths, forces can be different. Thus for exact
959 * continuation we should avoid extra log output.
960 * Note that the || bLastStep can result in non-exact continuation
961 * beyond the last step. But we don't consider that to be an issue.
963 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
964 do_verbose = bVerbose &&
965 (step % stepout == 0 || bFirstStep || bLastStep || bRerunMD);
967 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
969 if (bRerunMD)
971 bMasterState = TRUE;
973 else
975 bMasterState = FALSE;
976 /* Correct the new box if it is too skewed */
977 if (inputrecDynamicBox(ir))
979 if (correct_box(fplog, step, state->box, graph))
981 bMasterState = TRUE;
984 if (DOMAINDECOMP(cr) && bMasterState)
986 dd_collect_state(cr->dd, state, state_global);
990 if (DOMAINDECOMP(cr))
992 /* Repartition the domain decomposition */
993 dd_partition_system(fplog, step, cr,
994 bMasterState, nstglobalcomm,
995 state_global, top_global, ir,
996 state, &f, mdatoms, top, fr,
997 vsite, constr,
998 nrnb, wcycle,
999 do_verbose && !bPMETunePrinting);
1000 shouldCheckNumberOfBondedInteractions = true;
1001 update_realloc(upd, state->natoms);
1005 if (MASTER(cr) && do_log)
1007 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1010 if (ir->efep != efepNO)
1012 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1015 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1018 /* We need the kinetic energy at minus the half step for determining
1019 * the full step kinetic energy and possibly for T-coupling.*/
1020 /* This may not be quite working correctly yet . . . . */
1021 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1022 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1023 constr, &nullSignaller, state->box,
1024 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1025 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1026 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1027 top_global, top, state,
1028 &shouldCheckNumberOfBondedInteractions);
1030 clear_mat(force_vir);
1032 /* We write a checkpoint at this MD step when:
1033 * either at an NS step when we signalled through gs,
1034 * or at the last step (but not when we do not want confout),
1035 * but never at the first step or with rerun.
1037 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1038 (bLastStep && (Flags & MD_CONFOUT))) &&
1039 step > ir->init_step && !bRerunMD);
1040 if (bCPT)
1042 signals[eglsCHKPT].set = 0;
1045 /* Determine the energy and pressure:
1046 * at nstcalcenergy steps and at energy output steps (set below).
1048 if (EI_VV(ir->eI) && (!bInitStep))
1050 /* for vv, the first half of the integration actually corresponds
1051 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1052 but the virial needs to be calculated on both the current step and the 'next' step. Future
1053 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1055 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1056 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1057 bCalcVir = bCalcEnerStep ||
1058 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1060 else
1062 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1063 bCalcVir = bCalcEnerStep ||
1064 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1066 bCalcEner = bCalcEnerStep;
1068 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1070 if (do_ene || do_log || bDoReplEx)
1072 bCalcVir = TRUE;
1073 bCalcEner = TRUE;
1076 /* Do we need global communication ? */
1077 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1078 do_per_step(step, nstglobalcomm) ||
1079 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1081 force_flags = (GMX_FORCE_STATECHANGED |
1082 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1083 GMX_FORCE_ALLFORCES |
1084 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1085 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1086 (bDoFEP ? GMX_FORCE_DHDL : 0)
1089 if (shellfc)
1091 /* Now is the time to relax the shells */
1092 relax_shell_flexcon(fplog, cr, bVerbose, step,
1093 ir, bNS, force_flags, top,
1094 constr, enerd, fcd,
1095 state, &f, force_vir, mdatoms,
1096 nrnb, wcycle, graph, groups,
1097 shellfc, fr, bBornRadii, t, mu_tot,
1098 vsite);
1100 else
1102 /* The coordinates (x) are shifted (to get whole molecules)
1103 * in do_force.
1104 * This is parallellized as well, and does communication too.
1105 * Check comments in sim_util.c
1107 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1108 state->box, &state->x, &state->hist,
1109 &f, force_vir, mdatoms, enerd, fcd,
1110 state->lambda, graph,
1111 fr, vsite, mu_tot, t, ed, bBornRadii,
1112 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1115 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1116 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1118 rvec *vbuf = nullptr;
1120 wallcycle_start(wcycle, ewcUPDATE);
1121 if (ir->eI == eiVV && bInitStep)
1123 /* if using velocity verlet with full time step Ekin,
1124 * take the first half step only to compute the
1125 * virial for the first step. From there,
1126 * revert back to the initial coordinates
1127 * so that the input is actually the initial step.
1129 snew(vbuf, state->natoms);
1130 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1132 else
1134 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1135 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1138 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1139 ekind, M, upd, etrtVELOCITY1,
1140 cr, constr);
1142 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1144 wallcycle_stop(wcycle, ewcUPDATE);
1145 update_constraints(fplog, step, nullptr, ir, mdatoms,
1146 state, fr->bMolPBC, graph, &f,
1147 &top->idef, shake_vir,
1148 cr, nrnb, wcycle, upd, constr,
1149 TRUE, bCalcVir);
1150 wallcycle_start(wcycle, ewcUPDATE);
1152 else if (graph)
1154 /* Need to unshift here if a do_force has been
1155 called in the previous step */
1156 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1158 /* if VV, compute the pressure and constraints */
1159 /* For VV2, we strictly only need this if using pressure
1160 * control, but we really would like to have accurate pressures
1161 * printed out.
1162 * Think about ways around this in the future?
1163 * For now, keep this choice in comments.
1165 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1166 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1167 bPres = TRUE;
1168 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1169 if (bCalcEner && ir->eI == eiVVAK)
1171 bSumEkinhOld = TRUE;
1173 /* for vv, the first half of the integration actually corresponds to the previous step.
1174 So we need information from the last step in the first half of the integration */
1175 if (bGStat || do_per_step(step-1, nstglobalcomm))
1177 wallcycle_stop(wcycle, ewcUPDATE);
1178 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1179 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1180 constr, &nullSignaller, state->box,
1181 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1182 (bGStat ? CGLO_GSTAT : 0)
1183 | CGLO_ENERGY
1184 | (bTemp ? CGLO_TEMPERATURE : 0)
1185 | (bPres ? CGLO_PRESSURE : 0)
1186 | (bPres ? CGLO_CONSTRAINT : 0)
1187 | (bStopCM ? CGLO_STOPCM : 0)
1188 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1189 | CGLO_SCALEEKIN
1191 /* explanation of above:
1192 a) We compute Ekin at the full time step
1193 if 1) we are using the AveVel Ekin, and it's not the
1194 initial step, or 2) if we are using AveEkin, but need the full
1195 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1196 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1197 EkinAveVel because it's needed for the pressure */
1198 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1199 top_global, top, state,
1200 &shouldCheckNumberOfBondedInteractions);
1201 wallcycle_start(wcycle, ewcUPDATE);
1203 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1204 if (!bInitStep)
1206 if (bTrotter)
1208 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1209 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1211 /* TODO This is only needed when we're about to write
1212 * a checkpoint, because we use it after the restart
1213 * (in a kludge?). But what should we be doing if
1214 * startingFromCheckpoint or bInitStep are true? */
1215 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1217 copy_mat(shake_vir, state->svir_prev);
1218 copy_mat(force_vir, state->fvir_prev);
1220 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1222 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1223 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1224 enerd->term[F_EKIN] = trace(ekind->ekin);
1227 else if (bExchanged)
1229 wallcycle_stop(wcycle, ewcUPDATE);
1230 /* We need the kinetic energy at minus the half step for determining
1231 * the full step kinetic energy and possibly for T-coupling.*/
1232 /* This may not be quite working correctly yet . . . . */
1233 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1234 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1235 constr, &nullSignaller, state->box,
1236 nullptr, &bSumEkinhOld,
1237 CGLO_GSTAT | CGLO_TEMPERATURE);
1238 wallcycle_start(wcycle, ewcUPDATE);
1241 /* if it's the initial step, we performed this first step just to get the constraint virial */
1242 if (ir->eI == eiVV && bInitStep)
1244 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1245 sfree(vbuf);
1247 wallcycle_stop(wcycle, ewcUPDATE);
1250 /* compute the conserved quantity */
1251 if (EI_VV(ir->eI))
1253 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1254 if (ir->eI == eiVV)
1256 last_ekin = enerd->term[F_EKIN];
1258 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1260 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1262 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1263 if (ir->efep != efepNO && !bRerunMD)
1265 sum_dhdl(enerd, state->lambda, ir->fepvals);
1269 /* ######## END FIRST UPDATE STEP ############## */
1270 /* ######## If doing VV, we now have v(dt) ###### */
1271 if (bDoExpanded)
1273 /* perform extended ensemble sampling in lambda - we don't
1274 actually move to the new state before outputting
1275 statistics, but if performing simulated tempering, we
1276 do update the velocities and the tau_t. */
1278 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1279 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1280 copy_df_history(state_global->dfhist, state->dfhist);
1283 /* Now we have the energies and forces corresponding to the
1284 * coordinates at time t. We must output all of this before
1285 * the update.
1287 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1288 ir, state, state_global, observablesHistory,
1289 top_global, fr,
1290 outf, mdebin, ekind, &f,
1291 &nchkpt,
1292 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1293 bSumEkinhOld);
1294 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1295 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1297 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1298 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1300 copy_mat(state->svir_prev, shake_vir);
1301 copy_mat(state->fvir_prev, force_vir);
1304 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1306 /* Check whether everything is still allright */
1307 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1308 #if GMX_THREAD_MPI
1309 && MASTER(cr)
1310 #endif
1313 int nsteps_stop = -1;
1315 /* this just makes signals[].sig compatible with the hack
1316 of sending signals around by MPI_Reduce together with
1317 other floats */
1318 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1320 signals[eglsSTOPCOND].sig = 1;
1321 nsteps_stop = std::max(ir->nstlist, 2*nstglobalcomm);
1323 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1325 signals[eglsSTOPCOND].sig = -1;
1326 nsteps_stop = nstglobalcomm + 1;
1328 if (fplog)
1330 fprintf(fplog,
1331 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1332 gmx_get_signal_name(), nsteps_stop);
1333 fflush(fplog);
1335 fprintf(stderr,
1336 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1337 gmx_get_signal_name(), nsteps_stop);
1338 fflush(stderr);
1339 handled_stop_condition = (int)gmx_get_stop_condition();
1341 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1342 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1343 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1345 /* Signal to terminate the run */
1346 signals[eglsSTOPCOND].sig = 1;
1347 if (fplog)
1349 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1351 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1354 if (bResetCountersHalfMaxH && MASTER(cr) &&
1355 elapsed_time > max_hours*60.0*60.0*0.495)
1357 /* Set flag that will communicate the signal to all ranks in the simulation */
1358 signals[eglsRESETCOUNTERS].sig = 1;
1361 /* In parallel we only have to check for checkpointing in steps
1362 * where we do global communication,
1363 * otherwise the other nodes don't know.
1365 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1366 cpt_period >= 0 &&
1367 (cpt_period == 0 ||
1368 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1369 signals[eglsCHKPT].set == 0)
1371 signals[eglsCHKPT].sig = 1;
1374 /* ######### START SECOND UPDATE STEP ################# */
1376 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1377 in preprocessing */
1379 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1381 gmx_bool bIfRandomize;
1382 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1383 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1384 if (constr && bIfRandomize)
1386 update_constraints(fplog, step, nullptr, ir, mdatoms,
1387 state, fr->bMolPBC, graph, &f,
1388 &top->idef, tmp_vir,
1389 cr, nrnb, wcycle, upd, constr,
1390 TRUE, bCalcVir);
1393 /* Box is changed in update() when we do pressure coupling,
1394 * but we should still use the old box for energy corrections and when
1395 * writing it to the energy file, so it matches the trajectory files for
1396 * the same timestep above. Make a copy in a separate array.
1398 copy_mat(state->box, lastbox);
1400 dvdl_constr = 0;
1402 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1404 wallcycle_start(wcycle, ewcUPDATE);
1405 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1406 if (bTrotter)
1408 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1409 /* We can only do Berendsen coupling after we have summed
1410 * the kinetic energy or virial. Since the happens
1411 * in global_state after update, we should only do it at
1412 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1415 else
1417 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1418 update_pcouple_before_coordinates(fplog, step, ir, state,
1419 parrinellorahmanMu, M,
1420 bInitStep);
1423 if (EI_VV(ir->eI))
1425 /* velocity half-step update */
1426 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1427 ekind, M, upd, etrtVELOCITY2,
1428 cr, constr);
1431 /* Above, initialize just copies ekinh into ekin,
1432 * it doesn't copy position (for VV),
1433 * and entire integrator for MD.
1436 if (ir->eI == eiVVAK)
1438 /* We probably only need md->homenr, not state->natoms */
1439 if (state->natoms > cbuf_nalloc)
1441 cbuf_nalloc = state->natoms;
1442 srenew(cbuf, cbuf_nalloc);
1444 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1447 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1448 ekind, M, upd, etrtPOSITION, cr, constr);
1449 wallcycle_stop(wcycle, ewcUPDATE);
1451 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1452 fr->bMolPBC, graph, &f,
1453 &top->idef, shake_vir,
1454 cr, nrnb, wcycle, upd, constr,
1455 FALSE, bCalcVir);
1457 if (ir->eI == eiVVAK)
1459 /* erase F_EKIN and F_TEMP here? */
1460 /* just compute the kinetic energy at the half step to perform a trotter step */
1461 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1462 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1463 constr, &nullSignaller, lastbox,
1464 nullptr, &bSumEkinhOld,
1465 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1467 wallcycle_start(wcycle, ewcUPDATE);
1468 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1469 /* now we know the scaling, we can compute the positions again again */
1470 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1472 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1473 ekind, M, upd, etrtPOSITION, cr, constr);
1474 wallcycle_stop(wcycle, ewcUPDATE);
1476 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1477 /* are the small terms in the shake_vir here due
1478 * to numerical errors, or are they important
1479 * physically? I'm thinking they are just errors, but not completely sure.
1480 * For now, will call without actually constraining, constr=NULL*/
1481 update_constraints(fplog, step, nullptr, ir, mdatoms,
1482 state, fr->bMolPBC, graph, &f,
1483 &top->idef, tmp_vir,
1484 cr, nrnb, wcycle, upd, nullptr,
1485 FALSE, bCalcVir);
1487 if (EI_VV(ir->eI))
1489 /* this factor or 2 correction is necessary
1490 because half of the constraint force is removed
1491 in the vv step, so we have to double it. See
1492 the Redmine issue #1255. It is not yet clear
1493 if the factor of 2 is exact, or just a very
1494 good approximation, and this will be
1495 investigated. The next step is to see if this
1496 can be done adding a dhdl contribution from the
1497 rattle step, but this is somewhat more
1498 complicated with the current code. Will be
1499 investigated, hopefully for 4.6.3. However,
1500 this current solution is much better than
1501 having it completely wrong.
1503 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1505 else
1507 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1510 else if (graph)
1512 /* Need to unshift here */
1513 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1516 if (vsite != nullptr)
1518 wallcycle_start(wcycle, ewcVSITECONSTR);
1519 if (graph != nullptr)
1521 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1523 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1524 top->idef.iparams, top->idef.il,
1525 fr->ePBC, fr->bMolPBC, cr, state->box);
1527 if (graph != nullptr)
1529 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1531 wallcycle_stop(wcycle, ewcVSITECONSTR);
1534 /* ############## IF NOT VV, Calculate globals HERE ############ */
1535 /* With Leap-Frog we can skip compute_globals at
1536 * non-communication steps, but we need to calculate
1537 * the kinetic energy one step before communication.
1540 // Organize to do inter-simulation signalling on steps if
1541 // and when algorithms require it.
1542 bool doInterSimSignal = (!bFirstStep && bDoReplEx) || bUsingEnsembleRestraints;
1544 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1546 // Since we're already communicating at this step, we
1547 // can propagate intra-simulation signals. Note that
1548 // check_nstglobalcomm has the responsibility for
1549 // choosing the value of nstglobalcomm that is one way
1550 // bGStat becomes true, so we can't get into a
1551 // situation where e.g. checkpointing can't be
1552 // signalled.
1553 bool doIntraSimSignal = true;
1554 SimulationSignaller signaller(&signals, cr, doInterSimSignal, doIntraSimSignal);
1556 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1557 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1558 constr, &signaller,
1559 lastbox,
1560 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1561 (bGStat ? CGLO_GSTAT : 0)
1562 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1563 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1564 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1565 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1566 | CGLO_CONSTRAINT
1567 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1569 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1570 top_global, top, state,
1571 &shouldCheckNumberOfBondedInteractions);
1575 /* ############# END CALC EKIN AND PRESSURE ################# */
1577 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1578 the virial that should probably be addressed eventually. state->veta has better properies,
1579 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1580 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1582 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1584 /* Sum up the foreign energy and dhdl terms for md and sd.
1585 Currently done every step so that dhdl is correct in the .edr */
1586 sum_dhdl(enerd, state->lambda, ir->fepvals);
1589 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1590 pres, force_vir, shake_vir,
1591 parrinellorahmanMu,
1592 state, nrnb, upd);
1594 /* ################# END UPDATE STEP 2 ################# */
1595 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1597 /* The coordinates (x) were unshifted in update */
1598 if (!bGStat)
1600 /* We will not sum ekinh_old,
1601 * so signal that we still have to do it.
1603 bSumEkinhOld = TRUE;
1606 if (bCalcEner)
1608 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1610 /* use the directly determined last velocity, not actually the averaged half steps */
1611 if (bTrotter && ir->eI == eiVV)
1613 enerd->term[F_EKIN] = last_ekin;
1615 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1617 if (integratorHasConservedEnergyQuantity(ir))
1619 if (EI_VV(ir->eI))
1621 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1623 else
1625 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1628 /* ######### END PREPARING EDR OUTPUT ########### */
1631 /* Output stuff */
1632 if (MASTER(cr))
1634 if (fplog && do_log && bDoExpanded)
1636 /* only needed if doing expanded ensemble */
1637 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1638 state_global->dfhist, state->fep_state, ir->nstlog, step);
1640 if (bCalcEner)
1642 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1643 t, mdatoms->tmass, enerd, state,
1644 ir->fepvals, ir->expandedvals, lastbox,
1645 shake_vir, force_vir, total_vir, pres,
1646 ekind, mu_tot, constr);
1648 else
1650 upd_mdebin_step(mdebin);
1653 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1654 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1656 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1657 step, t,
1658 eprNORMAL, mdebin, fcd, groups, &(ir->opts));
1660 if (ir->bPull)
1662 pull_print_output(ir->pull_work, step, t);
1665 if (do_per_step(step, ir->nstlog))
1667 if (fflush(fplog) != 0)
1669 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1673 if (bDoExpanded)
1675 /* Have to do this part _after_ outputting the logfile and the edr file */
1676 /* Gets written into the state at the beginning of next loop*/
1677 state->fep_state = lamnew;
1679 /* Print the remaining wall clock time for the run */
1680 if (MULTIMASTER(cr) &&
1681 (do_verbose || gmx_got_usr_signal()) &&
1682 !bPMETunePrinting)
1684 if (shellfc)
1686 fprintf(stderr, "\n");
1688 print_time(stderr, walltime_accounting, step, ir, cr);
1691 /* Ion/water position swapping.
1692 * Not done in last step since trajectory writing happens before this call
1693 * in the MD loop and exchanges would be lost anyway. */
1694 bNeedRepartition = FALSE;
1695 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1696 do_per_step(step, ir->swap->nstswap))
1698 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1699 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1700 bRerunMD ? rerun_fr.box : state->box,
1701 MASTER(cr) && bVerbose, bRerunMD);
1703 if (bNeedRepartition && DOMAINDECOMP(cr))
1705 dd_collect_state(cr->dd, state, state_global);
1709 /* Replica exchange */
1710 bExchanged = FALSE;
1711 if (bDoReplEx)
1713 bExchanged = replica_exchange(fplog, cr, repl_ex,
1714 state_global, enerd,
1715 state, step, t);
1718 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1720 dd_partition_system(fplog, step, cr, TRUE, 1,
1721 state_global, top_global, ir,
1722 state, &f, mdatoms, top, fr,
1723 vsite, constr,
1724 nrnb, wcycle, FALSE);
1725 shouldCheckNumberOfBondedInteractions = true;
1726 update_realloc(upd, state->natoms);
1729 bFirstStep = FALSE;
1730 bInitStep = FALSE;
1731 startingFromCheckpoint = FALSE;
1733 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1734 /* With all integrators, except VV, we need to retain the pressure
1735 * at the current step for coupling at the next step.
1737 if ((state->flags & (1<<estPRES_PREV)) &&
1738 (bGStatEveryStep ||
1739 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1741 /* Store the pressure in t_state for pressure coupling
1742 * at the next MD step.
1744 copy_mat(pres, state->pres_prev);
1747 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1749 if ( (membed != nullptr) && (!bLastStep) )
1751 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1754 if (bRerunMD)
1756 if (MASTER(cr))
1758 /* read next frame from input trajectory */
1759 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1762 if (PAR(cr))
1764 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1768 cycles = wallcycle_stop(wcycle, ewcSTEP);
1769 if (DOMAINDECOMP(cr) && wcycle)
1771 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1774 if (!bRerunMD || !rerun_fr.bStep)
1776 /* increase the MD step number */
1777 step++;
1778 step_rel++;
1781 /* TODO make a counter-reset module */
1782 /* If it is time to reset counters, set a flag that remains
1783 true until counters actually get reset */
1784 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1785 signals[eglsRESETCOUNTERS].set != 0)
1787 if (pme_loadbal_is_active(pme_loadbal))
1789 /* Do not permit counter reset while PME load
1790 * balancing is active. The only purpose for resetting
1791 * counters is to measure reliable performance data,
1792 * and that can't be done before balancing
1793 * completes.
1795 * TODO consider fixing this by delaying the reset
1796 * until after load balancing completes,
1797 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1798 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1799 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1800 "resetting counters later in the run, e.g. with gmx "
1801 "mdrun -resetstep.", step);
1803 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1804 use_GPU(fr->nbv) ? fr->nbv : nullptr);
1805 wcycle_set_reset_counters(wcycle, -1);
1806 if (!(cr->duty & DUTY_PME))
1808 /* Tell our PME node to reset its counters */
1809 gmx_pme_send_resetcounters(cr, step);
1811 /* Correct max_hours for the elapsed time */
1812 max_hours -= elapsed_time/(60.0*60.0);
1813 /* If mdrun -maxh -resethway was active, it can only trigger once */
1814 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1815 /* Reset can only happen once, so clear the triggering flag. */
1816 signals[eglsRESETCOUNTERS].set = 0;
1819 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1820 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1823 /* End of main MD loop */
1825 /* Closing TNG files can include compressing data. Therefore it is good to do that
1826 * before stopping the time measurements. */
1827 mdoutf_tng_close(outf);
1829 /* Stop measuring walltime */
1830 walltime_accounting_end(walltime_accounting);
1832 if (bRerunMD && MASTER(cr))
1834 close_trx(status);
1837 if (!(cr->duty & DUTY_PME))
1839 /* Tell the PME only node to finish */
1840 gmx_pme_send_finish(cr);
1843 if (MASTER(cr))
1845 if (ir->nstcalcenergy > 0 && !bRerunMD)
1847 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1848 eprAVER, mdebin, fcd, groups, &(ir->opts));
1852 done_mdoutf(outf);
1854 if (bPMETune)
1856 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1859 done_shellfc(fplog, shellfc, step_rel);
1861 if (repl_ex_nst > 0 && MASTER(cr))
1863 print_replica_exchange_statistics(fplog, repl_ex);
1866 // Clean up swapcoords
1867 if (ir->eSwapCoords != eswapNO)
1869 finish_swapcoords(ir->swap);
1872 /* Do essential dynamics cleanup if needed. Close .edo file */
1873 done_ed(&ed);
1875 /* IMD cleanup, if bIMD is TRUE. */
1876 IMD_finalize(ir->bIMD, ir->imd);
1878 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1879 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1880 signals[eglsRESETCOUNTERS].set == 0 &&
1881 !bResetCountersHalfMaxH)
1883 walltime_accounting_set_valid_finish(walltime_accounting);
1886 return 0;