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, 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.
51 #include "thread_mpi/threads.h"
53 #include "gromacs/awh/awh.h"
54 #include "gromacs/commandline/filenm.h"
55 #include "gromacs/compat/make_unique.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/essentialdynamics/edsam.h"
60 #include "gromacs/ewald/pme.h"
61 #include "gromacs/ewald/pme-load-balancing.h"
62 #include "gromacs/fileio/trxio.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/gmxlib/nrnb.h"
65 #include "gromacs/gpu_utils/gpu_utils.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/listed-forces/manage-threading.h"
68 #include "gromacs/math/functions.h"
69 #include "gromacs/math/utilities.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/math/vectypes.h"
72 #include "gromacs/mdlib/compute_io.h"
73 #include "gromacs/mdlib/constr.h"
74 #include "gromacs/mdlib/ebin.h"
75 #include "gromacs/mdlib/expanded.h"
76 #include "gromacs/mdlib/force.h"
77 #include "gromacs/mdlib/force_flags.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/md_support.h"
80 #include "gromacs/mdlib/mdatoms.h"
81 #include "gromacs/mdlib/mdebin.h"
82 #include "gromacs/mdlib/mdoutf.h"
83 #include "gromacs/mdlib/mdrun.h"
84 #include "gromacs/mdlib/mdsetup.h"
85 #include "gromacs/mdlib/nb_verlet.h"
86 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
87 #include "gromacs/mdlib/ns.h"
88 #include "gromacs/mdlib/shellfc.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/sim_util.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/tgroup.h"
93 #include "gromacs/mdlib/trajectory_writing.h"
94 #include "gromacs/mdlib/update.h"
95 #include "gromacs/mdlib/vcm.h"
96 #include "gromacs/mdlib/vsite.h"
97 #include "gromacs/mdtypes/awh-history.h"
98 #include "gromacs/mdtypes/commrec.h"
99 #include "gromacs/mdtypes/df_history.h"
100 #include "gromacs/mdtypes/energyhistory.h"
101 #include "gromacs/mdtypes/fcdata.h"
102 #include "gromacs/mdtypes/forcerec.h"
103 #include "gromacs/mdtypes/group.h"
104 #include "gromacs/mdtypes/inputrec.h"
105 #include "gromacs/mdtypes/interaction_const.h"
106 #include "gromacs/mdtypes/md_enums.h"
107 #include "gromacs/mdtypes/mdatom.h"
108 #include "gromacs/mdtypes/observableshistory.h"
109 #include "gromacs/mdtypes/state.h"
110 #include "gromacs/pbcutil/mshift.h"
111 #include "gromacs/pbcutil/pbc.h"
112 #include "gromacs/pulling/pull.h"
113 #include "gromacs/swap/swapcoords.h"
114 #include "gromacs/timing/wallcycle.h"
115 #include "gromacs/timing/walltime_accounting.h"
116 #include "gromacs/topology/atoms.h"
117 #include "gromacs/topology/idef.h"
118 #include "gromacs/topology/mtop_util.h"
119 #include "gromacs/topology/topology.h"
120 #include "gromacs/trajectory/trajectoryframe.h"
121 #include "gromacs/utility/basedefinitions.h"
122 #include "gromacs/utility/cstringutil.h"
123 #include "gromacs/utility/fatalerror.h"
124 #include "gromacs/utility/logger.h"
125 #include "gromacs/utility/real.h"
126 #include "gromacs/utility/smalloc.h"
133 #include "corewrap.h"
136 using gmx::SimulationSignaller
;
138 /*! \brief Check whether bonded interactions are missing, if appropriate
140 * \param[in] fplog Log file pointer
141 * \param[in] cr Communication object
142 * \param[in] totalNumberOfBondedInteractions Result of the global reduction over the number of bonds treated in each domain
143 * \param[in] top_global Global topology for the error message
144 * \param[in] top_local Local topology for the error message
145 * \param[in] state Global state for the error message
146 * \param[inout] shouldCheckNumberOfBondedInteractions Whether we should do the check.
148 * \return Nothing, except that shouldCheckNumberOfBondedInteractions
149 * is always set to false after exit.
151 static void checkNumberOfBondedInteractions(FILE *fplog
, t_commrec
*cr
, int totalNumberOfBondedInteractions
,
152 gmx_mtop_t
*top_global
, gmx_localtop_t
*top_local
, t_state
*state
,
153 bool *shouldCheckNumberOfBondedInteractions
)
155 if (*shouldCheckNumberOfBondedInteractions
)
157 if (totalNumberOfBondedInteractions
!= cr
->dd
->nbonded_global
)
159 dd_print_missing_interactions(fplog
, cr
, totalNumberOfBondedInteractions
, top_global
, top_local
, state
); // Does not return
161 *shouldCheckNumberOfBondedInteractions
= false;
165 static void reset_all_counters(FILE *fplog
, const gmx::MDLogger
&mdlog
, t_commrec
*cr
,
167 gmx_int64_t
*step_rel
, t_inputrec
*ir
,
168 gmx_wallcycle_t wcycle
, t_nrnb
*nrnb
,
169 gmx_walltime_accounting_t walltime_accounting
,
170 struct nonbonded_verlet_t
*nbv
,
171 struct gmx_pme_t
*pme
)
173 char sbuf
[STEPSTRSIZE
];
175 /* Reset all the counters related to performance over the run */
176 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
177 "step %s: resetting all time and cycle counters",
178 gmx_step_str(step
, sbuf
));
182 nbnxn_gpu_reset_timings(nbv
);
185 if (pme_gpu_task_enabled(pme
))
187 pme_gpu_reset_timings(pme
);
190 if (use_GPU(nbv
) || pme_gpu_task_enabled(pme
))
195 wallcycle_stop(wcycle
, ewcRUN
);
196 wallcycle_reset_all(wcycle
);
197 if (DOMAINDECOMP(cr
))
199 reset_dd_statistics_counters(cr
->dd
);
202 ir
->init_step
+= *step_rel
;
203 ir
->nsteps
-= *step_rel
;
205 wallcycle_start(wcycle
, ewcRUN
);
206 walltime_accounting_start(walltime_accounting
);
207 print_date_and_time(fplog
, cr
->nodeid
, "Restarted time", gmx_gettime());
210 /*! \brief Copy the state from \p rerunFrame to \p globalState and, if requested, construct vsites
212 * \param[in] rerunFrame The trajectory frame to compute energy/forces for
213 * \param[in,out] globalState The global state container
214 * \param[in] constructVsites When true, vsite coordinates are constructed
215 * \param[in] vsite Vsite setup, can be nullptr when \p constructVsites = false
216 * \param[in] idef Topology parameters, used for constructing vsites
217 * \param[in] timeStep Time step, used for constructing vsites
218 * \param[in] forceRec Force record, used for constructing vsites
219 * \param[in,out] graph The molecular graph, used for constructing vsites when != nullptr
220 * \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
222 static void prepareRerunState(const t_trxframe
&rerunFrame
,
223 t_state
*globalState
,
224 bool constructVsites
,
225 const gmx_vsite_t
*vsite
,
228 const t_forcerec
&forceRec
,
230 gmx_bool
*warnWhenNoV
)
232 for (int i
= 0; i
< globalState
->natoms
; i
++)
234 copy_rvec(rerunFrame
.x
[i
], globalState
->x
[i
]);
238 for (int i
= 0; i
< globalState
->natoms
; i
++)
240 copy_rvec(rerunFrame
.v
[i
], globalState
->v
[i
]);
245 for (int i
= 0; i
< globalState
->natoms
; i
++)
247 clear_rvec(globalState
->v
[i
]);
251 fprintf(stderr
, "\nWARNING: Some frames do not contain velocities.\n"
252 " Ekin, temperature and pressure are incorrect,\n"
253 " the virial will be incorrect when constraints are present.\n"
255 *warnWhenNoV
= FALSE
;
258 copy_mat(rerunFrame
.box
, globalState
->box
);
262 GMX_ASSERT(vsite
, "Need valid vsite for constructing vsites");
266 /* Following is necessary because the graph may get out of sync
267 * with the coordinates if we only have every N'th coordinate set
269 mk_mshift(nullptr, graph
, forceRec
.ePBC
, globalState
->box
, as_rvec_array(globalState
->x
.data()));
270 shift_self(graph
, globalState
->box
, as_rvec_array(globalState
->x
.data()));
272 construct_vsites(vsite
, as_rvec_array(globalState
->x
.data()), timeStep
, as_rvec_array(globalState
->v
.data()),
273 idef
.iparams
, idef
.il
,
274 forceRec
.ePBC
, forceRec
.bMolPBC
, nullptr, globalState
->box
);
277 unshift_self(graph
, globalState
->box
, as_rvec_array(globalState
->x
.data()));
283 \copydoc integrator_t (FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
284 int nfile, const t_filenm fnm[],
285 const gmx_output_env_t *oenv,
286 const MdrunOptions &mdrunOptions,
287 gmx_vsite_t *vsite, gmx_constr_t constr,
288 gmx::IMDOutputProvider *outputProvider,
289 t_inputrec *inputrec,
290 gmx_mtop_t *top_global, t_fcdata *fcd,
291 t_state *state_global,
293 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
295 const ReplicaExchangeParameters &replExParams,
296 gmx_membed_t *membed,
297 gmx_walltime_accounting_t walltime_accounting)
299 double gmx::do_md(FILE *fplog
, t_commrec
*cr
, const gmx::MDLogger
&mdlog
,
300 int nfile
, const t_filenm fnm
[],
301 const gmx_output_env_t
*oenv
,
302 const MdrunOptions
&mdrunOptions
,
303 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
304 gmx::IMDOutputProvider
*outputProvider
,
306 gmx_mtop_t
*top_global
,
308 t_state
*state_global
,
309 ObservablesHistory
*observablesHistory
,
310 gmx::MDAtoms
*mdAtoms
,
311 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
313 const ReplicaExchangeParameters
&replExParams
,
314 gmx_membed_t
*membed
,
315 gmx_walltime_accounting_t walltime_accounting
)
317 gmx_mdoutf_t outf
= nullptr;
318 gmx_int64_t step
, step_rel
;
320 double t
, t0
, lam0
[efptNR
];
321 gmx_bool bGStatEveryStep
, bGStat
, bCalcVir
, bCalcEnerStep
, bCalcEner
;
322 gmx_bool bNS
, bNStList
, bSimAnn
, bStopCM
,
323 bFirstStep
, bInitStep
, bLastStep
= FALSE
,
324 bBornRadii
, bUsingEnsembleRestraints
;
325 gmx_bool bDoDHDL
= FALSE
, bDoFEP
= FALSE
, bDoExpanded
= FALSE
;
326 gmx_bool do_ene
, do_log
, do_verbose
, bRerunWarnNoV
= TRUE
,
327 bForceUpdate
= FALSE
, bCPT
;
328 gmx_bool bMasterState
;
329 int force_flags
, cglo_flags
;
330 tensor force_vir
, shake_vir
, total_vir
, tmp_vir
, pres
;
335 matrix parrinellorahmanMu
, M
;
337 gmx_repl_ex_t repl_ex
= nullptr;
340 t_mdebin
*mdebin
= nullptr;
341 gmx_enerdata_t
*enerd
;
342 PaddedRVecVector f
{};
343 gmx_global_stat_t gstat
;
344 gmx_update_t
*upd
= nullptr;
345 t_graph
*graph
= nullptr;
346 gmx_groups_t
*groups
;
347 gmx_ekindata_t
*ekind
;
348 gmx_shellfc_t
*shellfc
;
349 gmx_bool bSumEkinhOld
, bDoReplEx
, bExchanged
, bNeedRepartition
;
350 gmx_bool bResetCountersHalfMaxH
= FALSE
;
351 gmx_bool bTemp
, bPres
, bTrotter
;
353 rvec
*cbuf
= nullptr;
360 real saved_conserved_quantity
= 0;
364 char sbuf
[STEPSTRSIZE
], sbuf2
[STEPSTRSIZE
];
365 int handled_stop_condition
= gmx_stop_cond_none
; /* compare to get_stop_condition*/
368 /* PME load balancing data for GPU kernels */
369 pme_load_balancing_t
*pme_loadbal
= nullptr;
370 gmx_bool bPMETune
= FALSE
;
371 gmx_bool bPMETunePrinting
= FALSE
;
374 gmx_bool bIMDstep
= FALSE
;
377 /* Temporary addition for FAHCORE checkpointing */
380 /* Domain decomposition could incorrectly miss a bonded
381 interaction, but checking for that requires a global
382 communication stage, which does not otherwise happen in DD
383 code. So we do that alongside the first global energy reduction
384 after a new DD is made. These variables handle whether the
385 check happens, and the result it returns. */
386 bool shouldCheckNumberOfBondedInteractions
= false;
387 int totalNumberOfBondedInteractions
= -1;
389 SimulationSignals signals
;
390 // Most global communnication stages don't propagate mdrun
391 // signals, and will use this object to achieve that.
392 SimulationSignaller
nullSignaller(nullptr, nullptr, false, false);
394 if (mdrunOptions
.timingOptions
.resetHalfway
)
398 /* Signal to reset the counters half the simulation steps. */
399 wcycle_set_reset_counters(wcycle
, ir
->nsteps
/2);
401 /* Signal to reset the counters halfway the simulation time. */
402 bResetCountersHalfMaxH
= (mdrunOptions
.maximumHoursToRun
> 0);
405 /* md-vv uses averaged full step velocities for T-control
406 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
407 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
408 bTrotter
= (EI_VV(ir
->eI
) && (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
) || inputrecNvtTrotter(ir
)));
410 const gmx_bool bRerunMD
= mdrunOptions
.rerun
;
411 int nstglobalcomm
= mdrunOptions
.globalCommunicationInterval
;
415 /* Since we don't know if the frames read are related in any way,
416 * rebuild the neighborlist at every step.
419 ir
->nstcalcenergy
= 1;
423 nstglobalcomm
= check_nstglobalcomm(mdlog
, nstglobalcomm
, ir
);
424 bGStatEveryStep
= (nstglobalcomm
== 1);
428 ir
->nstxout_compressed
= 0;
430 groups
= &top_global
->groups
;
432 gmx_edsam
*ed
= nullptr;
433 if (opt2bSet("-ei", nfile
, fnm
) || observablesHistory
->edsamHistory
!= nullptr)
435 /* Initialize essential dynamics sampling */
436 ed
= init_edsam(opt2fn_null("-ei", nfile
, fnm
), opt2fn("-eo", nfile
, fnm
),
439 state_global
, observablesHistory
,
440 oenv
, mdrunOptions
.continuationOptions
.appendFiles
);
443 if (ir
->eSwapCoords
!= eswapNO
)
445 /* Initialize ion swapping code */
446 init_swapcoords(fplog
, ir
, opt2fn_master("-swap", nfile
, fnm
, cr
),
448 state_global
, observablesHistory
,
449 cr
, oenv
, mdrunOptions
);
453 init_md(fplog
, cr
, outputProvider
, ir
, oenv
, mdrunOptions
,
454 &t
, &t0
, state_global
, lam0
,
455 nrnb
, top_global
, &upd
,
456 nfile
, fnm
, &outf
, &mdebin
,
457 force_vir
, shake_vir
, mu_tot
, &bSimAnn
, &vcm
, wcycle
);
459 clear_mat(total_vir
);
461 /* Energy terms and groups */
463 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
466 /* Kinetic energy data */
468 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind
);
469 /* Copy the cos acceleration to the groups struct */
470 ekind
->cosacc
.cos_accel
= ir
->cos_accel
;
472 gstat
= global_stat_init(ir
);
474 /* Check for polarizable models and flexible constraints */
475 shellfc
= init_shell_flexcon(fplog
,
476 top_global
, n_flexible_constraints(constr
),
477 ir
->nstcalcenergy
, DOMAINDECOMP(cr
));
479 if (shellfc
&& ir
->nstcalcenergy
!= 1)
481 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
);
483 if (shellfc
&& DOMAINDECOMP(cr
))
485 gmx_fatal(FARGS
, "Shell particles are not implemented with domain decomposition, use a single rank");
487 if (shellfc
&& ir
->bDoAwh
)
489 gmx_fatal(FARGS
, "AWH biasing does not support shell particles.");
492 if (inputrecDeform(ir
))
494 tMPI_Thread_mutex_lock(&deform_init_box_mutex
);
495 set_deform_reference_box(upd
,
496 deform_init_init_step_tpx
,
497 deform_init_box_tpx
);
498 tMPI_Thread_mutex_unlock(&deform_init_box_mutex
);
502 double io
= compute_io(ir
, top_global
->natoms
, groups
, mdebin
->ebin
->nener
, 1);
503 if ((io
> 2000) && MASTER(cr
))
506 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
511 // Local state only becomes valid now.
512 std::unique_ptr
<t_state
> stateInstance
;
515 if (DOMAINDECOMP(cr
))
517 top
= dd_init_local_top(top_global
);
519 stateInstance
= gmx::compat::make_unique
<t_state
>();
520 state
= stateInstance
.get();
521 dd_init_local_state(cr
->dd
, state_global
, state
);
525 state_change_natoms(state_global
, state_global
->natoms
);
526 /* We need to allocate one element extra, since we might use
527 * (unaligned) 4-wide SIMD loads to access rvec entries.
529 f
.resize(gmx::paddedRVecVectorSize(state_global
->natoms
));
530 /* Copy the pointer to the global state */
531 state
= state_global
;
534 mdAlgorithmsSetupAtomData(cr
, ir
, top_global
, top
, fr
,
535 &graph
, mdAtoms
, vsite
, shellfc
);
537 update_realloc(upd
, state
->natoms
);
540 /* Set up interactive MD (IMD) */
541 init_IMD(ir
, cr
, top_global
, fplog
, ir
->nstcalcenergy
, MASTER(cr
) ? as_rvec_array(state_global
->x
.data()) : nullptr,
542 nfile
, fnm
, oenv
, mdrunOptions
);
544 if (DOMAINDECOMP(cr
))
546 /* Distribute the charge groups over the nodes from the master node */
547 dd_partition_system(fplog
, ir
->init_step
, cr
, TRUE
, 1,
548 state_global
, top_global
, ir
,
549 state
, &f
, mdAtoms
, top
, fr
,
551 nrnb
, nullptr, FALSE
);
552 shouldCheckNumberOfBondedInteractions
= true;
553 update_realloc(upd
, state
->natoms
);
556 auto mdatoms
= mdAtoms
->mdatoms();
558 // NOTE: The global state is no longer used at this point.
559 // But state_global is still used as temporary storage space for writing
560 // the global state to file and potentially for replica exchange.
561 // (Global topology should persist.)
563 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
565 const ContinuationOptions
&continuationOptions
= mdrunOptions
.continuationOptions
;
566 bool startingFromCheckpoint
= continuationOptions
.startedFromCheckpoint
;
570 init_expanded_ensemble(startingFromCheckpoint
, ir
, state
->dfhist
);
575 if (startingFromCheckpoint
)
577 /* Update mdebin with energy history if appending to output files */
578 if (continuationOptions
.appendFiles
)
580 restore_energyhistory_from_state(mdebin
, observablesHistory
->energyHistory
.get());
582 else if (observablesHistory
->energyHistory
.get() != nullptr)
584 /* We might have read an energy history from checkpoint.
585 * As we are not appending, we want to restart the statistics.
586 * Free the allocated memory and reset the counts.
588 observablesHistory
->energyHistory
= {};
591 if (observablesHistory
->energyHistory
.get() == nullptr)
593 observablesHistory
->energyHistory
= std::unique_ptr
<energyhistory_t
>(new energyhistory_t
{});
595 /* Set the initial energy history in state by updating once */
596 update_energyhistory(observablesHistory
->energyHistory
.get(), mdebin
);
599 /* Initialize constraints */
600 if (constr
&& !DOMAINDECOMP(cr
))
602 set_constraints(constr
, top
, ir
, mdatoms
, cr
);
605 /* Initialize AWH and restore state from history in checkpoint if needed. */
608 ir
->awh
= new gmx::Awh(fplog
, *ir
, cr
, *ir
->awhParams
, opt2fn("-awh", nfile
, fnm
), ir
->pull_work
);
610 if (startingFromCheckpoint
)
612 /* Restore the AWH history read from checkpoint */
613 ir
->awh
->restoreStateFromHistory(MASTER(cr
) ? state_global
->awhHistory
.get() : nullptr);
617 /* Initialize the AWH history here */
618 state_global
->awhHistory
= ir
->awh
->initHistoryFromState();
622 const bool useReplicaExchange
= (replExParams
.exchangeInterval
> 0);
623 if (useReplicaExchange
&& MASTER(cr
))
625 repl_ex
= init_replica_exchange(fplog
, cr
->ms
, top_global
->natoms
, ir
,
629 /* PME tuning is only supported in the Verlet scheme, with PME for
630 * Coulomb. It is not supported with only LJ PME, or for
632 bPMETune
= (mdrunOptions
.tunePme
&& EEL_PME(fr
->ic
->eeltype
) && !bRerunMD
&&
633 !mdrunOptions
.reproducible
&& ir
->cutoff_scheme
!= ecutsGROUP
);
636 pme_loadbal_init(&pme_loadbal
, cr
, mdlog
, ir
, state
->box
,
637 fr
->ic
, fr
->nbv
->listParams
.get(), fr
->pmedata
, use_GPU(fr
->nbv
),
641 if (!ir
->bContinuation
&& !bRerunMD
)
643 if (state
->flags
& (1 << estV
))
645 /* Set the velocities of vsites, shells and frozen atoms to zero */
646 for (i
= 0; i
< mdatoms
->homenr
; i
++)
648 if (mdatoms
->ptype
[i
] == eptVSite
||
649 mdatoms
->ptype
[i
] == eptShell
)
651 clear_rvec(state
->v
[i
]);
653 else if (mdatoms
->cFREEZE
)
655 for (m
= 0; m
< DIM
; m
++)
657 if (ir
->opts
.nFreeze
[mdatoms
->cFREEZE
[i
]][m
])
668 /* Constrain the initial coordinates and velocities */
669 do_constrain_first(fplog
, constr
, ir
, mdatoms
, state
,
674 /* Construct the virtual sites for the initial configuration */
675 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, nullptr,
676 top
->idef
.iparams
, top
->idef
.il
,
677 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
681 if (ir
->efep
!= efepNO
)
683 /* Set free energy calculation frequency as the greatest common
684 * denominator of nstdhdl and repl_ex_nst. */
685 nstfep
= ir
->fepvals
->nstdhdl
;
688 nstfep
= gmx_greatest_common_divisor(ir
->expandedvals
->nstexpanded
, nstfep
);
690 if (useReplicaExchange
)
692 nstfep
= gmx_greatest_common_divisor(replExParams
.exchangeInterval
, nstfep
);
696 /* Be REALLY careful about what flags you set here. You CANNOT assume
697 * this is the first step, since we might be restarting from a checkpoint,
698 * and in that case we should not do any modifications to the state.
700 bStopCM
= (ir
->comm_mode
!= ecmNO
&& !ir
->bContinuation
);
702 if (continuationOptions
.haveReadEkin
)
704 restore_ekinstate_from_state(cr
, ekind
, &state_global
->ekinstate
);
707 cglo_flags
= (CGLO_TEMPERATURE
| CGLO_GSTAT
708 | (EI_VV(ir
->eI
) ? CGLO_PRESSURE
: 0)
709 | (EI_VV(ir
->eI
) ? CGLO_CONSTRAINT
: 0)
710 | (continuationOptions
.haveReadEkin
? CGLO_READEKIN
: 0));
712 bSumEkinhOld
= FALSE
;
713 /* To minimize communication, compute_globals computes the COM velocity
714 * and the kinetic energy for the velocities without COM motion removed.
715 * Thus to get the kinetic energy without the COM contribution, we need
716 * to call compute_globals twice.
718 for (int cgloIteration
= 0; cgloIteration
< (bStopCM
? 2 : 1); cgloIteration
++)
720 int cglo_flags_iteration
= cglo_flags
;
721 if (bStopCM
&& cgloIteration
== 0)
723 cglo_flags_iteration
|= CGLO_STOPCM
;
725 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
726 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
727 constr
, &nullSignaller
, state
->box
,
728 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags_iteration
729 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
731 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
732 top_global
, top
, state
,
733 &shouldCheckNumberOfBondedInteractions
);
734 if (ir
->eI
== eiVVAK
)
736 /* a second call to get the half step temperature initialized as well */
737 /* we do the same call as above, but turn the pressure off -- internally to
738 compute_globals, this is recognized as a velocity verlet half-step
739 kinetic energy calculation. This minimized excess variables, but
740 perhaps loses some logic?*/
742 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
743 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
744 constr
, &nullSignaller
, state
->box
,
745 nullptr, &bSumEkinhOld
,
746 cglo_flags
& ~CGLO_PRESSURE
);
749 /* Calculate the initial half step temperature, and save the ekinh_old */
750 if (!continuationOptions
.startedFromCheckpoint
)
752 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
754 copy_mat(ekind
->tcstat
[i
].ekinh
, ekind
->tcstat
[i
].ekinh_old
);
758 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
759 temperature control */
760 trotter_seq
= init_npt_vars(ir
, state
, &MassQ
, bTrotter
);
764 if (!ir
->bContinuation
)
766 if (constr
&& ir
->eConstrAlg
== econtLINCS
)
769 "RMS relative constraint deviation after constraining: %.2e\n",
770 constr_rmsd(constr
));
772 if (EI_STATE_VELOCITY(ir
->eI
))
774 real temp
= enerd
->term
[F_TEMP
];
777 /* Result of Ekin averaged over velocities of -half
778 * and +half step, while we only have -half step here.
782 fprintf(fplog
, "Initial temperature: %g K\n", temp
);
788 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
789 " input trajectory '%s'\n\n",
790 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
791 if (mdrunOptions
.verbose
)
793 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
794 "run input file,\nwhich may not correspond to the time "
795 "needed to process input trajectory.\n\n");
801 fprintf(stderr
, "starting mdrun '%s'\n",
802 *(top_global
->name
));
805 sprintf(tbuf
, "%8.1f", (ir
->init_step
+ir
->nsteps
)*ir
->delta_t
);
809 sprintf(tbuf
, "%s", "infinite");
811 if (ir
->init_step
> 0)
813 fprintf(stderr
, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
814 gmx_step_str(ir
->init_step
+ir
->nsteps
, sbuf
), tbuf
,
815 gmx_step_str(ir
->init_step
, sbuf2
),
816 ir
->init_step
*ir
->delta_t
);
820 fprintf(stderr
, "%s steps, %s ps.\n",
821 gmx_step_str(ir
->nsteps
, sbuf
), tbuf
);
824 fprintf(fplog
, "\n");
827 walltime_accounting_start(walltime_accounting
);
828 wallcycle_start(wcycle
, ewcRUN
);
829 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
831 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
833 chkpt_ret
= fcCheckPointParallel( cr
->nodeid
,
837 gmx_fatal( 3, __FILE__
, __LINE__
, "Checkpoint error on step %d\n", 0 );
841 /***********************************************************
845 ************************************************************/
847 /* if rerunMD then read coordinates and velocities from input trajectory */
850 if (getenv("GMX_FORCE_UPDATE"))
858 bLastStep
= !read_first_frame(oenv
, &status
,
859 opt2fn("-rerun", nfile
, fnm
),
860 &rerun_fr
, TRX_NEED_X
| TRX_READ_V
);
861 if (rerun_fr
.natoms
!= top_global
->natoms
)
864 "Number of atoms in trajectory (%d) does not match the "
865 "run input file (%d)\n",
866 rerun_fr
.natoms
, top_global
->natoms
);
868 if (ir
->ePBC
!= epbcNONE
)
872 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
);
874 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < gmx::square(fr
->rlist
))
876 gmx_fatal(FARGS
, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
883 rerun_parallel_comm(cr
, &rerun_fr
, &bLastStep
);
886 if (ir
->ePBC
!= epbcNONE
)
888 /* Set the shift vectors.
889 * Necessary here when have a static box different from the tpr box.
891 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
895 /* Loop over MD steps or if rerunMD to end of input trajectory,
896 * or, if max_hours>0, until max_hours is reached.
898 real max_hours
= mdrunOptions
.maximumHoursToRun
;
900 /* Skip the first Nose-Hoover integration when we get the state from tpx */
901 bInitStep
= !startingFromCheckpoint
|| EI_VV(ir
->eI
);
902 bSumEkinhOld
= FALSE
;
904 bNeedRepartition
= FALSE
;
905 // TODO This implementation of ensemble orientation restraints is nasty because
906 // a user can't just do multi-sim with single-sim orientation restraints.
907 bUsingEnsembleRestraints
= (fcd
->disres
.nsystems
> 1) || (cr
->ms
&& fcd
->orires
.nr
);
910 // Replica exchange and ensemble restraints need all
911 // simulations to remain synchronized, so they need
912 // checkpoints and stop conditions to act on the same step, so
913 // the propagation of such signals must take place between
914 // simulations, not just within simulations.
915 bool checkpointIsLocal
= !useReplicaExchange
&& !bUsingEnsembleRestraints
;
916 bool stopConditionIsLocal
= !useReplicaExchange
&& !bUsingEnsembleRestraints
;
917 bool resetCountersIsLocal
= true;
918 signals
[eglsCHKPT
] = SimulationSignal(checkpointIsLocal
);
919 signals
[eglsSTOPCOND
] = SimulationSignal(stopConditionIsLocal
);
920 signals
[eglsRESETCOUNTERS
] = SimulationSignal(resetCountersIsLocal
);
923 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion
= (DOMAINDECOMP(cr
) ? DdOpenBalanceRegionBeforeForceComputation::yes
: DdOpenBalanceRegionBeforeForceComputation::no
);
924 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion
= (DOMAINDECOMP(cr
) ? DdCloseBalanceRegionAfterForceComputation::yes
: DdCloseBalanceRegionAfterForceComputation::no
);
926 step
= ir
->init_step
;
929 // TODO extract this to new multi-simulation module
930 if (MASTER(cr
) && MULTISIM(cr
) && !useReplicaExchange
)
932 if (!multisim_int_all_are_equal(cr
->ms
, ir
->nsteps
))
934 GMX_LOG(mdlog
.warning
).appendText(
935 "Note: The number of steps is not consistent across multi simulations,\n"
936 "but we are proceeding anyway!");
938 if (!multisim_int_all_are_equal(cr
->ms
, ir
->init_step
))
940 GMX_LOG(mdlog
.warning
).appendText(
941 "Note: The initial step is not consistent across multi simulations,\n"
942 "but we are proceeding anyway!");
946 /* and stop now if we should */
947 bLastStep
= (bLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
951 /* Determine if this is a neighbor search step */
952 bNStList
= (ir
->nstlist
> 0 && step
% ir
->nstlist
== 0);
954 if (bPMETune
&& bNStList
)
956 /* PME grid + cut-off optimization with GPUs or PME nodes */
957 pme_loadbal_do(pme_loadbal
, cr
,
958 (mdrunOptions
.verbose
&& MASTER(cr
)) ? stderr
: nullptr,
966 wallcycle_start(wcycle
, ewcSTEP
);
972 step
= rerun_fr
.step
;
973 step_rel
= step
- ir
->init_step
;
986 bLastStep
= (step_rel
== ir
->nsteps
);
987 t
= t0
+ step
*ir
->delta_t
;
990 // TODO Refactor this, so that nstfep does not need a default value of zero
991 if (ir
->efep
!= efepNO
|| ir
->bSimTemp
)
993 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
994 requiring different logic. */
999 setCurrentLambdasRerun(step
, ir
->fepvals
, &rerun_fr
, lam0
, state_global
);
1004 setCurrentLambdasLocal(step
, ir
->fepvals
, lam0
, state
);
1006 bDoDHDL
= do_per_step(step
, ir
->fepvals
->nstdhdl
);
1007 bDoFEP
= ((ir
->efep
!= efepNO
) && do_per_step(step
, nstfep
));
1008 bDoExpanded
= (do_per_step(step
, ir
->expandedvals
->nstexpanded
)
1009 && (ir
->bExpanded
) && (step
> 0) && (!startingFromCheckpoint
));
1012 bDoReplEx
= (useReplicaExchange
&& (step
> 0) && !bLastStep
&&
1013 do_per_step(step
, replExParams
.exchangeInterval
));
1017 update_annealing_target_temp(ir
, t
, upd
);
1020 if (bRerunMD
&& MASTER(cr
))
1022 const bool constructVsites
= (vsite
&& mdrunOptions
.rerunConstructVsites
);
1023 if (constructVsites
&& DOMAINDECOMP(cr
))
1025 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
1027 prepareRerunState(rerun_fr
, state_global
, constructVsites
, vsite
, top
->idef
, ir
->delta_t
, *fr
, graph
, &bRerunWarnNoV
);
1030 /* Stop Center of Mass motion */
1031 bStopCM
= (ir
->comm_mode
!= ecmNO
&& do_per_step(step
, ir
->nstcomm
));
1035 /* for rerun MD always do Neighbour Searching */
1036 bNS
= (bFirstStep
|| ir
->nstlist
!= 0);
1040 /* Determine whether or not to do Neighbour Searching */
1041 bNS
= (bFirstStep
|| bNStList
|| bExchanged
|| bNeedRepartition
);
1044 /* < 0 means stop at next step, > 0 means stop at next NS step */
1045 if ( (signals
[eglsSTOPCOND
].set
< 0) ||
1046 ( (signals
[eglsSTOPCOND
].set
> 0 ) && ( bNS
|| ir
->nstlist
== 0)))
1051 /* Determine whether or not to update the Born radii if doing GB */
1052 bBornRadii
= bFirstStep
;
1053 if (ir
->implicit_solvent
&& (step
% ir
->nstgbradii
== 0))
1058 /* do_log triggers energy and virial calculation. Because this leads
1059 * to different code paths, forces can be different. Thus for exact
1060 * continuation we should avoid extra log output.
1061 * Note that the || bLastStep can result in non-exact continuation
1062 * beyond the last step. But we don't consider that to be an issue.
1064 do_log
= do_per_step(step
, ir
->nstlog
) || (bFirstStep
&& !startingFromCheckpoint
) || bLastStep
|| bRerunMD
;
1065 do_verbose
= mdrunOptions
.verbose
&&
1066 (step
% mdrunOptions
.verboseStepPrintInterval
== 0 || bFirstStep
|| bLastStep
|| bRerunMD
);
1068 if (bNS
&& !(bFirstStep
&& ir
->bContinuation
&& !bRerunMD
))
1072 bMasterState
= TRUE
;
1076 bMasterState
= FALSE
;
1077 /* Correct the new box if it is too skewed */
1078 if (inputrecDynamicBox(ir
))
1080 if (correct_box(fplog
, step
, state
->box
, graph
))
1082 bMasterState
= TRUE
;
1085 if (DOMAINDECOMP(cr
) && bMasterState
)
1087 dd_collect_state(cr
->dd
, state
, state_global
);
1091 if (DOMAINDECOMP(cr
))
1093 /* Repartition the domain decomposition */
1094 dd_partition_system(fplog
, step
, cr
,
1095 bMasterState
, nstglobalcomm
,
1096 state_global
, top_global
, ir
,
1097 state
, &f
, mdAtoms
, top
, fr
,
1100 do_verbose
&& !bPMETunePrinting
);
1101 shouldCheckNumberOfBondedInteractions
= true;
1102 update_realloc(upd
, state
->natoms
);
1106 if (MASTER(cr
) && do_log
)
1108 print_ebin_header(fplog
, step
, t
); /* can we improve the information printed here? */
1111 if (ir
->efep
!= efepNO
)
1113 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
1116 if ((bRerunMD
&& rerun_fr
.bV
) || bExchanged
)
1119 /* We need the kinetic energy at minus the half step for determining
1120 * the full step kinetic energy and possibly for T-coupling.*/
1121 /* This may not be quite working correctly yet . . . . */
1122 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1123 wcycle
, enerd
, nullptr, nullptr, nullptr, nullptr, mu_tot
,
1124 constr
, &nullSignaller
, state
->box
,
1125 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1126 CGLO_GSTAT
| CGLO_TEMPERATURE
| CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
);
1127 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
1128 top_global
, top
, state
,
1129 &shouldCheckNumberOfBondedInteractions
);
1131 clear_mat(force_vir
);
1133 /* We write a checkpoint at this MD step when:
1134 * either at an NS step when we signalled through gs,
1135 * or at the last step (but not when we do not want confout),
1136 * but never at the first step or with rerun.
1138 bCPT
= (((signals
[eglsCHKPT
].set
&& (bNS
|| ir
->nstlist
== 0)) ||
1139 (bLastStep
&& mdrunOptions
.writeConfout
)) &&
1140 step
> ir
->init_step
&& !bRerunMD
);
1143 signals
[eglsCHKPT
].set
= 0;
1146 /* Determine the energy and pressure:
1147 * at nstcalcenergy steps and at energy output steps (set below).
1149 if (EI_VV(ir
->eI
) && (!bInitStep
))
1151 /* for vv, the first half of the integration actually corresponds
1152 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1153 but the virial needs to be calculated on both the current step and the 'next' step. Future
1154 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1156 /* TODO: This is probably not what we want, we will write to energy file one step after nstcalcenergy steps. */
1157 bCalcEnerStep
= do_per_step(step
- 1, ir
->nstcalcenergy
);
1158 bCalcVir
= bCalcEnerStep
||
1159 (ir
->epc
!= epcNO
&& (do_per_step(step
, ir
->nstpcouple
) || do_per_step(step
-1, ir
->nstpcouple
)));
1163 bCalcEnerStep
= do_per_step(step
, ir
->nstcalcenergy
);
1164 bCalcVir
= bCalcEnerStep
||
1165 (ir
->epc
!= epcNO
&& do_per_step(step
, ir
->nstpcouple
));
1167 bCalcEner
= bCalcEnerStep
;
1169 do_ene
= (do_per_step(step
, ir
->nstenergy
) || bLastStep
|| bRerunMD
);
1171 if (do_ene
|| do_log
|| bDoReplEx
)
1177 /* Do we need global communication ? */
1178 bGStat
= (bCalcVir
|| bCalcEner
|| bStopCM
||
1179 do_per_step(step
, nstglobalcomm
) ||
1180 (EI_VV(ir
->eI
) && inputrecNvtTrotter(ir
) && do_per_step(step
-1, nstglobalcomm
)));
1182 force_flags
= (GMX_FORCE_STATECHANGED
|
1183 ((inputrecDynamicBox(ir
) || bRerunMD
) ? GMX_FORCE_DYNAMICBOX
: 0) |
1184 GMX_FORCE_ALLFORCES
|
1185 (bCalcVir
? GMX_FORCE_VIRIAL
: 0) |
1186 (bCalcEner
? GMX_FORCE_ENERGY
: 0) |
1187 (bDoFEP
? GMX_FORCE_DHDL
: 0)
1192 /* Now is the time to relax the shells */
1193 relax_shell_flexcon(fplog
, cr
, mdrunOptions
.verbose
, step
,
1194 ir
, bNS
, force_flags
, top
,
1196 state
, &f
, force_vir
, mdatoms
,
1197 nrnb
, wcycle
, graph
, groups
,
1198 shellfc
, fr
, bBornRadii
, t
, mu_tot
,
1200 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
1204 /* The AWH history need to be saved _before_ doing force calculations where the AWH bias is updated
1205 (or the AWH update will be performed twice for one step when continuing). It would be best to
1206 call this update function from do_md_trajectory_writing but that would occur after do_force.
1207 One would have to divide the update_awh function into one function applying the AWH force
1208 and one doing the AWH bias update. The update AWH bias function could then be called after
1209 do_md_trajectory_writing (then containing update_awh_history).
1210 The checkpointing will in the future probably moved to the start of the md loop which will
1211 rid of this issue. */
1212 if (ir
->bDoAwh
&& bCPT
&& MASTER(cr
))
1214 ir
->awh
->updateHistory(state_global
->awhHistory
.get());
1217 /* The coordinates (x) are shifted (to get whole molecules)
1219 * This is parallellized as well, and does communication too.
1220 * Check comments in sim_util.c
1222 do_force(fplog
, cr
, ir
, step
, nrnb
, wcycle
, top
, groups
,
1223 state
->box
, state
->x
, &state
->hist
,
1224 f
, force_vir
, mdatoms
, enerd
, fcd
,
1225 state
->lambda
, graph
,
1226 fr
, vsite
, mu_tot
, t
, ed
, bBornRadii
,
1227 (bNS
? GMX_FORCE_NS
: 0) | force_flags
,
1228 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
1231 if (EI_VV(ir
->eI
) && !startingFromCheckpoint
&& !bRerunMD
)
1232 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1234 rvec
*vbuf
= nullptr;
1236 wallcycle_start(wcycle
, ewcUPDATE
);
1237 if (ir
->eI
== eiVV
&& bInitStep
)
1239 /* if using velocity verlet with full time step Ekin,
1240 * take the first half step only to compute the
1241 * virial for the first step. From there,
1242 * revert back to the initial coordinates
1243 * so that the input is actually the initial step.
1245 snew(vbuf
, state
->natoms
);
1246 copy_rvecn(as_rvec_array(state
->v
.data()), vbuf
, 0, state
->natoms
); /* should make this better for parallelizing? */
1250 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1251 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ1
);
1254 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1255 ekind
, M
, upd
, etrtVELOCITY1
,
1258 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
) /* Why is rerun_fr.bV here? Unclear. */
1260 wallcycle_stop(wcycle
, ewcUPDATE
);
1261 update_constraints(fplog
, step
, nullptr, ir
, mdatoms
,
1262 state
, fr
->bMolPBC
, graph
, f
,
1263 &top
->idef
, shake_vir
,
1264 cr
, nrnb
, wcycle
, upd
, constr
,
1266 wallcycle_start(wcycle
, ewcUPDATE
);
1270 /* Need to unshift here if a do_force has been
1271 called in the previous step */
1272 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1274 /* if VV, compute the pressure and constraints */
1275 /* For VV2, we strictly only need this if using pressure
1276 * control, but we really would like to have accurate pressures
1278 * Think about ways around this in the future?
1279 * For now, keep this choice in comments.
1281 /*bPres = (ir->eI==eiVV || inputrecNptTrotter(ir)); */
1282 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && inputrecNptTrotter(ir)));*/
1284 bTemp
= ((ir
->eI
== eiVV
&& (!bInitStep
)) || (ir
->eI
== eiVVAK
));
1285 if (bCalcEner
&& ir
->eI
== eiVVAK
)
1287 bSumEkinhOld
= TRUE
;
1289 /* for vv, the first half of the integration actually corresponds to the previous step.
1290 So we need information from the last step in the first half of the integration */
1291 if (bGStat
|| do_per_step(step
-1, nstglobalcomm
))
1293 wallcycle_stop(wcycle
, ewcUPDATE
);
1294 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1295 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1296 constr
, &nullSignaller
, state
->box
,
1297 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1298 (bGStat
? CGLO_GSTAT
: 0)
1300 | (bTemp
? CGLO_TEMPERATURE
: 0)
1301 | (bPres
? CGLO_PRESSURE
: 0)
1302 | (bPres
? CGLO_CONSTRAINT
: 0)
1303 | (bStopCM
? CGLO_STOPCM
: 0)
1304 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
1307 /* explanation of above:
1308 a) We compute Ekin at the full time step
1309 if 1) we are using the AveVel Ekin, and it's not the
1310 initial step, or 2) if we are using AveEkin, but need the full
1311 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1312 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1313 EkinAveVel because it's needed for the pressure */
1314 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
1315 top_global
, top
, state
,
1316 &shouldCheckNumberOfBondedInteractions
);
1317 wallcycle_start(wcycle
, ewcUPDATE
);
1319 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1324 m_add(force_vir
, shake_vir
, total_vir
); /* we need the un-dispersion corrected total vir here */
1325 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ2
);
1327 /* TODO This is only needed when we're about to write
1328 * a checkpoint, because we use it after the restart
1329 * (in a kludge?). But what should we be doing if
1330 * startingFromCheckpoint or bInitStep are true? */
1331 if (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
))
1333 copy_mat(shake_vir
, state
->svir_prev
);
1334 copy_mat(force_vir
, state
->fvir_prev
);
1336 if (inputrecNvtTrotter(ir
) && ir
->eI
== eiVV
)
1338 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1339 enerd
->term
[F_TEMP
] = sum_ekin(&(ir
->opts
), ekind
, nullptr, (ir
->eI
== eiVV
), FALSE
);
1340 enerd
->term
[F_EKIN
] = trace(ekind
->ekin
);
1343 else if (bExchanged
)
1345 wallcycle_stop(wcycle
, ewcUPDATE
);
1346 /* We need the kinetic energy at minus the half step for determining
1347 * the full step kinetic energy and possibly for T-coupling.*/
1348 /* This may not be quite working correctly yet . . . . */
1349 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1350 wcycle
, enerd
, nullptr, nullptr, nullptr, nullptr, mu_tot
,
1351 constr
, &nullSignaller
, state
->box
,
1352 nullptr, &bSumEkinhOld
,
1353 CGLO_GSTAT
| CGLO_TEMPERATURE
);
1354 wallcycle_start(wcycle
, ewcUPDATE
);
1357 /* if it's the initial step, we performed this first step just to get the constraint virial */
1358 if (ir
->eI
== eiVV
&& bInitStep
)
1360 copy_rvecn(vbuf
, as_rvec_array(state
->v
.data()), 0, state
->natoms
);
1363 wallcycle_stop(wcycle
, ewcUPDATE
);
1366 /* compute the conserved quantity */
1369 saved_conserved_quantity
= NPT_energy(ir
, state
, &MassQ
);
1372 last_ekin
= enerd
->term
[F_EKIN
];
1374 if ((ir
->eDispCorr
!= edispcEnerPres
) && (ir
->eDispCorr
!= edispcAllEnerPres
))
1376 saved_conserved_quantity
-= enerd
->term
[F_DISPCORR
];
1378 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1379 if (ir
->efep
!= efepNO
&& !bRerunMD
)
1381 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1385 /* ######## END FIRST UPDATE STEP ############## */
1386 /* ######## If doing VV, we now have v(dt) ###### */
1389 /* perform extended ensemble sampling in lambda - we don't
1390 actually move to the new state before outputting
1391 statistics, but if performing simulated tempering, we
1392 do update the velocities and the tau_t. */
1394 lamnew
= ExpandedEnsembleDynamics(fplog
, ir
, enerd
, state
, &MassQ
, state
->fep_state
, state
->dfhist
, step
, as_rvec_array(state
->v
.data()), mdatoms
);
1395 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1398 copy_df_history(state_global
->dfhist
, state
->dfhist
);
1402 /* Now we have the energies and forces corresponding to the
1403 * coordinates at time t. We must output all of this before
1406 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
1407 ir
, state
, state_global
, observablesHistory
,
1409 outf
, mdebin
, ekind
, f
,
1411 bCPT
, bRerunMD
, bLastStep
,
1412 mdrunOptions
.writeConfout
,
1414 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1415 bIMDstep
= do_IMD(ir
->bIMD
, step
, cr
, bNS
, state
->box
, as_rvec_array(state
->x
.data()), ir
, t
, wcycle
);
1417 /* kludge -- virial is lost with restart for MTTK NPT control. Must reload (saved earlier). */
1418 if (startingFromCheckpoint
&& (inputrecNptTrotter(ir
) || inputrecNphTrotter(ir
)))
1420 copy_mat(state
->svir_prev
, shake_vir
);
1421 copy_mat(state
->fvir_prev
, force_vir
);
1424 elapsed_time
= walltime_accounting_get_current_elapsed_time(walltime_accounting
);
1426 /* Check whether everything is still allright */
1427 if (((int)gmx_get_stop_condition() > handled_stop_condition
)
1433 int nsteps_stop
= -1;
1435 /* this just makes signals[].sig compatible with the hack
1436 of sending signals around by MPI_Reduce together with
1438 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns
)
1440 signals
[eglsSTOPCOND
].sig
= 1;
1441 nsteps_stop
= std::max(ir
->nstlist
, 2*nstglobalcomm
);
1443 if (gmx_get_stop_condition() == gmx_stop_cond_next
)
1445 signals
[eglsSTOPCOND
].sig
= -1;
1446 nsteps_stop
= nstglobalcomm
+ 1;
1451 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1452 gmx_get_signal_name(), nsteps_stop
);
1456 "\n\nReceived the %s signal, stopping within %d steps\n\n",
1457 gmx_get_signal_name(), nsteps_stop
);
1459 handled_stop_condition
= (int)gmx_get_stop_condition();
1461 else if (MASTER(cr
) && (bNS
|| ir
->nstlist
<= 0) &&
1462 (max_hours
> 0 && elapsed_time
> max_hours
*60.0*60.0*0.99) &&
1463 signals
[eglsSTOPCOND
].sig
== 0 && signals
[eglsSTOPCOND
].set
== 0)
1465 /* Signal to terminate the run */
1466 signals
[eglsSTOPCOND
].sig
= 1;
1469 fprintf(fplog
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1471 fprintf(stderr
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1474 if (bResetCountersHalfMaxH
&& MASTER(cr
) &&
1475 elapsed_time
> max_hours
*60.0*60.0*0.495)
1477 /* Set flag that will communicate the signal to all ranks in the simulation */
1478 signals
[eglsRESETCOUNTERS
].sig
= 1;
1481 /* In parallel we only have to check for checkpointing in steps
1482 * where we do global communication,
1483 * otherwise the other nodes don't know.
1485 const real cpt_period
= mdrunOptions
.checkpointOptions
.period
;
1486 if (MASTER(cr
) && ((bGStat
|| !PAR(cr
)) &&
1489 elapsed_time
>= nchkpt
*cpt_period
*60.0)) &&
1490 signals
[eglsCHKPT
].set
== 0)
1492 signals
[eglsCHKPT
].sig
= 1;
1495 /* ######### START SECOND UPDATE STEP ################# */
1497 /* at the start of step, randomize or scale the velocities ((if vv. Restriction of Andersen controlled
1500 if (ETC_ANDERSEN(ir
->etc
)) /* keep this outside of update_tcouple because of the extra info required to pass */
1502 gmx_bool bIfRandomize
;
1503 bIfRandomize
= update_randomize_velocities(ir
, step
, cr
, mdatoms
, state
, upd
, constr
);
1504 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1505 if (constr
&& bIfRandomize
)
1507 update_constraints(fplog
, step
, nullptr, ir
, mdatoms
,
1508 state
, fr
->bMolPBC
, graph
, f
,
1509 &top
->idef
, tmp_vir
,
1510 cr
, nrnb
, wcycle
, upd
, constr
,
1514 /* Box is changed in update() when we do pressure coupling,
1515 * but we should still use the old box for energy corrections and when
1516 * writing it to the energy file, so it matches the trajectory files for
1517 * the same timestep above. Make a copy in a separate array.
1519 copy_mat(state
->box
, lastbox
);
1523 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
)
1525 wallcycle_start(wcycle
, ewcUPDATE
);
1526 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1529 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ3
);
1530 /* We can only do Berendsen coupling after we have summed
1531 * the kinetic energy or virial. Since the happens
1532 * in global_state after update, we should only do it at
1533 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1538 update_tcouple(step
, ir
, state
, ekind
, &MassQ
, mdatoms
);
1539 update_pcouple_before_coordinates(fplog
, step
, ir
, state
,
1540 parrinellorahmanMu
, M
,
1546 /* velocity half-step update */
1547 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1548 ekind
, M
, upd
, etrtVELOCITY2
,
1552 /* Above, initialize just copies ekinh into ekin,
1553 * it doesn't copy position (for VV),
1554 * and entire integrator for MD.
1557 if (ir
->eI
== eiVVAK
)
1559 /* We probably only need md->homenr, not state->natoms */
1560 if (state
->natoms
> cbuf_nalloc
)
1562 cbuf_nalloc
= state
->natoms
;
1563 srenew(cbuf
, cbuf_nalloc
);
1565 copy_rvecn(as_rvec_array(state
->x
.data()), cbuf
, 0, state
->natoms
);
1568 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1569 ekind
, M
, upd
, etrtPOSITION
, cr
, constr
);
1570 wallcycle_stop(wcycle
, ewcUPDATE
);
1572 update_constraints(fplog
, step
, &dvdl_constr
, ir
, mdatoms
, state
,
1573 fr
->bMolPBC
, graph
, f
,
1574 &top
->idef
, shake_vir
,
1575 cr
, nrnb
, wcycle
, upd
, constr
,
1578 if (ir
->eI
== eiVVAK
)
1580 /* erase F_EKIN and F_TEMP here? */
1581 /* just compute the kinetic energy at the half step to perform a trotter step */
1582 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1583 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1584 constr
, &nullSignaller
, lastbox
,
1585 nullptr, &bSumEkinhOld
,
1586 (bGStat
? CGLO_GSTAT
: 0) | CGLO_TEMPERATURE
1588 wallcycle_start(wcycle
, ewcUPDATE
);
1589 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ4
);
1590 /* now we know the scaling, we can compute the positions again again */
1591 copy_rvecn(cbuf
, as_rvec_array(state
->x
.data()), 0, state
->natoms
);
1593 update_coords(fplog
, step
, ir
, mdatoms
, state
, f
, fcd
,
1594 ekind
, M
, upd
, etrtPOSITION
, cr
, constr
);
1595 wallcycle_stop(wcycle
, ewcUPDATE
);
1597 /* do we need an extra constraint here? just need to copy out of as_rvec_array(state->v.data()) to upd->xp? */
1598 /* are the small terms in the shake_vir here due
1599 * to numerical errors, or are they important
1600 * physically? I'm thinking they are just errors, but not completely sure.
1601 * For now, will call without actually constraining, constr=NULL*/
1602 update_constraints(fplog
, step
, nullptr, ir
, mdatoms
,
1603 state
, fr
->bMolPBC
, graph
, f
,
1604 &top
->idef
, tmp_vir
,
1605 cr
, nrnb
, wcycle
, upd
, nullptr,
1610 /* this factor or 2 correction is necessary
1611 because half of the constraint force is removed
1612 in the vv step, so we have to double it. See
1613 the Redmine issue #1255. It is not yet clear
1614 if the factor of 2 is exact, or just a very
1615 good approximation, and this will be
1616 investigated. The next step is to see if this
1617 can be done adding a dhdl contribution from the
1618 rattle step, but this is somewhat more
1619 complicated with the current code. Will be
1620 investigated, hopefully for 4.6.3. However,
1621 this current solution is much better than
1622 having it completely wrong.
1624 enerd
->term
[F_DVDL_CONSTR
] += 2*dvdl_constr
;
1628 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
1633 /* Need to unshift here */
1634 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1637 if (vsite
!= nullptr)
1639 wallcycle_start(wcycle
, ewcVSITECONSTR
);
1640 if (graph
!= nullptr)
1642 shift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1644 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, as_rvec_array(state
->v
.data()),
1645 top
->idef
.iparams
, top
->idef
.il
,
1646 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
1648 if (graph
!= nullptr)
1650 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
1652 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
1655 /* ############## IF NOT VV, Calculate globals HERE ############ */
1656 /* With Leap-Frog we can skip compute_globals at
1657 * non-communication steps, but we need to calculate
1658 * the kinetic energy one step before communication.
1661 // Organize to do inter-simulation signalling on steps if
1662 // and when algorithms require it.
1663 bool doInterSimSignal
= (!bFirstStep
&& bDoReplEx
) || bUsingEnsembleRestraints
;
1665 if (bGStat
|| (!EI_VV(ir
->eI
) && do_per_step(step
+1, nstglobalcomm
)) || doInterSimSignal
)
1667 // Since we're already communicating at this step, we
1668 // can propagate intra-simulation signals. Note that
1669 // check_nstglobalcomm has the responsibility for
1670 // choosing the value of nstglobalcomm that is one way
1671 // bGStat becomes true, so we can't get into a
1672 // situation where e.g. checkpointing can't be
1674 bool doIntraSimSignal
= true;
1675 SimulationSignaller
signaller(&signals
, cr
, doInterSimSignal
, doIntraSimSignal
);
1677 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
1678 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1681 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
1682 (bGStat
? CGLO_GSTAT
: 0)
1683 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_ENERGY
: 0)
1684 | (!EI_VV(ir
->eI
) && bStopCM
? CGLO_STOPCM
: 0)
1685 | (!EI_VV(ir
->eI
) ? CGLO_TEMPERATURE
: 0)
1686 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_PRESSURE
: 0)
1688 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
1690 checkNumberOfBondedInteractions(fplog
, cr
, totalNumberOfBondedInteractions
,
1691 top_global
, top
, state
,
1692 &shouldCheckNumberOfBondedInteractions
);
1696 /* ############# END CALC EKIN AND PRESSURE ################# */
1698 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1699 the virial that should probably be addressed eventually. state->veta has better properies,
1700 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1701 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1703 if (ir
->efep
!= efepNO
&& (!EI_VV(ir
->eI
) || bRerunMD
))
1705 /* Sum up the foreign energy and dhdl terms for md and sd.
1706 Currently done every step so that dhdl is correct in the .edr */
1707 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1710 update_pcouple_after_coordinates(fplog
, step
, ir
, mdatoms
,
1711 pres
, force_vir
, shake_vir
,
1715 /* ################# END UPDATE STEP 2 ################# */
1716 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1718 /* The coordinates (x) were unshifted in update */
1721 /* We will not sum ekinh_old,
1722 * so signal that we still have to do it.
1724 bSumEkinhOld
= TRUE
;
1729 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1731 /* use the directly determined last velocity, not actually the averaged half steps */
1732 if (bTrotter
&& ir
->eI
== eiVV
)
1734 enerd
->term
[F_EKIN
] = last_ekin
;
1736 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1738 if (integratorHasConservedEnergyQuantity(ir
))
1742 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + saved_conserved_quantity
;
1746 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + NPT_energy(ir
, state
, &MassQ
);
1749 /* ######### END PREPARING EDR OUTPUT ########### */
1755 if (fplog
&& do_log
&& bDoExpanded
)
1757 /* only needed if doing expanded ensemble */
1758 PrintFreeEnergyInfoToFile(fplog
, ir
->fepvals
, ir
->expandedvals
, ir
->bSimTemp
? ir
->simtempvals
: nullptr,
1759 state_global
->dfhist
, state
->fep_state
, ir
->nstlog
, step
);
1763 upd_mdebin(mdebin
, bDoDHDL
, bCalcEnerStep
,
1764 t
, mdatoms
->tmass
, enerd
, state
,
1765 ir
->fepvals
, ir
->expandedvals
, lastbox
,
1766 shake_vir
, force_vir
, total_vir
, pres
,
1767 ekind
, mu_tot
, constr
);
1771 upd_mdebin_step(mdebin
);
1774 gmx_bool do_dr
= do_per_step(step
, ir
->nstdisreout
);
1775 gmx_bool do_or
= do_per_step(step
, ir
->nstorireout
);
1777 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
, do_log
? fplog
: nullptr,
1779 eprNORMAL
, mdebin
, fcd
, groups
, &(ir
->opts
), ir
->awh
);
1783 pull_print_output(ir
->pull_work
, step
, t
);
1786 if (do_per_step(step
, ir
->nstlog
))
1788 if (fflush(fplog
) != 0)
1790 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
1796 /* Have to do this part _after_ outputting the logfile and the edr file */
1797 /* Gets written into the state at the beginning of next loop*/
1798 state
->fep_state
= lamnew
;
1800 /* Print the remaining wall clock time for the run */
1801 if (MULTIMASTER(cr
) &&
1802 (do_verbose
|| gmx_got_usr_signal()) &&
1807 fprintf(stderr
, "\n");
1809 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
1812 /* Ion/water position swapping.
1813 * Not done in last step since trajectory writing happens before this call
1814 * in the MD loop and exchanges would be lost anyway. */
1815 bNeedRepartition
= FALSE
;
1816 if ((ir
->eSwapCoords
!= eswapNO
) && (step
> 0) && !bLastStep
&&
1817 do_per_step(step
, ir
->swap
->nstswap
))
1819 bNeedRepartition
= do_swapcoords(cr
, step
, t
, ir
, wcycle
,
1820 bRerunMD
? rerun_fr
.x
: as_rvec_array(state
->x
.data()),
1821 bRerunMD
? rerun_fr
.box
: state
->box
,
1822 MASTER(cr
) && mdrunOptions
.verbose
,
1825 if (bNeedRepartition
&& DOMAINDECOMP(cr
))
1827 dd_collect_state(cr
->dd
, state
, state_global
);
1831 /* Replica exchange */
1835 bExchanged
= replica_exchange(fplog
, cr
, repl_ex
,
1836 state_global
, enerd
,
1840 if ( (bExchanged
|| bNeedRepartition
) && DOMAINDECOMP(cr
) )
1842 dd_partition_system(fplog
, step
, cr
, TRUE
, 1,
1843 state_global
, top_global
, ir
,
1844 state
, &f
, mdAtoms
, top
, fr
,
1846 nrnb
, wcycle
, FALSE
);
1847 shouldCheckNumberOfBondedInteractions
= true;
1848 update_realloc(upd
, state
->natoms
);
1853 startingFromCheckpoint
= false;
1855 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1856 /* With all integrators, except VV, we need to retain the pressure
1857 * at the current step for coupling at the next step.
1859 if ((state
->flags
& (1<<estPRES_PREV
)) &&
1861 (ir
->nstpcouple
> 0 && step
% ir
->nstpcouple
== 0)))
1863 /* Store the pressure in t_state for pressure coupling
1864 * at the next MD step.
1866 copy_mat(pres
, state
->pres_prev
);
1869 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1871 if ( (membed
!= nullptr) && (!bLastStep
) )
1873 rescale_membed(step_rel
, membed
, as_rvec_array(state_global
->x
.data()));
1880 /* read next frame from input trajectory */
1881 bLastStep
= !read_next_frame(oenv
, status
, &rerun_fr
);
1886 rerun_parallel_comm(cr
, &rerun_fr
, &bLastStep
);
1890 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
1891 if (DOMAINDECOMP(cr
) && wcycle
)
1893 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
1896 if (!bRerunMD
|| !rerun_fr
.bStep
)
1898 /* increase the MD step number */
1903 /* TODO make a counter-reset module */
1904 /* If it is time to reset counters, set a flag that remains
1905 true until counters actually get reset */
1906 if (step_rel
== wcycle_get_reset_counters(wcycle
) ||
1907 signals
[eglsRESETCOUNTERS
].set
!= 0)
1909 if (pme_loadbal_is_active(pme_loadbal
))
1911 /* Do not permit counter reset while PME load
1912 * balancing is active. The only purpose for resetting
1913 * counters is to measure reliable performance data,
1914 * and that can't be done before balancing
1917 * TODO consider fixing this by delaying the reset
1918 * until after load balancing completes,
1919 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1920 gmx_fatal(FARGS
, "PME tuning was still active when attempting to "
1921 "reset mdrun counters at step %" GMX_PRId64
". Try "
1922 "resetting counters later in the run, e.g. with gmx "
1923 "mdrun -resetstep.", step
);
1925 reset_all_counters(fplog
, mdlog
, cr
, step
, &step_rel
, ir
, wcycle
, nrnb
, walltime_accounting
,
1926 use_GPU(fr
->nbv
) ? fr
->nbv
: nullptr, fr
->pmedata
);
1927 wcycle_set_reset_counters(wcycle
, -1);
1928 if (!thisRankHasDuty(cr
, DUTY_PME
))
1930 /* Tell our PME node to reset its counters */
1931 gmx_pme_send_resetcounters(cr
, step
);
1933 /* Correct max_hours for the elapsed time */
1934 max_hours
-= elapsed_time
/(60.0*60.0);
1935 /* If mdrun -maxh -resethway was active, it can only trigger once */
1936 bResetCountersHalfMaxH
= FALSE
; /* TODO move this to where signals[eglsRESETCOUNTERS].sig is set */
1937 /* Reset can only happen once, so clear the triggering flag. */
1938 signals
[eglsRESETCOUNTERS
].set
= 0;
1941 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1942 IMD_prep_energies_send_positions(ir
->bIMD
&& MASTER(cr
), bIMDstep
, ir
->imd
, enerd
, step
, bCalcEner
, wcycle
);
1945 /* End of main MD loop */
1947 /* Closing TNG files can include compressing data. Therefore it is good to do that
1948 * before stopping the time measurements. */
1949 mdoutf_tng_close(outf
);
1951 /* Stop measuring walltime */
1952 walltime_accounting_end(walltime_accounting
);
1954 if (bRerunMD
&& MASTER(cr
))
1959 if (!thisRankHasDuty(cr
, DUTY_PME
))
1961 /* Tell the PME only node to finish */
1962 gmx_pme_send_finish(cr
);
1967 if (ir
->nstcalcenergy
> 0 && !bRerunMD
)
1969 print_ebin(mdoutf_get_fp_ene(outf
), FALSE
, FALSE
, FALSE
, fplog
, step
, t
,
1970 eprAVER
, mdebin
, fcd
, groups
, &(ir
->opts
), ir
->awh
);
1978 pme_loadbal_done(pme_loadbal
, fplog
, mdlog
, use_GPU(fr
->nbv
));
1981 done_shellfc(fplog
, shellfc
, step_rel
);
1983 if (useReplicaExchange
&& MASTER(cr
))
1985 print_replica_exchange_statistics(fplog
, repl_ex
);
1993 // Clean up swapcoords
1994 if (ir
->eSwapCoords
!= eswapNO
)
1996 finish_swapcoords(ir
->swap
);
1999 /* Do essential dynamics cleanup if needed. Close .edo file */
2002 /* IMD cleanup, if bIMD is TRUE. */
2003 IMD_finalize(ir
->bIMD
, ir
->imd
);
2005 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);
2006 if (step_rel
>= wcycle_get_reset_counters(wcycle
) &&
2007 signals
[eglsRESETCOUNTERS
].set
== 0 &&
2008 !bResetCountersHalfMaxH
)
2010 walltime_accounting_set_valid_finish(walltime_accounting
);