Fix membed initialization
[gromacs.git] / src / programs / mdrun / md.cpp
blobbf573e33a3f63f7fa2de49fb59e2e003acb633a8
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 "thread_mpi/threads.h"
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_network.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/ewald/pme.h"
54 #include "gromacs/ewald/pme-load-balancing.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/gmxlib/md_logging.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nrnb.h"
59 #include "gromacs/imd/imd.h"
60 #include "gromacs/listed-forces/manage-threading.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/mdlib/compute_io.h"
66 #include "gromacs/mdlib/constr.h"
67 #include "gromacs/mdlib/ebin.h"
68 #include "gromacs/mdlib/force.h"
69 #include "gromacs/mdlib/force_flags.h"
70 #include "gromacs/mdlib/forcerec.h"
71 #include "gromacs/mdlib/md_support.h"
72 #include "gromacs/mdlib/mdatoms.h"
73 #include "gromacs/mdlib/mdebin.h"
74 #include "gromacs/mdlib/mdoutf.h"
75 #include "gromacs/mdlib/mdrun.h"
76 #include "gromacs/mdlib/mdrun_signalling.h"
77 #include "gromacs/mdlib/nb_verlet.h"
78 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
79 #include "gromacs/mdlib/ns.h"
80 #include "gromacs/mdlib/shellfc.h"
81 #include "gromacs/mdlib/sighandler.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/tgroup.h"
84 #include "gromacs/mdlib/trajectory_writing.h"
85 #include "gromacs/mdlib/update.h"
86 #include "gromacs/mdlib/vcm.h"
87 #include "gromacs/mdlib/vsite.h"
88 #include "gromacs/mdtypes/commrec.h"
89 #include "gromacs/mdtypes/df_history.h"
90 #include "gromacs/mdtypes/energyhistory.h"
91 #include "gromacs/mdtypes/fcdata.h"
92 #include "gromacs/mdtypes/forcerec.h"
93 #include "gromacs/mdtypes/group.h"
94 #include "gromacs/mdtypes/inputrec.h"
95 #include "gromacs/mdtypes/interaction_const.h"
96 #include "gromacs/mdtypes/md_enums.h"
97 #include "gromacs/mdtypes/mdatom.h"
98 #include "gromacs/mdtypes/state.h"
99 #include "gromacs/pbcutil/mshift.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/pulling/pull_rotation.h"
103 #include "gromacs/swap/swapcoords.h"
104 #include "gromacs/timing/wallcycle.h"
105 #include "gromacs/timing/walltime_accounting.h"
106 #include "gromacs/topology/atoms.h"
107 #include "gromacs/topology/idef.h"
108 #include "gromacs/topology/mtop_util.h"
109 #include "gromacs/topology/topology.h"
110 #include "gromacs/trajectory/trajectoryframe.h"
111 #include "gromacs/utility/basedefinitions.h"
112 #include "gromacs/utility/cstringutil.h"
113 #include "gromacs/utility/fatalerror.h"
114 #include "gromacs/utility/real.h"
115 #include "gromacs/utility/smalloc.h"
117 #include "deform.h"
118 #include "membed.h"
119 #include "repl_ex.h"
121 #ifdef GMX_FAHCORE
122 #include "corewrap.h"
123 #endif
125 /*! \brief Check whether bonded interactions are missing, if appropriate
127 * \param[in] fplog Log file pointer
128 * \param[in] cr Communication object
129 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
130 * \param[in] top_global Global topology for the error message
131 * \param[in] top_local Local topology for the error message
132 * \param[in] state Global state for the error message
133 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
135 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
136 * is always set to false after exit.
138 static void checkNumberOfBondedInteractions(FILE *fplog, t_commrec *cr, int totalNumberOfBondedInteractions,
139 gmx_mtop_t *top_global, gmx_localtop_t *top_local, t_state *state,
140 bool *shouldCheckNumberOfBondedInteractions)
142 if (*shouldCheckNumberOfBondedInteractions)
144 if (totalNumberOfBondedInteractions != cr->dd->nbonded_global)
146 dd_print_missing_interactions(fplog, cr, totalNumberOfBondedInteractions, top_global, top_local, state); // Does not return
148 *shouldCheckNumberOfBondedInteractions = false;
152 static void reset_all_counters(FILE *fplog, t_commrec *cr,
153 gmx_int64_t step,
154 gmx_int64_t *step_rel, t_inputrec *ir,
155 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
156 gmx_walltime_accounting_t walltime_accounting,
157 struct nonbonded_verlet_t *nbv)
159 char sbuf[STEPSTRSIZE];
161 /* Reset all the counters related to performance over the run */
162 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
163 gmx_step_str(step, sbuf));
165 if (use_GPU(nbv))
167 nbnxn_gpu_reset_timings(nbv);
170 wallcycle_stop(wcycle, ewcRUN);
171 wallcycle_reset_all(wcycle);
172 if (DOMAINDECOMP(cr))
174 reset_dd_statistics_counters(cr->dd);
176 init_nrnb(nrnb);
177 ir->init_step += *step_rel;
178 ir->nsteps -= *step_rel;
179 *step_rel = 0;
180 wallcycle_start(wcycle, ewcRUN);
181 walltime_accounting_start(walltime_accounting);
182 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
185 /*! \libinternal
186 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
187 int nfile, const t_filenm fnm[],
188 const gmx_output_env_t *oenv, gmx_bool bVerbose,
189 int nstglobalcomm,
190 gmx_vsite_t *vsite, gmx_constr_t constr,
191 int stepout,
192 t_inputrec *inputrec,
193 gmx_mtop_t *top_global, t_fcdata *fcd,
194 t_state *state_global,
195 t_mdatoms *mdatoms,
196 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
197 gmx_edsam_t ed,
198 t_forcerec *fr,
199 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
200 real cpt_period, real max_hours,
201 int imdport,
202 unsigned long Flags,
203 gmx_walltime_accounting_t walltime_accounting)
205 double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
206 const gmx_output_env_t *oenv, gmx_bool bVerbose,
207 int nstglobalcomm,
208 gmx_vsite_t *vsite, gmx_constr_t constr,
209 int stepout, t_inputrec *ir,
210 gmx_mtop_t *top_global,
211 t_fcdata *fcd,
212 t_state *state_global,
213 t_mdatoms *mdatoms,
214 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
215 gmx_edsam_t ed, t_forcerec *fr,
216 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
217 real cpt_period, real max_hours,
218 int imdport,
219 unsigned long Flags,
220 gmx_walltime_accounting_t walltime_accounting)
222 gmx_mdoutf_t outf = NULL;
223 gmx_int64_t step, step_rel;
224 double elapsed_time;
225 double t, t0, lam0[efptNR];
226 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEnerStep, bCalcEner;
227 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD,
228 bFirstStep, startingFromCheckpoint, bInitStep, bLastStep = FALSE,
229 bBornRadii;
230 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
231 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
232 bForceUpdate = FALSE, bCPT;
233 gmx_bool bMasterState;
234 int force_flags, cglo_flags;
235 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
236 int i, m;
237 t_trxstatus *status;
238 rvec mu_tot;
239 t_vcm *vcm;
240 matrix pcoupl_mu, M;
241 t_trxframe rerun_fr;
242 gmx_repl_ex_t repl_ex = NULL;
243 int nchkpt = 1;
244 gmx_localtop_t *top;
245 t_mdebin *mdebin = NULL;
246 t_state *state = NULL;
247 gmx_enerdata_t *enerd;
248 rvec *f = NULL;
249 gmx_global_stat_t gstat;
250 gmx_update_t *upd = NULL;
251 t_graph *graph = NULL;
252 gmx_signalling_t gs;
253 gmx_groups_t *groups;
254 gmx_ekindata_t *ekind;
255 gmx_shellfc_t *shellfc;
256 gmx_bool bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
257 gmx_bool bResetCountersHalfMaxH = FALSE;
258 gmx_bool bTemp, bPres, bTrotter;
259 real dvdl_constr;
260 rvec *cbuf = NULL;
261 int cbuf_nalloc = 0;
262 matrix lastbox;
263 int lamnew = 0;
264 /* for FEP */
265 int nstfep = 0;
266 double cycles;
267 real saved_conserved_quantity = 0;
268 real last_ekin = 0;
269 t_extmass MassQ;
270 int **trotter_seq;
271 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
272 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
273 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
274 simulation stops. If equal to zero, don't
275 communicate any more between multisims.*/
276 /* PME load balancing data for GPU kernels */
277 pme_load_balancing_t *pme_loadbal = NULL;
278 gmx_bool bPMETune = FALSE;
279 gmx_bool bPMETunePrinting = FALSE;
281 /* Interactive MD */
282 gmx_bool bIMDstep = FALSE;
283 gmx_membed_t *membed = NULL;
285 #ifdef GMX_FAHCORE
286 /* Temporary addition for FAHCORE checkpointing */
287 int chkpt_ret;
288 #endif
289 /* Domain decomposition could incorrectly miss a bonded
290 interaction, but checking for that requires a global
291 communication stage, which does not otherwise happen in DD
292 code. So we do that alongside the first global energy reduction
293 after a new DD is made. These variables handle whether the
294 check happens, and the result it returns. */
295 bool shouldCheckNumberOfBondedInteractions = false;
296 int totalNumberOfBondedInteractions = -1;
298 /* Check for special mdrun options */
299 bRerunMD = (Flags & MD_RERUN);
300 if (Flags & MD_RESETCOUNTERSHALFWAY)
302 if (ir->nsteps > 0)
304 /* Signal to reset the counters half the simulation steps. */
305 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
307 /* Signal to reset the counters halfway the simulation time. */
308 bResetCountersHalfMaxH = (max_hours > 0);
311 /* md-vv uses averaged full step velocities for T-control
312 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
313 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
314 bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
316 if (bRerunMD)
318 /* Since we don't know if the frames read are related in any way,
319 * rebuild the neighborlist at every step.
321 ir->nstlist = 1;
322 ir->nstcalcenergy = 1;
323 nstglobalcomm = 1;
326 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
327 bGStatEveryStep = (nstglobalcomm == 1);
329 if (bRerunMD)
331 ir->nstxout_compressed = 0;
333 groups = &top_global->groups;
335 if (opt2bSet("-membed", nfile, fnm))
337 if (MASTER(cr))
339 fprintf(stderr, "Initializing membed");
341 /* Note that membed cannot work in parallel because mtop is
342 * changed here. Fix this if we ever want to make it run with
343 * multiple ranks. */
344 membed = init_membed(fplog, nfile, fnm, top_global, ir, state_global, cr, &cpt_period);
346 if (ir->bPull)
348 /* Initialize pull code */
349 ir->pull_work = init_pull(fplog, ir->pull, ir, nfile, fnm,
350 top_global, cr, oenv, ir->fepvals->init_lambda,
351 EI_DYNAMICS(ir->eI) && MASTER(cr), Flags);
354 if (ir->bRot)
356 /* Initialize enforced rotation code */
357 init_rot(fplog, ir, nfile, fnm, cr, state_global->x, state_global->box, top_global, oenv,
358 bVerbose, Flags);
361 if (ir->eSwapCoords != eswapNO)
363 /* Initialize ion swapping code */
364 init_swapcoords(fplog, bVerbose, ir, opt2fn_master("-swap", nfile, fnm, cr),
365 top_global, state_global->x, state_global->box, &state_global->swapstate, cr, oenv, Flags);
368 /* Initial values */
369 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
370 &(state_global->fep_state), lam0,
371 nrnb, top_global, &upd,
372 nfile, fnm, &outf, &mdebin,
373 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
375 clear_mat(total_vir);
376 clear_mat(pres);
377 /* Energy terms and groups */
378 snew(enerd, 1);
379 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
380 enerd);
381 if (DOMAINDECOMP(cr))
383 f = NULL;
385 else
387 snew(f, top_global->natoms);
390 /* Kinetic energy data */
391 snew(ekind, 1);
392 init_ekindata(fplog, top_global, &(ir->opts), ekind);
393 /* Copy the cos acceleration to the groups struct */
394 ekind->cosacc.cos_accel = ir->cos_accel;
396 gstat = global_stat_init(ir);
398 /* Check for polarizable models and flexible constraints */
399 shellfc = init_shell_flexcon(fplog,
400 top_global, n_flexible_constraints(constr),
401 ir->nstcalcenergy, DOMAINDECOMP(cr));
403 if (shellfc && ir->nstcalcenergy != 1)
405 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);
407 if (shellfc && DOMAINDECOMP(cr))
409 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
412 if (inputrecDeform(ir))
414 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
415 set_deform_reference_box(upd,
416 deform_init_init_step_tpx,
417 deform_init_box_tpx);
418 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
422 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
423 if ((io > 2000) && MASTER(cr))
425 fprintf(stderr,
426 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
427 io);
431 if (DOMAINDECOMP(cr))
433 top = dd_init_local_top(top_global);
435 snew(state, 1);
436 dd_init_local_state(cr->dd, state_global, state);
438 else
440 top = gmx_mtop_generate_local_top(top_global, ir->efep != efepNO);
442 forcerec_set_excl_load(fr, top);
444 state = serial_init_local_state(state_global);
446 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
448 if (vsite)
450 set_vsite_top(vsite, top, mdatoms, cr);
453 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
455 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
458 if (shellfc)
460 make_local_shells(cr, mdatoms, shellfc);
463 setup_bonded_threading(fr, &top->idef);
465 update_realloc(upd, state->nalloc);
468 /* Set up interactive MD (IMD) */
469 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
470 nfile, fnm, oenv, imdport, Flags);
472 if (DOMAINDECOMP(cr))
474 /* Distribute the charge groups over the nodes from the master node */
475 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
476 state_global, top_global, ir,
477 state, &f, mdatoms, top, fr,
478 vsite, constr,
479 nrnb, NULL, FALSE);
480 shouldCheckNumberOfBondedInteractions = true;
481 update_realloc(upd, state->nalloc);
484 update_mdatoms(mdatoms, state->lambda[efptMASS]);
486 startingFromCheckpoint = Flags & MD_STARTFROMCPT;
488 if (ir->bExpanded)
490 init_expanded_ensemble(startingFromCheckpoint, ir, &state->dfhist);
493 if (MASTER(cr))
495 if (startingFromCheckpoint)
497 /* Update mdebin with energy history if appending to output files */
498 if (Flags & MD_APPENDFILES)
500 restore_energyhistory_from_state(mdebin, state_global->enerhist);
502 else
504 /* We might have read an energy history from checkpoint,
505 * free the allocated memory and reset the counts.
507 done_energyhistory(state_global->enerhist);
508 init_energyhistory(state_global->enerhist);
511 /* Set the initial energy history in state by updating once */
512 update_energyhistory(state_global->enerhist, mdebin);
515 /* Initialize constraints */
516 if (constr && !DOMAINDECOMP(cr))
518 set_constraints(constr, top, ir, mdatoms, cr);
521 if (repl_ex_nst > 0 && MASTER(cr))
523 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
524 repl_ex_nst, repl_ex_nex, repl_ex_seed);
527 /* PME tuning is only supported with PME for Coulomb. Is is not supported
528 * with only LJ PME, or for reruns.
530 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
531 !(Flags & MD_REPRODUCIBLE));
532 if (bPMETune)
534 pme_loadbal_init(&pme_loadbal, cr, fplog, ir, state->box,
535 fr->ic, fr->pmedata, use_GPU(fr->nbv),
536 &bPMETunePrinting);
539 if (!ir->bContinuation && !bRerunMD)
541 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
543 /* Set the velocities of frozen particles to zero */
544 for (i = 0; i < mdatoms->homenr; i++)
546 for (m = 0; m < DIM; m++)
548 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
550 state->v[i][m] = 0;
556 if (constr)
558 /* Constrain the initial coordinates and velocities */
559 do_constrain_first(fplog, constr, ir, mdatoms, state,
560 cr, nrnb, fr, top);
562 if (vsite)
564 /* Construct the virtual sites for the initial configuration */
565 construct_vsites(vsite, state->x, ir->delta_t, NULL,
566 top->idef.iparams, top->idef.il,
567 fr->ePBC, fr->bMolPBC, cr, state->box);
571 if (ir->efep != efepNO)
573 /* Set free energy calculation frequency as the greatest common
574 * denominator of nstdhdl and repl_ex_nst. */
575 nstfep = ir->fepvals->nstdhdl;
576 if (ir->bExpanded)
578 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
580 if (repl_ex_nst > 0)
582 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
586 /* Be REALLY careful about what flags you set here. You CANNOT assume
587 * this is the first step, since we might be restarting from a checkpoint,
588 * and in that case we should not do any modifications to the state.
590 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
592 if (Flags & MD_READ_EKIN)
594 restore_ekinstate_from_state(cr, ekind, &state_global->ekinstate);
597 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
598 | (bStopCM ? CGLO_STOPCM : 0)
599 | (EI_VV(ir->eI) ? CGLO_PRESSURE : 0)
600 | (EI_VV(ir->eI) ? CGLO_CONSTRAINT : 0)
601 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
603 bSumEkinhOld = FALSE;
604 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
605 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
606 constr, NULL, FALSE, state->box,
607 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags
608 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
609 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
610 top_global, top, state,
611 &shouldCheckNumberOfBondedInteractions);
612 if (ir->eI == eiVVAK)
614 /* a second call to get the half step temperature initialized as well */
615 /* we do the same call as above, but turn the pressure off -- internally to
616 compute_globals, this is recognized as a velocity verlet half-step
617 kinetic energy calculation. This minimized excess variables, but
618 perhaps loses some logic?*/
620 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
621 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
622 constr, NULL, FALSE, state->box,
623 NULL, &bSumEkinhOld,
624 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
627 /* Calculate the initial half step temperature, and save the ekinh_old */
628 if (!(Flags & MD_STARTFROMCPT))
630 for (i = 0; (i < ir->opts.ngtc); i++)
632 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
635 if (ir->eI != eiVV)
637 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
638 and there is no previous step */
641 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
642 temperature control */
643 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
645 if (MASTER(cr))
647 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
649 fprintf(fplog,
650 "RMS relative constraint deviation after constraining: %.2e\n",
651 constr_rmsd(constr));
653 if (EI_STATE_VELOCITY(ir->eI))
655 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
657 if (bRerunMD)
659 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
660 " input trajectory '%s'\n\n",
661 *(top_global->name), opt2fn("-rerun", nfile, fnm));
662 if (bVerbose)
664 fprintf(stderr, "Calculated time to finish depends on nsteps from "
665 "run input file,\nwhich may not correspond to the time "
666 "needed to process input trajectory.\n\n");
669 else
671 char tbuf[20];
672 fprintf(stderr, "starting mdrun '%s'\n",
673 *(top_global->name));
674 if (ir->nsteps >= 0)
676 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
678 else
680 sprintf(tbuf, "%s", "infinite");
682 if (ir->init_step > 0)
684 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
685 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
686 gmx_step_str(ir->init_step, sbuf2),
687 ir->init_step*ir->delta_t);
689 else
691 fprintf(stderr, "%s steps, %s ps.\n",
692 gmx_step_str(ir->nsteps, sbuf), tbuf);
695 fprintf(fplog, "\n");
698 walltime_accounting_start(walltime_accounting);
699 wallcycle_start(wcycle, ewcRUN);
700 print_start(fplog, cr, walltime_accounting, "mdrun");
702 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
703 #ifdef GMX_FAHCORE
704 chkpt_ret = fcCheckPointParallel( cr->nodeid,
705 NULL, 0);
706 if (chkpt_ret == 0)
708 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
710 #endif
712 /***********************************************************
714 * Loop over MD steps
716 ************************************************************/
718 /* if rerunMD then read coordinates and velocities from input trajectory */
719 if (bRerunMD)
721 if (getenv("GMX_FORCE_UPDATE"))
723 bForceUpdate = TRUE;
726 rerun_fr.natoms = 0;
727 if (MASTER(cr))
729 bLastStep = !read_first_frame(oenv, &status,
730 opt2fn("-rerun", nfile, fnm),
731 &rerun_fr, TRX_NEED_X | TRX_READ_V);
732 if (rerun_fr.natoms != top_global->natoms)
734 gmx_fatal(FARGS,
735 "Number of atoms in trajectory (%d) does not match the "
736 "run input file (%d)\n",
737 rerun_fr.natoms, top_global->natoms);
739 if (ir->ePBC != epbcNONE)
741 if (!rerun_fr.bBox)
743 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);
745 if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist))
747 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
752 if (PAR(cr))
754 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
757 if (ir->ePBC != epbcNONE)
759 /* Set the shift vectors.
760 * Necessary here when have a static box different from the tpr box.
762 calc_shifts(rerun_fr.box, fr->shift_vec);
766 /* loop over MD steps or if rerunMD to end of input trajectory */
767 bFirstStep = TRUE;
768 /* Skip the first Nose-Hoover integration when we get the state from tpx */
769 bInitStep = !startingFromCheckpoint || EI_VV(ir->eI);
770 bSumEkinhOld = FALSE;
771 bExchanged = FALSE;
772 bNeedRepartition = FALSE;
774 init_global_signals(&gs, cr, ir, repl_ex_nst);
776 step = ir->init_step;
777 step_rel = 0;
779 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
781 /* check how many steps are left in other sims */
782 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
786 /* and stop now if we should */
787 bLastStep = (bLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
788 ((multisim_nsteps >= 0) && (step_rel >= multisim_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 : NULL,
800 fplog,
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 if (ir->efep != efepNO || ir->bSimTemp)
833 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
834 requiring different logic. */
836 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
837 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
838 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
839 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
840 && (ir->bExpanded) && (step > 0) && (!startingFromCheckpoint));
843 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
844 do_per_step(step, repl_ex_nst));
846 if (bSimAnn)
848 update_annealing_target_temp(ir, t, upd);
851 if (bRerunMD)
853 if (!DOMAINDECOMP(cr) || MASTER(cr))
855 for (i = 0; i < state_global->natoms; i++)
857 copy_rvec(rerun_fr.x[i], state_global->x[i]);
859 if (rerun_fr.bV)
861 for (i = 0; i < state_global->natoms; i++)
863 copy_rvec(rerun_fr.v[i], state_global->v[i]);
866 else
868 for (i = 0; i < state_global->natoms; i++)
870 clear_rvec(state_global->v[i]);
872 if (bRerunWarnNoV)
874 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
875 " Ekin, temperature and pressure are incorrect,\n"
876 " the virial will be incorrect when constraints are present.\n"
877 "\n");
878 bRerunWarnNoV = FALSE;
882 copy_mat(rerun_fr.box, state_global->box);
883 copy_mat(state_global->box, state->box);
885 if (vsite && (Flags & MD_RERUN_VSITE))
887 if (DOMAINDECOMP(cr))
889 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
891 if (graph)
893 /* Following is necessary because the graph may get out of sync
894 * with the coordinates if we only have every N'th coordinate set
896 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
897 shift_self(graph, state->box, state->x);
899 construct_vsites(vsite, state->x, ir->delta_t, state->v,
900 top->idef.iparams, top->idef.il,
901 fr->ePBC, fr->bMolPBC, cr, state->box);
902 if (graph)
904 unshift_self(graph, state->box, state->x);
909 /* Stop Center of Mass motion */
910 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
912 if (bRerunMD)
914 /* for rerun MD always do Neighbour Searching */
915 bNS = (bFirstStep || ir->nstlist != 0);
916 bNStList = bNS;
918 else
920 /* Determine whether or not to do Neighbour Searching */
921 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
924 /* check whether we should stop because another simulation has
925 stopped. */
926 if (MULTISIM(cr))
928 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
929 (multisim_nsteps != ir->nsteps) )
931 if (bNS)
933 if (MASTER(cr))
935 fprintf(stderr,
936 "Stopping simulation %d because another one has finished\n",
937 cr->ms->sim);
939 bLastStep = TRUE;
940 gs.sig[eglsCHKPT] = 1;
945 /* < 0 means stop at next step, > 0 means stop at next NS step */
946 if ( (gs.set[eglsSTOPCOND] < 0) ||
947 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
949 bLastStep = TRUE;
952 /* Determine whether or not to update the Born radii if doing GB */
953 bBornRadii = bFirstStep;
954 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
956 bBornRadii = TRUE;
959 /* do_log triggers energy and virial calculation. Because this leads
960 * to different code paths, forces can be different. Thus for exact
961 * continuation we should avoid extra log output.
962 * Note that the || bLastStep can result in non-exact continuation
963 * beyond the last step. But we don't consider that to be an issue.
965 do_log = do_per_step(step, ir->nstlog) || (bFirstStep && !startingFromCheckpoint) || bLastStep || bRerunMD;
966 do_verbose = bVerbose &&
967 (step % stepout == 0 || bFirstStep || bLastStep || bRerunMD);
969 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
971 if (bRerunMD)
973 bMasterState = TRUE;
975 else
977 bMasterState = FALSE;
978 /* Correct the new box if it is too skewed */
979 if (inputrecDynamicBox(ir))
981 if (correct_box(fplog, step, state->box, graph))
983 bMasterState = TRUE;
986 if (DOMAINDECOMP(cr) && bMasterState)
988 dd_collect_state(cr->dd, state, state_global);
992 if (DOMAINDECOMP(cr))
994 /* Repartition the domain decomposition */
995 dd_partition_system(fplog, step, cr,
996 bMasterState, nstglobalcomm,
997 state_global, top_global, ir,
998 state, &f, mdatoms, top, fr,
999 vsite, constr,
1000 nrnb, wcycle,
1001 do_verbose && !bPMETunePrinting);
1002 shouldCheckNumberOfBondedInteractions = true;
1003 update_realloc(upd, state->nalloc);
1007 if (MASTER(cr) && do_log)
1009 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
1012 if (ir->efep != efepNO)
1014 update_mdatoms(mdatoms, state->lambda[efptMASS]);
1017 if ((bRerunMD && rerun_fr.bV) || bExchanged)
1020 /* We need the kinetic energy at minus the half step for determining
1021 * the full step kinetic energy and possibly for T-coupling.*/
1022 /* This may not be quite working correctly yet . . . . */
1023 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1024 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1025 constr, NULL, FALSE, state->box,
1026 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1027 CGLO_GSTAT | CGLO_TEMPERATURE | CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS);
1028 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1029 top_global, top, state,
1030 &shouldCheckNumberOfBondedInteractions);
1032 clear_mat(force_vir);
1034 /* We write a checkpoint at this MD step when:
1035 * either at an NS step when we signalled through gs,
1036 * or at the last step (but not when we do not want confout),
1037 * but never at the first step or with rerun.
1039 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1040 (bLastStep && (Flags & MD_CONFOUT))) &&
1041 step > ir->init_step && !bRerunMD);
1042 if (bCPT)
1044 gs.set[eglsCHKPT] = 0;
1047 /* Determine the energy and pressure:
1048 * at nstcalcenergy steps and at energy output steps (set below).
1050 if (EI_VV(ir->eI) && (!bInitStep))
1052 /* for vv, the first half of the integration actually corresponds
1053 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1054 but the virial needs to be calculated on both the current step and the 'next' step. Future
1055 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1057 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1058 bCalcEnerStep = do_per_step(step - 1, ir->nstcalcenergy);
1059 bCalcVir = bCalcEnerStep ||
1060 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1062 else
1064 bCalcEnerStep = do_per_step(step, ir->nstcalcenergy);
1065 bCalcVir = bCalcEnerStep ||
1066 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1068 bCalcEner = bCalcEnerStep;
1070 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep || bRerunMD);
1072 if (do_ene || do_log || bDoReplEx)
1074 bCalcVir = TRUE;
1075 bCalcEner = TRUE;
1078 /* Do we need global communication ? */
1079 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1080 do_per_step(step, nstglobalcomm) ||
1081 (EI_VV(ir->eI) && inputrecNvtTrotter(ir) && do_per_step(step-1, nstglobalcomm)));
1083 force_flags = (GMX_FORCE_STATECHANGED |
1084 ((inputrecDynamicBox(ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1085 GMX_FORCE_ALLFORCES |
1086 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1087 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1088 (bDoFEP ? GMX_FORCE_DHDL : 0)
1091 if (shellfc)
1093 /* Now is the time to relax the shells */
1094 relax_shell_flexcon(fplog, cr, bVerbose, step,
1095 ir, bNS, force_flags, top,
1096 constr, enerd, fcd,
1097 state, f, force_vir, mdatoms,
1098 nrnb, wcycle, graph, groups,
1099 shellfc, fr, bBornRadii, t, mu_tot,
1100 vsite, mdoutf_get_fp_field(outf));
1102 else
1104 /* The coordinates (x) are shifted (to get whole molecules)
1105 * in do_force.
1106 * This is parallellized as well, and does communication too.
1107 * Check comments in sim_util.c
1109 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1110 state->box, state->x, &state->hist,
1111 f, force_vir, mdatoms, enerd, fcd,
1112 state->lambda, graph,
1113 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1114 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1117 if (EI_VV(ir->eI) && !startingFromCheckpoint && !bRerunMD)
1118 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1120 rvec *vbuf = NULL;
1122 wallcycle_start(wcycle, ewcUPDATE);
1123 if (ir->eI == eiVV && bInitStep)
1125 /* if using velocity verlet with full time step Ekin,
1126 * take the first half step only to compute the
1127 * virial for the first step. From there,
1128 * revert back to the initial coordinates
1129 * so that the input is actually the initial step.
1131 snew(vbuf, state->natoms);
1132 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1134 else
1136 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1137 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1140 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1141 ekind, M, upd, etrtVELOCITY1,
1142 cr, constr);
1144 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1146 wallcycle_stop(wcycle, ewcUPDATE);
1147 update_constraints(fplog, step, NULL, ir, mdatoms,
1148 state, fr->bMolPBC, graph, f,
1149 &top->idef, shake_vir,
1150 cr, nrnb, wcycle, upd, constr,
1151 TRUE, bCalcVir);
1152 wallcycle_start(wcycle, ewcUPDATE);
1154 else if (graph)
1156 /* Need to unshift here if a do_force has been
1157 called in the previous step */
1158 unshift_self(graph, state->box, state->x);
1160 /* if VV, compute the pressure and constraints */
1161 /* For VV2, we strictly only need this if using pressure
1162 * control, but we really would like to have accurate pressures
1163 * printed out.
1164 * Think about ways around this in the future?
1165 * For now, keep this choice in comments.
1167 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1168 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1169 bPres = TRUE;
1170 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1171 if (bCalcEner && ir->eI == eiVVAK)
1173 bSumEkinhOld = TRUE;
1175 /* for vv, the first half of the integration actually corresponds to the previous step.
1176 So we need information from the last step in the first half of the integration */
1177 if (bGStat || do_per_step(step-1, nstglobalcomm))
1179 wallcycle_stop(wcycle, ewcUPDATE);
1180 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1181 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1182 constr, NULL, FALSE, state->box,
1183 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1184 (bGStat ? CGLO_GSTAT : 0)
1185 | CGLO_ENERGY
1186 | (bTemp ? CGLO_TEMPERATURE : 0)
1187 | (bPres ? CGLO_PRESSURE : 0)
1188 | (bPres ? CGLO_CONSTRAINT : 0)
1189 | (bStopCM ? CGLO_STOPCM : 0)
1190 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1191 | CGLO_SCALEEKIN
1193 /* explanation of above:
1194 a) We compute Ekin at the full time step
1195 if 1) we are using the AveVel Ekin, and it's not the
1196 initial step, or 2) if we are using AveEkin, but need the full
1197 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1198 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1199 EkinAveVel because it's needed for the pressure */
1200 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1201 top_global, top, state,
1202 &shouldCheckNumberOfBondedInteractions);
1203 wallcycle_start(wcycle, ewcUPDATE);
1205 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1206 if (!bInitStep)
1208 if (bTrotter)
1210 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1211 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1213 copy_mat(shake_vir, state->svir_prev);
1214 copy_mat(force_vir, state->fvir_prev);
1215 if (inputrecNvtTrotter(ir) && ir->eI == eiVV)
1217 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1218 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1219 enerd->term[F_EKIN] = trace(ekind->ekin);
1222 else if (bExchanged)
1224 wallcycle_stop(wcycle, ewcUPDATE);
1225 /* We need the kinetic energy at minus the half step for determining
1226 * the full step kinetic energy and possibly for T-coupling.*/
1227 /* This may not be quite working correctly yet . . . . */
1228 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1229 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1230 constr, NULL, FALSE, state->box,
1231 NULL, &bSumEkinhOld,
1232 CGLO_GSTAT | CGLO_TEMPERATURE);
1233 wallcycle_start(wcycle, ewcUPDATE);
1236 /* if it's the initial step, we performed this first step just to get the constraint virial */
1237 if (ir->eI == eiVV && bInitStep)
1239 copy_rvecn(vbuf, state->v, 0, state->natoms);
1240 sfree(vbuf);
1242 wallcycle_stop(wcycle, ewcUPDATE);
1245 /* compute the conserved quantity */
1246 if (EI_VV(ir->eI))
1248 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1249 if (ir->eI == eiVV)
1251 last_ekin = enerd->term[F_EKIN];
1253 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1255 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1257 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1258 if (ir->efep != efepNO && !bRerunMD)
1260 sum_dhdl(enerd, state->lambda, ir->fepvals);
1264 /* ######## END FIRST UPDATE STEP ############## */
1265 /* ######## If doing VV, we now have v(dt) ###### */
1266 if (bDoExpanded)
1268 /* perform extended ensemble sampling in lambda - we don't
1269 actually move to the new state before outputting
1270 statistics, but if performing simulated tempering, we
1271 do update the velocities and the tau_t. */
1273 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1274 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1275 copy_df_history(&state_global->dfhist, &state->dfhist);
1278 /* Now we have the energies and forces corresponding to the
1279 * coordinates at time t. We must output all of this before
1280 * the update.
1282 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1283 ir, state, state_global, top_global, fr,
1284 outf, mdebin, ekind, f,
1285 &nchkpt,
1286 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1287 bSumEkinhOld);
1288 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1289 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1291 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1292 if (startingFromCheckpoint && bTrotter)
1294 copy_mat(state->svir_prev, shake_vir);
1295 copy_mat(state->fvir_prev, force_vir);
1298 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1300 /* Check whether everything is still allright */
1301 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1302 #if GMX_THREAD_MPI
1303 && MASTER(cr)
1304 #endif
1307 /* this is just make gs.sig compatible with the hack
1308 of sending signals around by MPI_Reduce with together with
1309 other floats */
1310 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1312 gs.sig[eglsSTOPCOND] = 1;
1314 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1316 gs.sig[eglsSTOPCOND] = -1;
1318 /* < 0 means stop at next step, > 0 means stop at next NS step */
1319 if (fplog)
1321 fprintf(fplog,
1322 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1323 gmx_get_signal_name(),
1324 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1325 fflush(fplog);
1327 fprintf(stderr,
1328 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1329 gmx_get_signal_name(),
1330 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1331 fflush(stderr);
1332 handled_stop_condition = (int)gmx_get_stop_condition();
1334 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1335 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1336 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1338 /* Signal to terminate the run */
1339 gs.sig[eglsSTOPCOND] = 1;
1340 if (fplog)
1342 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1344 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1347 if (bResetCountersHalfMaxH && MASTER(cr) &&
1348 elapsed_time > max_hours*60.0*60.0*0.495)
1350 /* Set flag that will communicate the signal to all ranks in the simulation */
1351 gs.sig[eglsRESETCOUNTERS] = 1;
1354 /* In parallel we only have to check for checkpointing in steps
1355 * where we do global communication,
1356 * otherwise the other nodes don't know.
1358 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1359 cpt_period >= 0 &&
1360 (cpt_period == 0 ||
1361 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1362 gs.set[eglsCHKPT] == 0)
1364 gs.sig[eglsCHKPT] = 1;
1367 /* ######### START SECOND UPDATE STEP ################# */
1369 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1370 in preprocessing */
1372 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1374 gmx_bool bIfRandomize;
1375 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1376 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1377 if (constr && bIfRandomize)
1379 update_constraints(fplog, step, NULL, ir, mdatoms,
1380 state, fr->bMolPBC, graph, f,
1381 &top->idef, tmp_vir,
1382 cr, nrnb, wcycle, upd, constr,
1383 TRUE, bCalcVir);
1386 /* Box is changed in update() when we do pressure coupling,
1387 * but we should still use the old box for energy corrections and when
1388 * writing it to the energy file, so it matches the trajectory files for
1389 * the same timestep above. Make a copy in a separate array.
1391 copy_mat(state->box, lastbox);
1393 dvdl_constr = 0;
1395 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1397 wallcycle_start(wcycle, ewcUPDATE);
1398 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1399 if (bTrotter)
1401 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1402 /* We can only do Berendsen coupling after we have summed
1403 * the kinetic energy or virial. Since the happens
1404 * in global_state after update, we should only do it at
1405 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1408 else
1410 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1411 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1414 if (EI_VV(ir->eI))
1416 /* velocity half-step update */
1417 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1418 ekind, M, upd, etrtVELOCITY2,
1419 cr, constr);
1422 /* Above, initialize just copies ekinh into ekin,
1423 * it doesn't copy position (for VV),
1424 * and entire integrator for MD.
1427 if (ir->eI == eiVVAK)
1429 /* We probably only need md->homenr, not state->natoms */
1430 if (state->natoms > cbuf_nalloc)
1432 cbuf_nalloc = state->natoms;
1433 srenew(cbuf, cbuf_nalloc);
1435 copy_rvecn(state->x, cbuf, 0, state->natoms);
1438 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1439 ekind, M, upd, etrtPOSITION, cr, constr);
1440 wallcycle_stop(wcycle, ewcUPDATE);
1442 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1443 fr->bMolPBC, graph, f,
1444 &top->idef, shake_vir,
1445 cr, nrnb, wcycle, upd, constr,
1446 FALSE, bCalcVir);
1448 if (ir->eI == eiVVAK)
1450 /* erase F_EKIN and F_TEMP here? */
1451 /* just compute the kinetic energy at the half step to perform a trotter step */
1452 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1453 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1454 constr, NULL, FALSE, lastbox,
1455 NULL, &bSumEkinhOld,
1456 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1458 wallcycle_start(wcycle, ewcUPDATE);
1459 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1460 /* now we know the scaling, we can compute the positions again again */
1461 copy_rvecn(cbuf, state->x, 0, state->natoms);
1463 update_coords(fplog, step, ir, mdatoms, state, f, fcd,
1464 ekind, M, upd, etrtPOSITION, cr, constr);
1465 wallcycle_stop(wcycle, ewcUPDATE);
1467 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1468 /* are the small terms in the shake_vir here due
1469 * to numerical errors, or are they important
1470 * physically? I'm thinking they are just errors, but not completely sure.
1471 * For now, will call without actually constraining, constr=NULL*/
1472 update_constraints(fplog, step, NULL, ir, mdatoms,
1473 state, fr->bMolPBC, graph, f,
1474 &top->idef, tmp_vir,
1475 cr, nrnb, wcycle, upd, NULL,
1476 FALSE, bCalcVir);
1478 if (EI_VV(ir->eI))
1480 /* this factor or 2 correction is necessary
1481 because half of the constraint force is removed
1482 in the vv step, so we have to double it. See
1483 the Redmine issue #1255. It is not yet clear
1484 if the factor of 2 is exact, or just a very
1485 good approximation, and this will be
1486 investigated. The next step is to see if this
1487 can be done adding a dhdl contribution from the
1488 rattle step, but this is somewhat more
1489 complicated with the current code. Will be
1490 investigated, hopefully for 4.6.3. However,
1491 this current solution is much better than
1492 having it completely wrong.
1494 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1496 else
1498 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1501 else if (graph)
1503 /* Need to unshift here */
1504 unshift_self(graph, state->box, state->x);
1507 if (vsite != NULL)
1509 wallcycle_start(wcycle, ewcVSITECONSTR);
1510 if (graph != NULL)
1512 shift_self(graph, state->box, state->x);
1514 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1515 top->idef.iparams, top->idef.il,
1516 fr->ePBC, fr->bMolPBC, cr, state->box);
1518 if (graph != NULL)
1520 unshift_self(graph, state->box, state->x);
1522 wallcycle_stop(wcycle, ewcVSITECONSTR);
1525 /* ############## IF NOT VV, Calculate globals HERE ############ */
1526 /* With Leap-Frog we can skip compute_globals at
1527 * non-communication steps, but we need to calculate
1528 * the kinetic energy one step before communication.
1530 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1532 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
1533 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1534 constr, &gs,
1535 (step_rel % gs.nstms == 0) &&
1536 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1537 lastbox,
1538 &totalNumberOfBondedInteractions, &bSumEkinhOld,
1539 (bGStat ? CGLO_GSTAT : 0)
1540 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1541 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1542 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1543 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1544 | CGLO_CONSTRAINT
1545 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
1547 checkNumberOfBondedInteractions(fplog, cr, totalNumberOfBondedInteractions,
1548 top_global, top, state,
1549 &shouldCheckNumberOfBondedInteractions);
1552 /* ############# END CALC EKIN AND PRESSURE ################# */
1554 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1555 the virial that should probably be addressed eventually. state->veta has better properies,
1556 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1557 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1559 if (ir->efep != efepNO && (!EI_VV(ir->eI) || bRerunMD))
1561 /* Sum up the foreign energy and dhdl terms for md and sd.
1562 Currently done every step so that dhdl is correct in the .edr */
1563 sum_dhdl(enerd, state->lambda, ir->fepvals);
1565 update_box(fplog, step, ir, mdatoms, state, f,
1566 pcoupl_mu, nrnb, upd);
1568 /* ################# END UPDATE STEP 2 ################# */
1569 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1571 /* The coordinates (x) were unshifted in update */
1572 if (!bGStat)
1574 /* We will not sum ekinh_old,
1575 * so signal that we still have to do it.
1577 bSumEkinhOld = TRUE;
1580 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1582 /* use the directly determined last velocity, not actually the averaged half steps */
1583 if (bTrotter && ir->eI == eiVV)
1585 enerd->term[F_EKIN] = last_ekin;
1587 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1589 if (EI_VV(ir->eI))
1591 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1593 else
1595 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1597 /* ######### END PREPARING EDR OUTPUT ########### */
1599 /* Output stuff */
1600 if (MASTER(cr))
1602 if (fplog && do_log && bDoExpanded)
1604 /* only needed if doing expanded ensemble */
1605 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1606 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1608 if (bCalcEner)
1610 upd_mdebin(mdebin, bDoDHDL, bCalcEnerStep,
1611 t, mdatoms->tmass, enerd, state,
1612 ir->fepvals, ir->expandedvals, lastbox,
1613 shake_vir, force_vir, total_vir, pres,
1614 ekind, mu_tot, constr);
1616 else
1618 upd_mdebin_step(mdebin);
1621 gmx_bool do_dr = do_per_step(step, ir->nstdisreout);
1622 gmx_bool do_or = do_per_step(step, ir->nstorireout);
1624 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1625 step, t,
1626 eprNORMAL, mdebin, fcd, groups, &(ir->opts));
1628 if (ir->bPull)
1630 pull_print_output(ir->pull_work, step, t);
1633 if (do_per_step(step, ir->nstlog))
1635 if (fflush(fplog) != 0)
1637 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1641 if (bDoExpanded)
1643 /* Have to do this part _after_ outputting the logfile and the edr file */
1644 /* Gets written into the state at the beginning of next loop*/
1645 state->fep_state = lamnew;
1647 /* Print the remaining wall clock time for the run */
1648 if (MULTIMASTER(cr) &&
1649 (do_verbose || gmx_got_usr_signal()) &&
1650 !bPMETunePrinting)
1652 if (shellfc)
1654 fprintf(stderr, "\n");
1656 print_time(stderr, walltime_accounting, step, ir, cr);
1659 /* Ion/water position swapping.
1660 * Not done in last step since trajectory writing happens before this call
1661 * in the MD loop and exchanges would be lost anyway. */
1662 bNeedRepartition = FALSE;
1663 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1664 do_per_step(step, ir->swap->nstswap))
1666 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1667 bRerunMD ? rerun_fr.x : state->x,
1668 bRerunMD ? rerun_fr.box : state->box,
1669 top_global, MASTER(cr) && bVerbose, bRerunMD);
1671 if (bNeedRepartition && DOMAINDECOMP(cr))
1673 dd_collect_state(cr->dd, state, state_global);
1677 /* Replica exchange */
1678 bExchanged = FALSE;
1679 if (bDoReplEx)
1681 bExchanged = replica_exchange(fplog, cr, repl_ex,
1682 state_global, enerd,
1683 state, step, t);
1686 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1688 dd_partition_system(fplog, step, cr, TRUE, 1,
1689 state_global, top_global, ir,
1690 state, &f, mdatoms, top, fr,
1691 vsite, constr,
1692 nrnb, wcycle, FALSE);
1693 shouldCheckNumberOfBondedInteractions = true;
1694 update_realloc(upd, state->nalloc);
1697 bFirstStep = FALSE;
1698 bInitStep = FALSE;
1699 startingFromCheckpoint = FALSE;
1701 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1702 /* With all integrators, except VV, we need to retain the pressure
1703 * at the current step for coupling at the next step.
1705 if ((state->flags & (1<<estPRES_PREV)) &&
1706 (bGStatEveryStep ||
1707 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1709 /* Store the pressure in t_state for pressure coupling
1710 * at the next MD step.
1712 copy_mat(pres, state->pres_prev);
1715 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1717 if ( (membed != NULL) && (!bLastStep) )
1719 rescale_membed(step_rel, membed, state_global->x);
1722 if (bRerunMD)
1724 if (MASTER(cr))
1726 /* read next frame from input trajectory */
1727 bLastStep = !read_next_frame(oenv, status, &rerun_fr);
1730 if (PAR(cr))
1732 rerun_parallel_comm(cr, &rerun_fr, &bLastStep);
1736 cycles = wallcycle_stop(wcycle, ewcSTEP);
1737 if (DOMAINDECOMP(cr) && wcycle)
1739 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1742 if (!bRerunMD || !rerun_fr.bStep)
1744 /* increase the MD step number */
1745 step++;
1746 step_rel++;
1749 /* TODO make a counter-reset module */
1750 /* If it is time to reset counters, set a flag that remains
1751 true until counters actually get reset */
1752 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1753 gs.set[eglsRESETCOUNTERS] != 0)
1755 if (pme_loadbal_is_active(pme_loadbal))
1757 /* Do not permit counter reset while PME load
1758 * balancing is active. The only purpose for resetting
1759 * counters is to measure reliable performance data,
1760 * and that can't be done before balancing
1761 * completes.
1763 * TODO consider fixing this by delaying the reset
1764 * until after load balancing completes,
1765 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1766 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1767 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1768 "resetting counters later in the run, e.g. with gmx "
1769 "mdrun -resetstep.", step);
1771 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1772 use_GPU(fr->nbv) ? fr->nbv : NULL);
1773 wcycle_set_reset_counters(wcycle, -1);
1774 if (!(cr->duty & DUTY_PME))
1776 /* Tell our PME node to reset its counters */
1777 gmx_pme_send_resetcounters(cr, step);
1779 /* Correct max_hours for the elapsed time */
1780 max_hours -= elapsed_time/(60.0*60.0);
1781 /* If mdrun -maxh -resethway was active, it can only trigger once */
1782 bResetCountersHalfMaxH = FALSE; /* TODO move this to where gs.sig[eglsRESETCOUNTERS] is set */
1783 /* Reset can only happen once, so clear the triggering flag. */
1784 gs.set[eglsRESETCOUNTERS] = 0;
1787 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1788 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1791 /* End of main MD loop */
1793 /* Closing TNG files can include compressing data. Therefore it is good to do that
1794 * before stopping the time measurements. */
1795 mdoutf_tng_close(outf);
1797 /* Stop measuring walltime */
1798 walltime_accounting_end(walltime_accounting);
1800 if (bRerunMD && MASTER(cr))
1802 close_trj(status);
1805 if (!(cr->duty & DUTY_PME))
1807 /* Tell the PME only node to finish */
1808 gmx_pme_send_finish(cr);
1811 if (MASTER(cr))
1813 if (ir->nstcalcenergy > 0 && !bRerunMD)
1815 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1816 eprAVER, mdebin, fcd, groups, &(ir->opts));
1820 done_mdoutf(outf);
1822 if (bPMETune)
1824 pme_loadbal_done(pme_loadbal, cr, fplog, use_GPU(fr->nbv));
1827 done_shellfc(fplog, shellfc, step_rel);
1829 if (repl_ex_nst > 0 && MASTER(cr))
1831 print_replica_exchange_statistics(fplog, repl_ex);
1834 // TODO clean up swapcoords
1836 if (ir->bRot)
1838 finish_rot(ir->rot);
1841 if (ir->bPull)
1843 finish_pull(ir->pull_work);
1846 if (membed != nullptr)
1848 free_membed(membed);
1851 /* IMD cleanup, if bIMD is TRUE. */
1852 IMD_finalize(ir->bIMD, ir->imd);
1854 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1856 return 0;