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,2017,2018, 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.
39 * \brief Implements the integrator for normal molecular dynamics simulations
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
56 #include "gromacs/awh/awh.h"
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/compat/make_unique.h"
59 #include "gromacs/domdec/collect.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_network.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/essentialdynamics/edsam.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme-load-balancing.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/gpu_utils/gpu_utils.h"
70 #include "gromacs/imd/imd.h"
71 #include "gromacs/listed-forces/manage-threading.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/math/vectypes.h"
76 #include "gromacs/mdlib/checkpointhandler.h"
77 #include "gromacs/mdlib/compute_io.h"
78 #include "gromacs/mdlib/constr.h"
79 #include "gromacs/mdlib/ebin.h"
80 #include "gromacs/mdlib/expanded.h"
81 #include "gromacs/mdlib/force.h"
82 #include "gromacs/mdlib/force_flags.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdebin.h"
87 #include "gromacs/mdlib/mdoutf.h"
88 #include "gromacs/mdlib/mdrun.h"
89 #include "gromacs/mdlib/mdsetup.h"
90 #include "gromacs/mdlib/membed.h"
91 #include "gromacs/mdlib/nb_verlet.h"
92 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
93 #include "gromacs/mdlib/ns.h"
94 #include "gromacs/mdlib/resethandler.h"
95 #include "gromacs/mdlib/shellfc.h"
96 #include "gromacs/mdlib/sighandler.h"
97 #include "gromacs/mdlib/sim_util.h"
98 #include "gromacs/mdlib/simulationsignal.h"
99 #include "gromacs/mdlib/tgroup.h"
100 #include "gromacs/mdlib/trajectory_writing.h"
101 #include "gromacs/mdlib/update.h"
102 #include "gromacs/mdlib/vcm.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdtypes/awh-history.h"
105 #include "gromacs/mdtypes/awh-params.h"
106 #include "gromacs/mdtypes/commrec.h"
107 #include "gromacs/mdtypes/df_history.h"
108 #include "gromacs/mdtypes/energyhistory.h"
109 #include "gromacs/mdtypes/fcdata.h"
110 #include "gromacs/mdtypes/forcerec.h"
111 #include "gromacs/mdtypes/group.h"
112 #include "gromacs/mdtypes/inputrec.h"
113 #include "gromacs/mdtypes/interaction_const.h"
114 #include "gromacs/mdtypes/md_enums.h"
115 #include "gromacs/mdtypes/mdatom.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/state.h"
118 #include "gromacs/pbcutil/mshift.h"
119 #include "gromacs/pbcutil/pbc.h"
120 #include "gromacs/pulling/pull.h"
121 #include "gromacs/swap/swapcoords.h"
122 #include "gromacs/timing/wallcycle.h"
123 #include "gromacs/timing/walltime_accounting.h"
124 #include "gromacs/topology/atoms.h"
125 #include "gromacs/topology/idef.h"
126 #include "gromacs/topology/mtop_util.h"
127 #include "gromacs/topology/topology.h"
128 #include "gromacs/trajectory/trajectoryframe.h"
129 #include "gromacs/utility/basedefinitions.h"
130 #include "gromacs/utility/cstringutil.h"
131 #include "gromacs/utility/fatalerror.h"
132 #include "gromacs/utility/logger.h"
133 #include "gromacs/utility/real.h"
134 #include "gromacs/utility/smalloc.h"
136 #include "integrator.h"
137 #include "replicaexchange.h"
140 #include "corewrap.h"
143 using gmx::SimulationSignaller
;
145 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
147 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
148 * \param[in,out] globalState The global state container
149 * \param[in] constructVsites When true, vsite coordinates are constructed
150 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
151 * \param[in] idef Topology parameters, used for constructing vsites
152 * \param[in] timeStep Time step, used for constructing vsites
153 * \param[in] forceRec Force record, used for constructing vsites
154 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
155 * \param[in,out] warnWhenNoV When true, issue a warning when no velocities are present in \p rerunFrame; is set to false when a warning was issued
157 static void prepareRerunState(const t_trxframe
&rerunFrame
,
158 t_state
*globalState
,
159 bool constructVsites
,
160 const gmx_vsite_t
*vsite
,
163 const t_forcerec
&forceRec
,
165 gmx_bool
*warnWhenNoV
)
167 for (int i
= 0; i
< globalState
->natoms
; i
++)
169 copy_rvec(rerunFrame
.x
[i
], globalState
->x
[i
]);
173 for (int i
= 0; i
< globalState
->natoms
; i
++)
175 copy_rvec(rerunFrame
.v
[i
], globalState
->v
[i
]);
180 for (int i
= 0; i
< globalState
->natoms
; i
++)
182 clear_rvec(globalState
->v
[i
]);
186 fprintf(stderr
, "\nWARNING: Some frames do not contain velocities.\n"
187 " Ekin, temperature and pressure are incorrect,\n"
188 " the virial will be incorrect when constraints are present.\n"
190 *warnWhenNoV
= FALSE
;
193 copy_mat(rerunFrame
.box
, globalState
->box
);
197 GMX_ASSERT(vsite
, "Need valid vsite for constructing vsites");
201 /* Following is necessary because the graph may get out of sync
202 * with the coordinates if we only have every N'th coordinate set
204 mk_mshift(nullptr, graph
, forceRec
.ePBC
, globalState
->box
, as_rvec_array(globalState
->x
.data()));
205 shift_self(graph
, globalState
->box
, as_rvec_array(globalState
->x
.data()));
207 construct_vsites(vsite
, as_rvec_array(globalState
->x
.data()), timeStep
, as_rvec_array(globalState
->v
.data()),
208 idef
.iparams
, idef
.il
,
209 forceRec
.ePBC
, forceRec
.bMolPBC
, nullptr, globalState
->box
);
212 unshift_self(graph
, globalState
->box
, as_rvec_array(globalState
->x
.data()));
217 void gmx::Integrator::do_md()
219 // TODO Historically, the EM and MD "integrators" used different
220 // names for the t_inputrec *parameter, but these must have the
221 // same name, now that it's a member of a struct. We use this ir
222 // alias to avoid a large ripple of nearly useless changes.
223 // t_inputrec is being replaced by IMdpOptionsProvider, so this
224 // will go away eventually.
225 t_inputrec
*ir
= inputrec
;
226 gmx_mdoutf
*outf
= nullptr;
227 int64_t step
, step_rel
;
228 double t
, t0
, lam0
[efptNR
];
229 gmx_bool bGStatEveryStep
, bGStat
, bCalcVir
, bCalcEnerStep
, bCalcEner
;
230 gmx_bool bNS
, bNStList
, bSimAnn
, bStopCM
,
231 bFirstStep
, bInitStep
, bLastStep
= FALSE
;
232 gmx_bool bDoDHDL
= FALSE
, bDoFEP
= FALSE
, bDoExpanded
= FALSE
;
233 gmx_bool do_ene
, do_log
, do_verbose
, bRerunWarnNoV
= TRUE
,
234 bForceUpdate
= FALSE
;
235 gmx_bool bMasterState
;
236 int force_flags
, cglo_flags
;
237 tensor force_vir
, shake_vir
, total_vir
, tmp_vir
, pres
;
242 matrix parrinellorahmanMu
, M
;
244 gmx_repl_ex_t repl_ex
= nullptr;
246 t_mdebin
*mdebin
= nullptr;
247 gmx_enerdata_t
*enerd
;
248 PaddedRVecVector f
{};
249 gmx_global_stat_t gstat
;
250 gmx_update_t
*upd
= nullptr;
251 t_graph
*graph
= nullptr;
252 gmx_groups_t
*groups
;
253 gmx_ekindata_t
*ekind
;
254 gmx_shellfc_t
*shellfc
;
255 gmx_bool bSumEkinhOld
, bDoReplEx
, bExchanged
, bNeedRepartition
;
256 gmx_bool bTemp
, bPres
, bTrotter
;
258 rvec
*cbuf
= nullptr;
265 real saved_conserved_quantity
= 0;
269 char sbuf
[STEPSTRSIZE
], sbuf2
[STEPSTRSIZE
];
270 int handled_stop_condition
= gmx_stop_cond_none
; /* compare to get_stop_condition*/
273 /* PME load balancing data for GPU kernels */
274 pme_load_balancing_t
*pme_loadbal
= nullptr;
275 gmx_bool bPMETune
= FALSE
;
276 gmx_bool bPMETunePrinting
= FALSE
;
279 gmx_bool bIMDstep
= FALSE
;
282 /* Temporary addition for FAHCORE checkpointing */
285 /* Domain decomposition could incorrectly miss a bonded
286 interaction, but checking for that requires a global
287 communication stage, which does not otherwise happen in DD
288 code. So we do that alongside the first global energy reduction
289 after a new DD is made. These variables handle whether the
290 check happens, and the result it returns. */
291 bool shouldCheckNumberOfBondedInteractions
= false;
292 int totalNumberOfBondedInteractions
= -1;
294 SimulationSignals signals
;
295 // Most global communnication stages don't propagate mdrun
296 // signals, and will use this object to achieve that.
297 SimulationSignaller
nullSignaller(nullptr, nullptr, nullptr, false, false);
299 if (!mdrunOptions
.writeConfout
)
301 // This is on by default, and the main known use case for
302 // turning it off is for convenience in benchmarking, which is
303 // something that should not show up in the general user
305 GMX_LOG(mdlog
.info
).asParagraph().
306 appendText("The -noconfout functionality is deprecated, and may be removed in a future version.");
309 /* md-vv uses averaged full step velocities for T-control
310 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
311 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
312 bTrotter
= (EI_VV(ir
->eI
) && (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
) || inputrecNvtTrotter(ir
)));
314 const gmx_bool bRerunMD
= mdrunOptions
.rerun
;
315 int nstglobalcomm
= mdrunOptions
.globalCommunicationInterval
;
319 /* Since we don't know if the frames read are related in any way,
320 * rebuild the neighborlist at every step.
323 ir
->nstcalcenergy
= 1;
327 nstglobalcomm
= check_nstglobalcomm(mdlog
, nstglobalcomm
, ir
);
328 bGStatEveryStep
= (nstglobalcomm
== 1);
332 ir
->nstxout_compressed
= 0;
334 groups
= &top_global
->groups
;
336 std::unique_ptr
<EssentialDynamics
> ed
= nullptr;
337 if (opt2bSet("-ei", nfile
, fnm
) || observablesHistory
->edsamHistory
!= nullptr)
339 /* Initialize essential dynamics sampling */
340 ed
= init_edsam(opt2fn_null("-ei", nfile
, fnm
), opt2fn("-eo", nfile
, fnm
),
343 state_global
, observablesHistory
,
344 oenv
, mdrunOptions
.continuationOptions
.appendFiles
);
348 init_md(fplog
, cr
, outputProvider
, ir
, oenv
, mdrunOptions
,
349 &t
, &t0
, state_global
, lam0
,
350 nrnb
, top_global
, &upd
, deform
,
351 nfile
, fnm
, &outf
, &mdebin
,
352 force_vir
, shake_vir
, total_vir
, pres
, mu_tot
, &bSimAnn
, &vcm
, wcycle
);
354 /* Energy terms and groups */
356 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
359 /* Kinetic energy data */
361 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind
);
362 /* Copy the cos acceleration to the groups struct */
363 ekind
->cosacc
.cos_accel
= ir
->cos_accel
;
365 gstat
= global_stat_init(ir
);
367 /* Check for polarizable models and flexible constraints */
368 shellfc
= init_shell_flexcon(fplog
,
369 top_global
, constr
? constr
->numFlexibleConstraints() : 0,
370 ir
->nstcalcenergy
, DOMAINDECOMP(cr
));
373 double io
= compute_io(ir
, top_global
->natoms
, groups
, mdebin
->ebin
->nener
, 1);
374 if ((io
> 2000) && MASTER(cr
))
377 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
382 /* Set up interactive MD (IMD) */
383 init_IMD(ir
, cr
, ms
, top_global
, fplog
, ir
->nstcalcenergy
,
384 MASTER(cr
) ? as_rvec_array(state_global
->x
.data()) : nullptr,
385 nfile
, fnm
, oenv
, mdrunOptions
);
387 // Local state only becomes valid now.
388 std::unique_ptr
<t_state
> stateInstance
;
391 if (DOMAINDECOMP(cr
))
393 top
= dd_init_local_top(top_global
);
395 stateInstance
= compat::make_unique
<t_state
>();
396 state
= stateInstance
.get();
397 dd_init_local_state(cr
->dd
, state_global
, state
);
399 /* Distribute the charge groups over the nodes from the master node */
400 dd_partition_system(fplog
, mdlog
, ir
->init_step
, cr
, TRUE
, 1,
401 state_global
, top_global
, ir
,
402 state
, &f
, mdAtoms
, top
, fr
,
404 nrnb
, nullptr, FALSE
);
405 shouldCheckNumberOfBondedInteractions
= true;
406 update_realloc(upd
, state
->natoms
);
410 state_change_natoms(state_global
, state_global
->natoms
);
411 /* We need to allocate one element extra, since we might use
412 * (unaligned) 4-wide SIMD loads to access rvec entries.
414 f
.resize(gmx::paddedRVecVectorSize(state_global
->natoms
));
415 /* Copy the pointer to the global state */
416 state
= state_global
;
419 mdAlgorithmsSetupAtomData(cr
, ir
, top_global
, top
, fr
,
420 &graph
, mdAtoms
, constr
, vsite
, shellfc
);
422 update_realloc(upd
, state
->natoms
);
425 auto mdatoms
= mdAtoms
->mdatoms();
427 // NOTE: The global state is no longer used at this point.
428 // But state_global is still used as temporary storage space for writing
429 // the global state to file and potentially for replica exchange.
430 // (Global topology should persist.)
432 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
434 const ContinuationOptions
&continuationOptions
= mdrunOptions
.continuationOptions
;
435 bool startingFromCheckpoint
= continuationOptions
.startedFromCheckpoint
;
439 init_expanded_ensemble(startingFromCheckpoint
, ir
, state
->dfhist
);
444 if (startingFromCheckpoint
)
446 /* Update mdebin with energy history if appending to output files */
447 if (continuationOptions
.appendFiles
)
449 restore_energyhistory_from_state(mdebin
, observablesHistory
->energyHistory
.get());
451 else if (observablesHistory
->energyHistory
!= nullptr)
453 /* We might have read an energy history from checkpoint.
454 * As we are not appending, we want to restart the statistics.
455 * Free the allocated memory and reset the counts.
457 observablesHistory
->energyHistory
= {};
460 if (observablesHistory
->energyHistory
== nullptr)
462 observablesHistory
->energyHistory
= compat::make_unique
<energyhistory_t
>();
464 /* Set the initial energy history in state by updating once */
465 update_energyhistory(observablesHistory
->energyHistory
.get(), mdebin
);
468 // TODO: Remove this by converting AWH into a ForceProvider
469 auto awh
= prepareAwhModule(fplog
, *ir
, state_global
, cr
, ms
, startingFromCheckpoint
,
471 opt2fn("-awh", nfile
, fnm
), ir
->pull_work
);
473 const bool useReplicaExchange
= (replExParams
.exchangeInterval
> 0);
474 if (useReplicaExchange
&& MASTER(cr
))
476 repl_ex
= init_replica_exchange(fplog
, ms
, top_global
->natoms
, ir
,
479 /* PME tuning is only supported in the Verlet scheme, with PME for
480 * Coulomb. It is not supported with only LJ PME, or for
482 bPMETune
= (mdrunOptions
.tunePme
&& EEL_PME(fr
->ic
->eeltype
) && !bRerunMD
&&
483 !mdrunOptions
.reproducible
&& ir
->cutoff_scheme
!= ecutsGROUP
);
486 pme_loadbal_init(&pme_loadbal
, cr
, mdlog
, *ir
, state
->box
,
487 *fr
->ic
, *fr
->nbv
->listParams
, fr
->pmedata
, use_GPU(fr
->nbv
),
491 if (!ir
->bContinuation
&& !bRerunMD
)
493 if (state
->flags
& (1 << estV
))
495 /* Set the velocities of vsites, shells and frozen atoms to zero */
496 for (i
= 0; i
< mdatoms
->homenr
; i
++)
498 if (mdatoms
->ptype
[i
] == eptVSite
||
499 mdatoms
->ptype
[i
] == eptShell
)
501 clear_rvec(state
->v
[i
]);
503 else if (mdatoms
->cFREEZE
)
505 for (m
= 0; m
< DIM
; m
++)
507 if (ir
->opts
.nFreeze
[mdatoms
->cFREEZE
[i
]][m
])
518 /* Constrain the initial coordinates and velocities */
519 do_constrain_first(fplog
, constr
, ir
, mdatoms
, state
);
523 /* Construct the virtual sites for the initial configuration */
524 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, nullptr,
525 top
->idef
.iparams
, top
->idef
.il
,
526 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
530 if (ir
->efep
!= efepNO
)
532 /* Set free energy calculation frequency as the greatest common
533 * denominator of nstdhdl and repl_ex_nst. */
534 nstfep
= ir
->fepvals
->nstdhdl
;
537 nstfep
= gmx_greatest_common_divisor(ir
->expandedvals
->nstexpanded
, nstfep
);
539 if (useReplicaExchange
)
541 nstfep
= gmx_greatest_common_divisor(replExParams
.exchangeInterval
, nstfep
);
545 /* Be REALLY careful about what flags you set here. You CANNOT assume
546 * this is the first step, since we might be restarting from a checkpoint,
547 * and in that case we should not do any modifications to the state.
549 bStopCM
= (ir
->comm_mode
!= ecmNO
&& !ir
->bContinuation
);
551 if (continuationOptions
.haveReadEkin
)
553 restore_ekinstate_from_state(cr
, ekind
, &state_global
->ekinstate
);
556 cglo_flags
= (CGLO_INITIALIZATION
| CGLO_TEMPERATURE
| CGLO_GSTAT
557 | (EI_VV(ir
->eI
) ? CGLO_PRESSURE
: 0)
558 | (EI_VV(ir
->eI
) ? CGLO_CONSTRAINT
: 0)
559 | (continuationOptions
.haveReadEkin
? CGLO_READEKIN
: 0));
561 bSumEkinhOld
= FALSE
;
562 /* To minimize communication, compute_globals computes the COM velocity
563 * and the kinetic energy for the velocities without COM motion removed.
564 * Thus to get the kinetic energy without the COM contribution, we need
565 * to call compute_globals twice.
567 for (int cgloIteration
= 0; cgloIteration
< (bStopCM
? 2 : 1); cgloIteration
++)
569 int cglo_flags_iteration
= cglo_flags
;
570 if (bStopCM
&& cgloIteration
== 0)
572 cglo_flags_iteration
|= CGLO_STOPCM
;
573 cglo_flags_iteration
&= ~CGLO_TEMPERATURE
;
575 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
576 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
577 constr
, &nullSignaller
, state
->box
,
578 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags_iteration
579 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
581 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
582 top_global
, top
, state
,
583 &shouldCheckNumberOfBondedInteractions
);
584 if (ir
->eI
== eiVVAK
)
586 /* a second call to get the half step temperature initialized as well */
587 /* we do the same call as above, but turn the pressure off -- internally to
588 compute_globals, this is recognized as a velocity verlet half-step
589 kinetic energy calculation. This minimized excess variables, but
590 perhaps loses some logic?*/
592 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
593 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
594 constr
, &nullSignaller
, state
->box
,
595 nullptr, &bSumEkinhOld
,
596 cglo_flags
& ~CGLO_PRESSURE
);
599 /* Calculate the initial half step temperature, and save the ekinh_old */
600 if (!continuationOptions
.startedFromCheckpoint
)
602 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
604 copy_mat(ekind
->tcstat
[i
].ekinh
, ekind
->tcstat
[i
].ekinh_old
);
608 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
609 temperature control */
610 trotter_seq
= init_npt_vars(ir
, state
, &MassQ
, bTrotter
);
614 if (!ir
->bContinuation
)
616 if (constr
&& ir
->eConstrAlg
== econtLINCS
)
619 "RMS relative constraint deviation after constraining: %.2e\n",
622 if (EI_STATE_VELOCITY(ir
->eI
))
624 real temp
= enerd
->term
[F_TEMP
];
627 /* Result of Ekin averaged over velocities of -half
628 * and +half step, while we only have -half step here.
632 fprintf(fplog
, "Initial temperature: %g K\n", temp
);
638 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
639 " input trajectory '%s'\n\n",
640 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
641 if (mdrunOptions
.verbose
)
643 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
644 "run input file,\nwhich may not correspond to the time "
645 "needed to process input trajectory.\n\n");
651 fprintf(stderr
, "starting mdrun '%s'\n",
652 *(top_global
->name
));
655 sprintf(tbuf
, "%8.1f", (ir
->init_step
+ir
->nsteps
)*ir
->delta_t
);
659 sprintf(tbuf
, "%s", "infinite");
661 if (ir
->init_step
> 0)
663 fprintf(stderr
, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
664 gmx_step_str(ir
->init_step
+ir
->nsteps
, sbuf
), tbuf
,
665 gmx_step_str(ir
->init_step
, sbuf2
),
666 ir
->init_step
*ir
->delta_t
);
670 fprintf(stderr
, "%s steps, %s ps.\n",
671 gmx_step_str(ir
->nsteps
, sbuf
), tbuf
);
674 fprintf(fplog
, "\n");
677 walltime_accounting_start_time(walltime_accounting
);
678 wallcycle_start(wcycle
, ewcRUN
);
679 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
681 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
683 chkpt_ret
= fcCheckPointParallel( cr
->nodeid
,
687 gmx_fatal( 3, __FILE__
, __LINE__
, "Checkpoint error on step %d\n", 0 );
691 /***********************************************************
695 ************************************************************/
697 /* if rerunMD then read coordinates and velocities from input trajectory */
700 if (getenv("GMX_FORCE_UPDATE"))
708 bLastStep
= !read_first_frame(oenv
, &status
,
709 opt2fn("-rerun", nfile
, fnm
),
710 &rerun_fr
, TRX_NEED_X
| TRX_READ_V
);
711 if (rerun_fr
.natoms
!= top_global
->natoms
)
714 "Number of atoms in trajectory (%d) does not match the "
715 "run input file (%d)\n",
716 rerun_fr
.natoms
, top_global
->natoms
);
718 if (ir
->ePBC
!= epbcNONE
)
722 gmx_fatal(FARGS
, "Rerun trajectory frame step %" PRId64
" time %f does not contain a box, while pbc is used", rerun_fr
.step
, rerun_fr
.time
);
724 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < gmx::square(fr
->rlist
))
726 gmx_fatal(FARGS
, "Rerun trajectory frame step %" PRId64
" time %f has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
733 rerun_parallel_comm(cr
, &rerun_fr
, &bLastStep
);
736 if (ir
->ePBC
!= epbcNONE
)
738 /* Set the shift vectors.
739 * Necessary here when have a static box different from the tpr box.
741 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
746 /* Skip the first Nose-Hoover integration when we get the state from tpx */
747 bInitStep
= !startingFromCheckpoint
|| EI_VV(ir
->eI
);
748 bSumEkinhOld
= FALSE
;
750 bNeedRepartition
= FALSE
;
752 bool simulationsShareState
= false;
753 int nstSignalComm
= nstglobalcomm
;
755 // TODO This implementation of ensemble orientation restraints is nasty because
756 // a user can't just do multi-sim with single-sim orientation restraints.
757 bool usingEnsembleRestraints
= (fcd
->disres
.nsystems
> 1) || ((ms
!= nullptr) && (fcd
->orires
.nr
!= 0));
758 bool awhUsesMultiSim
= (ir
->bDoAwh
&& ir
->awhParams
->shareBiasMultisim
&& (ms
!= nullptr));
760 // Replica exchange, ensemble restraints and AWH need all
761 // simulations to remain synchronized, so they need
762 // checkpoints and stop conditions to act on the same step, so
763 // the propagation of such signals must take place between
764 // simulations, not just within simulations.
765 // TODO: Make algorithm initializers set these flags.
766 simulationsShareState
= useReplicaExchange
|| usingEnsembleRestraints
|| awhUsesMultiSim
;
768 signals
[eglsSTOPCOND
] = SimulationSignal(!simulationsShareState
);
770 if (simulationsShareState
)
772 // Inter-simulation signal communication does not need to happen
773 // often, so we use a minimum of 200 steps to reduce overhead.
774 const int c_minimumInterSimulationSignallingInterval
= 200;
775 nstSignalComm
= ((c_minimumInterSimulationSignallingInterval
+ nstglobalcomm
- 1)/nstglobalcomm
)*nstglobalcomm
;
779 std::unique_ptr
<CheckpointHandler
> checkpointHandler
= nullptr;
783 checkpointHandler
= compat::make_unique
<CheckpointHandler
>(
784 compat::make_not_null
<SimulationSignal
*>(&signals
[eglsCHKPT
]),
785 simulationsShareState
, ir
->nstlist
== 0, MASTER(cr
),
786 mdrunOptions
.writeConfout
, mdrunOptions
.checkpointOptions
.period
);
789 const bool resetCountersIsLocal
= true;
790 auto resetHandler
= compat::make_unique
<ResetHandler
>(
791 compat::make_not_null
<SimulationSignal
*>(&signals
[eglsRESETCOUNTERS
]), !resetCountersIsLocal
,
792 ir
->nsteps
, MASTER(cr
), mdrunOptions
.timingOptions
.resetHalfway
,
793 mdrunOptions
.maximumHoursToRun
, mdlog
, wcycle
, walltime_accounting
);
795 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion
= (DOMAINDECOMP(cr
) ? DdOpenBalanceRegionBeforeForceComputation::yes
: DdOpenBalanceRegionBeforeForceComputation::no
);
796 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion
= (DOMAINDECOMP(cr
) ? DdCloseBalanceRegionAfterForceComputation::yes
: DdCloseBalanceRegionAfterForceComputation::no
);
798 step
= ir
->init_step
;
801 // TODO extract this to new multi-simulation module
802 if (MASTER(cr
) && isMultiSim(ms
) && !useReplicaExchange
)
804 if (!multisim_int_all_are_equal(ms
, ir
->nsteps
))
806 GMX_LOG(mdlog
.warning
).appendText(
807 "Note: The number of steps is not consistent across multi simulations,\n"
808 "but we are proceeding anyway!");
810 if (!multisim_int_all_are_equal(ms
, ir
->init_step
))
812 GMX_LOG(mdlog
.warning
).appendText(
813 "Note: The initial step is not consistent across multi simulations,\n"
814 "but we are proceeding anyway!");
818 /* and stop now if we should */
819 bLastStep
= (bLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
823 /* Determine if this is a neighbor search step */
824 bNStList
= (ir
->nstlist
> 0 && step
% ir
->nstlist
== 0);
826 if (bPMETune
&& bNStList
)
828 /* PME grid + cut-off optimization with GPUs or PME nodes */
829 pme_loadbal_do(pme_loadbal
, cr
,
830 (mdrunOptions
.verbose
&& MASTER(cr
)) ? stderr
: nullptr,
838 wallcycle_start(wcycle
, ewcSTEP
);
844 step
= rerun_fr
.step
;
845 step_rel
= step
- ir
->init_step
;
858 bLastStep
= (step_rel
== ir
->nsteps
);
859 t
= t0
+ step
*ir
->delta_t
;
862 // TODO Refactor this, so that nstfep does not need a default value of zero
863 if (ir
->efep
!= efepNO
|| ir
->bSimTemp
)
865 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
866 requiring different logic. */
871 setCurrentLambdasRerun(step
, ir
->fepvals
, &rerun_fr
, lam0
, state_global
);
876 setCurrentLambdasLocal(step
, ir
->fepvals
, lam0
, state
);
878 bDoDHDL
= do_per_step(step
, ir
->fepvals
->nstdhdl
);
879 bDoFEP
= ((ir
->efep
!= efepNO
) && do_per_step(step
, nstfep
));
880 bDoExpanded
= (do_per_step(step
, ir
->expandedvals
->nstexpanded
)
881 && (ir
->bExpanded
) && (step
> 0) && (!startingFromCheckpoint
));
884 bDoReplEx
= (useReplicaExchange
&& (step
> 0) && !bLastStep
&&
885 do_per_step(step
, replExParams
.exchangeInterval
));
889 update_annealing_target_temp(ir
, t
, upd
);
892 if (bRerunMD
&& MASTER(cr
))
894 const bool constructVsites
= ((vsite
!= nullptr) && mdrunOptions
.rerunConstructVsites
);
895 if (constructVsites
&& DOMAINDECOMP(cr
))
897 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
899 prepareRerunState(rerun_fr
, state_global
, constructVsites
, vsite
, top
->idef
, ir
->delta_t
, *fr
, graph
, &bRerunWarnNoV
);
902 /* Stop Center of Mass motion */
903 bStopCM
= (ir
->comm_mode
!= ecmNO
&& do_per_step(step
, ir
->nstcomm
));
907 /* for rerun MD always do Neighbour Searching */
908 bNS
= (bFirstStep
|| ir
->nstlist
!= 0);
912 /* Determine whether or not to do Neighbour Searching */
913 bNS
= (bFirstStep
|| bNStList
|| bExchanged
|| bNeedRepartition
);
916 /* < 0 means stop at next step, > 0 means stop at next NS step */
917 if ( (signals
[eglsSTOPCOND
].set
< 0) ||
918 ( (signals
[eglsSTOPCOND
].set
> 0 ) && ( bNS
|| ir
->nstlist
== 0)))
923 /* do_log triggers energy and virial calculation. Because this leads
924 * to different code paths, forces can be different. Thus for exact
925 * continuation we should avoid extra log output.
926 * Note that the || bLastStep can result in non-exact continuation
927 * beyond the last step. But we don't consider that to be an issue.
929 do_log
= do_per_step(step
, ir
->nstlog
) || (bFirstStep
&& !startingFromCheckpoint
) || bLastStep
|| bRerunMD
;
930 do_verbose
= mdrunOptions
.verbose
&&
931 (step
% mdrunOptions
.verboseStepPrintInterval
== 0 || bFirstStep
|| bLastStep
|| bRerunMD
);
933 if (bNS
&& !(bFirstStep
&& ir
->bContinuation
&& !bRerunMD
))
941 bMasterState
= FALSE
;
942 /* Correct the new box if it is too skewed */
943 if (inputrecDynamicBox(ir
))
945 if (correct_box(fplog
, step
, state
->box
, graph
))
950 if (DOMAINDECOMP(cr
) && bMasterState
)
952 dd_collect_state(cr
->dd
, state
, state_global
);
956 if (DOMAINDECOMP(cr
))
958 /* Repartition the domain decomposition */
959 dd_partition_system(fplog
, mdlog
, step
, cr
,
960 bMasterState
, nstglobalcomm
,
961 state_global
, top_global
, ir
,
962 state
, &f
, mdAtoms
, top
, fr
,
965 do_verbose
&& !bPMETunePrinting
);
966 shouldCheckNumberOfBondedInteractions
= true;
967 update_realloc(upd
, state
->natoms
);
971 if (MASTER(cr
) && do_log
)
973 print_ebin_header(fplog
, step
, t
); /* can we improve the information printed here? */
976 if (ir
->efep
!= efepNO
)
978 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
981 if ((bRerunMD
&& rerun_fr
.bV
) || bExchanged
)
984 /* We need the kinetic energy at minus the half step for determining
985 * the full step kinetic energy and possibly for T-coupling.*/
986 /* This may not be quite working correctly yet . . . . */
987 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
988 wcycle
, enerd
, nullptr, nullptr, nullptr, nullptr, mu_tot
,
989 constr
, &nullSignaller
, state
->box
,
990 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
991 CGLO_GSTAT
| CGLO_TEMPERATURE
| CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
);
992 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
993 top_global
, top
, state
,
994 &shouldCheckNumberOfBondedInteractions
);
996 clear_mat(force_vir
);
1000 checkpointHandler
->decideIfCheckpointingThisStep(bNS
, bFirstStep
, bLastStep
);
1003 /* Determine the energy and pressure:
1004 * at nstcalcenergy steps and at energy output steps (set below).
1006 if (EI_VV(ir
->eI
) && (!bInitStep
))
1008 /* for vv, the first half of the integration actually corresponds
1009 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1010 but the virial needs to be calculated on both the current step and the 'next' step. Future
1011 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1013 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1014 bCalcEnerStep
= do_per_step(step
- 1, ir
->nstcalcenergy
);
1015 bCalcVir
= bCalcEnerStep
||
1016 (ir
->epc
!= epcNO
&& (do_per_step(step
, ir
->nstpcouple
) || do_per_step(step
-1, ir
->nstpcouple
)));
1020 bCalcEnerStep
= do_per_step(step
, ir
->nstcalcenergy
);
1021 bCalcVir
= bCalcEnerStep
||
1022 (ir
->epc
!= epcNO
&& do_per_step(step
, ir
->nstpcouple
));
1024 bCalcEner
= bCalcEnerStep
;
1026 do_ene
= (do_per_step(step
, ir
->nstenergy
) || bLastStep
|| bRerunMD
);
1028 if (do_ene
|| do_log
|| bDoReplEx
)
1034 /* Do we need global communication ? */
1035 bGStat
= (bCalcVir
|| bCalcEner
|| bStopCM
||
1036 do_per_step(step
, nstglobalcomm
) ||
1037 (EI_VV(ir
->eI
) && inputrecNvtTrotter(ir
) && do_per_step(step
-1, nstglobalcomm
)));
1039 force_flags
= (GMX_FORCE_STATECHANGED
|
1040 ((inputrecDynamicBox(ir
) || bRerunMD
) ? GMX_FORCE_DYNAMICBOX
: 0) |
1041 GMX_FORCE_ALLFORCES
|
1042 (bCalcVir
? GMX_FORCE_VIRIAL
: 0) |
1043 (bCalcEner
? GMX_FORCE_ENERGY
: 0) |
1044 (bDoFEP
? GMX_FORCE_DHDL
: 0)
1049 /* Now is the time to relax the shells */
1050 relax_shell_flexcon(fplog
, cr
, ms
, mdrunOptions
.verbose
,
1051 enforcedRotation
, step
,
1052 ir
, bNS
, force_flags
, top
,
1054 state
, f
, force_vir
, mdatoms
,
1055 nrnb
, wcycle
, graph
, groups
,
1056 shellfc
, fr
, t
, mu_tot
,
1058 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
1062 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1063 (or the AWH update will be performed twice for one step when continuing). It would be best to
1064 call this update function from do_md_trajectory_writing but that would occur after do_force.
1065 One would have to divide the update_awh function into one function applying the AWH force
1066 and one doing the AWH bias update. The update AWH bias function could then be called after
1067 do_md_trajectory_writing (then containing update_awh_history).
1068 The checkpointing will in the future probably moved to the start of the md loop which will
1069 rid of this issue. */
1070 if (awh
&& !bRerunMD
&& checkpointHandler
->isCheckpointingStep() && MASTER(cr
))
1072 awh
->updateHistory(state_global
->awhHistory
.get());
1075 /* The coordinates (x) are shifted (to get whole molecules)
1077 * This is parallellized as well, and does communication too.
1078 * Check comments in sim_util.c
1080 do_force(fplog
, cr
, ms
, ir
, awh
.get(), enforcedRotation
,
1081 step
, nrnb
, wcycle
, top
, groups
,
1082 state
->box
, state
->x
, &state
->hist
,
1083 f
, force_vir
, mdatoms
, enerd
, fcd
,
1084 state
->lambda
, graph
,
1085 fr
, vsite
, mu_tot
, t
, ed
? ed
->getLegacyED() : nullptr,
1086 (bNS
? GMX_FORCE_NS
: 0) | force_flags
,
1087 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
1090 if (EI_VV(ir
->eI
) && !startingFromCheckpoint
&& !bRerunMD
)
1091 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1093 rvec
*vbuf
= nullptr;
1095 wallcycle_start(wcycle
, ewcUPDATE
);
1096 if (ir
->eI
== eiVV
&& bInitStep
)
1098 /* if using velocity verlet with full time step Ekin,
1099 * take the first half step only to compute the
1100 * virial for the first step. From there,
1101 * revert back to the initial coordinates
1102 * so that the input is actually the initial step.
1104 snew(vbuf
, state
->natoms
);
1105 copy_rvecn(as_rvec_array(state
->v
.data()), vbuf
, 0, state
->natoms
); /* should make this better for parallelizing? */
1109 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1110 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ1
);
1113 update_coords(step
, ir
, mdatoms
, state
, f
, fcd
,
1114 ekind
, M
, upd
, etrtVELOCITY1
,
1117 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
) /* Why is rerun_fr.bV here? Unclear. */
1119 wallcycle_stop(wcycle
, ewcUPDATE
);
1120 constrain_velocities(step
, nullptr,
1124 bCalcVir
, do_log
, do_ene
);
1125 wallcycle_start(wcycle
, ewcUPDATE
);
1129 /* Need to unshift here if a do_force has been
1130 called in the previous step */
1131 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1133 /* if VV, compute the pressure and constraints */
1134 /* For VV2, we strictly only need this if using pressure
1135 * control, but we really would like to have accurate pressures
1137 * Think about ways around this in the future?
1138 * For now, keep this choice in comments.
1140 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1141 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1143 bTemp
= ((ir
->eI
== eiVV
&& (!bInitStep
)) || (ir
->eI
== eiVVAK
));
1144 if (bCalcEner
&& ir
->eI
== eiVVAK
)
1146 bSumEkinhOld
= TRUE
;
1148 /* for vv, the first half of the integration actually corresponds to the previous step.
1149 So we need information from the last step in the first half of the integration */
1150 if (bGStat
|| do_per_step(step
-1, nstglobalcomm
))
1152 wallcycle_stop(wcycle
, ewcUPDATE
);
1153 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1154 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1155 constr
, &nullSignaller
, state
->box
,
1156 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1157 (bGStat
? CGLO_GSTAT
: 0)
1159 | (bTemp
? CGLO_TEMPERATURE
: 0)
1160 | (bPres
? CGLO_PRESSURE
: 0)
1161 | (bPres
? CGLO_CONSTRAINT
: 0)
1162 | (bStopCM
? CGLO_STOPCM
: 0)
1163 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
1166 /* explanation of above:
1167 a) We compute Ekin at the full time step
1168 if 1) we are using the AveVel Ekin, and it's not the
1169 initial step, or 2) if we are using AveEkin, but need the full
1170 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1171 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1172 EkinAveVel because it's needed for the pressure */
1173 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
1174 top_global
, top
, state
,
1175 &shouldCheckNumberOfBondedInteractions
);
1176 wallcycle_start(wcycle
, ewcUPDATE
);
1178 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1183 m_add(force_vir
, shake_vir
, total_vir
); /* we need the un-dispersion corrected total vir here */
1184 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ2
);
1186 /* TODO This is only needed when we're about to write
1187 * a checkpoint, because we use it after the restart
1188 * (in a kludge?). But what should we be doing if
1189 * startingFromCheckpoint or bInitStep are true? */
1190 if (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
))
1192 copy_mat(shake_vir
, state
->svir_prev
);
1193 copy_mat(force_vir
, state
->fvir_prev
);
1195 if (inputrecNvtTrotter(ir
) && ir
->eI
== eiVV
)
1197 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1198 enerd
->term
[F_TEMP
] = sum_ekin(&(ir
->opts
), ekind
, nullptr, (ir
->eI
== eiVV
), FALSE
);
1199 enerd
->term
[F_EKIN
] = trace(ekind
->ekin
);
1202 else if (bExchanged
)
1204 wallcycle_stop(wcycle
, ewcUPDATE
);
1205 /* We need the kinetic energy at minus the half step for determining
1206 * the full step kinetic energy and possibly for T-coupling.*/
1207 /* This may not be quite working correctly yet . . . . */
1208 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1209 wcycle
, enerd
, nullptr, nullptr, nullptr, nullptr, mu_tot
,
1210 constr
, &nullSignaller
, state
->box
,
1211 nullptr, &bSumEkinhOld
,
1212 CGLO_GSTAT
| CGLO_TEMPERATURE
);
1213 wallcycle_start(wcycle
, ewcUPDATE
);
1216 /* if it's the initial step, we performed this first step just to get the constraint virial */
1217 if (ir
->eI
== eiVV
&& bInitStep
)
1219 copy_rvecn(vbuf
, as_rvec_array(state
->v
.data()), 0, state
->natoms
);
1222 wallcycle_stop(wcycle
, ewcUPDATE
);
1225 /* compute the conserved quantity */
1228 saved_conserved_quantity
= NPT_energy(ir
, state
, &MassQ
);
1231 last_ekin
= enerd
->term
[F_EKIN
];
1233 if ((ir
->eDispCorr
!= edispcEnerPres
) && (ir
->eDispCorr
!= edispcAllEnerPres
))
1235 saved_conserved_quantity
-= enerd
->term
[F_DISPCORR
];
1237 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1238 if (ir
->efep
!= efepNO
&& !bRerunMD
)
1240 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1244 /* ######## END FIRST UPDATE STEP ############## */
1245 /* ######## If doing VV, we now have v(dt) ###### */
1248 /* perform extended ensemble sampling in lambda - we don't
1249 actually move to the new state before outputting
1250 statistics, but if performing simulated tempering, we
1251 do update the velocities and the tau_t. */
1253 lamnew
= ExpandedEnsembleDynamics(fplog
, ir
, enerd
, state
, &MassQ
, state
->fep_state
, state
->dfhist
, step
, as_rvec_array(state
->v
.data()), mdatoms
);
1254 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1257 copy_df_history(state_global
->dfhist
, state
->dfhist
);
1261 /* Now we have the energies and forces corresponding to the
1262 * coordinates at time t. We must output all of this before
1265 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
1266 ir
, state
, state_global
, observablesHistory
,
1268 outf
, mdebin
, ekind
, f
,
1269 !bRerunMD
&& checkpointHandler
->isCheckpointingStep(),
1270 bRerunMD
, bLastStep
,
1271 mdrunOptions
.writeConfout
,
1273 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1274 bIMDstep
= do_IMD(ir
->bIMD
, step
, cr
, bNS
, state
->box
, as_rvec_array(state
->x
.data()), ir
, t
, wcycle
);
1276 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1277 if (startingFromCheckpoint
&& (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
)))
1279 copy_mat(state
->svir_prev
, shake_vir
);
1280 copy_mat(state
->fvir_prev
, force_vir
);
1283 double secondsSinceStart
= walltime_accounting_get_time_since_start(walltime_accounting
);
1285 /* Check whether everything is still allright */
1286 if ((static_cast<int>(gmx_get_stop_condition()) > handled_stop_condition
)
1292 int nsteps_stop
= -1;
1294 /* this just makes signals[].sig compatible with the hack
1295 of sending signals around by MPI_Reduce together with
1297 if ((gmx_get_stop_condition() == gmx_stop_cond_next_ns
) ||
1298 (mdrunOptions
.reproducible
&&
1299 gmx_get_stop_condition() == gmx_stop_cond_next
))
1301 /* We need at least two global communication steps to pass
1302 * around the signal. We stop at a pair-list creation step
1303 * to allow for exact continuation, when possible.
1305 signals
[eglsSTOPCOND
].sig
= 1;
1306 nsteps_stop
= std::max(ir
->nstlist
, 2*nstSignalComm
);
1308 else if (gmx_get_stop_condition() == gmx_stop_cond_next
)
1310 /* Stop directly after the next global communication step.
1311 * This breaks exact continuation.
1313 signals
[eglsSTOPCOND
].sig
= -1;
1314 nsteps_stop
= nstSignalComm
+ 1;
1319 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1320 gmx_get_signal_name(), nsteps_stop
);
1324 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1325 gmx_get_signal_name(), nsteps_stop
);
1327 handled_stop_condition
= static_cast<int>(gmx_get_stop_condition());
1329 else if (MASTER(cr
) && (bNS
|| ir
->nstlist
<= 0) &&
1330 (mdrunOptions
.maximumHoursToRun
> 0 &&
1331 secondsSinceStart
> mdrunOptions
.maximumHoursToRun
*60.0*60.0*0.99) &&
1332 signals
[eglsSTOPCOND
].sig
== 0 && signals
[eglsSTOPCOND
].set
== 0)
1334 /* Signal to terminate the run */
1335 signals
[eglsSTOPCOND
].sig
= 1;
1338 fprintf(fplog
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",
1339 gmx_step_str(step
, sbuf
), mdrunOptions
.maximumHoursToRun
*0.99);
1341 fprintf(stderr
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",
1342 gmx_step_str(step
, sbuf
), mdrunOptions
.maximumHoursToRun
*0.99);
1345 resetHandler
->setSignal(walltime_accounting
);
1347 if ((bGStat
|| !PAR(cr
)) && !bRerunMD
)
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 checkpointHandler
->setSignal(walltime_accounting
);
1357 /* ######### START SECOND UPDATE STEP ################# */
1359 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1362 if (ETC_ANDERSEN(ir
->etc
)) /* keep this outside of update_tcouple because of the extra info required to pass */
1364 gmx_bool bIfRandomize
;
1365 bIfRandomize
= update_randomize_velocities(ir
, step
, cr
, mdatoms
, state
, upd
, constr
);
1366 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1367 if (constr
&& bIfRandomize
)
1369 constrain_velocities(step
, nullptr,
1373 bCalcVir
, do_log
, do_ene
);
1376 /* Box is changed in update() when we do pressure coupling,
1377 * but we should still use the old box for energy corrections and when
1378 * writing it to the energy file, so it matches the trajectory files for
1379 * the same timestep above. Make a copy in a separate array.
1381 copy_mat(state
->box
, lastbox
);
1385 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
)
1387 wallcycle_start(wcycle
, ewcUPDATE
);
1388 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1391 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ3
);
1392 /* We can only do Berendsen coupling after we have summed
1393 * the kinetic energy or virial. Since the happens
1394 * in global_state after update, we should only do it at
1395 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1400 update_tcouple(step
, ir
, state
, ekind
, &MassQ
, mdatoms
);
1401 update_pcouple_before_coordinates(fplog
, step
, ir
, state
,
1402 parrinellorahmanMu
, M
,
1408 /* velocity half-step update */
1409 update_coords(step
, ir
, mdatoms
, state
, f
, fcd
,
1410 ekind
, M
, upd
, etrtVELOCITY2
,
1414 /* Above, initialize just copies ekinh into ekin,
1415 * it doesn't copy position (for VV),
1416 * and entire integrator for MD.
1419 if (ir
->eI
== eiVVAK
)
1421 /* We probably only need md->homenr, not state->natoms */
1422 if (state
->natoms
> cbuf_nalloc
)
1424 cbuf_nalloc
= state
->natoms
;
1425 srenew(cbuf
, cbuf_nalloc
);
1427 copy_rvecn(as_rvec_array(state
->x
.data()), cbuf
, 0, state
->natoms
);
1430 update_coords(step
, ir
, mdatoms
, state
, f
, fcd
,
1431 ekind
, M
, upd
, etrtPOSITION
, cr
, constr
);
1432 wallcycle_stop(wcycle
, ewcUPDATE
);
1434 constrain_coordinates(step
, &dvdl_constr
, state
,
1436 wcycle
, upd
, constr
,
1437 bCalcVir
, do_log
, do_ene
);
1438 update_sd_second_half(step
, &dvdl_constr
, ir
, mdatoms
, state
,
1439 cr
, nrnb
, wcycle
, upd
, constr
, do_log
, do_ene
);
1440 finish_update(ir
, mdatoms
,
1442 nrnb
, wcycle
, upd
, constr
);
1444 if (ir
->eI
== eiVVAK
)
1446 /* erase F_EKIN and F_TEMP here? */
1447 /* just compute the kinetic energy at the half step to perform a trotter step */
1448 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1449 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1450 constr
, &nullSignaller
, lastbox
,
1451 nullptr, &bSumEkinhOld
,
1452 (bGStat
? CGLO_GSTAT
: 0) | CGLO_TEMPERATURE
1454 wallcycle_start(wcycle
, ewcUPDATE
);
1455 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ4
);
1456 /* now we know the scaling, we can compute the positions again again */
1457 copy_rvecn(cbuf
, as_rvec_array(state
->x
.data()), 0, state
->natoms
);
1459 update_coords(step
, ir
, mdatoms
, state
, f
, fcd
,
1460 ekind
, M
, upd
, etrtPOSITION
, cr
, constr
);
1461 wallcycle_stop(wcycle
, ewcUPDATE
);
1463 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1464 /* are the small terms in the shake_vir here due
1465 * to numerical errors, or are they important
1466 * physically? I'm thinking they are just errors, but not completely sure.
1467 * For now, will call without actually constraining, constr=NULL*/
1468 finish_update(ir
, mdatoms
,
1470 nrnb
, wcycle
, upd
, nullptr);
1474 /* this factor or 2 correction is necessary
1475 because half of the constraint force is removed
1476 in the vv step, so we have to double it. See
1477 the Redmine issue #1255. It is not yet clear
1478 if the factor of 2 is exact, or just a very
1479 good approximation, and this will be
1480 investigated. The next step is to see if this
1481 can be done adding a dhdl contribution from the
1482 rattle step, but this is somewhat more
1483 complicated with the current code. Will be
1484 investigated, hopefully for 4.6.3. However,
1485 this current solution is much better than
1486 having it completely wrong.
1488 enerd
->term
[F_DVDL_CONSTR
] += 2*dvdl_constr
;
1492 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
1497 /* Need to unshift here */
1498 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1501 if (vsite
!= nullptr)
1503 wallcycle_start(wcycle
, ewcVSITECONSTR
);
1504 if (graph
!= nullptr)
1506 shift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1508 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, as_rvec_array(state
->v
.data()),
1509 top
->idef
.iparams
, top
->idef
.il
,
1510 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
1512 if (graph
!= nullptr)
1514 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1516 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
1519 /* ############## IF NOT VV, Calculate globals HERE ############ */
1520 /* With Leap-Frog we can skip compute_globals at
1521 * non-communication steps, but we need to calculate
1522 * the kinetic energy one step before communication.
1525 // Organize to do inter-simulation signalling on steps if
1526 // and when algorithms require it.
1527 bool doInterSimSignal
= (simulationsShareState
&& do_per_step(step
, nstSignalComm
));
1529 if (bGStat
|| (!EI_VV(ir
->eI
) && do_per_step(step
+1, nstglobalcomm
)) || doInterSimSignal
)
1531 // Since we're already communicating at this step, we
1532 // can propagate intra-simulation signals. Note that
1533 // check_nstglobalcomm has the responsibility for
1534 // choosing the value of nstglobalcomm that is one way
1535 // bGStat becomes true, so we can't get into a
1536 // situation where e.g. checkpointing can't be
1538 bool doIntraSimSignal
= true;
1539 SimulationSignaller
signaller(&signals
, cr
, ms
, doInterSimSignal
, doIntraSimSignal
);
1541 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1542 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1545 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1546 (bGStat
? CGLO_GSTAT
: 0)
1547 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_ENERGY
: 0)
1548 | (!EI_VV(ir
->eI
) && bStopCM
? CGLO_STOPCM
: 0)
1549 | (!EI_VV(ir
->eI
) ? CGLO_TEMPERATURE
: 0)
1550 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_PRESSURE
: 0)
1552 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
1554 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
1555 top_global
, top
, state
,
1556 &shouldCheckNumberOfBondedInteractions
);
1560 /* ############# END CALC EKIN AND PRESSURE ################# */
1562 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1563 the virial that should probably be addressed eventually. state->veta has better properies,
1564 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1565 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1567 if (ir
->efep
!= efepNO
&& (!EI_VV(ir
->eI
) || bRerunMD
))
1569 /* Sum up the foreign energy and dhdl terms for md and sd.
1570 Currently done every step so that dhdl is correct in the .edr */
1571 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1574 update_pcouple_after_coordinates(fplog
, step
, ir
, mdatoms
,
1575 pres
, force_vir
, shake_vir
,
1579 /* ################# END UPDATE STEP 2 ################# */
1580 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1582 /* The coordinates (x) were unshifted in update */
1585 /* We will not sum ekinh_old,
1586 * so signal that we still have to do it.
1588 bSumEkinhOld
= TRUE
;
1593 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1595 /* use the directly determined last velocity, not actually the averaged half steps */
1596 if (bTrotter
&& ir
->eI
== eiVV
)
1598 enerd
->term
[F_EKIN
] = last_ekin
;
1600 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1602 if (integratorHasConservedEnergyQuantity(ir
))
1606 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + saved_conserved_quantity
;
1610 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + NPT_energy(ir
, state
, &MassQ
);
1613 /* ######### END PREPARING EDR OUTPUT ########### */
1619 if (fplog
&& do_log
&& bDoExpanded
)
1621 /* only needed if doing expanded ensemble */
1622 PrintFreeEnergyInfoToFile(fplog
, ir
->fepvals
, ir
->expandedvals
, ir
->bSimTemp
? ir
->simtempvals
: nullptr,
1623 state_global
->dfhist
, state
->fep_state
, ir
->nstlog
, step
);
1627 upd_mdebin(mdebin
, bDoDHDL
, bCalcEnerStep
,
1628 t
, mdatoms
->tmass
, enerd
, state
,
1629 ir
->fepvals
, ir
->expandedvals
, lastbox
,
1630 shake_vir
, force_vir
, total_vir
, pres
,
1631 ekind
, mu_tot
, constr
);
1635 upd_mdebin_step(mdebin
);
1638 gmx_bool do_dr
= do_per_step(step
, ir
->nstdisreout
);
1639 gmx_bool do_or
= do_per_step(step
, ir
->nstorireout
);
1641 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
, do_log
? fplog
: nullptr,
1643 eprNORMAL
, mdebin
, fcd
, groups
, &(ir
->opts
), awh
.get());
1647 pull_print_output(ir
->pull_work
, step
, t
);
1650 if (do_per_step(step
, ir
->nstlog
))
1652 if (fflush(fplog
) != 0)
1654 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
1660 /* Have to do this part _after_ outputting the logfile and the edr file */
1661 /* Gets written into the state at the beginning of next loop*/
1662 state
->fep_state
= lamnew
;
1664 /* Print the remaining wall clock time for the run */
1665 if (isMasterSimMasterRank(ms
, cr
) &&
1666 (do_verbose
|| gmx_got_usr_signal()) &&
1671 fprintf(stderr
, "\n");
1673 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
1676 /* Ion/water position swapping.
1677 * Not done in last step since trajectory writing happens before this call
1678 * in the MD loop and exchanges would be lost anyway. */
1679 bNeedRepartition
= FALSE
;
1680 if ((ir
->eSwapCoords
!= eswapNO
) && (step
> 0) && !bLastStep
&&
1681 do_per_step(step
, ir
->swap
->nstswap
))
1683 bNeedRepartition
= do_swapcoords(cr
, step
, t
, ir
, wcycle
,
1684 bRerunMD
? rerun_fr
.x
: as_rvec_array(state
->x
.data()),
1685 bRerunMD
? rerun_fr
.box
: state
->box
,
1686 MASTER(cr
) && mdrunOptions
.verbose
,
1689 if (bNeedRepartition
&& DOMAINDECOMP(cr
))
1691 dd_collect_state(cr
->dd
, state
, state_global
);
1695 /* Replica exchange */
1699 bExchanged
= replica_exchange(fplog
, cr
, ms
, repl_ex
,
1700 state_global
, enerd
,
1704 if ( (bExchanged
|| bNeedRepartition
) && DOMAINDECOMP(cr
) )
1706 dd_partition_system(fplog
, mdlog
, step
, cr
, TRUE
, 1,
1707 state_global
, top_global
, ir
,
1708 state
, &f
, mdAtoms
, top
, fr
,
1710 nrnb
, wcycle
, FALSE
);
1711 shouldCheckNumberOfBondedInteractions
= true;
1712 update_realloc(upd
, state
->natoms
);
1717 startingFromCheckpoint
= false;
1719 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1720 /* With all integrators, except VV, we need to retain the pressure
1721 * at the current step for coupling at the next step.
1723 if ((state
->flags
& (1<<estPRES_PREV
)) &&
1725 (ir
->nstpcouple
> 0 && step
% ir
->nstpcouple
== 0)))
1727 /* Store the pressure in t_state for pressure coupling
1728 * at the next MD step.
1730 copy_mat(pres
, state
->pres_prev
);
1733 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1735 if ( (membed
!= nullptr) && (!bLastStep
) )
1737 rescale_membed(step_rel
, membed
, as_rvec_array(state_global
->x
.data()));
1744 /* read next frame from input trajectory */
1745 bLastStep
= !read_next_frame(oenv
, status
, &rerun_fr
);
1750 rerun_parallel_comm(cr
, &rerun_fr
, &bLastStep
);
1754 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
1755 if (DOMAINDECOMP(cr
) && wcycle
)
1757 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
1760 if (!bRerunMD
|| !rerun_fr
.bStep
)
1762 /* increase the MD step number */
1767 resetHandler
->resetCounters(
1768 step
, step_rel
, mdlog
, fplog
, cr
, (use_GPU(fr
->nbv
) ? fr
->nbv
: nullptr),
1769 nrnb
, fr
->pmedata
, pme_loadbal
, wcycle
, walltime_accounting
);
1771 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1772 IMD_prep_energies_send_positions(ir
->bIMD
&& MASTER(cr
), bIMDstep
, ir
->imd
, enerd
, step
, bCalcEner
, wcycle
);
1775 /* End of main MD loop */
1777 /* Closing TNG files can include compressing data. Therefore it is good to do that
1778 * before stopping the time measurements. */
1779 mdoutf_tng_close(outf
);
1781 /* Stop measuring walltime */
1782 walltime_accounting_end_time(walltime_accounting
);
1784 if (bRerunMD
&& MASTER(cr
))
1789 if (!thisRankHasDuty(cr
, DUTY_PME
))
1791 /* Tell the PME only node to finish */
1792 gmx_pme_send_finish(cr
);
1797 if (ir
->nstcalcenergy
> 0 && !bRerunMD
)
1799 print_ebin(mdoutf_get_fp_ene(outf
), FALSE
, FALSE
, FALSE
, fplog
, step
, t
,
1800 eprAVER
, mdebin
, fcd
, groups
, &(ir
->opts
), awh
.get());
1803 done_mdebin(mdebin
);
1808 pme_loadbal_done(pme_loadbal
, fplog
, mdlog
, use_GPU(fr
->nbv
));
1811 done_shellfc(fplog
, shellfc
, step_rel
);
1813 if (useReplicaExchange
&& MASTER(cr
))
1815 print_replica_exchange_statistics(fplog
, repl_ex
);
1818 // Clean up swapcoords
1819 if (ir
->eSwapCoords
!= eswapNO
)
1821 finish_swapcoords(ir
->swap
);
1824 /* IMD cleanup, if bIMD is TRUE. */
1825 IMD_finalize(ir
->bIMD
, ir
->imd
);
1827 walltime_accounting_set_nsteps_done(walltime_accounting
, step
);
1829 destroy_enerdata(enerd
);