Add TNG writing and reading support
[gromacs.git] / src / programs / mdrun / md.c
blob6c0435ff591406796fa0938e3b72640416359d68
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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include "typedefs.h"
42 #include "smalloc.h"
43 #include "sysstuff.h"
44 #include "vec.h"
45 #include "vcm.h"
46 #include "mdebin.h"
47 #include "nrnb.h"
48 #include "calcmu.h"
49 #include "index.h"
50 #include "vsite.h"
51 #include "update.h"
52 #include "ns.h"
53 #include "mdrun.h"
54 #include "md_support.h"
55 #include "md_logging.h"
56 #include "network.h"
57 #include "xvgr.h"
58 #include "physics.h"
59 #include "names.h"
60 #include "force.h"
61 #include "disre.h"
62 #include "orires.h"
63 #include "pme.h"
64 #include "mdatoms.h"
65 #include "repl_ex.h"
66 #include "qmmm.h"
67 #include "domdec.h"
68 #include "domdec_network.h"
69 #include "partdec.h"
70 #include "coulomb.h"
71 #include "constr.h"
72 #include "shellfc.h"
73 #include "compute_io.h"
74 #include "mvdata.h"
75 #include "checkpoint.h"
76 #include "mtop_util.h"
77 #include "sighandler.h"
78 #include "txtdump.h"
79 #include "string2.h"
80 #include "pme_loadbal.h"
81 #include "bondf.h"
82 #include "membed.h"
83 #include "types/nlistheuristics.h"
84 #include "types/iteratedconstraints.h"
85 #include "nbnxn_cuda_data_mgmt.h"
87 #include "gromacs/utility/gmxmpi.h"
88 #include "gromacs/fileio/confio.h"
89 #include "gromacs/fileio/trajectory_writing.h"
90 #include "gromacs/fileio/trnio.h"
91 #include "gromacs/fileio/trxio.h"
92 #include "gromacs/fileio/xtcio.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/pulling/pull.h"
97 #ifdef GMX_FAHCORE
98 #include "corewrap.h"
99 #endif
101 static void reset_all_counters(FILE *fplog, t_commrec *cr,
102 gmx_int64_t step,
103 gmx_int64_t *step_rel, t_inputrec *ir,
104 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
105 gmx_walltime_accounting_t walltime_accounting,
106 nbnxn_cuda_ptr_t cu_nbv)
108 char sbuf[STEPSTRSIZE];
110 /* Reset all the counters related to performance over the run */
111 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
112 gmx_step_str(step, sbuf));
114 if (cu_nbv)
116 nbnxn_cuda_reset_timings(cu_nbv);
119 wallcycle_stop(wcycle, ewcRUN);
120 wallcycle_reset_all(wcycle);
121 if (DOMAINDECOMP(cr))
123 reset_dd_statistics_counters(cr->dd);
125 init_nrnb(nrnb);
126 ir->init_step += *step_rel;
127 ir->nsteps -= *step_rel;
128 *step_rel = 0;
129 wallcycle_start(wcycle, ewcRUN);
130 walltime_accounting_start(walltime_accounting);
131 print_date_and_time(fplog, cr->nodeid, "Restarted time", walltime_accounting);
134 double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
135 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
136 int nstglobalcomm,
137 gmx_vsite_t *vsite, gmx_constr_t constr,
138 int stepout, t_inputrec *ir,
139 gmx_mtop_t *top_global,
140 t_fcdata *fcd,
141 t_state *state_global,
142 t_mdatoms *mdatoms,
143 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
144 gmx_edsam_t ed, t_forcerec *fr,
145 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t membed,
146 real cpt_period, real max_hours,
147 const char gmx_unused *deviceOptions,
148 unsigned long Flags,
149 gmx_walltime_accounting_t walltime_accounting)
151 gmx_mdoutf_t outf = NULL;
152 gmx_int64_t step, step_rel;
153 double elapsed_time;
154 double t, t0, lam0[efptNR];
155 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
156 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
157 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
158 bBornRadii, bStartingFromCpt;
159 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
160 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
161 bForceUpdate = FALSE, bCPT;
162 gmx_bool bMasterState;
163 int force_flags, cglo_flags;
164 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
165 int i, m;
166 t_trxstatus *status;
167 rvec mu_tot;
168 t_vcm *vcm;
169 t_state *bufstate = NULL;
170 matrix *scale_tot, pcoupl_mu, M, ebox;
171 gmx_nlheur_t nlh;
172 t_trxframe rerun_fr;
173 gmx_repl_ex_t repl_ex = NULL;
174 int nchkpt = 1;
175 gmx_localtop_t *top;
176 t_mdebin *mdebin = NULL;
177 t_state *state = NULL;
178 rvec *f_global = NULL;
179 gmx_enerdata_t *enerd;
180 rvec *f = NULL;
181 gmx_global_stat_t gstat;
182 gmx_update_t upd = NULL;
183 t_graph *graph = NULL;
184 globsig_t gs;
185 gmx_rng_t mcrng = NULL;
186 gmx_groups_t *groups;
187 gmx_ekindata_t *ekind, *ekind_save;
188 gmx_shellfc_t shellfc;
189 int count, nconverged = 0;
190 real timestep = 0;
191 double tcount = 0;
192 gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged;
193 gmx_bool bAppend;
194 gmx_bool bResetCountersHalfMaxH = FALSE;
195 gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter;
196 gmx_bool bUpdateDoLR;
197 real dvdl_constr;
198 int a0, a1;
199 rvec *cbuf = NULL;
200 matrix lastbox;
201 real veta_save, scalevir, tracevir;
202 real vetanew = 0;
203 int lamnew = 0;
204 /* for FEP */
205 int nstfep;
206 double cycles;
207 real saved_conserved_quantity = 0;
208 real last_ekin = 0;
209 int iter_i;
210 t_extmass MassQ;
211 int **trotter_seq;
212 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
213 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
214 gmx_iterate_t iterate;
215 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
216 simulation stops. If equal to zero, don't
217 communicate any more between multisims.*/
218 /* PME load balancing data for GPU kernels */
219 pme_load_balancing_t pme_loadbal = NULL;
220 double cycles_pmes;
221 gmx_bool bPMETuneTry = FALSE, bPMETuneRunning = FALSE;
223 #ifdef GMX_FAHCORE
224 /* Temporary addition for FAHCORE checkpointing */
225 int chkpt_ret;
226 #endif
228 /* Check for special mdrun options */
229 bRerunMD = (Flags & MD_RERUN);
230 bAppend = (Flags & MD_APPENDFILES);
231 if (Flags & MD_RESETCOUNTERSHALFWAY)
233 if (ir->nsteps > 0)
235 /* Signal to reset the counters half the simulation steps. */
236 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
238 /* Signal to reset the counters halfway the simulation time. */
239 bResetCountersHalfMaxH = (max_hours > 0);
242 /* md-vv uses averaged full step velocities for T-control
243 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
244 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
245 bVV = EI_VV(ir->eI);
246 if (bVV) /* to store the initial velocities while computing virial */
248 snew(cbuf, top_global->natoms);
250 /* all the iteratative cases - only if there are constraints */
251 bIterativeCase = ((IR_NPH_TROTTER(ir) || IR_NPT_TROTTER(ir)) && (constr) && (!bRerunMD));
252 gmx_iterate_init(&iterate, FALSE); /* The default value of iterate->bIterationActive is set to
253 false in this step. The correct value, true or false,
254 is set at each step, as it depends on the frequency of temperature
255 and pressure control.*/
256 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
258 if (bRerunMD)
260 /* Since we don't know if the frames read are related in any way,
261 * rebuild the neighborlist at every step.
263 ir->nstlist = 1;
264 ir->nstcalcenergy = 1;
265 nstglobalcomm = 1;
268 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
270 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
271 bGStatEveryStep = (nstglobalcomm == 1);
273 if (!bGStatEveryStep && ir->nstlist == -1 && fplog != NULL)
275 fprintf(fplog,
276 "To reduce the energy communication with nstlist = -1\n"
277 "the neighbor list validity should not be checked at every step,\n"
278 "this means that exact integration is not guaranteed.\n"
279 "The neighbor list validity is checked after:\n"
280 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
281 "In most cases this will result in exact integration.\n"
282 "This reduces the energy communication by a factor of 2 to 3.\n"
283 "If you want less energy communication, set nstlist > 3.\n\n");
286 if (bRerunMD)
288 ir->nstxout_compressed = 0;
290 groups = &top_global->groups;
292 /* Initial values */
293 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
294 &(state_global->fep_state), lam0,
295 nrnb, top_global, &upd,
296 nfile, fnm, &outf, &mdebin,
297 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags);
299 clear_mat(total_vir);
300 clear_mat(pres);
301 /* Energy terms and groups */
302 snew(enerd, 1);
303 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
304 enerd);
305 if (DOMAINDECOMP(cr))
307 f = NULL;
309 else
311 snew(f, top_global->natoms);
314 /* Kinetic energy data */
315 snew(ekind, 1);
316 init_ekindata(fplog, top_global, &(ir->opts), ekind);
317 /* needed for iteration of constraints */
318 snew(ekind_save, 1);
319 init_ekindata(fplog, top_global, &(ir->opts), ekind_save);
320 /* Copy the cos acceleration to the groups struct */
321 ekind->cosacc.cos_accel = ir->cos_accel;
323 gstat = global_stat_init(ir);
324 debug_gmx();
326 /* Check for polarizable models and flexible constraints */
327 shellfc = init_shell_flexcon(fplog,
328 top_global, n_flexible_constraints(constr),
329 (ir->bContinuation ||
330 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
331 NULL : state_global->x);
333 if (DEFORM(*ir))
335 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
336 set_deform_reference_box(upd,
337 deform_init_init_step_tpx,
338 deform_init_box_tpx);
339 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
343 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
344 if ((io > 2000) && MASTER(cr))
346 fprintf(stderr,
347 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
348 io);
352 if (DOMAINDECOMP(cr))
354 top = dd_init_local_top(top_global);
356 snew(state, 1);
357 dd_init_local_state(cr->dd, state_global, state);
359 if (DDMASTER(cr->dd) && ir->nstfout)
361 snew(f_global, state_global->natoms);
364 else
366 if (PAR(cr))
368 /* Initialize the particle decomposition and split the topology */
369 top = split_system(fplog, top_global, ir, cr);
371 pd_cg_range(cr, &fr->cg0, &fr->hcg);
372 pd_at_range(cr, &a0, &a1);
374 else
376 top = gmx_mtop_generate_local_top(top_global, ir);
378 a0 = 0;
379 a1 = top_global->natoms;
382 forcerec_set_excl_load(fr, top, cr);
384 state = partdec_init_local_state(cr, state_global);
385 f_global = f;
387 atoms2md(top_global, ir, 0, NULL, a0, a1-a0, mdatoms);
389 if (vsite)
391 set_vsite_top(vsite, top, mdatoms, cr);
394 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
396 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
399 if (shellfc)
401 make_local_shells(cr, mdatoms, shellfc);
404 setup_bonded_threading(fr, &top->idef);
406 if (ir->pull && PAR(cr))
408 dd_make_local_pull_groups(NULL, ir->pull, mdatoms);
412 if (DOMAINDECOMP(cr))
414 /* Distribute the charge groups over the nodes from the master node */
415 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
416 state_global, top_global, ir,
417 state, &f, mdatoms, top, fr,
418 vsite, shellfc, constr,
419 nrnb, wcycle, FALSE);
423 update_mdatoms(mdatoms, state->lambda[efptMASS]);
425 if (opt2bSet("-cpi", nfile, fnm))
427 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
429 else
431 bStateFromCP = FALSE;
434 if (ir->bExpanded)
436 init_expanded_ensemble(bStateFromCP, ir, &mcrng, &state->dfhist);
439 if (MASTER(cr))
441 if (bStateFromCP)
443 /* Update mdebin with energy history if appending to output files */
444 if (Flags & MD_APPENDFILES)
446 restore_energyhistory_from_state(mdebin, &state_global->enerhist);
448 else
450 /* We might have read an energy history from checkpoint,
451 * free the allocated memory and reset the counts.
453 done_energyhistory(&state_global->enerhist);
454 init_energyhistory(&state_global->enerhist);
457 /* Set the initial energy history in state by updating once */
458 update_energyhistory(&state_global->enerhist, mdebin);
461 if ((state->flags & (1<<estLD_RNG)) && (Flags & MD_READ_RNG))
463 /* Set the random state if we read a checkpoint file */
464 set_stochd_state(upd, state);
467 if (state->flags & (1<<estMC_RNG))
469 set_mc_state(mcrng, state);
472 /* Initialize constraints */
473 if (constr)
475 if (!DOMAINDECOMP(cr))
477 set_constraints(constr, top, ir, mdatoms, cr);
481 if (repl_ex_nst > 0)
483 /* We need to be sure replica exchange can only occur
484 * when the energies are current */
485 check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy,
486 "repl_ex_nst", &repl_ex_nst);
487 /* This check needs to happen before inter-simulation
488 * signals are initialized, too */
490 if (repl_ex_nst > 0 && MASTER(cr))
492 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
493 repl_ex_nst, repl_ex_nex, repl_ex_seed);
496 /* PME tuning is only supported with GPUs or PME nodes and not with rerun or LJ-PME.
497 * With perturbed charges with soft-core we should not change the cut-off.
499 if ((Flags & MD_TUNEPME) &&
500 EEL_PME(fr->eeltype) &&
501 ( (fr->cutoff_scheme == ecutsVERLET && fr->nbv->bUseGPU) || !(cr->duty & DUTY_PME)) &&
502 !(ir->efep != efepNO && mdatoms->nChargePerturbed > 0 && ir->fepvals->bScCoul) &&
503 !bRerunMD && !EVDW_PME(fr->vdwtype))
505 pme_loadbal_init(&pme_loadbal, ir, state->box, fr->ic, fr->pmedata);
506 cycles_pmes = 0;
507 if (cr->duty & DUTY_PME)
509 /* Start tuning right away, as we can't measure the load */
510 bPMETuneRunning = TRUE;
512 else
514 /* Separate PME nodes, we can measure the PP/PME load balance */
515 bPMETuneTry = TRUE;
519 if (!ir->bContinuation && !bRerunMD)
521 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
523 /* Set the velocities of frozen particles to zero */
524 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
526 for (m = 0; m < DIM; m++)
528 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
530 state->v[i][m] = 0;
536 if (constr)
538 /* Constrain the initial coordinates and velocities */
539 do_constrain_first(fplog, constr, ir, mdatoms, state,
540 cr, nrnb, fr, top);
542 if (vsite)
544 /* Construct the virtual sites for the initial configuration */
545 construct_vsites(vsite, state->x, ir->delta_t, NULL,
546 top->idef.iparams, top->idef.il,
547 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
551 debug_gmx();
553 /* set free energy calculation frequency as the minimum
554 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
555 nstfep = ir->fepvals->nstdhdl;
556 if (ir->bExpanded)
558 nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl, nstfep);
560 if (repl_ex_nst > 0)
562 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
565 /* I'm assuming we need global communication the first time! MRS */
566 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
567 | ((ir->comm_mode != ecmNO) ? CGLO_STOPCM : 0)
568 | (bVV ? CGLO_PRESSURE : 0)
569 | (bVV ? CGLO_CONSTRAINT : 0)
570 | (bRerunMD ? CGLO_RERUNMD : 0)
571 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
573 bSumEkinhOld = FALSE;
574 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
575 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
576 constr, NULL, FALSE, state->box,
577 top_global, &bSumEkinhOld, cglo_flags);
578 if (ir->eI == eiVVAK)
580 /* a second call to get the half step temperature initialized as well */
581 /* we do the same call as above, but turn the pressure off -- internally to
582 compute_globals, this is recognized as a velocity verlet half-step
583 kinetic energy calculation. This minimized excess variables, but
584 perhaps loses some logic?*/
586 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
587 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
588 constr, NULL, FALSE, state->box,
589 top_global, &bSumEkinhOld,
590 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
593 /* Calculate the initial half step temperature, and save the ekinh_old */
594 if (!(Flags & MD_STARTFROMCPT))
596 for (i = 0; (i < ir->opts.ngtc); i++)
598 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
601 if (ir->eI != eiVV)
603 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
604 and there is no previous step */
607 /* if using an iterative algorithm, we need to create a working directory for the state. */
608 if (bIterativeCase)
610 bufstate = init_bufstate(state);
613 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
614 temperature control */
615 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
617 if (MASTER(cr))
619 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
621 fprintf(fplog,
622 "RMS relative constraint deviation after constraining: %.2e\n",
623 constr_rmsd(constr, FALSE));
625 if (EI_STATE_VELOCITY(ir->eI))
627 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
629 if (bRerunMD)
631 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
632 " input trajectory '%s'\n\n",
633 *(top_global->name), opt2fn("-rerun", nfile, fnm));
634 if (bVerbose)
636 fprintf(stderr, "Calculated time to finish depends on nsteps from "
637 "run input file,\nwhich may not correspond to the time "
638 "needed to process input trajectory.\n\n");
641 else
643 char tbuf[20];
644 fprintf(stderr, "starting mdrun '%s'\n",
645 *(top_global->name));
646 if (ir->nsteps >= 0)
648 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
650 else
652 sprintf(tbuf, "%s", "infinite");
654 if (ir->init_step > 0)
656 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
657 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
658 gmx_step_str(ir->init_step, sbuf2),
659 ir->init_step*ir->delta_t);
661 else
663 fprintf(stderr, "%s steps, %s ps.\n",
664 gmx_step_str(ir->nsteps, sbuf), tbuf);
667 fprintf(fplog, "\n");
670 /* Set and write start time */
671 walltime_accounting_start(walltime_accounting);
672 print_date_and_time(fplog, cr->nodeid, "Started mdrun", walltime_accounting);
673 wallcycle_start(wcycle, ewcRUN);
674 if (fplog)
676 fprintf(fplog, "\n");
679 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
680 #ifdef GMX_FAHCORE
681 chkpt_ret = fcCheckPointParallel( cr->nodeid,
682 NULL, 0);
683 if (chkpt_ret == 0)
685 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
687 #endif
689 debug_gmx();
690 /***********************************************************
692 * Loop over MD steps
694 ************************************************************/
696 /* if rerunMD then read coordinates and velocities from input trajectory */
697 if (bRerunMD)
699 if (getenv("GMX_FORCE_UPDATE"))
701 bForceUpdate = TRUE;
704 rerun_fr.natoms = 0;
705 if (MASTER(cr))
707 bNotLastFrame = read_first_frame(oenv, &status,
708 opt2fn("-rerun", nfile, fnm),
709 &rerun_fr, TRX_NEED_X | TRX_READ_V);
710 if (rerun_fr.natoms != top_global->natoms)
712 gmx_fatal(FARGS,
713 "Number of atoms in trajectory (%d) does not match the "
714 "run input file (%d)\n",
715 rerun_fr.natoms, top_global->natoms);
717 if (ir->ePBC != epbcNONE)
719 if (!rerun_fr.bBox)
721 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);
723 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
725 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
730 if (PAR(cr))
732 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
735 if (ir->ePBC != epbcNONE)
737 /* Set the shift vectors.
738 * Necessary here when have a static box different from the tpr box.
740 calc_shifts(rerun_fr.box, fr->shift_vec);
744 /* loop over MD steps or if rerunMD to end of input trajectory */
745 bFirstStep = TRUE;
746 /* Skip the first Nose-Hoover integration when we get the state from tpx */
747 bStateFromTPX = !bStateFromCP;
748 bInitStep = bFirstStep && (bStateFromTPX || bVV);
749 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
750 bLastStep = FALSE;
751 bSumEkinhOld = FALSE;
752 bExchanged = FALSE;
754 init_global_signals(&gs, cr, ir, repl_ex_nst);
756 step = ir->init_step;
757 step_rel = 0;
759 if (ir->nstlist == -1)
761 init_nlistheuristics(&nlh, bGStatEveryStep, step);
764 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
766 /* check how many steps are left in other sims */
767 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
771 /* and stop now if we should */
772 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
773 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
774 while (!bLastStep || (bRerunMD && bNotLastFrame))
777 wallcycle_start(wcycle, ewcSTEP);
779 if (bRerunMD)
781 if (rerun_fr.bStep)
783 step = rerun_fr.step;
784 step_rel = step - ir->init_step;
786 if (rerun_fr.bTime)
788 t = rerun_fr.time;
790 else
792 t = step;
795 else
797 bLastStep = (step_rel == ir->nsteps);
798 t = t0 + step*ir->delta_t;
801 if (ir->efep != efepNO || ir->bSimTemp)
803 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
804 requiring different logic. */
806 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
807 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
808 bDoFEP = (do_per_step(step, nstfep) && (ir->efep != efepNO));
809 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
810 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
813 if (bSimAnn)
815 update_annealing_target_temp(&(ir->opts), t);
818 if (bRerunMD)
820 if (!(DOMAINDECOMP(cr) && !MASTER(cr)))
822 for (i = 0; i < state_global->natoms; i++)
824 copy_rvec(rerun_fr.x[i], state_global->x[i]);
826 if (rerun_fr.bV)
828 for (i = 0; i < state_global->natoms; i++)
830 copy_rvec(rerun_fr.v[i], state_global->v[i]);
833 else
835 for (i = 0; i < state_global->natoms; i++)
837 clear_rvec(state_global->v[i]);
839 if (bRerunWarnNoV)
841 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
842 " Ekin, temperature and pressure are incorrect,\n"
843 " the virial will be incorrect when constraints are present.\n"
844 "\n");
845 bRerunWarnNoV = FALSE;
849 copy_mat(rerun_fr.box, state_global->box);
850 copy_mat(state_global->box, state->box);
852 if (vsite && (Flags & MD_RERUN_VSITE))
854 if (DOMAINDECOMP(cr))
856 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
858 if (graph)
860 /* Following is necessary because the graph may get out of sync
861 * with the coordinates if we only have every N'th coordinate set
863 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
864 shift_self(graph, state->box, state->x);
866 construct_vsites(vsite, state->x, ir->delta_t, state->v,
867 top->idef.iparams, top->idef.il,
868 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
869 if (graph)
871 unshift_self(graph, state->box, state->x);
876 /* Stop Center of Mass motion */
877 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
879 if (bRerunMD)
881 /* for rerun MD always do Neighbour Searching */
882 bNS = (bFirstStep || ir->nstlist != 0);
883 bNStList = bNS;
885 else
887 /* Determine whether or not to do Neighbour Searching and LR */
888 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
890 bNS = (bFirstStep || bExchanged || bNStList || bDoFEP ||
891 (ir->nstlist == -1 && nlh.nabnsb > 0));
893 if (bNS && ir->nstlist == -1)
895 set_nlistheuristics(&nlh, bFirstStep || bExchanged || bDoFEP, step);
899 /* check whether we should stop because another simulation has
900 stopped. */
901 if (MULTISIM(cr))
903 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
904 (multisim_nsteps != ir->nsteps) )
906 if (bNS)
908 if (MASTER(cr))
910 fprintf(stderr,
911 "Stopping simulation %d because another one has finished\n",
912 cr->ms->sim);
914 bLastStep = TRUE;
915 gs.sig[eglsCHKPT] = 1;
920 /* < 0 means stop at next step, > 0 means stop at next NS step */
921 if ( (gs.set[eglsSTOPCOND] < 0) ||
922 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
924 bLastStep = TRUE;
927 /* Determine whether or not to update the Born radii if doing GB */
928 bBornRadii = bFirstStep;
929 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
931 bBornRadii = TRUE;
934 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
935 do_verbose = bVerbose &&
936 (step % stepout == 0 || bFirstStep || bLastStep);
938 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
940 if (bRerunMD)
942 bMasterState = TRUE;
944 else
946 bMasterState = FALSE;
947 /* Correct the new box if it is too skewed */
948 if (DYNAMIC_BOX(*ir))
950 if (correct_box(fplog, step, state->box, graph))
952 bMasterState = TRUE;
955 if (DOMAINDECOMP(cr) && bMasterState)
957 dd_collect_state(cr->dd, state, state_global);
961 if (DOMAINDECOMP(cr))
963 /* Repartition the domain decomposition */
964 wallcycle_start(wcycle, ewcDOMDEC);
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 && !bPMETuneRunning);
972 wallcycle_stop(wcycle, ewcDOMDEC);
973 /* If using an iterative integrator, reallocate space to match the decomposition */
977 if (MASTER(cr) && do_log)
979 print_ebin_header(fplog, step, t, state->lambda[efptFEP]); /* can we improve the information printed here? */
982 if (ir->efep != efepNO)
984 update_mdatoms(mdatoms, state->lambda[efptMASS]);
987 if ((bRerunMD && rerun_fr.bV) || bExchanged)
990 /* We need the kinetic energy at minus the half step for determining
991 * the full step kinetic energy and possibly for T-coupling.*/
992 /* This may not be quite working correctly yet . . . . */
993 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
994 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
995 constr, NULL, FALSE, state->box,
996 top_global, &bSumEkinhOld,
997 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
999 clear_mat(force_vir);
1001 /* We write a checkpoint at this MD step when:
1002 * either at an NS step when we signalled through gs,
1003 * or at the last step (but not when we do not want confout),
1004 * but never at the first step or with rerun.
1006 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1007 (bLastStep && (Flags & MD_CONFOUT))) &&
1008 step > ir->init_step && !bRerunMD);
1009 if (bCPT)
1011 gs.set[eglsCHKPT] = 0;
1014 /* Determine the energy and pressure:
1015 * at nstcalcenergy steps and at energy output steps (set below).
1017 if (EI_VV(ir->eI) && (!bInitStep))
1019 /* for vv, the first half of the integration actually corresponds
1020 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1021 but the virial needs to be calculated on both the current step and the 'next' step. Future
1022 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1024 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1025 bCalcVir = bCalcEner ||
1026 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1028 else
1030 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1031 bCalcVir = bCalcEner ||
1032 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1035 /* Do we need global communication ? */
1036 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1037 do_per_step(step, nstglobalcomm) || (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)) ||
1038 (ir->nstlist == -1 && !bRerunMD && step >= nlh.step_nscheck));
1040 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1042 if (do_ene || do_log)
1044 bCalcVir = TRUE;
1045 bCalcEner = TRUE;
1046 bGStat = TRUE;
1049 /* these CGLO_ options remain the same throughout the iteration */
1050 cglo_flags = ((bRerunMD ? CGLO_RERUNMD : 0) |
1051 (bGStat ? CGLO_GSTAT : 0)
1054 force_flags = (GMX_FORCE_STATECHANGED |
1055 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1056 GMX_FORCE_ALLFORCES |
1057 GMX_FORCE_SEPLRF |
1058 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1059 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1060 (bDoFEP ? GMX_FORCE_DHDL : 0)
1063 if (fr->bTwinRange)
1065 if (do_per_step(step, ir->nstcalclr))
1067 force_flags |= GMX_FORCE_DO_LR;
1071 if (shellfc)
1073 /* Now is the time to relax the shells */
1074 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1075 ir, bNS, force_flags,
1076 top,
1077 constr, enerd, fcd,
1078 state, f, force_vir, mdatoms,
1079 nrnb, wcycle, graph, groups,
1080 shellfc, fr, bBornRadii, t, mu_tot,
1081 &bConverged, vsite,
1082 mdoutf_get_fp_field(outf));
1083 tcount += count;
1085 if (bConverged)
1087 nconverged++;
1090 else
1092 /* The coordinates (x) are shifted (to get whole molecules)
1093 * in do_force.
1094 * This is parallellized as well, and does communication too.
1095 * Check comments in sim_util.c
1097 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1098 state->box, state->x, &state->hist,
1099 f, force_vir, mdatoms, enerd, fcd,
1100 state->lambda, graph,
1101 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1102 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1105 if (bVV && !bStartingFromCpt && !bRerunMD)
1106 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1108 if (ir->eI == eiVV && bInitStep)
1110 /* if using velocity verlet with full time step Ekin,
1111 * take the first half step only to compute the
1112 * virial for the first step. From there,
1113 * revert back to the initial coordinates
1114 * so that the input is actually the initial step.
1116 copy_rvecn(state->v, cbuf, 0, state->natoms); /* should make this better for parallelizing? */
1118 else
1120 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1121 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1124 /* If we are using twin-range interactions where the long-range component
1125 * is only evaluated every nstcalclr>1 steps, we should do a special update
1126 * step to combine the long-range forces on these steps.
1127 * For nstcalclr=1 this is not done, since the forces would have been added
1128 * directly to the short-range forces already.
1130 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1132 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1133 f, bUpdateDoLR, fr->f_twin, fcd,
1134 ekind, M, upd, bInitStep, etrtVELOCITY1,
1135 cr, nrnb, constr, &top->idef);
1137 if (bIterativeCase && do_per_step(step-1, ir->nstpcouple) && !bInitStep)
1139 gmx_iterate_init(&iterate, TRUE);
1141 /* for iterations, we save these vectors, as we will be self-consistently iterating
1142 the calculations */
1144 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1146 /* save the state */
1147 if (iterate.bIterationActive)
1149 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1152 bFirstIterate = TRUE;
1153 while (bFirstIterate || iterate.bIterationActive)
1155 if (iterate.bIterationActive)
1157 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1158 if (bFirstIterate && bTrotter)
1160 /* The first time through, we need a decent first estimate
1161 of veta(t+dt) to compute the constraints. Do
1162 this by computing the box volume part of the
1163 trotter integration at this time. Nothing else
1164 should be changed by this routine here. If
1165 !(first time), we start with the previous value
1166 of veta. */
1168 veta_save = state->veta;
1169 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ0);
1170 vetanew = state->veta;
1171 state->veta = veta_save;
1175 bOK = TRUE;
1176 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1178 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1179 state, fr->bMolPBC, graph, f,
1180 &top->idef, shake_vir,
1181 cr, nrnb, wcycle, upd, constr,
1182 TRUE, bCalcVir, vetanew);
1184 if (!bOK)
1186 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1190 else if (graph)
1192 /* Need to unshift here if a do_force has been
1193 called in the previous step */
1194 unshift_self(graph, state->box, state->x);
1197 /* if VV, compute the pressure and constraints */
1198 /* For VV2, we strictly only need this if using pressure
1199 * control, but we really would like to have accurate pressures
1200 * printed out.
1201 * Think about ways around this in the future?
1202 * For now, keep this choice in comments.
1204 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1205 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1206 bPres = TRUE;
1207 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1208 if (bCalcEner && ir->eI == eiVVAK) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1210 bSumEkinhOld = TRUE;
1212 /* for vv, the first half of the integration actually corresponds to the previous step.
1213 So we need information from the last step in the first half of the integration */
1214 if (bGStat || do_per_step(step-1, nstglobalcomm))
1216 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1217 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1218 constr, NULL, FALSE, state->box,
1219 top_global, &bSumEkinhOld,
1220 cglo_flags
1221 | CGLO_ENERGY
1222 | (bTemp ? CGLO_TEMPERATURE : 0)
1223 | (bPres ? CGLO_PRESSURE : 0)
1224 | (bPres ? CGLO_CONSTRAINT : 0)
1225 | ((iterate.bIterationActive) ? CGLO_ITERATE : 0)
1226 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1227 | CGLO_SCALEEKIN
1229 /* explanation of above:
1230 a) We compute Ekin at the full time step
1231 if 1) we are using the AveVel Ekin, and it's not the
1232 initial step, or 2) if we are using AveEkin, but need the full
1233 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1234 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1235 EkinAveVel because it's needed for the pressure */
1237 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1238 if (!bInitStep)
1240 if (bTrotter)
1242 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1243 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1245 else
1247 if (bExchanged)
1250 /* We need the kinetic energy at minus the half step for determining
1251 * the full step kinetic energy and possibly for T-coupling.*/
1252 /* This may not be quite working correctly yet . . . . */
1253 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1254 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1255 constr, NULL, FALSE, state->box,
1256 top_global, &bSumEkinhOld,
1257 CGLO_RERUNMD | CGLO_GSTAT | CGLO_TEMPERATURE);
1262 if (iterate.bIterationActive &&
1263 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1264 state->veta, &vetanew))
1266 break;
1268 bFirstIterate = FALSE;
1271 if (bTrotter && !bInitStep)
1273 copy_mat(shake_vir, state->svir_prev);
1274 copy_mat(force_vir, state->fvir_prev);
1275 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1277 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1278 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1279 enerd->term[F_EKIN] = trace(ekind->ekin);
1282 /* if it's the initial step, we performed this first step just to get the constraint virial */
1283 if (bInitStep && ir->eI == eiVV)
1285 copy_rvecn(cbuf, state->v, 0, state->natoms);
1289 /* MRS -- now done iterating -- compute the conserved quantity */
1290 if (bVV)
1292 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1293 if (ir->eI == eiVV)
1295 last_ekin = enerd->term[F_EKIN];
1297 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1299 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1301 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1302 if (!bRerunMD)
1304 sum_dhdl(enerd, state->lambda, ir->fepvals);
1308 /* ######## END FIRST UPDATE STEP ############## */
1309 /* ######## If doing VV, we now have v(dt) ###### */
1310 if (bDoExpanded)
1312 /* perform extended ensemble sampling in lambda - we don't
1313 actually move to the new state before outputting
1314 statistics, but if performing simulated tempering, we
1315 do update the velocities and the tau_t. */
1317 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, mcrng, state->v, mdatoms);
1318 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1319 copy_df_history(&state_global->dfhist, &state->dfhist);
1322 /* Now we have the energies and forces corresponding to the
1323 * coordinates at time t. We must output all of this before
1324 * the update.
1326 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1327 ir, state, state_global, top_global, fr, upd,
1328 outf, mdebin, ekind, f, f_global,
1329 wcycle, mcrng, &nchkpt,
1330 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1331 bSumEkinhOld);
1333 /* kludge -- virial is lost with restart for NPT control. Must restart */
1334 if (bStartingFromCpt && bVV)
1336 copy_mat(state->svir_prev, shake_vir);
1337 copy_mat(state->fvir_prev, force_vir);
1340 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1342 /* Check whether everything is still allright */
1343 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1344 #ifdef GMX_THREAD_MPI
1345 && MASTER(cr)
1346 #endif
1349 /* this is just make gs.sig compatible with the hack
1350 of sending signals around by MPI_Reduce with together with
1351 other floats */
1352 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1354 gs.sig[eglsSTOPCOND] = 1;
1356 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1358 gs.sig[eglsSTOPCOND] = -1;
1360 /* < 0 means stop at next step, > 0 means stop at next NS step */
1361 if (fplog)
1363 fprintf(fplog,
1364 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1365 gmx_get_signal_name(),
1366 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1367 fflush(fplog);
1369 fprintf(stderr,
1370 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1371 gmx_get_signal_name(),
1372 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1373 fflush(stderr);
1374 handled_stop_condition = (int)gmx_get_stop_condition();
1376 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1377 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1378 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1380 /* Signal to terminate the run */
1381 gs.sig[eglsSTOPCOND] = 1;
1382 if (fplog)
1384 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1386 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1389 if (bResetCountersHalfMaxH && MASTER(cr) &&
1390 elapsed_time > max_hours*60.0*60.0*0.495)
1392 gs.sig[eglsRESETCOUNTERS] = 1;
1395 if (ir->nstlist == -1 && !bRerunMD)
1397 /* When bGStatEveryStep=FALSE, global_stat is only called
1398 * when we check the atom displacements, not at NS steps.
1399 * This means that also the bonded interaction count check is not
1400 * performed immediately after NS. Therefore a few MD steps could
1401 * be performed with missing interactions.
1402 * But wrong energies are never written to file,
1403 * since energies are only written after global_stat
1404 * has been called.
1406 if (step >= nlh.step_nscheck)
1408 nlh.nabnsb = natoms_beyond_ns_buffer(ir, fr, &top->cgs,
1409 nlh.scale_tot, state->x);
1411 else
1413 /* This is not necessarily true,
1414 * but step_nscheck is determined quite conservatively.
1416 nlh.nabnsb = 0;
1420 /* In parallel we only have to check for checkpointing in steps
1421 * where we do global communication,
1422 * otherwise the other nodes don't know.
1424 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1425 cpt_period >= 0 &&
1426 (cpt_period == 0 ||
1427 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1428 gs.set[eglsCHKPT] == 0)
1430 gs.sig[eglsCHKPT] = 1;
1433 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1434 if (EI_VV(ir->eI))
1436 if (!bInitStep)
1438 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1440 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1442 gmx_bool bIfRandomize;
1443 bIfRandomize = update_randomize_velocities(ir, step, mdatoms, state, upd, &top->idef, constr);
1444 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1445 if (constr && bIfRandomize)
1447 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1448 state, fr->bMolPBC, graph, f,
1449 &top->idef, tmp_vir,
1450 cr, nrnb, wcycle, upd, constr,
1451 TRUE, bCalcVir, vetanew);
1456 if (bIterativeCase && do_per_step(step, ir->nstpcouple))
1458 gmx_iterate_init(&iterate, TRUE);
1459 /* for iterations, we save these vectors, as we will be redoing the calculations */
1460 copy_coupling_state(state, bufstate, ekind, ekind_save, &(ir->opts));
1463 bFirstIterate = TRUE;
1464 while (bFirstIterate || iterate.bIterationActive)
1466 /* We now restore these vectors to redo the calculation with improved extended variables */
1467 if (iterate.bIterationActive)
1469 copy_coupling_state(bufstate, state, ekind_save, ekind, &(ir->opts));
1472 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1473 so scroll down for that logic */
1475 /* ######### START SECOND UPDATE STEP ################# */
1476 /* Box is changed in update() when we do pressure coupling,
1477 * but we should still use the old box for energy corrections and when
1478 * writing it to the energy file, so it matches the trajectory files for
1479 * the same timestep above. Make a copy in a separate array.
1481 copy_mat(state->box, lastbox);
1483 bOK = TRUE;
1484 dvdl_constr = 0;
1486 if (!(bRerunMD && !rerun_fr.bV && !bForceUpdate))
1488 wallcycle_start(wcycle, ewcUPDATE);
1489 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1490 if (bTrotter)
1492 if (iterate.bIterationActive)
1494 if (bFirstIterate)
1496 scalevir = 1;
1498 else
1500 /* we use a new value of scalevir to converge the iterations faster */
1501 scalevir = tracevir/trace(shake_vir);
1503 msmul(shake_vir, scalevir, shake_vir);
1504 m_add(force_vir, shake_vir, total_vir);
1505 clear_mat(shake_vir);
1507 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1508 /* We can only do Berendsen coupling after we have summed
1509 * the kinetic energy or virial. Since the happens
1510 * in global_state after update, we should only do it at
1511 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1514 else
1516 update_tcouple(step, ir, state, ekind, upd, &MassQ, mdatoms);
1517 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1520 if (bVV)
1522 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1524 /* velocity half-step update */
1525 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1526 bUpdateDoLR, fr->f_twin, fcd,
1527 ekind, M, upd, FALSE, etrtVELOCITY2,
1528 cr, nrnb, constr, &top->idef);
1531 /* Above, initialize just copies ekinh into ekin,
1532 * it doesn't copy position (for VV),
1533 * and entire integrator for MD.
1536 if (ir->eI == eiVVAK)
1538 copy_rvecn(state->x, cbuf, 0, state->natoms);
1540 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1542 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1543 bUpdateDoLR, fr->f_twin, fcd,
1544 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1545 wallcycle_stop(wcycle, ewcUPDATE);
1547 update_constraints(fplog, step, &dvdl_constr, ir, ekind, mdatoms, state,
1548 fr->bMolPBC, graph, f,
1549 &top->idef, shake_vir,
1550 cr, nrnb, wcycle, upd, constr,
1551 FALSE, bCalcVir, state->veta);
1553 if (ir->eI == eiVVAK)
1555 /* erase F_EKIN and F_TEMP here? */
1556 /* just compute the kinetic energy at the half step to perform a trotter step */
1557 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1558 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1559 constr, NULL, FALSE, lastbox,
1560 top_global, &bSumEkinhOld,
1561 cglo_flags | CGLO_TEMPERATURE
1563 wallcycle_start(wcycle, ewcUPDATE);
1564 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1565 /* now we know the scaling, we can compute the positions again again */
1566 copy_rvecn(cbuf, state->x, 0, state->natoms);
1568 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1570 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1571 bUpdateDoLR, fr->f_twin, fcd,
1572 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1573 wallcycle_stop(wcycle, ewcUPDATE);
1575 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1576 /* are the small terms in the shake_vir here due
1577 * to numerical errors, or are they important
1578 * physically? I'm thinking they are just errors, but not completely sure.
1579 * For now, will call without actually constraining, constr=NULL*/
1580 update_constraints(fplog, step, NULL, ir, ekind, mdatoms,
1581 state, fr->bMolPBC, graph, f,
1582 &top->idef, tmp_vir,
1583 cr, nrnb, wcycle, upd, NULL,
1584 FALSE, bCalcVir,
1585 state->veta);
1587 if (!bOK)
1589 gmx_fatal(FARGS, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1592 if (fr->bSepDVDL && fplog && do_log)
1594 gmx_print_sepdvdl(fplog, "Constraint dV/dl", 0.0, dvdl_constr);
1596 if (bVV)
1598 /* this factor or 2 correction is necessary
1599 because half of the constraint force is removed
1600 in the vv step, so we have to double it. See
1601 the Redmine issue #1255. It is not yet clear
1602 if the factor of 2 is exact, or just a very
1603 good approximation, and this will be
1604 investigated. The next step is to see if this
1605 can be done adding a dhdl contribution from the
1606 rattle step, but this is somewhat more
1607 complicated with the current code. Will be
1608 investigated, hopefully for 4.6.3. However,
1609 this current solution is much better than
1610 having it completely wrong.
1612 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1614 else
1616 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1619 else if (graph)
1621 /* Need to unshift here */
1622 unshift_self(graph, state->box, state->x);
1625 if (vsite != NULL)
1627 wallcycle_start(wcycle, ewcVSITECONSTR);
1628 if (graph != NULL)
1630 shift_self(graph, state->box, state->x);
1632 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1633 top->idef.iparams, top->idef.il,
1634 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1636 if (graph != NULL)
1638 unshift_self(graph, state->box, state->x);
1640 wallcycle_stop(wcycle, ewcVSITECONSTR);
1643 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1644 /* With Leap-Frog we can skip compute_globals at
1645 * non-communication steps, but we need to calculate
1646 * the kinetic energy one step before communication.
1648 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1650 if (ir->nstlist == -1 && bFirstIterate)
1652 gs.sig[eglsNABNSB] = nlh.nabnsb;
1654 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1655 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1656 constr,
1657 bFirstIterate ? &gs : NULL,
1658 (step_rel % gs.nstms == 0) &&
1659 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1660 lastbox,
1661 top_global, &bSumEkinhOld,
1662 cglo_flags
1663 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1664 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1665 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1666 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1667 | (iterate.bIterationActive ? CGLO_ITERATE : 0)
1668 | (bFirstIterate ? CGLO_FIRSTITERATE : 0)
1669 | CGLO_CONSTRAINT
1671 if (ir->nstlist == -1 && bFirstIterate)
1673 nlh.nabnsb = gs.set[eglsNABNSB];
1674 gs.set[eglsNABNSB] = 0;
1677 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1678 /* ############# END CALC EKIN AND PRESSURE ################# */
1680 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1681 the virial that should probably be addressed eventually. state->veta has better properies,
1682 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1683 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1685 if (iterate.bIterationActive &&
1686 done_iterating(cr, fplog, step, &iterate, bFirstIterate,
1687 trace(shake_vir), &tracevir))
1689 break;
1691 bFirstIterate = FALSE;
1694 if (!bVV || bRerunMD)
1696 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1697 sum_dhdl(enerd, state->lambda, ir->fepvals);
1699 update_box(fplog, step, ir, mdatoms, state, f,
1700 ir->nstlist == -1 ? &nlh.scale_tot : NULL, pcoupl_mu, nrnb, upd);
1702 /* ################# END UPDATE STEP 2 ################# */
1703 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1705 /* The coordinates (x) were unshifted in update */
1706 if (!bGStat)
1708 /* We will not sum ekinh_old,
1709 * so signal that we still have to do it.
1711 bSumEkinhOld = TRUE;
1714 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1716 /* use the directly determined last velocity, not actually the averaged half steps */
1717 if (bTrotter && ir->eI == eiVV)
1719 enerd->term[F_EKIN] = last_ekin;
1721 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1723 if (bVV)
1725 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1727 else
1729 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1731 /* ######### END PREPARING EDR OUTPUT ########### */
1733 /* Output stuff */
1734 if (MASTER(cr))
1736 gmx_bool do_dr, do_or;
1738 if (fplog && do_log && bDoExpanded)
1740 /* only needed if doing expanded ensemble */
1741 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1742 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1744 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1746 if (bCalcEner)
1748 upd_mdebin(mdebin, bDoDHDL, TRUE,
1749 t, mdatoms->tmass, enerd, state,
1750 ir->fepvals, ir->expandedvals, lastbox,
1751 shake_vir, force_vir, total_vir, pres,
1752 ekind, mu_tot, constr);
1754 else
1756 upd_mdebin_step(mdebin);
1759 do_dr = do_per_step(step, ir->nstdisreout);
1760 do_or = do_per_step(step, ir->nstorireout);
1762 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1763 step, t,
1764 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1766 if (ir->ePull != epullNO)
1768 pull_print_output(ir->pull, step, t);
1771 if (do_per_step(step, ir->nstlog))
1773 if (fflush(fplog) != 0)
1775 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1779 if (bDoExpanded)
1781 /* Have to do this part _after_ outputting the logfile and the edr file */
1782 /* Gets written into the state at the beginning of next loop*/
1783 state->fep_state = lamnew;
1785 /* Print the remaining wall clock time for the run */
1786 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal()) && !bPMETuneRunning)
1788 if (shellfc)
1790 fprintf(stderr, "\n");
1792 print_time(stderr, walltime_accounting, step, ir, cr);
1795 /* Replica exchange */
1796 bExchanged = FALSE;
1797 if ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
1798 do_per_step(step, repl_ex_nst))
1800 bExchanged = replica_exchange(fplog, cr, repl_ex,
1801 state_global, enerd,
1802 state, step, t);
1804 if (bExchanged && DOMAINDECOMP(cr))
1806 dd_partition_system(fplog, step, cr, TRUE, 1,
1807 state_global, top_global, ir,
1808 state, &f, mdatoms, top, fr,
1809 vsite, shellfc, constr,
1810 nrnb, wcycle, FALSE);
1814 bFirstStep = FALSE;
1815 bInitStep = FALSE;
1816 bStartingFromCpt = FALSE;
1818 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1819 /* With all integrators, except VV, we need to retain the pressure
1820 * at the current step for coupling at the next step.
1822 if ((state->flags & (1<<estPRES_PREV)) &&
1823 (bGStatEveryStep ||
1824 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1826 /* Store the pressure in t_state for pressure coupling
1827 * at the next MD step.
1829 copy_mat(pres, state->pres_prev);
1832 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1834 if ( (membed != NULL) && (!bLastStep) )
1836 rescale_membed(step_rel, membed, state_global->x);
1839 if (bRerunMD)
1841 if (MASTER(cr))
1843 /* read next frame from input trajectory */
1844 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1847 if (PAR(cr))
1849 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1853 if (!bRerunMD || !rerun_fr.bStep)
1855 /* increase the MD step number */
1856 step++;
1857 step_rel++;
1860 cycles = wallcycle_stop(wcycle, ewcSTEP);
1861 if (DOMAINDECOMP(cr) && wcycle)
1863 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1866 if (bPMETuneRunning || bPMETuneTry)
1868 /* PME grid + cut-off optimization with GPUs or PME nodes */
1870 /* Count the total cycles over the last steps */
1871 cycles_pmes += cycles;
1873 /* We can only switch cut-off at NS steps */
1874 if (step % ir->nstlist == 0)
1876 /* PME grid + cut-off optimization with GPUs or PME nodes */
1877 if (bPMETuneTry)
1879 if (DDMASTER(cr->dd))
1881 /* PME node load is too high, start tuning */
1882 bPMETuneRunning = (dd_pme_f_ratio(cr->dd) >= 1.05);
1884 dd_bcast(cr->dd, sizeof(gmx_bool), &bPMETuneRunning);
1886 if (bPMETuneRunning || step_rel > ir->nstlist*50)
1888 bPMETuneTry = FALSE;
1891 if (bPMETuneRunning)
1893 /* init_step might not be a multiple of nstlist,
1894 * but the first cycle is always skipped anyhow.
1896 bPMETuneRunning =
1897 pme_load_balance(pme_loadbal, cr,
1898 (bVerbose && MASTER(cr)) ? stderr : NULL,
1899 fplog,
1900 ir, state, cycles_pmes,
1901 fr->ic, fr->nbv, &fr->pmedata,
1902 step);
1904 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
1905 fr->ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1906 fr->rlist = fr->ic->rlist;
1907 fr->rlistlong = fr->ic->rlistlong;
1908 fr->rcoulomb = fr->ic->rcoulomb;
1909 fr->rvdw = fr->ic->rvdw;
1911 cycles_pmes = 0;
1915 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1916 gs.set[eglsRESETCOUNTERS] != 0)
1918 /* Reset all the counters related to performance over the run */
1919 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1920 fr->nbv != NULL && fr->nbv->bUseGPU ? fr->nbv->cu_nbv : NULL);
1921 wcycle_set_reset_counters(wcycle, -1);
1922 if (!(cr->duty & DUTY_PME))
1924 /* Tell our PME node to reset its counters */
1925 gmx_pme_send_resetcounters(cr, step);
1927 /* Correct max_hours for the elapsed time */
1928 max_hours -= elapsed_time/(60.0*60.0);
1929 bResetCountersHalfMaxH = FALSE;
1930 gs.set[eglsRESETCOUNTERS] = 0;
1934 /* End of main MD loop */
1935 debug_gmx();
1937 /* Stop measuring walltime */
1938 walltime_accounting_end(walltime_accounting);
1940 if (bRerunMD && MASTER(cr))
1942 close_trj(status);
1945 if (!(cr->duty & DUTY_PME))
1947 /* Tell the PME only node to finish */
1948 gmx_pme_send_finish(cr);
1951 if (MASTER(cr))
1953 if (ir->nstcalcenergy > 0 && !bRerunMD)
1955 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1956 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1960 done_mdoutf(outf);
1961 debug_gmx();
1963 if (ir->nstlist == -1 && nlh.nns > 0 && fplog)
1965 fprintf(fplog, "Average neighborlist lifetime: %.1f steps, std.dev.: %.1f steps\n", nlh.s1/nlh.nns, sqrt(nlh.s2/nlh.nns - sqr(nlh.s1/nlh.nns)));
1966 fprintf(fplog, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh.ab/nlh.nns);
1969 if (pme_loadbal != NULL)
1971 pme_loadbal_done(pme_loadbal, cr, fplog,
1972 fr->nbv != NULL && fr->nbv->bUseGPU);
1975 if (shellfc && fplog)
1977 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1978 (nconverged*100.0)/step_rel);
1979 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1980 tcount/step_rel);
1983 if (repl_ex_nst > 0 && MASTER(cr))
1985 print_replica_exchange_statistics(fplog, repl_ex);
1988 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1990 return 0;