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_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_stream_manager.h"
77 #include "gromacs/hardware/cpuinfo.h"
78 #include "gromacs/hardware/detecthardware.h"
79 #include "gromacs/hardware/device_management.h"
80 #include "gromacs/hardware/hardwaretopology.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/gpuforcereduction.h"
99 #include "gromacs/mdlib/makeconstraints.h"
100 #include "gromacs/mdlib/md_support.h"
101 #include "gromacs/mdlib/mdatoms.h"
102 #include "gromacs/mdlib/sighandler.h"
103 #include "gromacs/mdlib/stophandler.h"
104 #include "gromacs/mdlib/tgroup.h"
105 #include "gromacs/mdlib/updategroups.h"
106 #include "gromacs/mdlib/vsite.h"
107 #include "gromacs/mdrun/mdmodules.h"
108 #include "gromacs/mdrun/simulationcontext.h"
109 #include "gromacs/mdrun/simulationinput.h"
110 #include "gromacs/mdrun/simulationinputhandle.h"
111 #include "gromacs/mdrunutility/handlerestart.h"
112 #include "gromacs/mdrunutility/logging.h"
113 #include "gromacs/mdrunutility/multisim.h"
114 #include "gromacs/mdrunutility/printtime.h"
115 #include "gromacs/mdrunutility/threadaffinity.h"
116 #include "gromacs/mdtypes/checkpointdata.h"
117 #include "gromacs/mdtypes/commrec.h"
118 #include "gromacs/mdtypes/enerdata.h"
119 #include "gromacs/mdtypes/fcdata.h"
120 #include "gromacs/mdtypes/forcerec.h"
121 #include "gromacs/mdtypes/group.h"
122 #include "gromacs/mdtypes/inputrec.h"
123 #include "gromacs/mdtypes/interaction_const.h"
124 #include "gromacs/mdtypes/md_enums.h"
125 #include "gromacs/mdtypes/mdatom.h"
126 #include "gromacs/mdtypes/mdrunoptions.h"
127 #include "gromacs/mdtypes/observableshistory.h"
128 #include "gromacs/mdtypes/simulation_workload.h"
129 #include "gromacs/mdtypes/state.h"
130 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
131 #include "gromacs/modularsimulator/modularsimulator.h"
132 #include "gromacs/nbnxm/gpu_data_mgmt.h"
133 #include "gromacs/nbnxm/nbnxm.h"
134 #include "gromacs/nbnxm/pairlist_tuning.h"
135 #include "gromacs/pbcutil/pbc.h"
136 #include "gromacs/pulling/output.h"
137 #include "gromacs/pulling/pull.h"
138 #include "gromacs/pulling/pull_rotation.h"
139 #include "gromacs/restraint/manager.h"
140 #include "gromacs/restraint/restraintmdmodule.h"
141 #include "gromacs/restraint/restraintpotential.h"
142 #include "gromacs/swap/swapcoords.h"
143 #include "gromacs/taskassignment/decidegpuusage.h"
144 #include "gromacs/taskassignment/decidesimulationworkload.h"
145 #include "gromacs/taskassignment/resourcedivision.h"
146 #include "gromacs/taskassignment/taskassignment.h"
147 #include "gromacs/taskassignment/usergpuids.h"
148 #include "gromacs/timing/gpu_timing.h"
149 #include "gromacs/timing/wallcycle.h"
150 #include "gromacs/timing/wallcyclereporting.h"
151 #include "gromacs/topology/mtop_util.h"
152 #include "gromacs/trajectory/trajectoryframe.h"
153 #include "gromacs/utility/basenetwork.h"
154 #include "gromacs/utility/cstringutil.h"
155 #include "gromacs/utility/exceptions.h"
156 #include "gromacs/utility/fatalerror.h"
157 #include "gromacs/utility/filestream.h"
158 #include "gromacs/utility/gmxassert.h"
159 #include "gromacs/utility/gmxmpi.h"
160 #include "gromacs/utility/keyvaluetree.h"
161 #include "gromacs/utility/logger.h"
162 #include "gromacs/utility/loggerbuilder.h"
163 #include "gromacs/utility/mdmodulenotification.h"
164 #include "gromacs/utility/physicalnodecommunicator.h"
165 #include "gromacs/utility/pleasecite.h"
166 #include "gromacs/utility/programcontext.h"
167 #include "gromacs/utility/smalloc.h"
168 #include "gromacs/utility/stringutil.h"
170 #include "isimulator.h"
171 #include "membedholder.h"
172 #include "replicaexchange.h"
173 #include "simulatorbuilder.h"
179 /*! \brief Manage any development feature flag variables encountered
181 * The use of dev features indicated by environment variables is
182 * logged in order to ensure that runs with such features enabled can
183 * be identified from their log and standard output. Any cross
184 * dependencies are also checked, and if unsatisfied, a fatal error
187 * Note that some development features overrides are applied already here:
188 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
190 * \param[in] mdlog Logger object.
191 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
192 * \param[in] pmeRunMode The PME run mode for this run
193 * \returns The object populated with development feature flags.
195 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
196 const bool useGpuForNonbonded,
197 const PmeRunMode pmeRunMode)
199 DevelopmentFeatureFlags devFlags;
201 // Some builds of GCC 5 give false positive warnings that these
202 // getenv results are ignored when clearly they are used.
203 #pragma GCC diagnostic push
204 #pragma GCC diagnostic ignored "-Wunused-result"
206 devFlags.enableGpuBufferOps =
207 GMX_GPU_CUDA && useGpuForNonbonded && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
208 devFlags.enableGpuHaloExchange = GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_DD_COMMS") != nullptr;
209 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr) || GMX_FAHCORE;
210 devFlags.enableGpuPmePPComm =
211 GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
213 #pragma GCC diagnostic pop
215 if (devFlags.enableGpuBufferOps)
217 GMX_LOG(mdlog.warning)
219 .appendTextFormatted(
220 "This run uses the 'GPU buffer ops' feature, enabled by the "
221 "GMX_USE_GPU_BUFFER_OPS environment variable.");
224 if (devFlags.forceGpuUpdateDefault)
226 GMX_LOG(mdlog.warning)
228 .appendTextFormatted(
229 "This run will default to '-update gpu' as requested by the "
230 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
231 "decomposition lacks substantial testing and should be used with caution.");
234 if (devFlags.enableGpuHaloExchange)
236 if (useGpuForNonbonded)
238 if (!devFlags.enableGpuBufferOps)
240 GMX_LOG(mdlog.warning)
242 .appendTextFormatted(
243 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
244 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
245 devFlags.enableGpuBufferOps = true;
247 GMX_LOG(mdlog.warning)
249 .appendTextFormatted(
250 "This run has requested the 'GPU halo exchange' feature, enabled by "
252 "GMX_GPU_DD_COMMS environment variable.");
256 GMX_LOG(mdlog.warning)
258 .appendTextFormatted(
259 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
260 "halo exchange' feature will not be enabled as nonbonded interactions "
261 "are not offloaded.");
262 devFlags.enableGpuHaloExchange = false;
266 if (devFlags.enableGpuPmePPComm)
268 if (pmeRunMode == PmeRunMode::GPU)
270 if (!devFlags.enableGpuBufferOps)
272 GMX_LOG(mdlog.warning)
274 .appendTextFormatted(
275 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
276 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
277 devFlags.enableGpuBufferOps = true;
279 GMX_LOG(mdlog.warning)
281 .appendTextFormatted(
282 "This run uses the 'GPU PME-PP communications' feature, enabled "
283 "by the GMX_GPU_PME_PP_COMMS environment variable.");
287 std::string clarification;
288 if (pmeRunMode == PmeRunMode::Mixed)
291 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
296 clarification = "PME is not offloaded to the GPU.";
298 GMX_LOG(mdlog.warning)
301 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
302 "'GPU PME-PP communications' feature was not enabled as "
304 devFlags.enableGpuPmePPComm = false;
311 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
313 * Used to ensure that the master thread does not modify mdrunner during copy
314 * on the spawned threads. */
315 static void threadMpiMdrunnerAccessBarrier()
318 MPI_Barrier(MPI_COMM_WORLD);
322 Mdrunner Mdrunner::cloneOnSpawnedThread() const
324 auto newRunner = Mdrunner(std::make_unique<MDModules>());
326 // All runners in the same process share a restraint manager resource because it is
327 // part of the interface to the client code, which is associated only with the
328 // original thread. Handles to the same resources can be obtained by copy.
330 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
333 // Copy members of master runner.
334 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
335 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
336 newRunner.hw_opt = hw_opt;
337 newRunner.filenames = filenames;
339 newRunner.hwinfo_ = hwinfo_;
340 newRunner.oenv = oenv;
341 newRunner.mdrunOptions = mdrunOptions;
342 newRunner.domdecOptions = domdecOptions;
343 newRunner.nbpu_opt = nbpu_opt;
344 newRunner.pme_opt = pme_opt;
345 newRunner.pme_fft_opt = pme_fft_opt;
346 newRunner.bonded_opt = bonded_opt;
347 newRunner.update_opt = update_opt;
348 newRunner.nstlist_cmdline = nstlist_cmdline;
349 newRunner.replExParams = replExParams;
350 newRunner.pforce = pforce;
351 // Give the spawned thread the newly created valid communicator
352 // for the simulation.
353 newRunner.libraryWorldCommunicator = MPI_COMM_WORLD;
354 newRunner.simulationCommunicator = MPI_COMM_WORLD;
356 newRunner.startingBehavior = startingBehavior;
357 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
358 newRunner.inputHolder_ = inputHolder_;
360 threadMpiMdrunnerAccessBarrier();
365 /*! \brief The callback used for running on spawned threads.
367 * Obtains the pointer to the master mdrunner object from the one
368 * argument permitted to the thread-launch API call, copies it to make
369 * a new runner for this thread, reinitializes necessary data, and
370 * proceeds to the simulation. */
371 static void mdrunner_start_fn(const void* arg)
375 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
376 /* copy the arg list to make sure that it's thread-local. This
377 doesn't copy pointed-to items, of course; fnm, cr and fplog
378 are reset in the call below, all others should be const. */
379 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
382 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
386 void Mdrunner::spawnThreads(int numThreadsToLaunch)
389 /* now spawn new threads that start mdrunner_start_fn(), while
390 the main thread returns. Thread affinity is handled later. */
391 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
392 static_cast<const void*>(this))
395 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
398 // Give the master thread the newly created valid communicator for
400 libraryWorldCommunicator = MPI_COMM_WORLD;
401 simulationCommunicator = MPI_COMM_WORLD;
402 threadMpiMdrunnerAccessBarrier();
404 GMX_UNUSED_VALUE(numThreadsToLaunch);
405 GMX_UNUSED_VALUE(mdrunner_start_fn);
411 /*! \brief Initialize variables for Verlet scheme simulation */
412 static void prepare_verlet_scheme(FILE* fplog,
416 const gmx_mtop_t* mtop,
418 bool makeGpuPairList,
419 const gmx::CpuInfo& cpuinfo)
421 // We checked the cut-offs in grompp, but double-check here.
422 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
423 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
425 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
426 "With Verlet lists and PME we should have rcoulomb>=rvdw");
430 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
431 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
433 /* For NVE simulations, we will retain the initial list buffer */
434 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
436 /* Update the Verlet buffer size for the current run setup */
438 /* Here we assume SIMD-enabled kernels are being used. But as currently
439 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
440 * and 4x2 gives a larger buffer than 4x4, this is ok.
442 ListSetupType listType =
443 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
444 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
446 const real rlist_new =
447 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
449 if (rlist_new != ir->rlist)
451 if (fplog != nullptr)
454 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
455 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
457 ir->rlist = rlist_new;
461 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
463 gmx_fatal(FARGS, "Can not set nstlist without %s",
464 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
467 if (EI_DYNAMICS(ir->eI))
469 /* Set or try nstlist values */
470 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
474 /*! \brief Override the nslist value in inputrec
476 * with value passed on the command line (if any)
478 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
482 /* override with anything else than the default -2 */
483 if (nsteps_cmdline > -2)
485 char sbuf_steps[STEPSTRSIZE];
486 char sbuf_msg[STRLEN];
488 ir->nsteps = nsteps_cmdline;
489 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
492 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
493 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
497 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
498 gmx_step_str(nsteps_cmdline, sbuf_steps));
501 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
503 else if (nsteps_cmdline < -2)
505 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
507 /* Do nothing if nsteps_cmdline == -2 */
513 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
515 * If not, and if a warning may be issued, logs a warning about
516 * falling back to CPU code. With thread-MPI, only the first
517 * call to this function should have \c issueWarning true. */
518 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
520 bool gpuIsUseful = true;
523 if (ir.opts.ngener - ir.nwall > 1)
525 /* The GPU code does not support more than one energy group.
526 * If the user requested GPUs explicitly, a fatal error is given later.
530 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
531 "For better performance, run on the GPU without energy groups and then do "
532 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
538 warning = "TPI is not implemented for GPUs.";
541 if (!gpuIsUseful && issueWarning)
543 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
549 //! Initializes the logger for mdrun.
550 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
552 gmx::LoggerBuilder builder;
553 if (fplog != nullptr)
555 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
557 if (isSimulationMasterRank)
559 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
561 return builder.build();
564 //! Make a TaskTarget from an mdrun argument string.
565 static TaskTarget findTaskTarget(const char* optionString)
567 TaskTarget returnValue = TaskTarget::Auto;
569 if (strncmp(optionString, "auto", 3) == 0)
571 returnValue = TaskTarget::Auto;
573 else if (strncmp(optionString, "cpu", 3) == 0)
575 returnValue = TaskTarget::Cpu;
577 else if (strncmp(optionString, "gpu", 3) == 0)
579 returnValue = TaskTarget::Gpu;
583 GMX_ASSERT(false, "Option string should have been checked for sanity already");
589 //! Finish run, aggregate data to print performance info.
590 static void finish_run(FILE* fplog,
591 const gmx::MDLogger& mdlog,
593 const t_inputrec* inputrec,
595 gmx_wallcycle_t wcycle,
596 gmx_walltime_accounting_t walltime_accounting,
597 nonbonded_verlet_t* nbv,
598 const gmx_pme_t* pme,
602 double nbfs = 0, mflop = 0;
603 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
604 elapsed_time_over_all_threads_over_all_ranks;
605 /* Control whether it is valid to print a report. Only the
606 simulation master may print, but it should not do so if the run
607 terminated e.g. before a scheduled reset step. This is
608 complicated by the fact that PME ranks are unaware of the
609 reason why they were sent a pmerecvqxFINISH. To avoid
610 communication deadlocks, we always do the communication for the
611 report, even if we've decided not to write the report, because
612 how long it takes to finish the run is not important when we've
613 decided not to report on the simulation performance.
615 Further, we only report performance for dynamical integrators,
616 because those are the only ones for which we plan to
617 consider doing any optimizations. */
618 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
620 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
622 GMX_LOG(mdlog.warning)
624 .appendText("Simulation ended prematurely, no performance report will be written.");
629 std::unique_ptr<t_nrnb> nrnbTotalStorage;
632 nrnbTotalStorage = std::make_unique<t_nrnb>();
633 nrnb_tot = nrnbTotalStorage.get();
635 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
643 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
644 elapsed_time_over_all_threads =
645 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
649 /* reduce elapsed_time over all MPI ranks in the current simulation */
650 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
652 elapsed_time_over_all_ranks /= cr->nnodes;
653 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
654 * current simulation. */
655 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
656 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
661 elapsed_time_over_all_ranks = elapsed_time;
662 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
667 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
670 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
672 print_dd_statistics(cr, inputrec, fplog);
675 /* TODO Move the responsibility for any scaling by thread counts
676 * to the code that handled the thread region, so that there's a
677 * mechanism to keep cycle counting working during the transition
678 * to task parallelism. */
679 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
680 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
681 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
682 nthreads_pp, nthreads_pme);
683 auto cycle_sum(wallcycle_sum(cr, wcycle));
687 auto nbnxn_gpu_timings =
688 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
689 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
691 if (pme_gpu_task_enabled(pme))
693 pme_gpu_get_timings(pme, &pme_gpu_timings);
695 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
696 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
699 if (EI_DYNAMICS(inputrec->eI))
701 delta_t = inputrec->delta_t;
706 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
707 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
708 delta_t, nbfs, mflop);
712 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
713 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
714 delta_t, nbfs, mflop);
719 int Mdrunner::mdrunner()
722 t_forcerec* fr = nullptr;
723 real ewaldcoeff_q = 0;
724 real ewaldcoeff_lj = 0;
725 int nChargePerturbed = -1, nTypePerturbed = 0;
726 gmx_wallcycle_t wcycle;
727 gmx_walltime_accounting_t walltime_accounting = nullptr;
728 MembedHolder membedHolder(filenames.size(), filenames.data());
730 /* CAUTION: threads may be started later on in this function, so
731 cr doesn't reflect the final parallel state right now */
734 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
735 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
736 const bool doRerun = mdrunOptions.rerun;
738 // Handle task-assignment related user options.
739 EmulateGpuNonbonded emulateGpuNonbonded =
740 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
742 std::vector<int> userGpuTaskAssignment;
745 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
747 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
748 auto nonbondedTarget = findTaskTarget(nbpu_opt);
749 auto pmeTarget = findTaskTarget(pme_opt);
750 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
751 auto bondedTarget = findTaskTarget(bonded_opt);
752 auto updateTarget = findTaskTarget(update_opt);
754 FILE* fplog = nullptr;
755 // If we are appending, we don't write log output because we need
756 // to check that the old log file matches what the checkpoint file
757 // expects. Otherwise, we should start to write log output now if
758 // there is a file ready for it.
759 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
761 fplog = gmx_fio_getfp(logFileHandle);
763 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, simulationCommunicator);
764 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
765 gmx::MDLogger mdlog(logOwner.logger());
767 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo_);
769 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo_->deviceInfoList, hw_opt.gpuIdsAvailable);
770 const int numDevicesToUse = gmx::ssize(gpuIdsToUse);
772 // Print citation requests after all software/hardware printing
773 pleaseCiteGromacs(fplog);
775 // Note: legacy program logic relies on checking whether these pointers are assigned.
776 // Objects may or may not be allocated later.
777 std::unique_ptr<t_inputrec> inputrec;
778 std::unique_ptr<t_state> globalState;
780 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
782 if (isSimulationMasterRank)
784 // Allocate objects to be initialized by later function calls.
785 /* Only the master rank has the global state */
786 globalState = std::make_unique<t_state>();
787 inputrec = std::make_unique<t_inputrec>();
789 /* Read (nearly) all data required for the simulation
790 * and keep the partly serialized tpr contents to send to other ranks later
792 applyGlobalSimulationState(*inputHolder_.get(), partialDeserializedTpr.get(),
793 globalState.get(), inputrec.get(), &mtop);
796 /* Check and update the hardware options for internal consistency */
797 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
800 if (GMX_THREAD_MPI && isSimulationMasterRank)
802 bool useGpuForNonbonded = false;
803 bool useGpuForPme = false;
806 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
808 // If the user specified the number of ranks, then we must
809 // respect that, but in default mode, we need to allow for
810 // the number of GPUs to choose the number of ranks.
811 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
812 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
813 nonbondedTarget, numDevicesToUse, userGpuTaskAssignment, emulateGpuNonbonded,
814 canUseGpuForNonbonded,
815 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
816 hw_opt.nthreads_tmpi);
817 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
818 useGpuForNonbonded, pmeTarget, pmeFftTarget, numDevicesToUse, userGpuTaskAssignment,
819 *hwinfo_, *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
821 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
823 /* Determine how many thread-MPI ranks to start.
825 * TODO Over-writing the user-supplied value here does
826 * prevent any possible subsequent checks from working
828 hw_opt.nthreads_tmpi =
829 get_nthreads_mpi(hwinfo_, &hw_opt, numDevicesToUse, useGpuForNonbonded, useGpuForPme,
830 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
832 // Now start the threads for thread MPI.
833 spawnThreads(hw_opt.nthreads_tmpi);
834 // The spawned threads enter mdrunner() and execution of
835 // master and spawned threads joins at the end of this block.
838 GMX_RELEASE_ASSERT(ms || simulationCommunicator != MPI_COMM_NULL,
839 "Must have valid communicator unless running a multi-simulation");
840 CommrecHandle crHandle = init_commrec(simulationCommunicator);
841 t_commrec* cr = crHandle.get();
842 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
844 PhysicalNodeCommunicator physicalNodeComm(libraryWorldCommunicator, gmx_physicalnode_id_hash());
846 // If we detected the topology on this system, double-check that it makes sense
847 if (hwinfo_->hardwareTopology->isThisSystem())
849 hardwareTopologyDoubleCheckDetection(mdlog, *hwinfo_->hardwareTopology);
854 /* now broadcast everything to the non-master nodes/threads: */
855 if (!isSimulationMasterRank)
857 // Until now, only the master rank has a non-null pointer.
858 // On non-master ranks, allocate the object that will receive data in the following call.
859 inputrec = std::make_unique<t_inputrec>();
861 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec.get(), &mtop,
862 partialDeserializedTpr.get());
864 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
865 partialDeserializedTpr.reset(nullptr);
867 // Now the number of ranks is known to all ranks, and each knows
868 // the inputrec read by the master rank. The ranks can now all run
869 // the task-deciding functions and will agree on the result
870 // without needing to communicate.
871 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
873 // Note that these variables describe only their own node.
875 // Note that when bonded interactions run on a GPU they always run
876 // alongside a nonbonded task, so do not influence task assignment
877 // even though they affect the force calculation workload.
878 bool useGpuForNonbonded = false;
879 bool useGpuForPme = false;
880 bool useGpuForBonded = false;
881 bool useGpuForUpdate = false;
882 bool gpusWereDetected = hwinfo_->ngpu_compatible_tot > 0;
885 // It's possible that there are different numbers of GPUs on
886 // different nodes, which is the user's responsibility to
887 // handle. If unsuitable, we will notice that during task
889 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
890 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
891 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
892 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
893 useGpuForPme = decideWhetherToUseGpusForPme(
894 useGpuForNonbonded, pmeTarget, pmeFftTarget, userGpuTaskAssignment, *hwinfo_,
895 *inputrec, cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
896 useGpuForBonded = decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme,
897 bondedTarget, *inputrec, mtop,
898 domdecOptions.numPmeRanks, gpusWereDetected);
900 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
902 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
904 // Initialize development feature flags that enabled by environment variable
905 // and report those features that are enabled.
906 const DevelopmentFeatureFlags devFlags =
907 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
909 const bool useModularSimulator =
910 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
911 nullptr, doEssentialDynamics, membedHolder.doMembed());
914 // TODO: hide restraint implementation details from Mdrunner.
915 // There is nothing unique about restraints at this point as far as the
916 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
917 // factory functions from the SimulationContext on which to call mdModules_->add().
918 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
919 for (auto&& restraint : restraintManager_->getRestraints())
921 auto module = RestraintMDModule::create(restraint, restraint->sites());
922 mdModules_->add(std::move(module));
925 // TODO: Error handling
926 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
927 // now that the MdModules know their options, they know which callbacks to sign up to
928 mdModules_->subscribeToSimulationSetupNotifications();
929 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
931 if (inputrec->internalParameters != nullptr)
933 mdModulesNotifier.notify(*inputrec->internalParameters);
936 if (fplog != nullptr)
938 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
939 fprintf(fplog, "\n");
944 /* In rerun, set velocities to zero if present */
945 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
947 // rerun does not use velocities
951 "Rerun trajectory contains velocities. Rerun does only evaluate "
952 "potential energy and forces. The velocities will be ignored.");
953 for (int i = 0; i < globalState->natoms; i++)
955 clear_rvec(globalState->v[i]);
957 globalState->flags &= ~(1 << estV);
960 /* now make sure the state is initialized and propagated */
961 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
964 /* NM and TPI parallelize over force/energy calculations, not atoms,
965 * so we need to initialize and broadcast the global state.
967 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
971 globalState = std::make_unique<t_state>();
973 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
977 /* A parallel command line option consistency check that we can
978 only do after any threads have started. */
980 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
981 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
984 "The -dd or -npme option request a parallel simulation, "
986 "but %s was compiled without threads or MPI enabled",
987 output_env_get_program_display_name(oenv));
989 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
991 "but %s was not started through mpirun/mpiexec or only one rank was requested "
992 "through mpirun/mpiexec",
993 output_env_get_program_display_name(oenv));
997 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1000 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1001 "these are not compatible with mdrun -rerun");
1004 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1006 if (domdecOptions.numPmeRanks > 0)
1008 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1009 "PME-only ranks are requested, but the system does not use PME "
1010 "for electrostatics or LJ");
1013 domdecOptions.numPmeRanks = 0;
1016 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1018 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1019 * improve performance with many threads per GPU, since our OpenMP
1020 * scaling is bad, but it's difficult to automate the setup.
1022 domdecOptions.numPmeRanks = 0;
1026 if (domdecOptions.numPmeRanks < 0)
1028 domdecOptions.numPmeRanks = 0;
1029 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1033 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1034 "PME GPU decomposition is not supported");
1038 /* NMR restraints must be initialized before load_checkpoint,
1039 * since with time averaging the history is added to t_state.
1040 * For proper consistency check we therefore need to extend
1042 * So the PME-only nodes (if present) will also initialize
1043 * the distance restraints.
1046 /* This needs to be called before read_checkpoint to extend the state */
1047 t_disresdata* disresdata;
1048 snew(disresdata, 1);
1049 init_disres(fplog, &mtop, inputrec.get(), DisResRunMode::MDRun,
1050 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1051 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1052 globalState.get(), replExParams.exchangeInterval > 0);
1054 t_oriresdata* oriresdata;
1055 snew(oriresdata, 1);
1056 init_orires(fplog, &mtop, inputrec.get(), cr, ms, globalState.get(), oriresdata);
1058 auto deform = prepareBoxDeformation(
1059 globalState != nullptr ? globalState->box : box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1060 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mygroup, *inputrec);
1063 /* We have to remember the generation's first step before reading checkpoint.
1064 This way, we can report to the F@H core both the generation's first step
1065 and the restored first step, thus making it able to distinguish between
1066 an interruption/resume and start of the n-th generation simulation.
1067 Having this information, the F@H core can correctly calculate and report
1070 int gen_first_step = 0;
1073 gen_first_step = inputrec->init_step;
1077 ObservablesHistory observablesHistory = {};
1079 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1080 if (startingBehavior != StartingBehavior::NewSimulation)
1082 /* Check if checkpoint file exists before doing continuation.
1083 * This way we can use identical input options for the first and subsequent runs...
1085 if (mdrunOptions.numStepsCommandline > -2)
1087 /* Temporarily set the number of steps to unlimited to avoid
1088 * triggering the nsteps check in load_checkpoint().
1089 * This hack will go away soon when the -nsteps option is removed.
1091 inputrec->nsteps = -1;
1094 // Finish applying initial simulation state information from external sources on all ranks.
1095 // Reconcile checkpoint file data with Mdrunner state established up to this point.
1096 applyLocalState(*inputHolder_.get(), logFileHandle, cr, domdecOptions.numCells,
1097 inputrec.get(), globalState.get(), &observablesHistory,
1098 mdrunOptions.reproducible, mdModules_->notifier(),
1099 modularSimulatorCheckpointData.get(), useModularSimulator);
1100 // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1102 // on all code paths.
1103 // Write checkpoint or provide hook to update SimulationInput.
1104 // If there was a checkpoint file, SimulationInput contains more information
1105 // than if there wasn't. At this point, we have synchronized the in-memory
1106 // state with the filesystem state only for restarted simulations. We should
1107 // be calling applyLocalState unconditionally and expect that the completeness
1108 // of SimulationInput is not dependent on its creation method.
1110 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1112 // Now we can start normal logging to the truncated log file.
1113 fplog = gmx_fio_getfp(logFileHandle);
1114 prepareLogAppending(fplog);
1115 logOwner = buildLogger(fplog, MASTER(cr));
1116 mdlog = logOwner.logger();
1123 fcRegisterSteps(inputrec->nsteps + inputrec->init_step, gen_first_step);
1127 if (mdrunOptions.numStepsCommandline > -2)
1132 "The -nsteps functionality is deprecated, and may be removed in a future "
1134 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1137 /* override nsteps with value set on the commandline */
1138 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1140 if (isSimulationMasterRank)
1142 copy_mat(globalState->box, box);
1147 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1150 if (inputrec->cutoff_scheme != ecutsVERLET)
1153 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1154 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1156 /* Update rlist and nstlist. */
1157 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1158 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1159 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1161 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1162 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1165 // This builder is necessary while we have multi-part construction
1166 // of DD. Before DD is constructed, we use the existence of
1167 // the builder object to indicate that further construction of DD
1169 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1170 if (useDomainDecomposition)
1172 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1173 mdlog, cr, domdecOptions, mdrunOptions, mtop, *inputrec, box,
1174 positionsFromStatePointer(globalState.get()));
1178 /* PME, if used, is done on all nodes with 1D decomposition */
1179 cr->nnodes = cr->sizeOfDefaultCommunicator;
1180 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1181 cr->nodeid = cr->rankInDefaultCommunicator;
1183 cr->duty = (DUTY_PP | DUTY_PME);
1185 if (inputrec->pbcType == PbcType::Screw)
1187 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1191 // Produce the task assignment for this rank - done after DD is constructed
1192 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1193 gpuIdsToUse, userGpuTaskAssignment, *hwinfo_, simulationCommunicator, physicalNodeComm,
1194 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1195 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1196 // TODO cr->duty & DUTY_PME should imply that a PME
1197 // algorithm is active, but currently does not.
1198 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1200 // Get the device handles for the modules, nullptr when no task is assigned.
1202 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1204 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1205 bool useTiming = true;
1209 /* WARNING: CUDA timings are incorrect with multiple streams.
1210 * This is the main reason why they are disabled by default.
1212 // TODO: Consider turning on by default when we can detect nr of streams.
1213 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1215 else if (GMX_GPU_OPENCL)
1217 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1220 // TODO Currently this is always built, yet DD partition code
1221 // checks if it is built before using it. Probably it should
1222 // become an MDModule that is made only when another module
1223 // requires it (e.g. pull, CompEl, density fitting), so that we
1224 // don't update the local atom sets unilaterally every step.
1225 LocalAtomSetManager atomSets;
1228 // TODO Pass the GPU streams to ddBuilder to use in buffer
1229 // transfers (e.g. halo exchange)
1230 cr->dd = ddBuilder->build(&atomSets);
1231 // The builder's job is done, so destruct it
1232 ddBuilder.reset(nullptr);
1233 // Note that local state still does not exist yet.
1235 // Ensure that all atoms within the same update group are in the
1236 // same periodic image. Otherwise, a simulation that did not use
1237 // update groups (e.g. a single-rank simulation) cannot always be
1238 // correctly restarted in a way that does use update groups
1239 // (e.g. a multi-rank simulation).
1240 if (isSimulationMasterRank)
1242 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1243 if (useUpdateGroups)
1245 putUpdateGroupAtomsInSamePeriodicImage(*cr->dd, mtop, globalState->box, globalState->x);
1249 // The GPU update is decided here because we need to know whether the constraints or
1250 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1251 // defined). This is only known after DD is initialized, hence decision on using GPU
1252 // update is done so late.
1255 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1256 const bool haveFrozenAtoms = inputrecFrozenAtoms(inputrec.get());
1258 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1259 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1260 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1261 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1262 replExParams.exchangeInterval > 0, haveFrozenAtoms, doRerun, devFlags, mdlog);
1264 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1266 const bool printHostName = (cr->nnodes > 1);
1267 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1269 const bool disableNonbondedCalculation = (getenv("GMX_NO_NONBONDED") != nullptr);
1270 if (disableNonbondedCalculation)
1272 /* turn off non-bonded calculations */
1273 GMX_LOG(mdlog.warning)
1276 "Found environment variable GMX_NO_NONBONDED.\n"
1277 "Disabling nonbonded calculations.");
1280 MdrunScheduleWorkload runScheduleWork;
1282 bool useGpuDirectHalo = decideWhetherToUseGpuForHalo(
1283 devFlags, havePPDomainDecomposition(cr), useGpuForNonbonded, useModularSimulator,
1284 doRerun, EI_ENERGY_MINIMIZATION(inputrec->eI));
1286 // Also populates the simulation constant workload description.
1287 runScheduleWork.simulationWork = createSimulationWorkload(
1288 *inputrec, disableNonbondedCalculation, devFlags, useGpuForNonbonded, pmeRunMode,
1289 useGpuForBonded, useGpuForUpdate, useGpuDirectHalo);
1291 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1293 if (deviceInfo != nullptr)
1295 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1297 dd_setup_dlb_resource_sharing(cr, deviceId);
1299 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1300 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, useTiming);
1303 // If the user chose a task assignment, give them some hints
1304 // where appropriate.
1305 if (!userGpuTaskAssignment.empty())
1307 gpuTaskAssignments.logPerformanceHints(mdlog, numDevicesToUse);
1312 /* After possible communicator splitting in make_dd_communicators.
1313 * we can set up the intra/inter node communication.
1315 gmx_setup_nodecomm(fplog, cr);
1321 GMX_LOG(mdlog.warning)
1323 .appendTextFormatted(
1324 "This is simulation %d out of %d running as a composite GROMACS\n"
1325 "multi-simulation job. Setup for this simulation:\n",
1326 ms->simulationIndex_, ms->numSimulations_);
1328 GMX_LOG(mdlog.warning)
1329 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1331 cr->nnodes == 1 ? "thread" : "threads"
1333 cr->nnodes == 1 ? "process" : "processes"
1339 // If mdrun -pin auto honors any affinity setting that already
1340 // exists. If so, it is nice to provide feedback about whether
1341 // that existing affinity setting was from OpenMP or something
1342 // else, so we run this code both before and after we initialize
1343 // the OpenMP support.
1344 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, FALSE);
1345 /* Check and update the number of OpenMP threads requested */
1346 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo_, cr, ms, physicalNodeComm.size_,
1347 pmeRunMode, mtop, *inputrec);
1349 gmx_omp_nthreads_init(mdlog, cr, hwinfo_->nthreads_hw_avail, physicalNodeComm.size_,
1350 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1352 // Enable FP exception detection, but not in
1353 // Release mode and not for compilers with known buggy FP
1354 // exception support (clang with any optimization) or suspected
1355 // buggy FP exception support (gcc 7.* with optimization).
1356 #if !defined NDEBUG \
1357 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1358 && defined __OPTIMIZE__)
1359 const bool bEnableFPE = true;
1361 const bool bEnableFPE = false;
1363 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1366 gmx_feenableexcept();
1369 /* Now that we know the setup is consistent, check for efficiency */
1370 check_resource_division_efficiency(hwinfo_, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1371 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1373 /* getting number of PP/PME threads on this MPI / tMPI rank.
1374 PME: env variable should be read only on one node to make sure it is
1375 identical everywhere;
1377 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1378 : gmx_omp_nthreads_get(emntPME);
1379 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo_->hardwareTopology,
1380 physicalNodeComm, mdlog);
1382 // Enable Peer access between GPUs where available
1383 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1384 // any of the GPU communication features are active.
1385 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1386 && (runScheduleWork.simulationWork.useGpuHaloExchange
1387 || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1389 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1392 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1394 /* Before setting affinity, check whether the affinity has changed
1395 * - which indicates that probably the OpenMP library has changed it
1396 * since we first checked).
1398 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, TRUE);
1400 int numThreadsOnThisNode, intraNodeThreadOffset;
1401 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1402 &intraNodeThreadOffset);
1404 /* Set the CPU affinity */
1405 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo_->hardwareTopology, numThreadsOnThisRank,
1406 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1409 if (mdrunOptions.timingOptions.resetStep > -1)
1414 "The -resetstep functionality is deprecated, and may be removed in a "
1417 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1421 /* Master synchronizes its value of reset_counters with all nodes
1422 * including PME only nodes */
1423 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1424 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1425 wcycle_set_reset_counters(wcycle, reset_counters);
1428 // Membrane embedding must be initialized before we call init_forcerec()
1429 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1430 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1432 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1433 std::unique_ptr<MDAtoms> mdAtoms;
1434 std::unique_ptr<VirtualSitesHandler> vsite;
1435 std::unique_ptr<GpuBonded> gpuBonded;
1438 if (thisRankHasDuty(cr, DUTY_PP))
1440 mdModulesNotifier.notify(*cr);
1441 mdModulesNotifier.notify(&atomSets);
1442 mdModulesNotifier.notify(inputrec->pbcType);
1443 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1444 /* Initiate forcerecord */
1445 fr = new t_forcerec;
1446 fr->forceProviders = mdModules_->initForceProviders();
1447 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1448 opt2fn("-table", filenames.size(), filenames.data()),
1449 opt2fn("-tablep", filenames.size(), filenames.data()),
1450 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1451 // Dirty hack, for fixing disres and orires should be made mdmodules
1452 fr->fcdata->disres = disresdata;
1453 fr->fcdata->orires = oriresdata;
1455 // Save a handle to device stream manager to use elsewhere in the code
1456 // TODO: Forcerec is not a correct place to store it.
1457 fr->deviceStreamManager = deviceStreamManager.get();
1459 if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1462 deviceStreamManager != nullptr,
1463 "GPU device stream manager should be valid in order to use PME-PP direct "
1466 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1467 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1469 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1470 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1471 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1474 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo_,
1475 runScheduleWork.simulationWork.useGpuNonbonded,
1476 deviceStreamManager.get(), &mtop, box, wcycle);
1477 // TODO: Move the logic below to a GPU bonded builder
1478 if (runScheduleWork.simulationWork.useGpuBonded)
1480 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1481 "GPU device stream manager should be valid in order to use GPU "
1482 "version of bonded forces.");
1483 gpuBonded = std::make_unique<GpuBonded>(
1484 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1485 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1486 fr->gpuBonded = gpuBonded.get();
1489 /* Initialize the mdAtoms structure.
1490 * mdAtoms is not filled with atom data,
1491 * as this can not be done now with domain decomposition.
1493 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1494 if (globalState && thisRankHasPmeGpuTask)
1496 // The pinning of coordinates in the global state object works, because we only use
1497 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1498 // points to the global state object without DD.
1499 // FIXME: MD and EM separately set up the local state - this should happen in the same
1500 // function, which should also perform the pinning.
1501 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1504 /* Initialize the virtual site communication */
1505 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1507 calc_shifts(box, fr->shift_vec);
1509 /* With periodic molecules the charge groups should be whole at start up
1510 * and the virtual sites should not be far from their proper positions.
1512 if (!inputrec->bContinuation && MASTER(cr)
1513 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1515 /* Make molecules whole at start of run */
1516 if (fr->pbcType != PbcType::No)
1518 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1522 /* Correct initial vsite positions are required
1523 * for the initial distribution in the domain decomposition
1524 * and for the initial shell prediction.
1526 constructVirtualSitesGlobal(mtop, globalState->x);
1530 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1532 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1533 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1538 /* This is a PME only node */
1540 GMX_ASSERT(globalState == nullptr,
1541 "We don't need the state on a PME only rank and expect it to be unitialized");
1543 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1544 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1547 gmx_pme_t* sepPmeData = nullptr;
1548 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1549 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1550 "Double-checking that only PME-only ranks have no forcerec");
1551 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1553 // TODO should live in ewald module once its testing is improved
1555 // Later, this program could contain kernels that might be later
1556 // re-used as auto-tuning progresses, or subsequent simulations
1558 PmeGpuProgramStorage pmeGpuProgram;
1559 if (thisRankHasPmeGpuTask)
1562 (deviceStreamManager != nullptr),
1563 "GPU device stream manager should be initialized in order to use GPU for PME.");
1564 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1565 "GPU device should be initialized in order to use GPU for PME.");
1566 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1569 /* Initiate PME if necessary,
1570 * either on all nodes or on dedicated PME nodes only. */
1571 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1573 if (mdAtoms && mdAtoms->mdatoms())
1575 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1576 if (EVDW_PME(inputrec->vdwtype))
1578 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1581 if (cr->npmenodes > 0)
1583 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1584 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1585 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1588 if (thisRankHasDuty(cr, DUTY_PME))
1592 // TODO: This should be in the builder.
1593 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1594 || (deviceStreamManager != nullptr),
1595 "Device stream manager should be valid in order to use GPU "
1598 !runScheduleWork.simulationWork.useGpuPme
1599 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1600 "GPU PME stream should be valid in order to use GPU version of PME.");
1602 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1603 ? &deviceStreamManager->context()
1605 const DeviceStream* pmeStream =
1606 runScheduleWork.simulationWork.useGpuPme
1607 ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1610 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1611 nChargePerturbed != 0, nTypePerturbed != 0,
1612 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1613 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1614 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1616 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1621 if (EI_DYNAMICS(inputrec->eI))
1623 /* Turn on signal handling on all nodes */
1625 * (A user signal from the PME nodes (if any)
1626 * is communicated to the PP nodes.
1628 signal_handler_install();
1631 pull_t* pull_work = nullptr;
1632 if (thisRankHasDuty(cr, DUTY_PP))
1634 /* Assumes uniform use of the number of OpenMP threads */
1635 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1637 if (inputrec->bPull)
1639 /* Initialize pull code */
1640 pull_work = init_pull(fplog, inputrec->pull.get(), inputrec.get(), &mtop, cr, &atomSets,
1641 inputrec->fepvals->init_lambda);
1642 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1644 initPullHistory(pull_work, &observablesHistory);
1646 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1648 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1652 std::unique_ptr<EnforcedRotation> enforcedRotation;
1655 /* Initialize enforced rotation code */
1656 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1657 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1661 t_swap* swap = nullptr;
1662 if (inputrec->eSwapCoords != eswapNO)
1664 /* Initialize ion swapping code */
1665 swap = init_swapcoords(fplog, inputrec.get(),
1666 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1667 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1668 oenv, mdrunOptions, startingBehavior);
1671 /* Let makeConstraints know whether we have essential dynamics constraints. */
1672 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1673 ms, &nrnb, wcycle, fr->bMolPBC);
1675 /* Energy terms and groups */
1676 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1677 inputrec->fepvals->n_lambda);
1679 // cos acceleration is only supported by md, but older tpr
1680 // files might still combine it with other integrators
1681 GMX_RELEASE_ASSERT(inputrec->cos_accel == 0.0 || inputrec->eI == eiMD,
1682 "cos_acceleration is only supported by integrator=md");
1684 /* Kinetic energy data */
1685 gmx_ekindata_t ekind;
1686 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1688 /* Set up interactive MD (IMD) */
1690 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1691 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1692 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1694 if (DOMAINDECOMP(cr))
1696 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1697 /* This call is not included in init_domain_decomposition mainly
1698 * because fr->cginfo_mb is set later.
1700 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1701 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1704 if (runScheduleWork.simulationWork.useGpuBufferOps)
1706 fr->gpuForceReduction[gmx::AtomLocality::Local] = std::make_unique<gmx::GpuForceReduction>(
1707 deviceStreamManager->context(),
1708 deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedLocal), wcycle);
1709 fr->gpuForceReduction[gmx::AtomLocality::NonLocal] = std::make_unique<gmx::GpuForceReduction>(
1710 deviceStreamManager->context(),
1711 deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedNonLocal), wcycle);
1714 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1715 if (gpusWereDetected
1716 && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
1717 || runScheduleWork.simulationWork.useGpuBufferOps))
1719 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1720 ? GpuApiCallBehavior::Async
1721 : GpuApiCallBehavior::Sync;
1722 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1723 "GPU device stream manager should be initialized to use GPU.");
1724 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1725 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1726 fr->stateGpu = stateGpu.get();
1729 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1730 SimulatorBuilder simulatorBuilder;
1732 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1733 simulatorBuilder.add(std::move(membedHolder));
1734 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1735 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1738 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1739 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1740 simulatorBuilder.add(ConstraintsParam(
1741 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1743 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1744 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1745 inputrec.get(), fr));
1746 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1747 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1748 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1749 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1750 // Todo move to an MDModule
1751 simulatorBuilder.add(IonSwapping(swap));
1752 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1753 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1754 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
1756 // build and run simulator object based on user-input
1757 auto simulator = simulatorBuilder.build(useModularSimulator);
1760 if (fr->pmePpCommGpu)
1762 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1763 fr->pmePpCommGpu.reset();
1766 if (inputrec->bPull)
1768 finish_pull(pull_work);
1770 finish_swapcoords(swap);
1774 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1776 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1777 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1778 deviceStreamManager.get());
1781 wallcycle_stop(wcycle, ewcRUN);
1783 /* Finish up, write some stuff
1784 * if rerunMD, don't write last frame again
1786 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1787 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1789 // clean up cycle counter
1790 wallcycle_destroy(wcycle);
1792 deviceStreamManager.reset(nullptr);
1796 gmx_pme_destroy(pmedata);
1800 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1801 // before we destroy the GPU context(s)
1802 // Pinned buffers are associated with contexts in CUDA.
1803 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1804 mdAtoms.reset(nullptr);
1805 globalState.reset(nullptr);
1806 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1807 gpuBonded.reset(nullptr);
1808 /* Free pinned buffers in *fr */
1811 // TODO convert to C++ so we can get rid of these frees
1815 if (!hwinfo_->deviceInfoList.empty())
1817 /* stop the GPU profiler (only CUDA) */
1821 /* With tMPI we need to wait for all ranks to finish deallocation before
1822 * destroying the CUDA context as some tMPI ranks may be sharing
1825 * This is not a concern in OpenCL where we use one context per rank.
1827 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1828 * but it is easier and more futureproof to call it on the whole node.
1830 * Note that this function needs to be called even if GPUs are not used
1831 * in this run because the PME ranks have no knowledge of whether GPUs
1832 * are used or not, but all ranks need to enter the barrier below.
1833 * \todo Remove this physical node barrier after making sure
1834 * that it's not needed anymore (with a shared GPU run).
1838 physicalNodeComm.barrier();
1840 releaseDevice(deviceInfo);
1842 /* Does what it says */
1843 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1844 walltime_accounting_destroy(walltime_accounting);
1846 // Ensure log file content is written
1849 gmx_fio_flush(logFileHandle);
1852 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1853 * exceptions were enabled before function was called. */
1856 gmx_fedisableexcept();
1859 auto rc = static_cast<int>(gmx_get_stop_condition());
1862 /* we need to join all threads. The sub-threads join when they
1863 exit this function, but the master thread needs to be told to
1873 Mdrunner::~Mdrunner()
1875 // Clean up of the Manager.
1876 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1877 // but okay as long as threads synchronize some time before adding or accessing
1878 // a new set of restraints.
1879 if (restraintManager_)
1881 restraintManager_->clear();
1882 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1883 "restraints added during runner life time should be cleared at runner "
1888 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1890 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1891 // Not sure if this should be logged through the md logger or something else,
1892 // but it is helpful to have some sort of INFO level message sent somewhere.
1893 // std::cout << "Registering restraint named " << name << std::endl;
1895 // When multiple restraints are used, it may be wasteful to register them separately.
1896 // Maybe instead register an entire Restraint Manager as a force provider.
1897 restraintManager_->addToSpec(std::move(puller), name);
1900 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1902 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1904 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265 in CentOS 7
1905 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1907 class Mdrunner::BuilderImplementation
1910 BuilderImplementation() = delete;
1911 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1912 ~BuilderImplementation();
1914 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1915 real forceWarningThreshold,
1916 StartingBehavior startingBehavior);
1918 void addHardwareDetectionResult(const gmx_hw_info_t* hwinfo);
1920 void addDomdec(const DomdecOptions& options);
1922 void addInput(SimulationInputHandle inputHolder);
1924 void addVerletList(int nstlist);
1926 void addReplicaExchange(const ReplicaExchangeParameters& params);
1928 void addNonBonded(const char* nbpu_opt);
1930 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1932 void addBondedTaskAssignment(const char* bonded_opt);
1934 void addUpdateTaskAssignment(const char* update_opt);
1936 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1938 void addFilenames(ArrayRef<const t_filenm> filenames);
1940 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1942 void addLogFile(t_fileio* logFileHandle);
1944 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1949 // Default parameters copied from runner.h
1950 // \todo Clarify source(s) of default parameters.
1952 const char* nbpu_opt_ = nullptr;
1953 const char* pme_opt_ = nullptr;
1954 const char* pme_fft_opt_ = nullptr;
1955 const char* bonded_opt_ = nullptr;
1956 const char* update_opt_ = nullptr;
1958 MdrunOptions mdrunOptions_;
1960 DomdecOptions domdecOptions_;
1962 ReplicaExchangeParameters replicaExchangeParameters_;
1964 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1967 //! World communicator, used for hardware detection and task assignment
1968 MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
1970 //! Multisim communicator handle.
1971 gmx_multisim_t* multiSimulation_;
1973 //! mdrun communicator
1974 MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
1976 //! Print a warning if any force is larger than this (in kJ/mol nm).
1977 real forceWarningThreshold_ = -1;
1979 //! Whether the simulation will start afresh, or restart with/without appending.
1980 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1982 //! The modules that comprise the functionality of mdrun.
1983 std::unique_ptr<MDModules> mdModules_;
1985 //! Detected hardware.
1986 const gmx_hw_info_t* hwinfo_ = nullptr;
1988 //! \brief Parallelism information.
1989 gmx_hw_opt_t hardwareOptions_;
1991 //! filename options for simulation.
1992 ArrayRef<const t_filenm> filenames_;
1994 /*! \brief Handle to output environment.
1996 * \todo gmx_output_env_t needs lifetime management.
1998 gmx_output_env_t* outputEnvironment_ = nullptr;
2000 /*! \brief Non-owning handle to MD log file.
2002 * \todo Context should own output facilities for client.
2003 * \todo Improve log file handle management.
2005 * Code managing the FILE* relies on the ability to set it to
2006 * nullptr to check whether the filehandle is valid.
2008 t_fileio* logFileHandle_ = nullptr;
2011 * \brief Builder for simulation stop signal handler.
2013 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
2016 * \brief Sources for initial simulation state.
2018 * See issue #3652 for near-term refinements to the SimulationInput interface.
2020 * See issue #3379 for broader discussion on API aspects of simulation inputs and outputs.
2022 SimulationInputHandle inputHolder_;
2025 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
2026 compat::not_null<SimulationContext*> context) :
2027 mdModules_(std::move(mdModules))
2029 libraryWorldCommunicator_ = context->libraryWorldCommunicator_;
2030 simulationCommunicator_ = context->simulationCommunicator_;
2031 multiSimulation_ = context->multiSimulation_.get();
2034 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
2036 Mdrunner::BuilderImplementation&
2037 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
2038 const real forceWarningThreshold,
2039 const StartingBehavior startingBehavior)
2041 mdrunOptions_ = options;
2042 forceWarningThreshold_ = forceWarningThreshold;
2043 startingBehavior_ = startingBehavior;
2047 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
2049 domdecOptions_ = options;
2052 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
2057 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
2059 replicaExchangeParameters_ = params;
2062 Mdrunner Mdrunner::BuilderImplementation::build()
2064 auto newRunner = Mdrunner(std::move(mdModules_));
2066 newRunner.mdrunOptions = mdrunOptions_;
2067 newRunner.pforce = forceWarningThreshold_;
2068 newRunner.startingBehavior = startingBehavior_;
2069 newRunner.domdecOptions = domdecOptions_;
2071 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2072 newRunner.hw_opt = hardwareOptions_;
2074 // No invariant to check. This parameter exists to optionally override other behavior.
2075 newRunner.nstlist_cmdline = nstlist_;
2077 newRunner.replExParams = replicaExchangeParameters_;
2079 newRunner.filenames = filenames_;
2081 newRunner.libraryWorldCommunicator = libraryWorldCommunicator_;
2083 newRunner.simulationCommunicator = simulationCommunicator_;
2085 // nullptr is a valid value for the multisim handle
2086 newRunner.ms = multiSimulation_;
2090 newRunner.hwinfo_ = hwinfo_;
2094 GMX_THROW(gmx::APIError(
2095 "MdrunnerBuilder::addHardwareDetectionResult() is required before build()"));
2100 newRunner.inputHolder_ = std::move(inputHolder_);
2104 GMX_THROW(gmx::APIError("MdrunnerBuilder::addInput() is required before build()."));
2107 // \todo Clarify ownership and lifetime management for gmx_output_env_t
2108 // \todo Update sanity checking when output environment has clearly specified invariants.
2109 // Initialization and default values for oenv are not well specified in the current version.
2110 if (outputEnvironment_)
2112 newRunner.oenv = outputEnvironment_;
2116 GMX_THROW(gmx::APIError(
2117 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2120 newRunner.logFileHandle = logFileHandle_;
2124 newRunner.nbpu_opt = nbpu_opt_;
2128 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2131 if (pme_opt_ && pme_fft_opt_)
2133 newRunner.pme_opt = pme_opt_;
2134 newRunner.pme_fft_opt = pme_fft_opt_;
2138 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2143 newRunner.bonded_opt = bonded_opt_;
2147 GMX_THROW(gmx::APIError(
2148 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2153 newRunner.update_opt = update_opt_;
2157 GMX_THROW(gmx::APIError(
2158 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2162 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2164 if (stopHandlerBuilder_)
2166 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2170 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2176 void Mdrunner::BuilderImplementation::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2181 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2183 nbpu_opt_ = nbpu_opt;
2186 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2189 pme_fft_opt_ = pme_fft_opt;
2192 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2194 bonded_opt_ = bonded_opt;
2197 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2199 update_opt_ = update_opt;
2202 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2204 hardwareOptions_ = hardwareOptions;
2207 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2209 filenames_ = filenames;
2212 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2214 outputEnvironment_ = outputEnvironment;
2217 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2219 logFileHandle_ = logFileHandle;
2222 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2224 stopHandlerBuilder_ = std::move(builder);
2227 void Mdrunner::BuilderImplementation::addInput(SimulationInputHandle inputHolder)
2229 inputHolder_ = std::move(inputHolder);
2232 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2233 compat::not_null<SimulationContext*> context) :
2234 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2238 MdrunnerBuilder::~MdrunnerBuilder() = default;
2240 MdrunnerBuilder& MdrunnerBuilder::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2242 impl_->addHardwareDetectionResult(hwinfo);
2246 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2247 real forceWarningThreshold,
2248 const StartingBehavior startingBehavior)
2250 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2254 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2256 impl_->addDomdec(options);
2260 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2262 impl_->addVerletList(nstlist);
2266 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2268 impl_->addReplicaExchange(params);
2272 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2274 impl_->addNonBonded(nbpu_opt);
2278 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2280 // The builder method may become more general in the future, but in this version,
2281 // parameters for PME electrostatics are both required and the only parameters
2283 if (pme_opt && pme_fft_opt)
2285 impl_->addPME(pme_opt, pme_fft_opt);
2290 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2295 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2297 impl_->addBondedTaskAssignment(bonded_opt);
2301 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2303 impl_->addUpdateTaskAssignment(update_opt);
2307 Mdrunner MdrunnerBuilder::build()
2309 return impl_->build();
2312 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2314 impl_->addHardwareOptions(hardwareOptions);
2318 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2320 impl_->addFilenames(filenames);
2324 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2326 impl_->addOutputEnvironment(outputEnvironment);
2330 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2332 impl_->addLogFile(logFileHandle);
2336 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2338 impl_->addStopHandlerBuilder(std::move(builder));
2342 MdrunnerBuilder& MdrunnerBuilder::addInput(SimulationInputHandle input)
2344 impl_->addInput(std::move(input));
2348 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2350 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;