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/commrec.h"
114 #include "gromacs/mdtypes/enerdata.h"
115 #include "gromacs/mdtypes/fcdata.h"
116 #include "gromacs/mdtypes/forcerec.h"
117 #include "gromacs/mdtypes/group.h"
118 #include "gromacs/mdtypes/inputrec.h"
119 #include "gromacs/mdtypes/interaction_const.h"
120 #include "gromacs/mdtypes/md_enums.h"
121 #include "gromacs/mdtypes/mdatom.h"
122 #include "gromacs/mdtypes/mdrunoptions.h"
123 #include "gromacs/mdtypes/observableshistory.h"
124 #include "gromacs/mdtypes/simulation_workload.h"
125 #include "gromacs/mdtypes/state.h"
126 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
127 #include "gromacs/modularsimulator/modularsimulator.h"
128 #include "gromacs/nbnxm/gpu_data_mgmt.h"
129 #include "gromacs/nbnxm/nbnxm.h"
130 #include "gromacs/nbnxm/pairlist_tuning.h"
131 #include "gromacs/pbcutil/pbc.h"
132 #include "gromacs/pulling/output.h"
133 #include "gromacs/pulling/pull.h"
134 #include "gromacs/pulling/pull_rotation.h"
135 #include "gromacs/restraint/manager.h"
136 #include "gromacs/restraint/restraintmdmodule.h"
137 #include "gromacs/restraint/restraintpotential.h"
138 #include "gromacs/swap/swapcoords.h"
139 #include "gromacs/taskassignment/decidegpuusage.h"
140 #include "gromacs/taskassignment/decidesimulationworkload.h"
141 #include "gromacs/taskassignment/resourcedivision.h"
142 #include "gromacs/taskassignment/taskassignment.h"
143 #include "gromacs/taskassignment/usergpuids.h"
144 #include "gromacs/timing/gpu_timing.h"
145 #include "gromacs/timing/wallcycle.h"
146 #include "gromacs/timing/wallcyclereporting.h"
147 #include "gromacs/topology/mtop_util.h"
148 #include "gromacs/trajectory/trajectoryframe.h"
149 #include "gromacs/utility/basenetwork.h"
150 #include "gromacs/utility/cstringutil.h"
151 #include "gromacs/utility/exceptions.h"
152 #include "gromacs/utility/fatalerror.h"
153 #include "gromacs/utility/filestream.h"
154 #include "gromacs/utility/gmxassert.h"
155 #include "gromacs/utility/gmxmpi.h"
156 #include "gromacs/utility/keyvaluetree.h"
157 #include "gromacs/utility/logger.h"
158 #include "gromacs/utility/loggerbuilder.h"
159 #include "gromacs/utility/mdmodulenotification.h"
160 #include "gromacs/utility/physicalnodecommunicator.h"
161 #include "gromacs/utility/pleasecite.h"
162 #include "gromacs/utility/programcontext.h"
163 #include "gromacs/utility/smalloc.h"
164 #include "gromacs/utility/stringutil.h"
166 #include "isimulator.h"
167 #include "membedholder.h"
168 #include "replicaexchange.h"
169 #include "simulatorbuilder.h"
172 # include "corewrap.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.enableGpuPmePPComm =
210 GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
212 #pragma GCC diagnostic pop
214 if (devFlags.enableGpuBufferOps)
216 GMX_LOG(mdlog.warning)
218 .appendTextFormatted(
219 "This run uses the 'GPU buffer ops' feature, enabled by the "
220 "GMX_USE_GPU_BUFFER_OPS environment variable.");
223 if (devFlags.forceGpuUpdateDefault)
225 GMX_LOG(mdlog.warning)
227 .appendTextFormatted(
228 "This run will default to '-update gpu' as requested by the "
229 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
230 "decomposition lacks substantial testing and should be used with caution.");
233 if (devFlags.enableGpuHaloExchange)
235 if (useGpuForNonbonded)
237 if (!devFlags.enableGpuBufferOps)
239 GMX_LOG(mdlog.warning)
241 .appendTextFormatted(
242 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
243 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
244 devFlags.enableGpuBufferOps = true;
246 GMX_LOG(mdlog.warning)
248 .appendTextFormatted(
249 "This run has requested the 'GPU halo exchange' feature, enabled by "
251 "GMX_GPU_DD_COMMS environment variable.");
255 GMX_LOG(mdlog.warning)
257 .appendTextFormatted(
258 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
259 "halo exchange' feature will not be enabled as nonbonded interactions "
260 "are not offloaded.");
261 devFlags.enableGpuHaloExchange = false;
265 if (devFlags.enableGpuPmePPComm)
267 if (pmeRunMode == PmeRunMode::GPU)
269 GMX_LOG(mdlog.warning)
271 .appendTextFormatted(
272 "This run uses the 'GPU PME-PP communications' feature, enabled "
273 "by the GMX_GPU_PME_PP_COMMS environment variable.");
277 std::string clarification;
278 if (pmeRunMode == PmeRunMode::Mixed)
281 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
286 clarification = "PME is not offloaded to the GPU.";
288 GMX_LOG(mdlog.warning)
291 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
292 "'GPU PME-PP communications' feature was not enabled as "
294 devFlags.enableGpuPmePPComm = false;
301 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
303 * Used to ensure that the master thread does not modify mdrunner during copy
304 * on the spawned threads. */
305 static void threadMpiMdrunnerAccessBarrier()
308 MPI_Barrier(MPI_COMM_WORLD);
312 Mdrunner Mdrunner::cloneOnSpawnedThread() const
314 auto newRunner = Mdrunner(std::make_unique<MDModules>());
316 // All runners in the same process share a restraint manager resource because it is
317 // part of the interface to the client code, which is associated only with the
318 // original thread. Handles to the same resources can be obtained by copy.
320 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
323 // Copy members of master runner.
324 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
325 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
326 newRunner.hw_opt = hw_opt;
327 newRunner.filenames = filenames;
329 newRunner.oenv = oenv;
330 newRunner.mdrunOptions = mdrunOptions;
331 newRunner.domdecOptions = domdecOptions;
332 newRunner.nbpu_opt = nbpu_opt;
333 newRunner.pme_opt = pme_opt;
334 newRunner.pme_fft_opt = pme_fft_opt;
335 newRunner.bonded_opt = bonded_opt;
336 newRunner.update_opt = update_opt;
337 newRunner.nstlist_cmdline = nstlist_cmdline;
338 newRunner.replExParams = replExParams;
339 newRunner.pforce = pforce;
340 // Give the spawned thread the newly created valid communicator
341 // for the simulation.
342 newRunner.communicator = MPI_COMM_WORLD;
344 newRunner.startingBehavior = startingBehavior;
345 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
347 threadMpiMdrunnerAccessBarrier();
352 /*! \brief The callback used for running on spawned threads.
354 * Obtains the pointer to the master mdrunner object from the one
355 * argument permitted to the thread-launch API call, copies it to make
356 * a new runner for this thread, reinitializes necessary data, and
357 * proceeds to the simulation. */
358 static void mdrunner_start_fn(const void* arg)
362 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
363 /* copy the arg list to make sure that it's thread-local. This
364 doesn't copy pointed-to items, of course; fnm, cr and fplog
365 are reset in the call below, all others should be const. */
366 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
369 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
373 void Mdrunner::spawnThreads(int numThreadsToLaunch)
376 /* now spawn new threads that start mdrunner_start_fn(), while
377 the main thread returns. Thread affinity is handled later. */
378 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
379 static_cast<const void*>(this))
382 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
385 // Give the master thread the newly created valid communicator for
387 communicator = MPI_COMM_WORLD;
388 threadMpiMdrunnerAccessBarrier();
390 GMX_UNUSED_VALUE(numThreadsToLaunch);
391 GMX_UNUSED_VALUE(mdrunner_start_fn);
397 /*! \brief Initialize variables for Verlet scheme simulation */
398 static void prepare_verlet_scheme(FILE* fplog,
402 const gmx_mtop_t* mtop,
404 bool makeGpuPairList,
405 const gmx::CpuInfo& cpuinfo)
407 // We checked the cut-offs in grompp, but double-check here.
408 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
409 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
411 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
412 "With Verlet lists and PME we should have rcoulomb>=rvdw");
416 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
417 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
419 /* For NVE simulations, we will retain the initial list buffer */
420 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
422 /* Update the Verlet buffer size for the current run setup */
424 /* Here we assume SIMD-enabled kernels are being used. But as currently
425 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
426 * and 4x2 gives a larger buffer than 4x4, this is ok.
428 ListSetupType listType =
429 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
430 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
432 const real rlist_new =
433 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
435 if (rlist_new != ir->rlist)
437 if (fplog != nullptr)
440 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
441 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
443 ir->rlist = rlist_new;
447 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
449 gmx_fatal(FARGS, "Can not set nstlist without %s",
450 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
453 if (EI_DYNAMICS(ir->eI))
455 /* Set or try nstlist values */
456 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
460 /*! \brief Override the nslist value in inputrec
462 * with value passed on the command line (if any)
464 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
468 /* override with anything else than the default -2 */
469 if (nsteps_cmdline > -2)
471 char sbuf_steps[STEPSTRSIZE];
472 char sbuf_msg[STRLEN];
474 ir->nsteps = nsteps_cmdline;
475 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
478 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
479 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
483 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
484 gmx_step_str(nsteps_cmdline, sbuf_steps));
487 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
489 else if (nsteps_cmdline < -2)
491 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
493 /* Do nothing if nsteps_cmdline == -2 */
499 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
501 * If not, and if a warning may be issued, logs a warning about
502 * falling back to CPU code. With thread-MPI, only the first
503 * call to this function should have \c issueWarning true. */
504 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
506 bool gpuIsUseful = true;
509 if (ir.opts.ngener - ir.nwall > 1)
511 /* The GPU code does not support more than one energy group.
512 * If the user requested GPUs explicitly, a fatal error is given later.
516 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
517 "For better performance, run on the GPU without energy groups and then do "
518 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
524 warning = "TPI is not implemented for GPUs.";
527 if (!gpuIsUseful && issueWarning)
529 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
535 //! Initializes the logger for mdrun.
536 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
538 gmx::LoggerBuilder builder;
539 if (fplog != nullptr)
541 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
543 if (isSimulationMasterRank)
545 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
547 return builder.build();
550 //! Make a TaskTarget from an mdrun argument string.
551 static TaskTarget findTaskTarget(const char* optionString)
553 TaskTarget returnValue = TaskTarget::Auto;
555 if (strncmp(optionString, "auto", 3) == 0)
557 returnValue = TaskTarget::Auto;
559 else if (strncmp(optionString, "cpu", 3) == 0)
561 returnValue = TaskTarget::Cpu;
563 else if (strncmp(optionString, "gpu", 3) == 0)
565 returnValue = TaskTarget::Gpu;
569 GMX_ASSERT(false, "Option string should have been checked for sanity already");
575 //! Finish run, aggregate data to print performance info.
576 static void finish_run(FILE* fplog,
577 const gmx::MDLogger& mdlog,
579 const t_inputrec* inputrec,
581 gmx_wallcycle_t wcycle,
582 gmx_walltime_accounting_t walltime_accounting,
583 nonbonded_verlet_t* nbv,
584 const gmx_pme_t* pme,
588 double nbfs = 0, mflop = 0;
589 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
590 elapsed_time_over_all_threads_over_all_ranks;
591 /* Control whether it is valid to print a report. Only the
592 simulation master may print, but it should not do so if the run
593 terminated e.g. before a scheduled reset step. This is
594 complicated by the fact that PME ranks are unaware of the
595 reason why they were sent a pmerecvqxFINISH. To avoid
596 communication deadlocks, we always do the communication for the
597 report, even if we've decided not to write the report, because
598 how long it takes to finish the run is not important when we've
599 decided not to report on the simulation performance.
601 Further, we only report performance for dynamical integrators,
602 because those are the only ones for which we plan to
603 consider doing any optimizations. */
604 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
606 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
608 GMX_LOG(mdlog.warning)
610 .appendText("Simulation ended prematurely, no performance report will be written.");
615 std::unique_ptr<t_nrnb> nrnbTotalStorage;
618 nrnbTotalStorage = std::make_unique<t_nrnb>();
619 nrnb_tot = nrnbTotalStorage.get();
621 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
629 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
630 elapsed_time_over_all_threads =
631 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
635 /* reduce elapsed_time over all MPI ranks in the current simulation */
636 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
638 elapsed_time_over_all_ranks /= cr->nnodes;
639 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
640 * current simulation. */
641 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
642 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
647 elapsed_time_over_all_ranks = elapsed_time;
648 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
653 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
656 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
658 print_dd_statistics(cr, inputrec, fplog);
661 /* TODO Move the responsibility for any scaling by thread counts
662 * to the code that handled the thread region, so that there's a
663 * mechanism to keep cycle counting working during the transition
664 * to task parallelism. */
665 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
666 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
667 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
668 nthreads_pp, nthreads_pme);
669 auto cycle_sum(wallcycle_sum(cr, wcycle));
673 auto nbnxn_gpu_timings =
674 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
675 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
677 if (pme_gpu_task_enabled(pme))
679 pme_gpu_get_timings(pme, &pme_gpu_timings);
681 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
682 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
685 if (EI_DYNAMICS(inputrec->eI))
687 delta_t = inputrec->delta_t;
692 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
693 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
694 delta_t, nbfs, mflop);
698 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
699 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
700 delta_t, nbfs, mflop);
705 int Mdrunner::mdrunner()
708 t_forcerec* fr = nullptr;
709 real ewaldcoeff_q = 0;
710 real ewaldcoeff_lj = 0;
711 int nChargePerturbed = -1, nTypePerturbed = 0;
712 gmx_wallcycle_t wcycle;
713 gmx_walltime_accounting_t walltime_accounting = nullptr;
714 MembedHolder membedHolder(filenames.size(), filenames.data());
715 gmx_hw_info_t* hwinfo = nullptr;
717 /* CAUTION: threads may be started later on in this function, so
718 cr doesn't reflect the final parallel state right now */
721 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
722 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
723 const bool doRerun = mdrunOptions.rerun;
725 // Handle task-assignment related user options.
726 EmulateGpuNonbonded emulateGpuNonbonded =
727 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
729 std::vector<int> userGpuTaskAssignment;
732 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
734 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
735 auto nonbondedTarget = findTaskTarget(nbpu_opt);
736 auto pmeTarget = findTaskTarget(pme_opt);
737 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
738 auto bondedTarget = findTaskTarget(bonded_opt);
739 auto updateTarget = findTaskTarget(update_opt);
741 FILE* fplog = nullptr;
742 // If we are appending, we don't write log output because we need
743 // to check that the old log file matches what the checkpoint file
744 // expects. Otherwise, we should start to write log output now if
745 // there is a file ready for it.
746 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
748 fplog = gmx_fio_getfp(logFileHandle);
750 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
751 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
752 gmx::MDLogger mdlog(logOwner.logger());
754 // TODO The thread-MPI master rank makes a working
755 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
756 // after the threads have been launched. This works because no use
757 // is made of that communicator until after the execution paths
758 // have rejoined. But it is likely that we can improve the way
759 // this is expressed, e.g. by expressly running detection only the
760 // master rank for thread-MPI, rather than relying on the mutex
761 // and reference count.
762 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
763 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
765 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
767 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
769 // Print citation requests after all software/hardware printing
770 pleaseCiteGromacs(fplog);
772 // Note: legacy program logic relies on checking whether these pointers are assigned.
773 // Objects may or may not be allocated later.
774 std::unique_ptr<t_inputrec> inputrec;
775 std::unique_ptr<t_state> globalState;
777 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
779 if (isSimulationMasterRank)
781 // Allocate objects to be initialized by later function calls.
782 /* Only the master rank has the global state */
783 globalState = std::make_unique<t_state>();
784 inputrec = std::make_unique<t_inputrec>();
786 /* Read (nearly) all data required for the simulation
787 * and keep the partly serialized tpr contents to send to other ranks later
789 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
790 inputrec.get(), globalState.get(), &mtop);
793 /* Check and update the hardware options for internal consistency */
794 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
797 if (GMX_THREAD_MPI && isSimulationMasterRank)
799 bool useGpuForNonbonded = false;
800 bool useGpuForPme = false;
803 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
805 // If the user specified the number of ranks, then we must
806 // respect that, but in default mode, we need to allow for
807 // the number of GPUs to choose the number of ranks.
808 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
809 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
810 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
811 canUseGpuForNonbonded,
812 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
813 hw_opt.nthreads_tmpi);
814 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
815 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
816 *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
818 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
820 /* Determine how many thread-MPI ranks to start.
822 * TODO Over-writing the user-supplied value here does
823 * prevent any possible subsequent checks from working
825 hw_opt.nthreads_tmpi =
826 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
827 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
829 // Now start the threads for thread MPI.
830 spawnThreads(hw_opt.nthreads_tmpi);
831 // The spawned threads enter mdrunner() and execution of
832 // master and spawned threads joins at the end of this block.
833 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
836 GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
837 "Must have valid world communicator unless running a multi-simulation");
838 CommrecHandle crHandle = init_commrec(communicator);
839 t_commrec* cr = crHandle.get();
840 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
844 /* now broadcast everything to the non-master nodes/threads: */
845 if (!isSimulationMasterRank)
847 // Until now, only the master rank has a non-null pointer.
848 // On non-master ranks, allocate the object that will receive data in the following call.
849 inputrec = std::make_unique<t_inputrec>();
851 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec.get(), &mtop,
852 partialDeserializedTpr.get());
854 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
855 partialDeserializedTpr.reset(nullptr);
857 // Now the number of ranks is known to all ranks, and each knows
858 // the inputrec read by the master rank. The ranks can now all run
859 // the task-deciding functions and will agree on the result
860 // without needing to communicate.
861 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
863 // Note that these variables describe only their own node.
865 // Note that when bonded interactions run on a GPU they always run
866 // alongside a nonbonded task, so do not influence task assignment
867 // even though they affect the force calculation workload.
868 bool useGpuForNonbonded = false;
869 bool useGpuForPme = false;
870 bool useGpuForBonded = false;
871 bool useGpuForUpdate = false;
872 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
875 // It's possible that there are different numbers of GPUs on
876 // different nodes, which is the user's responsibility to
877 // handle. If unsuitable, we will notice that during task
879 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
880 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
881 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
882 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
883 useGpuForPme = decideWhetherToUseGpusForPme(
884 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec,
885 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
886 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
887 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
888 useGpuForBonded = decideWhetherToUseGpusForBonded(
889 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
890 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
891 domdecOptions.numPmeRanks, gpusWereDetected);
893 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
895 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
897 // Initialize development feature flags that enabled by environment variable
898 // and report those features that are enabled.
899 const DevelopmentFeatureFlags devFlags =
900 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
902 const bool useModularSimulator =
903 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
904 nullptr, doEssentialDynamics, membedHolder.doMembed());
907 // TODO: hide restraint implementation details from Mdrunner.
908 // There is nothing unique about restraints at this point as far as the
909 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
910 // factory functions from the SimulationContext on which to call mdModules_->add().
911 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
912 for (auto&& restraint : restraintManager_->getRestraints())
914 auto module = RestraintMDModule::create(restraint, restraint->sites());
915 mdModules_->add(std::move(module));
918 // TODO: Error handling
919 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
920 // now that the MdModules know their options, they know which callbacks to sign up to
921 mdModules_->subscribeToSimulationSetupNotifications();
922 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
924 if (inputrec->internalParameters != nullptr)
926 mdModulesNotifier.notify(*inputrec->internalParameters);
929 if (fplog != nullptr)
931 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
932 fprintf(fplog, "\n");
937 /* In rerun, set velocities to zero if present */
938 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
940 // rerun does not use velocities
944 "Rerun trajectory contains velocities. Rerun does only evaluate "
945 "potential energy and forces. The velocities will be ignored.");
946 for (int i = 0; i < globalState->natoms; i++)
948 clear_rvec(globalState->v[i]);
950 globalState->flags &= ~(1 << estV);
953 /* now make sure the state is initialized and propagated */
954 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
957 /* NM and TPI parallelize over force/energy calculations, not atoms,
958 * so we need to initialize and broadcast the global state.
960 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
964 globalState = std::make_unique<t_state>();
966 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
970 /* A parallel command line option consistency check that we can
971 only do after any threads have started. */
973 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
974 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
977 "The -dd or -npme option request a parallel simulation, "
979 "but %s was compiled without threads or MPI enabled",
980 output_env_get_program_display_name(oenv));
982 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
984 "but %s was not started through mpirun/mpiexec or only one rank was requested "
985 "through mpirun/mpiexec",
986 output_env_get_program_display_name(oenv));
990 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
993 "The .mdp file specified an energy mininization or normal mode algorithm, and "
994 "these are not compatible with mdrun -rerun");
997 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
999 if (domdecOptions.numPmeRanks > 0)
1001 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1002 "PME-only ranks are requested, but the system does not use PME "
1003 "for electrostatics or LJ");
1006 domdecOptions.numPmeRanks = 0;
1009 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1011 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1012 * improve performance with many threads per GPU, since our OpenMP
1013 * scaling is bad, but it's difficult to automate the setup.
1015 domdecOptions.numPmeRanks = 0;
1019 if (domdecOptions.numPmeRanks < 0)
1021 domdecOptions.numPmeRanks = 0;
1022 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1026 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1027 "PME GPU decomposition is not supported");
1034 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
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);
1062 ObservablesHistory observablesHistory = {};
1064 if (startingBehavior != StartingBehavior::NewSimulation)
1066 /* Check if checkpoint file exists before doing continuation.
1067 * This way we can use identical input options for the first and subsequent runs...
1069 if (mdrunOptions.numStepsCommandline > -2)
1071 /* Temporarily set the number of steps to unlimited to avoid
1072 * triggering the nsteps check in load_checkpoint().
1073 * This hack will go away soon when the -nsteps option is removed.
1075 inputrec->nsteps = -1;
1078 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1079 logFileHandle, cr, domdecOptions.numCells, inputrec.get(), globalState.get(),
1080 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1082 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1084 // Now we can start normal logging to the truncated log file.
1085 fplog = gmx_fio_getfp(logFileHandle);
1086 prepareLogAppending(fplog);
1087 logOwner = buildLogger(fplog, MASTER(cr));
1088 mdlog = logOwner.logger();
1092 if (mdrunOptions.numStepsCommandline > -2)
1097 "The -nsteps functionality is deprecated, and may be removed in a future "
1099 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1102 /* override nsteps with value set on the commandline */
1103 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1105 if (isSimulationMasterRank)
1107 copy_mat(globalState->box, box);
1112 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1115 if (inputrec->cutoff_scheme != ecutsVERLET)
1118 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1119 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1121 /* Update rlist and nstlist. */
1122 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1123 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1124 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1126 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1127 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1130 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1131 // This builder is necessary while we have multi-part construction
1132 // of DD. Before DD is constructed, we use the existence of
1133 // the builder object to indicate that further construction of DD
1135 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1136 if (useDomainDecomposition)
1138 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1139 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1140 positionsFromStatePointer(globalState.get()));
1144 /* PME, if used, is done on all nodes with 1D decomposition */
1145 cr->nnodes = cr->sizeOfDefaultCommunicator;
1146 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1147 cr->nodeid = cr->rankInDefaultCommunicator;
1149 cr->duty = (DUTY_PP | DUTY_PME);
1151 if (inputrec->pbcType == PbcType::Screw)
1153 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1157 // Produce the task assignment for this rank - done after DD is constructed
1158 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1159 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1160 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1161 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1162 // TODO cr->duty & DUTY_PME should imply that a PME
1163 // algorithm is active, but currently does not.
1164 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1166 // Get the device handles for the modules, nullptr when no task is assigned.
1168 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1170 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1171 bool useTiming = true;
1175 /* WARNING: CUDA timings are incorrect with multiple streams.
1176 * This is the main reason why they are disabled by default.
1178 // TODO: Consider turning on by default when we can detect nr of streams.
1179 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1181 else if (GMX_GPU_OPENCL)
1183 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1186 // TODO Currently this is always built, yet DD partition code
1187 // checks if it is built before using it. Probably it should
1188 // become an MDModule that is made only when another module
1189 // requires it (e.g. pull, CompEl, density fitting), so that we
1190 // don't update the local atom sets unilaterally every step.
1191 LocalAtomSetManager atomSets;
1194 // TODO Pass the GPU streams to ddBuilder to use in buffer
1195 // transfers (e.g. halo exchange)
1196 cr->dd = ddBuilder->build(&atomSets);
1197 // The builder's job is done, so destruct it
1198 ddBuilder.reset(nullptr);
1199 // Note that local state still does not exist yet.
1202 // The GPU update is decided here because we need to know whether the constraints or
1203 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1204 // defined). This is only known after DD is initialized, hence decision on using GPU
1205 // update is done so late.
1208 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1210 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1211 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1212 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1213 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1214 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1216 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1218 const bool printHostName = (cr->nnodes > 1);
1219 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1221 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1223 if (deviceInfo != nullptr)
1225 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1227 dd_setup_dlb_resource_sharing(cr, deviceId);
1229 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1230 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1231 useGpuForUpdate, useTiming);
1234 // If the user chose a task assignment, give them some hints
1235 // where appropriate.
1236 if (!userGpuTaskAssignment.empty())
1238 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1243 /* After possible communicator splitting in make_dd_communicators.
1244 * we can set up the intra/inter node communication.
1246 gmx_setup_nodecomm(fplog, cr);
1252 GMX_LOG(mdlog.warning)
1254 .appendTextFormatted(
1255 "This is simulation %d out of %d running as a composite GROMACS\n"
1256 "multi-simulation job. Setup for this simulation:\n",
1257 ms->simulationIndex_, ms->numSimulations_);
1259 GMX_LOG(mdlog.warning)
1260 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1262 cr->nnodes == 1 ? "thread" : "threads"
1264 cr->nnodes == 1 ? "process" : "processes"
1270 // If mdrun -pin auto honors any affinity setting that already
1271 // exists. If so, it is nice to provide feedback about whether
1272 // that existing affinity setting was from OpenMP or something
1273 // else, so we run this code both before and after we initialize
1274 // the OpenMP support.
1275 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1276 /* Check and update the number of OpenMP threads requested */
1277 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1278 pmeRunMode, mtop, *inputrec);
1280 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1281 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1283 // Enable FP exception detection, but not in
1284 // Release mode and not for compilers with known buggy FP
1285 // exception support (clang with any optimization) or suspected
1286 // buggy FP exception support (gcc 7.* with optimization).
1287 #if !defined NDEBUG \
1288 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1289 && defined __OPTIMIZE__)
1290 const bool bEnableFPE = true;
1292 const bool bEnableFPE = false;
1294 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1297 gmx_feenableexcept();
1300 /* Now that we know the setup is consistent, check for efficiency */
1301 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1302 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1304 /* getting number of PP/PME threads on this MPI / tMPI rank.
1305 PME: env variable should be read only on one node to make sure it is
1306 identical everywhere;
1308 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1309 : gmx_omp_nthreads_get(emntPME);
1310 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1311 physicalNodeComm, mdlog);
1313 // Enable Peer access between GPUs where available
1314 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1315 // any of the GPU communication features are active.
1316 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1317 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1319 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1322 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1324 /* Before setting affinity, check whether the affinity has changed
1325 * - which indicates that probably the OpenMP library has changed it
1326 * since we first checked).
1328 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1330 int numThreadsOnThisNode, intraNodeThreadOffset;
1331 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1332 &intraNodeThreadOffset);
1334 /* Set the CPU affinity */
1335 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1336 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1339 if (mdrunOptions.timingOptions.resetStep > -1)
1344 "The -resetstep functionality is deprecated, and may be removed in a "
1347 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1351 /* Master synchronizes its value of reset_counters with all nodes
1352 * including PME only nodes */
1353 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1354 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1355 wcycle_set_reset_counters(wcycle, reset_counters);
1358 // Membrane embedding must be initialized before we call init_forcerec()
1359 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1360 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1362 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1363 std::unique_ptr<MDAtoms> mdAtoms;
1364 std::unique_ptr<VirtualSitesHandler> vsite;
1365 std::unique_ptr<GpuBonded> gpuBonded;
1368 if (thisRankHasDuty(cr, DUTY_PP))
1370 mdModulesNotifier.notify(*cr);
1371 mdModulesNotifier.notify(&atomSets);
1372 mdModulesNotifier.notify(inputrec->pbcType);
1373 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1374 /* Initiate forcerecord */
1375 fr = new t_forcerec;
1376 fr->forceProviders = mdModules_->initForceProviders();
1377 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1378 opt2fn("-table", filenames.size(), filenames.data()),
1379 opt2fn("-tablep", filenames.size(), filenames.data()),
1380 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1381 // Dirty hack, for fixing disres and orires should be made mdmodules
1382 fr->listedForces->fcdata().disres = disresdata;
1383 fr->listedForces->fcdata().orires = oriresdata;
1385 // Save a handle to device stream manager to use elsewhere in the code
1386 // TODO: Forcerec is not a correct place to store it.
1387 fr->deviceStreamManager = deviceStreamManager.get();
1389 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1392 deviceStreamManager != nullptr,
1393 "GPU device stream manager should be valid in order to use PME-PP direct "
1396 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1397 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1399 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1400 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1401 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1404 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo, useGpuForNonbonded,
1405 deviceStreamManager.get(), &mtop, box, wcycle);
1406 // TODO: Move the logic below to a GPU bonded builder
1407 if (useGpuForBonded)
1409 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1410 "GPU device stream manager should be valid in order to use GPU "
1411 "version of bonded forces.");
1412 gpuBonded = std::make_unique<GpuBonded>(
1413 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1414 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1415 fr->gpuBonded = gpuBonded.get();
1418 /* Initialize the mdAtoms structure.
1419 * mdAtoms is not filled with atom data,
1420 * as this can not be done now with domain decomposition.
1422 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1423 if (globalState && thisRankHasPmeGpuTask)
1425 // The pinning of coordinates in the global state object works, because we only use
1426 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1427 // points to the global state object without DD.
1428 // FIXME: MD and EM separately set up the local state - this should happen in the same
1429 // function, which should also perform the pinning.
1430 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1433 /* Initialize the virtual site communication */
1434 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1436 calc_shifts(box, fr->shift_vec);
1438 /* With periodic molecules the charge groups should be whole at start up
1439 * and the virtual sites should not be far from their proper positions.
1441 if (!inputrec->bContinuation && MASTER(cr)
1442 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1444 /* Make molecules whole at start of run */
1445 if (fr->pbcType != PbcType::No)
1447 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1451 /* Correct initial vsite positions are required
1452 * for the initial distribution in the domain decomposition
1453 * and for the initial shell prediction.
1455 constructVirtualSitesGlobal(mtop, globalState->x);
1459 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1461 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1462 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1467 /* This is a PME only node */
1469 GMX_ASSERT(globalState == nullptr,
1470 "We don't need the state on a PME only rank and expect it to be unitialized");
1472 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1473 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1476 gmx_pme_t* sepPmeData = nullptr;
1477 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1478 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1479 "Double-checking that only PME-only ranks have no forcerec");
1480 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1482 // TODO should live in ewald module once its testing is improved
1484 // Later, this program could contain kernels that might be later
1485 // re-used as auto-tuning progresses, or subsequent simulations
1487 PmeGpuProgramStorage pmeGpuProgram;
1488 if (thisRankHasPmeGpuTask)
1491 (deviceStreamManager != nullptr),
1492 "GPU device stream manager should be initialized in order to use GPU for PME.");
1493 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1494 "GPU device should be initialized in order to use GPU for PME.");
1495 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1498 /* Initiate PME if necessary,
1499 * either on all nodes or on dedicated PME nodes only. */
1500 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1502 if (mdAtoms && mdAtoms->mdatoms())
1504 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1505 if (EVDW_PME(inputrec->vdwtype))
1507 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1510 if (cr->npmenodes > 0)
1512 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1513 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1514 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1517 if (thisRankHasDuty(cr, DUTY_PME))
1521 // TODO: This should be in the builder.
1522 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1523 "Device stream manager should be valid in order to use GPU "
1526 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1527 "GPU PME stream should be valid in order to use GPU version of PME.");
1529 const DeviceContext* deviceContext =
1530 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1531 const DeviceStream* pmeStream =
1532 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1534 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1535 nChargePerturbed != 0, nTypePerturbed != 0,
1536 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1537 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1538 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1540 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1545 if (EI_DYNAMICS(inputrec->eI))
1547 /* Turn on signal handling on all nodes */
1549 * (A user signal from the PME nodes (if any)
1550 * is communicated to the PP nodes.
1552 signal_handler_install();
1555 pull_t* pull_work = nullptr;
1556 if (thisRankHasDuty(cr, DUTY_PP))
1558 /* Assumes uniform use of the number of OpenMP threads */
1559 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1561 if (inputrec->bPull)
1563 /* Initialize pull code */
1564 pull_work = init_pull(fplog, inputrec->pull, inputrec.get(), &mtop, cr, &atomSets,
1565 inputrec->fepvals->init_lambda);
1566 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1568 initPullHistory(pull_work, &observablesHistory);
1570 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1572 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1576 std::unique_ptr<EnforcedRotation> enforcedRotation;
1579 /* Initialize enforced rotation code */
1580 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1581 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1585 t_swap* swap = nullptr;
1586 if (inputrec->eSwapCoords != eswapNO)
1588 /* Initialize ion swapping code */
1589 swap = init_swapcoords(fplog, inputrec.get(),
1590 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1591 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1592 oenv, mdrunOptions, startingBehavior);
1595 /* Let makeConstraints know whether we have essential dynamics constraints. */
1596 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1597 ms, &nrnb, wcycle, fr->bMolPBC);
1599 /* Energy terms and groups */
1600 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1601 inputrec->fepvals->n_lambda);
1603 /* Kinetic energy data */
1604 gmx_ekindata_t ekind;
1605 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1607 /* Set up interactive MD (IMD) */
1609 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1610 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1611 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1613 if (DOMAINDECOMP(cr))
1615 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1616 /* This call is not included in init_domain_decomposition mainly
1617 * because fr->cginfo_mb is set later.
1619 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1620 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1623 // TODO This is not the right place to manage the lifetime of
1624 // this data structure, but currently it's the easiest way to
1626 MdrunScheduleWorkload runScheduleWork;
1627 // Also populates the simulation constant workload description.
1628 runScheduleWork.simulationWork =
1629 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1630 useGpuForUpdate, devFlags.enableGpuBufferOps,
1631 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1633 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1634 if (gpusWereDetected
1635 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1636 || runScheduleWork.simulationWork.useGpuBufferOps))
1638 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1639 ? GpuApiCallBehavior::Async
1640 : GpuApiCallBehavior::Sync;
1641 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1642 "GPU device stream manager should be initialized to use GPU.");
1643 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1644 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1645 fr->stateGpu = stateGpu.get();
1648 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1649 SimulatorBuilder simulatorBuilder;
1651 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1652 simulatorBuilder.add(std::move(membedHolder));
1653 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1654 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1657 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1658 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1659 simulatorBuilder.add(ConstraintsParam(
1660 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1662 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1663 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1664 inputrec.get(), fr));
1665 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1666 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1667 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1668 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1669 // Todo move to an MDModule
1670 simulatorBuilder.add(IonSwapping(swap));
1671 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1672 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1674 // build and run simulator object based on user-input
1675 auto simulator = simulatorBuilder.build(useModularSimulator);
1678 if (fr->pmePpCommGpu)
1680 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1681 fr->pmePpCommGpu.reset();
1684 if (inputrec->bPull)
1686 finish_pull(pull_work);
1688 finish_swapcoords(swap);
1692 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1694 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1695 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1696 deviceStreamManager.get());
1699 wallcycle_stop(wcycle, ewcRUN);
1701 /* Finish up, write some stuff
1702 * if rerunMD, don't write last frame again
1704 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1705 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1707 // clean up cycle counter
1708 wallcycle_destroy(wcycle);
1710 deviceStreamManager.reset(nullptr);
1714 gmx_pme_destroy(pmedata);
1718 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1719 // before we destroy the GPU context(s) in free_gpu().
1720 // Pinned buffers are associated with contexts in CUDA.
1721 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1722 mdAtoms.reset(nullptr);
1723 globalState.reset(nullptr);
1724 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1725 gpuBonded.reset(nullptr);
1726 /* Free pinned buffers in *fr */
1729 // TODO convert to C++ so we can get rid of these frees
1733 if (hwinfo->gpu_info.n_dev > 0)
1735 /* stop the GPU profiler (only CUDA) */
1739 /* With tMPI we need to wait for all ranks to finish deallocation before
1740 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1743 * This is not a concern in OpenCL where we use one context per rank which
1744 * is freed in nbnxn_gpu_free().
1746 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1747 * but it is easier and more futureproof to call it on the whole node.
1749 * Note that this function needs to be called even if GPUs are not used
1750 * in this run because the PME ranks have no knowledge of whether GPUs
1751 * are used or not, but all ranks need to enter the barrier below.
1752 * \todo Remove this physical node barrier after making sure
1753 * that it's not needed anymore (with a shared GPU run).
1757 physicalNodeComm.barrier();
1760 free_gpu(deviceInfo);
1762 /* Does what it says */
1763 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1764 walltime_accounting_destroy(walltime_accounting);
1766 // Ensure log file content is written
1769 gmx_fio_flush(logFileHandle);
1772 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1773 * exceptions were enabled before function was called. */
1776 gmx_fedisableexcept();
1779 auto rc = static_cast<int>(gmx_get_stop_condition());
1782 /* we need to join all threads. The sub-threads join when they
1783 exit this function, but the master thread needs to be told to
1793 Mdrunner::~Mdrunner()
1795 // Clean up of the Manager.
1796 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1797 // but okay as long as threads synchronize some time before adding or accessing
1798 // a new set of restraints.
1799 if (restraintManager_)
1801 restraintManager_->clear();
1802 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1803 "restraints added during runner life time should be cleared at runner "
1808 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1810 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1811 // Not sure if this should be logged through the md logger or something else,
1812 // but it is helpful to have some sort of INFO level message sent somewhere.
1813 // std::cout << "Registering restraint named " << name << std::endl;
1815 // When multiple restraints are used, it may be wasteful to register them separately.
1816 // Maybe instead register an entire Restraint Manager as a force provider.
1817 restraintManager_->addToSpec(std::move(puller), name);
1820 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1822 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1824 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1826 class Mdrunner::BuilderImplementation
1829 BuilderImplementation() = delete;
1830 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1831 ~BuilderImplementation();
1833 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1834 real forceWarningThreshold,
1835 StartingBehavior startingBehavior);
1837 void addDomdec(const DomdecOptions& options);
1839 void addVerletList(int nstlist);
1841 void addReplicaExchange(const ReplicaExchangeParameters& params);
1843 void addNonBonded(const char* nbpu_opt);
1845 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1847 void addBondedTaskAssignment(const char* bonded_opt);
1849 void addUpdateTaskAssignment(const char* update_opt);
1851 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1853 void addFilenames(ArrayRef<const t_filenm> filenames);
1855 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1857 void addLogFile(t_fileio* logFileHandle);
1859 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1864 // Default parameters copied from runner.h
1865 // \todo Clarify source(s) of default parameters.
1867 const char* nbpu_opt_ = nullptr;
1868 const char* pme_opt_ = nullptr;
1869 const char* pme_fft_opt_ = nullptr;
1870 const char* bonded_opt_ = nullptr;
1871 const char* update_opt_ = nullptr;
1873 MdrunOptions mdrunOptions_;
1875 DomdecOptions domdecOptions_;
1877 ReplicaExchangeParameters replicaExchangeParameters_;
1879 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1882 //! Multisim communicator handle.
1883 gmx_multisim_t* multiSimulation_;
1885 //! mdrun communicator
1886 MPI_Comm communicator_ = MPI_COMM_NULL;
1888 //! Print a warning if any force is larger than this (in kJ/mol nm).
1889 real forceWarningThreshold_ = -1;
1891 //! Whether the simulation will start afresh, or restart with/without appending.
1892 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1894 //! The modules that comprise the functionality of mdrun.
1895 std::unique_ptr<MDModules> mdModules_;
1897 //! \brief Parallelism information.
1898 gmx_hw_opt_t hardwareOptions_;
1900 //! filename options for simulation.
1901 ArrayRef<const t_filenm> filenames_;
1903 /*! \brief Handle to output environment.
1905 * \todo gmx_output_env_t needs lifetime management.
1907 gmx_output_env_t* outputEnvironment_ = nullptr;
1909 /*! \brief Non-owning handle to MD log file.
1911 * \todo Context should own output facilities for client.
1912 * \todo Improve log file handle management.
1914 * Code managing the FILE* relies on the ability to set it to
1915 * nullptr to check whether the filehandle is valid.
1917 t_fileio* logFileHandle_ = nullptr;
1920 * \brief Builder for simulation stop signal handler.
1922 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1925 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1926 compat::not_null<SimulationContext*> context) :
1927 mdModules_(std::move(mdModules))
1929 communicator_ = context->communicator_;
1930 multiSimulation_ = context->multiSimulation_.get();
1933 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1935 Mdrunner::BuilderImplementation&
1936 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1937 const real forceWarningThreshold,
1938 const StartingBehavior startingBehavior)
1940 mdrunOptions_ = options;
1941 forceWarningThreshold_ = forceWarningThreshold;
1942 startingBehavior_ = startingBehavior;
1946 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1948 domdecOptions_ = options;
1951 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1956 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1958 replicaExchangeParameters_ = params;
1961 Mdrunner Mdrunner::BuilderImplementation::build()
1963 auto newRunner = Mdrunner(std::move(mdModules_));
1965 newRunner.mdrunOptions = mdrunOptions_;
1966 newRunner.pforce = forceWarningThreshold_;
1967 newRunner.startingBehavior = startingBehavior_;
1968 newRunner.domdecOptions = domdecOptions_;
1970 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1971 newRunner.hw_opt = hardwareOptions_;
1973 // No invariant to check. This parameter exists to optionally override other behavior.
1974 newRunner.nstlist_cmdline = nstlist_;
1976 newRunner.replExParams = replicaExchangeParameters_;
1978 newRunner.filenames = filenames_;
1980 newRunner.communicator = communicator_;
1982 // nullptr is a valid value for the multisim handle
1983 newRunner.ms = multiSimulation_;
1985 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1986 // \todo Update sanity checking when output environment has clearly specified invariants.
1987 // Initialization and default values for oenv are not well specified in the current version.
1988 if (outputEnvironment_)
1990 newRunner.oenv = outputEnvironment_;
1994 GMX_THROW(gmx::APIError(
1995 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1998 newRunner.logFileHandle = logFileHandle_;
2002 newRunner.nbpu_opt = nbpu_opt_;
2006 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2009 if (pme_opt_ && pme_fft_opt_)
2011 newRunner.pme_opt = pme_opt_;
2012 newRunner.pme_fft_opt = pme_fft_opt_;
2016 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2021 newRunner.bonded_opt = bonded_opt_;
2025 GMX_THROW(gmx::APIError(
2026 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2031 newRunner.update_opt = update_opt_;
2035 GMX_THROW(gmx::APIError(
2036 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2040 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2042 if (stopHandlerBuilder_)
2044 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2048 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2054 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2056 nbpu_opt_ = nbpu_opt;
2059 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2062 pme_fft_opt_ = pme_fft_opt;
2065 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2067 bonded_opt_ = bonded_opt;
2070 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2072 update_opt_ = update_opt;
2075 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2077 hardwareOptions_ = hardwareOptions;
2080 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2082 filenames_ = filenames;
2085 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2087 outputEnvironment_ = outputEnvironment;
2090 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2092 logFileHandle_ = logFileHandle;
2095 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2097 stopHandlerBuilder_ = std::move(builder);
2100 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2101 compat::not_null<SimulationContext*> context) :
2102 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2106 MdrunnerBuilder::~MdrunnerBuilder() = default;
2108 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2109 real forceWarningThreshold,
2110 const StartingBehavior startingBehavior)
2112 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2116 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2118 impl_->addDomdec(options);
2122 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2124 impl_->addVerletList(nstlist);
2128 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2130 impl_->addReplicaExchange(params);
2134 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2136 impl_->addNonBonded(nbpu_opt);
2140 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2142 // The builder method may become more general in the future, but in this version,
2143 // parameters for PME electrostatics are both required and the only parameters
2145 if (pme_opt && pme_fft_opt)
2147 impl_->addPME(pme_opt, pme_fft_opt);
2152 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2157 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2159 impl_->addBondedTaskAssignment(bonded_opt);
2163 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2165 impl_->addUpdateTaskAssignment(update_opt);
2169 Mdrunner MdrunnerBuilder::build()
2171 return impl_->build();
2174 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2176 impl_->addHardwareOptions(hardwareOptions);
2180 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2182 impl_->addFilenames(filenames);
2186 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2188 impl_->addOutputEnvironment(outputEnvironment);
2192 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2194 impl_->addLogFile(logFileHandle);
2198 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2200 impl_->addStopHandlerBuilder(std::move(builder));
2204 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2206 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;