Merge branch release-2018
[gromacs.git] / src / programs / mdrun / md.cpp
bloba8319e4bf8c89905242648c97b6766d92a4ec039
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 #include "gmxpre.h"
39 #include "md.h"
41 #include "config.h"
43 #include <stdio.h>
44 #include <stdlib.h>
46 #include <cmath>
48 #include <algorithm>
49 #include <memory>
51 #include "thread_mpi/threads.h"
53 #include "gromacs/awh/awh.h"
54 #include "gromacs/commandline/filenm.h"
55 #include "gromacs/compat/make_unique.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/essentialdynamics/edsam.h"
60 #include "gromacs/ewald/pme.h"
61 #include "gromacs/ewald/pme-load-balancing.h"
62 #include "gromacs/fileio/trxio.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/gmxlib/nrnb.h"
65 #include "gromacs/gpu_utils/gpu_utils.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/listed-forces/manage-threading.h"
68 #include "gromacs/math/functions.h"
69 #include "gromacs/math/utilities.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/math/vectypes.h"
72 #include "gromacs/mdlib/compute_io.h"
73 #include "gromacs/mdlib/constr.h"
74 #include "gromacs/mdlib/ebin.h"
75 #include "gromacs/mdlib/expanded.h"
76 #include "gromacs/mdlib/force.h"
77 #include "gromacs/mdlib/force_flags.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/md_support.h"
80 #include "gromacs/mdlib/mdatoms.h"
81 #include "gromacs/mdlib/mdebin.h"
82 #include "gromacs/mdlib/mdoutf.h"
83 #include "gromacs/mdlib/mdrun.h"
84 #include "gromacs/mdlib/mdsetup.h"
85 #include "gromacs/mdlib/nb_verlet.h"
86 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
87 #include "gromacs/mdlib/ns.h"
88 #include "gromacs/mdlib/shellfc.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/sim_util.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/tgroup.h"
93 #include "gromacs/mdlib/trajectory_writing.h"
94 #include "gromacs/mdlib/update.h"
95 #include "gromacs/mdlib/vcm.h"
96 #include "gromacs/mdlib/vsite.h"
97 #include "gromacs/mdtypes/awh-history.h"
98 #include "gromacs/mdtypes/awh-params.h"
99 #include "gromacs/mdtypes/commrec.h"
100 #include "gromacs/mdtypes/df_history.h"
101 #include "gromacs/mdtypes/energyhistory.h"
102 #include "gromacs/mdtypes/fcdata.h"
103 #include "gromacs/mdtypes/forcerec.h"
104 #include "gromacs/mdtypes/group.h"
105 #include "gromacs/mdtypes/inputrec.h"
106 #include "gromacs/mdtypes/interaction_const.h"
107 #include "gromacs/mdtypes/md_enums.h"
108 #include "gromacs/mdtypes/mdatom.h"
109 #include "gromacs/mdtypes/observableshistory.h"
110 #include "gromacs/mdtypes/state.h"
111 #include "gromacs/pbcutil/mshift.h"
112 #include "gromacs/pbcutil/pbc.h"
113 #include "gromacs/pulling/pull.h"
114 #include "gromacs/swap/swapcoords.h"
115 #include "gromacs/timing/wallcycle.h"
116 #include "gromacs/timing/walltime_accounting.h"
117 #include "gromacs/topology/atoms.h"
118 #include "gromacs/topology/idef.h"
119 #include "gromacs/topology/mtop_util.h"
120 #include "gromacs/topology/topology.h"
121 #include "gromacs/trajectory/trajectoryframe.h"
122 #include "gromacs/utility/basedefinitions.h"
123 #include "gromacs/utility/cstringutil.h"
124 #include "gromacs/utility/fatalerror.h"
125 #include "gromacs/utility/logger.h"
126 #include "gromacs/utility/real.h"
127 #include "gromacs/utility/smalloc.h"
129 #include "deform.h"
130 #include "membed.h"
131 #include "repl_ex.h"
133 #ifdef GMX_FAHCORE
134 #include "corewrap.h"
135 #endif
137 using gmx::SimulationSignaller;
139 /*! \brief Check whether bonded interactions are missing, if appropriate
141 * \param[in] fplog Log file pointer
142 * \param[in] cr Communication object
143 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
144 * \param[in] top_global Global topology for the error message
145 * \param[in] top_local Local topology for the error message
146 * \param[in] state Global state for the error message
147 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
149 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
150 * is always set to false after exit.
152 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
153 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
154 bool *shouldCheckNumberOfBondedInteractions)
156 if (*shouldCheckNumberOfBondedInteractions)
158 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
160 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
162 *shouldCheckNumberOfBondedInteractions = false;
166 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
167 gmx_int64_t step,
168 gmx_int64_t *step_rel, t_inputrec *ir,
169 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
170 gmx_walltime_accounting_t walltime_accounting,
171 struct nonbonded_verlet_t *nbv,
172 struct gmx_pme_t *pme)
174 char sbuf[STEPSTRSIZE];
176 /* Reset all the counters related to performance over the run */
177 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
178 "step %s: resetting all time and cycle counters",
179 gmx_step_str(step, sbuf));
181 if (use_GPU(nbv))
183 nbnxn_gpu_reset_timings(nbv);
186 if (pme_gpu_task_enabled(pme))
188 pme_gpu_reset_timings(pme);
191 if (use_GPU(nbv) || pme_gpu_task_enabled(pme))
193 resetGpuProfiler();
196 wallcycle_stop(wcycle, ewcRUN);
197 wallcycle_reset_all(wcycle);
198 if (DOMAINDECOMP(cr))
200 reset_dd_statistics_counters(cr->dd);
202 init_nrnb(nrnb);
203 ir->init_step += *step_rel;
204 ir->nsteps -= *step_rel;
205 *step_rel = 0;
206 wallcycle_start(wcycle, ewcRUN);
207 walltime_accounting_start(walltime_accounting);
208 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
211 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
213 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
214 * \param[in,out] globalState The global state container
215 * \param[in] constructVsites When true, vsite coordinates are constructed
216 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
217 * \param[in] idef Topology parameters, used for constructing vsites
218 * \param[in] timeStep Time step, used for constructing vsites
219 * \param[in] forceRec Force record, used for constructing vsites
220 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
221 * \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
223 static void prepareRerunState(const t_trxframe &rerunFrame,
224 t_state *globalState,
225 bool constructVsites,
226 const gmx_vsite_t *vsite,
227 const t_idef &idef,
228 double timeStep,
229 const t_forcerec &forceRec,
230 t_graph *graph,
231 gmx_bool *warnWhenNoV)
233 for (int i = 0; i < globalState->natoms; i++)
235 copy_rvec(rerunFrame.x[i], globalState->x[i]);
237 if (rerunFrame.bV)
239 for (int i = 0; i < globalState->natoms; i++)
241 copy_rvec(rerunFrame.v[i], globalState->v[i]);
244 else
246 for (int i = 0; i < globalState->natoms; i++)
248 clear_rvec(globalState->v[i]);
250 if (*warnWhenNoV)
252 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
253 " Ekin, temperature and pressure are incorrect,\n"
254 " the virial will be incorrect when constraints are present.\n"
255 "\n");
256 *warnWhenNoV = FALSE;
259 copy_mat(rerunFrame.box, globalState->box);
261 if (constructVsites)
263 GMX_ASSERT(vsite, "Need valid vsite for constructing vsites");
265 if (graph)
267 /* Following is necessary because the graph may get out of sync
268 * with the coordinates if we only have every N'th coordinate set
270 mk_mshift(nullptr, graph, forceRec.ePBC, globalState->box, as_rvec_array(globalState->x.data()));
271 shift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
273 construct_vsites(vsite, as_rvec_array(globalState->x.data()), timeStep, as_rvec_array(globalState->v.data()),
274 idef.iparams, idef.il,
275 forceRec.ePBC, forceRec.bMolPBC, nullptr, globalState->box);
276 if (graph)
278 unshift_self(graph, globalState->box, as_rvec_array(globalState->x.data()));
283 /*! \libinternal
284 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
285 const gmx_multisim_t *ms,
286 const gmx::MDLogger &mdlog,
287 int nfile, const t_filenm fnm[],
288 const gmx_output_env_t *oenv,
289 const MdrunOptions &mdrunOptions,
290 gmx_vsite_t *vsite, gmx_constr_t constr,
291 gmx::IMDOutputProvider *outputProvider,
292 t_inputrec *inputrec,
293 gmx_mtop_t *top_global, t_fcdata *fcd,
294 t_state *state_global,
295 t_mdatoms *mdatoms,
296 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
297 t_forcerec *fr,
298 const ReplicaExchangeParameters &replExParams,
299 gmx_membed_t *membed,
300 gmx_walltime_accounting_t walltime_accounting)
302 double gmx::do_md(FILE *fplog, t_commrec *cr,
303 const gmx_multisim_t *ms,
304 const gmx::MDLogger &mdlog,
305 int nfile, const t_filenm fnm[],
306 const gmx_output_env_t *oenv,
307 const MdrunOptions &mdrunOptions,
308 gmx_vsite_t *vsite, gmx_constr_t constr,
309 gmx::IMDOutputProvider *outputProvider,
310 t_inputrec *ir,
311 gmx_mtop_t *top_global,
312 t_fcdata *fcd,
313 t_state *state_global,
314 ObservablesHistory *observablesHistory,
315 gmx::MDAtoms *mdAtoms,
316 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
317 t_forcerec *fr,
318 const ReplicaExchangeParameters &replExParams,
319 gmx_membed_t *membed,
320 gmx_walltime_accounting_t walltime_accounting)
322 gmx_mdoutf_t outf = nullptr;
323 gmx_int64_t step, step_rel;
324 double elapsed_time;
325 double t, t0, lam0[efptNR];
326 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
327 gmx_bool bNS, bNStList, bSimAnn, bStopCM,
328 bFirstStep, bInitStep, bLastStep = FALSE;
329 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
330 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
331 bForceUpdate = FALSE, bCPT;
332 gmx_bool bMasterState;
333 int force_flags, cglo_flags;
334 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
335 int i, m;
336 t_trxstatus *status;
337 rvec mu_tot;
338 t_vcm *vcm;
339 matrix parrinellorahmanMu, M;
340 t_trxframe rerun_fr;
341 gmx_repl_ex_t repl_ex = nullptr;
342 int nchkpt = 1;
343 gmx_localtop_t *top;
344 t_mdebin *mdebin = nullptr;
345 gmx_enerdata_t *enerd;
346 PaddedRVecVector f {};
347 gmx_global_stat_t gstat;
348 gmx_update_t *upd = nullptr;
349 t_graph *graph = nullptr;
350 gmx_groups_t *groups;
351 gmx_ekindata_t *ekind;
352 gmx_shellfc_t *shellfc;
353 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
354 gmx_bool bResetCountersHalfMaxH = FALSE;
355 gmx_bool bTemp, bPres, bTrotter;
356 real dvdl_constr;
357 rvec *cbuf = nullptr;
358 int cbuf_nalloc = 0;
359 matrix lastbox;
360 int lamnew = 0;
361 /* for FEP */
362 int nstfep = 0;
363 double cycles;
364 real saved_conserved_quantity = 0;
365 real last_ekin = 0;
366 t_extmass MassQ;
367 int **trotter_seq;
368 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
369 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
372 /* PME load balancing data for GPU kernels */
373 pme_load_balancing_t *pme_loadbal = nullptr;
374 gmx_bool bPMETune = FALSE;
375 gmx_bool bPMETunePrinting = FALSE;
377 /* Interactive MD */
378 gmx_bool bIMDstep = FALSE;
380 #ifdef GMX_FAHCORE
381 /* Temporary addition for FAHCORE checkpointing */
382 int chkpt_ret;
383 #endif
384 /* Domain decomposition could incorrectly miss a bonded
385 interaction, but checking for that requires a global
386 communication stage, which does not otherwise happen in DD
387 code. So we do that alongside the first global energy reduction
388 after a new DD is made. These variables handle whether the
389 check happens, and the result it returns. */
390 bool shouldCheckNumberOfBondedInteractions = false;
391 int totalNumberOfBondedInteractions = -1;
393 SimulationSignals signals;
394 // Most global communnication stages don't propagate mdrun
395 // signals, and will use this object to achieve that.
396 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
398 if (!mdrunOptions.writeConfout)
400 // This is on by default, and the main known use case for
401 // turning it off is for convenience in benchmarking, which is
402 // something that should not show up in the general user
403 // interface.
404 GMX_LOG(mdlog.info).asParagraph().
405 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
407 if (mdrunOptions.timingOptions.resetHalfway)
409 GMX_LOG(mdlog.info).asParagraph().
410 appendText("The -resethway functionality is deprecated, and may be removed in a future version.");
411 if (ir->nsteps > 0)
413 /* Signal to reset the counters half the simulation steps. */
414 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
416 /* Signal to reset the counters halfway the simulation time. */
417 bResetCountersHalfMaxH = (mdrunOptions.maximumHoursToRun > 0);
420 /* md-vv uses averaged full step velocities for T-control
421 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
422 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
423 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
425 const gmx_bool bRerunMD = mdrunOptions.rerun;
426 int nstglobalcomm = mdrunOptions.globalCommunicationInterval;
428 if (bRerunMD)
430 /* Since we don't know if the frames read are related in any way,
431 * rebuild the neighborlist at every step.
433 ir->nstlist = 1;
434 ir->nstcalcenergy = 1;
435 nstglobalcomm = 1;
438 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
439 bGStatEveryStep = (nstglobalcomm == 1);
441 if (bRerunMD)
443 ir->nstxout_compressed = 0;
445 groups = &top_global->groups;
447 gmx_edsam *ed = nullptr;
448 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
450 /* Initialize essential dynamics sampling */
451 ed = init_edsam(opt2fn_null("-ei", nfile, fnm), opt2fn("-eo", nfile, fnm),
452 top_global,
453 ir, cr, constr,
454 state_global, observablesHistory,
455 oenv, mdrunOptions.continuationOptions.appendFiles);
458 if (ir->eSwapCoords != eswapNO)
460 /* Initialize ion swapping code */
461 init_swapcoords(fplog, ir, opt2fn_master("-swap", nfile, fnm, cr),
462 top_global,
463 state_global, observablesHistory,
464 cr, oenv, mdrunOptions);
467 /* Initial values */
468 init_md(fplog, cr, outputProvider, ir, oenv, mdrunOptions,
469 &t, &t0, state_global, lam0,
470 nrnb, top_global, &upd,
471 nfile, fnm, &outf, &mdebin,
472 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, wcycle);
474 clear_mat(total_vir);
475 clear_mat(pres);
476 /* Energy terms and groups */
477 snew(enerd, 1);
478 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
479 enerd);
481 /* Kinetic energy data */
482 snew(ekind, 1);
483 init_ekindata(fplog, top_global, &(ir->opts), ekind);
484 /* Copy the cos acceleration to the groups struct */
485 ekind->cosacc.cos_accel = ir->cos_accel;
487 gstat = global_stat_init(ir);
489 /* Check for polarizable models and flexible constraints */
490 shellfc = init_shell_flexcon(fplog,
491 top_global, n_flexible_constraints(constr),
492 ir->nstcalcenergy, DOMAINDECOMP(cr));
494 if (shellfc && ir->nstcalcenergy != 1)
496 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);
498 if (shellfc && DOMAINDECOMP(cr))
500 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
502 if (shellfc && ir->bDoAwh)
504 gmx_fatal(FARGS, "AWH biasing does not support shell particles.");
507 if (inputrecDeform(ir))
509 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
510 set_deform_reference_box(upd,
511 deform_init_init_step_tpx,
512 deform_init_box_tpx);
513 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
517 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
518 if ((io > 2000) && MASTER(cr))
520 fprintf(stderr,
521 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
522 io);
526 // Local state only becomes valid now.
527 std::unique_ptr<t_state> stateInstance;
528 t_state * state;
530 if (DOMAINDECOMP(cr))
532 top = dd_init_local_top(top_global);
534 stateInstance = gmx::compat::make_unique<t_state>();
535 state = stateInstance.get();
536 dd_init_local_state(cr->dd, state_global, state);
538 else
540 state_change_natoms(state_global, state_global->natoms);
541 /* We need to allocate one element extra, since we might use
542 * (unaligned) 4-wide SIMD loads to access rvec entries.
544 f.resize(gmx::paddedRVecVectorSize(state_global->natoms));
545 /* Copy the pointer to the global state */
546 state = state_global;
548 snew(top, 1);
549 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
550 &graph, mdAtoms, vsite, shellfc);
552 update_realloc(upd, state->natoms);
555 /* Set up interactive MD (IMD) */
556 init_IMD(ir, cr, ms, top_global, fplog, ir->nstcalcenergy,
557 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
558 nfile, fnm, oenv, mdrunOptions);
560 if (DOMAINDECOMP(cr))
562 /* Distribute the charge groups over the nodes from the master node */
563 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
564 state_global, top_global, ir,
565 state, &f, mdAtoms, top, fr,
566 vsite, constr,
567 nrnb, nullptr, FALSE);
568 shouldCheckNumberOfBondedInteractions = true;
569 update_realloc(upd, state->natoms);
572 auto mdatoms = mdAtoms->mdatoms();
574 // NOTE: The global state is no longer used at this point.
575 // But state_global is still used as temporary storage space for writing
576 // the global state to file and potentially for replica exchange.
577 // (Global topology should persist.)
579 update_mdatoms(mdatoms, state->lambda[efptMASS]);
581 const ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
582 bool startingFromCheckpoint = continuationOptions.startedFromCheckpoint;
584 if (ir->bExpanded)
586 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
589 if (MASTER(cr))
591 if (startingFromCheckpoint)
593 /* Update mdebin with energy history if appending to output files */
594 if (continuationOptions.appendFiles)
596 restore_energyhistory_from_state(mdebin, observablesHistory->energyHistory.get());
598 else if (observablesHistory->energyHistory.get() != nullptr)
600 /* We might have read an energy history from checkpoint.
601 * As we are not appending, we want to restart the statistics.
602 * Free the allocated memory and reset the counts.
604 observablesHistory->energyHistory = {};
607 if (observablesHistory->energyHistory.get() == nullptr)
609 observablesHistory->energyHistory = std::unique_ptr<energyhistory_t>(new energyhistory_t {});
611 /* Set the initial energy history in state by updating once */
612 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
615 /* Initialize constraints */
616 if (constr && !DOMAINDECOMP(cr))
618 set_constraints(constr, top, ir, mdatoms, cr);
621 /* Initialize AWH and restore state from history in checkpoint if needed. */
622 if (ir->bDoAwh)
624 ir->awh = new gmx::Awh(fplog, *ir, cr, ms, *ir->awhParams, opt2fn("-awh", nfile, fnm), ir->pull_work);
626 if (startingFromCheckpoint)
628 /* Restore the AWH history read from checkpoint */
629 ir->awh->restoreStateFromHistory(MASTER(cr) ? state_global->awhHistory.get() : nullptr);
631 else if (MASTER(cr))
633 /* Initialize the AWH history here */
634 state_global->awhHistory = ir->awh->initHistoryFromState();
638 const bool useReplicaExchange = (replExParams.exchangeInterval > 0);
639 if (useReplicaExchange && MASTER(cr))
641 repl_ex = init_replica_exchange(fplog, ms, top_global->natoms, ir,
642 replExParams);
644 /* PME tuning is only supported in the Verlet scheme, with PME for
645 * Coulomb. It is not supported with only LJ PME, or for
646 * reruns. */
647 bPMETune = (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !bRerunMD &&
648 !mdrunOptions.reproducible && ir->cutoff_scheme != ecutsGROUP);
649 if (bPMETune)
651 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
652 fr->ic, fr->nbv->listParams.get(), fr->pmedata, use_GPU(fr->nbv),
653 &bPMETunePrinting);
656 if (!ir->bContinuation && !bRerunMD)
658 if (state->flags & (1 << estV))
660 /* Set the velocities of vsites, shells and frozen atoms to zero */
661 for (i = 0; i < mdatoms->homenr; i++)
663 if (mdatoms->ptype[i] == eptVSite ||
664 mdatoms->ptype[i] == eptShell)
666 clear_rvec(state->v[i]);
668 else if (mdatoms->cFREEZE)
670 for (m = 0; m < DIM; m++)
672 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
674 state->v[i][m] = 0;
681 if (constr)
683 /* Constrain the initial coordinates and velocities */
684 do_constrain_first(fplog, constr, ir, mdatoms, state,
685 cr, ms, nrnb, fr, top);
687 if (vsite)
689 /* Construct the virtual sites for the initial configuration */
690 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
691 top->idef.iparams, top->idef.il,
692 fr->ePBC, fr->bMolPBC, cr, state->box);
696 if (ir->efep != efepNO)
698 /* Set free energy calculation frequency as the greatest common
699 * denominator of nstdhdl and repl_ex_nst. */
700 nstfep = ir->fepvals->nstdhdl;
701 if (ir->bExpanded)
703 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
705 if (useReplicaExchange)
707 nstfep = gmx_greatest_common_divisor(replExParams.exchangeInterval, nstfep);
711 /* Be REALLY careful about what flags you set here. You CANNOT assume
712 * this is the first step, since we might be restarting from a checkpoint,
713 * and in that case we should not do any modifications to the state.
715 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
717 if (continuationOptions.haveReadEkin)
719 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
722 cglo_flags = (CGLO_INITIALIZATION | CGLO_TEMPERATURE | CGLO_GSTAT
723 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
724 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
725 | (continuationOptions.haveReadEkin ? CGLO_READEKIN : 0));
727 bSumEkinhOld = FALSE;
728 /* To minimize communication, compute_globals computes the COM velocity
729 * and the kinetic energy for the velocities without COM motion removed.
730 * Thus to get the kinetic energy without the COM contribution, we need
731 * to call compute_globals twice.
733 for (int cgloIteration = 0; cgloIteration < (bStopCM ? 2 : 1); cgloIteration++)
735 int cglo_flags_iteration = cglo_flags;
736 if (bStopCM && cgloIteration == 0)
738 cglo_flags_iteration |= CGLO_STOPCM;
739 cglo_flags_iteration &= ~CGLO_TEMPERATURE;
741 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
742 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
743 constr, &nullSignaller, state->box,
744 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags_iteration
745 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
747 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
748 top_global, top, state,
749 &shouldCheckNumberOfBondedInteractions);
750 if (ir->eI == eiVVAK)
752 /* a second call to get the half step temperature initialized as well */
753 /* we do the same call as above, but turn the pressure off -- internally to
754 compute_globals, this is recognized as a velocity verlet half-step
755 kinetic energy calculation. This minimized excess variables, but
756 perhaps loses some logic?*/
758 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
759 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
760 constr, &nullSignaller, state->box,
761 nullptr, &bSumEkinhOld,
762 cglo_flags & ~CGLO_PRESSURE);
765 /* Calculate the initial half step temperature, and save the ekinh_old */
766 if (!continuationOptions.startedFromCheckpoint)
768 for (i = 0; (i < ir->opts.ngtc); i++)
770 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
774 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
775 temperature control */
776 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
778 if (MASTER(cr))
780 if (!ir->bContinuation)
782 if (constr && ir->eConstrAlg == econtLINCS)
784 fprintf(fplog,
785 "RMS relative constraint deviation after constraining: %.2e\n",
786 constr_rmsd(constr));
788 if (EI_STATE_VELOCITY(ir->eI))
790 real temp = enerd->term[F_TEMP];
791 if (ir->eI != eiVV)
793 /* Result of Ekin averaged over velocities of -half
794 * and +half step, while we only have -half step here.
796 temp *= 2;
798 fprintf(fplog, "Initial temperature: %g K\n", temp);
802 if (bRerunMD)
804 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
805 " input trajectory '%s'\n\n",
806 *(top_global->name), opt2fn("-rerun", nfile, fnm));
807 if (mdrunOptions.verbose)
809 fprintf(stderr, "Calculated time to finish depends on nsteps from "
810 "run input file,\nwhich may not correspond to the time "
811 "needed to process input trajectory.\n\n");
814 else
816 char tbuf[20];
817 fprintf(stderr, "starting mdrun '%s'\n",
818 *(top_global->name));
819 if (ir->nsteps >= 0)
821 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
823 else
825 sprintf(tbuf, "%s", "infinite");
827 if (ir->init_step > 0)
829 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
830 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
831 gmx_step_str(ir->init_step, sbuf2),
832 ir->init_step*ir->delta_t);
834 else
836 fprintf(stderr, "%s steps, %s ps.\n",
837 gmx_step_str(ir->nsteps, sbuf), tbuf);
840 fprintf(fplog, "\n");
843 walltime_accounting_start(walltime_accounting);
844 wallcycle_start(wcycle, ewcRUN);
845 print_start(fplog, cr, walltime_accounting, "mdrun");
847 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
848 #ifdef GMX_FAHCORE
849 chkpt_ret = fcCheckPointParallel( cr->nodeid,
850 NULL, 0);
851 if (chkpt_ret == 0)
853 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
855 #endif
857 /***********************************************************
859 * Loop over MD steps
861 ************************************************************/
863 /* if rerunMD then read coordinates and velocities from input trajectory */
864 if (bRerunMD)
866 if (getenv("GMX_FORCE_UPDATE"))
868 bForceUpdate = TRUE;
871 rerun_fr.natoms = 0;
872 if (MASTER(cr))
874 bLastStep = !read_first_frame(oenv, &status,
875 opt2fn("-rerun", nfile, fnm),
876 &rerun_fr, TRX_NEED_X | TRX_READ_V);
877 if (rerun_fr.natoms != top_global->natoms)
879 gmx_fatal(FARGS,
880 "Number of atoms in trajectory (%d) does not match the "
881 "run input file (%d)\n",
882 rerun_fr.natoms, top_global->natoms);
884 if (ir->ePBC != epbcNONE)
886 if (!rerun_fr.bBox)
888 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);
890 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
892 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
897 if (PAR(cr))
899 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
902 if (ir->ePBC != epbcNONE)
904 /* Set the shift vectors.
905 * Necessary here when have a static box different from the tpr box.
907 calc_shifts(rerun_fr.box, fr->shift_vec);
911 /* Loop over MD steps or if rerunMD to end of input trajectory,
912 * or, if max_hours>0, until max_hours is reached.
914 real max_hours = mdrunOptions.maximumHoursToRun;
915 bFirstStep = TRUE;
916 /* Skip the first Nose-Hoover integration when we get the state from tpx */
917 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
918 bSumEkinhOld = FALSE;
919 bExchanged = FALSE;
920 bNeedRepartition = FALSE;
922 bool simulationsShareState = false;
923 int nstSignalComm = nstglobalcomm;
925 // TODO This implementation of ensemble orientation restraints is nasty because
926 // a user can't just do multi-sim with single-sim orientation restraints.
927 bool usingEnsembleRestraints = (fcd->disres.nsystems > 1) || (ms && fcd->orires.nr);
928 bool awhUsesMultiSim = (ir->bDoAwh && ir->awhParams->shareBiasMultisim && ms);
930 // Replica exchange, ensemble restraints and AWH need all
931 // simulations to remain synchronized, so they need
932 // checkpoints and stop conditions to act on the same step, so
933 // the propagation of such signals must take place between
934 // simulations, not just within simulations.
935 // TODO: Make algorithm initializers set these flags.
936 simulationsShareState = useReplicaExchange || usingEnsembleRestraints || awhUsesMultiSim;
937 bool resetCountersIsLocal = true;
938 signals[eglsCHKPT] = SimulationSignal(!simulationsShareState);
939 signals[eglsSTOPCOND] = SimulationSignal(!simulationsShareState);
940 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
942 if (simulationsShareState)
944 // Inter-simulation signal communication does not need to happen
945 // often, so we use a minimum of 200 steps to reduce overhead.
946 const int c_minimumInterSimulationSignallingInterval = 200;
947 nstSignalComm = ((c_minimumInterSimulationSignallingInterval + nstglobalcomm - 1)/nstglobalcomm)*nstglobalcomm;
951 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion = (DOMAINDECOMP(cr) ? DdOpenBalanceRegionBeforeForceComputation::yes : DdOpenBalanceRegionBeforeForceComputation::no);
952 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion = (DOMAINDECOMP(cr) ? DdCloseBalanceRegionAfterForceComputation::yes : DdCloseBalanceRegionAfterForceComputation::no);
954 step = ir->init_step;
955 step_rel = 0;
957 // TODO extract this to new multi-simulation module
958 if (MASTER(cr) && isMultiSim(ms) && !useReplicaExchange)
960 if (!multisim_int_all_are_equal(ms, ir->nsteps))
962 GMX_LOG(mdlog.warning).appendText(
963 "Note: The number of steps is not consistent across multi simulations,\n"
964 "but we are proceeding anyway!");
966 if (!multisim_int_all_are_equal(ms, ir->init_step))
968 GMX_LOG(mdlog.warning).appendText(
969 "Note: The initial step is not consistent across multi simulations,\n"
970 "but we are proceeding anyway!");
974 /* and stop now if we should */
975 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
976 while (!bLastStep)
979 /* Determine if this is a neighbor search step */
980 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
982 if (bPMETune && bNStList)
984 /* PME grid + cut-off optimization with GPUs or PME nodes */
985 pme_loadbal_do(pme_loadbal, cr,
986 (mdrunOptions.verbose && MASTER(cr)) ? stderr : nullptr,
987 fplog, mdlog,
988 ir, fr, state,
989 wcycle,
990 step, step_rel,
991 &bPMETunePrinting);
994 wallcycle_start(wcycle, ewcSTEP);
996 if (bRerunMD)
998 if (rerun_fr.bStep)
1000 step = rerun_fr.step;
1001 step_rel = step - ir->init_step;
1003 if (rerun_fr.bTime)
1005 t = rerun_fr.time;
1007 else
1009 t = step;
1012 else
1014 bLastStep = (step_rel == ir->nsteps);
1015 t = t0 + step*ir->delta_t;
1018 // TODO Refactor this, so that nstfep does not need a default value of zero
1019 if (ir->efep != efepNO || ir->bSimTemp)
1021 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
1022 requiring different logic. */
1023 if (bRerunMD)
1025 if (MASTER(cr))
1027 setCurrentLambdasRerun(step, ir->fepvals, &rerun_fr, lam0, state_global);
1030 else
1032 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
1034 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
1035 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
1036 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
1037 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
1040 bDoReplEx = (useReplicaExchange && (step > 0) && !bLastStep &&
1041 do_per_step(step, replExParams.exchangeInterval));
1043 if (bSimAnn)
1045 update_annealing_target_temp(ir, t, upd);
1048 if (bRerunMD && MASTER(cr))
1050 const bool constructVsites = (vsite && mdrunOptions.rerunConstructVsites);
1051 if (constructVsites && DOMAINDECOMP(cr))
1053 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
1055 prepareRerunState(rerun_fr, state_global, constructVsites, vsite, top->idef, ir->delta_t, *fr, graph, &bRerunWarnNoV);
1058 /* Stop Center of Mass motion */
1059 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
1061 if (bRerunMD)
1063 /* for rerun MD always do Neighbour Searching */
1064 bNS = (bFirstStep || ir->nstlist != 0);
1066 else
1068 /* Determine whether or not to do Neighbour Searching */
1069 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
1072 /* < 0 means stop at next step, > 0 means stop at next NS step */
1073 if ( (signals[eglsSTOPCOND].set < 0) ||
1074 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
1076 bLastStep = TRUE;
1079 /* do_log triggers energy and virial calculation. Because this leads
1080 * to different code paths, forces can be different. Thus for exact
1081 * continuation we should avoid extra log output.
1082 * Note that the || bLastStep can result in non-exact continuation
1083 * beyond the last step. But we don't consider that to be an issue.
1085 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
1086 do_verbose = mdrunOptions.verbose &&
1087 (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep || bRerunMD);
1089 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
1091 if (bRerunMD)
1093 bMasterState = TRUE;
1095 else
1097 bMasterState = FALSE;
1098 /* Correct the new box if it is too skewed */
1099 if (inputrecDynamicBox(ir))
1101 if (correct_box(fplog, step, state->box, graph))
1103 bMasterState = TRUE;
1106 if (DOMAINDECOMP(cr) && bMasterState)
1108 dd_collect_state(cr->dd, state, state_global);
1112 if (DOMAINDECOMP(cr))
1114 /* Repartition the domain decomposition */
1115 dd_partition_system(fplog, step, cr,
1116 bMasterState, nstglobalcomm,
1117 state_global, top_global, ir,
1118 state, &f, mdAtoms, top, fr,
1119 vsite, constr,
1120 nrnb, wcycle,
1121 do_verbose && !bPMETunePrinting);
1122 shouldCheckNumberOfBondedInteractions = true;
1123 update_realloc(upd, state->natoms);
1127 if (MASTER(cr) && do_log)
1129 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1132 if (ir->efep != efepNO)
1134 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1137 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1140 /* We need the kinetic energy at minus the half step for determining
1141 * the full step kinetic energy and possibly for T-coupling.*/
1142 /* This may not be quite working correctly yet . . . . */
1143 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1144 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1145 constr, &nullSignaller, state->box,
1146 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1147 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1148 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1149 top_global, top, state,
1150 &shouldCheckNumberOfBondedInteractions);
1152 clear_mat(force_vir);
1154 /* We write a checkpoint at this MD step when:
1155 * either at an NS step when we signalled through gs,
1156 * or at the last step (but not when we do not want confout),
1157 * but never at the first step or with rerun.
1159 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1160 (bLastStep && mdrunOptions.writeConfout)) &&
1161 step > ir->init_step && !bRerunMD);
1162 if (bCPT)
1164 signals[eglsCHKPT].set = 0;
1167 /* Determine the energy and pressure:
1168 * at nstcalcenergy steps and at energy output steps (set below).
1170 if (EI_VV(ir->eI) && (!bInitStep))
1172 /* for vv, the first half of the integration actually corresponds
1173 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1174 but the virial needs to be calculated on both the current step and the 'next' step. Future
1175 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1177 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1178 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1179 bCalcVir = bCalcEnerStep ||
1180 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1182 else
1184 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1185 bCalcVir = bCalcEnerStep ||
1186 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1188 bCalcEner = bCalcEnerStep;
1190 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1192 if (do_ene || do_log || bDoReplEx)
1194 bCalcVir = TRUE;
1195 bCalcEner = TRUE;
1198 /* Do we need global communication ? */
1199 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1200 do_per_step(step, nstglobalcomm) ||
1201 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1203 force_flags = (GMX_FORCE_STATECHANGED |
1204 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1205 GMX_FORCE_ALLFORCES |
1206 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1207 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1208 (bDoFEP ? GMX_FORCE_DHDL : 0)
1211 if (shellfc)
1213 /* Now is the time to relax the shells */
1214 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, step,
1215 ir, bNS, force_flags, top,
1216 constr, enerd, fcd,
1217 state, &f, force_vir, mdatoms,
1218 nrnb, wcycle, graph, groups,
1219 shellfc, fr, t, mu_tot,
1220 vsite,
1221 ddOpenBalanceRegion, ddCloseBalanceRegion);
1223 else
1225 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1226 (or the AWH update will be performed twice for one step when continuing). It would be best to
1227 call this update function from do_md_trajectory_writing but that would occur after do_force.
1228 One would have to divide the update_awh function into one function applying the AWH force
1229 and one doing the AWH bias update. The update AWH bias function could then be called after
1230 do_md_trajectory_writing (then containing update_awh_history).
1231 The checkpointing will in the future probably moved to the start of the md loop which will
1232 rid of this issue. */
1233 if (ir->bDoAwh && bCPT && MASTER(cr))
1235 ir->awh->updateHistory(state_global->awhHistory.get());
1238 /* The coordinates (x) are shifted (to get whole molecules)
1239 * in do_force.
1240 * This is parallellized as well, and does communication too.
1241 * Check comments in sim_util.c
1243 do_force(fplog, cr, ms, ir, step, nrnb, wcycle, top, groups,
1244 state->box, state->x, &state->hist,
1245 f, force_vir, mdatoms, enerd, fcd,
1246 state->lambda, graph,
1247 fr, vsite, mu_tot, t, ed,
1248 (bNS ? GMX_FORCE_NS : 0) | force_flags,
1249 ddOpenBalanceRegion, ddCloseBalanceRegion);
1252 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1253 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1255 rvec *vbuf = nullptr;
1257 wallcycle_start(wcycle, ewcUPDATE);
1258 if (ir->eI == eiVV && bInitStep)
1260 /* if using velocity verlet with full time step Ekin,
1261 * take the first half step only to compute the
1262 * virial for the first step. From there,
1263 * revert back to the initial coordinates
1264 * so that the input is actually the initial step.
1266 snew(vbuf, state->natoms);
1267 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1269 else
1271 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1272 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1275 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1276 ekind, M, upd, etrtVELOCITY1,
1277 cr, constr);
1279 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1281 wallcycle_stop(wcycle, ewcUPDATE);
1282 update_constraints(fplog, step, nullptr, ir, mdatoms,
1283 state, fr->bMolPBC, graph, f,
1284 &top->idef, shake_vir,
1285 cr, ms, nrnb, wcycle, upd, constr,
1286 TRUE, bCalcVir);
1287 wallcycle_start(wcycle, ewcUPDATE);
1289 else if (graph)
1291 /* Need to unshift here if a do_force has been
1292 called in the previous step */
1293 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1295 /* if VV, compute the pressure and constraints */
1296 /* For VV2, we strictly only need this if using pressure
1297 * control, but we really would like to have accurate pressures
1298 * printed out.
1299 * Think about ways around this in the future?
1300 * For now, keep this choice in comments.
1302 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1303 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1304 bPres = TRUE;
1305 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1306 if (bCalcEner && ir->eI == eiVVAK)
1308 bSumEkinhOld = TRUE;
1310 /* for vv, the first half of the integration actually corresponds to the previous step.
1311 So we need information from the last step in the first half of the integration */
1312 if (bGStat || do_per_step(step-1, nstglobalcomm))
1314 wallcycle_stop(wcycle, ewcUPDATE);
1315 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1316 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1317 constr, &nullSignaller, state->box,
1318 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1319 (bGStat ? CGLO_GSTAT : 0)
1320 | CGLO_ENERGY
1321 | (bTemp ? CGLO_TEMPERATURE : 0)
1322 | (bPres ? CGLO_PRESSURE : 0)
1323 | (bPres ? CGLO_CONSTRAINT : 0)
1324 | (bStopCM ? CGLO_STOPCM : 0)
1325 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1326 | CGLO_SCALEEKIN
1328 /* explanation of above:
1329 a) We compute Ekin at the full time step
1330 if 1) we are using the AveVel Ekin, and it's not the
1331 initial step, or 2) if we are using AveEkin, but need the full
1332 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1333 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1334 EkinAveVel because it's needed for the pressure */
1335 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1336 top_global, top, state,
1337 &shouldCheckNumberOfBondedInteractions);
1338 wallcycle_start(wcycle, ewcUPDATE);
1340 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1341 if (!bInitStep)
1343 if (bTrotter)
1345 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1346 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1348 /* TODO This is only needed when we're about to write
1349 * a checkpoint, because we use it after the restart
1350 * (in a kludge?). But what should we be doing if
1351 * startingFromCheckpoint or bInitStep are true? */
1352 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir))
1354 copy_mat(shake_vir, state->svir_prev);
1355 copy_mat(force_vir, state->fvir_prev);
1357 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1359 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1360 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1361 enerd->term[F_EKIN] = trace(ekind->ekin);
1364 else if (bExchanged)
1366 wallcycle_stop(wcycle, ewcUPDATE);
1367 /* We need the kinetic energy at minus the half step for determining
1368 * the full step kinetic energy and possibly for T-coupling.*/
1369 /* This may not be quite working correctly yet . . . . */
1370 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1371 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1372 constr, &nullSignaller, state->box,
1373 nullptr, &bSumEkinhOld,
1374 CGLO_GSTAT | CGLO_TEMPERATURE);
1375 wallcycle_start(wcycle, ewcUPDATE);
1378 /* if it's the initial step, we performed this first step just to get the constraint virial */
1379 if (ir->eI == eiVV && bInitStep)
1381 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1382 sfree(vbuf);
1384 wallcycle_stop(wcycle, ewcUPDATE);
1387 /* compute the conserved quantity */
1388 if (EI_VV(ir->eI))
1390 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1391 if (ir->eI == eiVV)
1393 last_ekin = enerd->term[F_EKIN];
1395 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1397 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1399 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1400 if (ir->efep != efepNO && !bRerunMD)
1402 sum_dhdl(enerd, state->lambda, ir->fepvals);
1406 /* ######## END FIRST UPDATE STEP ############## */
1407 /* ######## If doing VV, we now have v(dt) ###### */
1408 if (bDoExpanded)
1410 /* perform extended ensemble sampling in lambda - we don't
1411 actually move to the new state before outputting
1412 statistics, but if performing simulated tempering, we
1413 do update the velocities and the tau_t. */
1415 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1416 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1417 if (MASTER(cr))
1419 copy_df_history(state_global->dfhist, state->dfhist);
1423 /* Now we have the energies and forces corresponding to the
1424 * coordinates at time t. We must output all of this before
1425 * the update.
1427 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1428 ir, state, state_global, observablesHistory,
1429 top_global, fr,
1430 outf, mdebin, ekind, f,
1431 &nchkpt,
1432 bCPT, bRerunMD, bLastStep,
1433 mdrunOptions.writeConfout,
1434 bSumEkinhOld);
1435 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1436 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1438 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1439 if (startingFromCheckpoint && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir)))
1441 copy_mat(state->svir_prev, shake_vir);
1442 copy_mat(state->fvir_prev, force_vir);
1445 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1447 /* Check whether everything is still allright */
1448 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1449 #if GMX_THREAD_MPI
1450 && MASTER(cr)
1451 #endif
1454 int nsteps_stop = -1;
1456 /* this just makes signals[].sig compatible with the hack
1457 of sending signals around by MPI_Reduce together with
1458 other floats */
1459 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns) ||
1460 (mdrunOptions.reproducible &&
1461 gmx_get_stop_condition() == gmx_stop_cond_next))
1463 /* We need at least two global communication steps to pass
1464 * around the signal. We stop at a pair-list creation step
1465 * to allow for exact continuation, when possible.
1467 signals[eglsSTOPCOND].sig = 1;
1468 nsteps_stop = std::max(ir->nstlist, 2*nstSignalComm);
1470 else if (gmx_get_stop_condition() == gmx_stop_cond_next)
1472 /* Stop directly after the next global communication step.
1473 * This breaks exact continuation.
1475 signals[eglsSTOPCOND].sig = -1;
1476 nsteps_stop = nstSignalComm + 1;
1478 if (fplog)
1480 fprintf(fplog,
1481 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1482 gmx_get_signal_name(), nsteps_stop);
1483 fflush(fplog);
1485 fprintf(stderr,
1486 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1487 gmx_get_signal_name(), nsteps_stop);
1488 fflush(stderr);
1489 handled_stop_condition = (int)gmx_get_stop_condition();
1491 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1492 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1493 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1495 /* Signal to terminate the run */
1496 signals[eglsSTOPCOND].sig = 1;
1497 if (fplog)
1499 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1501 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1504 if (bResetCountersHalfMaxH && MASTER(cr) &&
1505 elapsed_time > max_hours*60.0*60.0*0.495)
1507 /* Set flag that will communicate the signal to all ranks in the simulation */
1508 signals[eglsRESETCOUNTERS].sig = 1;
1511 /* In parallel we only have to check for checkpointing in steps
1512 * where we do global communication,
1513 * otherwise the other nodes don't know.
1515 const real cpt_period = mdrunOptions.checkpointOptions.period;
1516 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1517 cpt_period >= 0 &&
1518 (cpt_period == 0 ||
1519 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1520 signals[eglsCHKPT].set == 0)
1522 signals[eglsCHKPT].sig = 1;
1525 /* ######### START SECOND UPDATE STEP ################# */
1527 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1528 in preprocessing */
1530 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1532 gmx_bool bIfRandomize;
1533 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1534 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1535 if (constr && bIfRandomize)
1537 update_constraints(fplog, step, nullptr, ir, mdatoms,
1538 state, fr->bMolPBC, graph, f,
1539 &top->idef, tmp_vir,
1540 cr, ms, nrnb, wcycle, upd, constr,
1541 TRUE, bCalcVir);
1544 /* Box is changed in update() when we do pressure coupling,
1545 * but we should still use the old box for energy corrections and when
1546 * writing it to the energy file, so it matches the trajectory files for
1547 * the same timestep above. Make a copy in a separate array.
1549 copy_mat(state->box, lastbox);
1551 dvdl_constr = 0;
1553 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1555 wallcycle_start(wcycle, ewcUPDATE);
1556 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1557 if (bTrotter)
1559 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1560 /* We can only do Berendsen coupling after we have summed
1561 * the kinetic energy or virial. Since the happens
1562 * in global_state after update, we should only do it at
1563 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1566 else
1568 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1569 update_pcouple_before_coordinates(fplog, step, ir, state,
1570 parrinellorahmanMu, M,
1571 bInitStep);
1574 if (EI_VV(ir->eI))
1576 /* velocity half-step update */
1577 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1578 ekind, M, upd, etrtVELOCITY2,
1579 cr, constr);
1582 /* Above, initialize just copies ekinh into ekin,
1583 * it doesn't copy position (for VV),
1584 * and entire integrator for MD.
1587 if (ir->eI == eiVVAK)
1589 /* We probably only need md->homenr, not state->natoms */
1590 if (state->natoms > cbuf_nalloc)
1592 cbuf_nalloc = state->natoms;
1593 srenew(cbuf, cbuf_nalloc);
1595 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1598 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1599 ekind, M, upd, etrtPOSITION, cr, constr);
1600 wallcycle_stop(wcycle, ewcUPDATE);
1602 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1603 fr->bMolPBC, graph, f,
1604 &top->idef, shake_vir,
1605 cr, ms, nrnb, wcycle, upd, constr,
1606 FALSE, bCalcVir);
1608 if (ir->eI == eiVVAK)
1610 /* erase F_EKIN and F_TEMP here? */
1611 /* just compute the kinetic energy at the half step to perform a trotter step */
1612 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1613 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1614 constr, &nullSignaller, lastbox,
1615 nullptr, &bSumEkinhOld,
1616 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1618 wallcycle_start(wcycle, ewcUPDATE);
1619 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1620 /* now we know the scaling, we can compute the positions again again */
1621 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1623 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1624 ekind, M, upd, etrtPOSITION, cr, constr);
1625 wallcycle_stop(wcycle, ewcUPDATE);
1627 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1628 /* are the small terms in the shake_vir here due
1629 * to numerical errors, or are they important
1630 * physically? I'm thinking they are just errors, but not completely sure.
1631 * For now, will call without actually constraining, constr=NULL*/
1632 update_constraints(fplog, step, nullptr, ir, mdatoms,
1633 state, fr->bMolPBC, graph, f,
1634 &top->idef, tmp_vir,
1635 cr, ms, nrnb, wcycle, upd, nullptr,
1636 FALSE, bCalcVir);
1638 if (EI_VV(ir->eI))
1640 /* this factor or 2 correction is necessary
1641 because half of the constraint force is removed
1642 in the vv step, so we have to double it. See
1643 the Redmine issue #1255. It is not yet clear
1644 if the factor of 2 is exact, or just a very
1645 good approximation, and this will be
1646 investigated. The next step is to see if this
1647 can be done adding a dhdl contribution from the
1648 rattle step, but this is somewhat more
1649 complicated with the current code. Will be
1650 investigated, hopefully for 4.6.3. However,
1651 this current solution is much better than
1652 having it completely wrong.
1654 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1656 else
1658 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1661 else if (graph)
1663 /* Need to unshift here */
1664 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1667 if (vsite != nullptr)
1669 wallcycle_start(wcycle, ewcVSITECONSTR);
1670 if (graph != nullptr)
1672 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1674 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1675 top->idef.iparams, top->idef.il,
1676 fr->ePBC, fr->bMolPBC, cr, state->box);
1678 if (graph != nullptr)
1680 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1682 wallcycle_stop(wcycle, ewcVSITECONSTR);
1685 /* ############## IF NOT VV, Calculate globals HERE ############ */
1686 /* With Leap-Frog we can skip compute_globals at
1687 * non-communication steps, but we need to calculate
1688 * the kinetic energy one step before communication.
1691 // Organize to do inter-simulation signalling on steps if
1692 // and when algorithms require it.
1693 bool doInterSimSignal = (simulationsShareState && do_per_step(step, nstSignalComm));
1695 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1697 // Since we're already communicating at this step, we
1698 // can propagate intra-simulation signals. Note that
1699 // check_nstglobalcomm has the responsibility for
1700 // choosing the value of nstglobalcomm that is one way
1701 // bGStat becomes true, so we can't get into a
1702 // situation where e.g. checkpointing can't be
1703 // signalled.
1704 bool doIntraSimSignal = true;
1705 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
1707 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1708 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1709 constr, &signaller,
1710 lastbox,
1711 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1712 (bGStat ? CGLO_GSTAT : 0)
1713 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1714 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1715 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1716 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1717 | CGLO_CONSTRAINT
1718 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1720 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1721 top_global, top, state,
1722 &shouldCheckNumberOfBondedInteractions);
1726 /* ############# END CALC EKIN AND PRESSURE ################# */
1728 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1729 the virial that should probably be addressed eventually. state->veta has better properies,
1730 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1731 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1733 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1735 /* Sum up the foreign energy and dhdl terms for md and sd.
1736 Currently done every step so that dhdl is correct in the .edr */
1737 sum_dhdl(enerd, state->lambda, ir->fepvals);
1740 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1741 pres, force_vir, shake_vir,
1742 parrinellorahmanMu,
1743 state, nrnb, upd);
1745 /* ################# END UPDATE STEP 2 ################# */
1746 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1748 /* The coordinates (x) were unshifted in update */
1749 if (!bGStat)
1751 /* We will not sum ekinh_old,
1752 * so signal that we still have to do it.
1754 bSumEkinhOld = TRUE;
1757 if (bCalcEner)
1759 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1761 /* use the directly determined last velocity, not actually the averaged half steps */
1762 if (bTrotter && ir->eI == eiVV)
1764 enerd->term[F_EKIN] = last_ekin;
1766 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1768 if (integratorHasConservedEnergyQuantity(ir))
1770 if (EI_VV(ir->eI))
1772 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1774 else
1776 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1779 /* ######### END PREPARING EDR OUTPUT ########### */
1782 /* Output stuff */
1783 if (MASTER(cr))
1785 if (fplog && do_log && bDoExpanded)
1787 /* only needed if doing expanded ensemble */
1788 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1789 state_global->dfhist, state->fep_state, ir->nstlog, step);
1791 if (bCalcEner)
1793 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1794 t, mdatoms->tmass, enerd, state,
1795 ir->fepvals, ir->expandedvals, lastbox,
1796 shake_vir, force_vir, total_vir, pres,
1797 ekind, mu_tot, constr);
1799 else
1801 upd_mdebin_step(mdebin);
1804 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1805 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1807 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1808 step, t,
1809 eprNORMAL, mdebin, fcd, groups, &(ir->opts), ir->awh);
1811 if (ir->bPull)
1813 pull_print_output(ir->pull_work, step, t);
1816 if (do_per_step(step, ir->nstlog))
1818 if (fflush(fplog) != 0)
1820 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1824 if (bDoExpanded)
1826 /* Have to do this part _after_ outputting the logfile and the edr file */
1827 /* Gets written into the state at the beginning of next loop*/
1828 state->fep_state = lamnew;
1830 /* Print the remaining wall clock time for the run */
1831 if (isMasterSimMasterRank(ms, cr) &&
1832 (do_verbose || gmx_got_usr_signal()) &&
1833 !bPMETunePrinting)
1835 if (shellfc)
1837 fprintf(stderr, "\n");
1839 print_time(stderr, walltime_accounting, step, ir, cr);
1842 /* Ion/water position swapping.
1843 * Not done in last step since trajectory writing happens before this call
1844 * in the MD loop and exchanges would be lost anyway. */
1845 bNeedRepartition = FALSE;
1846 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1847 do_per_step(step, ir->swap->nstswap))
1849 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1850 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1851 bRerunMD ? rerun_fr.box : state->box,
1852 MASTER(cr) && mdrunOptions.verbose,
1853 bRerunMD);
1855 if (bNeedRepartition && DOMAINDECOMP(cr))
1857 dd_collect_state(cr->dd, state, state_global);
1861 /* Replica exchange */
1862 bExchanged = FALSE;
1863 if (bDoReplEx)
1865 bExchanged = replica_exchange(fplog, cr, ms, repl_ex,
1866 state_global, enerd,
1867 state, step, t);
1870 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1872 dd_partition_system(fplog, step, cr, TRUE, 1,
1873 state_global, top_global, ir,
1874 state, &f, mdAtoms, top, fr,
1875 vsite, constr,
1876 nrnb, wcycle, FALSE);
1877 shouldCheckNumberOfBondedInteractions = true;
1878 update_realloc(upd, state->natoms);
1881 bFirstStep = FALSE;
1882 bInitStep = FALSE;
1883 startingFromCheckpoint = false;
1885 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1886 /* With all integrators, except VV, we need to retain the pressure
1887 * at the current step for coupling at the next step.
1889 if ((state->flags & (1<<estPRES_PREV)) &&
1890 (bGStatEveryStep ||
1891 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1893 /* Store the pressure in t_state for pressure coupling
1894 * at the next MD step.
1896 copy_mat(pres, state->pres_prev);
1899 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1901 if ( (membed != nullptr) && (!bLastStep) )
1903 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1906 if (bRerunMD)
1908 if (MASTER(cr))
1910 /* read next frame from input trajectory */
1911 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1914 if (PAR(cr))
1916 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1920 cycles = wallcycle_stop(wcycle, ewcSTEP);
1921 if (DOMAINDECOMP(cr) && wcycle)
1923 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1926 if (!bRerunMD || !rerun_fr.bStep)
1928 /* increase the MD step number */
1929 step++;
1930 step_rel++;
1933 /* TODO make a counter-reset module */
1934 /* If it is time to reset counters, set a flag that remains
1935 true until counters actually get reset */
1936 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1937 signals[eglsRESETCOUNTERS].set != 0)
1939 if (pme_loadbal_is_active(pme_loadbal))
1941 /* Do not permit counter reset while PME load
1942 * balancing is active. The only purpose for resetting
1943 * counters is to measure reliable performance data,
1944 * and that can't be done before balancing
1945 * completes.
1947 * TODO consider fixing this by delaying the reset
1948 * until after load balancing completes,
1949 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1950 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1951 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1952 "resetting counters later in the run, e.g. with gmx "
1953 "mdrun -resetstep.", step);
1955 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1956 use_GPU(fr->nbv) ? fr->nbv : nullptr, fr->pmedata);
1957 wcycle_set_reset_counters(wcycle, -1);
1958 if (!thisRankHasDuty(cr, DUTY_PME))
1960 /* Tell our PME node to reset its counters */
1961 gmx_pme_send_resetcounters(cr, step);
1963 /* Correct max_hours for the elapsed time */
1964 max_hours -= elapsed_time/(60.0*60.0);
1965 /* If mdrun -maxh -resethway was active, it can only trigger once */
1966 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1967 /* Reset can only happen once, so clear the triggering flag. */
1968 signals[eglsRESETCOUNTERS].set = 0;
1971 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1972 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1975 /* End of main MD loop */
1977 /* Closing TNG files can include compressing data. Therefore it is good to do that
1978 * before stopping the time measurements. */
1979 mdoutf_tng_close(outf);
1981 /* Stop measuring walltime */
1982 walltime_accounting_end(walltime_accounting);
1984 if (bRerunMD && MASTER(cr))
1986 close_trx(status);
1989 if (!thisRankHasDuty(cr, DUTY_PME))
1991 /* Tell the PME only node to finish */
1992 gmx_pme_send_finish(cr);
1995 if (MASTER(cr))
1997 if (ir->nstcalcenergy > 0 && !bRerunMD)
1999 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
2000 eprAVER, mdebin, fcd, groups, &(ir->opts), ir->awh);
2003 done_ebin(mdebin->ebin);
2004 done_mdoutf(outf);
2006 if (bPMETune)
2008 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
2011 done_shellfc(fplog, shellfc, step_rel);
2013 if (useReplicaExchange && MASTER(cr))
2015 print_replica_exchange_statistics(fplog, repl_ex);
2018 if (ir->bDoAwh)
2020 delete ir->awh;
2023 // Clean up swapcoords
2024 if (ir->eSwapCoords != eswapNO)
2026 finish_swapcoords(ir->swap);
2029 /* Do essential dynamics cleanup if needed. Close .edo file */
2030 done_ed(&ed);
2032 /* IMD cleanup, if bIMD is TRUE. */
2033 IMD_finalize(ir->bIMD, ir->imd);
2035 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
2036 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
2037 signals[eglsRESETCOUNTERS].set == 0 &&
2038 !bResetCountersHalfMaxH)
2040 walltime_accounting_set_valid_finish(walltime_accounting);
2043 return 0;