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/gpu_utils.h"
77 #include "gromacs/hardware/cpuinfo.h"
78 #include "gromacs/hardware/detecthardware.h"
79 #include "gromacs/hardware/printhardware.h"
80 #include "gromacs/imd/imd.h"
81 #include "gromacs/listed_forces/disre.h"
82 #include "gromacs/listed_forces/gpubonded.h"
83 #include "gromacs/listed_forces/orires.h"
84 #include "gromacs/math/functions.h"
85 #include "gromacs/math/utilities.h"
86 #include "gromacs/math/vec.h"
87 #include "gromacs/mdlib/boxdeformation.h"
88 #include "gromacs/mdlib/broadcaststructs.h"
89 #include "gromacs/mdlib/calc_verletbuf.h"
90 #include "gromacs/mdlib/dispersioncorrection.h"
91 #include "gromacs/mdlib/enerdata_utils.h"
92 #include "gromacs/mdlib/force.h"
93 #include "gromacs/mdlib/forcerec.h"
94 #include "gromacs/mdlib/gmx_omp_nthreads.h"
95 #include "gromacs/mdlib/makeconstraints.h"
96 #include "gromacs/mdlib/md_support.h"
97 #include "gromacs/mdlib/mdatoms.h"
98 #include "gromacs/mdlib/membed.h"
99 #include "gromacs/mdlib/qmmm.h"
100 #include "gromacs/mdlib/sighandler.h"
101 #include "gromacs/mdlib/stophandler.h"
102 #include "gromacs/mdlib/updategroups.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdrun/mdmodules.h"
105 #include "gromacs/mdrun/simulationcontext.h"
106 #include "gromacs/mdrunutility/handlerestart.h"
107 #include "gromacs/mdrunutility/logging.h"
108 #include "gromacs/mdrunutility/multisim.h"
109 #include "gromacs/mdrunutility/printtime.h"
110 #include "gromacs/mdrunutility/threadaffinity.h"
111 #include "gromacs/mdtypes/commrec.h"
112 #include "gromacs/mdtypes/enerdata.h"
113 #include "gromacs/mdtypes/fcdata.h"
114 #include "gromacs/mdtypes/forcerec.h"
115 #include "gromacs/mdtypes/group.h"
116 #include "gromacs/mdtypes/inputrec.h"
117 #include "gromacs/mdtypes/interaction_const.h"
118 #include "gromacs/mdtypes/md_enums.h"
119 #include "gromacs/mdtypes/mdatom.h"
120 #include "gromacs/mdtypes/mdrunoptions.h"
121 #include "gromacs/mdtypes/observableshistory.h"
122 #include "gromacs/mdtypes/simulation_workload.h"
123 #include "gromacs/mdtypes/state.h"
124 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
125 #include "gromacs/nbnxm/gpu_data_mgmt.h"
126 #include "gromacs/nbnxm/nbnxm.h"
127 #include "gromacs/nbnxm/pairlist_tuning.h"
128 #include "gromacs/pbcutil/pbc.h"
129 #include "gromacs/pulling/output.h"
130 #include "gromacs/pulling/pull.h"
131 #include "gromacs/pulling/pull_rotation.h"
132 #include "gromacs/restraint/manager.h"
133 #include "gromacs/restraint/restraintmdmodule.h"
134 #include "gromacs/restraint/restraintpotential.h"
135 #include "gromacs/swap/swapcoords.h"
136 #include "gromacs/taskassignment/decidegpuusage.h"
137 #include "gromacs/taskassignment/decidesimulationworkload.h"
138 #include "gromacs/taskassignment/resourcedivision.h"
139 #include "gromacs/taskassignment/taskassignment.h"
140 #include "gromacs/taskassignment/usergpuids.h"
141 #include "gromacs/timing/gpu_timing.h"
142 #include "gromacs/timing/wallcycle.h"
143 #include "gromacs/timing/wallcyclereporting.h"
144 #include "gromacs/topology/mtop_util.h"
145 #include "gromacs/trajectory/trajectoryframe.h"
146 #include "gromacs/utility/basenetwork.h"
147 #include "gromacs/utility/cstringutil.h"
148 #include "gromacs/utility/exceptions.h"
149 #include "gromacs/utility/fatalerror.h"
150 #include "gromacs/utility/filestream.h"
151 #include "gromacs/utility/gmxassert.h"
152 #include "gromacs/utility/gmxmpi.h"
153 #include "gromacs/utility/keyvaluetree.h"
154 #include "gromacs/utility/logger.h"
155 #include "gromacs/utility/loggerbuilder.h"
156 #include "gromacs/utility/mdmodulenotification.h"
157 #include "gromacs/utility/physicalnodecommunicator.h"
158 #include "gromacs/utility/pleasecite.h"
159 #include "gromacs/utility/programcontext.h"
160 #include "gromacs/utility/smalloc.h"
161 #include "gromacs/utility/stringutil.h"
163 #include "isimulator.h"
164 #include "replicaexchange.h"
165 #include "simulatorbuilder.h"
168 # include "corewrap.h"
175 /*! \brief Manage any development feature flag variables encountered
177 * The use of dev features indicated by environment variables is
178 * logged in order to ensure that runs with such features enabled can
179 * be identified from their log and standard output. Any cross
180 * dependencies are also checked, and if unsatisfied, a fatal error
183 * Note that some development features overrides are applied already here:
184 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
186 * \param[in] mdlog Logger object.
187 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
188 * \param[in] pmeRunMode The PME run mode for this run
189 * \returns The object populated with development feature flags.
191 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
192 const bool useGpuForNonbonded,
193 const PmeRunMode pmeRunMode)
195 DevelopmentFeatureFlags devFlags;
197 // Some builds of GCC 5 give false positive warnings that these
198 // getenv results are ignored when clearly they are used.
199 #pragma GCC diagnostic push
200 #pragma GCC diagnostic ignored "-Wunused-result"
201 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
202 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
203 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
204 devFlags.enableGpuHaloExchange =
205 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
206 devFlags.enableGpuPmePPComm =
207 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
208 #pragma GCC diagnostic pop
210 if (devFlags.enableGpuBufferOps)
212 GMX_LOG(mdlog.warning)
214 .appendTextFormatted(
215 "This run uses the 'GPU buffer ops' feature, enabled by the "
216 "GMX_USE_GPU_BUFFER_OPS environment variable.");
219 if (devFlags.forceGpuUpdateDefault)
221 GMX_LOG(mdlog.warning)
223 .appendTextFormatted(
224 "This run will default to '-update gpu' as requested by the "
225 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
226 "decomposition lacks substantial testing and should be used with caution.");
229 if (devFlags.enableGpuHaloExchange)
231 if (useGpuForNonbonded)
233 if (!devFlags.enableGpuBufferOps)
235 GMX_LOG(mdlog.warning)
237 .appendTextFormatted(
238 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
239 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
240 devFlags.enableGpuBufferOps = true;
242 GMX_LOG(mdlog.warning)
244 .appendTextFormatted(
245 "This run has requested the 'GPU halo exchange' feature, enabled by "
247 "GMX_GPU_DD_COMMS environment variable.");
251 GMX_LOG(mdlog.warning)
253 .appendTextFormatted(
254 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
255 "halo exchange' feature will not be enabled as nonbonded interactions "
256 "are not offloaded.");
257 devFlags.enableGpuHaloExchange = false;
261 if (devFlags.enableGpuPmePPComm)
263 if (pmeRunMode == PmeRunMode::GPU)
265 GMX_LOG(mdlog.warning)
267 .appendTextFormatted(
268 "This run uses the 'GPU PME-PP communications' feature, enabled "
269 "by the GMX_GPU_PME_PP_COMMS environment variable.");
273 std::string clarification;
274 if (pmeRunMode == PmeRunMode::Mixed)
277 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
282 clarification = "PME is not offloaded to the GPU.";
284 GMX_LOG(mdlog.warning)
287 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
288 "'GPU PME-PP communications' feature was not enabled as "
290 devFlags.enableGpuPmePPComm = false;
297 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
299 * Used to ensure that the master thread does not modify mdrunner during copy
300 * on the spawned threads. */
301 static void threadMpiMdrunnerAccessBarrier()
304 MPI_Barrier(MPI_COMM_WORLD);
308 Mdrunner Mdrunner::cloneOnSpawnedThread() const
310 auto newRunner = Mdrunner(std::make_unique<MDModules>());
312 // All runners in the same process share a restraint manager resource because it is
313 // part of the interface to the client code, which is associated only with the
314 // original thread. Handles to the same resources can be obtained by copy.
316 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
319 // Copy members of master runner.
320 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
321 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
322 newRunner.hw_opt = hw_opt;
323 newRunner.filenames = filenames;
325 newRunner.oenv = oenv;
326 newRunner.mdrunOptions = mdrunOptions;
327 newRunner.domdecOptions = domdecOptions;
328 newRunner.nbpu_opt = nbpu_opt;
329 newRunner.pme_opt = pme_opt;
330 newRunner.pme_fft_opt = pme_fft_opt;
331 newRunner.bonded_opt = bonded_opt;
332 newRunner.update_opt = update_opt;
333 newRunner.nstlist_cmdline = nstlist_cmdline;
334 newRunner.replExParams = replExParams;
335 newRunner.pforce = pforce;
336 // Give the spawned thread the newly created valid communicator
337 // for the simulation.
338 newRunner.communicator = MPI_COMM_WORLD;
340 newRunner.startingBehavior = startingBehavior;
341 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
343 threadMpiMdrunnerAccessBarrier();
348 /*! \brief The callback used for running on spawned threads.
350 * Obtains the pointer to the master mdrunner object from the one
351 * argument permitted to the thread-launch API call, copies it to make
352 * a new runner for this thread, reinitializes necessary data, and
353 * proceeds to the simulation. */
354 static void mdrunner_start_fn(const void* arg)
358 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
359 /* copy the arg list to make sure that it's thread-local. This
360 doesn't copy pointed-to items, of course; fnm, cr and fplog
361 are reset in the call below, all others should be const. */
362 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
365 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
369 void Mdrunner::spawnThreads(int numThreadsToLaunch)
372 /* now spawn new threads that start mdrunner_start_fn(), while
373 the main thread returns. Thread affinity is handled later. */
374 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
375 static_cast<const void*>(this))
378 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
381 // Give the master thread the newly created valid communicator for
383 communicator = MPI_COMM_WORLD;
384 threadMpiMdrunnerAccessBarrier();
386 GMX_UNUSED_VALUE(numThreadsToLaunch);
387 GMX_UNUSED_VALUE(mdrunner_start_fn);
393 /*! \brief Initialize variables for Verlet scheme simulation */
394 static void prepare_verlet_scheme(FILE* fplog,
398 const gmx_mtop_t* mtop,
400 bool makeGpuPairList,
401 const gmx::CpuInfo& cpuinfo)
403 // We checked the cut-offs in grompp, but double-check here.
404 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
405 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
407 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
408 "With Verlet lists and PME we should have rcoulomb>=rvdw");
412 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
413 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
415 /* For NVE simulations, we will retain the initial list buffer */
416 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
418 /* Update the Verlet buffer size for the current run setup */
420 /* Here we assume SIMD-enabled kernels are being used. But as currently
421 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
422 * and 4x2 gives a larger buffer than 4x4, this is ok.
424 ListSetupType listType =
425 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
426 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
428 const real rlist_new =
429 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
431 if (rlist_new != ir->rlist)
433 if (fplog != nullptr)
436 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
437 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
439 ir->rlist = rlist_new;
443 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
445 gmx_fatal(FARGS, "Can not set nstlist without %s",
446 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
449 if (EI_DYNAMICS(ir->eI))
451 /* Set or try nstlist values */
452 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
456 /*! \brief Override the nslist value in inputrec
458 * with value passed on the command line (if any)
460 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
464 /* override with anything else than the default -2 */
465 if (nsteps_cmdline > -2)
467 char sbuf_steps[STEPSTRSIZE];
468 char sbuf_msg[STRLEN];
470 ir->nsteps = nsteps_cmdline;
471 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
474 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
475 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
479 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
480 gmx_step_str(nsteps_cmdline, sbuf_steps));
483 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
485 else if (nsteps_cmdline < -2)
487 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
489 /* Do nothing if nsteps_cmdline == -2 */
495 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
497 * If not, and if a warning may be issued, logs a warning about
498 * falling back to CPU code. With thread-MPI, only the first
499 * call to this function should have \c issueWarning true. */
500 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
502 bool gpuIsUseful = true;
505 if (ir.opts.ngener - ir.nwall > 1)
507 /* The GPU code does not support more than one energy group.
508 * If the user requested GPUs explicitly, a fatal error is given later.
512 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
513 "For better performance, run on the GPU without energy groups and then do "
514 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
520 warning = "TPI is not implemented for GPUs.";
523 if (!gpuIsUseful && issueWarning)
525 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
531 //! Initializes the logger for mdrun.
532 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
534 gmx::LoggerBuilder builder;
535 if (fplog != nullptr)
537 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
539 if (isSimulationMasterRank)
541 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
543 return builder.build();
546 //! Make a TaskTarget from an mdrun argument string.
547 static TaskTarget findTaskTarget(const char* optionString)
549 TaskTarget returnValue = TaskTarget::Auto;
551 if (strncmp(optionString, "auto", 3) == 0)
553 returnValue = TaskTarget::Auto;
555 else if (strncmp(optionString, "cpu", 3) == 0)
557 returnValue = TaskTarget::Cpu;
559 else if (strncmp(optionString, "gpu", 3) == 0)
561 returnValue = TaskTarget::Gpu;
565 GMX_ASSERT(false, "Option string should have been checked for sanity already");
571 //! Finish run, aggregate data to print performance info.
572 static void finish_run(FILE* fplog,
573 const gmx::MDLogger& mdlog,
575 const t_inputrec* inputrec,
577 gmx_wallcycle_t wcycle,
578 gmx_walltime_accounting_t walltime_accounting,
579 nonbonded_verlet_t* nbv,
580 const gmx_pme_t* pme,
584 double nbfs = 0, mflop = 0;
585 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
586 elapsed_time_over_all_threads_over_all_ranks;
587 /* Control whether it is valid to print a report. Only the
588 simulation master may print, but it should not do so if the run
589 terminated e.g. before a scheduled reset step. This is
590 complicated by the fact that PME ranks are unaware of the
591 reason why they were sent a pmerecvqxFINISH. To avoid
592 communication deadlocks, we always do the communication for the
593 report, even if we've decided not to write the report, because
594 how long it takes to finish the run is not important when we've
595 decided not to report on the simulation performance.
597 Further, we only report performance for dynamical integrators,
598 because those are the only ones for which we plan to
599 consider doing any optimizations. */
600 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
602 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
604 GMX_LOG(mdlog.warning)
606 .appendText("Simulation ended prematurely, no performance report will be written.");
611 std::unique_ptr<t_nrnb> nrnbTotalStorage;
614 nrnbTotalStorage = std::make_unique<t_nrnb>();
615 nrnb_tot = nrnbTotalStorage.get();
617 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
625 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
626 elapsed_time_over_all_threads =
627 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
631 /* reduce elapsed_time over all MPI ranks in the current simulation */
632 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
634 elapsed_time_over_all_ranks /= cr->nnodes;
635 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
636 * current simulation. */
637 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
638 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
643 elapsed_time_over_all_ranks = elapsed_time;
644 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
649 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
652 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
654 print_dd_statistics(cr, inputrec, fplog);
657 /* TODO Move the responsibility for any scaling by thread counts
658 * to the code that handled the thread region, so that there's a
659 * mechanism to keep cycle counting working during the transition
660 * to task parallelism. */
661 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
662 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
663 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
664 nthreads_pp, nthreads_pme);
665 auto cycle_sum(wallcycle_sum(cr, wcycle));
669 auto nbnxn_gpu_timings =
670 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
671 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
673 if (pme_gpu_task_enabled(pme))
675 pme_gpu_get_timings(pme, &pme_gpu_timings);
677 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
678 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
681 if (EI_DYNAMICS(inputrec->eI))
683 delta_t = inputrec->delta_t;
688 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
689 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
690 delta_t, nbfs, mflop);
694 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
695 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
696 delta_t, nbfs, mflop);
701 int Mdrunner::mdrunner()
704 t_forcerec* fr = nullptr;
705 t_fcdata* fcd = nullptr;
706 real ewaldcoeff_q = 0;
707 real ewaldcoeff_lj = 0;
708 int nChargePerturbed = -1, nTypePerturbed = 0;
709 gmx_wallcycle_t wcycle;
710 gmx_walltime_accounting_t walltime_accounting = nullptr;
711 gmx_membed_t* membed = nullptr;
712 gmx_hw_info_t* hwinfo = nullptr;
714 /* CAUTION: threads may be started later on in this function, so
715 cr doesn't reflect the final parallel state right now */
718 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
719 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
720 const bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
721 const bool doRerun = mdrunOptions.rerun;
723 // Handle task-assignment related user options.
724 EmulateGpuNonbonded emulateGpuNonbonded =
725 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
727 std::vector<int> userGpuTaskAssignment;
730 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
732 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
733 auto nonbondedTarget = findTaskTarget(nbpu_opt);
734 auto pmeTarget = findTaskTarget(pme_opt);
735 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
736 auto bondedTarget = findTaskTarget(bonded_opt);
737 auto updateTarget = findTaskTarget(update_opt);
739 FILE* fplog = nullptr;
740 // If we are appending, we don't write log output because we need
741 // to check that the old log file matches what the checkpoint file
742 // expects. Otherwise, we should start to write log output now if
743 // there is a file ready for it.
744 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
746 fplog = gmx_fio_getfp(logFileHandle);
748 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
749 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
750 gmx::MDLogger mdlog(logOwner.logger());
752 // TODO The thread-MPI master rank makes a working
753 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
754 // after the threads have been launched. This works because no use
755 // is made of that communicator until after the execution paths
756 // have rejoined. But it is likely that we can improve the way
757 // this is expressed, e.g. by expressly running detection only the
758 // master rank for thread-MPI, rather than relying on the mutex
759 // and reference count.
760 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
761 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
763 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
765 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
767 // Print citation requests after all software/hardware printing
768 pleaseCiteGromacs(fplog);
770 // TODO Replace this by unique_ptr once t_inputrec is C++
771 t_inputrec inputrecInstance;
772 t_inputrec* inputrec = nullptr;
773 std::unique_ptr<t_state> globalState;
775 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
777 if (isSimulationMasterRank)
779 /* Only the master rank has the global state */
780 globalState = std::make_unique<t_state>();
782 /* Read (nearly) all data required for the simulation
783 * and keep the partly serialized tpr contents to send to other ranks later
785 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
786 &inputrecInstance, globalState.get(), &mtop);
787 inputrec = &inputrecInstance;
790 /* Check and update the hardware options for internal consistency */
791 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
794 if (GMX_THREAD_MPI && isSimulationMasterRank)
796 bool useGpuForNonbonded = false;
797 bool useGpuForPme = false;
800 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
802 // If the user specified the number of ranks, then we must
803 // respect that, but in default mode, we need to allow for
804 // the number of GPUs to choose the number of ranks.
805 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
806 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
807 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
808 canUseGpuForNonbonded,
809 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
810 hw_opt.nthreads_tmpi);
811 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
812 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
813 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
815 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
817 /* Determine how many thread-MPI ranks to start.
819 * TODO Over-writing the user-supplied value here does
820 * prevent any possible subsequent checks from working
822 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
823 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
825 // Now start the threads for thread MPI.
826 spawnThreads(hw_opt.nthreads_tmpi);
827 // The spawned threads enter mdrunner() and execution of
828 // master and spawned threads joins at the end of this block.
829 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
832 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
833 CommrecHandle crHandle = init_commrec(communicator, ms);
834 t_commrec* cr = crHandle.get();
835 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
839 /* now broadcast everything to the non-master nodes/threads: */
840 if (!isSimulationMasterRank)
842 inputrec = &inputrecInstance;
844 init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
846 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
847 partialDeserializedTpr.reset(nullptr);
849 // Now the number of ranks is known to all ranks, and each knows
850 // the inputrec read by the master rank. The ranks can now all run
851 // the task-deciding functions and will agree on the result
852 // without needing to communicate.
853 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
855 // Note that these variables describe only their own node.
857 // Note that when bonded interactions run on a GPU they always run
858 // alongside a nonbonded task, so do not influence task assignment
859 // even though they affect the force calculation workload.
860 bool useGpuForNonbonded = false;
861 bool useGpuForPme = false;
862 bool useGpuForBonded = false;
863 bool useGpuForUpdate = false;
864 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
867 // It's possible that there are different numbers of GPUs on
868 // different nodes, which is the user's responsibility to
869 // handle. If unsuitable, we will notice that during task
871 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
872 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
873 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
874 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
875 useGpuForPme = decideWhetherToUseGpusForPme(
876 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
877 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
878 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
879 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
880 useGpuForBonded = decideWhetherToUseGpusForBonded(
881 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
882 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
883 domdecOptions.numPmeRanks, gpusWereDetected);
885 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
887 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
889 // Initialize development feature flags that enabled by environment variable
890 // and report those features that are enabled.
891 const DevelopmentFeatureFlags devFlags =
892 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
894 const bool useModularSimulator = checkUseModularSimulator(
895 false, inputrec, doRerun, mtop, ms, replExParams, nullptr, doEssentialDynamics, doMembed);
898 // TODO: hide restraint implementation details from Mdrunner.
899 // There is nothing unique about restraints at this point as far as the
900 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
901 // factory functions from the SimulationContext on which to call mdModules_->add().
902 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
903 for (auto&& restraint : restraintManager_->getRestraints())
905 auto module = RestraintMDModule::create(restraint, restraint->sites());
906 mdModules_->add(std::move(module));
909 // TODO: Error handling
910 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
911 // now that the MdModules know their options, they know which callbacks to sign up to
912 mdModules_->subscribeToSimulationSetupNotifications();
913 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
915 if (inputrec->internalParameters != nullptr)
917 mdModulesNotifier.notify(*inputrec->internalParameters);
920 if (fplog != nullptr)
922 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
923 fprintf(fplog, "\n");
928 /* In rerun, set velocities to zero if present */
929 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
931 // rerun does not use velocities
935 "Rerun trajectory contains velocities. Rerun does only evaluate "
936 "potential energy and forces. The velocities will be ignored.");
937 for (int i = 0; i < globalState->natoms; i++)
939 clear_rvec(globalState->v[i]);
941 globalState->flags &= ~(1 << estV);
944 /* now make sure the state is initialized and propagated */
945 set_state_entries(globalState.get(), inputrec, useModularSimulator);
948 /* NM and TPI parallelize over force/energy calculations, not atoms,
949 * so we need to initialize and broadcast the global state.
951 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
955 globalState = std::make_unique<t_state>();
957 broadcastStateWithoutDynamics(cr, globalState.get());
960 /* A parallel command line option consistency check that we can
961 only do after any threads have started. */
963 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
964 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
967 "The -dd or -npme option request a parallel simulation, "
969 "but %s was compiled without threads or MPI enabled",
970 output_env_get_program_display_name(oenv));
972 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
974 "but %s was not started through mpirun/mpiexec or only one rank was requested "
975 "through mpirun/mpiexec",
976 output_env_get_program_display_name(oenv));
980 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
983 "The .mdp file specified an energy mininization or normal mode algorithm, and "
984 "these are not compatible with mdrun -rerun");
987 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
989 if (domdecOptions.numPmeRanks > 0)
991 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
992 "PME-only ranks are requested, but the system does not use PME "
993 "for electrostatics or LJ");
996 domdecOptions.numPmeRanks = 0;
999 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1001 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1002 * improve performance with many threads per GPU, since our OpenMP
1003 * scaling is bad, but it's difficult to automate the setup.
1005 domdecOptions.numPmeRanks = 0;
1009 if (domdecOptions.numPmeRanks < 0)
1011 domdecOptions.numPmeRanks = 0;
1012 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1016 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1017 "PME GPU decomposition is not supported");
1024 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1028 /* NMR restraints must be initialized before load_checkpoint,
1029 * since with time averaging the history is added to t_state.
1030 * For proper consistency check we therefore need to extend
1032 * So the PME-only nodes (if present) will also initialize
1033 * the distance restraints.
1037 /* This needs to be called before read_checkpoint to extend the state */
1038 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1040 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1042 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1044 ObservablesHistory observablesHistory = {};
1046 if (startingBehavior != StartingBehavior::NewSimulation)
1048 /* Check if checkpoint file exists before doing continuation.
1049 * This way we can use identical input options for the first and subsequent runs...
1051 if (mdrunOptions.numStepsCommandline > -2)
1053 /* Temporarily set the number of steps to unlimited to avoid
1054 * triggering the nsteps check in load_checkpoint().
1055 * This hack will go away soon when the -nsteps option is removed.
1057 inputrec->nsteps = -1;
1060 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1061 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1062 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1064 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1066 // Now we can start normal logging to the truncated log file.
1067 fplog = gmx_fio_getfp(logFileHandle);
1068 prepareLogAppending(fplog);
1069 logOwner = buildLogger(fplog, MASTER(cr));
1070 mdlog = logOwner.logger();
1074 if (mdrunOptions.numStepsCommandline > -2)
1079 "The -nsteps functionality is deprecated, and may be removed in a future "
1081 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1084 /* override nsteps with value set on the commandline */
1085 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1089 copy_mat(globalState->box, box);
1094 gmx_bcast(sizeof(box), box, cr);
1097 if (inputrec->cutoff_scheme != ecutsVERLET)
1100 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1101 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1103 /* Update rlist and nstlist. */
1104 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1105 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1108 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1109 // This builder is necessary while we have multi-part construction
1110 // of DD. Before DD is constructed, we use the existence of
1111 // the builder object to indicate that further construction of DD
1113 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1114 if (useDomainDecomposition)
1116 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1117 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1118 positionsFromStatePointer(globalState.get()));
1122 /* PME, if used, is done on all nodes with 1D decomposition */
1124 cr->duty = (DUTY_PP | DUTY_PME);
1126 if (inputrec->pbcType == PbcType::Screw)
1128 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1132 // Produce the task assignment for this rank.
1133 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1134 GpuTaskAssignments gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1135 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1136 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1137 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1138 // TODO cr->duty & DUTY_PME should imply that a PME
1139 // algorithm is active, but currently does not.
1140 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1142 // Get the device handles for the modules, nullptr when no task is assigned.
1143 DeviceInformation* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1144 DeviceInformation* pmeDeviceInfo = gpuTaskAssignments.initPmeDevice();
1146 // TODO Initialize GPU streams here.
1148 // TODO Currently this is always built, yet DD partition code
1149 // checks if it is built before using it. Probably it should
1150 // become an MDModule that is made only when another module
1151 // requires it (e.g. pull, CompEl, density fitting), so that we
1152 // don't update the local atom sets unilaterally every step.
1153 LocalAtomSetManager atomSets;
1156 // TODO Pass the GPU streams to ddBuilder to use in buffer
1157 // transfers (e.g. halo exchange)
1158 cr->dd = ddBuilder->build(&atomSets);
1159 // The builder's job is done, so destruct it
1160 ddBuilder.reset(nullptr);
1161 // Note that local state still does not exist yet.
1164 // The GPU update is decided here because we need to know whether the constraints or
1165 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1166 // defined). This is only known after DD is initialized, hence decision on using GPU
1167 // update is done so late.
1170 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1172 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1173 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1174 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1175 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1176 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1178 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1180 const bool printHostName = (cr->nnodes > 1);
1181 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1183 // If the user chose a task assignment, give them some hints
1184 // where appropriate.
1185 if (!userGpuTaskAssignment.empty())
1187 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1192 /* After possible communicator splitting in make_dd_communicators.
1193 * we can set up the intra/inter node communication.
1195 gmx_setup_nodecomm(fplog, cr);
1201 GMX_LOG(mdlog.warning)
1203 .appendTextFormatted(
1204 "This is simulation %d out of %d running as a composite GROMACS\n"
1205 "multi-simulation job. Setup for this simulation:\n",
1208 GMX_LOG(mdlog.warning)
1209 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1211 cr->nnodes == 1 ? "thread" : "threads"
1213 cr->nnodes == 1 ? "process" : "processes"
1219 // If mdrun -pin auto honors any affinity setting that already
1220 // exists. If so, it is nice to provide feedback about whether
1221 // that existing affinity setting was from OpenMP or something
1222 // else, so we run this code both before and after we initialize
1223 // the OpenMP support.
1224 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1225 /* Check and update the number of OpenMP threads requested */
1226 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1227 pmeRunMode, mtop, *inputrec);
1229 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1230 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1232 // Enable FP exception detection, but not in
1233 // Release mode and not for compilers with known buggy FP
1234 // exception support (clang with any optimization) or suspected
1235 // buggy FP exception support (gcc 7.* with optimization).
1236 #if !defined NDEBUG \
1237 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1238 && defined __OPTIMIZE__)
1239 const bool bEnableFPE = true;
1241 const bool bEnableFPE = false;
1243 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1246 gmx_feenableexcept();
1249 /* Now that we know the setup is consistent, check for efficiency */
1250 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1251 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1253 /* getting number of PP/PME threads on this MPI / tMPI rank.
1254 PME: env variable should be read only on one node to make sure it is
1255 identical everywhere;
1257 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1258 : gmx_omp_nthreads_get(emntPME);
1259 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1260 physicalNodeComm, mdlog);
1262 // Enable Peer access between GPUs where available
1263 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1264 // any of the GPU communication features are active.
1265 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1266 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1268 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1271 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1273 /* Before setting affinity, check whether the affinity has changed
1274 * - which indicates that probably the OpenMP library has changed it
1275 * since we first checked).
1277 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1279 int numThreadsOnThisNode, intraNodeThreadOffset;
1280 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1281 &intraNodeThreadOffset);
1283 /* Set the CPU affinity */
1284 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1285 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1288 if (mdrunOptions.timingOptions.resetStep > -1)
1293 "The -resetstep functionality is deprecated, and may be removed in a "
1296 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1300 /* Master synchronizes its value of reset_counters with all nodes
1301 * including PME only nodes */
1302 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1303 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1304 wcycle_set_reset_counters(wcycle, reset_counters);
1307 // Membrane embedding must be initialized before we call init_forcerec()
1312 fprintf(stderr, "Initializing membed");
1314 /* Note that membed cannot work in parallel because mtop is
1315 * changed here. Fix this if we ever want to make it run with
1316 * multiple ranks. */
1317 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1318 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1321 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1322 std::unique_ptr<MDAtoms> mdAtoms;
1323 std::unique_ptr<gmx_vsite_t> vsite;
1324 std::unique_ptr<GpuBonded> gpuBonded;
1327 if (thisRankHasDuty(cr, DUTY_PP))
1329 mdModulesNotifier.notify(*cr);
1330 mdModulesNotifier.notify(&atomSets);
1331 mdModulesNotifier.notify(inputrec->pbcType);
1332 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1333 /* Initiate forcerecord */
1334 fr = new t_forcerec;
1335 fr->forceProviders = mdModules_->initForceProviders();
1336 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1337 opt2fn("-table", filenames.size(), filenames.data()),
1338 opt2fn("-tablep", filenames.size(), filenames.data()),
1339 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1341 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1343 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(cr->mpi_comm_mysim, cr->dd->pme_nodeid);
1346 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, nonbondedDeviceInfo,
1347 &mtop, box, wcycle);
1348 if (useGpuForBonded)
1350 auto stream = havePPDomainDecomposition(cr)
1351 ? Nbnxm::gpu_get_command_stream(
1352 fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1353 : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1354 gmx::InteractionLocality::Local);
1355 gpuBonded = std::make_unique<GpuBonded>(mtop.ffparams, stream, wcycle);
1356 fr->gpuBonded = gpuBonded.get();
1359 /* Initialize the mdAtoms structure.
1360 * mdAtoms is not filled with atom data,
1361 * as this can not be done now with domain decomposition.
1363 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1364 if (globalState && thisRankHasPmeGpuTask)
1366 // The pinning of coordinates in the global state object works, because we only use
1367 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1368 // points to the global state object without DD.
1369 // FIXME: MD and EM separately set up the local state - this should happen in the same
1370 // function, which should also perform the pinning.
1371 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1374 /* Initialize the virtual site communication */
1375 vsite = initVsite(mtop, cr);
1377 calc_shifts(box, fr->shift_vec);
1379 /* With periodic molecules the charge groups should be whole at start up
1380 * and the virtual sites should not be far from their proper positions.
1382 if (!inputrec->bContinuation && MASTER(cr)
1383 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1385 /* Make molecules whole at start of run */
1386 if (fr->pbcType != PbcType::No)
1388 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1392 /* Correct initial vsite positions are required
1393 * for the initial distribution in the domain decomposition
1394 * and for the initial shell prediction.
1396 constructVsitesGlobal(mtop, globalState->x);
1400 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1402 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1403 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1408 /* This is a PME only node */
1410 GMX_ASSERT(globalState == nullptr,
1411 "We don't need the state on a PME only rank and expect it to be unitialized");
1413 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1414 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1417 gmx_pme_t* sepPmeData = nullptr;
1418 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1419 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1420 "Double-checking that only PME-only ranks have no forcerec");
1421 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1423 // TODO should live in ewald module once its testing is improved
1425 // Later, this program could contain kernels that might be later
1426 // re-used as auto-tuning progresses, or subsequent simulations
1428 PmeGpuProgramStorage pmeGpuProgram;
1429 if (thisRankHasPmeGpuTask)
1431 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1434 /* Initiate PME if necessary,
1435 * either on all nodes or on dedicated PME nodes only. */
1436 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1438 if (mdAtoms && mdAtoms->mdatoms())
1440 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1441 if (EVDW_PME(inputrec->vdwtype))
1443 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1446 if (cr->npmenodes > 0)
1448 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1449 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1450 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1453 if (thisRankHasDuty(cr, DUTY_PME))
1457 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1458 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1459 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1460 nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1462 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1467 if (EI_DYNAMICS(inputrec->eI))
1469 /* Turn on signal handling on all nodes */
1471 * (A user signal from the PME nodes (if any)
1472 * is communicated to the PP nodes.
1474 signal_handler_install();
1477 pull_t* pull_work = nullptr;
1478 if (thisRankHasDuty(cr, DUTY_PP))
1480 /* Assumes uniform use of the number of OpenMP threads */
1481 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1483 if (inputrec->bPull)
1485 /* Initialize pull code */
1486 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1487 inputrec->fepvals->init_lambda);
1488 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1490 initPullHistory(pull_work, &observablesHistory);
1492 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1494 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1498 std::unique_ptr<EnforcedRotation> enforcedRotation;
1501 /* Initialize enforced rotation code */
1503 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1504 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1507 t_swap* swap = nullptr;
1508 if (inputrec->eSwapCoords != eswapNO)
1510 /* Initialize ion swapping code */
1511 swap = init_swapcoords(fplog, inputrec,
1512 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1513 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1514 oenv, mdrunOptions, startingBehavior);
1517 /* Let makeConstraints know whether we have essential dynamics constraints. */
1518 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1519 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1521 /* Energy terms and groups */
1522 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1523 inputrec->fepvals->n_lambda);
1525 /* Kinetic energy data */
1526 gmx_ekindata_t ekind;
1527 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1529 /* Set up interactive MD (IMD) */
1531 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1532 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1533 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1535 if (DOMAINDECOMP(cr))
1537 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1538 /* This call is not included in init_domain_decomposition mainly
1539 * because fr->cginfo_mb is set later.
1541 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1542 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1545 // TODO This is not the right place to manage the lifetime of
1546 // this data structure, but currently it's the easiest way to
1548 MdrunScheduleWorkload runScheduleWork;
1549 // Also populates the simulation constant workload description.
1550 runScheduleWork.simulationWork =
1551 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1552 useGpuForUpdate, devFlags.enableGpuBufferOps,
1553 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1555 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1556 if (gpusWereDetected
1557 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1558 || runScheduleWork.simulationWork.useGpuBufferOps))
1560 const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1561 const void* localStream =
1562 fr->nbv->gpu_nbv != nullptr
1563 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1565 const void* nonLocalStream =
1566 fr->nbv->gpu_nbv != nullptr
1567 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1569 const DeviceContext& deviceContext = *pme_gpu_get_device_context(fr->pmedata);
1570 const int paddingSize = pme_gpu_get_padding_size(fr->pmedata);
1571 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1572 ? GpuApiCallBehavior::Async
1573 : GpuApiCallBehavior::Sync;
1575 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1576 pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
1577 fr->stateGpu = stateGpu.get();
1580 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1581 SimulatorBuilder simulatorBuilder;
1583 // build and run simulator object based on user-input
1584 auto simulator = simulatorBuilder.build(
1585 useModularSimulator, fplog, cr, ms, mdlog, static_cast<int>(filenames.size()),
1586 filenames.data(), oenv, mdrunOptions, startingBehavior, vsite.get(), constr.get(),
1587 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1588 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1589 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1590 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1591 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1594 if (fr->pmePpCommGpu)
1596 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1597 fr->pmePpCommGpu.reset();
1600 if (inputrec->bPull)
1602 finish_pull(pull_work);
1604 finish_swapcoords(swap);
1608 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1610 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1611 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1614 wallcycle_stop(wcycle, ewcRUN);
1616 /* Finish up, write some stuff
1617 * if rerunMD, don't write last frame again
1619 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1620 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1622 // clean up cycle counter
1623 wallcycle_destroy(wcycle);
1628 gmx_pme_destroy(pmedata);
1632 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1633 // before we destroy the GPU context(s) in free_gpu().
1634 // Pinned buffers are associated with contexts in CUDA.
1635 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1636 mdAtoms.reset(nullptr);
1637 globalState.reset(nullptr);
1638 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1639 gpuBonded.reset(nullptr);
1640 /* Free pinned buffers in *fr */
1644 if (hwinfo->gpu_info.n_dev > 0)
1646 /* stop the GPU profiler (only CUDA) */
1650 /* With tMPI we need to wait for all ranks to finish deallocation before
1651 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1654 * This is not a concern in OpenCL where we use one context per rank which
1655 * is freed in nbnxn_gpu_free().
1657 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1658 * but it is easier and more futureproof to call it on the whole node.
1660 * Note that this function needs to be called even if GPUs are not used
1661 * in this run because the PME ranks have no knowledge of whether GPUs
1662 * are used or not, but all ranks need to enter the barrier below.
1663 * \todo Remove this physical node barrier after making sure
1664 * that it's not needed anymore (with a shared GPU run).
1668 physicalNodeComm.barrier();
1671 free_gpu(nonbondedDeviceInfo);
1672 free_gpu(pmeDeviceInfo);
1677 free_membed(membed);
1680 /* Does what it says */
1681 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1682 walltime_accounting_destroy(walltime_accounting);
1684 // Ensure log file content is written
1687 gmx_fio_flush(logFileHandle);
1690 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1691 * exceptions were enabled before function was called. */
1694 gmx_fedisableexcept();
1697 auto rc = static_cast<int>(gmx_get_stop_condition());
1700 /* we need to join all threads. The sub-threads join when they
1701 exit this function, but the master thread needs to be told to
1703 if (PAR(cr) && MASTER(cr))
1711 Mdrunner::~Mdrunner()
1713 // Clean up of the Manager.
1714 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1715 // but okay as long as threads synchronize some time before adding or accessing
1716 // a new set of restraints.
1717 if (restraintManager_)
1719 restraintManager_->clear();
1720 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1721 "restraints added during runner life time should be cleared at runner "
1726 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1728 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1729 // Not sure if this should be logged through the md logger or something else,
1730 // but it is helpful to have some sort of INFO level message sent somewhere.
1731 // std::cout << "Registering restraint named " << name << std::endl;
1733 // When multiple restraints are used, it may be wasteful to register them separately.
1734 // Maybe instead register an entire Restraint Manager as a force provider.
1735 restraintManager_->addToSpec(std::move(puller), name);
1738 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1740 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1742 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1743 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1745 class Mdrunner::BuilderImplementation
1748 BuilderImplementation() = delete;
1749 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1750 ~BuilderImplementation();
1752 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1753 real forceWarningThreshold,
1754 StartingBehavior startingBehavior);
1756 void addDomdec(const DomdecOptions& options);
1758 void addVerletList(int nstlist);
1760 void addReplicaExchange(const ReplicaExchangeParameters& params);
1762 void addNonBonded(const char* nbpu_opt);
1764 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1766 void addBondedTaskAssignment(const char* bonded_opt);
1768 void addUpdateTaskAssignment(const char* update_opt);
1770 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1772 void addFilenames(ArrayRef<const t_filenm> filenames);
1774 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1776 void addLogFile(t_fileio* logFileHandle);
1778 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1783 // Default parameters copied from runner.h
1784 // \todo Clarify source(s) of default parameters.
1786 const char* nbpu_opt_ = nullptr;
1787 const char* pme_opt_ = nullptr;
1788 const char* pme_fft_opt_ = nullptr;
1789 const char* bonded_opt_ = nullptr;
1790 const char* update_opt_ = nullptr;
1792 MdrunOptions mdrunOptions_;
1794 DomdecOptions domdecOptions_;
1796 ReplicaExchangeParameters replicaExchangeParameters_;
1798 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1801 //! Multisim communicator handle.
1802 gmx_multisim_t* multiSimulation_;
1804 //! mdrun communicator
1805 MPI_Comm communicator_ = MPI_COMM_NULL;
1807 //! Print a warning if any force is larger than this (in kJ/mol nm).
1808 real forceWarningThreshold_ = -1;
1810 //! Whether the simulation will start afresh, or restart with/without appending.
1811 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1813 //! The modules that comprise the functionality of mdrun.
1814 std::unique_ptr<MDModules> mdModules_;
1816 //! \brief Parallelism information.
1817 gmx_hw_opt_t hardwareOptions_;
1819 //! filename options for simulation.
1820 ArrayRef<const t_filenm> filenames_;
1822 /*! \brief Handle to output environment.
1824 * \todo gmx_output_env_t needs lifetime management.
1826 gmx_output_env_t* outputEnvironment_ = nullptr;
1828 /*! \brief Non-owning handle to MD log file.
1830 * \todo Context should own output facilities for client.
1831 * \todo Improve log file handle management.
1833 * Code managing the FILE* relies on the ability to set it to
1834 * nullptr to check whether the filehandle is valid.
1836 t_fileio* logFileHandle_ = nullptr;
1839 * \brief Builder for simulation stop signal handler.
1841 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1844 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1845 compat::not_null<SimulationContext*> context) :
1846 mdModules_(std::move(mdModules))
1848 communicator_ = context->communicator_;
1849 multiSimulation_ = context->multiSimulation_.get();
1852 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1854 Mdrunner::BuilderImplementation&
1855 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1856 const real forceWarningThreshold,
1857 const StartingBehavior startingBehavior)
1859 mdrunOptions_ = options;
1860 forceWarningThreshold_ = forceWarningThreshold;
1861 startingBehavior_ = startingBehavior;
1865 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1867 domdecOptions_ = options;
1870 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1875 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1877 replicaExchangeParameters_ = params;
1880 Mdrunner Mdrunner::BuilderImplementation::build()
1882 auto newRunner = Mdrunner(std::move(mdModules_));
1884 newRunner.mdrunOptions = mdrunOptions_;
1885 newRunner.pforce = forceWarningThreshold_;
1886 newRunner.startingBehavior = startingBehavior_;
1887 newRunner.domdecOptions = domdecOptions_;
1889 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1890 newRunner.hw_opt = hardwareOptions_;
1892 // No invariant to check. This parameter exists to optionally override other behavior.
1893 newRunner.nstlist_cmdline = nstlist_;
1895 newRunner.replExParams = replicaExchangeParameters_;
1897 newRunner.filenames = filenames_;
1899 newRunner.communicator = communicator_;
1901 // nullptr is a valid value for the multisim handle
1902 newRunner.ms = multiSimulation_;
1904 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1905 // \todo Update sanity checking when output environment has clearly specified invariants.
1906 // Initialization and default values for oenv are not well specified in the current version.
1907 if (outputEnvironment_)
1909 newRunner.oenv = outputEnvironment_;
1913 GMX_THROW(gmx::APIError(
1914 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1917 newRunner.logFileHandle = logFileHandle_;
1921 newRunner.nbpu_opt = nbpu_opt_;
1925 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1928 if (pme_opt_ && pme_fft_opt_)
1930 newRunner.pme_opt = pme_opt_;
1931 newRunner.pme_fft_opt = pme_fft_opt_;
1935 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1940 newRunner.bonded_opt = bonded_opt_;
1944 GMX_THROW(gmx::APIError(
1945 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1950 newRunner.update_opt = update_opt_;
1954 GMX_THROW(gmx::APIError(
1955 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1959 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1961 if (stopHandlerBuilder_)
1963 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1967 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1973 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1975 nbpu_opt_ = nbpu_opt;
1978 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
1981 pme_fft_opt_ = pme_fft_opt;
1984 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1986 bonded_opt_ = bonded_opt;
1989 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
1991 update_opt_ = update_opt;
1994 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
1996 hardwareOptions_ = hardwareOptions;
1999 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2001 filenames_ = filenames;
2004 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2006 outputEnvironment_ = outputEnvironment;
2009 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2011 logFileHandle_ = logFileHandle;
2014 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2016 stopHandlerBuilder_ = std::move(builder);
2019 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2020 compat::not_null<SimulationContext*> context) :
2021 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2025 MdrunnerBuilder::~MdrunnerBuilder() = default;
2027 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2028 real forceWarningThreshold,
2029 const StartingBehavior startingBehavior)
2031 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2035 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2037 impl_->addDomdec(options);
2041 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2043 impl_->addVerletList(nstlist);
2047 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2049 impl_->addReplicaExchange(params);
2053 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2055 impl_->addNonBonded(nbpu_opt);
2059 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2061 // The builder method may become more general in the future, but in this version,
2062 // parameters for PME electrostatics are both required and the only parameters
2064 if (pme_opt && pme_fft_opt)
2066 impl_->addPME(pme_opt, pme_fft_opt);
2071 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2076 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2078 impl_->addBondedTaskAssignment(bonded_opt);
2082 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2084 impl_->addUpdateTaskAssignment(update_opt);
2088 Mdrunner MdrunnerBuilder::build()
2090 return impl_->build();
2093 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2095 impl_->addHardwareOptions(hardwareOptions);
2099 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2101 impl_->addFilenames(filenames);
2105 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2107 impl_->addOutputEnvironment(outputEnvironment);
2111 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2113 impl_->addLogFile(logFileHandle);
2117 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2119 impl_->addStopHandlerBuilder(std::move(builder));
2123 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2125 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;