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-2019,2020, 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 MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/device_context.h"
77 #include "gromacs/gpu_utils/device_stream_manager.h"
78 #include "gromacs/hardware/cpuinfo.h"
79 #include "gromacs/hardware/detecthardware.h"
80 #include "gromacs/hardware/device_management.h"
81 #include "gromacs/hardware/printhardware.h"
82 #include "gromacs/imd/imd.h"
83 #include "gromacs/listed_forces/disre.h"
84 #include "gromacs/listed_forces/gpubonded.h"
85 #include "gromacs/listed_forces/listed_forces.h"
86 #include "gromacs/listed_forces/orires.h"
87 #include "gromacs/math/functions.h"
88 #include "gromacs/math/utilities.h"
89 #include "gromacs/math/vec.h"
90 #include "gromacs/mdlib/boxdeformation.h"
91 #include "gromacs/mdlib/broadcaststructs.h"
92 #include "gromacs/mdlib/calc_verletbuf.h"
93 #include "gromacs/mdlib/dispersioncorrection.h"
94 #include "gromacs/mdlib/enerdata_utils.h"
95 #include "gromacs/mdlib/force.h"
96 #include "gromacs/mdlib/forcerec.h"
97 #include "gromacs/mdlib/gmx_omp_nthreads.h"
98 #include "gromacs/mdlib/makeconstraints.h"
99 #include "gromacs/mdlib/md_support.h"
100 #include "gromacs/mdlib/mdatoms.h"
101 #include "gromacs/mdlib/sighandler.h"
102 #include "gromacs/mdlib/stophandler.h"
103 #include "gromacs/mdlib/tgroup.h"
104 #include "gromacs/mdlib/updategroups.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdrun/mdmodules.h"
107 #include "gromacs/mdrun/simulationcontext.h"
108 #include "gromacs/mdrunutility/handlerestart.h"
109 #include "gromacs/mdrunutility/logging.h"
110 #include "gromacs/mdrunutility/multisim.h"
111 #include "gromacs/mdrunutility/printtime.h"
112 #include "gromacs/mdrunutility/threadaffinity.h"
113 #include "gromacs/mdtypes/checkpointdata.h"
114 #include "gromacs/mdtypes/commrec.h"
115 #include "gromacs/mdtypes/enerdata.h"
116 #include "gromacs/mdtypes/fcdata.h"
117 #include "gromacs/mdtypes/forcerec.h"
118 #include "gromacs/mdtypes/group.h"
119 #include "gromacs/mdtypes/inputrec.h"
120 #include "gromacs/mdtypes/interaction_const.h"
121 #include "gromacs/mdtypes/md_enums.h"
122 #include "gromacs/mdtypes/mdatom.h"
123 #include "gromacs/mdtypes/mdrunoptions.h"
124 #include "gromacs/mdtypes/observableshistory.h"
125 #include "gromacs/mdtypes/simulation_workload.h"
126 #include "gromacs/mdtypes/state.h"
127 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
128 #include "gromacs/modularsimulator/modularsimulator.h"
129 #include "gromacs/nbnxm/gpu_data_mgmt.h"
130 #include "gromacs/nbnxm/nbnxm.h"
131 #include "gromacs/nbnxm/pairlist_tuning.h"
132 #include "gromacs/pbcutil/pbc.h"
133 #include "gromacs/pulling/output.h"
134 #include "gromacs/pulling/pull.h"
135 #include "gromacs/pulling/pull_rotation.h"
136 #include "gromacs/restraint/manager.h"
137 #include "gromacs/restraint/restraintmdmodule.h"
138 #include "gromacs/restraint/restraintpotential.h"
139 #include "gromacs/swap/swapcoords.h"
140 #include "gromacs/taskassignment/decidegpuusage.h"
141 #include "gromacs/taskassignment/decidesimulationworkload.h"
142 #include "gromacs/taskassignment/resourcedivision.h"
143 #include "gromacs/taskassignment/taskassignment.h"
144 #include "gromacs/taskassignment/usergpuids.h"
145 #include "gromacs/timing/gpu_timing.h"
146 #include "gromacs/timing/wallcycle.h"
147 #include "gromacs/timing/wallcyclereporting.h"
148 #include "gromacs/topology/mtop_util.h"
149 #include "gromacs/trajectory/trajectoryframe.h"
150 #include "gromacs/utility/basenetwork.h"
151 #include "gromacs/utility/cstringutil.h"
152 #include "gromacs/utility/exceptions.h"
153 #include "gromacs/utility/fatalerror.h"
154 #include "gromacs/utility/filestream.h"
155 #include "gromacs/utility/gmxassert.h"
156 #include "gromacs/utility/gmxmpi.h"
157 #include "gromacs/utility/keyvaluetree.h"
158 #include "gromacs/utility/logger.h"
159 #include "gromacs/utility/loggerbuilder.h"
160 #include "gromacs/utility/mdmodulenotification.h"
161 #include "gromacs/utility/physicalnodecommunicator.h"
162 #include "gromacs/utility/pleasecite.h"
163 #include "gromacs/utility/programcontext.h"
164 #include "gromacs/utility/smalloc.h"
165 #include "gromacs/utility/stringutil.h"
167 #include "isimulator.h"
168 #include "membedholder.h"
169 #include "replicaexchange.h"
170 #include "simulationinput.h"
171 #include "simulationinpututility.h"
172 #include "simulatorbuilder.h"
175 # include "corewrap.h"
182 /*! \brief Manage any development feature flag variables encountered
184 * The use of dev features indicated by environment variables is
185 * logged in order to ensure that runs with such features enabled can
186 * be identified from their log and standard output. Any cross
187 * dependencies are also checked, and if unsatisfied, a fatal error
190 * Note that some development features overrides are applied already here:
191 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
193 * \param[in] mdlog Logger object.
194 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
195 * \param[in] pmeRunMode The PME run mode for this run
196 * \returns The object populated with development feature flags.
198 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
199 const bool useGpuForNonbonded,
200 const PmeRunMode pmeRunMode)
202 DevelopmentFeatureFlags devFlags;
204 // Some builds of GCC 5 give false positive warnings that these
205 // getenv results are ignored when clearly they are used.
206 #pragma GCC diagnostic push
207 #pragma GCC diagnostic ignored "-Wunused-result"
209 devFlags.enableGpuBufferOps =
210 GMX_GPU_CUDA && useGpuForNonbonded && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
211 devFlags.enableGpuHaloExchange = GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_DD_COMMS") != nullptr;
212 devFlags.enableGpuPmePPComm =
213 GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
215 #pragma GCC diagnostic pop
217 if (devFlags.enableGpuBufferOps)
219 GMX_LOG(mdlog.warning)
221 .appendTextFormatted(
222 "This run uses the 'GPU buffer ops' feature, enabled by the "
223 "GMX_USE_GPU_BUFFER_OPS environment variable.");
226 if (devFlags.forceGpuUpdateDefault)
228 GMX_LOG(mdlog.warning)
230 .appendTextFormatted(
231 "This run will default to '-update gpu' as requested by the "
232 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
233 "decomposition lacks substantial testing and should be used with caution.");
236 if (devFlags.enableGpuHaloExchange)
238 if (useGpuForNonbonded)
240 if (!devFlags.enableGpuBufferOps)
242 GMX_LOG(mdlog.warning)
244 .appendTextFormatted(
245 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
246 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
247 devFlags.enableGpuBufferOps = true;
249 GMX_LOG(mdlog.warning)
251 .appendTextFormatted(
252 "This run has requested the 'GPU halo exchange' feature, enabled by "
254 "GMX_GPU_DD_COMMS environment variable.");
258 GMX_LOG(mdlog.warning)
260 .appendTextFormatted(
261 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
262 "halo exchange' feature will not be enabled as nonbonded interactions "
263 "are not offloaded.");
264 devFlags.enableGpuHaloExchange = false;
268 if (devFlags.enableGpuPmePPComm)
270 if (pmeRunMode == PmeRunMode::GPU)
272 if (!devFlags.enableGpuBufferOps)
274 GMX_LOG(mdlog.warning)
276 .appendTextFormatted(
277 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
278 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
279 devFlags.enableGpuBufferOps = true;
281 GMX_LOG(mdlog.warning)
283 .appendTextFormatted(
284 "This run uses the 'GPU PME-PP communications' feature, enabled "
285 "by the GMX_GPU_PME_PP_COMMS environment variable.");
289 std::string clarification;
290 if (pmeRunMode == PmeRunMode::Mixed)
293 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
298 clarification = "PME is not offloaded to the GPU.";
300 GMX_LOG(mdlog.warning)
303 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
304 "'GPU PME-PP communications' feature was not enabled as "
306 devFlags.enableGpuPmePPComm = false;
313 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
315 * Used to ensure that the master thread does not modify mdrunner during copy
316 * on the spawned threads. */
317 static void threadMpiMdrunnerAccessBarrier()
320 MPI_Barrier(MPI_COMM_WORLD);
324 Mdrunner Mdrunner::cloneOnSpawnedThread() const
326 auto newRunner = Mdrunner(std::make_unique<MDModules>());
328 // All runners in the same process share a restraint manager resource because it is
329 // part of the interface to the client code, which is associated only with the
330 // original thread. Handles to the same resources can be obtained by copy.
332 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
335 // Copy members of master runner.
336 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
337 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
338 newRunner.hw_opt = hw_opt;
339 newRunner.filenames = filenames;
341 newRunner.oenv = oenv;
342 newRunner.mdrunOptions = mdrunOptions;
343 newRunner.domdecOptions = domdecOptions;
344 newRunner.nbpu_opt = nbpu_opt;
345 newRunner.pme_opt = pme_opt;
346 newRunner.pme_fft_opt = pme_fft_opt;
347 newRunner.bonded_opt = bonded_opt;
348 newRunner.update_opt = update_opt;
349 newRunner.nstlist_cmdline = nstlist_cmdline;
350 newRunner.replExParams = replExParams;
351 newRunner.pforce = pforce;
352 // Give the spawned thread the newly created valid communicator
353 // for the simulation.
354 newRunner.communicator = MPI_COMM_WORLD;
356 newRunner.startingBehavior = startingBehavior;
357 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
359 threadMpiMdrunnerAccessBarrier();
364 /*! \brief The callback used for running on spawned threads.
366 * Obtains the pointer to the master mdrunner object from the one
367 * argument permitted to the thread-launch API call, copies it to make
368 * a new runner for this thread, reinitializes necessary data, and
369 * proceeds to the simulation. */
370 static void mdrunner_start_fn(const void* arg)
374 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
375 /* copy the arg list to make sure that it's thread-local. This
376 doesn't copy pointed-to items, of course; fnm, cr and fplog
377 are reset in the call below, all others should be const. */
378 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
381 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
385 void Mdrunner::spawnThreads(int numThreadsToLaunch)
388 /* now spawn new threads that start mdrunner_start_fn(), while
389 the main thread returns. Thread affinity is handled later. */
390 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
391 static_cast<const void*>(this))
394 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
397 // Give the master thread the newly created valid communicator for
399 communicator = MPI_COMM_WORLD;
400 threadMpiMdrunnerAccessBarrier();
402 GMX_UNUSED_VALUE(numThreadsToLaunch);
403 GMX_UNUSED_VALUE(mdrunner_start_fn);
409 /*! \brief Initialize variables for Verlet scheme simulation */
410 static void prepare_verlet_scheme(FILE* fplog,
414 const gmx_mtop_t* mtop,
416 bool makeGpuPairList,
417 const gmx::CpuInfo& cpuinfo)
419 // We checked the cut-offs in grompp, but double-check here.
420 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
421 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
423 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
424 "With Verlet lists and PME we should have rcoulomb>=rvdw");
428 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
429 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
431 /* For NVE simulations, we will retain the initial list buffer */
432 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
434 /* Update the Verlet buffer size for the current run setup */
436 /* Here we assume SIMD-enabled kernels are being used. But as currently
437 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
438 * and 4x2 gives a larger buffer than 4x4, this is ok.
440 ListSetupType listType =
441 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
442 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
444 const real rlist_new =
445 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
447 if (rlist_new != ir->rlist)
449 if (fplog != nullptr)
452 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
453 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
455 ir->rlist = rlist_new;
459 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
461 gmx_fatal(FARGS, "Can not set nstlist without %s",
462 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
465 if (EI_DYNAMICS(ir->eI))
467 /* Set or try nstlist values */
468 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
472 /*! \brief Override the nslist value in inputrec
474 * with value passed on the command line (if any)
476 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
480 /* override with anything else than the default -2 */
481 if (nsteps_cmdline > -2)
483 char sbuf_steps[STEPSTRSIZE];
484 char sbuf_msg[STRLEN];
486 ir->nsteps = nsteps_cmdline;
487 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
490 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
491 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
495 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
496 gmx_step_str(nsteps_cmdline, sbuf_steps));
499 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
501 else if (nsteps_cmdline < -2)
503 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
505 /* Do nothing if nsteps_cmdline == -2 */
511 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
513 * If not, and if a warning may be issued, logs a warning about
514 * falling back to CPU code. With thread-MPI, only the first
515 * call to this function should have \c issueWarning true. */
516 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
518 bool gpuIsUseful = true;
521 if (ir.opts.ngener - ir.nwall > 1)
523 /* The GPU code does not support more than one energy group.
524 * If the user requested GPUs explicitly, a fatal error is given later.
528 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
529 "For better performance, run on the GPU without energy groups and then do "
530 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
536 warning = "TPI is not implemented for GPUs.";
539 if (!gpuIsUseful && issueWarning)
541 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
547 //! Initializes the logger for mdrun.
548 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
550 gmx::LoggerBuilder builder;
551 if (fplog != nullptr)
553 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
555 if (isSimulationMasterRank)
557 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
559 return builder.build();
562 //! Make a TaskTarget from an mdrun argument string.
563 static TaskTarget findTaskTarget(const char* optionString)
565 TaskTarget returnValue = TaskTarget::Auto;
567 if (strncmp(optionString, "auto", 3) == 0)
569 returnValue = TaskTarget::Auto;
571 else if (strncmp(optionString, "cpu", 3) == 0)
573 returnValue = TaskTarget::Cpu;
575 else if (strncmp(optionString, "gpu", 3) == 0)
577 returnValue = TaskTarget::Gpu;
581 GMX_ASSERT(false, "Option string should have been checked for sanity already");
587 //! Finish run, aggregate data to print performance info.
588 static void finish_run(FILE* fplog,
589 const gmx::MDLogger& mdlog,
591 const t_inputrec* inputrec,
593 gmx_wallcycle_t wcycle,
594 gmx_walltime_accounting_t walltime_accounting,
595 nonbonded_verlet_t* nbv,
596 const gmx_pme_t* pme,
600 double nbfs = 0, mflop = 0;
601 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
602 elapsed_time_over_all_threads_over_all_ranks;
603 /* Control whether it is valid to print a report. Only the
604 simulation master may print, but it should not do so if the run
605 terminated e.g. before a scheduled reset step. This is
606 complicated by the fact that PME ranks are unaware of the
607 reason why they were sent a pmerecvqxFINISH. To avoid
608 communication deadlocks, we always do the communication for the
609 report, even if we've decided not to write the report, because
610 how long it takes to finish the run is not important when we've
611 decided not to report on the simulation performance.
613 Further, we only report performance for dynamical integrators,
614 because those are the only ones for which we plan to
615 consider doing any optimizations. */
616 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
618 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
620 GMX_LOG(mdlog.warning)
622 .appendText("Simulation ended prematurely, no performance report will be written.");
627 std::unique_ptr<t_nrnb> nrnbTotalStorage;
630 nrnbTotalStorage = std::make_unique<t_nrnb>();
631 nrnb_tot = nrnbTotalStorage.get();
633 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
641 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
642 elapsed_time_over_all_threads =
643 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
647 /* reduce elapsed_time over all MPI ranks in the current simulation */
648 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
650 elapsed_time_over_all_ranks /= cr->nnodes;
651 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
652 * current simulation. */
653 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
654 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
659 elapsed_time_over_all_ranks = elapsed_time;
660 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
665 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
668 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
670 print_dd_statistics(cr, inputrec, fplog);
673 /* TODO Move the responsibility for any scaling by thread counts
674 * to the code that handled the thread region, so that there's a
675 * mechanism to keep cycle counting working during the transition
676 * to task parallelism. */
677 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
678 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
679 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
680 nthreads_pp, nthreads_pme);
681 auto cycle_sum(wallcycle_sum(cr, wcycle));
685 auto nbnxn_gpu_timings =
686 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
687 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
689 if (pme_gpu_task_enabled(pme))
691 pme_gpu_get_timings(pme, &pme_gpu_timings);
693 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
694 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
697 if (EI_DYNAMICS(inputrec->eI))
699 delta_t = inputrec->delta_t;
704 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
705 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
706 delta_t, nbfs, mflop);
710 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
711 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
712 delta_t, nbfs, mflop);
717 int Mdrunner::mdrunner()
720 t_forcerec* fr = nullptr;
721 real ewaldcoeff_q = 0;
722 real ewaldcoeff_lj = 0;
723 int nChargePerturbed = -1, nTypePerturbed = 0;
724 gmx_wallcycle_t wcycle;
725 gmx_walltime_accounting_t walltime_accounting = nullptr;
726 MembedHolder membedHolder(filenames.size(), filenames.data());
727 gmx_hw_info_t* hwinfo = nullptr;
729 /* CAUTION: threads may be started later on in this function, so
730 cr doesn't reflect the final parallel state right now */
733 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
734 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
735 const bool doRerun = mdrunOptions.rerun;
737 // Handle task-assignment related user options.
738 EmulateGpuNonbonded emulateGpuNonbonded =
739 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
741 std::vector<int> userGpuTaskAssignment;
744 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
746 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
747 auto nonbondedTarget = findTaskTarget(nbpu_opt);
748 auto pmeTarget = findTaskTarget(pme_opt);
749 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
750 auto bondedTarget = findTaskTarget(bonded_opt);
751 auto updateTarget = findTaskTarget(update_opt);
753 FILE* fplog = nullptr;
754 // If we are appending, we don't write log output because we need
755 // to check that the old log file matches what the checkpoint file
756 // expects. Otherwise, we should start to write log output now if
757 // there is a file ready for it.
758 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
760 fplog = gmx_fio_getfp(logFileHandle);
762 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
763 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
764 gmx::MDLogger mdlog(logOwner.logger());
766 // TODO The thread-MPI master rank makes a working
767 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
768 // after the threads have been launched. This works because no use
769 // is made of that communicator until after the execution paths
770 // have rejoined. But it is likely that we can improve the way
771 // this is expressed, e.g. by expressly running detection only the
772 // master rank for thread-MPI, rather than relying on the mutex
773 // and reference count.
774 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
775 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
777 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
779 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->deviceInfoList, hw_opt.gpuIdsAvailable);
781 // Print citation requests after all software/hardware printing
782 pleaseCiteGromacs(fplog);
784 // TODO: Use SimulationInputHolder member to access SimulationInput. Issue #3374.
785 const auto* const tprFilename = ftp2fn(efTPR, filenames.size(), filenames.data());
786 const auto* const cpiFilename = opt2fn("-cpi", filenames.size(), filenames.data());
787 // Note that, as of this change, there is no public API for simulationInput or its creation.
788 // TODO: (#3374) Public API for providing simulationInput from client.
789 auto simulationInput = detail::makeSimulationInput(tprFilename, cpiFilename);
791 // Note: legacy program logic relies on checking whether these pointers are assigned.
792 // Objects may or may not be allocated later.
793 std::unique_ptr<t_inputrec> inputrec;
794 std::unique_ptr<t_state> globalState;
796 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
798 if (isSimulationMasterRank)
800 // Allocate objects to be initialized by later function calls.
801 /* Only the master rank has the global state */
802 globalState = std::make_unique<t_state>();
803 inputrec = std::make_unique<t_inputrec>();
805 /* Read (nearly) all data required for the simulation
806 * and keep the partly serialized tpr contents to send to other ranks later
808 applyGlobalSimulationState(*simulationInput.object_, partialDeserializedTpr.get(),
809 globalState.get(), inputrec.get(), &mtop);
812 /* Check and update the hardware options for internal consistency */
813 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
816 if (GMX_THREAD_MPI && isSimulationMasterRank)
818 bool useGpuForNonbonded = false;
819 bool useGpuForPme = false;
822 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
824 // If the user specified the number of ranks, then we must
825 // respect that, but in default mode, we need to allow for
826 // the number of GPUs to choose the number of ranks.
827 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
828 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
829 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
830 canUseGpuForNonbonded,
831 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
832 hw_opt.nthreads_tmpi);
833 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
834 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
835 *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
837 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
839 /* Determine how many thread-MPI ranks to start.
841 * TODO Over-writing the user-supplied value here does
842 * prevent any possible subsequent checks from working
844 hw_opt.nthreads_tmpi =
845 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
846 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
848 // Now start the threads for thread MPI.
849 spawnThreads(hw_opt.nthreads_tmpi);
850 // The spawned threads enter mdrunner() and execution of
851 // master and spawned threads joins at the end of this block.
852 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
855 GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
856 "Must have valid world communicator unless running a multi-simulation");
857 CommrecHandle crHandle = init_commrec(communicator);
858 t_commrec* cr = crHandle.get();
859 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
863 /* now broadcast everything to the non-master nodes/threads: */
864 if (!isSimulationMasterRank)
866 // Until now, only the master rank has a non-null pointer.
867 // On non-master ranks, allocate the object that will receive data in the following call.
868 inputrec = std::make_unique<t_inputrec>();
870 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec.get(), &mtop,
871 partialDeserializedTpr.get());
873 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
874 partialDeserializedTpr.reset(nullptr);
876 // Now the number of ranks is known to all ranks, and each knows
877 // the inputrec read by the master rank. The ranks can now all run
878 // the task-deciding functions and will agree on the result
879 // without needing to communicate.
880 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
882 // Note that these variables describe only their own node.
884 // Note that when bonded interactions run on a GPU they always run
885 // alongside a nonbonded task, so do not influence task assignment
886 // even though they affect the force calculation workload.
887 bool useGpuForNonbonded = false;
888 bool useGpuForPme = false;
889 bool useGpuForBonded = false;
890 bool useGpuForUpdate = false;
891 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
894 // It's possible that there are different numbers of GPUs on
895 // different nodes, which is the user's responsibility to
896 // handle. If unsuitable, we will notice that during task
898 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
899 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
900 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
901 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
902 useGpuForPme = decideWhetherToUseGpusForPme(
903 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec,
904 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
905 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
906 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
907 useGpuForBonded = decideWhetherToUseGpusForBonded(
908 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
909 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
910 domdecOptions.numPmeRanks, gpusWereDetected);
912 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
914 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
916 // Initialize development feature flags that enabled by environment variable
917 // and report those features that are enabled.
918 const DevelopmentFeatureFlags devFlags =
919 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
921 const bool useModularSimulator =
922 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
923 nullptr, doEssentialDynamics, membedHolder.doMembed());
926 // TODO: hide restraint implementation details from Mdrunner.
927 // There is nothing unique about restraints at this point as far as the
928 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
929 // factory functions from the SimulationContext on which to call mdModules_->add().
930 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
931 for (auto&& restraint : restraintManager_->getRestraints())
933 auto module = RestraintMDModule::create(restraint, restraint->sites());
934 mdModules_->add(std::move(module));
937 // TODO: Error handling
938 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
939 // now that the MdModules know their options, they know which callbacks to sign up to
940 mdModules_->subscribeToSimulationSetupNotifications();
941 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
943 if (inputrec->internalParameters != nullptr)
945 mdModulesNotifier.notify(*inputrec->internalParameters);
948 if (fplog != nullptr)
950 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
951 fprintf(fplog, "\n");
956 /* In rerun, set velocities to zero if present */
957 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
959 // rerun does not use velocities
963 "Rerun trajectory contains velocities. Rerun does only evaluate "
964 "potential energy and forces. The velocities will be ignored.");
965 for (int i = 0; i < globalState->natoms; i++)
967 clear_rvec(globalState->v[i]);
969 globalState->flags &= ~(1 << estV);
972 /* now make sure the state is initialized and propagated */
973 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
976 /* NM and TPI parallelize over force/energy calculations, not atoms,
977 * so we need to initialize and broadcast the global state.
979 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
983 globalState = std::make_unique<t_state>();
985 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
989 /* A parallel command line option consistency check that we can
990 only do after any threads have started. */
992 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
993 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
996 "The -dd or -npme option request a parallel simulation, "
998 "but %s was compiled without threads or MPI enabled",
999 output_env_get_program_display_name(oenv));
1000 #elif GMX_THREAD_MPI
1001 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
1003 "but %s was not started through mpirun/mpiexec or only one rank was requested "
1004 "through mpirun/mpiexec",
1005 output_env_get_program_display_name(oenv));
1009 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1012 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1013 "these are not compatible with mdrun -rerun");
1016 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1018 if (domdecOptions.numPmeRanks > 0)
1020 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1021 "PME-only ranks are requested, but the system does not use PME "
1022 "for electrostatics or LJ");
1025 domdecOptions.numPmeRanks = 0;
1028 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1030 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1031 * improve performance with many threads per GPU, since our OpenMP
1032 * scaling is bad, but it's difficult to automate the setup.
1034 domdecOptions.numPmeRanks = 0;
1038 if (domdecOptions.numPmeRanks < 0)
1040 domdecOptions.numPmeRanks = 0;
1041 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1045 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1046 "PME GPU decomposition is not supported");
1053 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1057 /* NMR restraints must be initialized before load_checkpoint,
1058 * since with time averaging the history is added to t_state.
1059 * For proper consistency check we therefore need to extend
1061 * So the PME-only nodes (if present) will also initialize
1062 * the distance restraints.
1065 /* This needs to be called before read_checkpoint to extend the state */
1066 t_disresdata* disresdata;
1067 snew(disresdata, 1);
1068 init_disres(fplog, &mtop, inputrec.get(), DisResRunMode::MDRun,
1069 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1070 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1071 globalState.get(), replExParams.exchangeInterval > 0);
1073 t_oriresdata* oriresdata;
1074 snew(oriresdata, 1);
1075 init_orires(fplog, &mtop, inputrec.get(), cr, ms, globalState.get(), oriresdata);
1077 auto deform = prepareBoxDeformation(
1078 globalState != nullptr ? globalState->box : box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1079 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mygroup, *inputrec);
1081 ObservablesHistory observablesHistory = {};
1083 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1084 if (startingBehavior != StartingBehavior::NewSimulation)
1086 /* Check if checkpoint file exists before doing continuation.
1087 * This way we can use identical input options for the first and subsequent runs...
1089 if (mdrunOptions.numStepsCommandline > -2)
1091 /* Temporarily set the number of steps to unlimited to avoid
1092 * triggering the nsteps check in load_checkpoint().
1093 * This hack will go away soon when the -nsteps option is removed.
1095 inputrec->nsteps = -1;
1098 // Finish applying initial simulation state information from external sources on all ranks.
1099 // Reconcile checkpoint file data with Mdrunner state established up to this point.
1100 applyLocalState(*simulationInput.object_, logFileHandle, cr, domdecOptions.numCells,
1101 inputrec.get(), globalState.get(), &observablesHistory,
1102 mdrunOptions.reproducible, mdModules_->notifier(),
1103 modularSimulatorCheckpointData.get(), useModularSimulator);
1104 // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1106 // on all code paths.
1107 // Write checkpoint or provide hook to update SimulationInput.
1108 // If there was a checkpoint file, SimulationInput contains more information
1109 // than if there wasn't. At this point, we have synchronized the in-memory
1110 // state with the filesystem state only for restarted simulations. We should
1111 // be calling applyLocalState unconditionally and expect that the completeness
1112 // of SimulationInput is not dependent on its creation method.
1114 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1116 // Now we can start normal logging to the truncated log file.
1117 fplog = gmx_fio_getfp(logFileHandle);
1118 prepareLogAppending(fplog);
1119 logOwner = buildLogger(fplog, MASTER(cr));
1120 mdlog = logOwner.logger();
1124 if (mdrunOptions.numStepsCommandline > -2)
1129 "The -nsteps functionality is deprecated, and may be removed in a future "
1131 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1134 /* override nsteps with value set on the commandline */
1135 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1137 if (isSimulationMasterRank)
1139 copy_mat(globalState->box, box);
1144 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1147 if (inputrec->cutoff_scheme != ecutsVERLET)
1150 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1151 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1153 /* Update rlist and nstlist. */
1154 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1155 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1156 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1158 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1159 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1162 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1163 // This builder is necessary while we have multi-part construction
1164 // of DD. Before DD is constructed, we use the existence of
1165 // the builder object to indicate that further construction of DD
1167 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1168 if (useDomainDecomposition)
1170 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1171 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1172 positionsFromStatePointer(globalState.get()));
1176 /* PME, if used, is done on all nodes with 1D decomposition */
1177 cr->nnodes = cr->sizeOfDefaultCommunicator;
1178 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1179 cr->nodeid = cr->rankInDefaultCommunicator;
1181 cr->duty = (DUTY_PP | DUTY_PME);
1183 if (inputrec->pbcType == PbcType::Screw)
1185 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1189 // Produce the task assignment for this rank - done after DD is constructed
1190 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1191 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1192 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1193 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1194 // TODO cr->duty & DUTY_PME should imply that a PME
1195 // algorithm is active, but currently does not.
1196 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1198 // Get the device handles for the modules, nullptr when no task is assigned.
1200 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1202 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1203 bool useTiming = true;
1207 /* WARNING: CUDA timings are incorrect with multiple streams.
1208 * This is the main reason why they are disabled by default.
1210 // TODO: Consider turning on by default when we can detect nr of streams.
1211 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1213 else if (GMX_GPU_OPENCL)
1215 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1218 // TODO Currently this is always built, yet DD partition code
1219 // checks if it is built before using it. Probably it should
1220 // become an MDModule that is made only when another module
1221 // requires it (e.g. pull, CompEl, density fitting), so that we
1222 // don't update the local atom sets unilaterally every step.
1223 LocalAtomSetManager atomSets;
1226 // TODO Pass the GPU streams to ddBuilder to use in buffer
1227 // transfers (e.g. halo exchange)
1228 cr->dd = ddBuilder->build(&atomSets);
1229 // The builder's job is done, so destruct it
1230 ddBuilder.reset(nullptr);
1231 // Note that local state still does not exist yet.
1234 // The GPU update is decided here because we need to know whether the constraints or
1235 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1236 // defined). This is only known after DD is initialized, hence decision on using GPU
1237 // update is done so late.
1240 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1242 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1243 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1244 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1245 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1246 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1248 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1250 const bool printHostName = (cr->nnodes > 1);
1251 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1253 MdrunScheduleWorkload runScheduleWork;
1254 // Also populates the simulation constant workload description.
1255 runScheduleWork.simulationWork = createSimulationWorkload(
1256 *inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded, useGpuForUpdate,
1257 devFlags.enableGpuBufferOps, devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1259 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1261 if (deviceInfo != nullptr)
1263 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1265 dd_setup_dlb_resource_sharing(cr, deviceId);
1267 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1268 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, useTiming);
1271 // If the user chose a task assignment, give them some hints
1272 // where appropriate.
1273 if (!userGpuTaskAssignment.empty())
1275 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1280 /* After possible communicator splitting in make_dd_communicators.
1281 * we can set up the intra/inter node communication.
1283 gmx_setup_nodecomm(fplog, cr);
1289 GMX_LOG(mdlog.warning)
1291 .appendTextFormatted(
1292 "This is simulation %d out of %d running as a composite GROMACS\n"
1293 "multi-simulation job. Setup for this simulation:\n",
1294 ms->simulationIndex_, ms->numSimulations_);
1296 GMX_LOG(mdlog.warning)
1297 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1299 cr->nnodes == 1 ? "thread" : "threads"
1301 cr->nnodes == 1 ? "process" : "processes"
1307 // If mdrun -pin auto honors any affinity setting that already
1308 // exists. If so, it is nice to provide feedback about whether
1309 // that existing affinity setting was from OpenMP or something
1310 // else, so we run this code both before and after we initialize
1311 // the OpenMP support.
1312 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1313 /* Check and update the number of OpenMP threads requested */
1314 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1315 pmeRunMode, mtop, *inputrec);
1317 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1318 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1320 // Enable FP exception detection, but not in
1321 // Release mode and not for compilers with known buggy FP
1322 // exception support (clang with any optimization) or suspected
1323 // buggy FP exception support (gcc 7.* with optimization).
1324 #if !defined NDEBUG \
1325 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1326 && defined __OPTIMIZE__)
1327 const bool bEnableFPE = true;
1329 const bool bEnableFPE = false;
1331 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1334 gmx_feenableexcept();
1337 /* Now that we know the setup is consistent, check for efficiency */
1338 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1339 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1341 /* getting number of PP/PME threads on this MPI / tMPI rank.
1342 PME: env variable should be read only on one node to make sure it is
1343 identical everywhere;
1345 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1346 : gmx_omp_nthreads_get(emntPME);
1347 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1348 physicalNodeComm, mdlog);
1350 // Enable Peer access between GPUs where available
1351 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1352 // any of the GPU communication features are active.
1353 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1354 && (runScheduleWork.simulationWork.useGpuHaloExchange
1355 || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1357 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1360 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1362 /* Before setting affinity, check whether the affinity has changed
1363 * - which indicates that probably the OpenMP library has changed it
1364 * since we first checked).
1366 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1368 int numThreadsOnThisNode, intraNodeThreadOffset;
1369 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1370 &intraNodeThreadOffset);
1372 /* Set the CPU affinity */
1373 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1374 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1377 if (mdrunOptions.timingOptions.resetStep > -1)
1382 "The -resetstep functionality is deprecated, and may be removed in a "
1385 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1389 /* Master synchronizes its value of reset_counters with all nodes
1390 * including PME only nodes */
1391 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1392 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1393 wcycle_set_reset_counters(wcycle, reset_counters);
1396 // Membrane embedding must be initialized before we call init_forcerec()
1397 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1398 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1400 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1401 std::unique_ptr<MDAtoms> mdAtoms;
1402 std::unique_ptr<VirtualSitesHandler> vsite;
1403 std::unique_ptr<GpuBonded> gpuBonded;
1406 if (thisRankHasDuty(cr, DUTY_PP))
1408 mdModulesNotifier.notify(*cr);
1409 mdModulesNotifier.notify(&atomSets);
1410 mdModulesNotifier.notify(inputrec->pbcType);
1411 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1412 /* Initiate forcerecord */
1413 fr = new t_forcerec;
1414 fr->forceProviders = mdModules_->initForceProviders();
1415 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1416 opt2fn("-table", filenames.size(), filenames.data()),
1417 opt2fn("-tablep", filenames.size(), filenames.data()),
1418 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1419 // Dirty hack, for fixing disres and orires should be made mdmodules
1420 fr->fcdata->disres = disresdata;
1421 fr->fcdata->orires = oriresdata;
1423 // Save a handle to device stream manager to use elsewhere in the code
1424 // TODO: Forcerec is not a correct place to store it.
1425 fr->deviceStreamManager = deviceStreamManager.get();
1427 if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1430 deviceStreamManager != nullptr,
1431 "GPU device stream manager should be valid in order to use PME-PP direct "
1434 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1435 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1437 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1438 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1439 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1442 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo,
1443 runScheduleWork.simulationWork.useGpuNonbonded,
1444 deviceStreamManager.get(), &mtop, box, wcycle);
1445 // TODO: Move the logic below to a GPU bonded builder
1446 if (runScheduleWork.simulationWork.useGpuBonded)
1448 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1449 "GPU device stream manager should be valid in order to use GPU "
1450 "version of bonded forces.");
1451 gpuBonded = std::make_unique<GpuBonded>(
1452 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1453 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1454 fr->gpuBonded = gpuBonded.get();
1457 /* Initialize the mdAtoms structure.
1458 * mdAtoms is not filled with atom data,
1459 * as this can not be done now with domain decomposition.
1461 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1462 if (globalState && thisRankHasPmeGpuTask)
1464 // The pinning of coordinates in the global state object works, because we only use
1465 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1466 // points to the global state object without DD.
1467 // FIXME: MD and EM separately set up the local state - this should happen in the same
1468 // function, which should also perform the pinning.
1469 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1472 /* Initialize the virtual site communication */
1473 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1475 calc_shifts(box, fr->shift_vec);
1477 /* With periodic molecules the charge groups should be whole at start up
1478 * and the virtual sites should not be far from their proper positions.
1480 if (!inputrec->bContinuation && MASTER(cr)
1481 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1483 /* Make molecules whole at start of run */
1484 if (fr->pbcType != PbcType::No)
1486 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1490 /* Correct initial vsite positions are required
1491 * for the initial distribution in the domain decomposition
1492 * and for the initial shell prediction.
1494 constructVirtualSitesGlobal(mtop, globalState->x);
1498 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1500 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1501 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1506 /* This is a PME only node */
1508 GMX_ASSERT(globalState == nullptr,
1509 "We don't need the state on a PME only rank and expect it to be unitialized");
1511 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1512 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1515 gmx_pme_t* sepPmeData = nullptr;
1516 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1517 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1518 "Double-checking that only PME-only ranks have no forcerec");
1519 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1521 // TODO should live in ewald module once its testing is improved
1523 // Later, this program could contain kernels that might be later
1524 // re-used as auto-tuning progresses, or subsequent simulations
1526 PmeGpuProgramStorage pmeGpuProgram;
1527 if (thisRankHasPmeGpuTask)
1530 (deviceStreamManager != nullptr),
1531 "GPU device stream manager should be initialized in order to use GPU for PME.");
1532 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1533 "GPU device should be initialized in order to use GPU for PME.");
1534 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1537 /* Initiate PME if necessary,
1538 * either on all nodes or on dedicated PME nodes only. */
1539 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1541 if (mdAtoms && mdAtoms->mdatoms())
1543 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1544 if (EVDW_PME(inputrec->vdwtype))
1546 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1549 if (cr->npmenodes > 0)
1551 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1552 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1553 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1556 if (thisRankHasDuty(cr, DUTY_PME))
1560 // TODO: This should be in the builder.
1561 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1562 || (deviceStreamManager != nullptr),
1563 "Device stream manager should be valid in order to use GPU "
1566 !runScheduleWork.simulationWork.useGpuPme
1567 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1568 "GPU PME stream should be valid in order to use GPU version of PME.");
1570 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1571 ? &deviceStreamManager->context()
1573 const DeviceStream* pmeStream =
1574 runScheduleWork.simulationWork.useGpuPme
1575 ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1578 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1579 nChargePerturbed != 0, nTypePerturbed != 0,
1580 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1581 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1582 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1584 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1589 if (EI_DYNAMICS(inputrec->eI))
1591 /* Turn on signal handling on all nodes */
1593 * (A user signal from the PME nodes (if any)
1594 * is communicated to the PP nodes.
1596 signal_handler_install();
1599 pull_t* pull_work = nullptr;
1600 if (thisRankHasDuty(cr, DUTY_PP))
1602 /* Assumes uniform use of the number of OpenMP threads */
1603 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1605 if (inputrec->bPull)
1607 /* Initialize pull code */
1608 pull_work = init_pull(fplog, inputrec->pull, inputrec.get(), &mtop, cr, &atomSets,
1609 inputrec->fepvals->init_lambda);
1610 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1612 initPullHistory(pull_work, &observablesHistory);
1614 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1616 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1620 std::unique_ptr<EnforcedRotation> enforcedRotation;
1623 /* Initialize enforced rotation code */
1624 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1625 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1629 t_swap* swap = nullptr;
1630 if (inputrec->eSwapCoords != eswapNO)
1632 /* Initialize ion swapping code */
1633 swap = init_swapcoords(fplog, inputrec.get(),
1634 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1635 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1636 oenv, mdrunOptions, startingBehavior);
1639 /* Let makeConstraints know whether we have essential dynamics constraints. */
1640 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1641 ms, &nrnb, wcycle, fr->bMolPBC);
1643 /* Energy terms and groups */
1644 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1645 inputrec->fepvals->n_lambda);
1647 /* Kinetic energy data */
1648 gmx_ekindata_t ekind;
1649 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1651 /* Set up interactive MD (IMD) */
1653 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1654 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1655 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1657 if (DOMAINDECOMP(cr))
1659 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1660 /* This call is not included in init_domain_decomposition mainly
1661 * because fr->cginfo_mb is set later.
1663 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1664 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1667 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1668 if (gpusWereDetected
1669 && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
1670 || runScheduleWork.simulationWork.useGpuBufferOps))
1672 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1673 ? GpuApiCallBehavior::Async
1674 : GpuApiCallBehavior::Sync;
1675 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1676 "GPU device stream manager should be initialized to use GPU.");
1677 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1678 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1679 fr->stateGpu = stateGpu.get();
1682 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1683 SimulatorBuilder simulatorBuilder;
1685 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1686 simulatorBuilder.add(std::move(membedHolder));
1687 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1688 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1691 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1692 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1693 simulatorBuilder.add(ConstraintsParam(
1694 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1696 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1697 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1698 inputrec.get(), fr));
1699 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1700 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1701 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1702 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1703 // Todo move to an MDModule
1704 simulatorBuilder.add(IonSwapping(swap));
1705 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1706 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1707 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
1709 // build and run simulator object based on user-input
1710 auto simulator = simulatorBuilder.build(useModularSimulator);
1713 if (fr->pmePpCommGpu)
1715 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1716 fr->pmePpCommGpu.reset();
1719 if (inputrec->bPull)
1721 finish_pull(pull_work);
1723 finish_swapcoords(swap);
1727 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1729 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1730 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1731 deviceStreamManager.get());
1734 wallcycle_stop(wcycle, ewcRUN);
1736 /* Finish up, write some stuff
1737 * if rerunMD, don't write last frame again
1739 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1740 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1742 // clean up cycle counter
1743 wallcycle_destroy(wcycle);
1745 deviceStreamManager.reset(nullptr);
1749 gmx_pme_destroy(pmedata);
1753 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1754 // before we destroy the GPU context(s)
1755 // Pinned buffers are associated with contexts in CUDA.
1756 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1757 mdAtoms.reset(nullptr);
1758 globalState.reset(nullptr);
1759 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1760 gpuBonded.reset(nullptr);
1761 /* Free pinned buffers in *fr */
1764 // TODO convert to C++ so we can get rid of these frees
1768 if (!hwinfo->deviceInfoList.empty())
1770 /* stop the GPU profiler (only CUDA) */
1774 /* With tMPI we need to wait for all ranks to finish deallocation before
1775 * destroying the CUDA context as some tMPI ranks may be sharing
1778 * This is not a concern in OpenCL where we use one context per rank.
1780 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1781 * but it is easier and more futureproof to call it on the whole node.
1783 * Note that this function needs to be called even if GPUs are not used
1784 * in this run because the PME ranks have no knowledge of whether GPUs
1785 * are used or not, but all ranks need to enter the barrier below.
1786 * \todo Remove this physical node barrier after making sure
1787 * that it's not needed anymore (with a shared GPU run).
1791 physicalNodeComm.barrier();
1793 releaseDevice(deviceInfo);
1795 /* Does what it says */
1796 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1797 walltime_accounting_destroy(walltime_accounting);
1799 // Ensure log file content is written
1802 gmx_fio_flush(logFileHandle);
1805 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1806 * exceptions were enabled before function was called. */
1809 gmx_fedisableexcept();
1812 auto rc = static_cast<int>(gmx_get_stop_condition());
1815 /* we need to join all threads. The sub-threads join when they
1816 exit this function, but the master thread needs to be told to
1826 Mdrunner::~Mdrunner()
1828 // Clean up of the Manager.
1829 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1830 // but okay as long as threads synchronize some time before adding or accessing
1831 // a new set of restraints.
1832 if (restraintManager_)
1834 restraintManager_->clear();
1835 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1836 "restraints added during runner life time should be cleared at runner "
1841 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1843 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1844 // Not sure if this should be logged through the md logger or something else,
1845 // but it is helpful to have some sort of INFO level message sent somewhere.
1846 // std::cout << "Registering restraint named " << name << std::endl;
1848 // When multiple restraints are used, it may be wasteful to register them separately.
1849 // Maybe instead register an entire Restraint Manager as a force provider.
1850 restraintManager_->addToSpec(std::move(puller), name);
1853 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1855 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1857 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1859 class Mdrunner::BuilderImplementation
1862 BuilderImplementation() = delete;
1863 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1864 ~BuilderImplementation();
1866 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1867 real forceWarningThreshold,
1868 StartingBehavior startingBehavior);
1870 void addDomdec(const DomdecOptions& options);
1872 void addVerletList(int nstlist);
1874 void addReplicaExchange(const ReplicaExchangeParameters& params);
1876 void addNonBonded(const char* nbpu_opt);
1878 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1880 void addBondedTaskAssignment(const char* bonded_opt);
1882 void addUpdateTaskAssignment(const char* update_opt);
1884 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1886 void addFilenames(ArrayRef<const t_filenm> filenames);
1888 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1890 void addLogFile(t_fileio* logFileHandle);
1892 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1897 // Default parameters copied from runner.h
1898 // \todo Clarify source(s) of default parameters.
1900 const char* nbpu_opt_ = nullptr;
1901 const char* pme_opt_ = nullptr;
1902 const char* pme_fft_opt_ = nullptr;
1903 const char* bonded_opt_ = nullptr;
1904 const char* update_opt_ = nullptr;
1906 MdrunOptions mdrunOptions_;
1908 DomdecOptions domdecOptions_;
1910 ReplicaExchangeParameters replicaExchangeParameters_;
1912 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1915 //! Multisim communicator handle.
1916 gmx_multisim_t* multiSimulation_;
1918 //! mdrun communicator
1919 MPI_Comm communicator_ = MPI_COMM_NULL;
1921 //! Print a warning if any force is larger than this (in kJ/mol nm).
1922 real forceWarningThreshold_ = -1;
1924 //! Whether the simulation will start afresh, or restart with/without appending.
1925 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1927 //! The modules that comprise the functionality of mdrun.
1928 std::unique_ptr<MDModules> mdModules_;
1930 //! \brief Parallelism information.
1931 gmx_hw_opt_t hardwareOptions_;
1933 //! filename options for simulation.
1934 ArrayRef<const t_filenm> filenames_;
1936 /*! \brief Handle to output environment.
1938 * \todo gmx_output_env_t needs lifetime management.
1940 gmx_output_env_t* outputEnvironment_ = nullptr;
1942 /*! \brief Non-owning handle to MD log file.
1944 * \todo Context should own output facilities for client.
1945 * \todo Improve log file handle management.
1947 * Code managing the FILE* relies on the ability to set it to
1948 * nullptr to check whether the filehandle is valid.
1950 t_fileio* logFileHandle_ = nullptr;
1953 * \brief Builder for simulation stop signal handler.
1955 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1958 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1959 compat::not_null<SimulationContext*> context) :
1960 mdModules_(std::move(mdModules))
1962 communicator_ = context->communicator_;
1963 multiSimulation_ = context->multiSimulation_.get();
1966 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1968 Mdrunner::BuilderImplementation&
1969 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1970 const real forceWarningThreshold,
1971 const StartingBehavior startingBehavior)
1973 mdrunOptions_ = options;
1974 forceWarningThreshold_ = forceWarningThreshold;
1975 startingBehavior_ = startingBehavior;
1979 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1981 domdecOptions_ = options;
1984 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1989 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1991 replicaExchangeParameters_ = params;
1994 Mdrunner Mdrunner::BuilderImplementation::build()
1996 auto newRunner = Mdrunner(std::move(mdModules_));
1998 newRunner.mdrunOptions = mdrunOptions_;
1999 newRunner.pforce = forceWarningThreshold_;
2000 newRunner.startingBehavior = startingBehavior_;
2001 newRunner.domdecOptions = domdecOptions_;
2003 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2004 newRunner.hw_opt = hardwareOptions_;
2006 // No invariant to check. This parameter exists to optionally override other behavior.
2007 newRunner.nstlist_cmdline = nstlist_;
2009 newRunner.replExParams = replicaExchangeParameters_;
2011 newRunner.filenames = filenames_;
2013 newRunner.communicator = communicator_;
2015 // nullptr is a valid value for the multisim handle
2016 newRunner.ms = multiSimulation_;
2018 // \todo Clarify ownership and lifetime management for gmx_output_env_t
2019 // \todo Update sanity checking when output environment has clearly specified invariants.
2020 // Initialization and default values for oenv are not well specified in the current version.
2021 if (outputEnvironment_)
2023 newRunner.oenv = outputEnvironment_;
2027 GMX_THROW(gmx::APIError(
2028 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2031 newRunner.logFileHandle = logFileHandle_;
2035 newRunner.nbpu_opt = nbpu_opt_;
2039 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2042 if (pme_opt_ && pme_fft_opt_)
2044 newRunner.pme_opt = pme_opt_;
2045 newRunner.pme_fft_opt = pme_fft_opt_;
2049 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2054 newRunner.bonded_opt = bonded_opt_;
2058 GMX_THROW(gmx::APIError(
2059 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2064 newRunner.update_opt = update_opt_;
2068 GMX_THROW(gmx::APIError(
2069 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2073 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2075 if (stopHandlerBuilder_)
2077 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2081 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2087 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2089 nbpu_opt_ = nbpu_opt;
2092 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2095 pme_fft_opt_ = pme_fft_opt;
2098 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2100 bonded_opt_ = bonded_opt;
2103 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2105 update_opt_ = update_opt;
2108 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2110 hardwareOptions_ = hardwareOptions;
2113 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2115 filenames_ = filenames;
2118 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2120 outputEnvironment_ = outputEnvironment;
2123 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2125 logFileHandle_ = logFileHandle;
2128 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2130 stopHandlerBuilder_ = std::move(builder);
2133 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2134 compat::not_null<SimulationContext*> context) :
2135 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2139 MdrunnerBuilder::~MdrunnerBuilder() = default;
2141 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2142 real forceWarningThreshold,
2143 const StartingBehavior startingBehavior)
2145 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2149 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2151 impl_->addDomdec(options);
2155 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2157 impl_->addVerletList(nstlist);
2161 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2163 impl_->addReplicaExchange(params);
2167 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2169 impl_->addNonBonded(nbpu_opt);
2173 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2175 // The builder method may become more general in the future, but in this version,
2176 // parameters for PME electrostatics are both required and the only parameters
2178 if (pme_opt && pme_fft_opt)
2180 impl_->addPME(pme_opt, pme_fft_opt);
2185 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2190 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2192 impl_->addBondedTaskAssignment(bonded_opt);
2196 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2198 impl_->addUpdateTaskAssignment(update_opt);
2202 Mdrunner MdrunnerBuilder::build()
2204 return impl_->build();
2207 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2209 impl_->addHardwareOptions(hardwareOptions);
2213 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2215 impl_->addFilenames(filenames);
2219 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2221 impl_->addOutputEnvironment(outputEnvironment);
2225 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2227 impl_->addLogFile(logFileHandle);
2231 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2233 impl_->addStopHandlerBuilder(std::move(builder));
2237 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2239 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;