2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015,2016, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
49 #include "thread_mpi/threads.h"
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_network.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/ewald/pme-load-balancing.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/gmxlib/md_logging.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/gpu_utils/gpu_utils.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed-forces/manage-threading.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/math/vectypes.h"
68 #include "gromacs/mdlib/compute_io.h"
69 #include "gromacs/mdlib/constr.h"
70 #include "gromacs/mdlib/ebin.h"
71 #include "gromacs/mdlib/force.h"
72 #include "gromacs/mdlib/force_flags.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/md_support.h"
75 #include "gromacs/mdlib/mdatoms.h"
76 #include "gromacs/mdlib/mdebin.h"
77 #include "gromacs/mdlib/mdoutf.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/nb_verlet.h"
80 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
81 #include "gromacs/mdlib/ns.h"
82 #include "gromacs/mdlib/shellfc.h"
83 #include "gromacs/mdlib/sighandler.h"
84 #include "gromacs/mdlib/sim_util.h"
85 #include "gromacs/mdlib/simulationsignal.h"
86 #include "gromacs/mdlib/tgroup.h"
87 #include "gromacs/mdlib/trajectory_writing.h"
88 #include "gromacs/mdlib/update.h"
89 #include "gromacs/mdlib/vcm.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdtypes/commrec.h"
92 #include "gromacs/mdtypes/df_history.h"
93 #include "gromacs/mdtypes/energyhistory.h"
94 #include "gromacs/mdtypes/fcdata.h"
95 #include "gromacs/mdtypes/forcerec.h"
96 #include "gromacs/mdtypes/group.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/interaction_const.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/mdtypes/mdatom.h"
101 #include "gromacs/mdtypes/state.h"
102 #include "gromacs/pbcutil/mshift.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/pulling/pull.h"
105 #include "gromacs/swap/swapcoords.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/timing/walltime_accounting.h"
108 #include "gromacs/topology/atoms.h"
109 #include "gromacs/topology/idef.h"
110 #include "gromacs/topology/mtop_util.h"
111 #include "gromacs/topology/topology.h"
112 #include "gromacs/trajectory/trajectoryframe.h"
113 #include "gromacs/utility/basedefinitions.h"
114 #include "gromacs/utility/cstringutil.h"
115 #include "gromacs/utility/fatalerror.h"
116 #include "gromacs/utility/real.h"
117 #include "gromacs/utility/smalloc.h"
124 #include "corewrap.h"
127 using gmx::SimulationSignaller
;
129 /*! \brief Check whether bonded interactions are missing, if appropriate
131 * \param[in] fplog Log file pointer
132 * \param[in] cr Communication object
133 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
134 * \param[in] top_global Global topology for the error message
135 * \param[in] top_local Local topology for the error message
136 * \param[in] state Global state for the error message
137 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
139 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
140 * is always set to false after exit.
142 static void checkNumberOfBondedInteractions(FILE *fplog
, t_commrec
*cr
, int totalNumberOfBondedInteractions
,
143 gmx_mtop_t
*top_global
, gmx_localtop_t
*top_local
, t_state
*state
,
144 bool *shouldCheckNumberOfBondedInteractions
)
146 if (*shouldCheckNumberOfBondedInteractions
)
148 if (totalNumberOfBondedInteractions
!= cr
->dd
->nbonded_global
)
150 dd_print_missing_interactions(fplog
, cr
, totalNumberOfBondedInteractions
, top_global
, top_local
, state
); // Does not return
152 *shouldCheckNumberOfBondedInteractions
= false;
156 static void reset_all_counters(FILE *fplog
, t_commrec
*cr
,
158 gmx_int64_t
*step_rel
, t_inputrec
*ir
,
159 gmx_wallcycle_t wcycle
, t_nrnb
*nrnb
,
160 gmx_walltime_accounting_t walltime_accounting
,
161 struct nonbonded_verlet_t
*nbv
)
163 char sbuf
[STEPSTRSIZE
];
165 /* Reset all the counters related to performance over the run */
166 md_print_warn(cr
, fplog
, "step %s: resetting all time and cycle counters\n",
167 gmx_step_str(step
, sbuf
));
171 nbnxn_gpu_reset_timings(nbv
);
175 wallcycle_stop(wcycle
, ewcRUN
);
176 wallcycle_reset_all(wcycle
);
177 if (DOMAINDECOMP(cr
))
179 reset_dd_statistics_counters(cr
->dd
);
182 ir
->init_step
+= *step_rel
;
183 ir
->nsteps
-= *step_rel
;
185 wallcycle_start(wcycle
, ewcRUN
);
186 walltime_accounting_start(walltime_accounting
);
187 print_date_and_time(fplog
, cr
->nodeid
, "Restarted time", gmx_gettime());
191 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
192 int nfile, const t_filenm fnm[],
193 const gmx_output_env_t *oenv, gmx_bool bVerbose,
195 gmx_vsite_t *vsite, gmx_constr_t constr,
197 t_inputrec *inputrec,
198 gmx_mtop_t *top_global, t_fcdata *fcd,
199 t_state *state_global,
201 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
204 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
205 real cpt_period, real max_hours,
208 gmx_walltime_accounting_t walltime_accounting)
210 double gmx::do_md(FILE *fplog
, t_commrec
*cr
, int nfile
, const t_filenm fnm
[],
211 const gmx_output_env_t
*oenv
, gmx_bool bVerbose
,
213 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
214 int stepout
, t_inputrec
*ir
,
215 gmx_mtop_t
*top_global
,
217 t_state
*state_global
,
219 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
220 gmx_edsam_t ed
, t_forcerec
*fr
,
221 int repl_ex_nst
, int repl_ex_nex
, int repl_ex_seed
,
222 real cpt_period
, real max_hours
,
225 gmx_walltime_accounting_t walltime_accounting
)
227 gmx_mdoutf_t outf
= NULL
;
228 gmx_int64_t step
, step_rel
;
230 double t
, t0
, lam0
[efptNR
];
231 gmx_bool bGStatEveryStep
, bGStat
, bCalcVir
, bCalcEnerStep
, bCalcEner
;
232 gmx_bool bNS
, bNStList
, bSimAnn
, bStopCM
, bRerunMD
,
233 bFirstStep
, startingFromCheckpoint
, bInitStep
, bLastStep
= FALSE
,
234 bBornRadii
, bUsingEnsembleRestraints
;
235 gmx_bool bDoDHDL
= FALSE
, bDoFEP
= FALSE
, bDoExpanded
= FALSE
;
236 gmx_bool do_ene
, do_log
, do_verbose
, bRerunWarnNoV
= TRUE
,
237 bForceUpdate
= FALSE
, bCPT
;
238 gmx_bool bMasterState
;
239 int force_flags
, cglo_flags
;
240 tensor force_vir
, shake_vir
, total_vir
, tmp_vir
, pres
;
247 gmx_repl_ex_t repl_ex
= NULL
;
250 t_mdebin
*mdebin
= NULL
;
251 t_state
*state
= NULL
;
252 gmx_enerdata_t
*enerd
;
254 gmx_global_stat_t gstat
;
255 gmx_update_t
*upd
= NULL
;
256 t_graph
*graph
= NULL
;
257 gmx_groups_t
*groups
;
258 gmx_ekindata_t
*ekind
;
259 gmx_shellfc_t
*shellfc
;
260 gmx_bool bSumEkinhOld
, bDoReplEx
, bExchanged
, bNeedRepartition
;
261 gmx_bool bResetCountersHalfMaxH
= FALSE
;
262 gmx_bool bTemp
, bPres
, bTrotter
;
271 real saved_conserved_quantity
= 0;
275 char sbuf
[STEPSTRSIZE
], sbuf2
[STEPSTRSIZE
];
276 int handled_stop_condition
= gmx_stop_cond_none
; /* compare to get_stop_condition*/
279 /* PME load balancing data for GPU kernels */
280 pme_load_balancing_t
*pme_loadbal
= NULL
;
281 gmx_bool bPMETune
= FALSE
;
282 gmx_bool bPMETunePrinting
= FALSE
;
285 gmx_bool bIMDstep
= FALSE
;
286 gmx_membed_t
*membed
= NULL
;
289 /* Temporary addition for FAHCORE checkpointing */
292 /* Domain decomposition could incorrectly miss a bonded
293 interaction, but checking for that requires a global
294 communication stage, which does not otherwise happen in DD
295 code. So we do that alongside the first global energy reduction
296 after a new DD is made. These variables handle whether the
297 check happens, and the result it returns. */
298 bool shouldCheckNumberOfBondedInteractions
= false;
299 int totalNumberOfBondedInteractions
= -1;
301 SimulationSignals signals
;
302 // Most global communnication stages don't propagate mdrun
303 // signals, and will use this object to achieve that.
304 SimulationSignaller
nullSignaller(nullptr, nullptr, false, false);
306 /* Check for special mdrun options */
307 bRerunMD
= (Flags
& MD_RERUN
);
308 if (Flags
& MD_RESETCOUNTERSHALFWAY
)
312 /* Signal to reset the counters half the simulation steps. */
313 wcycle_set_reset_counters(wcycle
, ir
->nsteps
/2);
315 /* Signal to reset the counters halfway the simulation time. */
316 bResetCountersHalfMaxH
= (max_hours
> 0);
319 /* md-vv uses averaged full step velocities for T-control
320 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
321 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
322 bTrotter
= (EI_VV(ir
->eI
) && (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
) || inputrecNvtTrotter(ir
)));
326 /* Since we don't know if the frames read are related in any way,
327 * rebuild the neighborlist at every step.
330 ir
->nstcalcenergy
= 1;
334 nstglobalcomm
= check_nstglobalcomm(fplog
, cr
, nstglobalcomm
, ir
);
335 bGStatEveryStep
= (nstglobalcomm
== 1);
339 ir
->nstxout_compressed
= 0;
341 groups
= &top_global
->groups
;
343 if (opt2bSet("-membed", nfile
, fnm
))
347 fprintf(stderr
, "Initializing membed");
349 /* Note that membed cannot work in parallel because mtop is
350 * changed here. Fix this if we ever want to make it run with
352 membed
= init_membed(fplog
, nfile
, fnm
, top_global
, ir
, state_global
, cr
, &cpt_period
);
355 if (ir
->eSwapCoords
!= eswapNO
)
357 /* Initialize ion swapping code */
358 init_swapcoords(fplog
, bVerbose
, ir
, opt2fn_master("-swap", nfile
, fnm
, cr
),
359 top_global
, state_global
->x
, state_global
->box
, &state_global
->swapstate
, cr
, oenv
, Flags
);
363 init_md(fplog
, cr
, ir
, oenv
, &t
, &t0
, state_global
->lambda
,
364 &(state_global
->fep_state
), lam0
,
365 nrnb
, top_global
, &upd
,
366 nfile
, fnm
, &outf
, &mdebin
,
367 force_vir
, shake_vir
, mu_tot
, &bSimAnn
, &vcm
, Flags
, wcycle
);
369 clear_mat(total_vir
);
371 /* Energy terms and groups */
373 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
375 if (DOMAINDECOMP(cr
))
381 snew(f
, top_global
->natoms
);
384 /* Kinetic energy data */
386 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind
);
387 /* Copy the cos acceleration to the groups struct */
388 ekind
->cosacc
.cos_accel
= ir
->cos_accel
;
390 gstat
= global_stat_init(ir
);
392 /* Check for polarizable models and flexible constraints */
393 shellfc
= init_shell_flexcon(fplog
,
394 top_global
, n_flexible_constraints(constr
),
395 ir
->nstcalcenergy
, DOMAINDECOMP(cr
));
397 if (shellfc
&& ir
->nstcalcenergy
!= 1)
399 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
);
401 if (shellfc
&& DOMAINDECOMP(cr
))
403 gmx_fatal(FARGS
, "Shell particles are not implemented with domain decomposition, use a single rank");
406 if (inputrecDeform(ir
))
408 tMPI_Thread_mutex_lock(&deform_init_box_mutex
);
409 set_deform_reference_box(upd
,
410 deform_init_init_step_tpx
,
411 deform_init_box_tpx
);
412 tMPI_Thread_mutex_unlock(&deform_init_box_mutex
);
416 double io
= compute_io(ir
, top_global
->natoms
, groups
, mdebin
->ebin
->nener
, 1);
417 if ((io
> 2000) && MASTER(cr
))
420 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
425 if (DOMAINDECOMP(cr
))
427 top
= dd_init_local_top(top_global
);
430 dd_init_local_state(cr
->dd
, state_global
, state
);
434 top
= gmx_mtop_generate_local_top(top_global
, ir
->efep
!= efepNO
);
436 forcerec_set_excl_load(fr
, top
);
438 state
= serial_init_local_state(state_global
);
440 atoms2md(top_global
, ir
, 0, NULL
, top_global
->natoms
, mdatoms
);
444 set_vsite_top(vsite
, top
, mdatoms
, cr
);
447 if (ir
->ePBC
!= epbcNONE
&& !fr
->bMolPBC
)
449 graph
= mk_graph(fplog
, &(top
->idef
), 0, top_global
->natoms
, FALSE
, FALSE
);
454 make_local_shells(cr
, mdatoms
, shellfc
);
457 setup_bonded_threading(fr
, &top
->idef
);
459 update_realloc(upd
, state
->nalloc
);
462 /* Set up interactive MD (IMD) */
463 init_IMD(ir
, cr
, top_global
, fplog
, ir
->nstcalcenergy
, state_global
->x
,
464 nfile
, fnm
, oenv
, imdport
, Flags
);
466 if (DOMAINDECOMP(cr
))
468 /* Distribute the charge groups over the nodes from the master node */
469 dd_partition_system(fplog
, ir
->init_step
, cr
, TRUE
, 1,
470 state_global
, top_global
, ir
,
471 state
, &f
, mdatoms
, top
, fr
,
474 shouldCheckNumberOfBondedInteractions
= true;
475 update_realloc(upd
, state
->nalloc
);
478 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
480 startingFromCheckpoint
= Flags
& MD_STARTFROMCPT
;
484 init_expanded_ensemble(startingFromCheckpoint
, ir
, &state
->dfhist
);
489 if (startingFromCheckpoint
)
491 /* Update mdebin with energy history if appending to output files */
492 if (Flags
& MD_APPENDFILES
)
494 restore_energyhistory_from_state(mdebin
, state_global
->enerhist
);
498 /* We might have read an energy history from checkpoint,
499 * free the allocated memory and reset the counts.
501 done_energyhistory(state_global
->enerhist
);
502 init_energyhistory(state_global
->enerhist
);
505 /* Set the initial energy history in state by updating once */
506 update_energyhistory(state_global
->enerhist
, mdebin
);
509 /* Initialize constraints */
510 if (constr
&& !DOMAINDECOMP(cr
))
512 set_constraints(constr
, top
, ir
, mdatoms
, cr
);
515 if (repl_ex_nst
> 0 && MASTER(cr
))
517 repl_ex
= init_replica_exchange(fplog
, cr
->ms
, state_global
, ir
,
518 repl_ex_nst
, repl_ex_nex
, repl_ex_seed
);
521 /* PME tuning is only supported with PME for Coulomb. Is is not supported
522 * with only LJ PME, or for reruns.
524 bPMETune
= ((Flags
& MD_TUNEPME
) && EEL_PME(fr
->eeltype
) && !bRerunMD
&&
525 !(Flags
& MD_REPRODUCIBLE
));
528 pme_loadbal_init(&pme_loadbal
, cr
, fplog
, ir
, state
->box
,
529 fr
->ic
, fr
->pmedata
, use_GPU(fr
->nbv
),
533 if (!ir
->bContinuation
&& !bRerunMD
)
535 if (mdatoms
->cFREEZE
&& (state
->flags
& (1<<estV
)))
537 /* Set the velocities of frozen particles to zero */
538 for (i
= 0; i
< mdatoms
->homenr
; i
++)
540 for (m
= 0; m
< DIM
; m
++)
542 if (ir
->opts
.nFreeze
[mdatoms
->cFREEZE
[i
]][m
])
552 /* Constrain the initial coordinates and velocities */
553 do_constrain_first(fplog
, constr
, ir
, mdatoms
, state
,
558 /* Construct the virtual sites for the initial configuration */
559 construct_vsites(vsite
, state
->x
, ir
->delta_t
, NULL
,
560 top
->idef
.iparams
, top
->idef
.il
,
561 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
565 if (ir
->efep
!= efepNO
)
567 /* Set free energy calculation frequency as the greatest common
568 * denominator of nstdhdl and repl_ex_nst. */
569 nstfep
= ir
->fepvals
->nstdhdl
;
572 nstfep
= gmx_greatest_common_divisor(ir
->expandedvals
->nstexpanded
, nstfep
);
576 nstfep
= gmx_greatest_common_divisor(repl_ex_nst
, nstfep
);
580 /* Be REALLY careful about what flags you set here. You CANNOT assume
581 * this is the first step, since we might be restarting from a checkpoint,
582 * and in that case we should not do any modifications to the state.
584 bStopCM
= (ir
->comm_mode
!= ecmNO
&& !ir
->bContinuation
);
586 if (Flags
& MD_READ_EKIN
)
588 restore_ekinstate_from_state(cr
, ekind
, &state_global
->ekinstate
);
591 cglo_flags
= (CGLO_TEMPERATURE
| CGLO_GSTAT
592 | (bStopCM
? CGLO_STOPCM
: 0)
593 | (EI_VV(ir
->eI
) ? CGLO_PRESSURE
: 0)
594 | (EI_VV(ir
->eI
) ? CGLO_CONSTRAINT
: 0)
595 | ((Flags
& MD_READ_EKIN
) ? CGLO_READEKIN
: 0));
597 bSumEkinhOld
= FALSE
;
598 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
599 NULL
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
600 constr
, &nullSignaller
, state
->box
,
601 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags
602 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
603 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
604 top_global
, top
, state
,
605 &shouldCheckNumberOfBondedInteractions
);
606 if (ir
->eI
== eiVVAK
)
608 /* a second call to get the half step temperature initialized as well */
609 /* we do the same call as above, but turn the pressure off -- internally to
610 compute_globals, this is recognized as a velocity verlet half-step
611 kinetic energy calculation. This minimized excess variables, but
612 perhaps loses some logic?*/
614 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
615 NULL
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
616 constr
, &nullSignaller
, state
->box
,
618 cglo_flags
&~(CGLO_STOPCM
| CGLO_PRESSURE
));
621 /* Calculate the initial half step temperature, and save the ekinh_old */
622 if (!(Flags
& MD_STARTFROMCPT
))
624 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
626 copy_mat(ekind
->tcstat
[i
].ekinh
, ekind
->tcstat
[i
].ekinh_old
);
631 enerd
->term
[F_TEMP
] *= 2; /* result of averages being done over previous and current step,
632 and there is no previous step */
635 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
636 temperature control */
637 trotter_seq
= init_npt_vars(ir
, state
, &MassQ
, bTrotter
);
641 if (constr
&& !ir
->bContinuation
&& ir
->eConstrAlg
== econtLINCS
)
644 "RMS relative constraint deviation after constraining: %.2e\n",
645 constr_rmsd(constr
));
647 if (EI_STATE_VELOCITY(ir
->eI
))
649 fprintf(fplog
, "Initial temperature: %g K\n", enerd
->term
[F_TEMP
]);
653 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
654 " input trajectory '%s'\n\n",
655 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
658 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
659 "run input file,\nwhich may not correspond to the time "
660 "needed to process input trajectory.\n\n");
666 fprintf(stderr
, "starting mdrun '%s'\n",
667 *(top_global
->name
));
670 sprintf(tbuf
, "%8.1f", (ir
->init_step
+ir
->nsteps
)*ir
->delta_t
);
674 sprintf(tbuf
, "%s", "infinite");
676 if (ir
->init_step
> 0)
678 fprintf(stderr
, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
679 gmx_step_str(ir
->init_step
+ir
->nsteps
, sbuf
), tbuf
,
680 gmx_step_str(ir
->init_step
, sbuf2
),
681 ir
->init_step
*ir
->delta_t
);
685 fprintf(stderr
, "%s steps, %s ps.\n",
686 gmx_step_str(ir
->nsteps
, sbuf
), tbuf
);
689 fprintf(fplog
, "\n");
692 walltime_accounting_start(walltime_accounting
);
693 wallcycle_start(wcycle
, ewcRUN
);
694 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
696 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
698 chkpt_ret
= fcCheckPointParallel( cr
->nodeid
,
702 gmx_fatal( 3, __FILE__
, __LINE__
, "Checkpoint error on step %d\n", 0 );
706 /***********************************************************
710 ************************************************************/
712 /* if rerunMD then read coordinates and velocities from input trajectory */
715 if (getenv("GMX_FORCE_UPDATE"))
723 bLastStep
= !read_first_frame(oenv
, &status
,
724 opt2fn("-rerun", nfile
, fnm
),
725 &rerun_fr
, TRX_NEED_X
| TRX_READ_V
);
726 if (rerun_fr
.natoms
!= top_global
->natoms
)
729 "Number of atoms in trajectory (%d) does not match the "
730 "run input file (%d)\n",
731 rerun_fr
.natoms
, top_global
->natoms
);
733 if (ir
->ePBC
!= epbcNONE
)
737 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
);
739 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < gmx::square(fr
->rlist
))
741 gmx_fatal(FARGS
, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
748 rerun_parallel_comm(cr
, &rerun_fr
, &bLastStep
);
751 if (ir
->ePBC
!= epbcNONE
)
753 /* Set the shift vectors.
754 * Necessary here when have a static box different from the tpr box.
756 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
760 /* loop over MD steps or if rerunMD to end of input trajectory */
762 /* Skip the first Nose-Hoover integration when we get the state from tpx */
763 bInitStep
= !startingFromCheckpoint
|| EI_VV(ir
->eI
);
764 bSumEkinhOld
= FALSE
;
766 bNeedRepartition
= FALSE
;
767 // TODO This implementation of ensemble orientation restraints is nasty because
768 // a user can't just do multi-sim with single-sim orientation restraints.
769 bUsingEnsembleRestraints
= (fcd
->disres
.nsystems
> 1) || (cr
->ms
&& fcd
->orires
.nr
);
772 // Replica exchange and ensemble restraints need all
773 // simulations to remain synchronized, so they need
774 // checkpoints and stop conditions to act on the same step, so
775 // the propagation of such signals must take place between
776 // simulations, not just within simulations.
777 bool checkpointIsLocal
= (repl_ex_nst
<= 0) && !bUsingEnsembleRestraints
;
778 bool stopConditionIsLocal
= (repl_ex_nst
<= 0) && !bUsingEnsembleRestraints
;
779 bool resetCountersIsLocal
= true;
780 signals
[eglsCHKPT
] = SimulationSignal(checkpointIsLocal
);
781 signals
[eglsSTOPCOND
] = SimulationSignal(stopConditionIsLocal
);
782 signals
[eglsRESETCOUNTERS
] = SimulationSignal(resetCountersIsLocal
);
785 step
= ir
->init_step
;
788 // TODO extract this to new multi-simulation module
789 if (MASTER(cr
) && MULTISIM(cr
) && (repl_ex_nst
<= 0 ))
791 if (!multisim_int_all_are_equal(cr
->ms
, ir
->nsteps
))
793 md_print_info(cr
, fplog
,
794 "Note: The number of steps is not consistent across multi simulations,\n"
795 "but we are proceeding anyway!\n");
797 if (!multisim_int_all_are_equal(cr
->ms
, ir
->init_step
))
799 md_print_info(cr
, fplog
,
800 "Note: The initial step is not consistent across multi simulations,\n"
801 "but we are proceeding anyway!\n");
805 /* and stop now if we should */
806 bLastStep
= (bLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
810 /* Determine if this is a neighbor search step */
811 bNStList
= (ir
->nstlist
> 0 && step
% ir
->nstlist
== 0);
813 if (bPMETune
&& bNStList
)
815 /* PME grid + cut-off optimization with GPUs or PME nodes */
816 pme_loadbal_do(pme_loadbal
, cr
,
817 (bVerbose
&& MASTER(cr
)) ? stderr
: NULL
,
825 wallcycle_start(wcycle
, ewcSTEP
);
831 step
= rerun_fr
.step
;
832 step_rel
= step
- ir
->init_step
;
845 bLastStep
= (step_rel
== ir
->nsteps
);
846 t
= t0
+ step
*ir
->delta_t
;
849 // TODO Refactor this, so that nstfep does not need a default value of zero
850 if (ir
->efep
!= efepNO
|| ir
->bSimTemp
)
852 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
853 requiring different logic. */
855 set_current_lambdas(step
, ir
->fepvals
, bRerunMD
, &rerun_fr
, state_global
, state
, lam0
);
856 bDoDHDL
= do_per_step(step
, ir
->fepvals
->nstdhdl
);
857 bDoFEP
= ((ir
->efep
!= efepNO
) && do_per_step(step
, nstfep
));
858 bDoExpanded
= (do_per_step(step
, ir
->expandedvals
->nstexpanded
)
859 && (ir
->bExpanded
) && (step
> 0) && (!startingFromCheckpoint
));
862 bDoReplEx
= ((repl_ex_nst
> 0) && (step
> 0) && !bLastStep
&&
863 do_per_step(step
, repl_ex_nst
));
867 update_annealing_target_temp(ir
, t
, upd
);
872 if (!DOMAINDECOMP(cr
) || MASTER(cr
))
874 for (i
= 0; i
< state_global
->natoms
; i
++)
876 copy_rvec(rerun_fr
.x
[i
], state_global
->x
[i
]);
880 for (i
= 0; i
< state_global
->natoms
; i
++)
882 copy_rvec(rerun_fr
.v
[i
], state_global
->v
[i
]);
887 for (i
= 0; i
< state_global
->natoms
; i
++)
889 clear_rvec(state_global
->v
[i
]);
893 fprintf(stderr
, "\nWARNING: Some frames do not contain velocities.\n"
894 " Ekin, temperature and pressure are incorrect,\n"
895 " the virial will be incorrect when constraints are present.\n"
897 bRerunWarnNoV
= FALSE
;
901 copy_mat(rerun_fr
.box
, state_global
->box
);
902 copy_mat(state_global
->box
, state
->box
);
904 if (vsite
&& (Flags
& MD_RERUN_VSITE
))
906 if (DOMAINDECOMP(cr
))
908 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
912 /* Following is necessary because the graph may get out of sync
913 * with the coordinates if we only have every N'th coordinate set
915 mk_mshift(fplog
, graph
, fr
->ePBC
, state
->box
, state
->x
);
916 shift_self(graph
, state
->box
, state
->x
);
918 construct_vsites(vsite
, state
->x
, ir
->delta_t
, state
->v
,
919 top
->idef
.iparams
, top
->idef
.il
,
920 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
923 unshift_self(graph
, state
->box
, state
->x
);
928 /* Stop Center of Mass motion */
929 bStopCM
= (ir
->comm_mode
!= ecmNO
&& do_per_step(step
, ir
->nstcomm
));
933 /* for rerun MD always do Neighbour Searching */
934 bNS
= (bFirstStep
|| ir
->nstlist
!= 0);
938 /* Determine whether or not to do Neighbour Searching */
939 bNS
= (bFirstStep
|| bNStList
|| bExchanged
|| bNeedRepartition
);
942 /* < 0 means stop at next step, > 0 means stop at next NS step */
943 if ( (signals
[eglsSTOPCOND
].set
< 0) ||
944 ( (signals
[eglsSTOPCOND
].set
> 0 ) && ( bNS
|| ir
->nstlist
== 0)))
949 /* Determine whether or not to update the Born radii if doing GB */
950 bBornRadii
= bFirstStep
;
951 if (ir
->implicit_solvent
&& (step
% ir
->nstgbradii
== 0))
956 /* do_log triggers energy and virial calculation. Because this leads
957 * to different code paths, forces can be different. Thus for exact
958 * continuation we should avoid extra log output.
959 * Note that the || bLastStep can result in non-exact continuation
960 * beyond the last step. But we don't consider that to be an issue.
962 do_log
= do_per_step(step
, ir
->nstlog
) || (bFirstStep
&& !startingFromCheckpoint
) || bLastStep
|| bRerunMD
;
963 do_verbose
= bVerbose
&&
964 (step
% stepout
== 0 || bFirstStep
|| bLastStep
|| bRerunMD
);
966 if (bNS
&& !(bFirstStep
&& ir
->bContinuation
&& !bRerunMD
))
974 bMasterState
= FALSE
;
975 /* Correct the new box if it is too skewed */
976 if (inputrecDynamicBox(ir
))
978 if (correct_box(fplog
, step
, state
->box
, graph
))
983 if (DOMAINDECOMP(cr
) && bMasterState
)
985 dd_collect_state(cr
->dd
, state
, state_global
);
989 if (DOMAINDECOMP(cr
))
991 /* Repartition the domain decomposition */
992 dd_partition_system(fplog
, step
, cr
,
993 bMasterState
, nstglobalcomm
,
994 state_global
, top_global
, ir
,
995 state
, &f
, mdatoms
, top
, fr
,
998 do_verbose
&& !bPMETunePrinting
);
999 shouldCheckNumberOfBondedInteractions
= true;
1000 update_realloc(upd
, state
->nalloc
);
1004 if (MASTER(cr
) && do_log
)
1006 print_ebin_header(fplog
, step
, t
); /* can we improve the information printed here? */
1009 if (ir
->efep
!= efepNO
)
1011 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
1014 if ((bRerunMD
&& rerun_fr
.bV
) || bExchanged
)
1017 /* We need the kinetic energy at minus the half step for determining
1018 * the full step kinetic energy and possibly for T-coupling.*/
1019 /* This may not be quite working correctly yet . . . . */
1020 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1021 wcycle
, enerd
, NULL
, NULL
, NULL
, NULL
, mu_tot
,
1022 constr
, &nullSignaller
, state
->box
,
1023 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1024 CGLO_GSTAT
| CGLO_TEMPERATURE
| CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
);
1025 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
1026 top_global
, top
, state
,
1027 &shouldCheckNumberOfBondedInteractions
);
1029 clear_mat(force_vir
);
1031 /* We write a checkpoint at this MD step when:
1032 * either at an NS step when we signalled through gs,
1033 * or at the last step (but not when we do not want confout),
1034 * but never at the first step or with rerun.
1036 bCPT
= (((signals
[eglsCHKPT
].set
&& (bNS
|| ir
->nstlist
== 0)) ||
1037 (bLastStep
&& (Flags
& MD_CONFOUT
))) &&
1038 step
> ir
->init_step
&& !bRerunMD
);
1041 signals
[eglsCHKPT
].set
= 0;
1044 /* Determine the energy and pressure:
1045 * at nstcalcenergy steps and at energy output steps (set below).
1047 if (EI_VV(ir
->eI
) && (!bInitStep
))
1049 /* for vv, the first half of the integration actually corresponds
1050 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1051 but the virial needs to be calculated on both the current step and the 'next' step. Future
1052 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1054 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1055 bCalcEnerStep
= do_per_step(step
- 1, ir
->nstcalcenergy
);
1056 bCalcVir
= bCalcEnerStep
||
1057 (ir
->epc
!= epcNO
&& (do_per_step(step
, ir
->nstpcouple
) || do_per_step(step
-1, ir
->nstpcouple
)));
1061 bCalcEnerStep
= do_per_step(step
, ir
->nstcalcenergy
);
1062 bCalcVir
= bCalcEnerStep
||
1063 (ir
->epc
!= epcNO
&& do_per_step(step
, ir
->nstpcouple
));
1065 bCalcEner
= bCalcEnerStep
;
1067 do_ene
= (do_per_step(step
, ir
->nstenergy
) || bLastStep
|| bRerunMD
);
1069 if (do_ene
|| do_log
|| bDoReplEx
)
1075 /* Do we need global communication ? */
1076 bGStat
= (bCalcVir
|| bCalcEner
|| bStopCM
||
1077 do_per_step(step
, nstglobalcomm
) ||
1078 (EI_VV(ir
->eI
) && inputrecNvtTrotter(ir
) && do_per_step(step
-1, nstglobalcomm
)));
1080 force_flags
= (GMX_FORCE_STATECHANGED
|
1081 ((inputrecDynamicBox(ir
) || bRerunMD
) ? GMX_FORCE_DYNAMICBOX
: 0) |
1082 GMX_FORCE_ALLFORCES
|
1083 (bCalcVir
? GMX_FORCE_VIRIAL
: 0) |
1084 (bCalcEner
? GMX_FORCE_ENERGY
: 0) |
1085 (bDoFEP
? GMX_FORCE_DHDL
: 0)
1090 /* Now is the time to relax the shells */
1091 relax_shell_flexcon(fplog
, cr
, bVerbose
, step
,
1092 ir
, bNS
, force_flags
, top
,
1094 state
, f
, force_vir
, mdatoms
,
1095 nrnb
, wcycle
, graph
, groups
,
1096 shellfc
, fr
, bBornRadii
, t
, mu_tot
,
1097 vsite
, mdoutf_get_fp_field(outf
));
1101 /* The coordinates (x) are shifted (to get whole molecules)
1103 * This is parallellized as well, and does communication too.
1104 * Check comments in sim_util.c
1106 do_force(fplog
, cr
, ir
, step
, nrnb
, wcycle
, top
, groups
,
1107 state
->box
, state
->x
, &state
->hist
,
1108 f
, force_vir
, mdatoms
, enerd
, fcd
,
1109 state
->lambda
, graph
,
1110 fr
, vsite
, mu_tot
, t
, mdoutf_get_fp_field(outf
), ed
, bBornRadii
,
1111 (bNS
? GMX_FORCE_NS
: 0) | force_flags
);
1114 if (EI_VV(ir
->eI
) && !startingFromCheckpoint
&& !bRerunMD
)
1115 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1119 wallcycle_start(wcycle
, ewcUPDATE
);
1120 if (ir
->eI
== eiVV
&& bInitStep
)
1122 /* if using velocity verlet with full time step Ekin,
1123 * take the first half step only to compute the
1124 * virial for the first step. From there,
1125 * revert back to the initial coordinates
1126 * so that the input is actually the initial step.
1128 snew(vbuf
, state
->natoms
);
1129 copy_rvecn(state
->v
, vbuf
, 0, state
->natoms
); /* should make this better for parallelizing? */
1133 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1134 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ1
);
1137 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1138 ekind
, M
, upd
, etrtVELOCITY1
,
1141 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
) /* Why is rerun_fr.bV here? Unclear. */
1143 wallcycle_stop(wcycle
, ewcUPDATE
);
1144 update_constraints(fplog
, step
, NULL
, ir
, mdatoms
,
1145 state
, fr
->bMolPBC
, graph
, f
,
1146 &top
->idef
, shake_vir
,
1147 cr
, nrnb
, wcycle
, upd
, constr
,
1149 wallcycle_start(wcycle
, ewcUPDATE
);
1153 /* Need to unshift here if a do_force has been
1154 called in the previous step */
1155 unshift_self(graph
, state
->box
, state
->x
);
1157 /* if VV, compute the pressure and constraints */
1158 /* For VV2, we strictly only need this if using pressure
1159 * control, but we really would like to have accurate pressures
1161 * Think about ways around this in the future?
1162 * For now, keep this choice in comments.
1164 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1165 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1167 bTemp
= ((ir
->eI
== eiVV
&& (!bInitStep
)) || (ir
->eI
== eiVVAK
));
1168 if (bCalcEner
&& ir
->eI
== eiVVAK
)
1170 bSumEkinhOld
= TRUE
;
1172 /* for vv, the first half of the integration actually corresponds to the previous step.
1173 So we need information from the last step in the first half of the integration */
1174 if (bGStat
|| do_per_step(step
-1, nstglobalcomm
))
1176 wallcycle_stop(wcycle
, ewcUPDATE
);
1177 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1178 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1179 constr
, &nullSignaller
, state
->box
,
1180 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1181 (bGStat
? CGLO_GSTAT
: 0)
1183 | (bTemp
? CGLO_TEMPERATURE
: 0)
1184 | (bPres
? CGLO_PRESSURE
: 0)
1185 | (bPres
? CGLO_CONSTRAINT
: 0)
1186 | (bStopCM
? CGLO_STOPCM
: 0)
1187 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
1190 /* explanation of above:
1191 a) We compute Ekin at the full time step
1192 if 1) we are using the AveVel Ekin, and it's not the
1193 initial step, or 2) if we are using AveEkin, but need the full
1194 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1195 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1196 EkinAveVel because it's needed for the pressure */
1197 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
1198 top_global
, top
, state
,
1199 &shouldCheckNumberOfBondedInteractions
);
1200 wallcycle_start(wcycle
, ewcUPDATE
);
1202 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1207 m_add(force_vir
, shake_vir
, total_vir
); /* we need the un-dispersion corrected total vir here */
1208 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ2
);
1210 copy_mat(shake_vir
, state
->svir_prev
);
1211 copy_mat(force_vir
, state
->fvir_prev
);
1212 if (inputrecNvtTrotter(ir
) && ir
->eI
== eiVV
)
1214 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1215 enerd
->term
[F_TEMP
] = sum_ekin(&(ir
->opts
), ekind
, NULL
, (ir
->eI
== eiVV
), FALSE
);
1216 enerd
->term
[F_EKIN
] = trace(ekind
->ekin
);
1219 else if (bExchanged
)
1221 wallcycle_stop(wcycle
, ewcUPDATE
);
1222 /* We need the kinetic energy at minus the half step for determining
1223 * the full step kinetic energy and possibly for T-coupling.*/
1224 /* This may not be quite working correctly yet . . . . */
1225 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1226 wcycle
, enerd
, NULL
, NULL
, NULL
, NULL
, mu_tot
,
1227 constr
, &nullSignaller
, state
->box
,
1228 NULL
, &bSumEkinhOld
,
1229 CGLO_GSTAT
| CGLO_TEMPERATURE
);
1230 wallcycle_start(wcycle
, ewcUPDATE
);
1233 /* if it's the initial step, we performed this first step just to get the constraint virial */
1234 if (ir
->eI
== eiVV
&& bInitStep
)
1236 copy_rvecn(vbuf
, state
->v
, 0, state
->natoms
);
1239 wallcycle_stop(wcycle
, ewcUPDATE
);
1242 /* compute the conserved quantity */
1245 saved_conserved_quantity
= compute_conserved_from_auxiliary(ir
, state
, &MassQ
);
1248 last_ekin
= enerd
->term
[F_EKIN
];
1250 if ((ir
->eDispCorr
!= edispcEnerPres
) && (ir
->eDispCorr
!= edispcAllEnerPres
))
1252 saved_conserved_quantity
-= enerd
->term
[F_DISPCORR
];
1254 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1255 if (ir
->efep
!= efepNO
&& !bRerunMD
)
1257 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1261 /* ######## END FIRST UPDATE STEP ############## */
1262 /* ######## If doing VV, we now have v(dt) ###### */
1265 /* perform extended ensemble sampling in lambda - we don't
1266 actually move to the new state before outputting
1267 statistics, but if performing simulated tempering, we
1268 do update the velocities and the tau_t. */
1270 lamnew
= ExpandedEnsembleDynamics(fplog
, ir
, enerd
, state
, &MassQ
, state
->fep_state
, &state
->dfhist
, step
, state
->v
, mdatoms
);
1271 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1272 copy_df_history(&state_global
->dfhist
, &state
->dfhist
);
1275 /* Now we have the energies and forces corresponding to the
1276 * coordinates at time t. We must output all of this before
1279 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
1280 ir
, state
, state_global
, top_global
, fr
,
1281 outf
, mdebin
, ekind
, f
,
1283 bCPT
, bRerunMD
, bLastStep
, (Flags
& MD_CONFOUT
),
1285 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1286 bIMDstep
= do_IMD(ir
->bIMD
, step
, cr
, bNS
, state
->box
, state
->x
, ir
, t
, wcycle
);
1288 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1289 if (startingFromCheckpoint
&& bTrotter
)
1291 copy_mat(state
->svir_prev
, shake_vir
);
1292 copy_mat(state
->fvir_prev
, force_vir
);
1295 elapsed_time
= walltime_accounting_get_current_elapsed_time(walltime_accounting
);
1297 /* Check whether everything is still allright */
1298 if (((int)gmx_get_stop_condition() > handled_stop_condition
)
1304 int nsteps_stop
= -1;
1306 /* this just makes signals[].sig compatible with the hack
1307 of sending signals around by MPI_Reduce together with
1309 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns
)
1311 signals
[eglsSTOPCOND
].sig
= 1;
1312 nsteps_stop
= std::max(ir
->nstlist
, 2*nstglobalcomm
);
1314 if (gmx_get_stop_condition() == gmx_stop_cond_next
)
1316 signals
[eglsSTOPCOND
].sig
= -1;
1317 nsteps_stop
= nstglobalcomm
+ 1;
1322 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1323 gmx_get_signal_name(), nsteps_stop
);
1327 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1328 gmx_get_signal_name(), nsteps_stop
);
1330 handled_stop_condition
= (int)gmx_get_stop_condition();
1332 else if (MASTER(cr
) && (bNS
|| ir
->nstlist
<= 0) &&
1333 (max_hours
> 0 && elapsed_time
> max_hours
*60.0*60.0*0.99) &&
1334 signals
[eglsSTOPCOND
].sig
== 0 && signals
[eglsSTOPCOND
].set
== 0)
1336 /* Signal to terminate the run */
1337 signals
[eglsSTOPCOND
].sig
= 1;
1340 fprintf(fplog
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1342 fprintf(stderr
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1345 if (bResetCountersHalfMaxH
&& MASTER(cr
) &&
1346 elapsed_time
> max_hours
*60.0*60.0*0.495)
1348 /* Set flag that will communicate the signal to all ranks in the simulation */
1349 signals
[eglsRESETCOUNTERS
].sig
= 1;
1352 /* In parallel we only have to check for checkpointing in steps
1353 * where we do global communication,
1354 * otherwise the other nodes don't know.
1356 if (MASTER(cr
) && ((bGStat
|| !PAR(cr
)) &&
1359 elapsed_time
>= nchkpt
*cpt_period
*60.0)) &&
1360 signals
[eglsCHKPT
].set
== 0)
1362 signals
[eglsCHKPT
].sig
= 1;
1365 /* ######### START SECOND UPDATE STEP ################# */
1367 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1370 if (ETC_ANDERSEN(ir
->etc
)) /* keep this outside of update_tcouple because of the extra info required to pass */
1372 gmx_bool bIfRandomize
;
1373 bIfRandomize
= update_randomize_velocities(ir
, step
, cr
, mdatoms
, state
, upd
, constr
);
1374 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1375 if (constr
&& bIfRandomize
)
1377 update_constraints(fplog
, step
, NULL
, ir
, mdatoms
,
1378 state
, fr
->bMolPBC
, graph
, f
,
1379 &top
->idef
, tmp_vir
,
1380 cr
, nrnb
, wcycle
, upd
, constr
,
1384 /* Box is changed in update() when we do pressure coupling,
1385 * but we should still use the old box for energy corrections and when
1386 * writing it to the energy file, so it matches the trajectory files for
1387 * the same timestep above. Make a copy in a separate array.
1389 copy_mat(state
->box
, lastbox
);
1393 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
)
1395 wallcycle_start(wcycle
, ewcUPDATE
);
1396 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1399 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ3
);
1400 /* We can only do Berendsen coupling after we have summed
1401 * the kinetic energy or virial. Since the happens
1402 * in global_state after update, we should only do it at
1403 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1408 update_tcouple(step
, ir
, state
, ekind
, &MassQ
, mdatoms
);
1409 update_pcouple(fplog
, step
, ir
, state
, pcoupl_mu
, M
, bInitStep
);
1414 /* velocity half-step update */
1415 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1416 ekind
, M
, upd
, etrtVELOCITY2
,
1420 /* Above, initialize just copies ekinh into ekin,
1421 * it doesn't copy position (for VV),
1422 * and entire integrator for MD.
1425 if (ir
->eI
== eiVVAK
)
1427 /* We probably only need md->homenr, not state->natoms */
1428 if (state
->natoms
> cbuf_nalloc
)
1430 cbuf_nalloc
= state
->natoms
;
1431 srenew(cbuf
, cbuf_nalloc
);
1433 copy_rvecn(state
->x
, cbuf
, 0, state
->natoms
);
1436 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1437 ekind
, M
, upd
, etrtPOSITION
, cr
, constr
);
1438 wallcycle_stop(wcycle
, ewcUPDATE
);
1440 update_constraints(fplog
, step
, &dvdl_constr
, ir
, mdatoms
, state
,
1441 fr
->bMolPBC
, graph
, f
,
1442 &top
->idef
, shake_vir
,
1443 cr
, nrnb
, wcycle
, upd
, constr
,
1446 if (ir
->eI
== eiVVAK
)
1448 /* erase F_EKIN and F_TEMP here? */
1449 /* just compute the kinetic energy at the half step to perform a trotter step */
1450 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1451 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1452 constr
, &nullSignaller
, lastbox
,
1453 NULL
, &bSumEkinhOld
,
1454 (bGStat
? CGLO_GSTAT
: 0) | CGLO_TEMPERATURE
1456 wallcycle_start(wcycle
, ewcUPDATE
);
1457 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ4
);
1458 /* now we know the scaling, we can compute the positions again again */
1459 copy_rvecn(cbuf
, state
->x
, 0, state
->natoms
);
1461 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1462 ekind
, M
, upd
, etrtPOSITION
, cr
, constr
);
1463 wallcycle_stop(wcycle
, ewcUPDATE
);
1465 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1466 /* are the small terms in the shake_vir here due
1467 * to numerical errors, or are they important
1468 * physically? I'm thinking they are just errors, but not completely sure.
1469 * For now, will call without actually constraining, constr=NULL*/
1470 update_constraints(fplog
, step
, NULL
, ir
, mdatoms
,
1471 state
, fr
->bMolPBC
, graph
, f
,
1472 &top
->idef
, tmp_vir
,
1473 cr
, nrnb
, wcycle
, upd
, NULL
,
1478 /* this factor or 2 correction is necessary
1479 because half of the constraint force is removed
1480 in the vv step, so we have to double it. See
1481 the Redmine issue #1255. It is not yet clear
1482 if the factor of 2 is exact, or just a very
1483 good approximation, and this will be
1484 investigated. The next step is to see if this
1485 can be done adding a dhdl contribution from the
1486 rattle step, but this is somewhat more
1487 complicated with the current code. Will be
1488 investigated, hopefully for 4.6.3. However,
1489 this current solution is much better than
1490 having it completely wrong.
1492 enerd
->term
[F_DVDL_CONSTR
] += 2*dvdl_constr
;
1496 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
1501 /* Need to unshift here */
1502 unshift_self(graph
, state
->box
, state
->x
);
1507 wallcycle_start(wcycle
, ewcVSITECONSTR
);
1510 shift_self(graph
, state
->box
, state
->x
);
1512 construct_vsites(vsite
, state
->x
, ir
->delta_t
, state
->v
,
1513 top
->idef
.iparams
, top
->idef
.il
,
1514 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
1518 unshift_self(graph
, state
->box
, state
->x
);
1520 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
1523 /* ############## IF NOT VV, Calculate globals HERE ############ */
1524 /* With Leap-Frog we can skip compute_globals at
1525 * non-communication steps, but we need to calculate
1526 * the kinetic energy one step before communication.
1529 // Organize to do inter-simulation signalling on steps if
1530 // and when algorithms require it.
1531 bool doInterSimSignal
= (!bFirstStep
&& bDoReplEx
) || bUsingEnsembleRestraints
;
1533 if (bGStat
|| (!EI_VV(ir
->eI
) && do_per_step(step
+1, nstglobalcomm
)) || doInterSimSignal
)
1535 // Since we're already communicating at this step, we
1536 // can propagate intra-simulation signals. Note that
1537 // check_nstglobalcomm has the responsibility for
1538 // choosing the value of nstglobalcomm that is one way
1539 // bGStat becomes true, so we can't get into a
1540 // situation where e.g. checkpointing can't be
1542 bool doIntraSimSignal
= true;
1543 SimulationSignaller
signaller(&signals
, cr
, doInterSimSignal
, doIntraSimSignal
);
1545 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1546 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1549 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1550 (bGStat
? CGLO_GSTAT
: 0)
1551 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_ENERGY
: 0)
1552 | (!EI_VV(ir
->eI
) && bStopCM
? CGLO_STOPCM
: 0)
1553 | (!EI_VV(ir
->eI
) ? CGLO_TEMPERATURE
: 0)
1554 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_PRESSURE
: 0)
1556 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
1558 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
1559 top_global
, top
, state
,
1560 &shouldCheckNumberOfBondedInteractions
);
1564 /* ############# END CALC EKIN AND PRESSURE ################# */
1566 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1567 the virial that should probably be addressed eventually. state->veta has better properies,
1568 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1569 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1571 if (ir
->efep
!= efepNO
&& (!EI_VV(ir
->eI
) || bRerunMD
))
1573 /* Sum up the foreign energy and dhdl terms for md and sd.
1574 Currently done every step so that dhdl is correct in the .edr */
1575 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1577 update_box(fplog
, step
, ir
, mdatoms
, state
, f
,
1578 pcoupl_mu
, nrnb
, upd
);
1580 /* ################# END UPDATE STEP 2 ################# */
1581 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1583 /* The coordinates (x) were unshifted in update */
1586 /* We will not sum ekinh_old,
1587 * so signal that we still have to do it.
1589 bSumEkinhOld
= TRUE
;
1592 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1594 /* use the directly determined last velocity, not actually the averaged half steps */
1595 if (bTrotter
&& ir
->eI
== eiVV
)
1597 enerd
->term
[F_EKIN
] = last_ekin
;
1599 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1603 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + saved_conserved_quantity
;
1607 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + compute_conserved_from_auxiliary(ir
, state
, &MassQ
);
1609 /* ######### END PREPARING EDR OUTPUT ########### */
1614 if (fplog
&& do_log
&& bDoExpanded
)
1616 /* only needed if doing expanded ensemble */
1617 PrintFreeEnergyInfoToFile(fplog
, ir
->fepvals
, ir
->expandedvals
, ir
->bSimTemp
? ir
->simtempvals
: NULL
,
1618 &state_global
->dfhist
, state
->fep_state
, ir
->nstlog
, step
);
1622 upd_mdebin(mdebin
, bDoDHDL
, bCalcEnerStep
,
1623 t
, mdatoms
->tmass
, enerd
, state
,
1624 ir
->fepvals
, ir
->expandedvals
, lastbox
,
1625 shake_vir
, force_vir
, total_vir
, pres
,
1626 ekind
, mu_tot
, constr
);
1630 upd_mdebin_step(mdebin
);
1633 gmx_bool do_dr
= do_per_step(step
, ir
->nstdisreout
);
1634 gmx_bool do_or
= do_per_step(step
, ir
->nstorireout
);
1636 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
, do_log
? fplog
: NULL
,
1638 eprNORMAL
, mdebin
, fcd
, groups
, &(ir
->opts
));
1642 pull_print_output(ir
->pull_work
, step
, t
);
1645 if (do_per_step(step
, ir
->nstlog
))
1647 if (fflush(fplog
) != 0)
1649 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
1655 /* Have to do this part _after_ outputting the logfile and the edr file */
1656 /* Gets written into the state at the beginning of next loop*/
1657 state
->fep_state
= lamnew
;
1659 /* Print the remaining wall clock time for the run */
1660 if (MULTIMASTER(cr
) &&
1661 (do_verbose
|| gmx_got_usr_signal()) &&
1666 fprintf(stderr
, "\n");
1668 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
1671 /* Ion/water position swapping.
1672 * Not done in last step since trajectory writing happens before this call
1673 * in the MD loop and exchanges would be lost anyway. */
1674 bNeedRepartition
= FALSE
;
1675 if ((ir
->eSwapCoords
!= eswapNO
) && (step
> 0) && !bLastStep
&&
1676 do_per_step(step
, ir
->swap
->nstswap
))
1678 bNeedRepartition
= do_swapcoords(cr
, step
, t
, ir
, wcycle
,
1679 bRerunMD
? rerun_fr
.x
: state
->x
,
1680 bRerunMD
? rerun_fr
.box
: state
->box
,
1681 top_global
, MASTER(cr
) && bVerbose
, bRerunMD
);
1683 if (bNeedRepartition
&& DOMAINDECOMP(cr
))
1685 dd_collect_state(cr
->dd
, state
, state_global
);
1689 /* Replica exchange */
1693 bExchanged
= replica_exchange(fplog
, cr
, repl_ex
,
1694 state_global
, enerd
,
1698 if ( (bExchanged
|| bNeedRepartition
) && DOMAINDECOMP(cr
) )
1700 dd_partition_system(fplog
, step
, cr
, TRUE
, 1,
1701 state_global
, top_global
, ir
,
1702 state
, &f
, mdatoms
, top
, fr
,
1704 nrnb
, wcycle
, FALSE
);
1705 shouldCheckNumberOfBondedInteractions
= true;
1706 update_realloc(upd
, state
->nalloc
);
1711 startingFromCheckpoint
= FALSE
;
1713 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1714 /* With all integrators, except VV, we need to retain the pressure
1715 * at the current step for coupling at the next step.
1717 if ((state
->flags
& (1<<estPRES_PREV
)) &&
1719 (ir
->nstpcouple
> 0 && step
% ir
->nstpcouple
== 0)))
1721 /* Store the pressure in t_state for pressure coupling
1722 * at the next MD step.
1724 copy_mat(pres
, state
->pres_prev
);
1727 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1729 if ( (membed
!= NULL
) && (!bLastStep
) )
1731 rescale_membed(step_rel
, membed
, state_global
->x
);
1738 /* read next frame from input trajectory */
1739 bLastStep
= !read_next_frame(oenv
, status
, &rerun_fr
);
1744 rerun_parallel_comm(cr
, &rerun_fr
, &bLastStep
);
1748 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
1749 if (DOMAINDECOMP(cr
) && wcycle
)
1751 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
1754 if (!bRerunMD
|| !rerun_fr
.bStep
)
1756 /* increase the MD step number */
1761 /* TODO make a counter-reset module */
1762 /* If it is time to reset counters, set a flag that remains
1763 true until counters actually get reset */
1764 if (step_rel
== wcycle_get_reset_counters(wcycle
) ||
1765 signals
[eglsRESETCOUNTERS
].set
!= 0)
1767 if (pme_loadbal_is_active(pme_loadbal
))
1769 /* Do not permit counter reset while PME load
1770 * balancing is active. The only purpose for resetting
1771 * counters is to measure reliable performance data,
1772 * and that can't be done before balancing
1775 * TODO consider fixing this by delaying the reset
1776 * until after load balancing completes,
1777 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1778 gmx_fatal(FARGS
, "PME tuning was still active when attempting to "
1779 "reset mdrun counters at step %" GMX_PRId64
". Try "
1780 "resetting counters later in the run, e.g. with gmx "
1781 "mdrun -resetstep.", step
);
1783 reset_all_counters(fplog
, cr
, step
, &step_rel
, ir
, wcycle
, nrnb
, walltime_accounting
,
1784 use_GPU(fr
->nbv
) ? fr
->nbv
: NULL
);
1785 wcycle_set_reset_counters(wcycle
, -1);
1786 if (!(cr
->duty
& DUTY_PME
))
1788 /* Tell our PME node to reset its counters */
1789 gmx_pme_send_resetcounters(cr
, step
);
1791 /* Correct max_hours for the elapsed time */
1792 max_hours
-= elapsed_time
/(60.0*60.0);
1793 /* If mdrun -maxh -resethway was active, it can only trigger once */
1794 bResetCountersHalfMaxH
= FALSE
; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1795 /* Reset can only happen once, so clear the triggering flag. */
1796 signals
[eglsRESETCOUNTERS
].set
= 0;
1799 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1800 IMD_prep_energies_send_positions(ir
->bIMD
&& MASTER(cr
), bIMDstep
, ir
->imd
, enerd
, step
, bCalcEner
, wcycle
);
1803 /* End of main MD loop */
1805 /* Closing TNG files can include compressing data. Therefore it is good to do that
1806 * before stopping the time measurements. */
1807 mdoutf_tng_close(outf
);
1809 /* Stop measuring walltime */
1810 walltime_accounting_end(walltime_accounting
);
1812 if (bRerunMD
&& MASTER(cr
))
1817 if (!(cr
->duty
& DUTY_PME
))
1819 /* Tell the PME only node to finish */
1820 gmx_pme_send_finish(cr
);
1825 if (ir
->nstcalcenergy
> 0 && !bRerunMD
)
1827 print_ebin(mdoutf_get_fp_ene(outf
), FALSE
, FALSE
, FALSE
, fplog
, step
, t
,
1828 eprAVER
, mdebin
, fcd
, groups
, &(ir
->opts
));
1836 pme_loadbal_done(pme_loadbal
, cr
, fplog
, use_GPU(fr
->nbv
));
1839 done_shellfc(fplog
, shellfc
, step_rel
);
1841 if (repl_ex_nst
> 0 && MASTER(cr
))
1843 print_replica_exchange_statistics(fplog
, repl_ex
);
1846 // Clean up swapcoords
1847 if (ir
->eSwapCoords
!= eswapNO
)
1849 finish_swapcoords(ir
->swap
);
1852 if (membed
!= nullptr)
1854 free_membed(membed
);
1857 /* IMD cleanup, if bIMD is TRUE. */
1858 IMD_finalize(ir
->bIMD
, ir
->imd
);
1860 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);