Remove state duplication in serial runs
[gromacs.git] / src / programs / mdrun / md.cpp
blob03ec3a0590401679596e35416f370cb29777b230
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, 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>
49 #include "thread_mpi/threads.h"
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_network.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/ewald/pme-load-balancing.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nrnb.h"
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/imd/imd.h"
62 #include "gromacs/listed-forces/manage-threading.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/utilities.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/math/vectypes.h"
67 #include "gromacs/mdlib/compute_io.h"
68 #include "gromacs/mdlib/constr.h"
69 #include "gromacs/mdlib/ebin.h"
70 #include "gromacs/mdlib/force.h"
71 #include "gromacs/mdlib/force_flags.h"
72 #include "gromacs/mdlib/forcerec.h"
73 #include "gromacs/mdlib/md_support.h"
74 #include "gromacs/mdlib/mdatoms.h"
75 #include "gromacs/mdlib/mdebin.h"
76 #include "gromacs/mdlib/mdoutf.h"
77 #include "gromacs/mdlib/mdrun.h"
78 #include "gromacs/mdlib/mdsetup.h"
79 #include "gromacs/mdlib/nb_verlet.h"
80 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
81 #include "gromacs/mdlib/ns.h"
82 #include "gromacs/mdlib/shellfc.h"
83 #include "gromacs/mdlib/sighandler.h"
84 #include "gromacs/mdlib/sim_util.h"
85 #include "gromacs/mdlib/simulationsignal.h"
86 #include "gromacs/mdlib/tgroup.h"
87 #include "gromacs/mdlib/trajectory_writing.h"
88 #include "gromacs/mdlib/update.h"
89 #include "gromacs/mdlib/vcm.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdtypes/commrec.h"
92 #include "gromacs/mdtypes/df_history.h"
93 #include "gromacs/mdtypes/energyhistory.h"
94 #include "gromacs/mdtypes/fcdata.h"
95 #include "gromacs/mdtypes/forcerec.h"
96 #include "gromacs/mdtypes/group.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/interaction_const.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/mdtypes/mdatom.h"
101 #include "gromacs/mdtypes/state.h"
102 #include "gromacs/pbcutil/mshift.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/pulling/pull.h"
105 #include "gromacs/swap/swapcoords.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/timing/walltime_accounting.h"
108 #include "gromacs/topology/atoms.h"
109 #include "gromacs/topology/idef.h"
110 #include "gromacs/topology/mtop_util.h"
111 #include "gromacs/topology/topology.h"
112 #include "gromacs/trajectory/trajectoryframe.h"
113 #include "gromacs/utility/basedefinitions.h"
114 #include "gromacs/utility/cstringutil.h"
115 #include "gromacs/utility/fatalerror.h"
116 #include "gromacs/utility/logger.h"
117 #include "gromacs/utility/real.h"
118 #include "gromacs/utility/smalloc.h"
120 #include "deform.h"
121 #include "membed.h"
122 #include "repl_ex.h"
124 #ifdef GMX_FAHCORE
125 #include "corewrap.h"
126 #endif
128 using gmx::SimulationSignaller;
130 /*! \brief Check whether bonded interactions are missing, if appropriate
132 * \param[in] fplog Log file pointer
133 * \param[in] cr Communication object
134 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
135 * \param[in] top_global Global topology for the error message
136 * \param[in] top_local Local topology for the error message
137 * \param[in] state Global state for the error message
138 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
140 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
141 * is always set to false after exit.
143 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
144 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
145 bool *shouldCheckNumberOfBondedInteractions)
147 if (*shouldCheckNumberOfBondedInteractions)
149 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
151 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
153 *shouldCheckNumberOfBondedInteractions = false;
157 static void reset_all_counters(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
158 gmx_int64_t step,
159 gmx_int64_t *step_rel, t_inputrec *ir,
160 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
161 gmx_walltime_accounting_t walltime_accounting,
162 struct nonbonded_verlet_t *nbv)
164 char sbuf[STEPSTRSIZE];
166 /* Reset all the counters related to performance over the run */
167 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
168 "step %s: resetting all time and cycle counters",
169 gmx_step_str(step, sbuf));
171 if (use_GPU(nbv))
173 nbnxn_gpu_reset_timings(nbv);
174 resetGpuProfiler();
177 wallcycle_stop(wcycle, ewcRUN);
178 wallcycle_reset_all(wcycle);
179 if (DOMAINDECOMP(cr))
181 reset_dd_statistics_counters(cr->dd);
183 init_nrnb(nrnb);
184 ir->init_step += *step_rel;
185 ir->nsteps -= *step_rel;
186 *step_rel = 0;
187 wallcycle_start(wcycle, ewcRUN);
188 walltime_accounting_start(walltime_accounting);
189 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
192 /*! \libinternal
193 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
194 int nfile, const t_filenm fnm[],
195 const gmx_output_env_t *oenv, gmx_bool bVerbose,
196 int nstglobalcomm,
197 gmx_vsite_t *vsite, gmx_constr_t constr,
198 int stepout,
199 t_inputrec *inputrec,
200 gmx_mtop_t *top_global, t_fcdata *fcd,
201 t_state *state_global,
202 t_mdatoms *mdatoms,
203 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
204 gmx_edsam_t ed,
205 t_forcerec *fr,
206 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
207 real cpt_period, real max_hours,
208 int imdport,
209 unsigned long Flags,
210 gmx_walltime_accounting_t walltime_accounting)
212 double gmx::do_md(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
213 int nfile, const t_filenm fnm[],
214 const gmx_output_env_t *oenv, gmx_bool bVerbose,
215 int nstglobalcomm,
216 gmx_vsite_t *vsite, gmx_constr_t constr,
217 int stepout, t_inputrec *ir,
218 gmx_mtop_t *top_global,
219 t_fcdata *fcd,
220 t_state *state_global,
221 t_mdatoms *mdatoms,
222 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
223 gmx_edsam_t ed, t_forcerec *fr,
224 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
225 gmx_membed_t *membed,
226 real cpt_period, real max_hours,
227 int imdport,
228 unsigned long Flags,
229 gmx_walltime_accounting_t walltime_accounting)
231 gmx_mdoutf_t outf = NULL;
232 gmx_int64_t step, step_rel;
233 double elapsed_time;
234 double t, t0, lam0[efptNR];
235 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
236 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD,
237 bFirstStep, startingFromCheckpoint, bInitStep, bLastStep = FALSE,
238 bBornRadii, bUsingEnsembleRestraints;
239 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
240 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
241 bForceUpdate = FALSE, bCPT;
242 gmx_bool bMasterState;
243 int force_flags, cglo_flags;
244 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
245 int i, m;
246 t_trxstatus *status;
247 rvec mu_tot;
248 t_vcm *vcm;
249 matrix pcoupl_mu, M;
250 t_trxframe rerun_fr;
251 gmx_repl_ex_t repl_ex = NULL;
252 int nchkpt = 1;
253 gmx_localtop_t *top;
254 t_mdebin *mdebin = NULL;
255 t_state *state = NULL;
256 gmx_enerdata_t *enerd;
257 rvec *f = NULL;
258 gmx_global_stat_t gstat;
259 gmx_update_t *upd = NULL;
260 t_graph *graph = NULL;
261 gmx_groups_t *groups;
262 gmx_ekindata_t *ekind;
263 gmx_shellfc_t *shellfc;
264 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
265 gmx_bool bResetCountersHalfMaxH = FALSE;
266 gmx_bool bTemp, bPres, bTrotter;
267 real dvdl_constr;
268 rvec *cbuf = NULL;
269 int cbuf_nalloc = 0;
270 matrix lastbox;
271 int lamnew = 0;
272 /* for FEP */
273 int nstfep = 0;
274 double cycles;
275 real saved_conserved_quantity = 0;
276 real last_ekin = 0;
277 t_extmass MassQ;
278 int **trotter_seq;
279 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
280 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
283 /* PME load balancing data for GPU kernels */
284 pme_load_balancing_t *pme_loadbal = NULL;
285 gmx_bool bPMETune = FALSE;
286 gmx_bool bPMETunePrinting = FALSE;
288 /* Interactive MD */
289 gmx_bool bIMDstep = FALSE;
291 #ifdef GMX_FAHCORE
292 /* Temporary addition for FAHCORE checkpointing */
293 int chkpt_ret;
294 #endif
295 /* Domain decomposition could incorrectly miss a bonded
296 interaction, but checking for that requires a global
297 communication stage, which does not otherwise happen in DD
298 code. So we do that alongside the first global energy reduction
299 after a new DD is made. These variables handle whether the
300 check happens, and the result it returns. */
301 bool shouldCheckNumberOfBondedInteractions = false;
302 int totalNumberOfBondedInteractions = -1;
304 SimulationSignals signals;
305 // Most global communnication stages don't propagate mdrun
306 // signals, and will use this object to achieve that.
307 SimulationSignaller nullSignaller(nullptr, nullptr, false, false);
309 /* Check for special mdrun options */
310 bRerunMD = (Flags & MD_RERUN);
311 if (Flags & MD_RESETCOUNTERSHALFWAY)
313 if (ir->nsteps > 0)
315 /* Signal to reset the counters half the simulation steps. */
316 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
318 /* Signal to reset the counters halfway the simulation time. */
319 bResetCountersHalfMaxH = (max_hours > 0);
322 /* md-vv uses averaged full step velocities for T-control
323 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
324 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
325 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
327 if (bRerunMD)
329 /* Since we don't know if the frames read are related in any way,
330 * rebuild the neighborlist at every step.
332 ir->nstlist = 1;
333 ir->nstcalcenergy = 1;
334 nstglobalcomm = 1;
337 nstglobalcomm = check_nstglobalcomm(mdlog, nstglobalcomm, ir);
338 bGStatEveryStep = (nstglobalcomm == 1);
340 if (bRerunMD)
342 ir->nstxout_compressed = 0;
344 groups = &top_global->groups;
346 if (ir->eSwapCoords != eswapNO)
348 /* Initialize ion swapping code */
349 init_swapcoords(fplog, bVerbose, ir, opt2fn_master("-swap", nfile, fnm, cr),
350 top_global, state_global->x, state_global->box, &state_global->swapstate, cr, oenv, Flags);
353 /* Initial values */
354 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
355 &(state_global->fep_state), lam0,
356 nrnb, top_global, &upd,
357 nfile, fnm, &outf, &mdebin,
358 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
360 clear_mat(total_vir);
361 clear_mat(pres);
362 /* Energy terms and groups */
363 snew(enerd, 1);
364 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
365 enerd);
366 if (DOMAINDECOMP(cr))
368 f = NULL;
370 else
372 snew(f, top_global->natoms);
375 /* Kinetic energy data */
376 snew(ekind, 1);
377 init_ekindata(fplog, top_global, &(ir->opts), ekind);
378 /* Copy the cos acceleration to the groups struct */
379 ekind->cosacc.cos_accel = ir->cos_accel;
381 gstat = global_stat_init(ir);
383 /* Check for polarizable models and flexible constraints */
384 shellfc = init_shell_flexcon(fplog,
385 top_global, n_flexible_constraints(constr),
386 ir->nstcalcenergy, DOMAINDECOMP(cr));
388 if (shellfc && ir->nstcalcenergy != 1)
390 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);
392 if (shellfc && DOMAINDECOMP(cr))
394 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
397 if (inputrecDeform(ir))
399 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
400 set_deform_reference_box(upd,
401 deform_init_init_step_tpx,
402 deform_init_box_tpx);
403 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
407 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
408 if ((io > 2000) && MASTER(cr))
410 fprintf(stderr,
411 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
412 io);
416 if (DOMAINDECOMP(cr))
418 top = dd_init_local_top(top_global);
420 snew(state, 1);
421 dd_init_local_state(cr->dd, state_global, state);
423 else
425 /* Copy the pointer to the global state */
426 state = state_global;
428 snew(top, 1);
429 mdAlgorithmsSetupAtomData(cr, ir, top_global, top, fr,
430 &graph, mdatoms, vsite, shellfc);
432 update_realloc(upd, state->nalloc);
435 /* Set up interactive MD (IMD) */
436 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
437 nfile, fnm, oenv, imdport, Flags);
439 if (DOMAINDECOMP(cr))
441 /* Distribute the charge groups over the nodes from the master node */
442 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
443 state_global, top_global, ir,
444 state, &f, mdatoms, top, fr,
445 vsite, constr,
446 nrnb, NULL, FALSE);
447 shouldCheckNumberOfBondedInteractions = true;
448 update_realloc(upd, state->nalloc);
451 update_mdatoms(mdatoms, state->lambda[efptMASS]);
453 startingFromCheckpoint = Flags & MD_STARTFROMCPT;
455 if (ir->bExpanded)
457 init_expanded_ensemble(startingFromCheckpoint, ir, state->dfhist);
460 if (MASTER(cr))
462 if (startingFromCheckpoint)
464 /* Update mdebin with energy history if appending to output files */
465 if (Flags & MD_APPENDFILES)
467 restore_energyhistory_from_state(mdebin, state_global->enerhist);
469 else
471 /* We might have read an energy history from checkpoint,
472 * free the allocated memory and reset the counts.
474 done_energyhistory(state_global->enerhist);
475 init_energyhistory(state_global->enerhist);
478 /* Set the initial energy history in state by updating once */
479 update_energyhistory(state_global->enerhist, mdebin);
482 /* Initialize constraints */
483 if (constr && !DOMAINDECOMP(cr))
485 set_constraints(constr, top, ir, mdatoms, cr);
488 if (repl_ex_nst > 0 && MASTER(cr))
490 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
491 repl_ex_nst, repl_ex_nex, repl_ex_seed);
494 /* PME tuning is only supported with PME for Coulomb. Is is not supported
495 * with only LJ PME, or for reruns.
497 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
498 !(Flags & MD_REPRODUCIBLE));
499 if (bPMETune)
501 pme_loadbal_init(&pme_loadbal, cr, mdlog, ir, state->box,
502 fr->ic, fr->pmedata, use_GPU(fr->nbv),
503 &bPMETunePrinting);
506 if (!ir->bContinuation && !bRerunMD)
508 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
510 /* Set the velocities of frozen particles to zero */
511 for (i = 0; i < mdatoms->homenr; i++)
513 for (m = 0; m < DIM; m++)
515 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
517 state->v[i][m] = 0;
523 if (constr)
525 /* Constrain the initial coordinates and velocities */
526 do_constrain_first(fplog, constr, ir, mdatoms, state,
527 cr, nrnb, fr, top);
529 if (vsite)
531 /* Construct the virtual sites for the initial configuration */
532 construct_vsites(vsite, state->x, ir->delta_t, NULL,
533 top->idef.iparams, top->idef.il,
534 fr->ePBC, fr->bMolPBC, cr, state->box);
538 if (ir->efep != efepNO)
540 /* Set free energy calculation frequency as the greatest common
541 * denominator of nstdhdl and repl_ex_nst. */
542 nstfep = ir->fepvals->nstdhdl;
543 if (ir->bExpanded)
545 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
547 if (repl_ex_nst > 0)
549 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
553 /* Be REALLY careful about what flags you set here. You CANNOT assume
554 * this is the first step, since we might be restarting from a checkpoint,
555 * and in that case we should not do any modifications to the state.
557 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
559 if (Flags & MD_READ_EKIN)
561 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
564 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
565 | (bStopCM ? CGLO_STOPCM : 0)
566 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
567 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
568 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
570 bSumEkinhOld = FALSE;
571 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
572 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
573 constr, &nullSignaller, state->box,
574 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags
575 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
576 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
577 top_global, top, state,
578 &shouldCheckNumberOfBondedInteractions);
579 if (ir->eI == eiVVAK)
581 /* a second call to get the half step temperature initialized as well */
582 /* we do the same call as above, but turn the pressure off -- internally to
583 compute_globals, this is recognized as a velocity verlet half-step
584 kinetic energy calculation. This minimized excess variables, but
585 perhaps loses some logic?*/
587 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
588 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
589 constr, &nullSignaller, state->box,
590 NULL, &bSumEkinhOld,
591 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
594 /* Calculate the initial half step temperature, and save the ekinh_old */
595 if (!(Flags & MD_STARTFROMCPT))
597 for (i = 0; (i < ir->opts.ngtc); i++)
599 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
602 if (ir->eI != eiVV)
604 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
605 and there is no previous step */
608 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
609 temperature control */
610 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
612 if (MASTER(cr))
614 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
616 fprintf(fplog,
617 "RMS relative constraint deviation after constraining: %.2e\n",
618 constr_rmsd(constr));
620 if (EI_STATE_VELOCITY(ir->eI))
622 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
624 if (bRerunMD)
626 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
627 " input trajectory '%s'\n\n",
628 *(top_global->name), opt2fn("-rerun", nfile, fnm));
629 if (bVerbose)
631 fprintf(stderr, "Calculated time to finish depends on nsteps from "
632 "run input file,\nwhich may not correspond to the time "
633 "needed to process input trajectory.\n\n");
636 else
638 char tbuf[20];
639 fprintf(stderr, "starting mdrun '%s'\n",
640 *(top_global->name));
641 if (ir->nsteps >= 0)
643 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
645 else
647 sprintf(tbuf, "%s", "infinite");
649 if (ir->init_step > 0)
651 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
652 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
653 gmx_step_str(ir->init_step, sbuf2),
654 ir->init_step*ir->delta_t);
656 else
658 fprintf(stderr, "%s steps, %s ps.\n",
659 gmx_step_str(ir->nsteps, sbuf), tbuf);
662 fprintf(fplog, "\n");
665 walltime_accounting_start(walltime_accounting);
666 wallcycle_start(wcycle, ewcRUN);
667 print_start(fplog, cr, walltime_accounting, "mdrun");
669 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
670 #ifdef GMX_FAHCORE
671 chkpt_ret = fcCheckPointParallel( cr->nodeid,
672 NULL, 0);
673 if (chkpt_ret == 0)
675 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
677 #endif
679 /***********************************************************
681 * Loop over MD steps
683 ************************************************************/
685 /* if rerunMD then read coordinates and velocities from input trajectory */
686 if (bRerunMD)
688 if (getenv("GMX_FORCE_UPDATE"))
690 bForceUpdate = TRUE;
693 rerun_fr.natoms = 0;
694 if (MASTER(cr))
696 bLastStep = !read_first_frame(oenv, &status,
697 opt2fn("-rerun", nfile, fnm),
698 &rerun_fr, TRX_NEED_X | TRX_READ_V);
699 if (rerun_fr.natoms != top_global->natoms)
701 gmx_fatal(FARGS,
702 "Number of atoms in trajectory (%d) does not match the "
703 "run input file (%d)\n",
704 rerun_fr.natoms, top_global->natoms);
706 if (ir->ePBC != epbcNONE)
708 if (!rerun_fr.bBox)
710 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);
712 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
714 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
719 if (PAR(cr))
721 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
724 if (ir->ePBC != epbcNONE)
726 /* Set the shift vectors.
727 * Necessary here when have a static box different from the tpr box.
729 calc_shifts(rerun_fr.box, fr->shift_vec);
733 /* loop over MD steps or if rerunMD to end of input trajectory */
734 bFirstStep = TRUE;
735 /* Skip the first Nose-Hoover integration when we get the state from tpx */
736 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
737 bSumEkinhOld = FALSE;
738 bExchanged = FALSE;
739 bNeedRepartition = FALSE;
740 // TODO This implementation of ensemble orientation restraints is nasty because
741 // a user can't just do multi-sim with single-sim orientation restraints.
742 bUsingEnsembleRestraints = (fcd->disres.nsystems > 1) || (cr->ms && fcd->orires.nr);
745 // Replica exchange and ensemble restraints need all
746 // simulations to remain synchronized, so they need
747 // checkpoints and stop conditions to act on the same step, so
748 // the propagation of such signals must take place between
749 // simulations, not just within simulations.
750 bool checkpointIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
751 bool stopConditionIsLocal = (repl_ex_nst <= 0) && !bUsingEnsembleRestraints;
752 bool resetCountersIsLocal = true;
753 signals[eglsCHKPT] = SimulationSignal(checkpointIsLocal);
754 signals[eglsSTOPCOND] = SimulationSignal(stopConditionIsLocal);
755 signals[eglsRESETCOUNTERS] = SimulationSignal(resetCountersIsLocal);
758 step = ir->init_step;
759 step_rel = 0;
761 // TODO extract this to new multi-simulation module
762 if (MASTER(cr) && MULTISIM(cr) && (repl_ex_nst <= 0 ))
764 if (!multisim_int_all_are_equal(cr->ms, ir->nsteps))
766 GMX_LOG(mdlog.warning).appendText(
767 "Note: The number of steps is not consistent across multi simulations,\n"
768 "but we are proceeding anyway!");
770 if (!multisim_int_all_are_equal(cr->ms, ir->init_step))
772 GMX_LOG(mdlog.warning).appendText(
773 "Note: The initial step is not consistent across multi simulations,\n"
774 "but we are proceeding anyway!");
778 /* and stop now if we should */
779 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
780 while (!bLastStep)
783 /* Determine if this is a neighbor search step */
784 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
786 if (bPMETune && bNStList)
788 /* PME grid + cut-off optimization with GPUs or PME nodes */
789 pme_loadbal_do(pme_loadbal, cr,
790 (bVerbose && MASTER(cr)) ? stderr : NULL,
791 fplog, mdlog,
792 ir, fr, state,
793 wcycle,
794 step, step_rel,
795 &bPMETunePrinting);
798 wallcycle_start(wcycle, ewcSTEP);
800 if (bRerunMD)
802 if (rerun_fr.bStep)
804 step = rerun_fr.step;
805 step_rel = step - ir->init_step;
807 if (rerun_fr.bTime)
809 t = rerun_fr.time;
811 else
813 t = step;
816 else
818 bLastStep = (step_rel == ir->nsteps);
819 t = t0 + step*ir->delta_t;
822 // TODO Refactor this, so that nstfep does not need a default value of zero
823 if (ir->efep != efepNO || ir->bSimTemp)
825 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
826 requiring different logic. */
828 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
829 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
830 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
831 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
832 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
835 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
836 do_per_step(step, repl_ex_nst));
838 if (bSimAnn)
840 update_annealing_target_temp(ir, t, upd);
843 if (bRerunMD)
845 if (!DOMAINDECOMP(cr) || MASTER(cr))
847 for (i = 0; i < state_global->natoms; i++)
849 copy_rvec(rerun_fr.x[i], state_global->x[i]);
851 if (rerun_fr.bV)
853 for (i = 0; i < state_global->natoms; i++)
855 copy_rvec(rerun_fr.v[i], state_global->v[i]);
858 else
860 for (i = 0; i < state_global->natoms; i++)
862 clear_rvec(state_global->v[i]);
864 if (bRerunWarnNoV)
866 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
867 " Ekin, temperature and pressure are incorrect,\n"
868 " the virial will be incorrect when constraints are present.\n"
869 "\n");
870 bRerunWarnNoV = FALSE;
874 copy_mat(rerun_fr.box, state_global->box);
875 copy_mat(state_global->box, state->box);
877 if (vsite && (Flags & MD_RERUN_VSITE))
879 if (DOMAINDECOMP(cr))
881 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
883 if (graph)
885 /* Following is necessary because the graph may get out of sync
886 * with the coordinates if we only have every N'th coordinate set
888 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
889 shift_self(graph, state->box, state->x);
891 construct_vsites(vsite, state->x, ir->delta_t, state->v,
892 top->idef.iparams, top->idef.il,
893 fr->ePBC, fr->bMolPBC, cr, state->box);
894 if (graph)
896 unshift_self(graph, state->box, state->x);
901 /* Stop Center of Mass motion */
902 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
904 if (bRerunMD)
906 /* for rerun MD always do Neighbour Searching */
907 bNS = (bFirstStep || ir->nstlist != 0);
909 else
911 /* Determine whether or not to do Neighbour Searching */
912 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
915 /* < 0 means stop at next step, > 0 means stop at next NS step */
916 if ( (signals[eglsSTOPCOND].set < 0) ||
917 ( (signals[eglsSTOPCOND].set > 0 ) && ( bNS || ir->nstlist == 0)))
919 bLastStep = TRUE;
922 /* Determine whether or not to update the Born radii if doing GB */
923 bBornRadii = bFirstStep;
924 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
926 bBornRadii = TRUE;
929 /* do_log triggers energy and virial calculation. Because this leads
930 * to different code paths, forces can be different. Thus for exact
931 * continuation we should avoid extra log output.
932 * Note that the || bLastStep can result in non-exact continuation
933 * beyond the last step. But we don't consider that to be an issue.
935 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
936 do_verbose = bVerbose &&
937 (step % stepout == 0 || bFirstStep || bLastStep || bRerunMD);
939 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
941 if (bRerunMD)
943 bMasterState = TRUE;
945 else
947 bMasterState = FALSE;
948 /* Correct the new box if it is too skewed */
949 if (inputrecDynamicBox(ir))
951 if (correct_box(fplog, step, state->box, graph))
953 bMasterState = TRUE;
956 if (DOMAINDECOMP(cr) && bMasterState)
958 dd_collect_state(cr->dd, state, state_global);
962 if (DOMAINDECOMP(cr))
964 /* Repartition the domain decomposition */
965 dd_partition_system(fplog, step, cr,
966 bMasterState, nstglobalcomm,
967 state_global, top_global, ir,
968 state, &f, mdatoms, top, fr,
969 vsite, constr,
970 nrnb, wcycle,
971 do_verbose && !bPMETunePrinting);
972 shouldCheckNumberOfBondedInteractions = true;
973 update_realloc(upd, state->nalloc);
977 if (MASTER(cr) && do_log)
979 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
982 if (ir->efep != efepNO)
984 update_mdatoms(mdatoms, state->lambda[efptMASS]);
987 if ((bRerunMD && rerun_fr.bV) || bExchanged)
990 /* We need the kinetic energy at minus the half step for determining
991 * the full step kinetic energy and possibly for T-coupling.*/
992 /* This may not be quite working correctly yet . . . . */
993 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
994 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
995 constr, &nullSignaller, state->box,
996 &totalNumberOfBondedInteractions, &bSumEkinhOld,
997 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
998 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
999 top_global, top, state,
1000 &shouldCheckNumberOfBondedInteractions);
1002 clear_mat(force_vir);
1004 /* We write a checkpoint at this MD step when:
1005 * either at an NS step when we signalled through gs,
1006 * or at the last step (but not when we do not want confout),
1007 * but never at the first step or with rerun.
1009 bCPT = (((signals[eglsCHKPT].set && (bNS || ir->nstlist == 0)) ||
1010 (bLastStep && (Flags & MD_CONFOUT))) &&
1011 step > ir->init_step && !bRerunMD);
1012 if (bCPT)
1014 signals[eglsCHKPT].set = 0;
1017 /* Determine the energy and pressure:
1018 * at nstcalcenergy steps and at energy output steps (set below).
1020 if (EI_VV(ir->eI) && (!bInitStep))
1022 /* for vv, the first half of the integration actually corresponds
1023 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1024 but the virial needs to be calculated on both the current step and the 'next' step. Future
1025 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1027 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1028 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1029 bCalcVir = bCalcEnerStep ||
1030 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1032 else
1034 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1035 bCalcVir = bCalcEnerStep ||
1036 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1038 bCalcEner = bCalcEnerStep;
1040 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1042 if (do_ene || do_log || bDoReplEx)
1044 bCalcVir = TRUE;
1045 bCalcEner = TRUE;
1048 /* Do we need global communication ? */
1049 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1050 do_per_step(step, nstglobalcomm) ||
1051 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1053 force_flags = (GMX_FORCE_STATECHANGED |
1054 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1055 GMX_FORCE_ALLFORCES |
1056 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1057 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1058 (bDoFEP ? GMX_FORCE_DHDL : 0)
1061 if (shellfc)
1063 /* Now is the time to relax the shells */
1064 relax_shell_flexcon(fplog, cr, bVerbose, step,
1065 ir, bNS, force_flags, top,
1066 constr, enerd, fcd,
1067 state, f, force_vir, mdatoms,
1068 nrnb, wcycle, graph, groups,
1069 shellfc, fr, bBornRadii, t, mu_tot,
1070 vsite, mdoutf_get_fp_field(outf));
1072 else
1074 /* The coordinates (x) are shifted (to get whole molecules)
1075 * in do_force.
1076 * This is parallellized as well, and does communication too.
1077 * Check comments in sim_util.c
1079 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1080 state->box, state->x, &state->hist,
1081 f, force_vir, mdatoms, enerd, fcd,
1082 state->lambda, graph,
1083 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1084 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1087 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1088 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1090 rvec *vbuf = NULL;
1092 wallcycle_start(wcycle, ewcUPDATE);
1093 if (ir->eI == eiVV && bInitStep)
1095 /* if using velocity verlet with full time step Ekin,
1096 * take the first half step only to compute the
1097 * virial for the first step. From there,
1098 * revert back to the initial coordinates
1099 * so that the input is actually the initial step.
1101 snew(vbuf, state->natoms);
1102 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1104 else
1106 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1107 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1110 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1111 ekind, M, upd, etrtVELOCITY1,
1112 cr, constr);
1114 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1116 wallcycle_stop(wcycle, ewcUPDATE);
1117 update_constraints(fplog, step, NULL, ir, mdatoms,
1118 state, fr->bMolPBC, graph, f,
1119 &top->idef, shake_vir,
1120 cr, nrnb, wcycle, upd, constr,
1121 TRUE, bCalcVir);
1122 wallcycle_start(wcycle, ewcUPDATE);
1124 else if (graph)
1126 /* Need to unshift here if a do_force has been
1127 called in the previous step */
1128 unshift_self(graph, state->box, state->x);
1130 /* if VV, compute the pressure and constraints */
1131 /* For VV2, we strictly only need this if using pressure
1132 * control, but we really would like to have accurate pressures
1133 * printed out.
1134 * Think about ways around this in the future?
1135 * For now, keep this choice in comments.
1137 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1138 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1139 bPres = TRUE;
1140 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1141 if (bCalcEner && ir->eI == eiVVAK)
1143 bSumEkinhOld = TRUE;
1145 /* for vv, the first half of the integration actually corresponds to the previous step.
1146 So we need information from the last step in the first half of the integration */
1147 if (bGStat || do_per_step(step-1, nstglobalcomm))
1149 wallcycle_stop(wcycle, ewcUPDATE);
1150 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1151 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1152 constr, &nullSignaller, state->box,
1153 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1154 (bGStat ? CGLO_GSTAT : 0)
1155 | CGLO_ENERGY
1156 | (bTemp ? CGLO_TEMPERATURE : 0)
1157 | (bPres ? CGLO_PRESSURE : 0)
1158 | (bPres ? CGLO_CONSTRAINT : 0)
1159 | (bStopCM ? CGLO_STOPCM : 0)
1160 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1161 | CGLO_SCALEEKIN
1163 /* explanation of above:
1164 a) We compute Ekin at the full time step
1165 if 1) we are using the AveVel Ekin, and it's not the
1166 initial step, or 2) if we are using AveEkin, but need the full
1167 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1168 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1169 EkinAveVel because it's needed for the pressure */
1170 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1171 top_global, top, state,
1172 &shouldCheckNumberOfBondedInteractions);
1173 wallcycle_start(wcycle, ewcUPDATE);
1175 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1176 if (!bInitStep)
1178 if (bTrotter)
1180 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1181 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1183 copy_mat(shake_vir, state->svir_prev);
1184 copy_mat(force_vir, state->fvir_prev);
1185 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1187 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1188 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1189 enerd->term[F_EKIN] = trace(ekind->ekin);
1192 else if (bExchanged)
1194 wallcycle_stop(wcycle, ewcUPDATE);
1195 /* We need the kinetic energy at minus the half step for determining
1196 * the full step kinetic energy and possibly for T-coupling.*/
1197 /* This may not be quite working correctly yet . . . . */
1198 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1199 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1200 constr, &nullSignaller, state->box,
1201 NULL, &bSumEkinhOld,
1202 CGLO_GSTAT | CGLO_TEMPERATURE);
1203 wallcycle_start(wcycle, ewcUPDATE);
1206 /* if it's the initial step, we performed this first step just to get the constraint virial */
1207 if (ir->eI == eiVV && bInitStep)
1209 copy_rvecn(vbuf, state->v, 0, state->natoms);
1210 sfree(vbuf);
1212 wallcycle_stop(wcycle, ewcUPDATE);
1215 /* compute the conserved quantity */
1216 if (EI_VV(ir->eI))
1218 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1219 if (ir->eI == eiVV)
1221 last_ekin = enerd->term[F_EKIN];
1223 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1225 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1227 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1228 if (ir->efep != efepNO && !bRerunMD)
1230 sum_dhdl(enerd, state->lambda, ir->fepvals);
1234 /* ######## END FIRST UPDATE STEP ############## */
1235 /* ######## If doing VV, we now have v(dt) ###### */
1236 if (bDoExpanded)
1238 /* perform extended ensemble sampling in lambda - we don't
1239 actually move to the new state before outputting
1240 statistics, but if performing simulated tempering, we
1241 do update the velocities and the tau_t. */
1243 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, state->dfhist, step, state->v, mdatoms);
1244 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1245 copy_df_history(state_global->dfhist, state->dfhist);
1248 /* Now we have the energies and forces corresponding to the
1249 * coordinates at time t. We must output all of this before
1250 * the update.
1252 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1253 ir, state, state_global, top_global, fr,
1254 outf, mdebin, ekind, f,
1255 &nchkpt,
1256 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1257 bSumEkinhOld);
1258 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1259 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1261 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1262 if (startingFromCheckpoint && bTrotter)
1264 copy_mat(state->svir_prev, shake_vir);
1265 copy_mat(state->fvir_prev, force_vir);
1268 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1270 /* Check whether everything is still allright */
1271 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1272 #if GMX_THREAD_MPI
1273 && MASTER(cr)
1274 #endif
1277 int nsteps_stop = -1;
1279 /* this just makes signals[].sig compatible with the hack
1280 of sending signals around by MPI_Reduce together with
1281 other floats */
1282 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1284 signals[eglsSTOPCOND].sig = 1;
1285 nsteps_stop = std::max(ir->nstlist, 2*nstglobalcomm);
1287 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1289 signals[eglsSTOPCOND].sig = -1;
1290 nsteps_stop = nstglobalcomm + 1;
1292 if (fplog)
1294 fprintf(fplog,
1295 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1296 gmx_get_signal_name(), nsteps_stop);
1297 fflush(fplog);
1299 fprintf(stderr,
1300 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1301 gmx_get_signal_name(), nsteps_stop);
1302 fflush(stderr);
1303 handled_stop_condition = (int)gmx_get_stop_condition();
1305 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1306 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1307 signals[eglsSTOPCOND].sig == 0 && signals[eglsSTOPCOND].set == 0)
1309 /* Signal to terminate the run */
1310 signals[eglsSTOPCOND].sig = 1;
1311 if (fplog)
1313 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1315 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1318 if (bResetCountersHalfMaxH && MASTER(cr) &&
1319 elapsed_time > max_hours*60.0*60.0*0.495)
1321 /* Set flag that will communicate the signal to all ranks in the simulation */
1322 signals[eglsRESETCOUNTERS].sig = 1;
1325 /* In parallel we only have to check for checkpointing in steps
1326 * where we do global communication,
1327 * otherwise the other nodes don't know.
1329 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1330 cpt_period >= 0 &&
1331 (cpt_period == 0 ||
1332 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1333 signals[eglsCHKPT].set == 0)
1335 signals[eglsCHKPT].sig = 1;
1338 /* ######### START SECOND UPDATE STEP ################# */
1340 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1341 in preprocessing */
1343 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1345 gmx_bool bIfRandomize;
1346 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1347 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1348 if (constr && bIfRandomize)
1350 update_constraints(fplog, step, NULL, ir, mdatoms,
1351 state, fr->bMolPBC, graph, f,
1352 &top->idef, tmp_vir,
1353 cr, nrnb, wcycle, upd, constr,
1354 TRUE, bCalcVir);
1357 /* Box is changed in update() when we do pressure coupling,
1358 * but we should still use the old box for energy corrections and when
1359 * writing it to the energy file, so it matches the trajectory files for
1360 * the same timestep above. Make a copy in a separate array.
1362 copy_mat(state->box, lastbox);
1364 dvdl_constr = 0;
1366 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1368 wallcycle_start(wcycle, ewcUPDATE);
1369 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1370 if (bTrotter)
1372 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1373 /* We can only do Berendsen coupling after we have summed
1374 * the kinetic energy or virial. Since the happens
1375 * in global_state after update, we should only do it at
1376 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1379 else
1381 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1382 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1385 if (EI_VV(ir->eI))
1387 /* velocity half-step update */
1388 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1389 ekind, M, upd, etrtVELOCITY2,
1390 cr, constr);
1393 /* Above, initialize just copies ekinh into ekin,
1394 * it doesn't copy position (for VV),
1395 * and entire integrator for MD.
1398 if (ir->eI == eiVVAK)
1400 /* We probably only need md->homenr, not state->natoms */
1401 if (state->natoms > cbuf_nalloc)
1403 cbuf_nalloc = state->natoms;
1404 srenew(cbuf, cbuf_nalloc);
1406 copy_rvecn(state->x, cbuf, 0, state->natoms);
1409 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1410 ekind, M, upd, etrtPOSITION, cr, constr);
1411 wallcycle_stop(wcycle, ewcUPDATE);
1413 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1414 fr->bMolPBC, graph, f,
1415 &top->idef, shake_vir,
1416 cr, nrnb, wcycle, upd, constr,
1417 FALSE, bCalcVir);
1419 if (ir->eI == eiVVAK)
1421 /* erase F_EKIN and F_TEMP here? */
1422 /* just compute the kinetic energy at the half step to perform a trotter step */
1423 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1424 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1425 constr, &nullSignaller, lastbox,
1426 NULL, &bSumEkinhOld,
1427 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1429 wallcycle_start(wcycle, ewcUPDATE);
1430 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1431 /* now we know the scaling, we can compute the positions again again */
1432 copy_rvecn(cbuf, state->x, 0, state->natoms);
1434 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1435 ekind, M, upd, etrtPOSITION, cr, constr);
1436 wallcycle_stop(wcycle, ewcUPDATE);
1438 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1439 /* are the small terms in the shake_vir here due
1440 * to numerical errors, or are they important
1441 * physically? I'm thinking they are just errors, but not completely sure.
1442 * For now, will call without actually constraining, constr=NULL*/
1443 update_constraints(fplog, step, NULL, ir, mdatoms,
1444 state, fr->bMolPBC, graph, f,
1445 &top->idef, tmp_vir,
1446 cr, nrnb, wcycle, upd, NULL,
1447 FALSE, bCalcVir);
1449 if (EI_VV(ir->eI))
1451 /* this factor or 2 correction is necessary
1452 because half of the constraint force is removed
1453 in the vv step, so we have to double it. See
1454 the Redmine issue #1255. It is not yet clear
1455 if the factor of 2 is exact, or just a very
1456 good approximation, and this will be
1457 investigated. The next step is to see if this
1458 can be done adding a dhdl contribution from the
1459 rattle step, but this is somewhat more
1460 complicated with the current code. Will be
1461 investigated, hopefully for 4.6.3. However,
1462 this current solution is much better than
1463 having it completely wrong.
1465 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1467 else
1469 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1472 else if (graph)
1474 /* Need to unshift here */
1475 unshift_self(graph, state->box, state->x);
1478 if (vsite != NULL)
1480 wallcycle_start(wcycle, ewcVSITECONSTR);
1481 if (graph != NULL)
1483 shift_self(graph, state->box, state->x);
1485 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1486 top->idef.iparams, top->idef.il,
1487 fr->ePBC, fr->bMolPBC, cr, state->box);
1489 if (graph != NULL)
1491 unshift_self(graph, state->box, state->x);
1493 wallcycle_stop(wcycle, ewcVSITECONSTR);
1496 /* ############## IF NOT VV, Calculate globals HERE ############ */
1497 /* With Leap-Frog we can skip compute_globals at
1498 * non-communication steps, but we need to calculate
1499 * the kinetic energy one step before communication.
1502 // Organize to do inter-simulation signalling on steps if
1503 // and when algorithms require it.
1504 bool doInterSimSignal = (!bFirstStep && bDoReplEx) || bUsingEnsembleRestraints;
1506 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)) || doInterSimSignal)
1508 // Since we're already communicating at this step, we
1509 // can propagate intra-simulation signals. Note that
1510 // check_nstglobalcomm has the responsibility for
1511 // choosing the value of nstglobalcomm that is one way
1512 // bGStat becomes true, so we can't get into a
1513 // situation where e.g. checkpointing can't be
1514 // signalled.
1515 bool doIntraSimSignal = true;
1516 SimulationSignaller signaller(&signals, cr, doInterSimSignal, doIntraSimSignal);
1518 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1519 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1520 constr, &signaller,
1521 lastbox,
1522 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1523 (bGStat ? CGLO_GSTAT : 0)
1524 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1525 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1526 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1527 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1528 | CGLO_CONSTRAINT
1529 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1531 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1532 top_global, top, state,
1533 &shouldCheckNumberOfBondedInteractions);
1537 /* ############# END CALC EKIN AND PRESSURE ################# */
1539 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1540 the virial that should probably be addressed eventually. state->veta has better properies,
1541 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1542 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1544 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1546 /* Sum up the foreign energy and dhdl terms for md and sd.
1547 Currently done every step so that dhdl is correct in the .edr */
1548 sum_dhdl(enerd, state->lambda, ir->fepvals);
1550 update_box(fplog, step, ir, mdatoms, state, f,
1551 pcoupl_mu, nrnb, upd);
1553 /* ################# END UPDATE STEP 2 ################# */
1554 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1556 /* The coordinates (x) were unshifted in update */
1557 if (!bGStat)
1559 /* We will not sum ekinh_old,
1560 * so signal that we still have to do it.
1562 bSumEkinhOld = TRUE;
1565 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1567 /* use the directly determined last velocity, not actually the averaged half steps */
1568 if (bTrotter && ir->eI == eiVV)
1570 enerd->term[F_EKIN] = last_ekin;
1572 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1574 if (EI_VV(ir->eI))
1576 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1578 else
1580 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1582 /* ######### END PREPARING EDR OUTPUT ########### */
1584 /* Output stuff */
1585 if (MASTER(cr))
1587 if (fplog && do_log && bDoExpanded)
1589 /* only needed if doing expanded ensemble */
1590 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1591 state_global->dfhist, state->fep_state, ir->nstlog, step);
1593 if (bCalcEner)
1595 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1596 t, mdatoms->tmass, enerd, state,
1597 ir->fepvals, ir->expandedvals, lastbox,
1598 shake_vir, force_vir, total_vir, pres,
1599 ekind, mu_tot, constr);
1601 else
1603 upd_mdebin_step(mdebin);
1606 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1607 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1609 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1610 step, t,
1611 eprNORMAL, mdebin, fcd, groups, &(ir->opts));
1613 if (ir->bPull)
1615 pull_print_output(ir->pull_work, step, t);
1618 if (do_per_step(step, ir->nstlog))
1620 if (fflush(fplog) != 0)
1622 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1626 if (bDoExpanded)
1628 /* Have to do this part _after_ outputting the logfile and the edr file */
1629 /* Gets written into the state at the beginning of next loop*/
1630 state->fep_state = lamnew;
1632 /* Print the remaining wall clock time for the run */
1633 if (MULTIMASTER(cr) &&
1634 (do_verbose || gmx_got_usr_signal()) &&
1635 !bPMETunePrinting)
1637 if (shellfc)
1639 fprintf(stderr, "\n");
1641 print_time(stderr, walltime_accounting, step, ir, cr);
1644 /* Ion/water position swapping.
1645 * Not done in last step since trajectory writing happens before this call
1646 * in the MD loop and exchanges would be lost anyway. */
1647 bNeedRepartition = FALSE;
1648 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1649 do_per_step(step, ir->swap->nstswap))
1651 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1652 bRerunMD ? rerun_fr.x : state->x,
1653 bRerunMD ? rerun_fr.box : state->box,
1654 top_global, MASTER(cr) && bVerbose, bRerunMD);
1656 if (bNeedRepartition && DOMAINDECOMP(cr))
1658 dd_collect_state(cr->dd, state, state_global);
1662 /* Replica exchange */
1663 bExchanged = FALSE;
1664 if (bDoReplEx)
1666 bExchanged = replica_exchange(fplog, cr, repl_ex,
1667 state_global, enerd,
1668 state, step, t);
1671 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1673 dd_partition_system(fplog, step, cr, TRUE, 1,
1674 state_global, top_global, ir,
1675 state, &f, mdatoms, top, fr,
1676 vsite, constr,
1677 nrnb, wcycle, FALSE);
1678 shouldCheckNumberOfBondedInteractions = true;
1679 update_realloc(upd, state->nalloc);
1682 bFirstStep = FALSE;
1683 bInitStep = FALSE;
1684 startingFromCheckpoint = FALSE;
1686 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1687 /* With all integrators, except VV, we need to retain the pressure
1688 * at the current step for coupling at the next step.
1690 if ((state->flags & (1<<estPRES_PREV)) &&
1691 (bGStatEveryStep ||
1692 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1694 /* Store the pressure in t_state for pressure coupling
1695 * at the next MD step.
1697 copy_mat(pres, state->pres_prev);
1700 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1702 if ( (membed != NULL) && (!bLastStep) )
1704 rescale_membed(step_rel, membed, state_global->x);
1707 if (bRerunMD)
1709 if (MASTER(cr))
1711 /* read next frame from input trajectory */
1712 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1715 if (PAR(cr))
1717 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1721 cycles = wallcycle_stop(wcycle, ewcSTEP);
1722 if (DOMAINDECOMP(cr) && wcycle)
1724 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1727 if (!bRerunMD || !rerun_fr.bStep)
1729 /* increase the MD step number */
1730 step++;
1731 step_rel++;
1734 /* TODO make a counter-reset module */
1735 /* If it is time to reset counters, set a flag that remains
1736 true until counters actually get reset */
1737 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1738 signals[eglsRESETCOUNTERS].set != 0)
1740 if (pme_loadbal_is_active(pme_loadbal))
1742 /* Do not permit counter reset while PME load
1743 * balancing is active. The only purpose for resetting
1744 * counters is to measure reliable performance data,
1745 * and that can't be done before balancing
1746 * completes.
1748 * TODO consider fixing this by delaying the reset
1749 * until after load balancing completes,
1750 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1751 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1752 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1753 "resetting counters later in the run, e.g. with gmx "
1754 "mdrun -resetstep.", step);
1756 reset_all_counters(fplog, mdlog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1757 use_GPU(fr->nbv) ? fr->nbv : NULL);
1758 wcycle_set_reset_counters(wcycle, -1);
1759 if (!(cr->duty & DUTY_PME))
1761 /* Tell our PME node to reset its counters */
1762 gmx_pme_send_resetcounters(cr, step);
1764 /* Correct max_hours for the elapsed time */
1765 max_hours -= elapsed_time/(60.0*60.0);
1766 /* If mdrun -maxh -resethway was active, it can only trigger once */
1767 bResetCountersHalfMaxH = FALSE; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1768 /* Reset can only happen once, so clear the triggering flag. */
1769 signals[eglsRESETCOUNTERS].set = 0;
1772 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1773 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1776 /* End of main MD loop */
1778 /* Closing TNG files can include compressing data. Therefore it is good to do that
1779 * before stopping the time measurements. */
1780 mdoutf_tng_close(outf);
1782 /* Stop measuring walltime */
1783 walltime_accounting_end(walltime_accounting);
1785 if (bRerunMD && MASTER(cr))
1787 close_trj(status);
1790 if (!(cr->duty & DUTY_PME))
1792 /* Tell the PME only node to finish */
1793 gmx_pme_send_finish(cr);
1796 if (MASTER(cr))
1798 if (ir->nstcalcenergy > 0 && !bRerunMD)
1800 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1801 eprAVER, mdebin, fcd, groups, &(ir->opts));
1805 done_mdoutf(outf);
1807 if (bPMETune)
1809 pme_loadbal_done(pme_loadbal, fplog, mdlog, use_GPU(fr->nbv));
1812 done_shellfc(fplog, shellfc, step_rel);
1814 if (repl_ex_nst > 0 && MASTER(cr))
1816 print_replica_exchange_statistics(fplog, repl_ex);
1819 // Clean up swapcoords
1820 if (ir->eSwapCoords != eswapNO)
1822 finish_swapcoords(ir->swap);
1825 /* IMD cleanup, if bIMD is TRUE. */
1826 IMD_finalize(ir->bIMD, ir->imd);
1828 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1830 return 0;