Merged legacyheaders/types/inputrec.h into mdtypes/inputrec.h
[gromacs.git] / src / programs / mdrun / md.cpp
blobf37310c306d847170f715b9c3b53f9dd447c17e5
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/md_logging.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/sighandler.h"
59 #include "gromacs/imd/imd.h"
60 #include "gromacs/legacyheaders/nrnb.h"
61 #include "gromacs/legacyheaders/types/commrec.h"
62 #include "gromacs/legacyheaders/types/enums.h"
63 #include "gromacs/legacyheaders/types/fcdata.h"
64 #include "gromacs/legacyheaders/types/force_flags.h"
65 #include "gromacs/legacyheaders/types/forcerec.h"
66 #include "gromacs/legacyheaders/types/group.h"
67 #include "gromacs/legacyheaders/types/interaction_const.h"
68 #include "gromacs/legacyheaders/types/mdatom.h"
69 #include "gromacs/legacyheaders/types/nrnb.h"
70 #include "gromacs/listed-forces/manage-threading.h"
71 #include "gromacs/math/utilities.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vectypes.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/md_support.h"
80 #include "gromacs/mdlib/mdatoms.h"
81 #include "gromacs/mdlib/mdebin.h"
82 #include "gromacs/mdlib/mdoutf.h"
83 #include "gromacs/mdlib/mdrun.h"
84 #include "gromacs/mdlib/mdrun_signalling.h"
85 #include "gromacs/mdlib/nb_verlet.h"
86 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
87 #include "gromacs/mdlib/ns.h"
88 #include "gromacs/mdlib/shellfc.h"
89 #include "gromacs/mdlib/sim_util.h"
90 #include "gromacs/mdlib/tgroup.h"
91 #include "gromacs/mdlib/trajectory_writing.h"
92 #include "gromacs/mdlib/update.h"
93 #include "gromacs/mdlib/vcm.h"
94 #include "gromacs/mdlib/vsite.h"
95 #include "gromacs/mdtypes/df_history.h"
96 #include "gromacs/mdtypes/energyhistory.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/state.h"
99 #include "gromacs/pbcutil/mshift.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/swap/swapcoords.h"
103 #include "gromacs/timing/wallcycle.h"
104 #include "gromacs/timing/walltime_accounting.h"
105 #include "gromacs/topology/atoms.h"
106 #include "gromacs/topology/idef.h"
107 #include "gromacs/topology/mtop_util.h"
108 #include "gromacs/topology/topology.h"
109 #include "gromacs/utility/basedefinitions.h"
110 #include "gromacs/utility/cstringutil.h"
111 #include "gromacs/utility/fatalerror.h"
112 #include "gromacs/utility/real.h"
113 #include "gromacs/utility/smalloc.h"
115 #include "deform.h"
116 #include "membed.h"
117 #include "repl_ex.h"
119 #ifdef GMX_FAHCORE
120 #include "corewrap.h"
121 #endif
123 static void reset_all_counters(FILE *fplog, t_commrec *cr,
124 gmx_int64_t step,
125 gmx_int64_t *step_rel, t_inputrec *ir,
126 gmx_wallcycle_t wcycle, t_nrnb *nrnb,
127 gmx_walltime_accounting_t walltime_accounting,
128 struct nonbonded_verlet_t *nbv)
130 char sbuf[STEPSTRSIZE];
132 /* Reset all the counters related to performance over the run */
133 md_print_warn(cr, fplog, "step %s: resetting all time and cycle counters\n",
134 gmx_step_str(step, sbuf));
136 if (use_GPU(nbv))
138 nbnxn_gpu_reset_timings(nbv);
141 wallcycle_stop(wcycle, ewcRUN);
142 wallcycle_reset_all(wcycle);
143 if (DOMAINDECOMP(cr))
145 reset_dd_statistics_counters(cr->dd);
147 init_nrnb(nrnb);
148 ir->init_step += *step_rel;
149 ir->nsteps -= *step_rel;
150 *step_rel = 0;
151 wallcycle_start(wcycle, ewcRUN);
152 walltime_accounting_start(walltime_accounting);
153 print_date_and_time(fplog, cr->nodeid, "Restarted time", gmx_gettime());
156 /*! \libinternal
157 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
158 int nfile, const t_filenm fnm[],
159 const gmx_output_env_t *oenv, gmx_bool bVerbose,
160 gmx_bool bCompact, int nstglobalcomm,
161 gmx_vsite_t *vsite, gmx_constr_t constr,
162 int stepout,
163 t_inputrec *inputrec,
164 gmx_mtop_t *top_global, t_fcdata *fcd,
165 t_state *state_global,
166 t_mdatoms *mdatoms,
167 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
168 gmx_edsam_t ed,
169 t_forcerec *fr,
170 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
171 gmx_membed_t *membed,
172 real cpt_period, real max_hours,
173 int imdport,
174 unsigned long Flags,
175 gmx_walltime_accounting_t walltime_accounting)
177 double gmx::do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
178 const gmx_output_env_t *oenv, gmx_bool bVerbose, gmx_bool bCompact,
179 int nstglobalcomm,
180 gmx_vsite_t *vsite, gmx_constr_t constr,
181 int stepout, t_inputrec *ir,
182 gmx_mtop_t *top_global,
183 t_fcdata *fcd,
184 t_state *state_global,
185 t_mdatoms *mdatoms,
186 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
187 gmx_edsam_t ed, t_forcerec *fr,
188 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed, gmx_membed_t *membed,
189 real cpt_period, real max_hours,
190 int imdport,
191 unsigned long Flags,
192 gmx_walltime_accounting_t walltime_accounting)
194 gmx_mdoutf_t outf = NULL;
195 gmx_int64_t step, step_rel;
196 double elapsed_time;
197 double t, t0, lam0[efptNR];
198 gmx_bool bGStatEveryStep, bGStat, bCalcVir, bCalcEner;
199 gmx_bool bNS, bNStList, bSimAnn, bStopCM, bRerunMD, bNotLastFrame = FALSE,
200 bFirstStep, bStateFromCP, bStateFromTPX, bInitStep, bLastStep,
201 bBornRadii, bStartingFromCpt;
202 gmx_bool bDoDHDL = FALSE, bDoFEP = FALSE, bDoExpanded = FALSE;
203 gmx_bool do_ene, do_log, do_verbose, bRerunWarnNoV = TRUE,
204 bForceUpdate = FALSE, bCPT;
205 gmx_bool bMasterState;
206 int force_flags, cglo_flags;
207 tensor force_vir, shake_vir, total_vir, tmp_vir, pres;
208 int i, m;
209 t_trxstatus *status;
210 rvec mu_tot;
211 t_vcm *vcm;
212 matrix pcoupl_mu, M;
213 t_trxframe rerun_fr;
214 gmx_repl_ex_t repl_ex = NULL;
215 int nchkpt = 1;
216 gmx_localtop_t *top;
217 t_mdebin *mdebin = NULL;
218 t_state *state = NULL;
219 rvec *f_global = NULL;
220 gmx_enerdata_t *enerd;
221 rvec *f = NULL;
222 gmx_global_stat_t gstat;
223 gmx_update_t upd = NULL;
224 t_graph *graph = NULL;
225 gmx_signalling_t gs;
226 gmx_groups_t *groups;
227 gmx_ekindata_t *ekind;
228 gmx_shellfc_t *shellfc;
229 int count, nconverged = 0;
230 double tcount = 0;
231 gmx_bool bConverged = TRUE, bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition;
232 gmx_bool bResetCountersHalfMaxH = FALSE;
233 gmx_bool bVV, bTemp, bPres, bTrotter;
234 gmx_bool bUpdateDoLR;
235 real dvdl_constr;
236 rvec *cbuf = NULL;
237 int cbuf_nalloc = 0;
238 matrix lastbox;
239 int lamnew = 0;
240 /* for FEP */
241 int nstfep = 0;
242 double cycles;
243 real saved_conserved_quantity = 0;
244 real last_ekin = 0;
245 t_extmass MassQ;
246 int **trotter_seq;
247 char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
248 int handled_stop_condition = gmx_stop_cond_none; /* compare to get_stop_condition*/
249 gmx_int64_t multisim_nsteps = -1; /* number of steps to do before first multisim
250 simulation stops. If equal to zero, don't
251 communicate any more between multisims.*/
252 /* PME load balancing data for GPU kernels */
253 pme_load_balancing_t *pme_loadbal = NULL;
254 gmx_bool bPMETune = FALSE;
255 gmx_bool bPMETunePrinting = FALSE;
257 /* Interactive MD */
258 gmx_bool bIMDstep = FALSE;
260 #ifdef GMX_FAHCORE
261 /* Temporary addition for FAHCORE checkpointing */
262 int chkpt_ret;
263 #endif
265 /* Check for special mdrun options */
266 bRerunMD = (Flags & MD_RERUN);
267 if (Flags & MD_RESETCOUNTERSHALFWAY)
269 if (ir->nsteps > 0)
271 /* Signal to reset the counters half the simulation steps. */
272 wcycle_set_reset_counters(wcycle, ir->nsteps/2);
274 /* Signal to reset the counters halfway the simulation time. */
275 bResetCountersHalfMaxH = (max_hours > 0);
278 /* md-vv uses averaged full step velocities for T-control
279 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
280 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
281 bVV = EI_VV(ir->eI);
282 bTrotter = (bVV && (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)));
284 if (bRerunMD)
286 /* Since we don't know if the frames read are related in any way,
287 * rebuild the neighborlist at every step.
289 ir->nstlist = 1;
290 ir->nstcalcenergy = 1;
291 nstglobalcomm = 1;
294 check_ir_old_tpx_versions(cr, fplog, ir, top_global);
296 nstglobalcomm = check_nstglobalcomm(fplog, cr, nstglobalcomm, ir);
297 bGStatEveryStep = (nstglobalcomm == 1);
299 if (bRerunMD)
301 ir->nstxout_compressed = 0;
303 groups = &top_global->groups;
305 /* Initial values */
306 init_md(fplog, cr, ir, oenv, &t, &t0, state_global->lambda,
307 &(state_global->fep_state), lam0,
308 nrnb, top_global, &upd,
309 nfile, fnm, &outf, &mdebin,
310 force_vir, shake_vir, mu_tot, &bSimAnn, &vcm, Flags, wcycle);
312 clear_mat(total_vir);
313 clear_mat(pres);
314 /* Energy terms and groups */
315 snew(enerd, 1);
316 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
317 enerd);
318 if (DOMAINDECOMP(cr))
320 f = NULL;
322 else
324 snew(f, top_global->natoms);
327 /* Kinetic energy data */
328 snew(ekind, 1);
329 init_ekindata(fplog, top_global, &(ir->opts), ekind);
330 /* Copy the cos acceleration to the groups struct */
331 ekind->cosacc.cos_accel = ir->cos_accel;
333 gstat = global_stat_init(ir);
335 /* Check for polarizable models and flexible constraints */
336 shellfc = init_shell_flexcon(fplog,
337 top_global, n_flexible_constraints(constr),
338 (ir->bContinuation ||
339 (DOMAINDECOMP(cr) && !MASTER(cr))) ?
340 NULL : state_global->x);
341 if (shellfc && ir->nstcalcenergy != 1)
343 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);
345 if (shellfc && DOMAINDECOMP(cr))
347 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
349 if (shellfc && ir->eI == eiNM)
351 /* Currently shells don't work with Normal Modes */
352 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");
355 if (vsite && ir->eI == eiNM)
357 /* Currently virtual sites don't work with Normal Modes */
358 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");
361 if (DEFORM(*ir))
363 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
364 set_deform_reference_box(upd,
365 deform_init_init_step_tpx,
366 deform_init_box_tpx);
367 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
371 double io = compute_io(ir, top_global->natoms, groups, mdebin->ebin->nener, 1);
372 if ((io > 2000) && MASTER(cr))
374 fprintf(stderr,
375 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
376 io);
380 if (DOMAINDECOMP(cr))
382 top = dd_init_local_top(top_global);
384 snew(state, 1);
385 dd_init_local_state(cr->dd, state_global, state);
387 if (DDMASTER(cr->dd) && ir->nstfout)
389 snew(f_global, state_global->natoms);
392 else
394 top = gmx_mtop_generate_local_top(top_global, ir);
396 forcerec_set_excl_load(fr, top);
398 state = serial_init_local_state(state_global);
399 f_global = f;
401 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
403 if (vsite)
405 set_vsite_top(vsite, top, mdatoms, cr);
408 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
410 graph = mk_graph(fplog, &(top->idef), 0, top_global->natoms, FALSE, FALSE);
413 if (shellfc)
415 make_local_shells(cr, mdatoms, shellfc);
418 setup_bonded_threading(fr, &top->idef);
421 /* Set up interactive MD (IMD) */
422 init_IMD(ir, cr, top_global, fplog, ir->nstcalcenergy, state_global->x,
423 nfile, fnm, oenv, imdport, Flags);
425 if (DOMAINDECOMP(cr))
427 /* Distribute the charge groups over the nodes from the master node */
428 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
429 state_global, top_global, ir,
430 state, &f, mdatoms, top, fr,
431 vsite, shellfc, constr,
432 nrnb, NULL, FALSE);
436 update_mdatoms(mdatoms, state->lambda[efptMASS]);
438 if (opt2bSet("-cpi", nfile, fnm))
440 bStateFromCP = gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr);
442 else
444 bStateFromCP = FALSE;
447 if (ir->bExpanded)
449 init_expanded_ensemble(bStateFromCP, ir, &state->dfhist);
452 if (MASTER(cr))
454 if (bStateFromCP)
456 /* Update mdebin with energy history if appending to output files */
457 if (Flags & MD_APPENDFILES)
459 restore_energyhistory_from_state(mdebin, state_global->enerhist);
461 else
463 /* We might have read an energy history from checkpoint,
464 * free the allocated memory and reset the counts.
466 done_energyhistory(state_global->enerhist);
467 init_energyhistory(state_global->enerhist);
470 /* Set the initial energy history in state by updating once */
471 update_energyhistory(state_global->enerhist, mdebin);
474 /* Initialize constraints */
475 if (constr && !DOMAINDECOMP(cr))
477 set_constraints(constr, top, ir, mdatoms, cr);
480 if (repl_ex_nst > 0 && MASTER(cr))
482 repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir,
483 repl_ex_nst, repl_ex_nex, repl_ex_seed);
486 /* PME tuning is only supported with PME for Coulomb. Is is not supported
487 * with only LJ PME, or for reruns.
489 bPMETune = ((Flags & MD_TUNEPME) && EEL_PME(fr->eeltype) && !bRerunMD &&
490 !(Flags & MD_REPRODUCIBLE));
491 if (bPMETune)
493 pme_loadbal_init(&pme_loadbal, cr, fplog, ir, state->box,
494 fr->ic, fr->pmedata, use_GPU(fr->nbv),
495 &bPMETunePrinting);
498 if (!ir->bContinuation && !bRerunMD)
500 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
502 /* Set the velocities of frozen particles to zero */
503 for (i = 0; i < mdatoms->homenr; i++)
505 for (m = 0; m < DIM; m++)
507 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
509 state->v[i][m] = 0;
515 if (constr)
517 /* Constrain the initial coordinates and velocities */
518 do_constrain_first(fplog, constr, ir, mdatoms, state,
519 cr, nrnb, fr, top);
521 if (vsite)
523 /* Construct the virtual sites for the initial configuration */
524 construct_vsites(vsite, state->x, ir->delta_t, NULL,
525 top->idef.iparams, top->idef.il,
526 fr->ePBC, fr->bMolPBC, cr, state->box);
530 if (IR_TWINRANGE(*ir) && repl_ex_nst % ir->nstcalclr != 0)
532 /* We should exchange at nstcalclr steps to get correct integration */
533 gmx_fatal(FARGS, "The replica exchange period (%d) is not divisible by nstcalclr (%d)", repl_ex_nst, ir->nstcalclr);
536 if (ir->efep != efepNO)
538 /* Set free energy calculation frequency as the greatest common
539 * denominator of nstdhdl and repl_ex_nst.
540 * Check for nstcalclr with twin-range, since we need the long-range
541 * contribution to the free-energy at the correct (nstcalclr) steps.
543 nstfep = ir->fepvals->nstdhdl;
544 if (ir->bExpanded)
546 if (IR_TWINRANGE(*ir) &&
547 ir->expandedvals->nstexpanded % ir->nstcalclr != 0)
549 gmx_fatal(FARGS, "nstexpanded should be divisible by nstcalclr");
551 nstfep = gmx_greatest_common_divisor(ir->expandedvals->nstexpanded, nstfep);
553 if (repl_ex_nst > 0)
555 nstfep = gmx_greatest_common_divisor(repl_ex_nst, nstfep);
557 /* We checked divisibility of repl_ex_nst and nstcalclr above */
558 if (IR_TWINRANGE(*ir) && nstfep % ir->nstcalclr != 0)
560 gmx_incons("nstfep not divisible by nstcalclr");
564 /* Be REALLY careful about what flags you set here. You CANNOT assume
565 * this is the first step, since we might be restarting from a checkpoint,
566 * and in that case we should not do any modifications to the state.
568 bStopCM = (ir->comm_mode != ecmNO && !ir->bContinuation);
570 cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
571 | (bStopCM ? CGLO_STOPCM : 0)
572 | (bVV ? CGLO_PRESSURE : 0)
573 | (bVV ? CGLO_CONSTRAINT : 0)
574 | ((Flags & MD_READ_EKIN) ? CGLO_READEKIN : 0));
576 bSumEkinhOld = FALSE;
577 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
578 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
579 constr, NULL, FALSE, state->box,
580 top_global, &bSumEkinhOld, cglo_flags);
581 if (ir->eI == eiVVAK)
583 /* a second call to get the half step temperature initialized as well */
584 /* we do the same call as above, but turn the pressure off -- internally to
585 compute_globals, this is recognized as a velocity verlet half-step
586 kinetic energy calculation. This minimized excess variables, but
587 perhaps loses some logic?*/
589 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
590 NULL, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
591 constr, NULL, FALSE, state->box,
592 top_global, &bSumEkinhOld,
593 cglo_flags &~(CGLO_STOPCM | CGLO_PRESSURE));
596 /* Calculate the initial half step temperature, and save the ekinh_old */
597 if (!(Flags & MD_STARTFROMCPT))
599 for (i = 0; (i < ir->opts.ngtc); i++)
601 copy_mat(ekind->tcstat[i].ekinh, ekind->tcstat[i].ekinh_old);
604 if (ir->eI != eiVV)
606 enerd->term[F_TEMP] *= 2; /* result of averages being done over previous and current step,
607 and there is no previous step */
610 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
611 temperature control */
612 trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter);
614 if (MASTER(cr))
616 if (constr && !ir->bContinuation && ir->eConstrAlg == econtLINCS)
618 fprintf(fplog,
619 "RMS relative constraint deviation after constraining: %.2e\n",
620 constr_rmsd(constr, FALSE));
622 if (EI_STATE_VELOCITY(ir->eI))
624 fprintf(fplog, "Initial temperature: %g K\n", enerd->term[F_TEMP]);
626 if (bRerunMD)
628 fprintf(stderr, "starting md rerun '%s', reading coordinates from"
629 " input trajectory '%s'\n\n",
630 *(top_global->name), opt2fn("-rerun", nfile, fnm));
631 if (bVerbose)
633 fprintf(stderr, "Calculated time to finish depends on nsteps from "
634 "run input file,\nwhich may not correspond to the time "
635 "needed to process input trajectory.\n\n");
638 else
640 char tbuf[20];
641 fprintf(stderr, "starting mdrun '%s'\n",
642 *(top_global->name));
643 if (ir->nsteps >= 0)
645 sprintf(tbuf, "%8.1f", (ir->init_step+ir->nsteps)*ir->delta_t);
647 else
649 sprintf(tbuf, "%s", "infinite");
651 if (ir->init_step > 0)
653 fprintf(stderr, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
654 gmx_step_str(ir->init_step+ir->nsteps, sbuf), tbuf,
655 gmx_step_str(ir->init_step, sbuf2),
656 ir->init_step*ir->delta_t);
658 else
660 fprintf(stderr, "%s steps, %s ps.\n",
661 gmx_step_str(ir->nsteps, sbuf), tbuf);
664 fprintf(fplog, "\n");
667 walltime_accounting_start(walltime_accounting);
668 wallcycle_start(wcycle, ewcRUN);
669 print_start(fplog, cr, walltime_accounting, "mdrun");
671 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
672 #ifdef GMX_FAHCORE
673 chkpt_ret = fcCheckPointParallel( cr->nodeid,
674 NULL, 0);
675 if (chkpt_ret == 0)
677 gmx_fatal( 3, __FILE__, __LINE__, "Checkpoint error on step %d\n", 0 );
679 #endif
681 /***********************************************************
683 * Loop over MD steps
685 ************************************************************/
687 /* if rerunMD then read coordinates and velocities from input trajectory */
688 if (bRerunMD)
690 if (getenv("GMX_FORCE_UPDATE"))
692 bForceUpdate = TRUE;
695 rerun_fr.natoms = 0;
696 if (MASTER(cr))
698 bNotLastFrame = read_first_frame(oenv, &status,
699 opt2fn("-rerun", nfile, fnm),
700 &rerun_fr, TRX_NEED_X | TRX_READ_V);
701 if (rerun_fr.natoms != top_global->natoms)
703 gmx_fatal(FARGS,
704 "Number of atoms in trajectory (%d) does not match the "
705 "run input file (%d)\n",
706 rerun_fr.natoms, top_global->natoms);
708 if (ir->ePBC != epbcNONE)
710 if (!rerun_fr.bBox)
712 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);
714 if (max_cutoff2(ir->ePBC, rerun_fr.box) < sqr(fr->rlistlong))
716 gmx_fatal(FARGS, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr.step, rerun_fr.time);
721 if (PAR(cr))
723 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
726 if (ir->ePBC != epbcNONE)
728 /* Set the shift vectors.
729 * Necessary here when have a static box different from the tpr box.
731 calc_shifts(rerun_fr.box, fr->shift_vec);
735 /* loop over MD steps or if rerunMD to end of input trajectory */
736 bFirstStep = TRUE;
737 /* Skip the first Nose-Hoover integration when we get the state from tpx */
738 bStateFromTPX = !bStateFromCP;
739 bInitStep = bFirstStep && (bStateFromTPX || bVV);
740 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
741 bSumEkinhOld = FALSE;
742 bExchanged = FALSE;
743 bNeedRepartition = FALSE;
745 init_global_signals(&gs, cr, ir, repl_ex_nst);
747 step = ir->init_step;
748 step_rel = 0;
750 if (MULTISIM(cr) && (repl_ex_nst <= 0 ))
752 /* check how many steps are left in other sims */
753 multisim_nsteps = get_multisim_nsteps(cr, ir->nsteps);
757 /* and stop now if we should */
758 bLastStep = (bRerunMD || (ir->nsteps >= 0 && step_rel > ir->nsteps) ||
759 ((multisim_nsteps >= 0) && (step_rel >= multisim_nsteps )));
760 while (!bLastStep || (bRerunMD && bNotLastFrame))
763 /* Determine if this is a neighbor search step */
764 bNStList = (ir->nstlist > 0 && step % ir->nstlist == 0);
766 if (bPMETune && bNStList)
768 /* PME grid + cut-off optimization with GPUs or PME nodes */
769 pme_loadbal_do(pme_loadbal, cr,
770 (bVerbose && MASTER(cr)) ? stderr : NULL,
771 fplog,
772 ir, fr, state,
773 wcycle,
774 step, step_rel,
775 &bPMETunePrinting);
778 wallcycle_start(wcycle, ewcSTEP);
780 if (bRerunMD)
782 if (rerun_fr.bStep)
784 step = rerun_fr.step;
785 step_rel = step - ir->init_step;
787 if (rerun_fr.bTime)
789 t = rerun_fr.time;
791 else
793 t = step;
796 else
798 bLastStep = (step_rel == ir->nsteps);
799 t = t0 + step*ir->delta_t;
802 if (ir->efep != efepNO || ir->bSimTemp)
804 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
805 requiring different logic. */
807 set_current_lambdas(step, ir->fepvals, bRerunMD, &rerun_fr, state_global, state, lam0);
808 bDoDHDL = do_per_step(step, ir->fepvals->nstdhdl);
809 bDoFEP = ((ir->efep != efepNO) && do_per_step(step, nstfep));
810 bDoExpanded = (do_per_step(step, ir->expandedvals->nstexpanded)
811 && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt));
814 bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep &&
815 do_per_step(step, repl_ex_nst));
817 if (bSimAnn)
819 update_annealing_target_temp(&(ir->opts), t);
822 if (bRerunMD)
824 if (!DOMAINDECOMP(cr) || MASTER(cr))
826 for (i = 0; i < state_global->natoms; i++)
828 copy_rvec(rerun_fr.x[i], state_global->x[i]);
830 if (rerun_fr.bV)
832 for (i = 0; i < state_global->natoms; i++)
834 copy_rvec(rerun_fr.v[i], state_global->v[i]);
837 else
839 for (i = 0; i < state_global->natoms; i++)
841 clear_rvec(state_global->v[i]);
843 if (bRerunWarnNoV)
845 fprintf(stderr, "\nWARNING: Some frames do not contain velocities.\n"
846 " Ekin, temperature and pressure are incorrect,\n"
847 " the virial will be incorrect when constraints are present.\n"
848 "\n");
849 bRerunWarnNoV = FALSE;
853 copy_mat(rerun_fr.box, state_global->box);
854 copy_mat(state_global->box, state->box);
856 if (vsite && (Flags & MD_RERUN_VSITE))
858 if (DOMAINDECOMP(cr))
860 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
862 if (graph)
864 /* Following is necessary because the graph may get out of sync
865 * with the coordinates if we only have every N'th coordinate set
867 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
868 shift_self(graph, state->box, state->x);
870 construct_vsites(vsite, state->x, ir->delta_t, state->v,
871 top->idef.iparams, top->idef.il,
872 fr->ePBC, fr->bMolPBC, cr, state->box);
873 if (graph)
875 unshift_self(graph, state->box, state->x);
880 /* Stop Center of Mass motion */
881 bStopCM = (ir->comm_mode != ecmNO && do_per_step(step, ir->nstcomm));
883 if (bRerunMD)
885 /* for rerun MD always do Neighbour Searching */
886 bNS = (bFirstStep || ir->nstlist != 0);
887 bNStList = bNS;
889 else
891 /* Determine whether or not to do Neighbour Searching and LR */
892 bNS = (bFirstStep || bNStList || bExchanged || bNeedRepartition);
895 /* check whether we should stop because another simulation has
896 stopped. */
897 if (MULTISIM(cr))
899 if ( (multisim_nsteps >= 0) && (step_rel >= multisim_nsteps) &&
900 (multisim_nsteps != ir->nsteps) )
902 if (bNS)
904 if (MASTER(cr))
906 fprintf(stderr,
907 "Stopping simulation %d because another one has finished\n",
908 cr->ms->sim);
910 bLastStep = TRUE;
911 gs.sig[eglsCHKPT] = 1;
916 /* < 0 means stop at next step, > 0 means stop at next NS step */
917 if ( (gs.set[eglsSTOPCOND] < 0) ||
918 ( (gs.set[eglsSTOPCOND] > 0) && (bNStList || ir->nstlist == 0) ) )
920 bLastStep = TRUE;
923 /* Determine whether or not to update the Born radii if doing GB */
924 bBornRadii = bFirstStep;
925 if (ir->implicit_solvent && (step % ir->nstgbradii == 0))
927 bBornRadii = TRUE;
930 do_log = do_per_step(step, ir->nstlog) || bFirstStep || bLastStep;
931 do_verbose = bVerbose &&
932 (step % stepout == 0 || bFirstStep || bLastStep);
934 if (bNS && !(bFirstStep && ir->bContinuation && !bRerunMD))
936 if (bRerunMD)
938 bMasterState = TRUE;
940 else
942 bMasterState = FALSE;
943 /* Correct the new box if it is too skewed */
944 if (DYNAMIC_BOX(*ir))
946 if (correct_box(fplog, step, state->box, graph))
948 bMasterState = TRUE;
951 if (DOMAINDECOMP(cr) && bMasterState)
953 dd_collect_state(cr->dd, state, state_global);
957 if (DOMAINDECOMP(cr))
959 /* Repartition the domain decomposition */
960 dd_partition_system(fplog, step, cr,
961 bMasterState, nstglobalcomm,
962 state_global, top_global, ir,
963 state, &f, mdatoms, top, fr,
964 vsite, shellfc, constr,
965 nrnb, wcycle,
966 do_verbose && !bPMETunePrinting);
970 if (MASTER(cr) && do_log)
972 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
975 if (ir->efep != efepNO)
977 update_mdatoms(mdatoms, state->lambda[efptMASS]);
980 if ((bRerunMD && rerun_fr.bV) || bExchanged)
983 /* We need the kinetic energy at minus the half step for determining
984 * the full step kinetic energy and possibly for T-coupling.*/
985 /* This may not be quite working correctly yet . . . . */
986 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
987 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
988 constr, NULL, FALSE, state->box,
989 top_global, &bSumEkinhOld,
990 CGLO_GSTAT | CGLO_TEMPERATURE);
992 clear_mat(force_vir);
994 /* We write a checkpoint at this MD step when:
995 * either at an NS step when we signalled through gs,
996 * or at the last step (but not when we do not want confout),
997 * but never at the first step or with rerun.
999 bCPT = (((gs.set[eglsCHKPT] && (bNS || ir->nstlist == 0)) ||
1000 (bLastStep && (Flags & MD_CONFOUT))) &&
1001 step > ir->init_step && !bRerunMD);
1002 if (bCPT)
1004 gs.set[eglsCHKPT] = 0;
1007 /* Determine the energy and pressure:
1008 * at nstcalcenergy steps and at energy output steps (set below).
1010 if (EI_VV(ir->eI) && (!bInitStep))
1012 /* for vv, the first half of the integration actually corresponds
1013 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1014 but the virial needs to be calculated on both the current step and the 'next' step. Future
1015 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1017 bCalcEner = do_per_step(step-1, ir->nstcalcenergy);
1018 bCalcVir = bCalcEner ||
1019 (ir->epc != epcNO && (do_per_step(step, ir->nstpcouple) || do_per_step(step-1, ir->nstpcouple)));
1021 else
1023 bCalcEner = do_per_step(step, ir->nstcalcenergy);
1024 bCalcVir = bCalcEner ||
1025 (ir->epc != epcNO && do_per_step(step, ir->nstpcouple));
1028 /* Do we need global communication ? */
1029 bGStat = (bCalcVir || bCalcEner || bStopCM ||
1030 do_per_step(step, nstglobalcomm) ||
1031 (bVV && IR_NVT_TROTTER(ir) && do_per_step(step-1, nstglobalcomm)));
1033 do_ene = (do_per_step(step, ir->nstenergy) || bLastStep);
1035 if (do_ene || do_log || bDoReplEx)
1037 bCalcVir = TRUE;
1038 bCalcEner = TRUE;
1039 bGStat = TRUE;
1042 force_flags = (GMX_FORCE_STATECHANGED |
1043 ((DYNAMIC_BOX(*ir) || bRerunMD) ? GMX_FORCE_DYNAMICBOX : 0) |
1044 GMX_FORCE_ALLFORCES |
1045 GMX_FORCE_SEPLRF |
1046 (bCalcVir ? GMX_FORCE_VIRIAL : 0) |
1047 (bCalcEner ? GMX_FORCE_ENERGY : 0) |
1048 (bDoFEP ? GMX_FORCE_DHDL : 0)
1051 if (fr->bTwinRange)
1053 if (do_per_step(step, ir->nstcalclr))
1055 force_flags |= GMX_FORCE_DO_LR;
1059 if (shellfc)
1061 /* Now is the time to relax the shells */
1062 count = relax_shell_flexcon(fplog, cr, bVerbose, step,
1063 ir, bNS, force_flags,
1064 top,
1065 constr, enerd, fcd,
1066 state, f, force_vir, mdatoms,
1067 nrnb, wcycle, graph, groups,
1068 shellfc, fr, bBornRadii, t, mu_tot,
1069 &bConverged, vsite,
1070 mdoutf_get_fp_field(outf));
1071 tcount += count;
1073 if (bConverged)
1075 nconverged++;
1078 else
1080 /* The coordinates (x) are shifted (to get whole molecules)
1081 * in do_force.
1082 * This is parallellized as well, and does communication too.
1083 * Check comments in sim_util.c
1085 do_force(fplog, cr, ir, step, nrnb, wcycle, top, groups,
1086 state->box, state->x, &state->hist,
1087 f, force_vir, mdatoms, enerd, fcd,
1088 state->lambda, graph,
1089 fr, vsite, mu_tot, t, mdoutf_get_fp_field(outf), ed, bBornRadii,
1090 (bNS ? GMX_FORCE_NS : 0) | force_flags);
1093 if (bVV && !bStartingFromCpt && !bRerunMD)
1094 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1096 rvec *vbuf = NULL;
1098 wallcycle_start(wcycle, ewcUPDATE);
1099 if (ir->eI == eiVV && bInitStep)
1101 /* if using velocity verlet with full time step Ekin,
1102 * take the first half step only to compute the
1103 * virial for the first step. From there,
1104 * revert back to the initial coordinates
1105 * so that the input is actually the initial step.
1107 snew(vbuf, state->natoms);
1108 copy_rvecn(state->v, vbuf, 0, state->natoms); /* should make this better for parallelizing? */
1110 else
1112 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1113 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ1);
1116 /* If we are using twin-range interactions where the long-range component
1117 * is only evaluated every nstcalclr>1 steps, we should do a special update
1118 * step to combine the long-range forces on these steps.
1119 * For nstcalclr=1 this is not done, since the forces would have been added
1120 * directly to the short-range forces already.
1122 * TODO Remove various aspects of VV+twin-range in master
1123 * branch, because VV integrators did not ever support
1124 * twin-range multiple time stepping with constraints.
1126 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1128 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC,
1129 f, bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1130 ekind, M, upd, bInitStep, etrtVELOCITY1,
1131 cr, nrnb, constr, &top->idef);
1133 if (!bRerunMD || rerun_fr.bV || bForceUpdate) /* Why is rerun_fr.bV here? Unclear. */
1135 wallcycle_stop(wcycle, ewcUPDATE);
1136 update_constraints(fplog, step, NULL, ir, mdatoms,
1137 state, fr->bMolPBC, graph, f,
1138 &top->idef, shake_vir,
1139 cr, nrnb, wcycle, upd, constr,
1140 TRUE, bCalcVir);
1141 wallcycle_start(wcycle, ewcUPDATE);
1142 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1144 /* Correct the virial for multiple time stepping */
1145 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1148 else if (graph)
1150 /* Need to unshift here if a do_force has been
1151 called in the previous step */
1152 unshift_self(graph, state->box, state->x);
1154 /* if VV, compute the pressure and constraints */
1155 /* For VV2, we strictly only need this if using pressure
1156 * control, but we really would like to have accurate pressures
1157 * printed out.
1158 * Think about ways around this in the future?
1159 * For now, keep this choice in comments.
1161 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1162 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1163 bPres = TRUE;
1164 bTemp = ((ir->eI == eiVV && (!bInitStep)) || (ir->eI == eiVVAK));
1165 if (bCalcEner && ir->eI == eiVVAK)
1167 bSumEkinhOld = TRUE;
1169 /* for vv, the first half of the integration actually corresponds to the previous step.
1170 So we need information from the last step in the first half of the integration */
1171 if (bGStat || do_per_step(step-1, nstglobalcomm))
1173 wallcycle_stop(wcycle, ewcUPDATE);
1174 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1175 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1176 constr, NULL, FALSE, state->box,
1177 top_global, &bSumEkinhOld,
1178 (bGStat ? CGLO_GSTAT : 0)
1179 | CGLO_ENERGY
1180 | (bTemp ? CGLO_TEMPERATURE : 0)
1181 | (bPres ? CGLO_PRESSURE : 0)
1182 | (bPres ? CGLO_CONSTRAINT : 0)
1183 | (bStopCM ? CGLO_STOPCM : 0)
1184 | CGLO_SCALEEKIN
1186 /* explanation of above:
1187 a) We compute Ekin at the full time step
1188 if 1) we are using the AveVel Ekin, and it's not the
1189 initial step, or 2) if we are using AveEkin, but need the full
1190 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1191 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1192 EkinAveVel because it's needed for the pressure */
1193 wallcycle_start(wcycle, ewcUPDATE);
1195 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1196 if (!bInitStep)
1198 if (bTrotter)
1200 m_add(force_vir, shake_vir, total_vir); /* we need the un-dispersion corrected total vir here */
1201 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ2);
1203 else
1205 if (bExchanged)
1207 wallcycle_stop(wcycle, ewcUPDATE);
1208 /* We need the kinetic energy at minus the half step for determining
1209 * the full step kinetic energy and possibly for T-coupling.*/
1210 /* This may not be quite working correctly yet . . . . */
1211 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1212 wcycle, enerd, NULL, NULL, NULL, NULL, mu_tot,
1213 constr, NULL, FALSE, state->box,
1214 top_global, &bSumEkinhOld,
1215 CGLO_GSTAT | CGLO_TEMPERATURE);
1216 wallcycle_start(wcycle, ewcUPDATE);
1220 if (bTrotter && !bInitStep)
1222 copy_mat(shake_vir, state->svir_prev);
1223 copy_mat(force_vir, state->fvir_prev);
1224 if (IR_NVT_TROTTER(ir) && ir->eI == eiVV)
1226 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1227 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, NULL, (ir->eI == eiVV), FALSE);
1228 enerd->term[F_EKIN] = trace(ekind->ekin);
1231 /* if it's the initial step, we performed this first step just to get the constraint virial */
1232 if (ir->eI == eiVV && bInitStep)
1234 copy_rvecn(vbuf, state->v, 0, state->natoms);
1235 sfree(vbuf);
1237 wallcycle_stop(wcycle, ewcUPDATE);
1240 /* compute the conserved quantity */
1241 if (bVV)
1243 saved_conserved_quantity = compute_conserved_from_auxiliary(ir, state, &MassQ);
1244 if (ir->eI == eiVV)
1246 last_ekin = enerd->term[F_EKIN];
1248 if ((ir->eDispCorr != edispcEnerPres) && (ir->eDispCorr != edispcAllEnerPres))
1250 saved_conserved_quantity -= enerd->term[F_DISPCORR];
1252 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1253 if (ir->efep != efepNO && !bRerunMD)
1255 sum_dhdl(enerd, state->lambda, ir->fepvals);
1259 /* ######## END FIRST UPDATE STEP ############## */
1260 /* ######## If doing VV, we now have v(dt) ###### */
1261 if (bDoExpanded)
1263 /* perform extended ensemble sampling in lambda - we don't
1264 actually move to the new state before outputting
1265 statistics, but if performing simulated tempering, we
1266 do update the velocities and the tau_t. */
1268 lamnew = ExpandedEnsembleDynamics(fplog, ir, enerd, state, &MassQ, state->fep_state, &state->dfhist, step, state->v, mdatoms);
1269 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1270 copy_df_history(&state_global->dfhist, &state->dfhist);
1273 /* Now we have the energies and forces corresponding to the
1274 * coordinates at time t. We must output all of this before
1275 * the update.
1277 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
1278 ir, state, state_global, top_global, fr,
1279 outf, mdebin, ekind, f, f_global,
1280 &nchkpt,
1281 bCPT, bRerunMD, bLastStep, (Flags & MD_CONFOUT),
1282 bSumEkinhOld);
1283 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1284 bIMDstep = do_IMD(ir->bIMD, step, cr, bNS, state->box, state->x, ir, t, wcycle);
1286 /* kludge -- virial is lost with restart for NPT control. Must restart */
1287 if (bStartingFromCpt && bVV)
1289 copy_mat(state->svir_prev, shake_vir);
1290 copy_mat(state->fvir_prev, force_vir);
1293 elapsed_time = walltime_accounting_get_current_elapsed_time(walltime_accounting);
1295 /* Check whether everything is still allright */
1296 if (((int)gmx_get_stop_condition() > handled_stop_condition)
1297 #ifdef GMX_THREAD_MPI
1298 && MASTER(cr)
1299 #endif
1302 /* this is just make gs.sig compatible with the hack
1303 of sending signals around by MPI_Reduce with together with
1304 other floats */
1305 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns)
1307 gs.sig[eglsSTOPCOND] = 1;
1309 if (gmx_get_stop_condition() == gmx_stop_cond_next)
1311 gs.sig[eglsSTOPCOND] = -1;
1313 /* < 0 means stop at next step, > 0 means stop at next NS step */
1314 if (fplog)
1316 fprintf(fplog,
1317 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1318 gmx_get_signal_name(),
1319 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1320 fflush(fplog);
1322 fprintf(stderr,
1323 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1324 gmx_get_signal_name(),
1325 gs.sig[eglsSTOPCOND] == 1 ? "NS " : "");
1326 fflush(stderr);
1327 handled_stop_condition = (int)gmx_get_stop_condition();
1329 else if (MASTER(cr) && (bNS || ir->nstlist <= 0) &&
1330 (max_hours > 0 && elapsed_time > max_hours*60.0*60.0*0.99) &&
1331 gs.sig[eglsSTOPCOND] == 0 && gs.set[eglsSTOPCOND] == 0)
1333 /* Signal to terminate the run */
1334 gs.sig[eglsSTOPCOND] = 1;
1335 if (fplog)
1337 fprintf(fplog, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1339 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step, sbuf), max_hours*0.99);
1342 if (bResetCountersHalfMaxH && MASTER(cr) &&
1343 elapsed_time > max_hours*60.0*60.0*0.495)
1345 /* Set flag that will communicate the signal to all ranks in the simulation */
1346 gs.sig[eglsRESETCOUNTERS] = 1;
1349 /* In parallel we only have to check for checkpointing in steps
1350 * where we do global communication,
1351 * otherwise the other nodes don't know.
1353 if (MASTER(cr) && ((bGStat || !PAR(cr)) &&
1354 cpt_period >= 0 &&
1355 (cpt_period == 0 ||
1356 elapsed_time >= nchkpt*cpt_period*60.0)) &&
1357 gs.set[eglsCHKPT] == 0)
1359 gs.sig[eglsCHKPT] = 1;
1362 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1363 if (EI_VV(ir->eI))
1365 if (!bInitStep)
1367 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1369 if (ETC_ANDERSEN(ir->etc)) /* keep this outside of update_tcouple because of the extra info required to pass */
1371 gmx_bool bIfRandomize;
1372 bIfRandomize = update_randomize_velocities(ir, step, cr, mdatoms, state, upd, constr);
1373 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1374 if (constr && bIfRandomize)
1376 update_constraints(fplog, step, NULL, ir, mdatoms,
1377 state, fr->bMolPBC, graph, f,
1378 &top->idef, tmp_vir,
1379 cr, nrnb, wcycle, upd, constr,
1380 TRUE, bCalcVir);
1384 /* ######### START SECOND UPDATE STEP ################# */
1385 /* Box is changed in update() when we do pressure coupling,
1386 * but we should still use the old box for energy corrections and when
1387 * writing it to the energy file, so it matches the trajectory files for
1388 * the same timestep above. Make a copy in a separate array.
1390 copy_mat(state->box, lastbox);
1392 dvdl_constr = 0;
1394 if (!bRerunMD || rerun_fr.bV || bForceUpdate)
1396 wallcycle_start(wcycle, ewcUPDATE);
1397 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1398 if (bTrotter)
1400 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ3);
1401 /* We can only do Berendsen coupling after we have summed
1402 * the kinetic energy or virial. Since the happens
1403 * in global_state after update, we should only do it at
1404 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1407 else
1409 update_tcouple(step, ir, state, ekind, &MassQ, mdatoms);
1410 update_pcouple(fplog, step, ir, state, pcoupl_mu, M, bInitStep);
1413 if (bVV)
1415 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1417 /* velocity half-step update */
1418 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1419 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1420 ekind, M, upd, FALSE, etrtVELOCITY2,
1421 cr, nrnb, constr, &top->idef);
1424 /* Above, initialize just copies ekinh into ekin,
1425 * it doesn't copy position (for VV),
1426 * and entire integrator for MD.
1429 if (ir->eI == eiVVAK)
1431 /* We probably only need md->homenr, not state->natoms */
1432 if (state->natoms > cbuf_nalloc)
1434 cbuf_nalloc = state->natoms;
1435 srenew(cbuf, cbuf_nalloc);
1437 copy_rvecn(state->x, cbuf, 0, state->natoms);
1439 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1441 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1442 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1443 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1444 wallcycle_stop(wcycle, ewcUPDATE);
1446 update_constraints(fplog, step, &dvdl_constr, ir, mdatoms, state,
1447 fr->bMolPBC, graph, f,
1448 &top->idef, shake_vir,
1449 cr, nrnb, wcycle, upd, constr,
1450 FALSE, bCalcVir);
1452 if (bCalcVir && bUpdateDoLR && ir->nstcalclr > 1)
1454 /* Correct the virial for multiple time stepping */
1455 m_sub(shake_vir, fr->vir_twin_constr, shake_vir);
1458 if (ir->eI == eiVVAK)
1460 /* erase F_EKIN and F_TEMP here? */
1461 /* just compute the kinetic energy at the half step to perform a trotter step */
1462 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1463 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1464 constr, NULL, FALSE, lastbox,
1465 top_global, &bSumEkinhOld,
1466 (bGStat ? CGLO_GSTAT : 0) | CGLO_TEMPERATURE
1468 wallcycle_start(wcycle, ewcUPDATE);
1469 trotter_update(ir, step, ekind, enerd, state, total_vir, mdatoms, &MassQ, trotter_seq, ettTSEQ4);
1470 /* now we know the scaling, we can compute the positions again again */
1471 copy_rvecn(cbuf, state->x, 0, state->natoms);
1473 bUpdateDoLR = (fr->bTwinRange && do_per_step(step, ir->nstcalclr));
1475 update_coords(fplog, step, ir, mdatoms, state, fr->bMolPBC, f,
1476 bUpdateDoLR, fr->f_twin, bCalcVir ? &fr->vir_twin_constr : NULL, fcd,
1477 ekind, M, upd, bInitStep, etrtPOSITION, cr, nrnb, constr, &top->idef);
1478 wallcycle_stop(wcycle, ewcUPDATE);
1480 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1481 /* are the small terms in the shake_vir here due
1482 * to numerical errors, or are they important
1483 * physically? I'm thinking they are just errors, but not completely sure.
1484 * For now, will call without actually constraining, constr=NULL*/
1485 update_constraints(fplog, step, NULL, ir, mdatoms,
1486 state, fr->bMolPBC, graph, f,
1487 &top->idef, tmp_vir,
1488 cr, nrnb, wcycle, upd, NULL,
1489 FALSE, bCalcVir);
1491 if (bVV)
1493 /* this factor or 2 correction is necessary
1494 because half of the constraint force is removed
1495 in the vv step, so we have to double it. See
1496 the Redmine issue #1255. It is not yet clear
1497 if the factor of 2 is exact, or just a very
1498 good approximation, and this will be
1499 investigated. The next step is to see if this
1500 can be done adding a dhdl contribution from the
1501 rattle step, but this is somewhat more
1502 complicated with the current code. Will be
1503 investigated, hopefully for 4.6.3. However,
1504 this current solution is much better than
1505 having it completely wrong.
1507 enerd->term[F_DVDL_CONSTR] += 2*dvdl_constr;
1509 else
1511 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1514 else if (graph)
1516 /* Need to unshift here */
1517 unshift_self(graph, state->box, state->x);
1520 if (vsite != NULL)
1522 wallcycle_start(wcycle, ewcVSITECONSTR);
1523 if (graph != NULL)
1525 shift_self(graph, state->box, state->x);
1527 construct_vsites(vsite, state->x, ir->delta_t, state->v,
1528 top->idef.iparams, top->idef.il,
1529 fr->ePBC, fr->bMolPBC, cr, state->box);
1531 if (graph != NULL)
1533 unshift_self(graph, state->box, state->x);
1535 wallcycle_stop(wcycle, ewcVSITECONSTR);
1538 /* ############## IF NOT VV, Calculate globals HERE ############ */
1539 /* With Leap-Frog we can skip compute_globals at
1540 * non-communication steps, but we need to calculate
1541 * the kinetic energy one step before communication.
1543 if (bGStat || (!EI_VV(ir->eI) && do_per_step(step+1, nstglobalcomm)))
1545 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, state_global, mdatoms, nrnb, vcm,
1546 wcycle, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
1547 constr, &gs,
1548 (step_rel % gs.nstms == 0) &&
1549 (multisim_nsteps < 0 || (step_rel < multisim_nsteps)),
1550 lastbox,
1551 top_global, &bSumEkinhOld,
1552 (bGStat ? CGLO_GSTAT : 0)
1553 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_ENERGY : 0)
1554 | (!EI_VV(ir->eI) && bStopCM ? CGLO_STOPCM : 0)
1555 | (!EI_VV(ir->eI) ? CGLO_TEMPERATURE : 0)
1556 | (!EI_VV(ir->eI) || bRerunMD ? CGLO_PRESSURE : 0)
1557 | CGLO_CONSTRAINT
1561 /* ############# END CALC EKIN AND PRESSURE ################# */
1563 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1564 the virial that should probably be addressed eventually. state->veta has better properies,
1565 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1566 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1568 if (ir->efep != efepNO && (!bVV || bRerunMD))
1570 /* Sum up the foreign energy and dhdl terms for md and sd.
1571 Currently done every step so that dhdl is correct in the .edr */
1572 sum_dhdl(enerd, state->lambda, ir->fepvals);
1574 update_box(fplog, step, ir, mdatoms, state, f,
1575 pcoupl_mu, nrnb, upd);
1577 /* ################# END UPDATE STEP 2 ################# */
1578 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1580 /* The coordinates (x) were unshifted in update */
1581 if (!bGStat)
1583 /* We will not sum ekinh_old,
1584 * so signal that we still have to do it.
1586 bSumEkinhOld = TRUE;
1589 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1591 /* use the directly determined last velocity, not actually the averaged half steps */
1592 if (bTrotter && ir->eI == eiVV)
1594 enerd->term[F_EKIN] = last_ekin;
1596 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1598 if (bVV)
1600 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + saved_conserved_quantity;
1602 else
1604 enerd->term[F_ECONSERVED] = enerd->term[F_ETOT] + compute_conserved_from_auxiliary(ir, state, &MassQ);
1606 /* ######### END PREPARING EDR OUTPUT ########### */
1608 /* Output stuff */
1609 if (MASTER(cr))
1611 gmx_bool do_dr, do_or;
1613 if (fplog && do_log && bDoExpanded)
1615 /* only needed if doing expanded ensemble */
1616 PrintFreeEnergyInfoToFile(fplog, ir->fepvals, ir->expandedvals, ir->bSimTemp ? ir->simtempvals : NULL,
1617 &state_global->dfhist, state->fep_state, ir->nstlog, step);
1619 if (!(bStartingFromCpt && (EI_VV(ir->eI))))
1621 if (bCalcEner)
1623 upd_mdebin(mdebin, bDoDHDL, TRUE,
1624 t, mdatoms->tmass, enerd, state,
1625 ir->fepvals, ir->expandedvals, lastbox,
1626 shake_vir, force_vir, total_vir, pres,
1627 ekind, mu_tot, constr);
1629 else
1631 upd_mdebin_step(mdebin);
1634 do_dr = do_per_step(step, ir->nstdisreout);
1635 do_or = do_per_step(step, ir->nstorireout);
1637 print_ebin(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or, do_log ? fplog : NULL,
1638 step, t,
1639 eprNORMAL, bCompact, mdebin, fcd, groups, &(ir->opts));
1641 if (ir->bPull)
1643 pull_print_output(ir->pull_work, step, t);
1646 if (do_per_step(step, ir->nstlog))
1648 if (fflush(fplog) != 0)
1650 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
1654 if (bDoExpanded)
1656 /* Have to do this part _after_ outputting the logfile and the edr file */
1657 /* Gets written into the state at the beginning of next loop*/
1658 state->fep_state = lamnew;
1660 /* Print the remaining wall clock time for the run */
1661 if (MULTIMASTER(cr) &&
1662 (do_verbose || gmx_got_usr_signal()) &&
1663 !bPMETunePrinting)
1665 if (shellfc)
1667 fprintf(stderr, "\n");
1669 print_time(stderr, walltime_accounting, step, ir, cr);
1672 /* Ion/water position swapping.
1673 * Not done in last step since trajectory writing happens before this call
1674 * in the MD loop and exchanges would be lost anyway. */
1675 bNeedRepartition = FALSE;
1676 if ((ir->eSwapCoords != eswapNO) && (step > 0) && !bLastStep &&
1677 do_per_step(step, ir->swap->nstswap))
1679 bNeedRepartition = do_swapcoords(cr, step, t, ir, wcycle,
1680 bRerunMD ? rerun_fr.x : state->x,
1681 bRerunMD ? rerun_fr.box : state->box,
1682 top_global, MASTER(cr) && bVerbose, bRerunMD);
1684 if (bNeedRepartition && DOMAINDECOMP(cr))
1686 dd_collect_state(cr->dd, state, state_global);
1690 /* Replica exchange */
1691 bExchanged = FALSE;
1692 if (bDoReplEx)
1694 bExchanged = replica_exchange(fplog, cr, repl_ex,
1695 state_global, enerd,
1696 state, step, t);
1699 if ( (bExchanged || bNeedRepartition) && DOMAINDECOMP(cr) )
1701 dd_partition_system(fplog, step, cr, TRUE, 1,
1702 state_global, top_global, ir,
1703 state, &f, mdatoms, top, fr,
1704 vsite, shellfc, constr,
1705 nrnb, wcycle, FALSE);
1708 bFirstStep = FALSE;
1709 bInitStep = FALSE;
1710 bStartingFromCpt = FALSE;
1712 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1713 /* With all integrators, except VV, we need to retain the pressure
1714 * at the current step for coupling at the next step.
1716 if ((state->flags & (1<<estPRES_PREV)) &&
1717 (bGStatEveryStep ||
1718 (ir->nstpcouple > 0 && step % ir->nstpcouple == 0)))
1720 /* Store the pressure in t_state for pressure coupling
1721 * at the next MD step.
1723 copy_mat(pres, state->pres_prev);
1726 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1728 if ( (membed != NULL) && (!bLastStep) )
1730 rescale_membed(step_rel, membed, state_global->x);
1733 if (bRerunMD)
1735 if (MASTER(cr))
1737 /* read next frame from input trajectory */
1738 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
1741 if (PAR(cr))
1743 rerun_parallel_comm(cr, &rerun_fr, &bNotLastFrame);
1747 cycles = wallcycle_stop(wcycle, ewcSTEP);
1748 if (DOMAINDECOMP(cr) && wcycle)
1750 dd_cycles_add(cr->dd, cycles, ddCyclStep);
1753 if (!bRerunMD || !rerun_fr.bStep)
1755 /* increase the MD step number */
1756 step++;
1757 step_rel++;
1760 /* TODO make a counter-reset module */
1761 /* If it is time to reset counters, set a flag that remains
1762 true until counters actually get reset */
1763 if (step_rel == wcycle_get_reset_counters(wcycle) ||
1764 gs.set[eglsRESETCOUNTERS] != 0)
1766 if (pme_loadbal_is_active(pme_loadbal))
1768 /* Do not permit counter reset while PME load
1769 * balancing is active. The only purpose for resetting
1770 * counters is to measure reliable performance data,
1771 * and that can't be done before balancing
1772 * completes.
1774 * TODO consider fixing this by delaying the reset
1775 * until after load balancing completes,
1776 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1777 gmx_fatal(FARGS, "PME tuning was still active when attempting to "
1778 "reset mdrun counters at step %" GMX_PRId64 ". Try "
1779 "resetting counters later in the run, e.g. with gmx "
1780 "mdrun -resetstep.", step);
1782 reset_all_counters(fplog, cr, step, &step_rel, ir, wcycle, nrnb, walltime_accounting,
1783 use_GPU(fr->nbv) ? fr->nbv : NULL);
1784 wcycle_set_reset_counters(wcycle, -1);
1785 if (!(cr->duty & DUTY_PME))
1787 /* Tell our PME node to reset its counters */
1788 gmx_pme_send_resetcounters(cr, step);
1790 /* Correct max_hours for the elapsed time */
1791 max_hours -= elapsed_time/(60.0*60.0);
1792 /* If mdrun -maxh -resethway was active, it can only trigger once */
1793 bResetCountersHalfMaxH = FALSE; /* TODO move this to where gs.sig[eglsRESETCOUNTERS] is set */
1794 /* Reset can only happen once, so clear the triggering flag. */
1795 gs.set[eglsRESETCOUNTERS] = 0;
1798 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1799 IMD_prep_energies_send_positions(ir->bIMD && MASTER(cr), bIMDstep, ir->imd, enerd, step, bCalcEner, wcycle);
1802 /* End of main MD loop */
1804 /* Closing TNG files can include compressing data. Therefore it is good to do that
1805 * before stopping the time measurements. */
1806 mdoutf_tng_close(outf);
1808 /* Stop measuring walltime */
1809 walltime_accounting_end(walltime_accounting);
1811 if (bRerunMD && MASTER(cr))
1813 close_trj(status);
1816 if (!(cr->duty & DUTY_PME))
1818 /* Tell the PME only node to finish */
1819 gmx_pme_send_finish(cr);
1822 if (MASTER(cr))
1824 if (ir->nstcalcenergy > 0 && !bRerunMD)
1826 print_ebin(mdoutf_get_fp_ene(outf), FALSE, FALSE, FALSE, fplog, step, t,
1827 eprAVER, FALSE, mdebin, fcd, groups, &(ir->opts));
1831 done_mdoutf(outf);
1833 if (bPMETune)
1835 pme_loadbal_done(pme_loadbal, cr, fplog, use_GPU(fr->nbv));
1838 if (shellfc && fplog)
1840 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1841 (nconverged*100.0)/step_rel);
1842 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1843 tcount/step_rel);
1846 if (repl_ex_nst > 0 && MASTER(cr))
1848 print_replica_exchange_statistics(fplog, repl_ex);
1851 /* IMD cleanup, if bIMD is TRUE. */
1852 IMD_finalize(ir->bIMD, ir->imd);
1854 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
1856 return 0;