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,2021, 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.h"
68 #include "gromacs/ewald/pme_gpu_program.h"
69 #include "gromacs/ewald/pme_only.h"
70 #include "gromacs/ewald/pme_pp_comm_gpu.h"
71 #include "gromacs/fileio/checkpoint.h"
72 #include "gromacs/fileio/gmxfio.h"
73 #include "gromacs/fileio/oenv.h"
74 #include "gromacs/fileio/tpxio.h"
75 #include "gromacs/gmxlib/network.h"
76 #include "gromacs/gmxlib/nrnb.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/hardwaretopology.h"
82 #include "gromacs/hardware/printhardware.h"
83 #include "gromacs/imd/imd.h"
84 #include "gromacs/listed_forces/disre.h"
85 #include "gromacs/listed_forces/gpubonded.h"
86 #include "gromacs/listed_forces/listed_forces.h"
87 #include "gromacs/listed_forces/orires.h"
88 #include "gromacs/math/functions.h"
89 #include "gromacs/math/utilities.h"
90 #include "gromacs/math/vec.h"
91 #include "gromacs/mdlib/boxdeformation.h"
92 #include "gromacs/mdlib/broadcaststructs.h"
93 #include "gromacs/mdlib/calc_verletbuf.h"
94 #include "gromacs/mdlib/dispersioncorrection.h"
95 #include "gromacs/mdlib/enerdata_utils.h"
96 #include "gromacs/mdlib/force.h"
97 #include "gromacs/mdlib/forcerec.h"
98 #include "gromacs/mdlib/gmx_omp_nthreads.h"
99 #include "gromacs/mdlib/gpuforcereduction.h"
100 #include "gromacs/mdlib/makeconstraints.h"
101 #include "gromacs/mdlib/md_support.h"
102 #include "gromacs/mdlib/mdatoms.h"
103 #include "gromacs/mdlib/sighandler.h"
104 #include "gromacs/mdlib/stophandler.h"
105 #include "gromacs/mdlib/tgroup.h"
106 #include "gromacs/mdlib/updategroups.h"
107 #include "gromacs/mdlib/vsite.h"
108 #include "gromacs/mdrun/mdmodules.h"
109 #include "gromacs/mdrun/simulationcontext.h"
110 #include "gromacs/mdrun/simulationinput.h"
111 #include "gromacs/mdrun/simulationinputhandle.h"
112 #include "gromacs/mdrunutility/handlerestart.h"
113 #include "gromacs/mdrunutility/logging.h"
114 #include "gromacs/mdrunutility/multisim.h"
115 #include "gromacs/mdrunutility/printtime.h"
116 #include "gromacs/mdrunutility/threadaffinity.h"
117 #include "gromacs/mdtypes/checkpointdata.h"
118 #include "gromacs/mdtypes/commrec.h"
119 #include "gromacs/mdtypes/enerdata.h"
120 #include "gromacs/mdtypes/fcdata.h"
121 #include "gromacs/mdtypes/forcerec.h"
122 #include "gromacs/mdtypes/group.h"
123 #include "gromacs/mdtypes/inputrec.h"
124 #include "gromacs/mdtypes/interaction_const.h"
125 #include "gromacs/mdtypes/md_enums.h"
126 #include "gromacs/mdtypes/mdatom.h"
127 #include "gromacs/mdtypes/mdrunoptions.h"
128 #include "gromacs/mdtypes/observableshistory.h"
129 #include "gromacs/mdtypes/simulation_workload.h"
130 #include "gromacs/mdtypes/state.h"
131 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
132 #include "gromacs/modularsimulator/modularsimulator.h"
133 #include "gromacs/nbnxm/gpu_data_mgmt.h"
134 #include "gromacs/nbnxm/nbnxm.h"
135 #include "gromacs/nbnxm/pairlist_tuning.h"
136 #include "gromacs/pbcutil/pbc.h"
137 #include "gromacs/pulling/output.h"
138 #include "gromacs/pulling/pull.h"
139 #include "gromacs/pulling/pull_rotation.h"
140 #include "gromacs/restraint/manager.h"
141 #include "gromacs/restraint/restraintmdmodule.h"
142 #include "gromacs/restraint/restraintpotential.h"
143 #include "gromacs/swap/swapcoords.h"
144 #include "gromacs/taskassignment/decidegpuusage.h"
145 #include "gromacs/taskassignment/decidesimulationworkload.h"
146 #include "gromacs/taskassignment/resourcedivision.h"
147 #include "gromacs/taskassignment/taskassignment.h"
148 #include "gromacs/taskassignment/usergpuids.h"
149 #include "gromacs/timing/gpu_timing.h"
150 #include "gromacs/timing/wallcycle.h"
151 #include "gromacs/timing/wallcyclereporting.h"
152 #include "gromacs/topology/mtop_util.h"
153 #include "gromacs/trajectory/trajectoryframe.h"
154 #include "gromacs/utility/basenetwork.h"
155 #include "gromacs/utility/cstringutil.h"
156 #include "gromacs/utility/exceptions.h"
157 #include "gromacs/utility/fatalerror.h"
158 #include "gromacs/utility/filestream.h"
159 #include "gromacs/utility/gmxassert.h"
160 #include "gromacs/utility/gmxmpi.h"
161 #include "gromacs/utility/keyvaluetree.h"
162 #include "gromacs/utility/logger.h"
163 #include "gromacs/utility/loggerbuilder.h"
164 #include "gromacs/utility/mdmodulesnotifiers.h"
165 #include "gromacs/utility/physicalnodecommunicator.h"
166 #include "gromacs/utility/pleasecite.h"
167 #include "gromacs/utility/programcontext.h"
168 #include "gromacs/utility/smalloc.h"
169 #include "gromacs/utility/stringutil.h"
170 #include "gromacs/utility/mpiinfo.h"
172 #include "isimulator.h"
173 #include "membedholder.h"
174 #include "replicaexchange.h"
175 #include "simulatorbuilder.h"
181 /*! \brief Manage any development feature flag variables encountered
183 * The use of dev features indicated by environment variables is
184 * logged in order to ensure that runs with such features enabled can
185 * be identified from their log and standard output. Any cross
186 * dependencies are also checked, and if unsatisfied, a fatal error
189 * Note that some development features overrides are applied already here:
190 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
192 * \param[in] mdlog Logger object.
193 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
194 * \param[in] pmeRunMode The PME run mode for this run
195 * \returns The object populated with development feature flags.
197 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
198 const bool useGpuForNonbonded,
199 const PmeRunMode pmeRunMode)
201 DevelopmentFeatureFlags devFlags;
203 // Some builds of GCC 5 give false positive warnings that these
204 // getenv results are ignored when clearly they are used.
205 #pragma GCC diagnostic push
206 #pragma GCC diagnostic ignored "-Wunused-result"
208 devFlags.enableGpuBufferOps =
209 GMX_GPU_CUDA && useGpuForNonbonded && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
210 devFlags.enableGpuHaloExchange = GMX_GPU_CUDA && getenv("GMX_GPU_DD_COMMS") != nullptr;
211 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr) || GMX_FAHCORE;
212 devFlags.enableGpuPmePPComm = GMX_GPU_CUDA && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
214 #pragma GCC diagnostic pop
216 // Direct GPU comm path is being used with CUDA_AWARE_MPI
217 // make sure underlying MPI implementation is CUDA-aware
218 if (!GMX_THREAD_MPI && (devFlags.enableGpuPmePPComm || devFlags.enableGpuHaloExchange))
220 const bool haveDetectedCudaAwareMpi =
221 (checkMpiCudaAwareSupport() == CudaAwareMpiStatus::Supported);
222 const bool forceCudaAwareMpi = (getenv("GMX_FORCE_CUDA_AWARE_MPI") != nullptr);
224 if (!haveDetectedCudaAwareMpi && forceCudaAwareMpi)
226 // CUDA-aware support not detected in MPI library but, user has forced it's use
227 GMX_LOG(mdlog.warning)
229 .appendTextFormatted(
230 "This run has forced use of 'CUDA-aware MPI'. "
231 "But, GROMACS cannot determine if underlying MPI "
232 "is CUDA-aware. GROMACS recommends use of latest openMPI version "
233 "for CUDA-aware support. "
234 "If you observe failures at runtime, try unsetting "
235 "GMX_FORCE_CUDA_AWARE_MPI environment variable.");
238 if (haveDetectedCudaAwareMpi || forceCudaAwareMpi)
240 devFlags.usingCudaAwareMpi = true;
241 GMX_LOG(mdlog.warning)
243 .appendTextFormatted(
244 "Using CUDA-aware MPI for 'GPU halo exchange' or 'GPU PME-PP "
245 "communications' feature.");
249 if (devFlags.enableGpuHaloExchange)
251 GMX_LOG(mdlog.warning)
253 .appendTextFormatted(
254 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
255 "halo exchange' feature will not be enabled as GROMACS couldn't "
256 "detect CUDA_aware support in underlying MPI implementation.");
257 devFlags.enableGpuHaloExchange = false;
259 if (devFlags.enableGpuPmePPComm)
261 GMX_LOG(mdlog.warning)
264 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
265 "'GPU PME-PP communications' feature will not be enabled as "
267 "detect CUDA_aware support in underlying MPI implementation.");
268 devFlags.enableGpuPmePPComm = false;
271 GMX_LOG(mdlog.warning)
273 .appendTextFormatted(
274 "GROMACS recommends use of latest OpenMPI version for CUDA-aware "
276 "If you are certain about CUDA-aware support in your MPI library, "
277 "you can force it's use by setting environment variable "
278 " GMX_FORCE_CUDA_AWARE_MPI.");
282 if (devFlags.enableGpuBufferOps)
284 GMX_LOG(mdlog.warning)
286 .appendTextFormatted(
287 "This run uses the 'GPU buffer ops' feature, enabled by the "
288 "GMX_USE_GPU_BUFFER_OPS environment variable.");
291 if (devFlags.forceGpuUpdateDefault)
293 GMX_LOG(mdlog.warning)
295 .appendTextFormatted(
296 "This run will default to '-update gpu' as requested by the "
297 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
298 "decomposition lacks substantial testing and should be used with caution.");
301 if (devFlags.enableGpuHaloExchange)
303 if (useGpuForNonbonded)
305 if (!devFlags.enableGpuBufferOps)
307 GMX_LOG(mdlog.warning)
309 .appendTextFormatted(
310 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
311 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
312 devFlags.enableGpuBufferOps = true;
314 GMX_LOG(mdlog.warning)
316 .appendTextFormatted(
317 "This run has requested the 'GPU halo exchange' feature, enabled by "
319 "GMX_GPU_DD_COMMS environment variable.");
323 GMX_LOG(mdlog.warning)
325 .appendTextFormatted(
326 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
327 "halo exchange' feature will not be enabled as nonbonded interactions "
328 "are not offloaded.");
329 devFlags.enableGpuHaloExchange = false;
333 if (devFlags.enableGpuPmePPComm)
335 if (pmeRunMode == PmeRunMode::GPU)
337 if (!devFlags.enableGpuBufferOps)
339 GMX_LOG(mdlog.warning)
341 .appendTextFormatted(
342 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
343 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
344 devFlags.enableGpuBufferOps = true;
346 GMX_LOG(mdlog.warning)
348 .appendTextFormatted(
349 "This run uses the 'GPU PME-PP communications' feature, enabled "
350 "by the GMX_GPU_PME_PP_COMMS environment variable.");
354 std::string clarification;
355 if (pmeRunMode == PmeRunMode::Mixed)
358 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
363 clarification = "PME is not offloaded to the GPU.";
365 GMX_LOG(mdlog.warning)
368 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
369 "'GPU PME-PP communications' feature was not enabled as "
371 devFlags.enableGpuPmePPComm = false;
378 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
380 * Used to ensure that the master thread does not modify mdrunner during copy
381 * on the spawned threads. */
382 static void threadMpiMdrunnerAccessBarrier()
385 MPI_Barrier(MPI_COMM_WORLD);
389 Mdrunner Mdrunner::cloneOnSpawnedThread() const
391 auto newRunner = Mdrunner(std::make_unique<MDModules>());
393 // All runners in the same process share a restraint manager resource because it is
394 // part of the interface to the client code, which is associated only with the
395 // original thread. Handles to the same resources can be obtained by copy.
397 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
400 // Copy members of master runner.
401 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
402 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
403 newRunner.hw_opt = hw_opt;
404 newRunner.filenames = filenames;
406 newRunner.hwinfo_ = hwinfo_;
407 newRunner.oenv = oenv;
408 newRunner.mdrunOptions = mdrunOptions;
409 newRunner.domdecOptions = domdecOptions;
410 newRunner.nbpu_opt = nbpu_opt;
411 newRunner.pme_opt = pme_opt;
412 newRunner.pme_fft_opt = pme_fft_opt;
413 newRunner.bonded_opt = bonded_opt;
414 newRunner.update_opt = update_opt;
415 newRunner.nstlist_cmdline = nstlist_cmdline;
416 newRunner.replExParams = replExParams;
417 newRunner.pforce = pforce;
418 // Give the spawned thread the newly created valid communicator
419 // for the simulation.
420 newRunner.libraryWorldCommunicator = MPI_COMM_WORLD;
421 newRunner.simulationCommunicator = MPI_COMM_WORLD;
423 newRunner.startingBehavior = startingBehavior;
424 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
425 newRunner.inputHolder_ = inputHolder_;
427 threadMpiMdrunnerAccessBarrier();
432 /*! \brief The callback used for running on spawned threads.
434 * Obtains the pointer to the master mdrunner object from the one
435 * argument permitted to the thread-launch API call, copies it to make
436 * a new runner for this thread, reinitializes necessary data, and
437 * proceeds to the simulation. */
438 static void mdrunner_start_fn(const void* arg)
442 const auto* masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
443 /* copy the arg list to make sure that it's thread-local. This
444 doesn't copy pointed-to items, of course; fnm, cr and fplog
445 are reset in the call below, all others should be const. */
446 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
449 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
453 void Mdrunner::spawnThreads(int numThreadsToLaunch)
456 /* now spawn new threads that start mdrunner_start_fn(), while
457 the main thread returns. Thread affinity is handled later. */
458 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn, static_cast<const void*>(this))
461 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
464 // Give the master thread the newly created valid communicator for
466 libraryWorldCommunicator = MPI_COMM_WORLD;
467 simulationCommunicator = MPI_COMM_WORLD;
468 threadMpiMdrunnerAccessBarrier();
470 GMX_UNUSED_VALUE(numThreadsToLaunch);
471 GMX_UNUSED_VALUE(mdrunner_start_fn);
477 /*! \brief Initialize variables for Verlet scheme simulation */
478 static void prepare_verlet_scheme(FILE* fplog,
482 const gmx_mtop_t& mtop,
484 bool makeGpuPairList,
485 const gmx::CpuInfo& cpuinfo)
487 // We checked the cut-offs in grompp, but double-check here.
488 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
489 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == VanDerWaalsType::Cut)
491 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
492 "With Verlet lists and PME we should have rcoulomb>=rvdw");
496 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
497 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
499 /* For NVE simulations, we will retain the initial list buffer */
500 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0
501 && !(EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No))
503 /* Update the Verlet buffer size for the current run setup */
505 /* Here we assume SIMD-enabled kernels are being used. But as currently
506 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
507 * and 4x2 gives a larger buffer than 4x4, this is ok.
509 ListSetupType listType =
510 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
511 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
513 const real rlist_new =
514 calcVerletBufferSize(mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
516 if (rlist_new != ir->rlist)
518 if (fplog != nullptr)
521 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
524 listSetup.cluster_size_i,
525 listSetup.cluster_size_j);
527 ir->rlist = rlist_new;
531 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
534 "Can not set nstlist without %s",
535 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
538 if (EI_DYNAMICS(ir->eI))
540 /* Set or try nstlist values */
541 increaseNstlist(fplog, cr, ir, nstlist_cmdline, &mtop, box, makeGpuPairList, cpuinfo);
545 /*! \brief Override the nslist value in inputrec
547 * with value passed on the command line (if any)
549 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
553 /* override with anything else than the default -2 */
554 if (nsteps_cmdline > -2)
556 char sbuf_steps[STEPSTRSIZE];
557 char sbuf_msg[STRLEN];
559 ir->nsteps = nsteps_cmdline;
560 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
563 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
564 gmx_step_str(nsteps_cmdline, sbuf_steps),
565 fabs(nsteps_cmdline * ir->delta_t));
570 "Overriding nsteps with value passed on the command line: %s steps",
571 gmx_step_str(nsteps_cmdline, sbuf_steps));
574 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
576 else if (nsteps_cmdline < -2)
578 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
580 /* Do nothing if nsteps_cmdline == -2 */
586 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
588 * If not, and if a warning may be issued, logs a warning about
589 * falling back to CPU code. With thread-MPI, only the first
590 * call to this function should have \c issueWarning true. */
591 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
593 bool gpuIsUseful = true;
596 if (ir.opts.ngener - ir.nwall > 1)
598 /* The GPU code does not support more than one energy group.
599 * If the user requested GPUs explicitly, a fatal error is given later.
603 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
604 "For better performance, run on the GPU without energy groups and then do "
605 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
611 warning = "TPI is not implemented for GPUs.";
614 if (!gpuIsUseful && issueWarning)
616 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
622 //! Initializes the logger for mdrun.
623 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
625 gmx::LoggerBuilder builder;
626 if (fplog != nullptr)
628 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
630 if (isSimulationMasterRank)
632 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
634 return builder.build();
637 //! Make a TaskTarget from an mdrun argument string.
638 static TaskTarget findTaskTarget(const char* optionString)
640 TaskTarget returnValue = TaskTarget::Auto;
642 if (strncmp(optionString, "auto", 3) == 0)
644 returnValue = TaskTarget::Auto;
646 else if (strncmp(optionString, "cpu", 3) == 0)
648 returnValue = TaskTarget::Cpu;
650 else if (strncmp(optionString, "gpu", 3) == 0)
652 returnValue = TaskTarget::Gpu;
656 GMX_ASSERT(false, "Option string should have been checked for sanity already");
662 //! Finish run, aggregate data to print performance info.
663 static void finish_run(FILE* fplog,
664 const gmx::MDLogger& mdlog,
666 const t_inputrec& inputrec,
668 gmx_wallcycle* wcycle,
669 gmx_walltime_accounting_t walltime_accounting,
670 nonbonded_verlet_t* nbv,
671 const gmx_pme_t* pme,
675 double nbfs = 0, mflop = 0;
676 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
677 elapsed_time_over_all_threads_over_all_ranks;
678 /* Control whether it is valid to print a report. Only the
679 simulation master may print, but it should not do so if the run
680 terminated e.g. before a scheduled reset step. This is
681 complicated by the fact that PME ranks are unaware of the
682 reason why they were sent a pmerecvqxFINISH. To avoid
683 communication deadlocks, we always do the communication for the
684 report, even if we've decided not to write the report, because
685 how long it takes to finish the run is not important when we've
686 decided not to report on the simulation performance.
688 Further, we only report performance for dynamical integrators,
689 because those are the only ones for which we plan to
690 consider doing any optimizations. */
691 bool printReport = EI_DYNAMICS(inputrec.eI) && SIMMASTER(cr);
693 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
695 GMX_LOG(mdlog.warning)
697 .appendText("Simulation ended prematurely, no performance report will be written.");
702 std::unique_ptr<t_nrnb> nrnbTotalStorage;
705 nrnbTotalStorage = std::make_unique<t_nrnb>();
706 nrnb_tot = nrnbTotalStorage.get();
708 MPI_Allreduce(nrnb->n.data(), nrnb_tot->n.data(), eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
716 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
717 elapsed_time_over_all_threads =
718 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
722 /* reduce elapsed_time over all MPI ranks in the current simulation */
723 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
724 elapsed_time_over_all_ranks /= cr->nnodes;
725 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
726 * current simulation. */
727 MPI_Allreduce(&elapsed_time_over_all_threads,
728 &elapsed_time_over_all_threads_over_all_ranks,
737 elapsed_time_over_all_ranks = elapsed_time;
738 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
743 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
746 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
748 print_dd_statistics(cr, inputrec, fplog);
751 /* TODO Move the responsibility for any scaling by thread counts
752 * to the code that handled the thread region, so that there's a
753 * mechanism to keep cycle counting working during the transition
754 * to task parallelism. */
755 int nthreads_pp = gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded);
756 int nthreads_pme = gmx_omp_nthreads_get(ModuleMultiThread::Pme);
757 wallcycle_scale_by_num_threads(
758 wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
759 auto cycle_sum(wallcycle_sum(cr, wcycle));
763 auto* nbnxn_gpu_timings =
764 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
765 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
767 if (pme_gpu_task_enabled(pme))
769 pme_gpu_get_timings(pme, &pme_gpu_timings);
771 wallcycle_print(fplog,
777 elapsed_time_over_all_ranks,
783 if (EI_DYNAMICS(inputrec.eI))
785 delta_t = inputrec.delta_t;
791 elapsed_time_over_all_threads_over_all_ranks,
792 elapsed_time_over_all_ranks,
793 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
801 elapsed_time_over_all_threads_over_all_ranks,
802 elapsed_time_over_all_ranks,
803 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
811 int Mdrunner::mdrunner()
814 std::unique_ptr<t_forcerec> fr;
815 real ewaldcoeff_q = 0;
816 real ewaldcoeff_lj = 0;
817 int nChargePerturbed = -1, nTypePerturbed = 0;
818 gmx_walltime_accounting_t walltime_accounting = nullptr;
819 MembedHolder membedHolder(filenames.size(), filenames.data());
821 /* CAUTION: threads may be started later on in this function, so
822 cr doesn't reflect the final parallel state right now */
825 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
826 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
827 const bool doRerun = mdrunOptions.rerun;
829 // Handle task-assignment related user options.
830 EmulateGpuNonbonded emulateGpuNonbonded =
831 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
833 std::vector<int> userGpuTaskAssignment;
836 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
838 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
839 auto nonbondedTarget = findTaskTarget(nbpu_opt);
840 auto pmeTarget = findTaskTarget(pme_opt);
841 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
842 auto bondedTarget = findTaskTarget(bonded_opt);
843 auto updateTarget = findTaskTarget(update_opt);
845 FILE* fplog = nullptr;
846 // If we are appending, we don't write log output because we need
847 // to check that the old log file matches what the checkpoint file
848 // expects. Otherwise, we should start to write log output now if
849 // there is a file ready for it.
850 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
852 fplog = gmx_fio_getfp(logFileHandle);
854 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, simulationCommunicator);
855 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
856 gmx::MDLogger mdlog(logOwner.logger());
858 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo_);
860 std::vector<int> availableDevices =
861 makeListOfAvailableDevices(hwinfo_->deviceInfoList, hw_opt.devicesSelectedByUser);
862 const int numAvailableDevices = gmx::ssize(availableDevices);
864 // Print citation requests after all software/hardware printing
865 pleaseCiteGromacs(fplog);
867 // Note: legacy program logic relies on checking whether these pointers are assigned.
868 // Objects may or may not be allocated later.
869 std::unique_ptr<t_inputrec> inputrec;
870 std::unique_ptr<t_state> globalState;
872 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
874 if (isSimulationMasterRank)
876 // Allocate objects to be initialized by later function calls.
877 /* Only the master rank has the global state */
878 globalState = std::make_unique<t_state>();
879 inputrec = std::make_unique<t_inputrec>();
881 /* Read (nearly) all data required for the simulation
882 * and keep the partly serialized tpr contents to send to other ranks later
884 applyGlobalSimulationState(
885 *inputHolder_.get(), partialDeserializedTpr.get(), globalState.get(), inputrec.get(), &mtop);
888 /* Check and update the hardware options for internal consistency */
889 checkAndUpdateHardwareOptions(
890 mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks, inputrec.get());
892 if (GMX_THREAD_MPI && isSimulationMasterRank)
894 bool useGpuForNonbonded = false;
895 bool useGpuForPme = false;
898 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
900 // If the user specified the number of ranks, then we must
901 // respect that, but in default mode, we need to allow for
902 // the number of GPUs to choose the number of ranks.
903 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
904 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
906 numAvailableDevices > 0,
907 userGpuTaskAssignment,
909 canUseGpuForNonbonded,
910 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
911 hw_opt.nthreads_tmpi);
912 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(useGpuForNonbonded,
915 userGpuTaskAssignment,
918 hw_opt.nthreads_tmpi,
919 domdecOptions.numPmeRanks);
921 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
923 /* Determine how many thread-MPI ranks to start.
925 * TODO Over-writing the user-supplied value here does
926 * prevent any possible subsequent checks from working
928 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo_,
936 membedHolder.doMembed());
938 // Now start the threads for thread MPI.
939 spawnThreads(hw_opt.nthreads_tmpi);
940 // The spawned threads enter mdrunner() and execution of
941 // master and spawned threads joins at the end of this block.
944 GMX_RELEASE_ASSERT(ms || simulationCommunicator != MPI_COMM_NULL,
945 "Must have valid communicator unless running a multi-simulation");
946 CommrecHandle crHandle = init_commrec(simulationCommunicator);
947 t_commrec* cr = crHandle.get();
948 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
950 PhysicalNodeCommunicator physicalNodeComm(libraryWorldCommunicator, gmx_physicalnode_id_hash());
952 // If we detected the topology on this system, double-check that it makes sense
953 if (hwinfo_->hardwareTopology->isThisSystem())
955 hardwareTopologyDoubleCheckDetection(mdlog, *hwinfo_->hardwareTopology);
960 /* now broadcast everything to the non-master nodes/threads: */
961 if (!isSimulationMasterRank)
963 // Until now, only the master rank has a non-null pointer.
964 // On non-master ranks, allocate the object that will receive data in the following call.
965 inputrec = std::make_unique<t_inputrec>();
967 init_parallel(cr->mpiDefaultCommunicator,
971 partialDeserializedTpr.get());
973 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
974 partialDeserializedTpr.reset(nullptr);
977 !inputrec->useConstantAcceleration,
978 "Linear acceleration has been removed in GROMACS 2022, and was broken for many years "
979 "before that. Use GROMACS 4.5 or earlier if you need this feature.");
981 // Now the number of ranks is known to all ranks, and each knows
982 // the inputrec read by the master rank. The ranks can now all run
983 // the task-deciding functions and will agree on the result
984 // without needing to communicate.
985 const bool useDomainDecomposition =
986 (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == IntegrationAlgorithm::NM));
988 // Note that these variables describe only their own node.
990 // Note that when bonded interactions run on a GPU they always run
991 // alongside a nonbonded task, so do not influence task assignment
992 // even though they affect the force calculation workload.
993 bool useGpuForNonbonded = false;
994 bool useGpuForPme = false;
995 bool useGpuForBonded = false;
996 bool useGpuForUpdate = false;
997 bool gpusWereDetected = hwinfo_->ngpu_compatible_tot > 0;
1000 // It's possible that there are different numbers of GPUs on
1001 // different nodes, which is the user's responsibility to
1002 // handle. If unsuitable, we will notice that during task
1004 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
1005 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
1007 userGpuTaskAssignment,
1008 emulateGpuNonbonded,
1009 canUseGpuForNonbonded,
1010 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI),
1012 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded,
1014 userGpuTaskAssignment,
1017 cr->sizeOfDefaultCommunicator,
1018 domdecOptions.numPmeRanks,
1020 useGpuForBonded = decideWhetherToUseGpusForBonded(
1021 useGpuForNonbonded, useGpuForPme, bondedTarget, *inputrec, mtop, domdecOptions.numPmeRanks, gpusWereDetected);
1023 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1025 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
1027 // Initialize development feature flags that enabled by environment variable
1028 // and report those features that are enabled.
1029 const DevelopmentFeatureFlags devFlags =
1030 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
1032 const bool useModularSimulator = checkUseModularSimulator(false,
1039 doEssentialDynamics,
1040 membedHolder.doMembed());
1042 // Build restraints.
1043 // TODO: hide restraint implementation details from Mdrunner.
1044 // There is nothing unique about restraints at this point as far as the
1045 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
1046 // factory functions from the SimulationContext on which to call mdModules_->add().
1047 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
1048 for (auto&& restraint : restraintManager_->getRestraints())
1050 auto module = RestraintMDModule::create(restraint, restraint->sites());
1051 mdModules_->add(std::move(module));
1054 // TODO: Error handling
1055 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
1056 // now that the MDModules know their options, they know which callbacks to sign up to
1057 mdModules_->subscribeToSimulationSetupNotifications();
1058 const auto& setupNotifier = mdModules_->notifiers().simulationSetupNotifier_;
1060 if (inputrec->internalParameters != nullptr)
1062 setupNotifier.notify(*inputrec->internalParameters);
1065 if (fplog != nullptr)
1067 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
1068 fprintf(fplog, "\n");
1073 /* In rerun, set velocities to zero if present */
1074 if (doRerun && ((globalState->flags & enumValueToBitMask(StateEntry::V)) != 0))
1076 // rerun does not use velocities
1080 "Rerun trajectory contains velocities. Rerun does only evaluate "
1081 "potential energy and forces. The velocities will be ignored.");
1082 for (int i = 0; i < globalState->natoms; i++)
1084 clear_rvec(globalState->v[i]);
1086 globalState->flags &= ~enumValueToBitMask(StateEntry::V);
1089 /* now make sure the state is initialized and propagated */
1090 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
1093 /* NM and TPI parallelize over force/energy calculations, not atoms,
1094 * so we need to initialize and broadcast the global state.
1096 if (inputrec->eI == IntegrationAlgorithm::NM || inputrec->eI == IntegrationAlgorithm::TPI)
1100 globalState = std::make_unique<t_state>();
1102 broadcastStateWithoutDynamics(
1103 cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr), globalState.get());
1106 /* A parallel command line option consistency check that we can
1107 only do after any threads have started. */
1109 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
1110 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
1113 "The -dd or -npme option request a parallel simulation, "
1115 "but %s was compiled without threads or MPI enabled",
1116 output_env_get_program_display_name(oenv));
1117 #elif GMX_THREAD_MPI
1118 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
1120 "but %s was not started through mpirun/mpiexec or only one rank was requested "
1121 "through mpirun/mpiexec",
1122 output_env_get_program_display_name(oenv));
1126 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || IntegrationAlgorithm::NM == inputrec->eI))
1129 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1130 "these are not compatible with mdrun -rerun");
1133 /* Object for collecting reasons for not using PME-only ranks */
1134 SeparatePmeRanksPermitted separatePmeRanksPermitted;
1136 /* Permit MDModules to notify whether they want to use PME-only ranks */
1137 setupNotifier.notify(&separatePmeRanksPermitted);
1139 /* If simulation is not using PME then disable PME-only ranks */
1140 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1142 separatePmeRanksPermitted.disablePmeRanks(
1143 "PME-only ranks are requested, but the system does not use PME "
1144 "for electrostatics or LJ");
1147 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1148 * improve performance with many threads per GPU, since our OpenMP
1149 * scaling is bad, but it's difficult to automate the setup.
1151 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1153 separatePmeRanksPermitted.disablePmeRanks(
1154 "PME-only CPU ranks are not automatically used when "
1155 "non-bonded interactions are computed on GPUs");
1158 /* If GPU is used for PME then only 1 PME rank is permitted */
1159 if (useGpuForPme && (domdecOptions.numPmeRanks < 0 || domdecOptions.numPmeRanks > 1))
1161 separatePmeRanksPermitted.disablePmeRanks(
1162 "PME GPU decomposition is not supported. Only one separate PME-only GPU rank "
1166 /* Disable PME-only ranks if some parts of the code requested so and it's up to GROMACS to decide */
1167 if (!separatePmeRanksPermitted.permitSeparatePmeRanks() && domdecOptions.numPmeRanks < 0)
1169 domdecOptions.numPmeRanks = 0;
1172 .appendText("Simulation will not use PME-only ranks because: "
1173 + separatePmeRanksPermitted.reasonsWhyDisabled());
1176 /* If some parts of the code could not use PME-only ranks and
1177 * user explicitly used mdrun -npme option then throw an error */
1178 if (!separatePmeRanksPermitted.permitSeparatePmeRanks() && domdecOptions.numPmeRanks > 0)
1180 gmx_fatal_collective(FARGS,
1181 cr->mpiDefaultCommunicator,
1183 "Requested -npme %d option is not viable because: %s",
1184 domdecOptions.numPmeRanks,
1185 separatePmeRanksPermitted.reasonsWhyDisabled().c_str());
1188 /* NMR restraints must be initialized before load_checkpoint,
1189 * since with time averaging the history is added to t_state.
1190 * For proper consistency check we therefore need to extend
1192 * So the PME-only nodes (if present) will also initialize
1193 * the distance restraints.
1196 /* This needs to be called before read_checkpoint to extend the state */
1197 t_disresdata* disresdata;
1198 snew(disresdata, 1);
1202 DisResRunMode::MDRun,
1203 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1204 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1209 replExParams.exchangeInterval > 0);
1211 std::unique_ptr<t_oriresdata> oriresData;
1212 if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1214 oriresData = std::make_unique<t_oriresdata>(fplog, mtop, *inputrec, cr, ms, globalState.get());
1217 auto deform = prepareBoxDeformation(globalState != nullptr ? globalState->box : box,
1218 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1219 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1220 cr->mpi_comm_mygroup,
1224 /* We have to remember the generation's first step before reading checkpoint.
1225 This way, we can report to the F@H core both the generation's first step
1226 and the restored first step, thus making it able to distinguish between
1227 an interruption/resume and start of the n-th generation simulation.
1228 Having this information, the F@H core can correctly calculate and report
1231 int gen_first_step = 0;
1234 gen_first_step = inputrec->init_step;
1238 ObservablesHistory observablesHistory = {};
1240 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1241 if (startingBehavior != StartingBehavior::NewSimulation)
1243 /* Check if checkpoint file exists before doing continuation.
1244 * This way we can use identical input options for the first and subsequent runs...
1246 if (mdrunOptions.numStepsCommandline > -2)
1248 /* Temporarily set the number of steps to unlimited to avoid
1249 * triggering the nsteps check in load_checkpoint().
1250 * This hack will go away soon when the -nsteps option is removed.
1252 inputrec->nsteps = -1;
1255 // Finish applying initial simulation state information from external sources on all ranks.
1256 // Reconcile checkpoint file data with Mdrunner state established up to this point.
1257 applyLocalState(*inputHolder_.get(),
1260 domdecOptions.numCells,
1263 &observablesHistory,
1264 mdrunOptions.reproducible,
1265 mdModules_->notifiers(),
1266 modularSimulatorCheckpointData.get(),
1267 useModularSimulator);
1268 // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1270 // on all code paths.
1271 // Write checkpoint or provide hook to update SimulationInput.
1272 // If there was a checkpoint file, SimulationInput contains more information
1273 // than if there wasn't. At this point, we have synchronized the in-memory
1274 // state with the filesystem state only for restarted simulations. We should
1275 // be calling applyLocalState unconditionally and expect that the completeness
1276 // of SimulationInput is not dependent on its creation method.
1278 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1280 // Now we can start normal logging to the truncated log file.
1281 fplog = gmx_fio_getfp(logFileHandle);
1282 prepareLogAppending(fplog);
1283 logOwner = buildLogger(fplog, MASTER(cr));
1284 mdlog = logOwner.logger();
1291 fcRegisterSteps(inputrec->nsteps + inputrec->init_step, gen_first_step);
1295 if (mdrunOptions.numStepsCommandline > -2)
1300 "The -nsteps functionality is deprecated, and may be removed in a future "
1302 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1305 /* override nsteps with value set on the commandline */
1306 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1308 if (isSimulationMasterRank)
1310 copy_mat(globalState->box, box);
1315 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1318 if (inputrec->cutoff_scheme != CutoffScheme::Verlet)
1321 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1322 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1324 /* Update rlist and nstlist. */
1325 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1326 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1327 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1329 prepare_verlet_scheme(fplog,
1335 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1338 // This builder is necessary while we have multi-part construction
1339 // of DD. Before DD is constructed, we use the existence of
1340 // the builder object to indicate that further construction of DD
1342 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1343 if (useDomainDecomposition)
1345 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1353 positionsFromStatePointer(globalState.get()));
1357 /* PME, if used, is done on all nodes with 1D decomposition */
1358 cr->nnodes = cr->sizeOfDefaultCommunicator;
1359 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1360 cr->nodeid = cr->rankInDefaultCommunicator;
1362 cr->duty = (DUTY_PP | DUTY_PME);
1364 if (inputrec->pbcType == PbcType::Screw)
1366 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1370 // Produce the task assignment for this rank - done after DD is constructed
1371 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1373 userGpuTaskAssignment,
1375 simulationCommunicator,
1383 thisRankHasDuty(cr, DUTY_PP),
1384 // TODO cr->duty & DUTY_PME should imply that a PME
1385 // algorithm is active, but currently does not.
1386 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1388 // Get the device handles for the modules, nullptr when no task is assigned.
1390 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1392 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1393 bool useTiming = true;
1397 /* WARNING: CUDA timings are incorrect with multiple streams.
1398 * This is the main reason why they are disabled by default.
1400 // TODO: Consider turning on by default when we can detect nr of streams.
1401 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1403 else if (GMX_GPU_OPENCL)
1405 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1408 // TODO Currently this is always built, yet DD partition code
1409 // checks if it is built before using it. Probably it should
1410 // become an MDModule that is made only when another module
1411 // requires it (e.g. pull, CompEl, density fitting), so that we
1412 // don't update the local atom sets unilaterally every step.
1413 LocalAtomSetManager atomSets;
1416 // TODO Pass the GPU streams to ddBuilder to use in buffer
1417 // transfers (e.g. halo exchange)
1418 cr->dd = ddBuilder->build(&atomSets);
1419 // The builder's job is done, so destruct it
1420 ddBuilder.reset(nullptr);
1421 // Note that local state still does not exist yet.
1424 // The GPU update is decided here because we need to know whether the constraints or
1425 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1426 // defined). This is only known after DD is initialized, hence decision on using GPU
1427 // update is done so late.
1430 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1431 const bool haveFrozenAtoms = inputrecFrozenAtoms(inputrec.get());
1433 useGpuForUpdate = decideWhetherToUseGpuForUpdate(useDomainDecomposition,
1436 domdecOptions.numPmeRanks > 0,
1442 doEssentialDynamics,
1443 gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1444 replExParams.exchangeInterval > 0,
1450 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1452 const bool printHostName = (cr->nnodes > 1);
1453 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1455 const bool disableNonbondedCalculation = (getenv("GMX_NO_NONBONDED") != nullptr);
1456 if (disableNonbondedCalculation)
1458 /* turn off non-bonded calculations */
1459 GMX_LOG(mdlog.warning)
1462 "Found environment variable GMX_NO_NONBONDED.\n"
1463 "Disabling nonbonded calculations.");
1466 MdrunScheduleWorkload runScheduleWork;
1468 bool useGpuDirectHalo = decideWhetherToUseGpuForHalo(devFlags,
1469 havePPDomainDecomposition(cr),
1471 useModularSimulator,
1473 EI_ENERGY_MINIMIZATION(inputrec->eI));
1475 // Also populates the simulation constant workload description.
1476 runScheduleWork.simulationWork = createSimulationWorkload(*inputrec,
1477 disableNonbondedCalculation,
1485 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1487 if (deviceInfo != nullptr)
1489 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1491 dd_setup_dlb_resource_sharing(cr, deviceId);
1493 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1494 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, useTiming);
1497 // If the user chose a task assignment, give them some hints
1498 // where appropriate.
1499 if (!userGpuTaskAssignment.empty())
1501 gpuTaskAssignments.logPerformanceHints(mdlog, numAvailableDevices);
1506 /* After possible communicator splitting in make_dd_communicators.
1507 * we can set up the intra/inter node communication.
1509 gmx_setup_nodecomm(fplog, cr);
1515 GMX_LOG(mdlog.warning)
1517 .appendTextFormatted(
1518 "This is simulation %d out of %d running as a composite GROMACS\n"
1519 "multi-simulation job. Setup for this simulation:\n",
1520 ms->simulationIndex_,
1521 ms->numSimulations_);
1523 GMX_LOG(mdlog.warning)
1524 .appendTextFormatted("Using %d MPI %s\n",
1527 cr->nnodes == 1 ? "thread" : "threads"
1529 cr->nnodes == 1 ? "process" : "processes"
1535 // If mdrun -pin auto honors any affinity setting that already
1536 // exists. If so, it is nice to provide feedback about whether
1537 // that existing affinity setting was from OpenMP or something
1538 // else, so we run this code both before and after we initialize
1539 // the OpenMP support.
1540 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, FALSE);
1541 /* Check and update the number of OpenMP threads requested */
1542 checkAndUpdateRequestedNumOpenmpThreads(
1543 &hw_opt, *hwinfo_, cr, ms, physicalNodeComm.size_, pmeRunMode, mtop, *inputrec);
1545 gmx_omp_nthreads_init(mdlog,
1547 hwinfo_->nthreads_hw_avail,
1548 physicalNodeComm.size_,
1549 hw_opt.nthreads_omp,
1550 hw_opt.nthreads_omp_pme,
1551 !thisRankHasDuty(cr, DUTY_PP));
1553 const bool bEnableFPE = gmxShouldEnableFPExceptions();
1554 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1557 gmx_feenableexcept();
1560 /* Now that we know the setup is consistent, check for efficiency */
1561 check_resource_division_efficiency(
1562 hwinfo_, gpuTaskAssignments.thisRankHasAnyGpuTask(), mdrunOptions.ntompOptionIsSet, cr, mdlog);
1564 /* getting number of PP/PME threads on this MPI / tMPI rank.
1565 PME: env variable should be read only on one node to make sure it is
1566 identical everywhere;
1568 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP)
1569 ? gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded)
1570 : gmx_omp_nthreads_get(ModuleMultiThread::Pme);
1571 checkHardwareOversubscription(
1572 numThreadsOnThisRank, cr->nodeid, *hwinfo_->hardwareTopology, physicalNodeComm, mdlog);
1574 // Enable Peer access between GPUs where available
1575 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1576 // any of the GPU communication features are active.
1577 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1578 && (runScheduleWork.simulationWork.useGpuHaloExchange
1579 || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1581 setupGpuDevicePeerAccess(gpuTaskAssignments.deviceIdsAssigned(), mdlog);
1584 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1586 /* Before setting affinity, check whether the affinity has changed
1587 * - which indicates that probably the OpenMP library has changed it
1588 * since we first checked).
1590 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, TRUE);
1592 int numThreadsOnThisNode, intraNodeThreadOffset;
1593 analyzeThreadsOnThisNode(
1594 physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode, &intraNodeThreadOffset);
1596 /* Set the CPU affinity */
1597 gmx_set_thread_affinity(mdlog,
1600 *hwinfo_->hardwareTopology,
1601 numThreadsOnThisRank,
1602 numThreadsOnThisNode,
1603 intraNodeThreadOffset,
1607 if (mdrunOptions.timingOptions.resetStep > -1)
1612 "The -resetstep functionality is deprecated, and may be removed in a "
1615 std::unique_ptr<gmx_wallcycle> wcycle =
1616 wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1620 /* Master synchronizes its value of reset_counters with all nodes
1621 * including PME only nodes */
1622 int64_t reset_counters = wcycle_get_reset_counters(wcycle.get());
1623 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1624 wcycle_set_reset_counters(wcycle.get(), reset_counters);
1627 // Membrane embedding must be initialized before we call init_forcerec()
1628 membedHolder.initializeMembed(fplog,
1635 &mdrunOptions.checkpointOptions.period);
1637 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1638 std::unique_ptr<MDAtoms> mdAtoms;
1639 std::unique_ptr<VirtualSitesHandler> vsite;
1640 std::unique_ptr<GpuBonded> gpuBonded;
1643 if (thisRankHasDuty(cr, DUTY_PP))
1645 setupNotifier.notify(*cr);
1646 setupNotifier.notify(&atomSets);
1647 setupNotifier.notify(mtop);
1648 setupNotifier.notify(inputrec->pbcType);
1649 setupNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1650 /* Initiate forcerecord */
1651 fr = std::make_unique<t_forcerec>();
1652 fr->forceProviders = mdModules_->initForceProviders();
1653 init_forcerec(fplog,
1660 opt2fn("-table", filenames.size(), filenames.data()),
1661 opt2fn("-tablep", filenames.size(), filenames.data()),
1662 opt2fns("-tableb", filenames.size(), filenames.data()),
1664 // Dirty hack, for fixing disres and orires should be made mdmodules
1665 fr->fcdata->disres = disresdata;
1666 fr->fcdata->orires.swap(oriresData);
1668 // Save a handle to device stream manager to use elsewhere in the code
1669 // TODO: Forcerec is not a correct place to store it.
1670 fr->deviceStreamManager = deviceStreamManager.get();
1672 if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1675 deviceStreamManager != nullptr,
1676 "GPU device stream manager should be valid in order to use PME-PP direct "
1679 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1680 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1682 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1685 deviceStreamManager->context(),
1686 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1689 fr->nbv = Nbnxm::init_nb_verlet(mdlog,
1694 runScheduleWork.simulationWork.useGpuNonbonded,
1695 deviceStreamManager.get(),
1699 // TODO: Move the logic below to a GPU bonded builder
1700 if (runScheduleWork.simulationWork.useGpuBonded)
1702 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1703 "GPU device stream manager should be valid in order to use GPU "
1704 "version of bonded forces.");
1705 gpuBonded = std::make_unique<GpuBonded>(
1707 fr->ic->epsfac * fr->fudgeQQ,
1708 deviceStreamManager->context(),
1709 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)),
1711 fr->gpuBonded = gpuBonded.get();
1714 /* Initialize the mdAtoms structure.
1715 * mdAtoms is not filled with atom data,
1716 * as this can not be done now with domain decomposition.
1718 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1719 if (globalState && thisRankHasPmeGpuTask)
1721 // The pinning of coordinates in the global state object works, because we only use
1722 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1723 // points to the global state object without DD.
1724 // FIXME: MD and EM separately set up the local state - this should happen in the same
1725 // function, which should also perform the pinning.
1726 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1729 /* Initialize the virtual site communication */
1730 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingsPerMoleculeType;
1731 if (DOMAINDECOMP(cr))
1733 updateGroupingsPerMoleculeType = getUpdateGroupingsPerMoleculeType(*cr->dd);
1735 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType, updateGroupingsPerMoleculeType);
1737 calc_shifts(box, fr->shift_vec);
1739 /* With periodic molecules the charge groups should be whole at start up
1740 * and the virtual sites should not be far from their proper positions.
1742 if (!inputrec->bContinuation && MASTER(cr)
1743 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1745 /* Make molecules whole at start of run */
1746 if (fr->pbcType != PbcType::No)
1748 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1752 /* Correct initial vsite positions are required
1753 * for the initial distribution in the domain decomposition
1754 * and for the initial shell prediction.
1756 constructVirtualSitesGlobal(mtop, globalState->x);
1759 // Make the DD reverse topology, now that any vsites that are present are available
1760 if (DOMAINDECOMP(cr))
1762 dd_make_reverse_top(fplog, cr->dd, mtop, vsite.get(), *inputrec, domdecOptions.ddBondedChecking);
1765 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1767 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1768 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1773 /* This is a PME only node */
1775 GMX_ASSERT(globalState == nullptr,
1776 "We don't need the state on a PME only rank and expect it to be unitialized");
1778 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1779 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1782 gmx_pme_t* sepPmeData = nullptr;
1783 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1784 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1785 "Double-checking that only PME-only ranks have no forcerec");
1786 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1788 // TODO should live in ewald module once its testing is improved
1790 // Later, this program could contain kernels that might be later
1791 // re-used as auto-tuning progresses, or subsequent simulations
1793 PmeGpuProgramStorage pmeGpuProgram;
1794 if (thisRankHasPmeGpuTask)
1797 (deviceStreamManager != nullptr),
1798 "GPU device stream manager should be initialized in order to use GPU for PME.");
1799 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1800 "GPU device should be initialized in order to use GPU for PME.");
1801 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1804 /* Initiate PME if necessary,
1805 * either on all nodes or on dedicated PME nodes only. */
1806 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1808 if (mdAtoms && mdAtoms->mdatoms())
1810 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1811 if (EVDW_PME(inputrec->vdwtype))
1813 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1816 if (cr->npmenodes > 0)
1818 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1819 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1820 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1823 if (thisRankHasDuty(cr, DUTY_PME))
1827 // TODO: This should be in the builder.
1828 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1829 || (deviceStreamManager != nullptr),
1830 "Device stream manager should be valid in order to use GPU "
1833 !runScheduleWork.simulationWork.useGpuPme
1834 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1835 "GPU PME stream should be valid in order to use GPU version of PME.");
1837 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1838 ? &deviceStreamManager->context()
1840 const DeviceStream* pmeStream =
1841 runScheduleWork.simulationWork.useGpuPme
1842 ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1845 pmedata = gmx_pme_init(cr,
1846 getNumPmeDomains(cr->dd),
1848 nChargePerturbed != 0,
1849 nTypePerturbed != 0,
1850 mdrunOptions.reproducible,
1853 gmx_omp_nthreads_get(ModuleMultiThread::Pme),
1858 pmeGpuProgram.get(),
1861 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1866 if (EI_DYNAMICS(inputrec->eI))
1868 /* Turn on signal handling on all nodes */
1870 * (A user signal from the PME nodes (if any)
1871 * is communicated to the PP nodes.
1873 signal_handler_install();
1876 pull_t* pull_work = nullptr;
1877 if (thisRankHasDuty(cr, DUTY_PP))
1879 /* Assumes uniform use of the number of OpenMP threads */
1880 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(ModuleMultiThread::Default));
1882 if (inputrec->bPull)
1884 /* Initialize pull code */
1885 pull_work = init_pull(fplog,
1886 inputrec->pull.get(),
1891 inputrec->fepvals->init_lambda);
1892 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1894 initPullHistory(pull_work, &observablesHistory);
1896 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1898 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1902 std::unique_ptr<EnforcedRotation> enforcedRotation;
1905 /* Initialize enforced rotation code */
1906 enforcedRotation = init_rot(fplog,
1919 t_swap* swap = nullptr;
1920 if (inputrec->eSwapCoords != SwapType::No)
1922 /* Initialize ion swapping code */
1923 swap = init_swapcoords(fplog,
1925 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1928 &observablesHistory,
1936 /* Let makeConstraints know whether we have essential dynamics constraints. */
1937 auto constr = makeConstraints(mtop,
1940 doEssentialDynamics,
1943 DOMAINDECOMP(cr) && ddMayHaveSplitConstraints(*cr->dd),
1949 /* Energy terms and groups */
1950 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1951 inputrec->fepvals->n_lambda);
1953 // cos acceleration is only supported by md, but older tpr
1954 // files might still combine it with other integrators
1955 GMX_RELEASE_ASSERT(inputrec->cos_accel == 0.0 || inputrec->eI == IntegrationAlgorithm::MD,
1956 "cos_acceleration is only supported by integrator=md");
1958 /* Kinetic energy data */
1959 gmx_ekindata_t ekind(inputrec->opts.ngtc,
1960 inputrec->cos_accel,
1961 gmx_omp_nthreads_get(ModuleMultiThread::Update));
1963 /* Set up interactive MD (IMD) */
1964 auto imdSession = makeImdSession(inputrec.get(),
1971 MASTER(cr) ? globalState->x : gmx::ArrayRef<gmx::RVec>(),
1975 mdrunOptions.imdOptions,
1978 if (DOMAINDECOMP(cr))
1980 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1981 /* This call is not included in init_domain_decomposition
1982 * because fr->cginfo_mb is set later.
1984 makeBondedLinks(cr->dd, mtop, fr->cginfo_mb);
1987 if (runScheduleWork.simulationWork.useGpuBufferOps)
1989 fr->gpuForceReduction[gmx::AtomLocality::Local] = std::make_unique<gmx::GpuForceReduction>(
1990 deviceStreamManager->context(),
1991 deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedLocal),
1993 fr->gpuForceReduction[gmx::AtomLocality::NonLocal] = std::make_unique<gmx::GpuForceReduction>(
1994 deviceStreamManager->context(),
1995 deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedNonLocal),
1999 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
2000 if (gpusWereDetected
2001 && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
2002 || runScheduleWork.simulationWork.useGpuBufferOps))
2004 GpuApiCallBehavior transferKind =
2005 (inputrec->eI == IntegrationAlgorithm::MD && !doRerun && !useModularSimulator)
2006 ? GpuApiCallBehavior::Async
2007 : GpuApiCallBehavior::Sync;
2008 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
2009 "GPU device stream manager should be initialized to use GPU.");
2010 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
2011 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle.get());
2012 fr->stateGpu = stateGpu.get();
2015 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
2016 SimulatorBuilder simulatorBuilder;
2018 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
2019 simulatorBuilder.add(std::move(membedHolder));
2020 simulatorBuilder.add(std::move(stopHandlerBuilder_));
2021 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
2024 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
2025 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle.get()));
2026 simulatorBuilder.add(ConstraintsParam(
2027 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, vsite.get()));
2028 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
2029 simulatorBuilder.add(LegacyInput(
2030 static_cast<int>(filenames.size()), filenames.data(), inputrec.get(), fr.get()));
2031 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
2032 simulatorBuilder.add(InteractiveMD(imdSession.get()));
2033 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifiers()));
2034 simulatorBuilder.add(CenterOfMassPulling(pull_work));
2035 // Todo move to an MDModule
2036 simulatorBuilder.add(IonSwapping(swap));
2037 simulatorBuilder.add(TopologyData(mtop, mdAtoms.get()));
2038 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
2039 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
2041 // build and run simulator object based on user-input
2042 auto simulator = simulatorBuilder.build(useModularSimulator);
2045 if (fr->pmePpCommGpu)
2047 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
2048 fr->pmePpCommGpu.reset();
2051 if (inputrec->bPull)
2053 finish_pull(pull_work);
2055 finish_swapcoords(swap);
2059 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
2061 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(ModuleMultiThread::Pme));
2062 gmx_pmeonly(pmedata,
2066 walltime_accounting,
2069 runScheduleWork.simulationWork.useGpuPmePpCommunication,
2070 deviceStreamManager.get());
2073 wallcycle_stop(wcycle.get(), WallCycleCounter::Run);
2075 /* Finish up, write some stuff
2076 * if rerunMD, don't write last frame again
2084 walltime_accounting,
2085 fr ? fr->nbv.get() : nullptr,
2087 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
2090 deviceStreamManager.reset(nullptr);
2094 gmx_pme_destroy(pmedata);
2098 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
2099 // before we destroy the GPU context(s)
2100 // Pinned buffers are associated with contexts in CUDA.
2101 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
2102 mdAtoms.reset(nullptr);
2103 globalState.reset(nullptr);
2104 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
2105 gpuBonded.reset(nullptr);
2106 fr.reset(nullptr); // destruct forcerec before gpu
2107 // TODO convert to C++ so we can get rid of these frees
2110 if (!hwinfo_->deviceInfoList.empty())
2112 /* stop the GPU profiler (only CUDA) */
2116 /* With tMPI we need to wait for all ranks to finish deallocation before
2117 * destroying the CUDA context as some tMPI ranks may be sharing
2120 * This is not a concern in OpenCL where we use one context per rank.
2122 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
2123 * but it is easier and more futureproof to call it on the whole node.
2125 * Note that this function needs to be called even if GPUs are not used
2126 * in this run because the PME ranks have no knowledge of whether GPUs
2127 * are used or not, but all ranks need to enter the barrier below.
2128 * \todo Remove this physical node barrier after making sure
2129 * that it's not needed anymore (with a shared GPU run).
2133 physicalNodeComm.barrier();
2136 if (!devFlags.usingCudaAwareMpi)
2138 // Don't reset GPU in case of CUDA-AWARE MPI
2139 // UCX creates CUDA buffers which are cleaned-up as part of MPI_Finalize()
2140 // resetting the device before MPI_Finalize() results in crashes inside UCX
2141 releaseDevice(deviceInfo);
2144 /* Does what it says */
2145 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
2146 walltime_accounting_destroy(walltime_accounting);
2148 // Ensure log file content is written
2151 gmx_fio_flush(logFileHandle);
2154 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
2155 * exceptions were enabled before function was called. */
2158 gmx_fedisableexcept();
2161 auto rc = static_cast<int>(gmx_get_stop_condition());
2164 /* we need to join all threads. The sub-threads join when they
2165 exit this function, but the master thread needs to be told to
2175 Mdrunner::~Mdrunner()
2177 // Clean up of the Manager.
2178 // This will end up getting called on every thread-MPI rank, which is unnecessary,
2179 // but okay as long as threads synchronize some time before adding or accessing
2180 // a new set of restraints.
2181 if (restraintManager_)
2183 restraintManager_->clear();
2184 GMX_ASSERT(restraintManager_->countRestraints() == 0,
2185 "restraints added during runner life time should be cleared at runner "
2190 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
2192 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
2193 // Not sure if this should be logged through the md logger or something else,
2194 // but it is helpful to have some sort of INFO level message sent somewhere.
2195 // std::cout << "Registering restraint named " << name << std::endl;
2197 // When multiple restraints are used, it may be wasteful to register them separately.
2198 // Maybe instead register an entire Restraint Manager as a force provider.
2199 restraintManager_->addToSpec(std::move(puller), name);
2202 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
2204 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
2206 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265 in CentOS 7
2207 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
2209 class Mdrunner::BuilderImplementation
2212 BuilderImplementation() = delete;
2213 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
2214 ~BuilderImplementation();
2216 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
2217 real forceWarningThreshold,
2218 StartingBehavior startingBehavior);
2220 void addHardwareDetectionResult(const gmx_hw_info_t* hwinfo);
2222 void addDomdec(const DomdecOptions& options);
2224 void addInput(SimulationInputHandle inputHolder);
2226 void addVerletList(int nstlist);
2228 void addReplicaExchange(const ReplicaExchangeParameters& params);
2230 void addNonBonded(const char* nbpu_opt);
2232 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
2234 void addBondedTaskAssignment(const char* bonded_opt);
2236 void addUpdateTaskAssignment(const char* update_opt);
2238 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
2240 void addFilenames(ArrayRef<const t_filenm> filenames);
2242 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
2244 void addLogFile(t_fileio* logFileHandle);
2246 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
2251 // Default parameters copied from runner.h
2252 // \todo Clarify source(s) of default parameters.
2254 const char* nbpu_opt_ = nullptr;
2255 const char* pme_opt_ = nullptr;
2256 const char* pme_fft_opt_ = nullptr;
2257 const char* bonded_opt_ = nullptr;
2258 const char* update_opt_ = nullptr;
2260 MdrunOptions mdrunOptions_;
2262 DomdecOptions domdecOptions_;
2264 ReplicaExchangeParameters replicaExchangeParameters_;
2266 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
2269 //! World communicator, used for hardware detection and task assignment
2270 MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
2272 //! Multisim communicator handle.
2273 gmx_multisim_t* multiSimulation_;
2275 //! mdrun communicator
2276 MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
2278 //! Print a warning if any force is larger than this (in kJ/mol nm).
2279 real forceWarningThreshold_ = -1;
2281 //! Whether the simulation will start afresh, or restart with/without appending.
2282 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
2284 //! The modules that comprise the functionality of mdrun.
2285 std::unique_ptr<MDModules> mdModules_;
2287 //! Detected hardware.
2288 const gmx_hw_info_t* hwinfo_ = nullptr;
2290 //! \brief Parallelism information.
2291 gmx_hw_opt_t hardwareOptions_;
2293 //! filename options for simulation.
2294 ArrayRef<const t_filenm> filenames_;
2296 /*! \brief Handle to output environment.
2298 * \todo gmx_output_env_t needs lifetime management.
2300 gmx_output_env_t* outputEnvironment_ = nullptr;
2302 /*! \brief Non-owning handle to MD log file.
2304 * \todo Context should own output facilities for client.
2305 * \todo Improve log file handle management.
2307 * Code managing the FILE* relies on the ability to set it to
2308 * nullptr to check whether the filehandle is valid.
2310 t_fileio* logFileHandle_ = nullptr;
2313 * \brief Builder for simulation stop signal handler.
2315 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
2318 * \brief Sources for initial simulation state.
2320 * See issue #3652 for near-term refinements to the SimulationInput interface.
2322 * See issue #3379 for broader discussion on API aspects of simulation inputs and outputs.
2324 SimulationInputHandle inputHolder_;
2327 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
2328 compat::not_null<SimulationContext*> context) :
2329 mdModules_(std::move(mdModules))
2331 libraryWorldCommunicator_ = context->libraryWorldCommunicator_;
2332 simulationCommunicator_ = context->simulationCommunicator_;
2333 multiSimulation_ = context->multiSimulation_.get();
2336 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
2338 Mdrunner::BuilderImplementation&
2339 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
2340 const real forceWarningThreshold,
2341 const StartingBehavior startingBehavior)
2343 mdrunOptions_ = options;
2344 forceWarningThreshold_ = forceWarningThreshold;
2345 startingBehavior_ = startingBehavior;
2349 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
2351 domdecOptions_ = options;
2354 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
2359 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
2361 replicaExchangeParameters_ = params;
2364 Mdrunner Mdrunner::BuilderImplementation::build()
2366 auto newRunner = Mdrunner(std::move(mdModules_));
2368 newRunner.mdrunOptions = mdrunOptions_;
2369 newRunner.pforce = forceWarningThreshold_;
2370 newRunner.startingBehavior = startingBehavior_;
2371 newRunner.domdecOptions = domdecOptions_;
2373 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2374 newRunner.hw_opt = hardwareOptions_;
2376 // No invariant to check. This parameter exists to optionally override other behavior.
2377 newRunner.nstlist_cmdline = nstlist_;
2379 newRunner.replExParams = replicaExchangeParameters_;
2381 newRunner.filenames = filenames_;
2383 newRunner.libraryWorldCommunicator = libraryWorldCommunicator_;
2385 newRunner.simulationCommunicator = simulationCommunicator_;
2387 // nullptr is a valid value for the multisim handle
2388 newRunner.ms = multiSimulation_;
2392 newRunner.hwinfo_ = hwinfo_;
2396 GMX_THROW(gmx::APIError(
2397 "MdrunnerBuilder::addHardwareDetectionResult() is required before build()"));
2402 newRunner.inputHolder_ = std::move(inputHolder_);
2406 GMX_THROW(gmx::APIError("MdrunnerBuilder::addInput() is required before build()."));
2409 // \todo Clarify ownership and lifetime management for gmx_output_env_t
2410 // \todo Update sanity checking when output environment has clearly specified invariants.
2411 // Initialization and default values for oenv are not well specified in the current version.
2412 if (outputEnvironment_)
2414 newRunner.oenv = outputEnvironment_;
2418 GMX_THROW(gmx::APIError(
2419 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2422 newRunner.logFileHandle = logFileHandle_;
2426 newRunner.nbpu_opt = nbpu_opt_;
2430 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2433 if (pme_opt_ && pme_fft_opt_)
2435 newRunner.pme_opt = pme_opt_;
2436 newRunner.pme_fft_opt = pme_fft_opt_;
2440 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2445 newRunner.bonded_opt = bonded_opt_;
2449 GMX_THROW(gmx::APIError(
2450 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2455 newRunner.update_opt = update_opt_;
2459 GMX_THROW(gmx::APIError(
2460 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2464 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2466 if (stopHandlerBuilder_)
2468 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2472 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2478 void Mdrunner::BuilderImplementation::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2483 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2485 nbpu_opt_ = nbpu_opt;
2488 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2491 pme_fft_opt_ = pme_fft_opt;
2494 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2496 bonded_opt_ = bonded_opt;
2499 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2501 update_opt_ = update_opt;
2504 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2506 hardwareOptions_ = hardwareOptions;
2509 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2511 filenames_ = filenames;
2514 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2516 outputEnvironment_ = outputEnvironment;
2519 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2521 logFileHandle_ = logFileHandle;
2524 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2526 stopHandlerBuilder_ = std::move(builder);
2529 void Mdrunner::BuilderImplementation::addInput(SimulationInputHandle inputHolder)
2531 inputHolder_ = std::move(inputHolder);
2534 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2535 compat::not_null<SimulationContext*> context) :
2536 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2540 MdrunnerBuilder::~MdrunnerBuilder() = default;
2542 MdrunnerBuilder& MdrunnerBuilder::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2544 impl_->addHardwareDetectionResult(hwinfo);
2548 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2549 real forceWarningThreshold,
2550 const StartingBehavior startingBehavior)
2552 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2556 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2558 impl_->addDomdec(options);
2562 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2564 impl_->addVerletList(nstlist);
2568 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2570 impl_->addReplicaExchange(params);
2574 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2576 impl_->addNonBonded(nbpu_opt);
2580 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2582 // The builder method may become more general in the future, but in this version,
2583 // parameters for PME electrostatics are both required and the only parameters
2585 if (pme_opt && pme_fft_opt)
2587 impl_->addPME(pme_opt, pme_fft_opt);
2592 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2597 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2599 impl_->addBondedTaskAssignment(bonded_opt);
2603 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2605 impl_->addUpdateTaskAssignment(update_opt);
2609 Mdrunner MdrunnerBuilder::build()
2611 return impl_->build();
2614 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2616 impl_->addHardwareOptions(hardwareOptions);
2620 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2622 impl_->addFilenames(filenames);
2626 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2628 impl_->addOutputEnvironment(outputEnvironment);
2632 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2634 impl_->addLogFile(logFileHandle);
2638 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2640 impl_->addStopHandlerBuilder(std::move(builder));
2644 MdrunnerBuilder& MdrunnerBuilder::addInput(SimulationInputHandle input)
2646 impl_->addInput(std::move(input));
2650 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2652 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;