Moved df_history utility routines from typedefs to new file.
[gromacs.git] / src / programs / mdrun / md.cpp
blob4b96aec3400898a1292fb1a4455785c7066504f6
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, 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/domdec/domdec.h"
50 #include "gromacs/domdec/domdec_network.h"
51 #include "gromacs/ewald/pme.h"
52 #include "gromacs/ewald/pme-load-balancing.h"
53 #include "gromacs/fileio/filenm.h"
54 #include "gromacs/fileio/trx.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/gmxlib/df_history.h"
57 #include "gromacs/gmxlib/energyhistory.h"
58 #include "gromacs/gmxlib/md_logging.h"
59 #include "gromacs/gmxlib/sighandler.h"
60 #include "gromacs/imd/imd.h"
61 #include "gromacs/legacyheaders/force.h"
62 #include "gromacs/legacyheaders/network.h"
63 #include "gromacs/legacyheaders/nrnb.h"
64 #include "gromacs/legacyheaders/ns.h"
65 #include "gromacs/legacyheaders/typedefs.h"
66 #include "gromacs/legacyheaders/vsite.h"
67 #include "gromacs/legacyheaders/types/commrec.h"
68 #include "gromacs/legacyheaders/types/enums.h"
69 #include "gromacs/legacyheaders/types/fcdata.h"
70 #include "gromacs/legacyheaders/types/force_flags.h"
71 #include "gromacs/legacyheaders/types/forcerec.h"
72 #include "gromacs/legacyheaders/types/group.h"
73 #include "gromacs/legacyheaders/types/inputrec.h"
74 #include "gromacs/legacyheaders/types/interaction_const.h"
75 #include "gromacs/legacyheaders/types/mdatom.h"
76 #include "gromacs/legacyheaders/types/nrnb.h"
77 #include "gromacs/legacyheaders/types/state.h"
78 #include "gromacs/listed-forces/manage-threading.h"
79 #include "gromacs/math/utilities.h"
80 #include "gromacs/math/vec.h"
81 #include "gromacs/math/vectypes.h"
82 #include "gromacs/mdlib/compute_io.h"
83 #include "gromacs/mdlib/constr.h"
84 #include "gromacs/mdlib/ebin.h"
85 #include "gromacs/mdlib/forcerec.h"
86 #include "gromacs/mdlib/md_support.h"
87 #include "gromacs/mdlib/mdatoms.h"
88 #include "gromacs/mdlib/mdebin.h"
89 #include "gromacs/mdlib/mdoutf.h"
90 #include "gromacs/mdlib/mdrun.h"
91 #include "gromacs/mdlib/mdrun_signalling.h"
92 #include "gromacs/mdlib/nb_verlet.h"
93 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
94 #include "gromacs/mdlib/shellfc.h"
95 #include "gromacs/mdlib/sim_util.h"
96 #include "gromacs/mdlib/tgroup.h"
97 #include "gromacs/mdlib/trajectory_writing.h"
98 #include "gromacs/mdlib/update.h"
99 #include "gromacs/mdlib/vcm.h"
100 #include "gromacs/pbcutil/mshift.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/pulling/pull.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/utility/basedefinitions.h"
111 #include "gromacs/utility/cstringutil.h"
112 #include "gromacs/utility/fatalerror.h"
113 #include "gromacs/utility/real.h"
114 #include "gromacs/utility/smalloc.h"
116 #include "deform.h"
117 #include "membed.h"
118 #include "repl_ex.h"
120 #ifdef GMX_FAHCORE
121 #include "corewrap.h"
122 #endif
124 static void reset_all_counters(FILE *fplog, t_commrec *cr,
125 gmx_int64_t step,
126 gmx_int64_t *step_rel, t_inputrec *ir,
127 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
128 gmx_walltime_accounting_t walltime_accounting,
129 struct nonbonded_verlet_t *nbv)
131 char sbuf[STEPSTRSIZE];
133 /* Reset all the counters related to performance over the run */
134 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
135 gmx_step_str(step, sbuf));
137 if (use_GPU(nbv))
139 nbnxn_gpu_reset_timings(nbv);
142 wallcycle_stop(wcycle, ewcRUN);
143 wallcycle_reset_all(wcycle);
144 if (DOMAINDECOMP(cr))
146 reset_dd_statistics_counters(cr->dd);
148 init_nrnb(nrnb);
149 ir->init_step += *step_rel;
150 ir->nsteps -= *step_rel;
151 *step_rel = 0;
152 wallcycle_start(wcycle, ewcRUN);
153 walltime_accounting_start(walltime_accounting);
154 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
157 /*! \libinternal
158 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
159 int nfile, const t_filenm fnm[],
160 const gmx_output_env_t *oenv, gmx_bool bVerbose,
161 gmx_bool bCompact, int nstglobalcomm,
162 gmx_vsite_t *vsite, gmx_constr_t constr,
163 int stepout,
164 t_inputrec *inputrec,
165 gmx_mtop_t *top_global, t_fcdata *fcd,
166 t_state *state_global,
167 t_mdatoms *mdatoms,
168 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
169 gmx_edsam_t ed,
170 t_forcerec *fr,
171 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
172 gmx_membed_t *membed,
173 real cpt_period, real max_hours,
174 int imdport,
175 unsigned long Flags,
176 gmx_walltime_accounting_t walltime_accounting)
178 double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
179 const gmx_output_env_t *oenv, gmx_bool bVerbose, gmx_bool bCompact,
180 int nstglobalcomm,
181 gmx_vsite_t *vsite, gmx_constr_t constr,
182 int stepout, t_inputrec *ir,
183 gmx_mtop_t *top_global,
184 t_fcdata *fcd,
185 t_state *state_global,
186 t_mdatoms *mdatoms,
187 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
188 gmx_edsam_t ed, t_forcerec *fr,
189 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t *membed,
190 real cpt_period, real max_hours,
191 int imdport,
192 unsigned long Flags,
193 gmx_walltime_accounting_t walltime_accounting)
195 gmx_mdoutf_t outf = NULL;
196 gmx_int64_t step, step_rel;
197 double elapsed_time;
198 double t, t0, lam0[efptNR];
199 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
200 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
201 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
202 bBornRadii, bStartingFromCpt;
203 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
204 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
205 bForceUpdate = FALSE, bCPT;
206 gmx_bool bMasterState;
207 int force_flags, cglo_flags;
208 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
209 int i, m;
210 t_trxstatus *status;
211 rvec mu_tot;
212 t_vcm *vcm;
213 matrix pcoupl_mu, M;
214 t_trxframe rerun_fr;
215 gmx_repl_ex_t repl_ex = NULL;
216 int nchkpt = 1;
217 gmx_localtop_t *top;
218 t_mdebin *mdebin = NULL;
219 t_state *state = NULL;
220 rvec *f_global = NULL;
221 gmx_enerdata_t *enerd;
222 rvec *f = NULL;
223 gmx_global_stat_t gstat;
224 gmx_update_t upd = NULL;
225 t_graph *graph = NULL;
226 gmx_signalling_t gs;
227 gmx_groups_t *groups;
228 gmx_ekindata_t *ekind;
229 gmx_shellfc_t *shellfc;
230 int count, nconverged = 0;
231 double tcount = 0;
232 gmx_bool bConverged = TRUE, bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
233 gmx_bool bResetCountersHalfMaxH = FALSE;
234 gmx_bool bVV, bTemp, bPres, bTrotter;
235 gmx_bool bUpdateDoLR;
236 real dvdl_constr;
237 rvec *cbuf = NULL;
238 int cbuf_nalloc = 0;
239 matrix lastbox;
240 int lamnew = 0;
241 /* for FEP */
242 int nstfep = 0;
243 double cycles;
244 real saved_conserved_quantity = 0;
245 real last_ekin = 0;
246 t_extmass MassQ;
247 int **trotter_seq;
248 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
249 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
250 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
251 simulation stops. If equal to zero, don't
252 communicate any more between multisims.*/
253 /* PME load balancing data for GPU kernels */
254 pme_load_balancing_t *pme_loadbal = NULL;
255 gmx_bool bPMETune = FALSE;
256 gmx_bool bPMETunePrinting = FALSE;
258 /* Interactive MD */
259 gmx_bool bIMDstep = FALSE;
261 #ifdef GMX_FAHCORE
262 /* Temporary addition for FAHCORE checkpointing */
263 int chkpt_ret;
264 #endif
266 /* Check for special mdrun options */
267 bRerunMD = (Flags & MD_RERUN);
268 if (Flags & MD_RESETCOUNTERSHALFWAY)
270 if (ir->nsteps > 0)
272 /* Signal to reset the counters half the simulation steps. */
273 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
275 /* Signal to reset the counters halfway the simulation time. */
276 bResetCountersHalfMaxH = (max_hours > 0);
279 /* md-vv uses averaged full step velocities for T-control
280 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
281 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
282 bVV = EI_VV(ir->eI);
283 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
285 if (bRerunMD)
287 /* Since we don't know if the frames read are related in any way,
288 * rebuild the neighborlist at every step.
290 ir->nstlist = 1;
291 ir->nstcalcenergy = 1;
292 nstglobalcomm = 1;
295 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
297 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
298 bGStatEveryStep = (nstglobalcomm == 1);
300 if (bRerunMD)
302 ir->nstxout_compressed = 0;
304 groups = &top_global->groups;
306 /* Initial values */
307 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
308 &(state_global->fep_state), lam0,
309 nrnb, top_global, &upd,
310 nfile, fnm, &outf, &mdebin,
311 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
313 clear_mat(total_vir);
314 clear_mat(pres);
315 /* Energy terms and groups */
316 snew(enerd, 1);
317 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
318 enerd);
319 if (DOMAINDECOMP(cr))
321 f = NULL;
323 else
325 snew(f, top_global->natoms);
328 /* Kinetic energy data */
329 snew(ekind, 1);
330 init_ekindata(fplog, top_global, &(ir->opts), ekind);
331 /* Copy the cos acceleration to the groups struct */
332 ekind->cosacc.cos_accel = ir->cos_accel;
334 gstat = global_stat_init(ir);
335 debug_gmx();
337 /* Check for polarizable models and flexible constraints */
338 shellfc = init_shell_flexcon(fplog,
339 top_global, n_flexible_constraints(constr),
340 (ir->bContinuation ||
341 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
342 NULL : state_global->x);
343 if (shellfc && ir->nstcalcenergy != 1)
345 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);
347 if (shellfc && DOMAINDECOMP(cr))
349 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
351 if (shellfc && ir->eI == eiNM)
353 /* Currently shells don't work with Normal Modes */
354 gmx_fatal(FARGS, "Normal Mode analysis is not supported with shells.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
357 if (vsite && ir->eI == eiNM)
359 /* Currently virtual sites don't work with Normal Modes */
360 gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
363 if (DEFORM(*ir))
365 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
366 set_deform_reference_box(upd,
367 deform_init_init_step_tpx,
368 deform_init_box_tpx);
369 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
373 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
374 if ((io > 2000) && MASTER(cr))
376 fprintf(stderr,
377 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
378 io);
382 if (DOMAINDECOMP(cr))
384 top = dd_init_local_top(top_global);
386 snew(state, 1);
387 dd_init_local_state(cr->dd, state_global, state);
389 if (DDMASTER(cr->dd) && ir->nstfout)
391 snew(f_global, state_global->natoms);
394 else
396 top = gmx_mtop_generate_local_top(top_global, ir);
398 forcerec_set_excl_load(fr, top);
400 state = serial_init_local_state(state_global);
401 f_global = f;
403 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
405 if (vsite)
407 set_vsite_top(vsite, top, mdatoms, cr);
410 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
412 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
415 if (shellfc)
417 make_local_shells(cr, mdatoms, shellfc);
420 setup_bonded_threading(fr, &top->idef);
423 /* Set up interactive MD (IMD) */
424 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
425 nfile, fnm, oenv, imdport, Flags);
427 if (DOMAINDECOMP(cr))
429 /* Distribute the charge groups over the nodes from the master node */
430 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
431 state_global, top_global, ir,
432 state, &f, mdatoms, top, fr,
433 vsite, shellfc, constr,
434 nrnb, NULL, FALSE);
438 update_mdatoms(mdatoms, state->lambda[efptMASS]);
440 if (opt2bSet("-cpi", nfile, fnm))
442 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
444 else
446 bStateFromCP = FALSE;
449 if (ir->bExpanded)
451 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
454 if (MASTER(cr))
456 if (bStateFromCP)
458 /* Update mdebin with energy history if appending to output files */
459 if (Flags & MD_APPENDFILES)
461 restore_energyhistory_from_state(mdebin, state_global->enerhist);
463 else
465 /* We might have read an energy history from checkpoint,
466 * free the allocated memory and reset the counts.
468 done_energyhistory(state_global->enerhist);
469 init_energyhistory(state_global->enerhist);
472 /* Set the initial energy history in state by updating once */
473 update_energyhistory(state_global->enerhist, mdebin);
476 /* Initialize constraints */
477 if (constr && !DOMAINDECOMP(cr))
479 set_constraints(constr, top, ir, mdatoms, cr);
482 if (repl_ex_nst > 0 && MASTER(cr))
484 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
485 repl_ex_nst, repl_ex_nex, repl_ex_seed);
488 /* PME tuning is only supported with PME for Coulomb. Is is not supported
489 * with only LJ PME, or for reruns.
491 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
492 !(Flags & MD_REPRODUCIBLE));
493 if (bPMETune)
495 pme_loadbal_init(&pme_loadbal, cr, fplog, ir, state->box,
496 fr->ic, fr->pmedata, use_GPU(fr->nbv),
497 &bPMETunePrinting);
500 if (!ir->bContinuation && !bRerunMD)
502 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
504 /* Set the velocities of frozen particles to zero */
505 for (i = 0; i < mdatoms->homenr; i++)
507 for (m = 0; m < DIM; m++)
509 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
511 state->v[i][m] = 0;
517 if (constr)
519 /* Constrain the initial coordinates and velocities */
520 do_constrain_first(fplog, constr, ir, mdatoms, state,
521 cr, nrnb, fr, top);
523 if (vsite)
525 /* Construct the virtual sites for the initial configuration */
526 construct_vsites(vsite, state->x, ir->delta_t, NULL,
527 top->idef.iparams, top->idef.il,
528 fr->ePBC, fr->bMolPBC, cr, state->box);
532 debug_gmx();
534 if (IR_TWINRANGE(*ir) && repl_ex_nst % ir->nstcalclr != 0)
536 /* We should exchange at nstcalclr steps to get correct integration */
537 gmx_fatal(FARGS, "The replica exchange period (%d) is not divisible by nstcalclr (%d)", repl_ex_nst, ir->nstcalclr);
540 if (ir->efep != efepNO)
542 /* Set free energy calculation frequency as the greatest common
543 * denominator of nstdhdl and repl_ex_nst.
544 * Check for nstcalclr with twin-range, since we need the long-range
545 * contribution to the free-energy at the correct (nstcalclr) steps.
547 nstfep = ir->fepvals->nstdhdl;
548 if (ir->bExpanded)
550 if (IR_TWINRANGE(*ir) &&
551 ir->expandedvals->nstexpanded % ir->nstcalclr != 0)
553 gmx_fatal(FARGS, "nstexpanded should be divisible by nstcalclr");
555 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
557 if (repl_ex_nst > 0)
559 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
561 /* We checked divisibility of repl_ex_nst and nstcalclr above */
562 if (IR_TWINRANGE(*ir) && nstfep % ir->nstcalclr != 0)
564 gmx_incons("nstfep not divisible by nstcalclr");
568 /* Be REALLY careful about what flags you set here. You CANNOT assume
569 * this is the first step, since we might be restarting from a checkpoint,
570 * and in that case we should not do any modifications to the state.
572 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
574 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
575 | (bStopCM ? CGLO_STOPCM : 0)
576 | (bVV ? CGLO_PRESSURE : 0)
577 | (bVV ? CGLO_CONSTRAINT : 0)
578 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
580 bSumEkinhOld = FALSE;
581 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
582 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
583 constr, NULL, FALSE, state->box,
584 top_global, &bSumEkinhOld, cglo_flags);
585 if (ir->eI == eiVVAK)
587 /* a second call to get the half step temperature initialized as well */
588 /* we do the same call as above, but turn the pressure off -- internally to
589 compute_globals, this is recognized as a velocity verlet half-step
590 kinetic energy calculation. This minimized excess variables, but
591 perhaps loses some logic?*/
593 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
594 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
595 constr, NULL, FALSE, state->box,
596 top_global, &bSumEkinhOld,
597 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
600 /* Calculate the initial half step temperature, and save the ekinh_old */
601 if (!(Flags & MD_STARTFROMCPT))
603 for (i = 0; (i < ir->opts.ngtc); i++)
605 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
608 if (ir->eI != eiVV)
610 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
611 and there is no previous step */
614 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
615 temperature control */
616 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
618 if (MASTER(cr))
620 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
622 fprintf(fplog,
623 "RMS relative constraint deviation after constraining: %.2e\n",
624 constr_rmsd(constr, FALSE));
626 if (EI_STATE_VELOCITY(ir->eI))
628 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
630 if (bRerunMD)
632 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
633 " input trajectory '%s'\n\n",
634 *(top_global->name), opt2fn("-rerun", nfile, fnm));
635 if (bVerbose)
637 fprintf(stderr, "Calculated time to finish depends on nsteps from "
638 "run input file,\nwhich may not correspond to the time "
639 "needed to process input trajectory.\n\n");
642 else
644 char tbuf[20];
645 fprintf(stderr, "starting mdrun '%s'\n",
646 *(top_global->name));
647 if (ir->nsteps >= 0)
649 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
651 else
653 sprintf(tbuf, "%s", "infinite");
655 if (ir->init_step > 0)
657 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
658 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
659 gmx_step_str(ir->init_step, sbuf2),
660 ir->init_step*ir->delta_t);
662 else
664 fprintf(stderr, "%s steps, %s ps.\n",
665 gmx_step_str(ir->nsteps, sbuf), tbuf);
668 fprintf(fplog, "\n");
671 walltime_accounting_start(walltime_accounting);
672 wallcycle_start(wcycle, ewcRUN);
673 print_start(fplog, cr, walltime_accounting, "mdrun");
675 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
676 #ifdef GMX_FAHCORE
677 chkpt_ret = fcCheckPointParallel( cr->nodeid,
678 NULL, 0);
679 if (chkpt_ret == 0)
681 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
683 #endif
685 debug_gmx();
686 /***********************************************************
688 * Loop over MD steps
690 ************************************************************/
692 /* if rerunMD then read coordinates and velocities from input trajectory */
693 if (bRerunMD)
695 if (getenv("GMX_FORCE_UPDATE"))
697 bForceUpdate = TRUE;
700 rerun_fr.natoms = 0;
701 if (MASTER(cr))
703 bNotLastFrame = read_first_frame(oenv, &status,
704 opt2fn("-rerun", nfile, fnm),
705 &rerun_fr, TRX_NEED_X | TRX_READ_V);
706 if (rerun_fr.natoms != top_global->natoms)
708 gmx_fatal(FARGS,
709 "Number of atoms in trajectory (%d) does not match the "
710 "run input file (%d)\n",
711 rerun_fr.natoms, top_global->natoms);
713 if (ir->ePBC != epbcNONE)
715 if (!rerun_fr.bBox)
717 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);
719 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
721 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
726 if (PAR(cr))
728 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
731 if (ir->ePBC != epbcNONE)
733 /* Set the shift vectors.
734 * Necessary here when have a static box different from the tpr box.
736 calc_shifts(rerun_fr.box, fr->shift_vec);
740 /* loop over MD steps or if rerunMD to end of input trajectory */
741 bFirstStep = TRUE;
742 /* Skip the first Nose-Hoover integration when we get the state from tpx */
743 bStateFromTPX = !bStateFromCP;
744 bInitStep = bFirstStep && (bStateFromTPX || bVV);
745 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
746 bSumEkinhOld = FALSE;
747 bExchanged = FALSE;
748 bNeedRepartition = FALSE;
750 init_global_signals(&gs, cr, ir, repl_ex_nst);
752 step = ir->init_step;
753 step_rel = 0;
755 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
757 /* check how many steps are left in other sims */
758 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
762 /* and stop now if we should */
763 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
764 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
765 while (!bLastStep || (bRerunMD && bNotLastFrame))
768 /* Determine if this is a neighbor search step */
769 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
771 if (bPMETune && bNStList)
773 /* PME grid + cut-off optimization with GPUs or PME nodes */
774 pme_loadbal_do(pme_loadbal, cr,
775 (bVerbose && MASTER(cr)) ? stderr : NULL,
776 fplog,
777 ir, fr, state,
778 wcycle,
779 step, step_rel,
780 &bPMETunePrinting);
783 wallcycle_start(wcycle, ewcSTEP);
785 if (bRerunMD)
787 if (rerun_fr.bStep)
789 step = rerun_fr.step;
790 step_rel = step - ir->init_step;
792 if (rerun_fr.bTime)
794 t = rerun_fr.time;
796 else
798 t = step;
801 else
803 bLastStep = (step_rel == ir->nsteps);
804 t = t0 + step*ir->delta_t;
807 if (ir->efep != efepNO || ir->bSimTemp)
809 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
810 requiring different logic. */
812 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
813 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
814 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
815 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
816 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
819 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
820 do_per_step(step, repl_ex_nst));
822 if (bSimAnn)
824 update_annealing_target_temp(&(ir->opts), t);
827 if (bRerunMD)
829 if (!DOMAINDECOMP(cr) || MASTER(cr))
831 for (i = 0; i < state_global->natoms; i++)
833 copy_rvec(rerun_fr.x[i], state_global->x[i]);
835 if (rerun_fr.bV)
837 for (i = 0; i < state_global->natoms; i++)
839 copy_rvec(rerun_fr.v[i], state_global->v[i]);
842 else
844 for (i = 0; i < state_global->natoms; i++)
846 clear_rvec(state_global->v[i]);
848 if (bRerunWarnNoV)
850 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
851 " Ekin, temperature and pressure are incorrect,\n"
852 " the virial will be incorrect when constraints are present.\n"
853 "\n");
854 bRerunWarnNoV = FALSE;
858 copy_mat(rerun_fr.box, state_global->box);
859 copy_mat(state_global->box, state->box);
861 if (vsite && (Flags & MD_RERUN_VSITE))
863 if (DOMAINDECOMP(cr))
865 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
867 if (graph)
869 /* Following is necessary because the graph may get out of sync
870 * with the coordinates if we only have every N'th coordinate set
872 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
873 shift_self(graph, state->box, state->x);
875 construct_vsites(vsite, state->x, ir->delta_t, state->v,
876 top->idef.iparams, top->idef.il,
877 fr->ePBC, fr->bMolPBC, cr, state->box);
878 if (graph)
880 unshift_self(graph, state->box, state->x);
885 /* Stop Center of Mass motion */
886 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
888 if (bRerunMD)
890 /* for rerun MD always do Neighbour Searching */
891 bNS = (bFirstStep || ir->nstlist != 0);
892 bNStList = bNS;
894 else
896 /* Determine whether or not to do Neighbour Searching and LR */
897 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
900 /* check whether we should stop because another simulation has
901 stopped. */
902 if (MULTISIM(cr))
904 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
905 (multisim_nsteps != ir->nsteps) )
907 if (bNS)
909 if (MASTER(cr))
911 fprintf(stderr,
912 "Stopping simulation %d because another one has finished\n",
913 cr->ms->sim);
915 bLastStep = TRUE;
916 gs.sig[eglsCHKPT] = 1;
921 /* < 0 means stop at next step, > 0 means stop at next NS step */
922 if ( (gs.set[eglsSTOPCOND] < 0) ||
923 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
925 bLastStep = TRUE;
928 /* Determine whether or not to update the Born radii if doing GB */
929 bBornRadii = bFirstStep;
930 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
932 bBornRadii = TRUE;
935 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
936 do_verbose = bVerbose &&
937 (step % stepout == 0 || bFirstStep || bLastStep);
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 (DYNAMIC_BOX(*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, shellfc, constr,
970 nrnb, wcycle,
971 do_verbose && !bPMETunePrinting);
975 if (MASTER(cr) && do_log)
977 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
980 if (ir->efep != efepNO)
982 update_mdatoms(mdatoms, state->lambda[efptMASS]);
985 if ((bRerunMD && rerun_fr.bV) || bExchanged)
988 /* We need the kinetic energy at minus the half step for determining
989 * the full step kinetic energy and possibly for T-coupling.*/
990 /* This may not be quite working correctly yet . . . . */
991 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
992 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
993 constr, NULL, FALSE, state->box,
994 top_global, &bSumEkinhOld,
995 CGLO_GSTAT | CGLO_TEMPERATURE);
997 clear_mat(force_vir);
999 /* We write a checkpoint at this MD step when:
1000 * either at an NS step when we signalled through gs,
1001 * or at the last step (but not when we do not want confout),
1002 * but never at the first step or with rerun.
1004 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1005 (bLastStep && (Flags & MD_CONFOUT))) &&
1006 step > ir->init_step && !bRerunMD);
1007 if (bCPT)
1009 gs.set[eglsCHKPT] = 0;
1012 /* Determine the energy and pressure:
1013 * at nstcalcenergy steps and at energy output steps (set below).
1015 if (EI_VV(ir->eI) && (!bInitStep))
1017 /* for vv, the first half of the integration actually corresponds
1018 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1019 but the virial needs to be calculated on both the current step and the 'next' step. Future
1020 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1022 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1023 bCalcVir = bCalcEner ||
1024 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1026 else
1028 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1029 bCalcVir = bCalcEner ||
1030 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1033 /* Do we need global communication ? */
1034 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1035 do_per_step(step, nstglobalcomm) ||
1036 (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)));
1038 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1040 if (do_ene || do_log || bDoReplEx)
1042 bCalcVir = TRUE;
1043 bCalcEner = TRUE;
1044 bGStat = TRUE;
1047 force_flags = (GMX_FORCE_STATECHANGED |
1048 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1049 GMX_FORCE_ALLFORCES |
1050 GMX_FORCE_SEPLRF |
1051 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1052 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1053 (bDoFEP ? GMX_FORCE_DHDL : 0)
1056 if (fr->bTwinRange)
1058 if (do_per_step(step, ir->nstcalclr))
1060 force_flags |= GMX_FORCE_DO_LR;
1064 if (shellfc)
1066 /* Now is the time to relax the shells */
1067 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1068 ir, bNS, force_flags,
1069 top,
1070 constr, enerd, fcd,
1071 state, f, force_vir, mdatoms,
1072 nrnb, wcycle, graph, groups,
1073 shellfc, fr, bBornRadii, t, mu_tot,
1074 &bConverged, vsite,
1075 mdoutf_get_fp_field(outf));
1076 tcount += count;
1078 if (bConverged)
1080 nconverged++;
1083 else
1085 /* The coordinates (x) are shifted (to get whole molecules)
1086 * in do_force.
1087 * This is parallellized as well, and does communication too.
1088 * Check comments in sim_util.c
1090 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1091 state->box, state->x, &state->hist,
1092 f, force_vir, mdatoms, enerd, fcd,
1093 state->lambda, graph,
1094 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1095 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1098 if (bVV && !bStartingFromCpt && !bRerunMD)
1099 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1101 rvec *vbuf = NULL;
1103 wallcycle_start(wcycle, ewcUPDATE);
1104 if (ir->eI == eiVV && bInitStep)
1106 /* if using velocity verlet with full time step Ekin,
1107 * take the first half step only to compute the
1108 * virial for the first step. From there,
1109 * revert back to the initial coordinates
1110 * so that the input is actually the initial step.
1112 snew(vbuf, state->natoms);
1113 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1115 else
1117 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1118 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1121 /* If we are using twin-range interactions where the long-range component
1122 * is only evaluated every nstcalclr>1 steps, we should do a special update
1123 * step to combine the long-range forces on these steps.
1124 * For nstcalclr=1 this is not done, since the forces would have been added
1125 * directly to the short-range forces already.
1127 * TODO Remove various aspects of VV+twin-range in master
1128 * branch, because VV integrators did not ever support
1129 * twin-range multiple time stepping with constraints.
1131 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1133 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1134 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1135 ekind, M, upd, bInitStep, etrtVELOCITY1,
1136 cr, nrnb, constr, &top->idef);
1138 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1140 wallcycle_stop(wcycle, ewcUPDATE);
1141 update_constraints(fplog, step, NULL, ir, mdatoms,
1142 state, fr->bMolPBC, graph, f,
1143 &top->idef, shake_vir,
1144 cr, nrnb, wcycle, upd, constr,
1145 TRUE, bCalcVir);
1146 wallcycle_start(wcycle, ewcUPDATE);
1147 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1149 /* Correct the virial for multiple time stepping */
1150 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1153 else if (graph)
1155 /* Need to unshift here if a do_force has been
1156 called in the previous step */
1157 unshift_self(graph, state->box, state->x);
1159 /* if VV, compute the pressure and constraints */
1160 /* For VV2, we strictly only need this if using pressure
1161 * control, but we really would like to have accurate pressures
1162 * printed out.
1163 * Think about ways around this in the future?
1164 * For now, keep this choice in comments.
1166 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1167 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1168 bPres = TRUE;
1169 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1170 if (bCalcEner && ir->eI == eiVVAK)
1172 bSumEkinhOld = TRUE;
1174 /* for vv, the first half of the integration actually corresponds to the previous step.
1175 So we need information from the last step in the first half of the integration */
1176 if (bGStat || do_per_step(step-1, nstglobalcomm))
1178 wallcycle_stop(wcycle, ewcUPDATE);
1179 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1180 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1181 constr, NULL, FALSE, state->box,
1182 top_global, &bSumEkinhOld,
1183 (bGStat ? CGLO_GSTAT : 0)
1184 | CGLO_ENERGY
1185 | (bTemp ? CGLO_TEMPERATURE : 0)
1186 | (bPres ? CGLO_PRESSURE : 0)
1187 | (bPres ? CGLO_CONSTRAINT : 0)
1188 | (bStopCM ? CGLO_STOPCM : 0)
1189 | CGLO_SCALEEKIN
1191 /* explanation of above:
1192 a) We compute Ekin at the full time step
1193 if 1) we are using the AveVel Ekin, and it's not the
1194 initial step, or 2) if we are using AveEkin, but need the full
1195 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1196 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1197 EkinAveVel because it's needed for the pressure */
1198 wallcycle_start(wcycle, ewcUPDATE);
1200 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1201 if (!bInitStep)
1203 if (bTrotter)
1205 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1206 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1208 else
1210 if (bExchanged)
1212 wallcycle_stop(wcycle, ewcUPDATE);
1213 /* We need the kinetic energy at minus the half step for determining
1214 * the full step kinetic energy and possibly for T-coupling.*/
1215 /* This may not be quite working correctly yet . . . . */
1216 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1217 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1218 constr, NULL, FALSE, state->box,
1219 top_global, &bSumEkinhOld,
1220 CGLO_GSTAT | CGLO_TEMPERATURE);
1221 wallcycle_start(wcycle, ewcUPDATE);
1225 if (bTrotter && !bInitStep)
1227 copy_mat(shake_vir, state->svir_prev);
1228 copy_mat(force_vir, state->fvir_prev);
1229 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1231 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1232 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1233 enerd->term[F_EKIN] = trace(ekind->ekin);
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 (bVV)
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, f_global,
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 NPT control. Must restart */
1292 if (bStartingFromCpt && bVV)
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 #ifdef 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 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1368 if (EI_VV(ir->eI))
1370 if (!bInitStep)
1372 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1374 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1376 gmx_bool bIfRandomize;
1377 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1378 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1379 if (constr && bIfRandomize)
1381 update_constraints(fplog, step, NULL, ir, mdatoms,
1382 state, fr->bMolPBC, graph, f,
1383 &top->idef, tmp_vir,
1384 cr, nrnb, wcycle, upd, constr,
1385 TRUE, bCalcVir);
1389 /* ######### START SECOND UPDATE STEP ################# */
1390 /* Box is changed in update() when we do pressure coupling,
1391 * but we should still use the old box for energy corrections and when
1392 * writing it to the energy file, so it matches the trajectory files for
1393 * the same timestep above. Make a copy in a separate array.
1395 copy_mat(state->box, lastbox);
1397 dvdl_constr = 0;
1399 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1401 wallcycle_start(wcycle, ewcUPDATE);
1402 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1403 if (bTrotter)
1405 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1406 /* We can only do Berendsen coupling after we have summed
1407 * the kinetic energy or virial. Since the happens
1408 * in global_state after update, we should only do it at
1409 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1412 else
1414 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1415 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1418 if (bVV)
1420 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1422 /* velocity half-step update */
1423 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1424 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1425 ekind, M, upd, FALSE, etrtVELOCITY2,
1426 cr, nrnb, constr, &top->idef);
1429 /* Above, initialize just copies ekinh into ekin,
1430 * it doesn't copy position (for VV),
1431 * and entire integrator for MD.
1434 if (ir->eI == eiVVAK)
1436 /* We probably only need md->homenr, not state->natoms */
1437 if (state->natoms > cbuf_nalloc)
1439 cbuf_nalloc = state->natoms;
1440 srenew(cbuf, cbuf_nalloc);
1442 copy_rvecn(state->x, cbuf, 0, state->natoms);
1444 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1446 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1447 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1448 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1449 wallcycle_stop(wcycle, ewcUPDATE);
1451 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1452 fr->bMolPBC, graph, f,
1453 &top->idef, shake_vir,
1454 cr, nrnb, wcycle, upd, constr,
1455 FALSE, bCalcVir);
1457 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1459 /* Correct the virial for multiple time stepping */
1460 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1463 if (ir->eI == eiVVAK)
1465 /* erase F_EKIN and F_TEMP here? */
1466 /* just compute the kinetic energy at the half step to perform a trotter step */
1467 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1468 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1469 constr, NULL, FALSE, lastbox,
1470 top_global, &bSumEkinhOld,
1471 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1473 wallcycle_start(wcycle, ewcUPDATE);
1474 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1475 /* now we know the scaling, we can compute the positions again again */
1476 copy_rvecn(cbuf, state->x, 0, state->natoms);
1478 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1480 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1481 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1482 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1483 wallcycle_stop(wcycle, ewcUPDATE);
1485 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1486 /* are the small terms in the shake_vir here due
1487 * to numerical errors, or are they important
1488 * physically? I'm thinking they are just errors, but not completely sure.
1489 * For now, will call without actually constraining, constr=NULL*/
1490 update_constraints(fplog, step, NULL, ir, mdatoms,
1491 state, fr->bMolPBC, graph, f,
1492 &top->idef, tmp_vir,
1493 cr, nrnb, wcycle, upd, NULL,
1494 FALSE, bCalcVir);
1496 if (bVV)
1498 /* this factor or 2 correction is necessary
1499 because half of the constraint force is removed
1500 in the vv step, so we have to double it. See
1501 the Redmine issue #1255. It is not yet clear
1502 if the factor of 2 is exact, or just a very
1503 good approximation, and this will be
1504 investigated. The next step is to see if this
1505 can be done adding a dhdl contribution from the
1506 rattle step, but this is somewhat more
1507 complicated with the current code. Will be
1508 investigated, hopefully for 4.6.3. However,
1509 this current solution is much better than
1510 having it completely wrong.
1512 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1514 else
1516 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1519 else if (graph)
1521 /* Need to unshift here */
1522 unshift_self(graph, state->box, state->x);
1525 if (vsite != NULL)
1527 wallcycle_start(wcycle, ewcVSITECONSTR);
1528 if (graph != NULL)
1530 shift_self(graph, state->box, state->x);
1532 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1533 top->idef.iparams, top->idef.il,
1534 fr->ePBC, fr->bMolPBC, cr, state->box);
1536 if (graph != NULL)
1538 unshift_self(graph, state->box, state->x);
1540 wallcycle_stop(wcycle, ewcVSITECONSTR);
1543 /* ############## IF NOT VV, Calculate globals HERE ############ */
1544 /* With Leap-Frog we can skip compute_globals at
1545 * non-communication steps, but we need to calculate
1546 * the kinetic energy one step before communication.
1548 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1550 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1551 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1552 constr, &gs,
1553 (step_rel % gs.nstms == 0) &&
1554 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1555 lastbox,
1556 top_global, &bSumEkinhOld,
1557 (bGStat ? CGLO_GSTAT : 0)
1558 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1559 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1560 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1561 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1562 | CGLO_CONSTRAINT
1566 /* ############# END CALC EKIN AND PRESSURE ################# */
1568 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1569 the virial that should probably be addressed eventually. state->veta has better properies,
1570 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1571 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1573 if (ir->efep != efepNO && (!bVV || bRerunMD))
1575 /* Sum up the foreign energy and dhdl terms for md and sd.
1576 Currently done every step so that dhdl is correct in the .edr */
1577 sum_dhdl(enerd, state->lambda, ir->fepvals);
1579 update_box(fplog, step, ir, mdatoms, state, f,
1580 pcoupl_mu, nrnb, upd);
1582 /* ################# END UPDATE STEP 2 ################# */
1583 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1585 /* The coordinates (x) were unshifted in update */
1586 if (!bGStat)
1588 /* We will not sum ekinh_old,
1589 * so signal that we still have to do it.
1591 bSumEkinhOld = TRUE;
1594 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1596 /* use the directly determined last velocity, not actually the averaged half steps */
1597 if (bTrotter && ir->eI == eiVV)
1599 enerd->term[F_EKIN] = last_ekin;
1601 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1603 if (bVV)
1605 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1607 else
1609 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1611 /* ######### END PREPARING EDR OUTPUT ########### */
1613 /* Output stuff */
1614 if (MASTER(cr))
1616 gmx_bool do_dr, do_or;
1618 if (fplog && do_log && bDoExpanded)
1620 /* only needed if doing expanded ensemble */
1621 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1622 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1624 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1626 if (bCalcEner)
1628 upd_mdebin(mdebin, bDoDHDL, TRUE,
1629 t, mdatoms->tmass, enerd, state,
1630 ir->fepvals, ir->expandedvals, lastbox,
1631 shake_vir, force_vir, total_vir, pres,
1632 ekind, mu_tot, constr);
1634 else
1636 upd_mdebin_step(mdebin);
1639 do_dr = do_per_step(step, ir->nstdisreout);
1640 do_or = do_per_step(step, ir->nstorireout);
1642 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1643 step, t,
1644 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1646 if (ir->bPull)
1648 pull_print_output(ir->pull_work, step, t);
1651 if (do_per_step(step, ir->nstlog))
1653 if (fflush(fplog) != 0)
1655 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1659 if (bDoExpanded)
1661 /* Have to do this part _after_ outputting the logfile and the edr file */
1662 /* Gets written into the state at the beginning of next loop*/
1663 state->fep_state = lamnew;
1665 /* Print the remaining wall clock time for the run */
1666 if (MULTIMASTER(cr) &&
1667 (do_verbose || gmx_got_usr_signal()) &&
1668 !bPMETunePrinting)
1670 if (shellfc)
1672 fprintf(stderr, "\n");
1674 print_time(stderr, walltime_accounting, step, ir, cr);
1677 /* Ion/water position swapping.
1678 * Not done in last step since trajectory writing happens before this call
1679 * in the MD loop and exchanges would be lost anyway. */
1680 bNeedRepartition = FALSE;
1681 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1682 do_per_step(step, ir->swap->nstswap))
1684 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1685 bRerunMD ? rerun_fr.x : state->x,
1686 bRerunMD ? rerun_fr.box : state->box,
1687 top_global, MASTER(cr) && bVerbose, bRerunMD);
1689 if (bNeedRepartition && DOMAINDECOMP(cr))
1691 dd_collect_state(cr->dd, state, state_global);
1695 /* Replica exchange */
1696 bExchanged = FALSE;
1697 if (bDoReplEx)
1699 bExchanged = replica_exchange(fplog, cr, repl_ex,
1700 state_global, enerd,
1701 state, step, t);
1704 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1706 dd_partition_system(fplog, step, cr, TRUE, 1,
1707 state_global, top_global, ir,
1708 state, &f, mdatoms, top, fr,
1709 vsite, shellfc, constr,
1710 nrnb, wcycle, FALSE);
1713 bFirstStep = FALSE;
1714 bInitStep = FALSE;
1715 bStartingFromCpt = FALSE;
1717 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1718 /* With all integrators, except VV, we need to retain the pressure
1719 * at the current step for coupling at the next step.
1721 if ((state->flags & (1<<estPRES_PREV)) &&
1722 (bGStatEveryStep ||
1723 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1725 /* Store the pressure in t_state for pressure coupling
1726 * at the next MD step.
1728 copy_mat(pres, state->pres_prev);
1731 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1733 if ( (membed != NULL) && (!bLastStep) )
1735 rescale_membed(step_rel, membed, state_global->x);
1738 if (bRerunMD)
1740 if (MASTER(cr))
1742 /* read next frame from input trajectory */
1743 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1746 if (PAR(cr))
1748 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1752 cycles = wallcycle_stop(wcycle, ewcSTEP);
1753 if (DOMAINDECOMP(cr) && wcycle)
1755 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1758 if (!bRerunMD || !rerun_fr.bStep)
1760 /* increase the MD step number */
1761 step++;
1762 step_rel++;
1765 /* TODO make a counter-reset module */
1766 /* If it is time to reset counters, set a flag that remains
1767 true until counters actually get reset */
1768 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1769 gs.set[eglsRESETCOUNTERS] != 0)
1771 if (pme_loadbal_is_active(pme_loadbal))
1773 /* Do not permit counter reset while PME load
1774 * balancing is active. The only purpose for resetting
1775 * counters is to measure reliable performance data,
1776 * and that can't be done before balancing
1777 * completes.
1779 * TODO consider fixing this by delaying the reset
1780 * until after load balancing completes,
1781 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1782 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1783 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1784 "resetting counters later in the run, e.g. with gmx "
1785 "mdrun -resetstep.", step);
1787 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1788 use_GPU(fr->nbv) ? fr->nbv : NULL);
1789 wcycle_set_reset_counters(wcycle, -1);
1790 if (!(cr->duty & DUTY_PME))
1792 /* Tell our PME node to reset its counters */
1793 gmx_pme_send_resetcounters(cr, step);
1795 /* Correct max_hours for the elapsed time */
1796 max_hours -= elapsed_time/(60.0*60.0);
1797 /* If mdrun -maxh -resethway was active, it can only trigger once */
1798 bResetCountersHalfMaxH = FALSE; /* TODO move this to where gs.sig[eglsRESETCOUNTERS] is set */
1799 /* Reset can only happen once, so clear the triggering flag. */
1800 gs.set[eglsRESETCOUNTERS] = 0;
1803 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1804 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1807 /* End of main MD loop */
1808 debug_gmx();
1810 /* Closing TNG files can include compressing data. Therefore it is good to do that
1811 * before stopping the time measurements. */
1812 mdoutf_tng_close(outf);
1814 /* Stop measuring walltime */
1815 walltime_accounting_end(walltime_accounting);
1817 if (bRerunMD && MASTER(cr))
1819 close_trj(status);
1822 if (!(cr->duty & DUTY_PME))
1824 /* Tell the PME only node to finish */
1825 gmx_pme_send_finish(cr);
1828 if (MASTER(cr))
1830 if (ir->nstcalcenergy > 0 && !bRerunMD)
1832 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1833 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1837 done_mdoutf(outf);
1838 debug_gmx();
1840 if (bPMETune)
1842 pme_loadbal_done(pme_loadbal, cr, fplog, use_GPU(fr->nbv));
1845 if (shellfc && fplog)
1847 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1848 (nconverged*100.0)/step_rel);
1849 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1850 tcount/step_rel);
1853 if (repl_ex_nst > 0 && MASTER(cr))
1855 print_replica_exchange_statistics(fplog, repl_ex);
1858 /* IMD cleanup, if bIMD is TRUE. */
1859 IMD_finalize(ir->bIMD, ir->imd);
1861 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1863 return 0;