Update Awh initialization and lifetime management
[gromacs.git] / src / gromacs / mdrun / md.cpp
blob2c1ff9962ebd2585c11be16a88532e6436553c14
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,2018, 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 /*! \internal \file
39 * \brief Implements the integrator for normal molecular dynamics simulations
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
44 #include "gmxpre.h"
46 #include "config.h"
48 #include <stdio.h>
49 #include <stdlib.h>
51 #include <cmath>
53 #include <algorithm>
54 #include <memory>
56 #include "thread_mpi/threads.h"
58 #include "gromacs/awh/awh.h"
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/compat/make_unique.h"
61 #include "gromacs/domdec/collect.h"
62 #include "gromacs/domdec/domdec.h"
63 #include "gromacs/domdec/domdec_network.h"
64 #include "gromacs/domdec/domdec_struct.h"
65 #include "gromacs/essentialdynamics/edsam.h"
66 #include "gromacs/ewald/pme.h"
67 #include "gromacs/ewald/pme-load-balancing.h"
68 #include "gromacs/fileio/trxio.h"
69 #include "gromacs/gmxlib/network.h"
70 #include "gromacs/gmxlib/nrnb.h"
71 #include "gromacs/gpu_utils/gpu_utils.h"
72 #include "gromacs/imd/imd.h"
73 #include "gromacs/listed-forces/manage-threading.h"
74 #include "gromacs/math/functions.h"
75 #include "gromacs/math/utilities.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/math/vectypes.h"
78 #include "gromacs/mdlib/compute_io.h"
79 #include "gromacs/mdlib/constr.h"
80 #include "gromacs/mdlib/deform.h"
81 #include "gromacs/mdlib/ebin.h"
82 #include "gromacs/mdlib/expanded.h"
83 #include "gromacs/mdlib/force.h"
84 #include "gromacs/mdlib/force_flags.h"
85 #include "gromacs/mdlib/forcerec.h"
86 #include "gromacs/mdlib/md_support.h"
87 #include "gromacs/mdlib/mdatoms.h"
88 #include "gromacs/mdlib/mdebin.h"
89 #include "gromacs/mdlib/mdoutf.h"
90 #include "gromacs/mdlib/mdrun.h"
91 #include "gromacs/mdlib/mdsetup.h"
92 #include "gromacs/mdlib/membed.h"
93 #include "gromacs/mdlib/nb_verlet.h"
94 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
95 #include "gromacs/mdlib/ns.h"
96 #include "gromacs/mdlib/repl_ex.h"
97 #include "gromacs/mdlib/shellfc.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/sim_util.h"
100 #include "gromacs/mdlib/simulationsignal.h"
101 #include "gromacs/mdlib/tgroup.h"
102 #include "gromacs/mdlib/trajectory_writing.h"
103 #include "gromacs/mdlib/update.h"
104 #include "gromacs/mdlib/vcm.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdtypes/awh-history.h"
107 #include "gromacs/mdtypes/awh-params.h"
108 #include "gromacs/mdtypes/commrec.h"
109 #include "gromacs/mdtypes/df_history.h"
110 #include "gromacs/mdtypes/energyhistory.h"
111 #include "gromacs/mdtypes/fcdata.h"
112 #include "gromacs/mdtypes/forcerec.h"
113 #include "gromacs/mdtypes/group.h"
114 #include "gromacs/mdtypes/inputrec.h"
115 #include "gromacs/mdtypes/interaction_const.h"
116 #include "gromacs/mdtypes/md_enums.h"
117 #include "gromacs/mdtypes/mdatom.h"
118 #include "gromacs/mdtypes/observableshistory.h"
119 #include "gromacs/mdtypes/state.h"
120 #include "gromacs/pbcutil/mshift.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/swap/swapcoords.h"
124 #include "gromacs/timing/wallcycle.h"
125 #include "gromacs/timing/walltime_accounting.h"
126 #include "gromacs/topology/atoms.h"
127 #include "gromacs/topology/idef.h"
128 #include "gromacs/topology/mtop_util.h"
129 #include "gromacs/topology/topology.h"
130 #include "gromacs/trajectory/trajectoryframe.h"
131 #include "gromacs/utility/basedefinitions.h"
132 #include "gromacs/utility/cstringutil.h"
133 #include "gromacs/utility/fatalerror.h"
134 #include "gromacs/utility/logger.h"
135 #include "gromacs/utility/real.h"
136 #include "gromacs/utility/smalloc.h"
138 #include "integrator.h"
140 #ifdef GMX_FAHCORE
141 #include "corewrap.h"
142 #endif
144 using gmx::SimulationSignaller;
146 //! Resets all the counters.
147 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
148 gmx_int64_t step,
149 gmx_int64_t *step_rel, t_inputrec *ir,
150 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
151 gmx_walltime_accounting_t walltime_accounting,
152 struct nonbonded_verlet_t *nbv,
153 struct gmx_pme_t *pme)
155 char sbuf[STEPSTRSIZE];
157 /* Reset all the counters related to performance over the run */
158 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
159 "step %s: resetting all time and cycle counters",
160 gmx_step_str(step, sbuf));
162 if (use_GPU(nbv))
164 nbnxn_gpu_reset_timings(nbv);
167 if (pme_gpu_task_enabled(pme))
169 pme_gpu_reset_timings(pme);
172 if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
174 resetGpuProfiler();
177 wallcycle_stop(wcycle, ewcRUN);
178 wallcycle_reset_all(wcycle);
179 if (DOMAINDECOMP(cr))
181 reset_dd_statistics_counters(cr->dd);
183 init_nrnb(nrnb);
184 ir->init_step += *step_rel;
185 ir->nsteps -= *step_rel;
186 *step_rel = 0;
187 wallcycle_start(wcycle, ewcRUN);
188 walltime_accounting_start(walltime_accounting);
189 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
192 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
194 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
195 * \param[in,out] globalState The global state container
196 * \param[in] constructVsites When true, vsite coordinates are constructed
197 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
198 * \param[in] idef Topology parameters, used for constructing vsites
199 * \param[in] timeStep Time step, used for constructing vsites
200 * \param[in] forceRec Force record, used for constructing vsites
201 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
202 * \param[in,out] warnWhenNoV When true, issue a warning when no velocities are present in \p rerunFrame; is set to false when a warning was issued
204 static void prepareRerunState(const t_trxframe &rerunFrame,
205 t_state *globalState,
206 bool constructVsites,
207 const gmx_vsite_t *vsite,
208 const t_idef &idef,
209 double timeStep,
210 const t_forcerec &forceRec,
211 t_graph *graph,
212 gmx_bool *warnWhenNoV)
214 for (int i = 0; i < globalState->natoms; i++)
216 copy_rvec(rerunFrame.x[i], globalState->x[i]);
218 if (rerunFrame.bV)
220 for (int i = 0; i < globalState->natoms; i++)
222 copy_rvec(rerunFrame.v[i], globalState->v[i]);
225 else
227 for (int i = 0; i < globalState->natoms; i++)
229 clear_rvec(globalState->v[i]);
231 if (*warnWhenNoV)
233 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
234 " Ekin, temperature and pressure are incorrect,\n"
235 " the virial will be incorrect when constraints are present.\n"
236 "\n");
237 *warnWhenNoV = FALSE;
240 copy_mat(rerunFrame.box, globalState->box);
242 if (constructVsites)
244 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
246 if (graph)
248 /* Following is necessary because the graph may get out of sync
249 * with the coordinates if we only have every N'th coordinate set
251 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
252 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
254 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
255 idef.iparams, idef.il,
256 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
257 if (graph)
259 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
264 void gmx::Integrator::do_md()
266 // TODO Historically, the EM and MD "integrators" used different
267 // names for the t_inputrec *parameter, but these must have the
268 // same name, now that it's a member of a struct. We use this ir
269 // alias to avoid a large ripple of nearly useless changes.
270 // t_inputrec is being replaced by IMdpOptionsProvider, so this
271 // will go away eventually.
272 t_inputrec *ir = inputrec;
273 gmx_mdoutf *outf = nullptr;
274 gmx_int64_t step, step_rel;
275 double elapsed_time;
276 double t, t0, lam0[efptNR];
277 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
278 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
279 bFirstStep, bInitStep, bLastStep = FALSE;
280 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
281 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
282 bForceUpdate = FALSE, bCPT;
283 gmx_bool bMasterState;
284 int force_flags, cglo_flags;
285 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
286 int i, m;
287 t_trxstatus *status;
288 rvec mu_tot;
289 t_vcm *vcm;
290 matrix parrinellorahmanMu, M;
291 t_trxframe rerun_fr;
292 gmx_repl_ex_t repl_ex = nullptr;
293 int nchkpt = 1;
294 gmx_localtop_t *top;
295 t_mdebin *mdebin = nullptr;
296 gmx_enerdata_t *enerd;
297 PaddedRVecVector f {};
298 gmx_global_stat_t gstat;
299 gmx_update_t *upd = nullptr;
300 t_graph *graph = nullptr;
301 gmx_groups_t *groups;
302 gmx_ekindata_t *ekind;
303 gmx_shellfc_t *shellfc;
304 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
305 gmx_bool bResetCountersHalfMaxH = FALSE;
306 gmx_bool bTemp, bPres, bTrotter;
307 real dvdl_constr;
308 rvec *cbuf = nullptr;
309 int cbuf_nalloc = 0;
310 matrix lastbox;
311 int lamnew = 0;
312 /* for FEP */
313 int nstfep = 0;
314 double cycles;
315 real saved_conserved_quantity = 0;
316 real last_ekin = 0;
317 t_extmass MassQ;
318 int **trotter_seq;
319 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
320 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
323 /* PME load balancing data for GPU kernels */
324 pme_load_balancing_t *pme_loadbal = nullptr;
325 gmx_bool bPMETune = FALSE;
326 gmx_bool bPMETunePrinting = FALSE;
328 /* Interactive MD */
329 gmx_bool bIMDstep = FALSE;
331 #ifdef GMX_FAHCORE
332 /* Temporary addition for FAHCORE checkpointing */
333 int chkpt_ret;
334 #endif
335 /* Domain decomposition could incorrectly miss a bonded
336 interaction, but checking for that requires a global
337 communication stage, which does not otherwise happen in DD
338 code. So we do that alongside the first global energy reduction
339 after a new DD is made. These variables handle whether the
340 check happens, and the result it returns. */
341 bool shouldCheckNumberOfBondedInteractions = false;
342 int totalNumberOfBondedInteractions = -1;
344 SimulationSignals signals;
345 // Most global communnication stages don't propagate mdrun
346 // signals, and will use this object to achieve that.
347 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
349 if (!mdrunOptions.writeConfout)
351 // This is on by default, and the main known use case for
352 // turning it off is for convenience in benchmarking, which is
353 // something that should not show up in the general user
354 // interface.
355 GMX_LOG(mdlog.info).asParagraph().
356 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
358 if (mdrunOptions.timingOptions.resetHalfway)
360 GMX_LOG(mdlog.info).asParagraph().
361 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
362 if (ir->nsteps > 0)
364 /* Signal to reset the counters half the simulation steps. */
365 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
367 /* Signal to reset the counters halfway the simulation time. */
368 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
371 /* md-vv uses averaged full step velocities for T-control
372 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
373 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
374 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
376 const gmx_bool bRerunMD = mdrunOptions.rerun;
377 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
379 if (bRerunMD)
381 /* Since we don't know if the frames read are related in any way,
382 * rebuild the neighborlist at every step.
384 ir->nstlist = 1;
385 ir->nstcalcenergy = 1;
386 nstglobalcomm = 1;
389 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
390 bGStatEveryStep = (nstglobalcomm == 1);
392 if (bRerunMD)
394 ir->nstxout_compressed = 0;
396 groups = &top_global->groups;
398 gmx_edsam *ed = nullptr;
399 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
401 /* Initialize essential dynamics sampling */
402 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
403 top_global,
404 ir, cr, constr,
405 state_global, observablesHistory,
406 oenv, mdrunOptions.continuationOptions.appendFiles);
409 if (ir->eSwapCoords != eswapNO)
411 /* Initialize ion swapping code */
412 init_swapcoords(fplog, ir, opt2fn_master("-swap", nfile, fnm, cr),
413 top_global,
414 state_global, observablesHistory,
415 cr, oenv, mdrunOptions);
418 /* Initial values */
419 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
420 &t, &t0, state_global, lam0,
421 nrnb, top_global, &upd,
422 nfile, fnm, &outf, &mdebin,
423 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, wcycle);
425 clear_mat(total_vir);
426 clear_mat(pres);
427 /* Energy terms and groups */
428 snew(enerd, 1);
429 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
430 enerd);
432 /* Kinetic energy data */
433 snew(ekind, 1);
434 init_ekindata(fplog, top_global, &(ir->opts), ekind);
435 /* Copy the cos acceleration to the groups struct */
436 ekind->cosacc.cos_accel = ir->cos_accel;
438 gstat = global_stat_init(ir);
440 /* Check for polarizable models and flexible constraints */
441 shellfc = init_shell_flexcon(fplog,
442 top_global, n_flexible_constraints(constr),
443 ir->nstcalcenergy, DOMAINDECOMP(cr));
445 if (inputrecDeform(ir))
447 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
448 set_deform_reference_box(upd,
449 deform_init_init_step_tpx,
450 deform_init_box_tpx);
451 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
455 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
456 if ((io > 2000) && MASTER(cr))
458 fprintf(stderr,
459 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
460 io);
464 // Local state only becomes valid now.
465 std::unique_ptr<t_state> stateInstance;
466 t_state * state;
468 if (DOMAINDECOMP(cr))
470 top = dd_init_local_top(top_global);
472 stateInstance = gmx::compat::make_unique<t_state>();
473 state = stateInstance.get();
474 dd_init_local_state(cr->dd, state_global, state);
476 else
478 state_change_natoms(state_global, state_global->natoms);
479 /* We need to allocate one element extra, since we might use
480 * (unaligned) 4-wide SIMD loads to access rvec entries.
482 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
483 /* Copy the pointer to the global state */
484 state = state_global;
486 snew(top, 1);
487 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
488 &graph, mdAtoms, vsite, shellfc);
490 update_realloc(upd, state->natoms);
493 /* Set up interactive MD (IMD) */
494 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
495 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
496 nfile, fnm, oenv, mdrunOptions);
498 if (DOMAINDECOMP(cr))
500 /* Distribute the charge groups over the nodes from the master node */
501 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
502 state_global, top_global, ir,
503 state, &f, mdAtoms, top, fr,
504 vsite, constr,
505 nrnb, nullptr, FALSE);
506 shouldCheckNumberOfBondedInteractions = true;
507 update_realloc(upd, state->natoms);
510 auto mdatoms = mdAtoms->mdatoms();
512 // NOTE: The global state is no longer used at this point.
513 // But state_global is still used as temporary storage space for writing
514 // the global state to file and potentially for replica exchange.
515 // (Global topology should persist.)
517 update_mdatoms(mdatoms, state->lambda[efptMASS]);
519 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
520 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
522 if (ir->bExpanded)
524 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
527 if (MASTER(cr))
529 if (startingFromCheckpoint)
531 /* Update mdebin with energy history if appending to output files */
532 if (continuationOptions.appendFiles)
534 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
536 else if (observablesHistory->energyHistory.get() != nullptr)
538 /* We might have read an energy history from checkpoint.
539 * As we are not appending, we want to restart the statistics.
540 * Free the allocated memory and reset the counts.
542 observablesHistory->energyHistory = {};
545 if (observablesHistory->energyHistory.get() == nullptr)
547 observablesHistory->energyHistory = std::unique_ptr<energyhistory_t>(new energyhistory_t {});
549 /* Set the initial energy history in state by updating once */
550 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
553 /* Initialize constraints */
554 if (constr && !DOMAINDECOMP(cr))
556 set_constraints(constr, top, ir, mdatoms, cr);
559 // TODO: Remove this by converting AWH into a ForceProvider
560 auto awh = prepareAwhModule(fplog, *ir, state_global, cr, ms, startingFromCheckpoint,
561 shellfc != nullptr,
562 opt2fn("-awh", nfile, fnm), ir->pull_work);
564 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
565 if (useReplicaExchange && MASTER(cr))
567 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
568 replExParams);
570 /* PME tuning is only supported in the Verlet scheme, with PME for
571 * Coulomb. It is not supported with only LJ PME, or for
572 * reruns. */
573 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
574 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
575 if (bPMETune)
577 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
578 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
579 &bPMETunePrinting);
582 if (!ir->bContinuation && !bRerunMD)
584 if (state->flags & (1 << estV))
586 /* Set the velocities of vsites, shells and frozen atoms to zero */
587 for (i = 0; i < mdatoms->homenr; i++)
589 if (mdatoms->ptype[i] == eptVSite ||
590 mdatoms->ptype[i] == eptShell)
592 clear_rvec(state->v[i]);
594 else if (mdatoms->cFREEZE)
596 for (m = 0; m < DIM; m++)
598 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
600 state->v[i][m] = 0;
607 if (constr)
609 /* Constrain the initial coordinates and velocities */
610 do_constrain_first(fplog, constr, ir, mdatoms, state,
611 cr, ms, nrnb, fr, top);
613 if (vsite)
615 /* Construct the virtual sites for the initial configuration */
616 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
617 top->idef.iparams, top->idef.il,
618 fr->ePBC, fr->bMolPBC, cr, state->box);
622 if (ir->efep != efepNO)
624 /* Set free energy calculation frequency as the greatest common
625 * denominator of nstdhdl and repl_ex_nst. */
626 nstfep = ir->fepvals->nstdhdl;
627 if (ir->bExpanded)
629 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
631 if (useReplicaExchange)
633 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
637 /* Be REALLY careful about what flags you set here. You CANNOT assume
638 * this is the first step, since we might be restarting from a checkpoint,
639 * and in that case we should not do any modifications to the state.
641 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
643 if (continuationOptions.haveReadEkin)
645 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
648 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
649 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
650 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
651 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
653 bSumEkinhOld = FALSE;
654 /* To minimize communication, compute_globals computes the COM velocity
655 * and the kinetic energy for the velocities without COM motion removed.
656 * Thus to get the kinetic energy without the COM contribution, we need
657 * to call compute_globals twice.
659 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
661 int cglo_flags_iteration = cglo_flags;
662 if (bStopCM && cgloIteration == 0)
664 cglo_flags_iteration |= CGLO_STOPCM;
665 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
667 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
668 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
669 constr, &nullSignaller, state->box,
670 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
671 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
673 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
674 top_global, top, state,
675 &shouldCheckNumberOfBondedInteractions);
676 if (ir->eI == eiVVAK)
678 /* a second call to get the half step temperature initialized as well */
679 /* we do the same call as above, but turn the pressure off -- internally to
680 compute_globals, this is recognized as a velocity verlet half-step
681 kinetic energy calculation. This minimized excess variables, but
682 perhaps loses some logic?*/
684 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
685 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
686 constr, &nullSignaller, state->box,
687 nullptr, &bSumEkinhOld,
688 cglo_flags & ~CGLO_PRESSURE);
691 /* Calculate the initial half step temperature, and save the ekinh_old */
692 if (!continuationOptions.startedFromCheckpoint)
694 for (i = 0; (i < ir->opts.ngtc); i++)
696 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
700 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
701 temperature control */
702 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
704 if (MASTER(cr))
706 if (!ir->bContinuation)
708 if (constr && ir->eConstrAlg == econtLINCS)
710 fprintf(fplog,
711 "RMS relative constraint deviation after constraining: %.2e\n",
712 constr_rmsd(constr));
714 if (EI_STATE_VELOCITY(ir->eI))
716 real temp = enerd->term[F_TEMP];
717 if (ir->eI != eiVV)
719 /* Result of Ekin averaged over velocities of -half
720 * and +half step, while we only have -half step here.
722 temp *= 2;
724 fprintf(fplog, "Initial temperature: %g K\n", temp);
728 if (bRerunMD)
730 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
731 " input trajectory '%s'\n\n",
732 *(top_global->name), opt2fn("-rerun", nfile, fnm));
733 if (mdrunOptions.verbose)
735 fprintf(stderr, "Calculated time to finish depends on nsteps from "
736 "run input file,\nwhich may not correspond to the time "
737 "needed to process input trajectory.\n\n");
740 else
742 char tbuf[20];
743 fprintf(stderr, "starting mdrun '%s'\n",
744 *(top_global->name));
745 if (ir->nsteps >= 0)
747 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
749 else
751 sprintf(tbuf, "%s", "infinite");
753 if (ir->init_step > 0)
755 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
756 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
757 gmx_step_str(ir->init_step, sbuf2),
758 ir->init_step*ir->delta_t);
760 else
762 fprintf(stderr, "%s steps, %s ps.\n",
763 gmx_step_str(ir->nsteps, sbuf), tbuf);
766 fprintf(fplog, "\n");
769 walltime_accounting_start(walltime_accounting);
770 wallcycle_start(wcycle, ewcRUN);
771 print_start(fplog, cr, walltime_accounting, "mdrun");
773 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
774 #ifdef GMX_FAHCORE
775 chkpt_ret = fcCheckPointParallel( cr->nodeid,
776 NULL, 0);
777 if (chkpt_ret == 0)
779 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
781 #endif
783 /***********************************************************
785 * Loop over MD steps
787 ************************************************************/
789 /* if rerunMD then read coordinates and velocities from input trajectory */
790 if (bRerunMD)
792 if (getenv("GMX_FORCE_UPDATE"))
794 bForceUpdate = TRUE;
797 rerun_fr.natoms = 0;
798 if (MASTER(cr))
800 bLastStep = !read_first_frame(oenv, &status,
801 opt2fn("-rerun", nfile, fnm),
802 &rerun_fr, TRX_NEED_X | TRX_READ_V);
803 if (rerun_fr.natoms != top_global->natoms)
805 gmx_fatal(FARGS,
806 "Number of atoms in trajectory (%d) does not match the "
807 "run input file (%d)\n",
808 rerun_fr.natoms, top_global->natoms);
810 if (ir->ePBC != epbcNONE)
812 if (!rerun_fr.bBox)
814 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);
816 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
818 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
823 if (PAR(cr))
825 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
828 if (ir->ePBC != epbcNONE)
830 /* Set the shift vectors.
831 * Necessary here when have a static box different from the tpr box.
833 calc_shifts(rerun_fr.box, fr->shift_vec);
837 /* Loop over MD steps or if rerunMD to end of input trajectory,
838 * or, if max_hours>0, until max_hours is reached.
840 real max_hours = mdrunOptions.maximumHoursToRun;
841 bFirstStep = TRUE;
842 /* Skip the first Nose-Hoover integration when we get the state from tpx */
843 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
844 bSumEkinhOld = FALSE;
845 bExchanged = FALSE;
846 bNeedRepartition = FALSE;
848 bool simulationsShareState = false;
849 int nstSignalComm = nstglobalcomm;
851 // TODO This implementation of ensemble orientation restraints is nasty because
852 // a user can't just do multi-sim with single-sim orientation restraints.
853 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || (ms && fcd->orires.nr);
854 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && ms);
856 // Replica exchange, ensemble restraints and AWH need all
857 // simulations to remain synchronized, so they need
858 // checkpoints and stop conditions to act on the same step, so
859 // the propagation of such signals must take place between
860 // simulations, not just within simulations.
861 // TODO: Make algorithm initializers set these flags.
862 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
863 bool resetCountersIsLocal = true;
864 signals[eglsCHKPT] = SimulationSignal(!simulationsShareState);
865 signals[eglsSTOPCOND] = SimulationSignal(!simulationsShareState);
866 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
868 if (simulationsShareState)
870 // Inter-simulation signal communication does not need to happen
871 // often, so we use a minimum of 200 steps to reduce overhead.
872 const int c_minimumInterSimulationSignallingInterval = 200;
873 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
877 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
878 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
880 step = ir->init_step;
881 step_rel = 0;
883 // TODO extract this to new multi-simulation module
884 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
886 if (!multisim_int_all_are_equal(ms, ir->nsteps))
888 GMX_LOG(mdlog.warning).appendText(
889 "Note: The number of steps is not consistent across multi simulations,\n"
890 "but we are proceeding anyway!");
892 if (!multisim_int_all_are_equal(ms, ir->init_step))
894 GMX_LOG(mdlog.warning).appendText(
895 "Note: The initial step is not consistent across multi simulations,\n"
896 "but we are proceeding anyway!");
900 /* and stop now if we should */
901 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
902 while (!bLastStep)
905 /* Determine if this is a neighbor search step */
906 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
908 if (bPMETune && bNStList)
910 /* PME grid + cut-off optimization with GPUs or PME nodes */
911 pme_loadbal_do(pme_loadbal, cr,
912 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
913 fplog, mdlog,
914 ir, fr, state,
915 wcycle,
916 step, step_rel,
917 &bPMETunePrinting);
920 wallcycle_start(wcycle, ewcSTEP);
922 if (bRerunMD)
924 if (rerun_fr.bStep)
926 step = rerun_fr.step;
927 step_rel = step - ir->init_step;
929 if (rerun_fr.bTime)
931 t = rerun_fr.time;
933 else
935 t = step;
938 else
940 bLastStep = (step_rel == ir->nsteps);
941 t = t0 + step*ir->delta_t;
944 // TODO Refactor this, so that nstfep does not need a default value of zero
945 if (ir->efep != efepNO || ir->bSimTemp)
947 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
948 requiring different logic. */
949 if (bRerunMD)
951 if (MASTER(cr))
953 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
956 else
958 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
960 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
961 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
962 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
963 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
966 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
967 do_per_step(step, replExParams.exchangeInterval));
969 if (bSimAnn)
971 update_annealing_target_temp(ir, t, upd);
974 if (bRerunMD && MASTER(cr))
976 const bool constructVsites = (vsite && mdrunOptions.rerunConstructVsites);
977 if (constructVsites && DOMAINDECOMP(cr))
979 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
981 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
984 /* Stop Center of Mass motion */
985 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
987 if (bRerunMD)
989 /* for rerun MD always do Neighbour Searching */
990 bNS = (bFirstStep || ir->nstlist != 0);
992 else
994 /* Determine whether or not to do Neighbour Searching */
995 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
998 /* < 0 means stop at next step, > 0 means stop at next NS step */
999 if ( (signals[eglsSTOPCOND].set < 0) ||
1000 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
1002 bLastStep = TRUE;
1005 /* do_log triggers energy and virial calculation. Because this leads
1006 * to different code paths, forces can be different. Thus for exact
1007 * continuation we should avoid extra log output.
1008 * Note that the || bLastStep can result in non-exact continuation
1009 * beyond the last step. But we don't consider that to be an issue.
1011 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
1012 do_verbose = mdrunOptions.verbose &&
1013 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
1015 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
1017 if (bRerunMD)
1019 bMasterState = TRUE;
1021 else
1023 bMasterState = FALSE;
1024 /* Correct the new box if it is too skewed */
1025 if (inputrecDynamicBox(ir))
1027 if (correct_box(fplog, step, state->box, graph))
1029 bMasterState = TRUE;
1032 if (DOMAINDECOMP(cr) && bMasterState)
1034 dd_collect_state(cr->dd, state, state_global);
1038 if (DOMAINDECOMP(cr))
1040 /* Repartition the domain decomposition */
1041 dd_partition_system(fplog, step, cr,
1042 bMasterState, nstglobalcomm,
1043 state_global, top_global, ir,
1044 state, &f, mdAtoms, top, fr,
1045 vsite, constr,
1046 nrnb, wcycle,
1047 do_verbose && !bPMETunePrinting);
1048 shouldCheckNumberOfBondedInteractions = true;
1049 update_realloc(upd, state->natoms);
1053 if (MASTER(cr) && do_log)
1055 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1058 if (ir->efep != efepNO)
1060 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1063 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1066 /* We need the kinetic energy at minus the half step for determining
1067 * the full step kinetic energy and possibly for T-coupling.*/
1068 /* This may not be quite working correctly yet . . . . */
1069 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1070 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1071 constr, &nullSignaller, state->box,
1072 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1073 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1074 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1075 top_global, top, state,
1076 &shouldCheckNumberOfBondedInteractions);
1078 clear_mat(force_vir);
1080 /* We write a checkpoint at this MD step when:
1081 * either at an NS step when we signalled through gs,
1082 * or at the last step (but not when we do not want confout),
1083 * but never at the first step or with rerun.
1085 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1086 (bLastStep && mdrunOptions.writeConfout)) &&
1087 step > ir->init_step && !bRerunMD);
1088 if (bCPT)
1090 signals[eglsCHKPT].set = 0;
1093 /* Determine the energy and pressure:
1094 * at nstcalcenergy steps and at energy output steps (set below).
1096 if (EI_VV(ir->eI) && (!bInitStep))
1098 /* for vv, the first half of the integration actually corresponds
1099 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1100 but the virial needs to be calculated on both the current step and the 'next' step. Future
1101 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1103 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1104 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1105 bCalcVir = bCalcEnerStep ||
1106 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1108 else
1110 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1111 bCalcVir = bCalcEnerStep ||
1112 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1114 bCalcEner = bCalcEnerStep;
1116 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1118 if (do_ene || do_log || bDoReplEx)
1120 bCalcVir = TRUE;
1121 bCalcEner = TRUE;
1124 /* Do we need global communication ? */
1125 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1126 do_per_step(step, nstglobalcomm) ||
1127 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1129 force_flags = (GMX_FORCE_STATECHANGED |
1130 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1131 GMX_FORCE_ALLFORCES |
1132 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1133 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1134 (bDoFEP ? GMX_FORCE_DHDL : 0)
1137 if (shellfc)
1139 /* Now is the time to relax the shells */
1140 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, step,
1141 ir, bNS, force_flags, top,
1142 constr, enerd, fcd,
1143 state, f, force_vir, mdatoms,
1144 nrnb, wcycle, graph, groups,
1145 shellfc, fr, t, mu_tot,
1146 vsite,
1147 ddOpenBalanceRegion, ddCloseBalanceRegion);
1149 else
1151 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1152 (or the AWH update will be performed twice for one step when continuing). It would be best to
1153 call this update function from do_md_trajectory_writing but that would occur after do_force.
1154 One would have to divide the update_awh function into one function applying the AWH force
1155 and one doing the AWH bias update. The update AWH bias function could then be called after
1156 do_md_trajectory_writing (then containing update_awh_history).
1157 The checkpointing will in the future probably moved to the start of the md loop which will
1158 rid of this issue. */
1159 if (awh && bCPT && MASTER(cr))
1161 awh->updateHistory(state_global->awhHistory.get());
1164 /* The coordinates (x) are shifted (to get whole molecules)
1165 * in do_force.
1166 * This is parallellized as well, and does communication too.
1167 * Check comments in sim_util.c
1169 do_force(fplog, cr, ms, ir, awh.get(),
1170 step, nrnb, wcycle, top, groups,
1171 state->box, state->x, &state->hist,
1172 f, force_vir, mdatoms, enerd, fcd,
1173 state->lambda, graph,
1174 fr, vsite, mu_tot, t, ed,
1175 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1176 ddOpenBalanceRegion, ddCloseBalanceRegion);
1179 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1180 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1182 rvec *vbuf = nullptr;
1184 wallcycle_start(wcycle, ewcUPDATE);
1185 if (ir->eI == eiVV && bInitStep)
1187 /* if using velocity verlet with full time step Ekin,
1188 * take the first half step only to compute the
1189 * virial for the first step. From there,
1190 * revert back to the initial coordinates
1191 * so that the input is actually the initial step.
1193 snew(vbuf, state->natoms);
1194 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1196 else
1198 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1199 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1202 update_coords(step, ir, mdatoms, state, f, fcd,
1203 ekind, M, upd, etrtVELOCITY1,
1204 cr, constr);
1206 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1208 wallcycle_stop(wcycle, ewcUPDATE);
1209 constrain_velocities(step, nullptr, ir, mdatoms,
1210 state, fr->bMolPBC,
1211 &top->idef, shake_vir,
1212 cr, ms, nrnb, wcycle, constr,
1213 bCalcVir, do_log, do_ene);
1214 wallcycle_start(wcycle, ewcUPDATE);
1216 else if (graph)
1218 /* Need to unshift here if a do_force has been
1219 called in the previous step */
1220 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1222 /* if VV, compute the pressure and constraints */
1223 /* For VV2, we strictly only need this if using pressure
1224 * control, but we really would like to have accurate pressures
1225 * printed out.
1226 * Think about ways around this in the future?
1227 * For now, keep this choice in comments.
1229 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1230 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1231 bPres = TRUE;
1232 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1233 if (bCalcEner && ir->eI == eiVVAK)
1235 bSumEkinhOld = TRUE;
1237 /* for vv, the first half of the integration actually corresponds to the previous step.
1238 So we need information from the last step in the first half of the integration */
1239 if (bGStat || do_per_step(step-1, nstglobalcomm))
1241 wallcycle_stop(wcycle, ewcUPDATE);
1242 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1243 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1244 constr, &nullSignaller, state->box,
1245 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1246 (bGStat ? CGLO_GSTAT : 0)
1247 | CGLO_ENERGY
1248 | (bTemp ? CGLO_TEMPERATURE : 0)
1249 | (bPres ? CGLO_PRESSURE : 0)
1250 | (bPres ? CGLO_CONSTRAINT : 0)
1251 | (bStopCM ? CGLO_STOPCM : 0)
1252 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1253 | CGLO_SCALEEKIN
1255 /* explanation of above:
1256 a) We compute Ekin at the full time step
1257 if 1) we are using the AveVel Ekin, and it's not the
1258 initial step, or 2) if we are using AveEkin, but need the full
1259 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1260 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1261 EkinAveVel because it's needed for the pressure */
1262 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1263 top_global, top, state,
1264 &shouldCheckNumberOfBondedInteractions);
1265 wallcycle_start(wcycle, ewcUPDATE);
1267 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1268 if (!bInitStep)
1270 if (bTrotter)
1272 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1273 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1275 /* TODO This is only needed when we're about to write
1276 * a checkpoint, because we use it after the restart
1277 * (in a kludge?). But what should we be doing if
1278 * startingFromCheckpoint or bInitStep are true? */
1279 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1281 copy_mat(shake_vir, state->svir_prev);
1282 copy_mat(force_vir, state->fvir_prev);
1284 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1286 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1287 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1288 enerd->term[F_EKIN] = trace(ekind->ekin);
1291 else if (bExchanged)
1293 wallcycle_stop(wcycle, ewcUPDATE);
1294 /* We need the kinetic energy at minus the half step for determining
1295 * the full step kinetic energy and possibly for T-coupling.*/
1296 /* This may not be quite working correctly yet . . . . */
1297 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1298 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1299 constr, &nullSignaller, state->box,
1300 nullptr, &bSumEkinhOld,
1301 CGLO_GSTAT | CGLO_TEMPERATURE);
1302 wallcycle_start(wcycle, ewcUPDATE);
1305 /* if it's the initial step, we performed this first step just to get the constraint virial */
1306 if (ir->eI == eiVV && bInitStep)
1308 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1309 sfree(vbuf);
1311 wallcycle_stop(wcycle, ewcUPDATE);
1314 /* compute the conserved quantity */
1315 if (EI_VV(ir->eI))
1317 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1318 if (ir->eI == eiVV)
1320 last_ekin = enerd->term[F_EKIN];
1322 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1324 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1326 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1327 if (ir->efep != efepNO && !bRerunMD)
1329 sum_dhdl(enerd, state->lambda, ir->fepvals);
1333 /* ######## END FIRST UPDATE STEP ############## */
1334 /* ######## If doing VV, we now have v(dt) ###### */
1335 if (bDoExpanded)
1337 /* perform extended ensemble sampling in lambda - we don't
1338 actually move to the new state before outputting
1339 statistics, but if performing simulated tempering, we
1340 do update the velocities and the tau_t. */
1342 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1343 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1344 if (MASTER(cr))
1346 copy_df_history(state_global->dfhist, state->dfhist);
1350 /* Now we have the energies and forces corresponding to the
1351 * coordinates at time t. We must output all of this before
1352 * the update.
1354 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1355 ir, state, state_global, observablesHistory,
1356 top_global, fr,
1357 outf, mdebin, ekind, f,
1358 &nchkpt,
1359 bCPT, bRerunMD, bLastStep,
1360 mdrunOptions.writeConfout,
1361 bSumEkinhOld);
1362 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1363 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1365 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1366 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1368 copy_mat(state->svir_prev, shake_vir);
1369 copy_mat(state->fvir_prev, force_vir);
1372 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1374 /* Check whether everything is still allright */
1375 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1376 #if GMX_THREAD_MPI
1377 && MASTER(cr)
1378 #endif
1381 int nsteps_stop = -1;
1383 /* this just makes signals[].sig compatible with the hack
1384 of sending signals around by MPI_Reduce together with
1385 other floats */
1386 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1387 (mdrunOptions.reproducible &&
1388 gmx_get_stop_condition() == gmx_stop_cond_next))
1390 /* We need at least two global communication steps to pass
1391 * around the signal. We stop at a pair-list creation step
1392 * to allow for exact continuation, when possible.
1394 signals[eglsSTOPCOND].sig = 1;
1395 nsteps_stop = std::max(ir->nstlist, 2*nstSignalComm);
1397 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1399 /* Stop directly after the next global communication step.
1400 * This breaks exact continuation.
1402 signals[eglsSTOPCOND].sig = -1;
1403 nsteps_stop = nstSignalComm + 1;
1405 if (fplog)
1407 fprintf(fplog,
1408 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1409 gmx_get_signal_name(), nsteps_stop);
1410 fflush(fplog);
1412 fprintf(stderr,
1413 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1414 gmx_get_signal_name(), nsteps_stop);
1415 fflush(stderr);
1416 handled_stop_condition = (int)gmx_get_stop_condition();
1418 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1419 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1420 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1422 /* Signal to terminate the run */
1423 signals[eglsSTOPCOND].sig = 1;
1424 if (fplog)
1426 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1428 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1431 if (bResetCountersHalfMaxH && MASTER(cr) &&
1432 elapsed_time > max_hours*60.0*60.0*0.495)
1434 /* Set flag that will communicate the signal to all ranks in the simulation */
1435 signals[eglsRESETCOUNTERS].sig = 1;
1438 /* In parallel we only have to check for checkpointing in steps
1439 * where we do global communication,
1440 * otherwise the other nodes don't know.
1442 const real cpt_period = mdrunOptions.checkpointOptions.period;
1443 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1444 cpt_period >= 0 &&
1445 (cpt_period == 0 ||
1446 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1447 signals[eglsCHKPT].set == 0)
1449 signals[eglsCHKPT].sig = 1;
1452 /* ######### START SECOND UPDATE STEP ################# */
1454 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1455 in preprocessing */
1457 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1459 gmx_bool bIfRandomize;
1460 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1461 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1462 if (constr && bIfRandomize)
1464 constrain_velocities(step, nullptr, ir, mdatoms,
1465 state, fr->bMolPBC,
1466 &top->idef, tmp_vir,
1467 cr, ms, nrnb, wcycle, constr,
1468 bCalcVir, do_log, do_ene);
1471 /* Box is changed in update() when we do pressure coupling,
1472 * but we should still use the old box for energy corrections and when
1473 * writing it to the energy file, so it matches the trajectory files for
1474 * the same timestep above. Make a copy in a separate array.
1476 copy_mat(state->box, lastbox);
1478 dvdl_constr = 0;
1480 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1482 wallcycle_start(wcycle, ewcUPDATE);
1483 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1484 if (bTrotter)
1486 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1487 /* We can only do Berendsen coupling after we have summed
1488 * the kinetic energy or virial. Since the happens
1489 * in global_state after update, we should only do it at
1490 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1493 else
1495 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1496 update_pcouple_before_coordinates(fplog, step, ir, state,
1497 parrinellorahmanMu, M,
1498 bInitStep);
1501 if (EI_VV(ir->eI))
1503 /* velocity half-step update */
1504 update_coords(step, ir, mdatoms, state, f, fcd,
1505 ekind, M, upd, etrtVELOCITY2,
1506 cr, constr);
1509 /* Above, initialize just copies ekinh into ekin,
1510 * it doesn't copy position (for VV),
1511 * and entire integrator for MD.
1514 if (ir->eI == eiVVAK)
1516 /* We probably only need md->homenr, not state->natoms */
1517 if (state->natoms > cbuf_nalloc)
1519 cbuf_nalloc = state->natoms;
1520 srenew(cbuf, cbuf_nalloc);
1522 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1525 update_coords(step, ir, mdatoms, state, f, fcd,
1526 ekind, M, upd, etrtPOSITION, cr, constr);
1527 wallcycle_stop(wcycle, ewcUPDATE);
1529 constrain_coordinates(step, &dvdl_constr, ir, mdatoms, state,
1530 fr->bMolPBC,
1531 &top->idef, shake_vir,
1532 cr, ms, nrnb, wcycle, upd, constr,
1533 bCalcVir, do_log, do_ene);
1534 update_sd_second_half(step, &dvdl_constr, ir, mdatoms, state,
1535 fr->bMolPBC, &top->idef, cr, ms,
1536 nrnb, wcycle, upd, constr, do_log, do_ene);
1537 finish_update(ir, mdatoms,
1538 state, graph,
1539 nrnb, wcycle, upd, constr);
1541 if (ir->eI == eiVVAK)
1543 /* erase F_EKIN and F_TEMP here? */
1544 /* just compute the kinetic energy at the half step to perform a trotter step */
1545 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1546 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1547 constr, &nullSignaller, lastbox,
1548 nullptr, &bSumEkinhOld,
1549 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1551 wallcycle_start(wcycle, ewcUPDATE);
1552 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1553 /* now we know the scaling, we can compute the positions again again */
1554 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1556 update_coords(step, ir, mdatoms, state, f, fcd,
1557 ekind, M, upd, etrtPOSITION, cr, constr);
1558 wallcycle_stop(wcycle, ewcUPDATE);
1560 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1561 /* are the small terms in the shake_vir here due
1562 * to numerical errors, or are they important
1563 * physically? I'm thinking they are just errors, but not completely sure.
1564 * For now, will call without actually constraining, constr=NULL*/
1565 finish_update(ir, mdatoms,
1566 state, graph,
1567 nrnb, wcycle, upd, nullptr);
1569 if (EI_VV(ir->eI))
1571 /* this factor or 2 correction is necessary
1572 because half of the constraint force is removed
1573 in the vv step, so we have to double it. See
1574 the Redmine issue #1255. It is not yet clear
1575 if the factor of 2 is exact, or just a very
1576 good approximation, and this will be
1577 investigated. The next step is to see if this
1578 can be done adding a dhdl contribution from the
1579 rattle step, but this is somewhat more
1580 complicated with the current code. Will be
1581 investigated, hopefully for 4.6.3. However,
1582 this current solution is much better than
1583 having it completely wrong.
1585 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1587 else
1589 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1592 else if (graph)
1594 /* Need to unshift here */
1595 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1598 if (vsite != nullptr)
1600 wallcycle_start(wcycle, ewcVSITECONSTR);
1601 if (graph != nullptr)
1603 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1605 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1606 top->idef.iparams, top->idef.il,
1607 fr->ePBC, fr->bMolPBC, cr, state->box);
1609 if (graph != nullptr)
1611 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1613 wallcycle_stop(wcycle, ewcVSITECONSTR);
1616 /* ############## IF NOT VV, Calculate globals HERE ############ */
1617 /* With Leap-Frog we can skip compute_globals at
1618 * non-communication steps, but we need to calculate
1619 * the kinetic energy one step before communication.
1622 // Organize to do inter-simulation signalling on steps if
1623 // and when algorithms require it.
1624 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1626 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1628 // Since we're already communicating at this step, we
1629 // can propagate intra-simulation signals. Note that
1630 // check_nstglobalcomm has the responsibility for
1631 // choosing the value of nstglobalcomm that is one way
1632 // bGStat becomes true, so we can't get into a
1633 // situation where e.g. checkpointing can't be
1634 // signalled.
1635 bool doIntraSimSignal = true;
1636 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1638 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1639 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1640 constr, &signaller,
1641 lastbox,
1642 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1643 (bGStat ? CGLO_GSTAT : 0)
1644 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1645 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1646 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1647 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1648 | CGLO_CONSTRAINT
1649 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1651 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1652 top_global, top, state,
1653 &shouldCheckNumberOfBondedInteractions);
1657 /* ############# END CALC EKIN AND PRESSURE ################# */
1659 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1660 the virial that should probably be addressed eventually. state->veta has better properies,
1661 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1662 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1664 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1666 /* Sum up the foreign energy and dhdl terms for md and sd.
1667 Currently done every step so that dhdl is correct in the .edr */
1668 sum_dhdl(enerd, state->lambda, ir->fepvals);
1671 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1672 pres, force_vir, shake_vir,
1673 parrinellorahmanMu,
1674 state, nrnb, upd);
1676 /* ################# END UPDATE STEP 2 ################# */
1677 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1679 /* The coordinates (x) were unshifted in update */
1680 if (!bGStat)
1682 /* We will not sum ekinh_old,
1683 * so signal that we still have to do it.
1685 bSumEkinhOld = TRUE;
1688 if (bCalcEner)
1690 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1692 /* use the directly determined last velocity, not actually the averaged half steps */
1693 if (bTrotter && ir->eI == eiVV)
1695 enerd->term[F_EKIN] = last_ekin;
1697 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1699 if (integratorHasConservedEnergyQuantity(ir))
1701 if (EI_VV(ir->eI))
1703 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1705 else
1707 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1710 /* ######### END PREPARING EDR OUTPUT ########### */
1713 /* Output stuff */
1714 if (MASTER(cr))
1716 if (fplog && do_log && bDoExpanded)
1718 /* only needed if doing expanded ensemble */
1719 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1720 state_global->dfhist, state->fep_state, ir->nstlog, step);
1722 if (bCalcEner)
1724 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1725 t, mdatoms->tmass, enerd, state,
1726 ir->fepvals, ir->expandedvals, lastbox,
1727 shake_vir, force_vir, total_vir, pres,
1728 ekind, mu_tot, constr);
1730 else
1732 upd_mdebin_step(mdebin);
1735 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1736 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1738 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1739 step, t,
1740 eprNORMAL, mdebin, fcd, groups, &(ir->opts), awh.get());
1742 if (ir->bPull)
1744 pull_print_output(ir->pull_work, step, t);
1747 if (do_per_step(step, ir->nstlog))
1749 if (fflush(fplog) != 0)
1751 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1755 if (bDoExpanded)
1757 /* Have to do this part _after_ outputting the logfile and the edr file */
1758 /* Gets written into the state at the beginning of next loop*/
1759 state->fep_state = lamnew;
1761 /* Print the remaining wall clock time for the run */
1762 if (isMasterSimMasterRank(ms, cr) &&
1763 (do_verbose || gmx_got_usr_signal()) &&
1764 !bPMETunePrinting)
1766 if (shellfc)
1768 fprintf(stderr, "\n");
1770 print_time(stderr, walltime_accounting, step, ir, cr);
1773 /* Ion/water position swapping.
1774 * Not done in last step since trajectory writing happens before this call
1775 * in the MD loop and exchanges would be lost anyway. */
1776 bNeedRepartition = FALSE;
1777 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1778 do_per_step(step, ir->swap->nstswap))
1780 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1781 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1782 bRerunMD ? rerun_fr.box : state->box,
1783 MASTER(cr) && mdrunOptions.verbose,
1784 bRerunMD);
1786 if (bNeedRepartition && DOMAINDECOMP(cr))
1788 dd_collect_state(cr->dd, state, state_global);
1792 /* Replica exchange */
1793 bExchanged = FALSE;
1794 if (bDoReplEx)
1796 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1797 state_global, enerd,
1798 state, step, t);
1801 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1803 dd_partition_system(fplog, step, cr, TRUE, 1,
1804 state_global, top_global, ir,
1805 state, &f, mdAtoms, top, fr,
1806 vsite, constr,
1807 nrnb, wcycle, FALSE);
1808 shouldCheckNumberOfBondedInteractions = true;
1809 update_realloc(upd, state->natoms);
1812 bFirstStep = FALSE;
1813 bInitStep = FALSE;
1814 startingFromCheckpoint = false;
1816 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1817 /* With all integrators, except VV, we need to retain the pressure
1818 * at the current step for coupling at the next step.
1820 if ((state->flags & (1<<estPRES_PREV)) &&
1821 (bGStatEveryStep ||
1822 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1824 /* Store the pressure in t_state for pressure coupling
1825 * at the next MD step.
1827 copy_mat(pres, state->pres_prev);
1830 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1832 if ( (membed != nullptr) && (!bLastStep) )
1834 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1837 if (bRerunMD)
1839 if (MASTER(cr))
1841 /* read next frame from input trajectory */
1842 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1845 if (PAR(cr))
1847 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1851 cycles = wallcycle_stop(wcycle, ewcSTEP);
1852 if (DOMAINDECOMP(cr) && wcycle)
1854 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1857 if (!bRerunMD || !rerun_fr.bStep)
1859 /* increase the MD step number */
1860 step++;
1861 step_rel++;
1864 /* TODO make a counter-reset module */
1865 /* If it is time to reset counters, set a flag that remains
1866 true until counters actually get reset */
1867 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1868 signals[eglsRESETCOUNTERS].set != 0)
1870 if (pme_loadbal_is_active(pme_loadbal))
1872 /* Do not permit counter reset while PME load
1873 * balancing is active. The only purpose for resetting
1874 * counters is to measure reliable performance data,
1875 * and that can't be done before balancing
1876 * completes.
1878 * TODO consider fixing this by delaying the reset
1879 * until after load balancing completes,
1880 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1881 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1882 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1883 "resetting counters later in the run, e.g. with gmx "
1884 "mdrun -resetstep.", step);
1886 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1887 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1888 wcycle_set_reset_counters(wcycle, -1);
1889 if (!thisRankHasDuty(cr, DUTY_PME))
1891 /* Tell our PME node to reset its counters */
1892 gmx_pme_send_resetcounters(cr, step);
1894 /* Correct max_hours for the elapsed time */
1895 max_hours -= elapsed_time/(60.0*60.0);
1896 /* If mdrun -maxh -resethway was active, it can only trigger once */
1897 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1898 /* Reset can only happen once, so clear the triggering flag. */
1899 signals[eglsRESETCOUNTERS].set = 0;
1902 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1903 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1906 /* End of main MD loop */
1908 /* Closing TNG files can include compressing data. Therefore it is good to do that
1909 * before stopping the time measurements. */
1910 mdoutf_tng_close(outf);
1912 /* Stop measuring walltime */
1913 walltime_accounting_end(walltime_accounting);
1915 if (bRerunMD && MASTER(cr))
1917 close_trx(status);
1920 if (!thisRankHasDuty(cr, DUTY_PME))
1922 /* Tell the PME only node to finish */
1923 gmx_pme_send_finish(cr);
1926 if (MASTER(cr))
1928 if (ir->nstcalcenergy > 0 && !bRerunMD)
1930 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1931 eprAVER, mdebin, fcd, groups, &(ir->opts), awh.get());
1934 done_mdebin(mdebin);
1935 done_mdoutf(outf);
1937 if (bPMETune)
1939 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1942 done_shellfc(fplog, shellfc, step_rel);
1944 if (useReplicaExchange && MASTER(cr))
1946 print_replica_exchange_statistics(fplog, repl_ex);
1949 // Clean up swapcoords
1950 if (ir->eSwapCoords != eswapNO)
1952 finish_swapcoords(ir->swap);
1955 /* Do essential dynamics cleanup if needed. Close .edo file */
1956 done_ed(&ed);
1958 /* IMD cleanup, if bIMD is TRUE. */
1959 IMD_finalize(ir->bIMD, ir->imd);
1961 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1962 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1963 signals[eglsRESETCOUNTERS].set == 0 &&
1964 !bResetCountersHalfMaxH)
1966 walltime_accounting_set_valid_finish(walltime_accounting);
1969 destroy_enerdata(enerd);
1970 sfree(enerd);
1971 sfree(top);