2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/device_context.h"
77 #include "gromacs/gpu_utils/device_stream_manager.h"
78 #include "gromacs/hardware/cpuinfo.h"
79 #include "gromacs/hardware/detecthardware.h"
80 #include "gromacs/hardware/device_management.h"
81 #include "gromacs/hardware/printhardware.h"
82 #include "gromacs/imd/imd.h"
83 #include "gromacs/listed_forces/disre.h"
84 #include "gromacs/listed_forces/gpubonded.h"
85 #include "gromacs/listed_forces/listed_forces.h"
86 #include "gromacs/listed_forces/orires.h"
87 #include "gromacs/math/functions.h"
88 #include "gromacs/math/utilities.h"
89 #include "gromacs/math/vec.h"
90 #include "gromacs/mdlib/boxdeformation.h"
91 #include "gromacs/mdlib/broadcaststructs.h"
92 #include "gromacs/mdlib/calc_verletbuf.h"
93 #include "gromacs/mdlib/dispersioncorrection.h"
94 #include "gromacs/mdlib/enerdata_utils.h"
95 #include "gromacs/mdlib/force.h"
96 #include "gromacs/mdlib/forcerec.h"
97 #include "gromacs/mdlib/gmx_omp_nthreads.h"
98 #include "gromacs/mdlib/makeconstraints.h"
99 #include "gromacs/mdlib/md_support.h"
100 #include "gromacs/mdlib/mdatoms.h"
101 #include "gromacs/mdlib/sighandler.h"
102 #include "gromacs/mdlib/stophandler.h"
103 #include "gromacs/mdlib/tgroup.h"
104 #include "gromacs/mdlib/updategroups.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdrun/mdmodules.h"
107 #include "gromacs/mdrun/simulationcontext.h"
108 #include "gromacs/mdrunutility/handlerestart.h"
109 #include "gromacs/mdrunutility/logging.h"
110 #include "gromacs/mdrunutility/multisim.h"
111 #include "gromacs/mdrunutility/printtime.h"
112 #include "gromacs/mdrunutility/threadaffinity.h"
113 #include "gromacs/mdtypes/checkpointdata.h"
114 #include "gromacs/mdtypes/commrec.h"
115 #include "gromacs/mdtypes/enerdata.h"
116 #include "gromacs/mdtypes/fcdata.h"
117 #include "gromacs/mdtypes/forcerec.h"
118 #include "gromacs/mdtypes/group.h"
119 #include "gromacs/mdtypes/inputrec.h"
120 #include "gromacs/mdtypes/interaction_const.h"
121 #include "gromacs/mdtypes/md_enums.h"
122 #include "gromacs/mdtypes/mdatom.h"
123 #include "gromacs/mdtypes/mdrunoptions.h"
124 #include "gromacs/mdtypes/observableshistory.h"
125 #include "gromacs/mdtypes/simulation_workload.h"
126 #include "gromacs/mdtypes/state.h"
127 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
128 #include "gromacs/modularsimulator/modularsimulator.h"
129 #include "gromacs/nbnxm/gpu_data_mgmt.h"
130 #include "gromacs/nbnxm/nbnxm.h"
131 #include "gromacs/nbnxm/pairlist_tuning.h"
132 #include "gromacs/pbcutil/pbc.h"
133 #include "gromacs/pulling/output.h"
134 #include "gromacs/pulling/pull.h"
135 #include "gromacs/pulling/pull_rotation.h"
136 #include "gromacs/restraint/manager.h"
137 #include "gromacs/restraint/restraintmdmodule.h"
138 #include "gromacs/restraint/restraintpotential.h"
139 #include "gromacs/swap/swapcoords.h"
140 #include "gromacs/taskassignment/decidegpuusage.h"
141 #include "gromacs/taskassignment/decidesimulationworkload.h"
142 #include "gromacs/taskassignment/resourcedivision.h"
143 #include "gromacs/taskassignment/taskassignment.h"
144 #include "gromacs/taskassignment/usergpuids.h"
145 #include "gromacs/timing/gpu_timing.h"
146 #include "gromacs/timing/wallcycle.h"
147 #include "gromacs/timing/wallcyclereporting.h"
148 #include "gromacs/topology/mtop_util.h"
149 #include "gromacs/trajectory/trajectoryframe.h"
150 #include "gromacs/utility/basenetwork.h"
151 #include "gromacs/utility/cstringutil.h"
152 #include "gromacs/utility/exceptions.h"
153 #include "gromacs/utility/fatalerror.h"
154 #include "gromacs/utility/filestream.h"
155 #include "gromacs/utility/gmxassert.h"
156 #include "gromacs/utility/gmxmpi.h"
157 #include "gromacs/utility/keyvaluetree.h"
158 #include "gromacs/utility/logger.h"
159 #include "gromacs/utility/loggerbuilder.h"
160 #include "gromacs/utility/mdmodulenotification.h"
161 #include "gromacs/utility/physicalnodecommunicator.h"
162 #include "gromacs/utility/pleasecite.h"
163 #include "gromacs/utility/programcontext.h"
164 #include "gromacs/utility/smalloc.h"
165 #include "gromacs/utility/stringutil.h"
167 #include "isimulator.h"
168 #include "membedholder.h"
169 #include "replicaexchange.h"
170 #include "simulatorbuilder.h"
173 # include "corewrap.h"
180 /*! \brief Manage any development feature flag variables encountered
182 * The use of dev features indicated by environment variables is
183 * logged in order to ensure that runs with such features enabled can
184 * be identified from their log and standard output. Any cross
185 * dependencies are also checked, and if unsatisfied, a fatal error
188 * Note that some development features overrides are applied already here:
189 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
191 * \param[in] mdlog Logger object.
192 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
193 * \param[in] pmeRunMode The PME run mode for this run
194 * \returns The object populated with development feature flags.
196 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
197 const bool useGpuForNonbonded,
198 const PmeRunMode pmeRunMode)
200 DevelopmentFeatureFlags devFlags;
202 // Some builds of GCC 5 give false positive warnings that these
203 // getenv results are ignored when clearly they are used.
204 #pragma GCC diagnostic push
205 #pragma GCC diagnostic ignored "-Wunused-result"
207 devFlags.enableGpuBufferOps =
208 GMX_GPU_CUDA && useGpuForNonbonded && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
209 devFlags.enableGpuHaloExchange = GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_DD_COMMS") != nullptr;
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.oenv = oenv;
340 newRunner.mdrunOptions = mdrunOptions;
341 newRunner.domdecOptions = domdecOptions;
342 newRunner.nbpu_opt = nbpu_opt;
343 newRunner.pme_opt = pme_opt;
344 newRunner.pme_fft_opt = pme_fft_opt;
345 newRunner.bonded_opt = bonded_opt;
346 newRunner.update_opt = update_opt;
347 newRunner.nstlist_cmdline = nstlist_cmdline;
348 newRunner.replExParams = replExParams;
349 newRunner.pforce = pforce;
350 // Give the spawned thread the newly created valid communicator
351 // for the simulation.
352 newRunner.communicator = MPI_COMM_WORLD;
354 newRunner.startingBehavior = startingBehavior;
355 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
357 threadMpiMdrunnerAccessBarrier();
362 /*! \brief The callback used for running on spawned threads.
364 * Obtains the pointer to the master mdrunner object from the one
365 * argument permitted to the thread-launch API call, copies it to make
366 * a new runner for this thread, reinitializes necessary data, and
367 * proceeds to the simulation. */
368 static void mdrunner_start_fn(const void* arg)
372 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
373 /* copy the arg list to make sure that it's thread-local. This
374 doesn't copy pointed-to items, of course; fnm, cr and fplog
375 are reset in the call below, all others should be const. */
376 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
379 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
383 void Mdrunner::spawnThreads(int numThreadsToLaunch)
386 /* now spawn new threads that start mdrunner_start_fn(), while
387 the main thread returns. Thread affinity is handled later. */
388 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
389 static_cast<const void*>(this))
392 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
395 // Give the master thread the newly created valid communicator for
397 communicator = MPI_COMM_WORLD;
398 threadMpiMdrunnerAccessBarrier();
400 GMX_UNUSED_VALUE(numThreadsToLaunch);
401 GMX_UNUSED_VALUE(mdrunner_start_fn);
407 /*! \brief Initialize variables for Verlet scheme simulation */
408 static void prepare_verlet_scheme(FILE* fplog,
412 const gmx_mtop_t* mtop,
414 bool makeGpuPairList,
415 const gmx::CpuInfo& cpuinfo)
417 // We checked the cut-offs in grompp, but double-check here.
418 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
419 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
421 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
422 "With Verlet lists and PME we should have rcoulomb>=rvdw");
426 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
427 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
429 /* For NVE simulations, we will retain the initial list buffer */
430 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
432 /* Update the Verlet buffer size for the current run setup */
434 /* Here we assume SIMD-enabled kernels are being used. But as currently
435 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
436 * and 4x2 gives a larger buffer than 4x4, this is ok.
438 ListSetupType listType =
439 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
440 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
442 const real rlist_new =
443 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
445 if (rlist_new != ir->rlist)
447 if (fplog != nullptr)
450 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
451 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
453 ir->rlist = rlist_new;
457 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
459 gmx_fatal(FARGS, "Can not set nstlist without %s",
460 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
463 if (EI_DYNAMICS(ir->eI))
465 /* Set or try nstlist values */
466 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
470 /*! \brief Override the nslist value in inputrec
472 * with value passed on the command line (if any)
474 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
478 /* override with anything else than the default -2 */
479 if (nsteps_cmdline > -2)
481 char sbuf_steps[STEPSTRSIZE];
482 char sbuf_msg[STRLEN];
484 ir->nsteps = nsteps_cmdline;
485 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
488 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
489 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
493 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
494 gmx_step_str(nsteps_cmdline, sbuf_steps));
497 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
499 else if (nsteps_cmdline < -2)
501 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
503 /* Do nothing if nsteps_cmdline == -2 */
509 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
511 * If not, and if a warning may be issued, logs a warning about
512 * falling back to CPU code. With thread-MPI, only the first
513 * call to this function should have \c issueWarning true. */
514 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
516 bool gpuIsUseful = true;
519 if (ir.opts.ngener - ir.nwall > 1)
521 /* The GPU code does not support more than one energy group.
522 * If the user requested GPUs explicitly, a fatal error is given later.
526 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
527 "For better performance, run on the GPU without energy groups and then do "
528 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
534 warning = "TPI is not implemented for GPUs.";
537 if (!gpuIsUseful && issueWarning)
539 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
545 //! Initializes the logger for mdrun.
546 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
548 gmx::LoggerBuilder builder;
549 if (fplog != nullptr)
551 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
553 if (isSimulationMasterRank)
555 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
557 return builder.build();
560 //! Make a TaskTarget from an mdrun argument string.
561 static TaskTarget findTaskTarget(const char* optionString)
563 TaskTarget returnValue = TaskTarget::Auto;
565 if (strncmp(optionString, "auto", 3) == 0)
567 returnValue = TaskTarget::Auto;
569 else if (strncmp(optionString, "cpu", 3) == 0)
571 returnValue = TaskTarget::Cpu;
573 else if (strncmp(optionString, "gpu", 3) == 0)
575 returnValue = TaskTarget::Gpu;
579 GMX_ASSERT(false, "Option string should have been checked for sanity already");
585 //! Finish run, aggregate data to print performance info.
586 static void finish_run(FILE* fplog,
587 const gmx::MDLogger& mdlog,
589 const t_inputrec* inputrec,
591 gmx_wallcycle_t wcycle,
592 gmx_walltime_accounting_t walltime_accounting,
593 nonbonded_verlet_t* nbv,
594 const gmx_pme_t* pme,
598 double nbfs = 0, mflop = 0;
599 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
600 elapsed_time_over_all_threads_over_all_ranks;
601 /* Control whether it is valid to print a report. Only the
602 simulation master may print, but it should not do so if the run
603 terminated e.g. before a scheduled reset step. This is
604 complicated by the fact that PME ranks are unaware of the
605 reason why they were sent a pmerecvqxFINISH. To avoid
606 communication deadlocks, we always do the communication for the
607 report, even if we've decided not to write the report, because
608 how long it takes to finish the run is not important when we've
609 decided not to report on the simulation performance.
611 Further, we only report performance for dynamical integrators,
612 because those are the only ones for which we plan to
613 consider doing any optimizations. */
614 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
616 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
618 GMX_LOG(mdlog.warning)
620 .appendText("Simulation ended prematurely, no performance report will be written.");
625 std::unique_ptr<t_nrnb> nrnbTotalStorage;
628 nrnbTotalStorage = std::make_unique<t_nrnb>();
629 nrnb_tot = nrnbTotalStorage.get();
631 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
639 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
640 elapsed_time_over_all_threads =
641 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
645 /* reduce elapsed_time over all MPI ranks in the current simulation */
646 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
648 elapsed_time_over_all_ranks /= cr->nnodes;
649 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
650 * current simulation. */
651 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
652 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
657 elapsed_time_over_all_ranks = elapsed_time;
658 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
663 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
666 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
668 print_dd_statistics(cr, inputrec, fplog);
671 /* TODO Move the responsibility for any scaling by thread counts
672 * to the code that handled the thread region, so that there's a
673 * mechanism to keep cycle counting working during the transition
674 * to task parallelism. */
675 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
676 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
677 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
678 nthreads_pp, nthreads_pme);
679 auto cycle_sum(wallcycle_sum(cr, wcycle));
683 auto nbnxn_gpu_timings =
684 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
685 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
687 if (pme_gpu_task_enabled(pme))
689 pme_gpu_get_timings(pme, &pme_gpu_timings);
691 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
692 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
695 if (EI_DYNAMICS(inputrec->eI))
697 delta_t = inputrec->delta_t;
702 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
703 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
704 delta_t, nbfs, mflop);
708 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
709 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
710 delta_t, nbfs, mflop);
715 int Mdrunner::mdrunner()
718 t_forcerec* fr = nullptr;
719 real ewaldcoeff_q = 0;
720 real ewaldcoeff_lj = 0;
721 int nChargePerturbed = -1, nTypePerturbed = 0;
722 gmx_wallcycle_t wcycle;
723 gmx_walltime_accounting_t walltime_accounting = nullptr;
724 MembedHolder membedHolder(filenames.size(), filenames.data());
725 gmx_hw_info_t* hwinfo = nullptr;
727 /* CAUTION: threads may be started later on in this function, so
728 cr doesn't reflect the final parallel state right now */
731 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
732 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
733 const bool doRerun = mdrunOptions.rerun;
735 // Handle task-assignment related user options.
736 EmulateGpuNonbonded emulateGpuNonbonded =
737 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
739 std::vector<int> userGpuTaskAssignment;
742 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
744 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
745 auto nonbondedTarget = findTaskTarget(nbpu_opt);
746 auto pmeTarget = findTaskTarget(pme_opt);
747 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
748 auto bondedTarget = findTaskTarget(bonded_opt);
749 auto updateTarget = findTaskTarget(update_opt);
751 FILE* fplog = nullptr;
752 // If we are appending, we don't write log output because we need
753 // to check that the old log file matches what the checkpoint file
754 // expects. Otherwise, we should start to write log output now if
755 // there is a file ready for it.
756 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
758 fplog = gmx_fio_getfp(logFileHandle);
760 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
761 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
762 gmx::MDLogger mdlog(logOwner.logger());
764 // TODO The thread-MPI master rank makes a working
765 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
766 // after the threads have been launched. This works because no use
767 // is made of that communicator until after the execution paths
768 // have rejoined. But it is likely that we can improve the way
769 // this is expressed, e.g. by expressly running detection only the
770 // master rank for thread-MPI, rather than relying on the mutex
771 // and reference count.
772 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
773 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
775 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
777 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->deviceInfoList, hw_opt.gpuIdsAvailable);
779 // Print citation requests after all software/hardware printing
780 pleaseCiteGromacs(fplog);
782 // Note: legacy program logic relies on checking whether these pointers are assigned.
783 // Objects may or may not be allocated later.
784 std::unique_ptr<t_inputrec> inputrec;
785 std::unique_ptr<t_state> globalState;
787 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
789 if (isSimulationMasterRank)
791 // Allocate objects to be initialized by later function calls.
792 /* Only the master rank has the global state */
793 globalState = std::make_unique<t_state>();
794 inputrec = std::make_unique<t_inputrec>();
796 /* Read (nearly) all data required for the simulation
797 * and keep the partly serialized tpr contents to send to other ranks later
799 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
800 inputrec.get(), globalState.get(), &mtop);
803 /* Check and update the hardware options for internal consistency */
804 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
807 if (GMX_THREAD_MPI && isSimulationMasterRank)
809 bool useGpuForNonbonded = false;
810 bool useGpuForPme = false;
813 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
815 // If the user specified the number of ranks, then we must
816 // respect that, but in default mode, we need to allow for
817 // the number of GPUs to choose the number of ranks.
818 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
819 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
820 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
821 canUseGpuForNonbonded,
822 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
823 hw_opt.nthreads_tmpi);
824 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
825 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
826 *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
828 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
830 /* Determine how many thread-MPI ranks to start.
832 * TODO Over-writing the user-supplied value here does
833 * prevent any possible subsequent checks from working
835 hw_opt.nthreads_tmpi =
836 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
837 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
839 // Now start the threads for thread MPI.
840 spawnThreads(hw_opt.nthreads_tmpi);
841 // The spawned threads enter mdrunner() and execution of
842 // master and spawned threads joins at the end of this block.
843 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
846 GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
847 "Must have valid world communicator unless running a multi-simulation");
848 CommrecHandle crHandle = init_commrec(communicator);
849 t_commrec* cr = crHandle.get();
850 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
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, userGpuTaskAssignment, *hwinfo, *inputrec,
895 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
896 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
897 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
898 useGpuForBonded = decideWhetherToUseGpusForBonded(
899 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
900 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
901 domdecOptions.numPmeRanks, gpusWereDetected);
903 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
905 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
907 // Initialize development feature flags that enabled by environment variable
908 // and report those features that are enabled.
909 const DevelopmentFeatureFlags devFlags =
910 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
912 const bool useModularSimulator =
913 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
914 nullptr, doEssentialDynamics, membedHolder.doMembed());
917 // TODO: hide restraint implementation details from Mdrunner.
918 // There is nothing unique about restraints at this point as far as the
919 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
920 // factory functions from the SimulationContext on which to call mdModules_->add().
921 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
922 for (auto&& restraint : restraintManager_->getRestraints())
924 auto module = RestraintMDModule::create(restraint, restraint->sites());
925 mdModules_->add(std::move(module));
928 // TODO: Error handling
929 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
930 // now that the MdModules know their options, they know which callbacks to sign up to
931 mdModules_->subscribeToSimulationSetupNotifications();
932 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
934 if (inputrec->internalParameters != nullptr)
936 mdModulesNotifier.notify(*inputrec->internalParameters);
939 if (fplog != nullptr)
941 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
942 fprintf(fplog, "\n");
947 /* In rerun, set velocities to zero if present */
948 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
950 // rerun does not use velocities
954 "Rerun trajectory contains velocities. Rerun does only evaluate "
955 "potential energy and forces. The velocities will be ignored.");
956 for (int i = 0; i < globalState->natoms; i++)
958 clear_rvec(globalState->v[i]);
960 globalState->flags &= ~(1 << estV);
963 /* now make sure the state is initialized and propagated */
964 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
967 /* NM and TPI parallelize over force/energy calculations, not atoms,
968 * so we need to initialize and broadcast the global state.
970 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
974 globalState = std::make_unique<t_state>();
976 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
980 /* A parallel command line option consistency check that we can
981 only do after any threads have started. */
983 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
984 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
987 "The -dd or -npme option request a parallel simulation, "
989 "but %s was compiled without threads or MPI enabled",
990 output_env_get_program_display_name(oenv));
992 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
994 "but %s was not started through mpirun/mpiexec or only one rank was requested "
995 "through mpirun/mpiexec",
996 output_env_get_program_display_name(oenv));
1000 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1003 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1004 "these are not compatible with mdrun -rerun");
1007 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1009 if (domdecOptions.numPmeRanks > 0)
1011 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1012 "PME-only ranks are requested, but the system does not use PME "
1013 "for electrostatics or LJ");
1016 domdecOptions.numPmeRanks = 0;
1019 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1021 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1022 * improve performance with many threads per GPU, since our OpenMP
1023 * scaling is bad, but it's difficult to automate the setup.
1025 domdecOptions.numPmeRanks = 0;
1029 if (domdecOptions.numPmeRanks < 0)
1031 domdecOptions.numPmeRanks = 0;
1032 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1036 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1037 "PME GPU decomposition is not supported");
1044 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1048 /* NMR restraints must be initialized before load_checkpoint,
1049 * since with time averaging the history is added to t_state.
1050 * For proper consistency check we therefore need to extend
1052 * So the PME-only nodes (if present) will also initialize
1053 * the distance restraints.
1056 /* This needs to be called before read_checkpoint to extend the state */
1057 t_disresdata* disresdata;
1058 snew(disresdata, 1);
1059 init_disres(fplog, &mtop, inputrec.get(), DisResRunMode::MDRun,
1060 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1061 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1062 globalState.get(), replExParams.exchangeInterval > 0);
1064 t_oriresdata* oriresdata;
1065 snew(oriresdata, 1);
1066 init_orires(fplog, &mtop, inputrec.get(), cr, ms, globalState.get(), oriresdata);
1068 auto deform = prepareBoxDeformation(
1069 globalState != nullptr ? globalState->box : box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1070 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mygroup, *inputrec);
1072 ObservablesHistory observablesHistory = {};
1074 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1075 if (startingBehavior != StartingBehavior::NewSimulation)
1077 /* Check if checkpoint file exists before doing continuation.
1078 * This way we can use identical input options for the first and subsequent runs...
1080 if (mdrunOptions.numStepsCommandline > -2)
1082 /* Temporarily set the number of steps to unlimited to avoid
1083 * triggering the nsteps check in load_checkpoint().
1084 * This hack will go away soon when the -nsteps option is removed.
1086 inputrec->nsteps = -1;
1089 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1090 logFileHandle, cr, domdecOptions.numCells, inputrec.get(), globalState.get(),
1091 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier(),
1092 modularSimulatorCheckpointData.get(), useModularSimulator);
1094 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1096 // Now we can start normal logging to the truncated log file.
1097 fplog = gmx_fio_getfp(logFileHandle);
1098 prepareLogAppending(fplog);
1099 logOwner = buildLogger(fplog, MASTER(cr));
1100 mdlog = logOwner.logger();
1104 if (mdrunOptions.numStepsCommandline > -2)
1109 "The -nsteps functionality is deprecated, and may be removed in a future "
1111 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1114 /* override nsteps with value set on the commandline */
1115 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1117 if (isSimulationMasterRank)
1119 copy_mat(globalState->box, box);
1124 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1127 if (inputrec->cutoff_scheme != ecutsVERLET)
1130 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1131 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1133 /* Update rlist and nstlist. */
1134 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1135 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1136 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1138 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1139 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1142 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1143 // This builder is necessary while we have multi-part construction
1144 // of DD. Before DD is constructed, we use the existence of
1145 // the builder object to indicate that further construction of DD
1147 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1148 if (useDomainDecomposition)
1150 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1151 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1152 positionsFromStatePointer(globalState.get()));
1156 /* PME, if used, is done on all nodes with 1D decomposition */
1157 cr->nnodes = cr->sizeOfDefaultCommunicator;
1158 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1159 cr->nodeid = cr->rankInDefaultCommunicator;
1161 cr->duty = (DUTY_PP | DUTY_PME);
1163 if (inputrec->pbcType == PbcType::Screw)
1165 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1169 // Produce the task assignment for this rank - done after DD is constructed
1170 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1171 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1172 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1173 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1174 // TODO cr->duty & DUTY_PME should imply that a PME
1175 // algorithm is active, but currently does not.
1176 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1178 // Get the device handles for the modules, nullptr when no task is assigned.
1180 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1182 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1183 bool useTiming = true;
1187 /* WARNING: CUDA timings are incorrect with multiple streams.
1188 * This is the main reason why they are disabled by default.
1190 // TODO: Consider turning on by default when we can detect nr of streams.
1191 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1193 else if (GMX_GPU_OPENCL)
1195 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1198 // TODO Currently this is always built, yet DD partition code
1199 // checks if it is built before using it. Probably it should
1200 // become an MDModule that is made only when another module
1201 // requires it (e.g. pull, CompEl, density fitting), so that we
1202 // don't update the local atom sets unilaterally every step.
1203 LocalAtomSetManager atomSets;
1206 // TODO Pass the GPU streams to ddBuilder to use in buffer
1207 // transfers (e.g. halo exchange)
1208 cr->dd = ddBuilder->build(&atomSets);
1209 // The builder's job is done, so destruct it
1210 ddBuilder.reset(nullptr);
1211 // Note that local state still does not exist yet.
1214 // The GPU update is decided here because we need to know whether the constraints or
1215 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1216 // defined). This is only known after DD is initialized, hence decision on using GPU
1217 // update is done so late.
1220 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1222 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1223 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1224 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1225 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1226 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1228 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1230 const bool printHostName = (cr->nnodes > 1);
1231 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1233 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1235 if (deviceInfo != nullptr)
1237 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1239 dd_setup_dlb_resource_sharing(cr, deviceId);
1241 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1242 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1243 useGpuForUpdate, useTiming);
1246 // If the user chose a task assignment, give them some hints
1247 // where appropriate.
1248 if (!userGpuTaskAssignment.empty())
1250 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1255 /* After possible communicator splitting in make_dd_communicators.
1256 * we can set up the intra/inter node communication.
1258 gmx_setup_nodecomm(fplog, cr);
1264 GMX_LOG(mdlog.warning)
1266 .appendTextFormatted(
1267 "This is simulation %d out of %d running as a composite GROMACS\n"
1268 "multi-simulation job. Setup for this simulation:\n",
1269 ms->simulationIndex_, ms->numSimulations_);
1271 GMX_LOG(mdlog.warning)
1272 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1274 cr->nnodes == 1 ? "thread" : "threads"
1276 cr->nnodes == 1 ? "process" : "processes"
1282 // If mdrun -pin auto honors any affinity setting that already
1283 // exists. If so, it is nice to provide feedback about whether
1284 // that existing affinity setting was from OpenMP or something
1285 // else, so we run this code both before and after we initialize
1286 // the OpenMP support.
1287 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1288 /* Check and update the number of OpenMP threads requested */
1289 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1290 pmeRunMode, mtop, *inputrec);
1292 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1293 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1295 // Enable FP exception detection, but not in
1296 // Release mode and not for compilers with known buggy FP
1297 // exception support (clang with any optimization) or suspected
1298 // buggy FP exception support (gcc 7.* with optimization).
1299 #if !defined NDEBUG \
1300 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1301 && defined __OPTIMIZE__)
1302 const bool bEnableFPE = true;
1304 const bool bEnableFPE = false;
1306 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1309 gmx_feenableexcept();
1312 /* Now that we know the setup is consistent, check for efficiency */
1313 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1314 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1316 /* getting number of PP/PME threads on this MPI / tMPI rank.
1317 PME: env variable should be read only on one node to make sure it is
1318 identical everywhere;
1320 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1321 : gmx_omp_nthreads_get(emntPME);
1322 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1323 physicalNodeComm, mdlog);
1325 // Enable Peer access between GPUs where available
1326 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1327 // any of the GPU communication features are active.
1328 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1329 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1331 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1334 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1336 /* Before setting affinity, check whether the affinity has changed
1337 * - which indicates that probably the OpenMP library has changed it
1338 * since we first checked).
1340 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1342 int numThreadsOnThisNode, intraNodeThreadOffset;
1343 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1344 &intraNodeThreadOffset);
1346 /* Set the CPU affinity */
1347 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1348 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1351 if (mdrunOptions.timingOptions.resetStep > -1)
1356 "The -resetstep functionality is deprecated, and may be removed in a "
1359 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1363 /* Master synchronizes its value of reset_counters with all nodes
1364 * including PME only nodes */
1365 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1366 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1367 wcycle_set_reset_counters(wcycle, reset_counters);
1370 // Membrane embedding must be initialized before we call init_forcerec()
1371 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1372 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1374 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1375 std::unique_ptr<MDAtoms> mdAtoms;
1376 std::unique_ptr<VirtualSitesHandler> vsite;
1377 std::unique_ptr<GpuBonded> gpuBonded;
1380 if (thisRankHasDuty(cr, DUTY_PP))
1382 mdModulesNotifier.notify(*cr);
1383 mdModulesNotifier.notify(&atomSets);
1384 mdModulesNotifier.notify(inputrec->pbcType);
1385 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1386 /* Initiate forcerecord */
1387 fr = new t_forcerec;
1388 fr->forceProviders = mdModules_->initForceProviders();
1389 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1390 opt2fn("-table", filenames.size(), filenames.data()),
1391 opt2fn("-tablep", filenames.size(), filenames.data()),
1392 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1393 // Dirty hack, for fixing disres and orires should be made mdmodules
1394 fr->fcdata->disres = disresdata;
1395 fr->fcdata->orires = oriresdata;
1397 // Save a handle to device stream manager to use elsewhere in the code
1398 // TODO: Forcerec is not a correct place to store it.
1399 fr->deviceStreamManager = deviceStreamManager.get();
1401 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1404 deviceStreamManager != nullptr,
1405 "GPU device stream manager should be valid in order to use PME-PP direct "
1408 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1409 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1411 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1412 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1413 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1416 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo, useGpuForNonbonded,
1417 deviceStreamManager.get(), &mtop, box, wcycle);
1418 // TODO: Move the logic below to a GPU bonded builder
1419 if (useGpuForBonded)
1421 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1422 "GPU device stream manager should be valid in order to use GPU "
1423 "version of bonded forces.");
1424 gpuBonded = std::make_unique<GpuBonded>(
1425 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1426 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1427 fr->gpuBonded = gpuBonded.get();
1430 /* Initialize the mdAtoms structure.
1431 * mdAtoms is not filled with atom data,
1432 * as this can not be done now with domain decomposition.
1434 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1435 if (globalState && thisRankHasPmeGpuTask)
1437 // The pinning of coordinates in the global state object works, because we only use
1438 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1439 // points to the global state object without DD.
1440 // FIXME: MD and EM separately set up the local state - this should happen in the same
1441 // function, which should also perform the pinning.
1442 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1445 /* Initialize the virtual site communication */
1446 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1448 calc_shifts(box, fr->shift_vec);
1450 /* With periodic molecules the charge groups should be whole at start up
1451 * and the virtual sites should not be far from their proper positions.
1453 if (!inputrec->bContinuation && MASTER(cr)
1454 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1456 /* Make molecules whole at start of run */
1457 if (fr->pbcType != PbcType::No)
1459 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1463 /* Correct initial vsite positions are required
1464 * for the initial distribution in the domain decomposition
1465 * and for the initial shell prediction.
1467 constructVirtualSitesGlobal(mtop, globalState->x);
1471 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1473 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1474 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1479 /* This is a PME only node */
1481 GMX_ASSERT(globalState == nullptr,
1482 "We don't need the state on a PME only rank and expect it to be unitialized");
1484 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1485 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1488 gmx_pme_t* sepPmeData = nullptr;
1489 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1490 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1491 "Double-checking that only PME-only ranks have no forcerec");
1492 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1494 // TODO should live in ewald module once its testing is improved
1496 // Later, this program could contain kernels that might be later
1497 // re-used as auto-tuning progresses, or subsequent simulations
1499 PmeGpuProgramStorage pmeGpuProgram;
1500 if (thisRankHasPmeGpuTask)
1503 (deviceStreamManager != nullptr),
1504 "GPU device stream manager should be initialized in order to use GPU for PME.");
1505 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1506 "GPU device should be initialized in order to use GPU for PME.");
1507 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1510 /* Initiate PME if necessary,
1511 * either on all nodes or on dedicated PME nodes only. */
1512 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1514 if (mdAtoms && mdAtoms->mdatoms())
1516 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1517 if (EVDW_PME(inputrec->vdwtype))
1519 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1522 if (cr->npmenodes > 0)
1524 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1525 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1526 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1529 if (thisRankHasDuty(cr, DUTY_PME))
1533 // TODO: This should be in the builder.
1534 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1535 "Device stream manager should be valid in order to use GPU "
1538 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1539 "GPU PME stream should be valid in order to use GPU version of PME.");
1541 const DeviceContext* deviceContext =
1542 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1543 const DeviceStream* pmeStream =
1544 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1546 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1547 nChargePerturbed != 0, nTypePerturbed != 0,
1548 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1549 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1550 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1552 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1557 if (EI_DYNAMICS(inputrec->eI))
1559 /* Turn on signal handling on all nodes */
1561 * (A user signal from the PME nodes (if any)
1562 * is communicated to the PP nodes.
1564 signal_handler_install();
1567 pull_t* pull_work = nullptr;
1568 if (thisRankHasDuty(cr, DUTY_PP))
1570 /* Assumes uniform use of the number of OpenMP threads */
1571 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1573 if (inputrec->bPull)
1575 /* Initialize pull code */
1576 pull_work = init_pull(fplog, inputrec->pull, inputrec.get(), &mtop, cr, &atomSets,
1577 inputrec->fepvals->init_lambda);
1578 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1580 initPullHistory(pull_work, &observablesHistory);
1582 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1584 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1588 std::unique_ptr<EnforcedRotation> enforcedRotation;
1591 /* Initialize enforced rotation code */
1592 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1593 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1597 t_swap* swap = nullptr;
1598 if (inputrec->eSwapCoords != eswapNO)
1600 /* Initialize ion swapping code */
1601 swap = init_swapcoords(fplog, inputrec.get(),
1602 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1603 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1604 oenv, mdrunOptions, startingBehavior);
1607 /* Let makeConstraints know whether we have essential dynamics constraints. */
1608 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1609 ms, &nrnb, wcycle, fr->bMolPBC);
1611 /* Energy terms and groups */
1612 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1613 inputrec->fepvals->n_lambda);
1615 /* Kinetic energy data */
1616 gmx_ekindata_t ekind;
1617 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1619 /* Set up interactive MD (IMD) */
1621 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1622 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1623 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1625 if (DOMAINDECOMP(cr))
1627 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1628 /* This call is not included in init_domain_decomposition mainly
1629 * because fr->cginfo_mb is set later.
1631 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1632 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1635 // TODO This is not the right place to manage the lifetime of
1636 // this data structure, but currently it's the easiest way to
1638 MdrunScheduleWorkload runScheduleWork;
1639 // Also populates the simulation constant workload description.
1640 runScheduleWork.simulationWork =
1641 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1642 useGpuForUpdate, devFlags.enableGpuBufferOps,
1643 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1645 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1646 if (gpusWereDetected
1647 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1648 || runScheduleWork.simulationWork.useGpuBufferOps))
1650 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1651 ? GpuApiCallBehavior::Async
1652 : GpuApiCallBehavior::Sync;
1653 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1654 "GPU device stream manager should be initialized to use GPU.");
1655 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1656 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1657 fr->stateGpu = stateGpu.get();
1660 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1661 SimulatorBuilder simulatorBuilder;
1663 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1664 simulatorBuilder.add(std::move(membedHolder));
1665 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1666 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1669 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1670 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1671 simulatorBuilder.add(ConstraintsParam(
1672 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1674 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1675 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1676 inputrec.get(), fr));
1677 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1678 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1679 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1680 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1681 // Todo move to an MDModule
1682 simulatorBuilder.add(IonSwapping(swap));
1683 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1684 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1685 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
1687 // build and run simulator object based on user-input
1688 auto simulator = simulatorBuilder.build(useModularSimulator);
1691 if (fr->pmePpCommGpu)
1693 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1694 fr->pmePpCommGpu.reset();
1697 if (inputrec->bPull)
1699 finish_pull(pull_work);
1701 finish_swapcoords(swap);
1705 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1707 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1708 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1709 deviceStreamManager.get());
1712 wallcycle_stop(wcycle, ewcRUN);
1714 /* Finish up, write some stuff
1715 * if rerunMD, don't write last frame again
1717 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1718 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1720 // clean up cycle counter
1721 wallcycle_destroy(wcycle);
1723 deviceStreamManager.reset(nullptr);
1727 gmx_pme_destroy(pmedata);
1731 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1732 // before we destroy the GPU context(s)
1733 // Pinned buffers are associated with contexts in CUDA.
1734 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1735 mdAtoms.reset(nullptr);
1736 globalState.reset(nullptr);
1737 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1738 gpuBonded.reset(nullptr);
1739 /* Free pinned buffers in *fr */
1742 // TODO convert to C++ so we can get rid of these frees
1746 if (!hwinfo->deviceInfoList.empty())
1748 /* stop the GPU profiler (only CUDA) */
1752 /* With tMPI we need to wait for all ranks to finish deallocation before
1753 * destroying the CUDA context as some tMPI ranks may be sharing
1756 * This is not a concern in OpenCL where we use one context per rank.
1758 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1759 * but it is easier and more futureproof to call it on the whole node.
1761 * Note that this function needs to be called even if GPUs are not used
1762 * in this run because the PME ranks have no knowledge of whether GPUs
1763 * are used or not, but all ranks need to enter the barrier below.
1764 * \todo Remove this physical node barrier after making sure
1765 * that it's not needed anymore (with a shared GPU run).
1769 physicalNodeComm.barrier();
1771 releaseDevice(deviceInfo);
1773 /* Does what it says */
1774 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1775 walltime_accounting_destroy(walltime_accounting);
1777 // Ensure log file content is written
1780 gmx_fio_flush(logFileHandle);
1783 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1784 * exceptions were enabled before function was called. */
1787 gmx_fedisableexcept();
1790 auto rc = static_cast<int>(gmx_get_stop_condition());
1793 /* we need to join all threads. The sub-threads join when they
1794 exit this function, but the master thread needs to be told to
1804 Mdrunner::~Mdrunner()
1806 // Clean up of the Manager.
1807 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1808 // but okay as long as threads synchronize some time before adding or accessing
1809 // a new set of restraints.
1810 if (restraintManager_)
1812 restraintManager_->clear();
1813 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1814 "restraints added during runner life time should be cleared at runner "
1819 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1821 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1822 // Not sure if this should be logged through the md logger or something else,
1823 // but it is helpful to have some sort of INFO level message sent somewhere.
1824 // std::cout << "Registering restraint named " << name << std::endl;
1826 // When multiple restraints are used, it may be wasteful to register them separately.
1827 // Maybe instead register an entire Restraint Manager as a force provider.
1828 restraintManager_->addToSpec(std::move(puller), name);
1831 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1833 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1835 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1837 class Mdrunner::BuilderImplementation
1840 BuilderImplementation() = delete;
1841 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1842 ~BuilderImplementation();
1844 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1845 real forceWarningThreshold,
1846 StartingBehavior startingBehavior);
1848 void addDomdec(const DomdecOptions& options);
1850 void addVerletList(int nstlist);
1852 void addReplicaExchange(const ReplicaExchangeParameters& params);
1854 void addNonBonded(const char* nbpu_opt);
1856 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1858 void addBondedTaskAssignment(const char* bonded_opt);
1860 void addUpdateTaskAssignment(const char* update_opt);
1862 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1864 void addFilenames(ArrayRef<const t_filenm> filenames);
1866 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1868 void addLogFile(t_fileio* logFileHandle);
1870 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1875 // Default parameters copied from runner.h
1876 // \todo Clarify source(s) of default parameters.
1878 const char* nbpu_opt_ = nullptr;
1879 const char* pme_opt_ = nullptr;
1880 const char* pme_fft_opt_ = nullptr;
1881 const char* bonded_opt_ = nullptr;
1882 const char* update_opt_ = nullptr;
1884 MdrunOptions mdrunOptions_;
1886 DomdecOptions domdecOptions_;
1888 ReplicaExchangeParameters replicaExchangeParameters_;
1890 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1893 //! Multisim communicator handle.
1894 gmx_multisim_t* multiSimulation_;
1896 //! mdrun communicator
1897 MPI_Comm communicator_ = MPI_COMM_NULL;
1899 //! Print a warning if any force is larger than this (in kJ/mol nm).
1900 real forceWarningThreshold_ = -1;
1902 //! Whether the simulation will start afresh, or restart with/without appending.
1903 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1905 //! The modules that comprise the functionality of mdrun.
1906 std::unique_ptr<MDModules> mdModules_;
1908 //! \brief Parallelism information.
1909 gmx_hw_opt_t hardwareOptions_;
1911 //! filename options for simulation.
1912 ArrayRef<const t_filenm> filenames_;
1914 /*! \brief Handle to output environment.
1916 * \todo gmx_output_env_t needs lifetime management.
1918 gmx_output_env_t* outputEnvironment_ = nullptr;
1920 /*! \brief Non-owning handle to MD log file.
1922 * \todo Context should own output facilities for client.
1923 * \todo Improve log file handle management.
1925 * Code managing the FILE* relies on the ability to set it to
1926 * nullptr to check whether the filehandle is valid.
1928 t_fileio* logFileHandle_ = nullptr;
1931 * \brief Builder for simulation stop signal handler.
1933 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1936 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1937 compat::not_null<SimulationContext*> context) :
1938 mdModules_(std::move(mdModules))
1940 communicator_ = context->communicator_;
1941 multiSimulation_ = context->multiSimulation_.get();
1944 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1946 Mdrunner::BuilderImplementation&
1947 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1948 const real forceWarningThreshold,
1949 const StartingBehavior startingBehavior)
1951 mdrunOptions_ = options;
1952 forceWarningThreshold_ = forceWarningThreshold;
1953 startingBehavior_ = startingBehavior;
1957 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1959 domdecOptions_ = options;
1962 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1967 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1969 replicaExchangeParameters_ = params;
1972 Mdrunner Mdrunner::BuilderImplementation::build()
1974 auto newRunner = Mdrunner(std::move(mdModules_));
1976 newRunner.mdrunOptions = mdrunOptions_;
1977 newRunner.pforce = forceWarningThreshold_;
1978 newRunner.startingBehavior = startingBehavior_;
1979 newRunner.domdecOptions = domdecOptions_;
1981 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1982 newRunner.hw_opt = hardwareOptions_;
1984 // No invariant to check. This parameter exists to optionally override other behavior.
1985 newRunner.nstlist_cmdline = nstlist_;
1987 newRunner.replExParams = replicaExchangeParameters_;
1989 newRunner.filenames = filenames_;
1991 newRunner.communicator = communicator_;
1993 // nullptr is a valid value for the multisim handle
1994 newRunner.ms = multiSimulation_;
1996 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1997 // \todo Update sanity checking when output environment has clearly specified invariants.
1998 // Initialization and default values for oenv are not well specified in the current version.
1999 if (outputEnvironment_)
2001 newRunner.oenv = outputEnvironment_;
2005 GMX_THROW(gmx::APIError(
2006 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2009 newRunner.logFileHandle = logFileHandle_;
2013 newRunner.nbpu_opt = nbpu_opt_;
2017 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2020 if (pme_opt_ && pme_fft_opt_)
2022 newRunner.pme_opt = pme_opt_;
2023 newRunner.pme_fft_opt = pme_fft_opt_;
2027 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2032 newRunner.bonded_opt = bonded_opt_;
2036 GMX_THROW(gmx::APIError(
2037 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2042 newRunner.update_opt = update_opt_;
2046 GMX_THROW(gmx::APIError(
2047 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2051 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2053 if (stopHandlerBuilder_)
2055 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2059 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2065 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2067 nbpu_opt_ = nbpu_opt;
2070 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2073 pme_fft_opt_ = pme_fft_opt;
2076 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2078 bonded_opt_ = bonded_opt;
2081 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2083 update_opt_ = update_opt;
2086 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2088 hardwareOptions_ = hardwareOptions;
2091 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2093 filenames_ = filenames;
2096 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2098 outputEnvironment_ = outputEnvironment;
2101 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2103 logFileHandle_ = logFileHandle;
2106 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2108 stopHandlerBuilder_ = std::move(builder);
2111 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2112 compat::not_null<SimulationContext*> context) :
2113 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2117 MdrunnerBuilder::~MdrunnerBuilder() = default;
2119 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2120 real forceWarningThreshold,
2121 const StartingBehavior startingBehavior)
2123 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2127 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2129 impl_->addDomdec(options);
2133 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2135 impl_->addVerletList(nstlist);
2139 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2141 impl_->addReplicaExchange(params);
2145 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2147 impl_->addNonBonded(nbpu_opt);
2151 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2153 // The builder method may become more general in the future, but in this version,
2154 // parameters for PME electrostatics are both required and the only parameters
2156 if (pme_opt && pme_fft_opt)
2158 impl_->addPME(pme_opt, pme_fft_opt);
2163 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2168 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2170 impl_->addBondedTaskAssignment(bonded_opt);
2174 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2176 impl_->addUpdateTaskAssignment(update_opt);
2180 Mdrunner MdrunnerBuilder::build()
2182 return impl_->build();
2185 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2187 impl_->addHardwareOptions(hardwareOptions);
2191 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2193 impl_->addFilenames(filenames);
2197 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2199 impl_->addOutputEnvironment(outputEnvironment);
2203 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2205 impl_->addLogFile(logFileHandle);
2209 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2211 impl_->addStopHandlerBuilder(std::move(builder));
2215 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2217 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;