Improved state_change_natoms and used it more
[gromacs.git] / src / programs / mdrun / md.cpp
blobc01e69fc97dbb58555c5b497dc881751a79370c5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "md.h"
41 #include "config.h"
43 #include <math.h>
44 #include <stdio.h>
45 #include <stdlib.h>
47 #include <algorithm>
48 #include <memory>
50 #include "thread_mpi/threads.h"
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_network.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme-load-balancing.h"
58 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/gpu_utils/gpu_utils.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed-forces/manage-threading.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/math/vectypes.h"
68 #include "gromacs/mdlib/compute_io.h"
69 #include "gromacs/mdlib/constr.h"
70 #include "gromacs/mdlib/ebin.h"
71 #include "gromacs/mdlib/force.h"
72 #include "gromacs/mdlib/force_flags.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/md_support.h"
75 #include "gromacs/mdlib/mdatoms.h"
76 #include "gromacs/mdlib/mdebin.h"
77 #include "gromacs/mdlib/mdoutf.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/nb_verlet.h"
81 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
82 #include "gromacs/mdlib/ns.h"
83 #include "gromacs/mdlib/shellfc.h"
84 #include "gromacs/mdlib/sighandler.h"
85 #include "gromacs/mdlib/sim_util.h"
86 #include "gromacs/mdlib/simulationsignal.h"
87 #include "gromacs/mdlib/tgroup.h"
88 #include "gromacs/mdlib/trajectory_writing.h"
89 #include "gromacs/mdlib/update.h"
90 #include "gromacs/mdlib/vcm.h"
91 #include "gromacs/mdlib/vsite.h"
92 #include "gromacs/mdtypes/commrec.h"
93 #include "gromacs/mdtypes/df_history.h"
94 #include "gromacs/mdtypes/energyhistory.h"
95 #include "gromacs/mdtypes/fcdata.h"
96 #include "gromacs/mdtypes/forcerec.h"
97 #include "gromacs/mdtypes/group.h"
98 #include "gromacs/mdtypes/inputrec.h"
99 #include "gromacs/mdtypes/interaction_const.h"
100 #include "gromacs/mdtypes/md_enums.h"
101 #include "gromacs/mdtypes/mdatom.h"
102 #include "gromacs/mdtypes/state.h"
103 #include "gromacs/pbcutil/mshift.h"
104 #include "gromacs/pbcutil/pbc.h"
105 #include "gromacs/pulling/pull.h"
106 #include "gromacs/swap/swapcoords.h"
107 #include "gromacs/timing/wallcycle.h"
108 #include "gromacs/timing/walltime_accounting.h"
109 #include "gromacs/topology/atoms.h"
110 #include "gromacs/topology/idef.h"
111 #include "gromacs/topology/mtop_util.h"
112 #include "gromacs/topology/topology.h"
113 #include "gromacs/trajectory/trajectoryframe.h"
114 #include "gromacs/utility/basedefinitions.h"
115 #include "gromacs/utility/cstringutil.h"
116 #include "gromacs/utility/fatalerror.h"
117 #include "gromacs/utility/logger.h"
118 #include "gromacs/utility/real.h"
119 #include "gromacs/utility/smalloc.h"
121 #include "deform.h"
122 #include "membed.h"
123 #include "repl_ex.h"
125 #ifdef GMX_FAHCORE
126 #include "corewrap.h"
127 #endif
129 using gmx::SimulationSignaller;
131 /*! \brief Check whether bonded interactions are missing, if appropriate
133 * \param[in] fplog Log file pointer
134 * \param[in] cr Communication object
135 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
136 * \param[in] top_global Global topology for the error message
137 * \param[in] top_local Local topology for the error message
138 * \param[in] state Global state for the error message
139 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
141 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
142 * is always set to false after exit.
144 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
145 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
146 bool *shouldCheckNumberOfBondedInteractions)
148 if (*shouldCheckNumberOfBondedInteractions)
150 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
152 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
154 *shouldCheckNumberOfBondedInteractions = false;
158 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
159 gmx_int64_t step,
160 gmx_int64_t *step_rel, t_inputrec *ir,
161 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
162 gmx_walltime_accounting_t walltime_accounting,
163 struct nonbonded_verlet_t *nbv)
165 char sbuf[STEPSTRSIZE];
167 /* Reset all the counters related to performance over the run */
168 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
169 "step %s: resetting all time and cycle counters",
170 gmx_step_str(step, sbuf));
172 if (use_GPU(nbv))
174 nbnxn_gpu_reset_timings(nbv);
175 resetGpuProfiler();
178 wallcycle_stop(wcycle, ewcRUN);
179 wallcycle_reset_all(wcycle);
180 if (DOMAINDECOMP(cr))
182 reset_dd_statistics_counters(cr->dd);
184 init_nrnb(nrnb);
185 ir->init_step += *step_rel;
186 ir->nsteps -= *step_rel;
187 *step_rel = 0;
188 wallcycle_start(wcycle, ewcRUN);
189 walltime_accounting_start(walltime_accounting);
190 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
193 /*! \libinternal
194 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
195 int nfile, const t_filenm fnm[],
196 const gmx_output_env_t *oenv, gmx_bool bVerbose,
197 int nstglobalcomm,
198 gmx_vsite_t *vsite, gmx_constr_t constr,
199 int stepout,
200 t_inputrec *inputrec,
201 gmx_mtop_t *top_global, t_fcdata *fcd,
202 t_state *state_global,
203 t_mdatoms *mdatoms,
204 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
205 gmx_edsam_t ed,
206 t_forcerec *fr,
207 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
208 real cpt_period, real max_hours,
209 int imdport,
210 unsigned long Flags,
211 gmx_walltime_accounting_t walltime_accounting)
213 double gmx::do_md(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
214 int nfile, const t_filenm fnm[],
215 const gmx_output_env_t *oenv, gmx_bool bVerbose,
216 int nstglobalcomm,
217 gmx_vsite_t *vsite, gmx_constr_t constr,
218 int stepout, t_inputrec *ir,
219 gmx_mtop_t *top_global,
220 t_fcdata *fcd,
221 t_state *state_global,
222 energyhistory_t *energyHistory,
223 t_mdatoms *mdatoms,
224 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
225 gmx_edsam_t ed, t_forcerec *fr,
226 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
227 gmx_membed_t *membed,
228 real cpt_period, real max_hours,
229 int imdport,
230 unsigned long Flags,
231 gmx_walltime_accounting_t walltime_accounting)
233 gmx_mdoutf_t outf = nullptr;
234 gmx_int64_t step, step_rel;
235 double elapsed_time;
236 double t, t0, lam0[efptNR];
237 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
238 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD,
239 bFirstStep, startingFromCheckpoint, bInitStep, bLastStep = FALSE,
240 bBornRadii, bUsingEnsembleRestraints;
241 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
242 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
243 bForceUpdate = FALSE, bCPT;
244 gmx_bool bMasterState;
245 int force_flags, cglo_flags;
246 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
247 int i, m;
248 t_trxstatus *status;
249 rvec mu_tot;
250 t_vcm *vcm;
251 matrix parrinellorahmanMu, M;
252 t_trxframe rerun_fr;
253 gmx_repl_ex_t repl_ex = nullptr;
254 int nchkpt = 1;
255 gmx_localtop_t *top;
256 t_mdebin *mdebin = nullptr;
257 gmx_enerdata_t *enerd;
258 PaddedRVecVector f {};
259 gmx_global_stat_t gstat;
260 gmx_update_t *upd = nullptr;
261 t_graph *graph = nullptr;
262 gmx_groups_t *groups;
263 gmx_ekindata_t *ekind;
264 gmx_shellfc_t *shellfc;
265 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
266 gmx_bool bResetCountersHalfMaxH = FALSE;
267 gmx_bool bTemp, bPres, bTrotter;
268 real dvdl_constr;
269 rvec *cbuf = nullptr;
270 int cbuf_nalloc = 0;
271 matrix lastbox;
272 int lamnew = 0;
273 /* for FEP */
274 int nstfep = 0;
275 double cycles;
276 real saved_conserved_quantity = 0;
277 real last_ekin = 0;
278 t_extmass MassQ;
279 int **trotter_seq;
280 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
281 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
284 /* PME load balancing data for GPU kernels */
285 pme_load_balancing_t *pme_loadbal = nullptr;
286 gmx_bool bPMETune = FALSE;
287 gmx_bool bPMETunePrinting = FALSE;
289 /* Interactive MD */
290 gmx_bool bIMDstep = FALSE;
292 #ifdef GMX_FAHCORE
293 /* Temporary addition for FAHCORE checkpointing */
294 int chkpt_ret;
295 #endif
296 /* Domain decomposition could incorrectly miss a bonded
297 interaction, but checking for that requires a global
298 communication stage, which does not otherwise happen in DD
299 code. So we do that alongside the first global energy reduction
300 after a new DD is made. These variables handle whether the
301 check happens, and the result it returns. */
302 bool shouldCheckNumberOfBondedInteractions = false;
303 int totalNumberOfBondedInteractions = -1;
305 SimulationSignals signals;
306 // Most global communnication stages don't propagate mdrun
307 // signals, and will use this object to achieve that.
308 SimulationSignaller nullSignaller(nullptr, nullptr, false, false);
310 /* Check for special mdrun options */
311 bRerunMD = (Flags & MD_RERUN);
312 if (Flags & MD_RESETCOUNTERSHALFWAY)
314 if (ir->nsteps > 0)
316 /* Signal to reset the counters half the simulation steps. */
317 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
319 /* Signal to reset the counters halfway the simulation time. */
320 bResetCountersHalfMaxH = (max_hours > 0);
323 /* md-vv uses averaged full step velocities for T-control
324 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
325 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
326 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
328 if (bRerunMD)
330 /* Since we don't know if the frames read are related in any way,
331 * rebuild the neighborlist at every step.
333 ir->nstlist = 1;
334 ir->nstcalcenergy = 1;
335 nstglobalcomm = 1;
338 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
339 bGStatEveryStep = (nstglobalcomm == 1);
341 if (bRerunMD)
343 ir->nstxout_compressed = 0;
345 groups = &top_global->groups;
347 if (ir->eSwapCoords != eswapNO)
349 /* Initialize ion swapping code */
350 init_swapcoords(fplog, bVerbose, ir, opt2fn_master("-swap", nfile, fnm, cr),
351 top_global, as_rvec_array(state_global->x.data()), state_global->box, &state_global->swapstate, cr, oenv, Flags);
354 /* Initial values */
355 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
356 &(state_global->fep_state), lam0,
357 nrnb, top_global, &upd,
358 nfile, fnm, &outf, &mdebin,
359 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
361 clear_mat(total_vir);
362 clear_mat(pres);
363 /* Energy terms and groups */
364 snew(enerd, 1);
365 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
366 enerd);
368 /* Kinetic energy data */
369 snew(ekind, 1);
370 init_ekindata(fplog, top_global, &(ir->opts), ekind);
371 /* Copy the cos acceleration to the groups struct */
372 ekind->cosacc.cos_accel = ir->cos_accel;
374 gstat = global_stat_init(ir);
376 /* Check for polarizable models and flexible constraints */
377 shellfc = init_shell_flexcon(fplog,
378 top_global, n_flexible_constraints(constr),
379 ir->nstcalcenergy, DOMAINDECOMP(cr));
381 if (shellfc && ir->nstcalcenergy != 1)
383 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);
385 if (shellfc && DOMAINDECOMP(cr))
387 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
390 if (inputrecDeform(ir))
392 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
393 set_deform_reference_box(upd,
394 deform_init_init_step_tpx,
395 deform_init_box_tpx);
396 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
400 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
401 if ((io > 2000) && MASTER(cr))
403 fprintf(stderr,
404 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
405 io);
409 std::unique_ptr<t_state> stateInstance;
410 t_state * state;
412 if (DOMAINDECOMP(cr))
414 top = dd_init_local_top(top_global);
416 stateInstance = std::unique_ptr<t_state>(new t_state);
417 state = stateInstance.get();
418 dd_init_local_state(cr->dd, state_global, state);
420 else
422 state_change_natoms(state_global, state_global->natoms);
423 /* We need to allocate one element extra, since we might use
424 * (unaligned) 4-wide SIMD loads to access rvec entries.
426 f.resize(state_global->natoms + 1);
427 /* Copy the pointer to the global state */
428 state = state_global;
430 snew(top, 1);
431 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
432 &graph, mdatoms, vsite, shellfc);
434 update_realloc(upd, state->natoms);
437 /* Set up interactive MD (IMD) */
438 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, as_rvec_array(state_global->x.data()),
439 nfile, fnm, oenv, imdport, Flags);
441 if (DOMAINDECOMP(cr))
443 /* Distribute the charge groups over the nodes from the master node */
444 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
445 state_global, top_global, ir,
446 state, &f, mdatoms, top, fr,
447 vsite, constr,
448 nrnb, nullptr, FALSE);
449 shouldCheckNumberOfBondedInteractions = true;
450 update_realloc(upd, state->natoms);
453 update_mdatoms(mdatoms, state->lambda[efptMASS]);
455 startingFromCheckpoint = Flags & MD_STARTFROMCPT;
457 if (ir->bExpanded)
459 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
462 if (MASTER(cr))
464 if (startingFromCheckpoint)
466 /* Update mdebin with energy history if appending to output files */
467 if (Flags & MD_APPENDFILES)
469 restore_energyhistory_from_state(mdebin, energyHistory);
471 else
473 /* We might have read an energy history from checkpoint,
474 * free the allocated memory and reset the counts.
476 *energyHistory = {};
479 /* Set the initial energy history in state by updating once */
480 update_energyhistory(energyHistory, mdebin);
483 /* Initialize constraints */
484 if (constr && !DOMAINDECOMP(cr))
486 set_constraints(constr, top, ir, mdatoms, cr);
489 if (repl_ex_nst > 0 && MASTER(cr))
491 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
492 repl_ex_nst, repl_ex_nex, repl_ex_seed);
495 /* PME tuning is only supported with PME for Coulomb. Is is not supported
496 * with only LJ PME, or for reruns.
498 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
499 !(Flags & MD_REPRODUCIBLE));
500 if (bPMETune)
502 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
503 fr->ic, fr->pmedata, use_GPU(fr->nbv),
504 &bPMETunePrinting);
507 if (!ir->bContinuation && !bRerunMD)
509 if (state->flags & (1 << estV))
511 /* Set the velocities of vsites, shells and frozen atoms to zero */
512 for (i = 0; i < mdatoms->homenr; i++)
514 if (mdatoms->ptype[i] == eptVSite ||
515 mdatoms->ptype[i] == eptShell)
517 clear_rvec(state->v[i]);
519 else if (mdatoms->cFREEZE)
521 for (m = 0; m < DIM; m++)
523 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
525 state->v[i][m] = 0;
532 if (constr)
534 /* Constrain the initial coordinates and velocities */
535 do_constrain_first(fplog, constr, ir, mdatoms, state,
536 cr, nrnb, fr, top);
538 if (vsite)
540 /* Construct the virtual sites for the initial configuration */
541 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, nullptr,
542 top->idef.iparams, top->idef.il,
543 fr->ePBC, fr->bMolPBC, cr, state->box);
547 if (ir->efep != efepNO)
549 /* Set free energy calculation frequency as the greatest common
550 * denominator of nstdhdl and repl_ex_nst. */
551 nstfep = ir->fepvals->nstdhdl;
552 if (ir->bExpanded)
554 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
556 if (repl_ex_nst > 0)
558 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
562 /* Be REALLY careful about what flags you set here. You CANNOT assume
563 * this is the first step, since we might be restarting from a checkpoint,
564 * and in that case we should not do any modifications to the state.
566 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
568 if (Flags & MD_READ_EKIN)
570 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
573 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
574 | (bStopCM ? CGLO_STOPCM : 0)
575 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
576 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
577 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
579 bSumEkinhOld = FALSE;
580 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
581 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
582 constr, &nullSignaller, state->box,
583 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags
584 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
585 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
586 top_global, top, state,
587 &shouldCheckNumberOfBondedInteractions);
588 if (ir->eI == eiVVAK)
590 /* a second call to get the half step temperature initialized as well */
591 /* we do the same call as above, but turn the pressure off -- internally to
592 compute_globals, this is recognized as a velocity verlet half-step
593 kinetic energy calculation. This minimized excess variables, but
594 perhaps loses some logic?*/
596 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
597 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
598 constr, &nullSignaller, state->box,
599 nullptr, &bSumEkinhOld,
600 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
603 /* Calculate the initial half step temperature, and save the ekinh_old */
604 if (!(Flags & MD_STARTFROMCPT))
606 for (i = 0; (i < ir->opts.ngtc); i++)
608 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
611 if (ir->eI != eiVV)
613 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
614 and there is no previous step */
617 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
618 temperature control */
619 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
621 if (MASTER(cr))
623 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
625 fprintf(fplog,
626 "RMS relative constraint deviation after constraining: %.2e\n",
627 constr_rmsd(constr));
629 if (EI_STATE_VELOCITY(ir->eI))
631 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
633 if (bRerunMD)
635 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
636 " input trajectory '%s'\n\n",
637 *(top_global->name), opt2fn("-rerun", nfile, fnm));
638 if (bVerbose)
640 fprintf(stderr, "Calculated time to finish depends on nsteps from "
641 "run input file,\nwhich may not correspond to the time "
642 "needed to process input trajectory.\n\n");
645 else
647 char tbuf[20];
648 fprintf(stderr, "starting mdrun '%s'\n",
649 *(top_global->name));
650 if (ir->nsteps >= 0)
652 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
654 else
656 sprintf(tbuf, "%s", "infinite");
658 if (ir->init_step > 0)
660 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
661 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
662 gmx_step_str(ir->init_step, sbuf2),
663 ir->init_step*ir->delta_t);
665 else
667 fprintf(stderr, "%s steps, %s ps.\n",
668 gmx_step_str(ir->nsteps, sbuf), tbuf);
671 fprintf(fplog, "\n");
674 walltime_accounting_start(walltime_accounting);
675 wallcycle_start(wcycle, ewcRUN);
676 print_start(fplog, cr, walltime_accounting, "mdrun");
678 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
679 #ifdef GMX_FAHCORE
680 chkpt_ret = fcCheckPointParallel( cr->nodeid,
681 NULL, 0);
682 if (chkpt_ret == 0)
684 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
686 #endif
688 /***********************************************************
690 * Loop over MD steps
692 ************************************************************/
694 /* if rerunMD then read coordinates and velocities from input trajectory */
695 if (bRerunMD)
697 if (getenv("GMX_FORCE_UPDATE"))
699 bForceUpdate = TRUE;
702 rerun_fr.natoms = 0;
703 if (MASTER(cr))
705 bLastStep = !read_first_frame(oenv, &status,
706 opt2fn("-rerun", nfile, fnm),
707 &rerun_fr, TRX_NEED_X | TRX_READ_V);
708 if (rerun_fr.natoms != top_global->natoms)
710 gmx_fatal(FARGS,
711 "Number of atoms in trajectory (%d) does not match the "
712 "run input file (%d)\n",
713 rerun_fr.natoms, top_global->natoms);
715 if (ir->ePBC != epbcNONE)
717 if (!rerun_fr.bBox)
719 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);
721 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
723 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
728 if (PAR(cr))
730 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
733 if (ir->ePBC != epbcNONE)
735 /* Set the shift vectors.
736 * Necessary here when have a static box different from the tpr box.
738 calc_shifts(rerun_fr.box, fr->shift_vec);
742 /* loop over MD steps or if rerunMD to end of input trajectory */
743 bFirstStep = TRUE;
744 /* Skip the first Nose-Hoover integration when we get the state from tpx */
745 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
746 bSumEkinhOld = FALSE;
747 bExchanged = FALSE;
748 bNeedRepartition = FALSE;
749 // TODO This implementation of ensemble orientation restraints is nasty because
750 // a user can't just do multi-sim with single-sim orientation restraints.
751 bUsingEnsembleRestraints = (fcd->disres.nsystems > 1) || (cr->ms && fcd->orires.nr);
754 // Replica exchange and ensemble restraints need all
755 // simulations to remain synchronized, so they need
756 // checkpoints and stop conditions to act on the same step, so
757 // the propagation of such signals must take place between
758 // simulations, not just within simulations.
759 bool checkpointIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
760 bool stopConditionIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
761 bool resetCountersIsLocal = true;
762 signals[eglsCHKPT] = SimulationSignal(checkpointIsLocal);
763 signals[eglsSTOPCOND] = SimulationSignal(stopConditionIsLocal);
764 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
767 step = ir->init_step;
768 step_rel = 0;
770 // TODO extract this to new multi-simulation module
771 if (MASTER(cr) && MULTISIM(cr) && (repl_ex_nst <= 0 ))
773 if (!multisim_int_all_are_equal(cr->ms, ir->nsteps))
775 GMX_LOG(mdlog.warning).appendText(
776 "Note: The number of steps is not consistent across multi simulations,\n"
777 "but we are proceeding anyway!");
779 if (!multisim_int_all_are_equal(cr->ms, ir->init_step))
781 GMX_LOG(mdlog.warning).appendText(
782 "Note: The initial step is not consistent across multi simulations,\n"
783 "but we are proceeding anyway!");
787 /* and stop now if we should */
788 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
789 while (!bLastStep)
792 /* Determine if this is a neighbor search step */
793 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
795 if (bPMETune && bNStList)
797 /* PME grid + cut-off optimization with GPUs or PME nodes */
798 pme_loadbal_do(pme_loadbal, cr,
799 (bVerbose && MASTER(cr)) ? stderr : nullptr,
800 fplog, mdlog,
801 ir, fr, state,
802 wcycle,
803 step, step_rel,
804 &bPMETunePrinting);
807 wallcycle_start(wcycle, ewcSTEP);
809 if (bRerunMD)
811 if (rerun_fr.bStep)
813 step = rerun_fr.step;
814 step_rel = step - ir->init_step;
816 if (rerun_fr.bTime)
818 t = rerun_fr.time;
820 else
822 t = step;
825 else
827 bLastStep = (step_rel == ir->nsteps);
828 t = t0 + step*ir->delta_t;
831 // TODO Refactor this, so that nstfep does not need a default value of zero
832 if (ir->efep != efepNO || ir->bSimTemp)
834 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
835 requiring different logic. */
837 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
838 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
839 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
840 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
841 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
844 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
845 do_per_step(step, repl_ex_nst));
847 if (bSimAnn)
849 update_annealing_target_temp(ir, t, upd);
852 if (bRerunMD)
854 if (!DOMAINDECOMP(cr) || MASTER(cr))
856 for (i = 0; i < state_global->natoms; i++)
858 copy_rvec(rerun_fr.x[i], state_global->x[i]);
860 if (rerun_fr.bV)
862 for (i = 0; i < state_global->natoms; i++)
864 copy_rvec(rerun_fr.v[i], state_global->v[i]);
867 else
869 for (i = 0; i < state_global->natoms; i++)
871 clear_rvec(state_global->v[i]);
873 if (bRerunWarnNoV)
875 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
876 " Ekin, temperature and pressure are incorrect,\n"
877 " the virial will be incorrect when constraints are present.\n"
878 "\n");
879 bRerunWarnNoV = FALSE;
883 copy_mat(rerun_fr.box, state_global->box);
884 copy_mat(state_global->box, state->box);
886 if (vsite && (Flags & MD_RERUN_VSITE))
888 if (DOMAINDECOMP(cr))
890 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
892 if (graph)
894 /* Following is necessary because the graph may get out of sync
895 * with the coordinates if we only have every N'th coordinate set
897 mk_mshift(fplog, graph, fr->ePBC, state->box, as_rvec_array(state->x.data()));
898 shift_self(graph, state->box, as_rvec_array(state->x.data()));
900 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
901 top->idef.iparams, top->idef.il,
902 fr->ePBC, fr->bMolPBC, cr, state->box);
903 if (graph)
905 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
910 /* Stop Center of Mass motion */
911 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
913 if (bRerunMD)
915 /* for rerun MD always do Neighbour Searching */
916 bNS = (bFirstStep || ir->nstlist != 0);
918 else
920 /* Determine whether or not to do Neighbour Searching */
921 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
924 /* < 0 means stop at next step, > 0 means stop at next NS step */
925 if ( (signals[eglsSTOPCOND].set < 0) ||
926 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
928 bLastStep = TRUE;
931 /* Determine whether or not to update the Born radii if doing GB */
932 bBornRadii = bFirstStep;
933 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
935 bBornRadii = TRUE;
938 /* do_log triggers energy and virial calculation. Because this leads
939 * to different code paths, forces can be different. Thus for exact
940 * continuation we should avoid extra log output.
941 * Note that the || bLastStep can result in non-exact continuation
942 * beyond the last step. But we don't consider that to be an issue.
944 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
945 do_verbose = bVerbose &&
946 (step % stepout == 0 || bFirstStep || bLastStep || bRerunMD);
948 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
950 if (bRerunMD)
952 bMasterState = TRUE;
954 else
956 bMasterState = FALSE;
957 /* Correct the new box if it is too skewed */
958 if (inputrecDynamicBox(ir))
960 if (correct_box(fplog, step, state->box, graph))
962 bMasterState = TRUE;
965 if (DOMAINDECOMP(cr) && bMasterState)
967 dd_collect_state(cr->dd, state, state_global);
971 if (DOMAINDECOMP(cr))
973 /* Repartition the domain decomposition */
974 dd_partition_system(fplog, step, cr,
975 bMasterState, nstglobalcomm,
976 state_global, top_global, ir,
977 state, &f, mdatoms, top, fr,
978 vsite, constr,
979 nrnb, wcycle,
980 do_verbose && !bPMETunePrinting);
981 shouldCheckNumberOfBondedInteractions = true;
982 update_realloc(upd, state->natoms);
986 if (MASTER(cr) && do_log)
988 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
991 if (ir->efep != efepNO)
993 update_mdatoms(mdatoms, state->lambda[efptMASS]);
996 if ((bRerunMD && rerun_fr.bV) || bExchanged)
999 /* We need the kinetic energy at minus the half step for determining
1000 * the full step kinetic energy and possibly for T-coupling.*/
1001 /* This may not be quite working correctly yet . . . . */
1002 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1003 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1004 constr, &nullSignaller, state->box,
1005 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1006 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1007 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1008 top_global, top, state,
1009 &shouldCheckNumberOfBondedInteractions);
1011 clear_mat(force_vir);
1013 /* We write a checkpoint at this MD step when:
1014 * either at an NS step when we signalled through gs,
1015 * or at the last step (but not when we do not want confout),
1016 * but never at the first step or with rerun.
1018 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1019 (bLastStep && (Flags & MD_CONFOUT))) &&
1020 step > ir->init_step && !bRerunMD);
1021 if (bCPT)
1023 signals[eglsCHKPT].set = 0;
1026 /* Determine the energy and pressure:
1027 * at nstcalcenergy steps and at energy output steps (set below).
1029 if (EI_VV(ir->eI) && (!bInitStep))
1031 /* for vv, the first half of the integration actually corresponds
1032 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1033 but the virial needs to be calculated on both the current step and the 'next' step. Future
1034 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1036 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1037 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1038 bCalcVir = bCalcEnerStep ||
1039 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1041 else
1043 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1044 bCalcVir = bCalcEnerStep ||
1045 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1047 bCalcEner = bCalcEnerStep;
1049 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1051 if (do_ene || do_log || bDoReplEx)
1053 bCalcVir = TRUE;
1054 bCalcEner = TRUE;
1057 /* Do we need global communication ? */
1058 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1059 do_per_step(step, nstglobalcomm) ||
1060 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1062 force_flags = (GMX_FORCE_STATECHANGED |
1063 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1064 GMX_FORCE_ALLFORCES |
1065 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1066 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1067 (bDoFEP ? GMX_FORCE_DHDL : 0)
1070 if (shellfc)
1072 /* Now is the time to relax the shells */
1073 relax_shell_flexcon(fplog, cr, bVerbose, step,
1074 ir, bNS, force_flags, top,
1075 constr, enerd, fcd,
1076 state, &f, force_vir, mdatoms,
1077 nrnb, wcycle, graph, groups,
1078 shellfc, fr, bBornRadii, t, mu_tot,
1079 vsite);
1081 else
1083 /* The coordinates (x) are shifted (to get whole molecules)
1084 * in do_force.
1085 * This is parallellized as well, and does communication too.
1086 * Check comments in sim_util.c
1088 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1089 state->box, &state->x, &state->hist,
1090 &f, force_vir, mdatoms, enerd, fcd,
1091 state->lambda, graph,
1092 fr, vsite, mu_tot, t, ed, bBornRadii,
1093 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1096 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1097 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1099 rvec *vbuf = nullptr;
1101 wallcycle_start(wcycle, ewcUPDATE);
1102 if (ir->eI == eiVV && bInitStep)
1104 /* if using velocity verlet with full time step Ekin,
1105 * take the first half step only to compute the
1106 * virial for the first step. From there,
1107 * revert back to the initial coordinates
1108 * so that the input is actually the initial step.
1110 snew(vbuf, state->natoms);
1111 copy_rvecn(as_rvec_array(state->v.data()), vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1113 else
1115 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1116 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1119 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1120 ekind, M, upd, etrtVELOCITY1,
1121 cr, constr);
1123 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1125 wallcycle_stop(wcycle, ewcUPDATE);
1126 update_constraints(fplog, step, nullptr, ir, mdatoms,
1127 state, fr->bMolPBC, graph, &f,
1128 &top->idef, shake_vir,
1129 cr, nrnb, wcycle, upd, constr,
1130 TRUE, bCalcVir);
1131 wallcycle_start(wcycle, ewcUPDATE);
1133 else if (graph)
1135 /* Need to unshift here if a do_force has been
1136 called in the previous step */
1137 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1139 /* if VV, compute the pressure and constraints */
1140 /* For VV2, we strictly only need this if using pressure
1141 * control, but we really would like to have accurate pressures
1142 * printed out.
1143 * Think about ways around this in the future?
1144 * For now, keep this choice in comments.
1146 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1147 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1148 bPres = TRUE;
1149 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1150 if (bCalcEner && ir->eI == eiVVAK)
1152 bSumEkinhOld = TRUE;
1154 /* for vv, the first half of the integration actually corresponds to the previous step.
1155 So we need information from the last step in the first half of the integration */
1156 if (bGStat || do_per_step(step-1, nstglobalcomm))
1158 wallcycle_stop(wcycle, ewcUPDATE);
1159 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1160 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1161 constr, &nullSignaller, state->box,
1162 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1163 (bGStat ? CGLO_GSTAT : 0)
1164 | CGLO_ENERGY
1165 | (bTemp ? CGLO_TEMPERATURE : 0)
1166 | (bPres ? CGLO_PRESSURE : 0)
1167 | (bPres ? CGLO_CONSTRAINT : 0)
1168 | (bStopCM ? CGLO_STOPCM : 0)
1169 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1170 | CGLO_SCALEEKIN
1172 /* explanation of above:
1173 a) We compute Ekin at the full time step
1174 if 1) we are using the AveVel Ekin, and it's not the
1175 initial step, or 2) if we are using AveEkin, but need the full
1176 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1177 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1178 EkinAveVel because it's needed for the pressure */
1179 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1180 top_global, top, state,
1181 &shouldCheckNumberOfBondedInteractions);
1182 wallcycle_start(wcycle, ewcUPDATE);
1184 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1185 if (!bInitStep)
1187 if (bTrotter)
1189 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1190 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1192 copy_mat(shake_vir, state->svir_prev);
1193 copy_mat(force_vir, state->fvir_prev);
1194 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1196 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1197 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, nullptr, (ir->eI == eiVV), FALSE);
1198 enerd->term[F_EKIN] = trace(ekind->ekin);
1201 else if (bExchanged)
1203 wallcycle_stop(wcycle, ewcUPDATE);
1204 /* We need the kinetic energy at minus the half step for determining
1205 * the full step kinetic energy and possibly for T-coupling.*/
1206 /* This may not be quite working correctly yet . . . . */
1207 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1208 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
1209 constr, &nullSignaller, state->box,
1210 nullptr, &bSumEkinhOld,
1211 CGLO_GSTAT | CGLO_TEMPERATURE);
1212 wallcycle_start(wcycle, ewcUPDATE);
1215 /* if it's the initial step, we performed this first step just to get the constraint virial */
1216 if (ir->eI == eiVV && bInitStep)
1218 copy_rvecn(vbuf, as_rvec_array(state->v.data()), 0, state->natoms);
1219 sfree(vbuf);
1221 wallcycle_stop(wcycle, ewcUPDATE);
1224 /* compute the conserved quantity */
1225 if (EI_VV(ir->eI))
1227 saved_conserved_quantity = NPT_energy(ir, state, &MassQ);
1228 if (ir->eI == eiVV)
1230 last_ekin = enerd->term[F_EKIN];
1232 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1234 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1236 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1237 if (ir->efep != efepNO && !bRerunMD)
1239 sum_dhdl(enerd, state->lambda, ir->fepvals);
1243 /* ######## END FIRST UPDATE STEP ############## */
1244 /* ######## If doing VV, we now have v(dt) ###### */
1245 if (bDoExpanded)
1247 /* perform extended ensemble sampling in lambda - we don't
1248 actually move to the new state before outputting
1249 statistics, but if performing simulated tempering, we
1250 do update the velocities and the tau_t. */
1252 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, as_rvec_array(state->v.data()), mdatoms);
1253 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1254 copy_df_history(state_global->dfhist, state->dfhist);
1257 /* Now we have the energies and forces corresponding to the
1258 * coordinates at time t. We must output all of this before
1259 * the update.
1261 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1262 ir, state, state_global, energyHistory,
1263 top_global, fr,
1264 outf, mdebin, ekind, &f,
1265 &nchkpt,
1266 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1267 bSumEkinhOld);
1268 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1269 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, as_rvec_array(state->x.data()), ir, t, wcycle);
1271 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1272 if (startingFromCheckpoint && bTrotter)
1274 copy_mat(state->svir_prev, shake_vir);
1275 copy_mat(state->fvir_prev, force_vir);
1278 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1280 /* Check whether everything is still allright */
1281 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1282 #if GMX_THREAD_MPI
1283 && MASTER(cr)
1284 #endif
1287 int nsteps_stop = -1;
1289 /* this just makes signals[].sig compatible with the hack
1290 of sending signals around by MPI_Reduce together with
1291 other floats */
1292 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1294 signals[eglsSTOPCOND].sig = 1;
1295 nsteps_stop = std::max(ir->nstlist, 2*nstglobalcomm);
1297 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1299 signals[eglsSTOPCOND].sig = -1;
1300 nsteps_stop = nstglobalcomm + 1;
1302 if (fplog)
1304 fprintf(fplog,
1305 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1306 gmx_get_signal_name(), nsteps_stop);
1307 fflush(fplog);
1309 fprintf(stderr,
1310 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1311 gmx_get_signal_name(), nsteps_stop);
1312 fflush(stderr);
1313 handled_stop_condition = (int)gmx_get_stop_condition();
1315 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1316 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1317 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1319 /* Signal to terminate the run */
1320 signals[eglsSTOPCOND].sig = 1;
1321 if (fplog)
1323 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1325 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1328 if (bResetCountersHalfMaxH && MASTER(cr) &&
1329 elapsed_time > max_hours*60.0*60.0*0.495)
1331 /* Set flag that will communicate the signal to all ranks in the simulation */
1332 signals[eglsRESETCOUNTERS].sig = 1;
1335 /* In parallel we only have to check for checkpointing in steps
1336 * where we do global communication,
1337 * otherwise the other nodes don't know.
1339 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1340 cpt_period >= 0 &&
1341 (cpt_period == 0 ||
1342 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1343 signals[eglsCHKPT].set == 0)
1345 signals[eglsCHKPT].sig = 1;
1348 /* ######### START SECOND UPDATE STEP ################# */
1350 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1351 in preprocessing */
1353 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1355 gmx_bool bIfRandomize;
1356 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1357 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1358 if (constr && bIfRandomize)
1360 update_constraints(fplog, step, nullptr, ir, mdatoms,
1361 state, fr->bMolPBC, graph, &f,
1362 &top->idef, tmp_vir,
1363 cr, nrnb, wcycle, upd, constr,
1364 TRUE, bCalcVir);
1367 /* Box is changed in update() when we do pressure coupling,
1368 * but we should still use the old box for energy corrections and when
1369 * writing it to the energy file, so it matches the trajectory files for
1370 * the same timestep above. Make a copy in a separate array.
1372 copy_mat(state->box, lastbox);
1374 dvdl_constr = 0;
1376 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1378 wallcycle_start(wcycle, ewcUPDATE);
1379 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1380 if (bTrotter)
1382 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1383 /* We can only do Berendsen coupling after we have summed
1384 * the kinetic energy or virial. Since the happens
1385 * in global_state after update, we should only do it at
1386 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1389 else
1391 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1392 update_pcouple_before_coordinates(fplog, step, ir, state,
1393 parrinellorahmanMu, M,
1394 bInitStep);
1397 if (EI_VV(ir->eI))
1399 /* velocity half-step update */
1400 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1401 ekind, M, upd, etrtVELOCITY2,
1402 cr, constr);
1405 /* Above, initialize just copies ekinh into ekin,
1406 * it doesn't copy position (for VV),
1407 * and entire integrator for MD.
1410 if (ir->eI == eiVVAK)
1412 /* We probably only need md->homenr, not state->natoms */
1413 if (state->natoms > cbuf_nalloc)
1415 cbuf_nalloc = state->natoms;
1416 srenew(cbuf, cbuf_nalloc);
1418 copy_rvecn(as_rvec_array(state->x.data()), cbuf, 0, state->natoms);
1421 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1422 ekind, M, upd, etrtPOSITION, cr, constr);
1423 wallcycle_stop(wcycle, ewcUPDATE);
1425 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1426 fr->bMolPBC, graph, &f,
1427 &top->idef, shake_vir,
1428 cr, nrnb, wcycle, upd, constr,
1429 FALSE, bCalcVir);
1431 if (ir->eI == eiVVAK)
1433 /* erase F_EKIN and F_TEMP here? */
1434 /* just compute the kinetic energy at the half step to perform a trotter step */
1435 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1436 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1437 constr, &nullSignaller, lastbox,
1438 nullptr, &bSumEkinhOld,
1439 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1441 wallcycle_start(wcycle, ewcUPDATE);
1442 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1443 /* now we know the scaling, we can compute the positions again again */
1444 copy_rvecn(cbuf, as_rvec_array(state->x.data()), 0, state->natoms);
1446 update_coords(fplog, step, ir, mdatoms, state, &f, fcd,
1447 ekind, M, upd, etrtPOSITION, cr, constr);
1448 wallcycle_stop(wcycle, ewcUPDATE);
1450 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1451 /* are the small terms in the shake_vir here due
1452 * to numerical errors, or are they important
1453 * physically? I'm thinking they are just errors, but not completely sure.
1454 * For now, will call without actually constraining, constr=NULL*/
1455 update_constraints(fplog, step, nullptr, ir, mdatoms,
1456 state, fr->bMolPBC, graph, &f,
1457 &top->idef, tmp_vir,
1458 cr, nrnb, wcycle, upd, nullptr,
1459 FALSE, bCalcVir);
1461 if (EI_VV(ir->eI))
1463 /* this factor or 2 correction is necessary
1464 because half of the constraint force is removed
1465 in the vv step, so we have to double it. See
1466 the Redmine issue #1255. It is not yet clear
1467 if the factor of 2 is exact, or just a very
1468 good approximation, and this will be
1469 investigated. The next step is to see if this
1470 can be done adding a dhdl contribution from the
1471 rattle step, but this is somewhat more
1472 complicated with the current code. Will be
1473 investigated, hopefully for 4.6.3. However,
1474 this current solution is much better than
1475 having it completely wrong.
1477 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1479 else
1481 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1484 else if (graph)
1486 /* Need to unshift here */
1487 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1490 if (vsite != nullptr)
1492 wallcycle_start(wcycle, ewcVSITECONSTR);
1493 if (graph != nullptr)
1495 shift_self(graph, state->box, as_rvec_array(state->x.data()));
1497 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
1498 top->idef.iparams, top->idef.il,
1499 fr->ePBC, fr->bMolPBC, cr, state->box);
1501 if (graph != nullptr)
1503 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
1505 wallcycle_stop(wcycle, ewcVSITECONSTR);
1508 /* ############## IF NOT VV, Calculate globals HERE ############ */
1509 /* With Leap-Frog we can skip compute_globals at
1510 * non-communication steps, but we need to calculate
1511 * the kinetic energy one step before communication.
1514 // Organize to do inter-simulation signalling on steps if
1515 // and when algorithms require it.
1516 bool doInterSimSignal = (!bFirstStep && bDoReplEx) || bUsingEnsembleRestraints;
1518 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1520 // Since we're already communicating at this step, we
1521 // can propagate intra-simulation signals. Note that
1522 // check_nstglobalcomm has the responsibility for
1523 // choosing the value of nstglobalcomm that is one way
1524 // bGStat becomes true, so we can't get into a
1525 // situation where e.g. checkpointing can't be
1526 // signalled.
1527 bool doIntraSimSignal = true;
1528 SimulationSignaller signaller(&signals, cr, doInterSimSignal, doIntraSimSignal);
1530 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1531 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1532 constr, &signaller,
1533 lastbox,
1534 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1535 (bGStat ? CGLO_GSTAT : 0)
1536 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1537 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1538 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1539 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1540 | CGLO_CONSTRAINT
1541 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1543 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1544 top_global, top, state,
1545 &shouldCheckNumberOfBondedInteractions);
1549 /* ############# END CALC EKIN AND PRESSURE ################# */
1551 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1552 the virial that should probably be addressed eventually. state->veta has better properies,
1553 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1554 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1556 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1558 /* Sum up the foreign energy and dhdl terms for md and sd.
1559 Currently done every step so that dhdl is correct in the .edr */
1560 sum_dhdl(enerd, state->lambda, ir->fepvals);
1563 update_pcouple_after_coordinates(fplog, step, ir, mdatoms,
1564 pres, parrinellorahmanMu,
1565 state, nrnb, upd);
1567 /* ################# END UPDATE STEP 2 ################# */
1568 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1570 /* The coordinates (x) were unshifted in update */
1571 if (!bGStat)
1573 /* We will not sum ekinh_old,
1574 * so signal that we still have to do it.
1576 bSumEkinhOld = TRUE;
1579 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1581 /* use the directly determined last velocity, not actually the averaged half steps */
1582 if (bTrotter && ir->eI == eiVV)
1584 enerd->term[F_EKIN] = last_ekin;
1586 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1588 if (EI_VV(ir->eI))
1590 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1592 else
1594 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + NPT_energy(ir, state, &MassQ);
1596 /* ######### END PREPARING EDR OUTPUT ########### */
1598 /* Output stuff */
1599 if (MASTER(cr))
1601 if (fplog && do_log && bDoExpanded)
1603 /* only needed if doing expanded ensemble */
1604 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : nullptr,
1605 state_global->dfhist, state->fep_state, ir->nstlog, step);
1607 if (bCalcEner)
1609 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1610 t, mdatoms->tmass, enerd, state,
1611 ir->fepvals, ir->expandedvals, lastbox,
1612 shake_vir, force_vir, total_vir, pres,
1613 ekind, mu_tot, constr);
1615 else
1617 upd_mdebin_step(mdebin);
1620 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1621 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1623 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : nullptr,
1624 step, t,
1625 eprNORMAL, mdebin, fcd, groups, &(ir->opts));
1627 if (ir->bPull)
1629 pull_print_output(ir->pull_work, step, t);
1632 if (do_per_step(step, ir->nstlog))
1634 if (fflush(fplog) != 0)
1636 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1640 if (bDoExpanded)
1642 /* Have to do this part _after_ outputting the logfile and the edr file */
1643 /* Gets written into the state at the beginning of next loop*/
1644 state->fep_state = lamnew;
1646 /* Print the remaining wall clock time for the run */
1647 if (MULTIMASTER(cr) &&
1648 (do_verbose || gmx_got_usr_signal()) &&
1649 !bPMETunePrinting)
1651 if (shellfc)
1653 fprintf(stderr, "\n");
1655 print_time(stderr, walltime_accounting, step, ir, cr);
1658 /* Ion/water position swapping.
1659 * Not done in last step since trajectory writing happens before this call
1660 * in the MD loop and exchanges would be lost anyway. */
1661 bNeedRepartition = FALSE;
1662 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1663 do_per_step(step, ir->swap->nstswap))
1665 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1666 bRerunMD ? rerun_fr.x : as_rvec_array(state->x.data()),
1667 bRerunMD ? rerun_fr.box : state->box,
1668 MASTER(cr) && bVerbose, bRerunMD);
1670 if (bNeedRepartition && DOMAINDECOMP(cr))
1672 dd_collect_state(cr->dd, state, state_global);
1676 /* Replica exchange */
1677 bExchanged = FALSE;
1678 if (bDoReplEx)
1680 bExchanged = replica_exchange(fplog, cr, repl_ex,
1681 state_global, enerd,
1682 state, step, t);
1685 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1687 dd_partition_system(fplog, step, cr, TRUE, 1,
1688 state_global, top_global, ir,
1689 state, &f, mdatoms, top, fr,
1690 vsite, constr,
1691 nrnb, wcycle, FALSE);
1692 shouldCheckNumberOfBondedInteractions = true;
1693 update_realloc(upd, state->natoms);
1696 bFirstStep = FALSE;
1697 bInitStep = FALSE;
1698 startingFromCheckpoint = FALSE;
1700 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1701 /* With all integrators, except VV, we need to retain the pressure
1702 * at the current step for coupling at the next step.
1704 if ((state->flags & (1<<estPRES_PREV)) &&
1705 (bGStatEveryStep ||
1706 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1708 /* Store the pressure in t_state for pressure coupling
1709 * at the next MD step.
1711 copy_mat(pres, state->pres_prev);
1714 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1716 if ( (membed != nullptr) && (!bLastStep) )
1718 rescale_membed(step_rel, membed, as_rvec_array(state_global->x.data()));
1721 if (bRerunMD)
1723 if (MASTER(cr))
1725 /* read next frame from input trajectory */
1726 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1729 if (PAR(cr))
1731 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1735 cycles = wallcycle_stop(wcycle, ewcSTEP);
1736 if (DOMAINDECOMP(cr) && wcycle)
1738 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1741 if (!bRerunMD || !rerun_fr.bStep)
1743 /* increase the MD step number */
1744 step++;
1745 step_rel++;
1748 /* TODO make a counter-reset module */
1749 /* If it is time to reset counters, set a flag that remains
1750 true until counters actually get reset */
1751 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1752 signals[eglsRESETCOUNTERS].set != 0)
1754 if (pme_loadbal_is_active(pme_loadbal))
1756 /* Do not permit counter reset while PME load
1757 * balancing is active. The only purpose for resetting
1758 * counters is to measure reliable performance data,
1759 * and that can't be done before balancing
1760 * completes.
1762 * TODO consider fixing this by delaying the reset
1763 * until after load balancing completes,
1764 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1765 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1766 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1767 "resetting counters later in the run, e.g. with gmx "
1768 "mdrun -resetstep.", step);
1770 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1771 use_GPU(fr->nbv) ? fr->nbv : nullptr);
1772 wcycle_set_reset_counters(wcycle, -1);
1773 if (!(cr->duty & DUTY_PME))
1775 /* Tell our PME node to reset its counters */
1776 gmx_pme_send_resetcounters(cr, step);
1778 /* Correct max_hours for the elapsed time */
1779 max_hours -= elapsed_time/(60.0*60.0);
1780 /* If mdrun -maxh -resethway was active, it can only trigger once */
1781 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1782 /* Reset can only happen once, so clear the triggering flag. */
1783 signals[eglsRESETCOUNTERS].set = 0;
1786 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1787 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1790 /* End of main MD loop */
1792 /* Closing TNG files can include compressing data. Therefore it is good to do that
1793 * before stopping the time measurements. */
1794 mdoutf_tng_close(outf);
1796 /* Stop measuring walltime */
1797 walltime_accounting_end(walltime_accounting);
1799 if (bRerunMD && MASTER(cr))
1801 close_trj(status);
1804 if (!(cr->duty & DUTY_PME))
1806 /* Tell the PME only node to finish */
1807 gmx_pme_send_finish(cr);
1810 if (MASTER(cr))
1812 if (ir->nstcalcenergy > 0 && !bRerunMD)
1814 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1815 eprAVER, mdebin, fcd, groups, &(ir->opts));
1819 done_mdoutf(outf, ir);
1821 if (bPMETune)
1823 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1826 done_shellfc(fplog, shellfc, step_rel);
1828 if (repl_ex_nst > 0 && MASTER(cr))
1830 print_replica_exchange_statistics(fplog, repl_ex);
1833 // Clean up swapcoords
1834 if (ir->eSwapCoords != eswapNO)
1836 finish_swapcoords(ir->swap);
1839 /* IMD cleanup, if bIMD is TRUE. */
1840 IMD_finalize(ir->bIMD, ir->imd);
1842 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1843 if (step_rel >= wcycle_get_reset_counters(wcycle) &&
1844 signals[eglsRESETCOUNTERS].set == 0 &&
1845 !bResetCountersHalfMaxH)
1847 walltime_accounting_set_valid_finish(walltime_accounting);
1850 return 0;