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/gpu_utils.h"
78 #include "gromacs/hardware/cpuinfo.h"
79 #include "gromacs/hardware/detecthardware.h"
80 #include "gromacs/hardware/printhardware.h"
81 #include "gromacs/imd/imd.h"
82 #include "gromacs/listed_forces/disre.h"
83 #include "gromacs/listed_forces/gpubonded.h"
84 #include "gromacs/listed_forces/orires.h"
85 #include "gromacs/math/functions.h"
86 #include "gromacs/math/utilities.h"
87 #include "gromacs/math/vec.h"
88 #include "gromacs/mdlib/boxdeformation.h"
89 #include "gromacs/mdlib/broadcaststructs.h"
90 #include "gromacs/mdlib/calc_verletbuf.h"
91 #include "gromacs/mdlib/dispersioncorrection.h"
92 #include "gromacs/mdlib/enerdata_utils.h"
93 #include "gromacs/mdlib/force.h"
94 #include "gromacs/mdlib/forcerec.h"
95 #include "gromacs/mdlib/gmx_omp_nthreads.h"
96 #include "gromacs/mdlib/makeconstraints.h"
97 #include "gromacs/mdlib/md_support.h"
98 #include "gromacs/mdlib/mdatoms.h"
99 #include "gromacs/mdlib/membed.h"
100 #include "gromacs/mdlib/qmmm.h"
101 #include "gromacs/mdlib/sighandler.h"
102 #include "gromacs/mdlib/stophandler.h"
103 #include "gromacs/mdlib/updategroups.h"
104 #include "gromacs/mdlib/vsite.h"
105 #include "gromacs/mdrun/mdmodules.h"
106 #include "gromacs/mdrun/simulationcontext.h"
107 #include "gromacs/mdrunutility/handlerestart.h"
108 #include "gromacs/mdrunutility/logging.h"
109 #include "gromacs/mdrunutility/multisim.h"
110 #include "gromacs/mdrunutility/printtime.h"
111 #include "gromacs/mdrunutility/threadaffinity.h"
112 #include "gromacs/mdtypes/commrec.h"
113 #include "gromacs/mdtypes/enerdata.h"
114 #include "gromacs/mdtypes/fcdata.h"
115 #include "gromacs/mdtypes/forcerec.h"
116 #include "gromacs/mdtypes/group.h"
117 #include "gromacs/mdtypes/inputrec.h"
118 #include "gromacs/mdtypes/interaction_const.h"
119 #include "gromacs/mdtypes/md_enums.h"
120 #include "gromacs/mdtypes/mdatom.h"
121 #include "gromacs/mdtypes/mdrunoptions.h"
122 #include "gromacs/mdtypes/observableshistory.h"
123 #include "gromacs/mdtypes/simulation_workload.h"
124 #include "gromacs/mdtypes/state.h"
125 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
126 #include "gromacs/nbnxm/gpu_data_mgmt.h"
127 #include "gromacs/nbnxm/nbnxm.h"
128 #include "gromacs/nbnxm/pairlist_tuning.h"
129 #include "gromacs/pbcutil/pbc.h"
130 #include "gromacs/pulling/output.h"
131 #include "gromacs/pulling/pull.h"
132 #include "gromacs/pulling/pull_rotation.h"
133 #include "gromacs/restraint/manager.h"
134 #include "gromacs/restraint/restraintmdmodule.h"
135 #include "gromacs/restraint/restraintpotential.h"
136 #include "gromacs/swap/swapcoords.h"
137 #include "gromacs/taskassignment/decidegpuusage.h"
138 #include "gromacs/taskassignment/decidesimulationworkload.h"
139 #include "gromacs/taskassignment/resourcedivision.h"
140 #include "gromacs/taskassignment/taskassignment.h"
141 #include "gromacs/taskassignment/usergpuids.h"
142 #include "gromacs/timing/gpu_timing.h"
143 #include "gromacs/timing/wallcycle.h"
144 #include "gromacs/timing/wallcyclereporting.h"
145 #include "gromacs/topology/mtop_util.h"
146 #include "gromacs/trajectory/trajectoryframe.h"
147 #include "gromacs/utility/basenetwork.h"
148 #include "gromacs/utility/cstringutil.h"
149 #include "gromacs/utility/exceptions.h"
150 #include "gromacs/utility/fatalerror.h"
151 #include "gromacs/utility/filestream.h"
152 #include "gromacs/utility/gmxassert.h"
153 #include "gromacs/utility/gmxmpi.h"
154 #include "gromacs/utility/keyvaluetree.h"
155 #include "gromacs/utility/logger.h"
156 #include "gromacs/utility/loggerbuilder.h"
157 #include "gromacs/utility/mdmodulenotification.h"
158 #include "gromacs/utility/physicalnodecommunicator.h"
159 #include "gromacs/utility/pleasecite.h"
160 #include "gromacs/utility/programcontext.h"
161 #include "gromacs/utility/smalloc.h"
162 #include "gromacs/utility/stringutil.h"
164 #include "isimulator.h"
165 #include "replicaexchange.h"
166 #include "simulatorbuilder.h"
169 # include "corewrap.h"
176 /*! \brief Manage any development feature flag variables encountered
178 * The use of dev features indicated by environment variables is
179 * logged in order to ensure that runs with such features enabled can
180 * be identified from their log and standard output. Any cross
181 * dependencies are also checked, and if unsatisfied, a fatal error
184 * Note that some development features overrides are applied already here:
185 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
187 * \param[in] mdlog Logger object.
188 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
189 * \param[in] pmeRunMode The PME run mode for this run
190 * \returns The object populated with development feature flags.
192 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
193 const bool useGpuForNonbonded,
194 const PmeRunMode pmeRunMode)
196 DevelopmentFeatureFlags devFlags;
198 // Some builds of GCC 5 give false positive warnings that these
199 // getenv results are ignored when clearly they are used.
200 #pragma GCC diagnostic push
201 #pragma GCC diagnostic ignored "-Wunused-result"
202 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
203 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
204 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
205 devFlags.enableGpuHaloExchange =
206 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
207 devFlags.enableGpuPmePPComm =
208 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
209 #pragma GCC diagnostic pop
211 if (devFlags.enableGpuBufferOps)
213 GMX_LOG(mdlog.warning)
215 .appendTextFormatted(
216 "This run uses the 'GPU buffer ops' feature, enabled by the "
217 "GMX_USE_GPU_BUFFER_OPS environment variable.");
220 if (devFlags.forceGpuUpdateDefault)
222 GMX_LOG(mdlog.warning)
224 .appendTextFormatted(
225 "This run will default to '-update gpu' as requested by the "
226 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
227 "decomposition lacks substantial testing and should be used with caution.");
230 if (devFlags.enableGpuHaloExchange)
232 if (useGpuForNonbonded)
234 if (!devFlags.enableGpuBufferOps)
236 GMX_LOG(mdlog.warning)
238 .appendTextFormatted(
239 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
240 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
241 devFlags.enableGpuBufferOps = true;
243 GMX_LOG(mdlog.warning)
245 .appendTextFormatted(
246 "This run has requested the 'GPU halo exchange' feature, enabled by "
248 "GMX_GPU_DD_COMMS environment variable.");
252 GMX_LOG(mdlog.warning)
254 .appendTextFormatted(
255 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
256 "halo exchange' feature will not be enabled as nonbonded interactions "
257 "are not offloaded.");
258 devFlags.enableGpuHaloExchange = false;
262 if (devFlags.enableGpuPmePPComm)
264 if (pmeRunMode == PmeRunMode::GPU)
266 GMX_LOG(mdlog.warning)
268 .appendTextFormatted(
269 "This run uses the 'GPU PME-PP communications' feature, enabled "
270 "by the GMX_GPU_PME_PP_COMMS environment variable.");
274 std::string clarification;
275 if (pmeRunMode == PmeRunMode::Mixed)
278 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
283 clarification = "PME is not offloaded to the GPU.";
285 GMX_LOG(mdlog.warning)
288 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
289 "'GPU PME-PP communications' feature was not enabled as "
291 devFlags.enableGpuPmePPComm = false;
298 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
300 * Used to ensure that the master thread does not modify mdrunner during copy
301 * on the spawned threads. */
302 static void threadMpiMdrunnerAccessBarrier()
305 MPI_Barrier(MPI_COMM_WORLD);
309 Mdrunner Mdrunner::cloneOnSpawnedThread() const
311 auto newRunner = Mdrunner(std::make_unique<MDModules>());
313 // All runners in the same process share a restraint manager resource because it is
314 // part of the interface to the client code, which is associated only with the
315 // original thread. Handles to the same resources can be obtained by copy.
317 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
320 // Copy members of master runner.
321 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
322 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
323 newRunner.hw_opt = hw_opt;
324 newRunner.filenames = filenames;
326 newRunner.oenv = oenv;
327 newRunner.mdrunOptions = mdrunOptions;
328 newRunner.domdecOptions = domdecOptions;
329 newRunner.nbpu_opt = nbpu_opt;
330 newRunner.pme_opt = pme_opt;
331 newRunner.pme_fft_opt = pme_fft_opt;
332 newRunner.bonded_opt = bonded_opt;
333 newRunner.update_opt = update_opt;
334 newRunner.nstlist_cmdline = nstlist_cmdline;
335 newRunner.replExParams = replExParams;
336 newRunner.pforce = pforce;
337 // Give the spawned thread the newly created valid communicator
338 // for the simulation.
339 newRunner.communicator = MPI_COMM_WORLD;
341 newRunner.startingBehavior = startingBehavior;
342 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
344 threadMpiMdrunnerAccessBarrier();
349 /*! \brief The callback used for running on spawned threads.
351 * Obtains the pointer to the master mdrunner object from the one
352 * argument permitted to the thread-launch API call, copies it to make
353 * a new runner for this thread, reinitializes necessary data, and
354 * proceeds to the simulation. */
355 static void mdrunner_start_fn(const void* arg)
359 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
360 /* copy the arg list to make sure that it's thread-local. This
361 doesn't copy pointed-to items, of course; fnm, cr and fplog
362 are reset in the call below, all others should be const. */
363 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
366 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
370 void Mdrunner::spawnThreads(int numThreadsToLaunch)
373 /* now spawn new threads that start mdrunner_start_fn(), while
374 the main thread returns. Thread affinity is handled later. */
375 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
376 static_cast<const void*>(this))
379 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
382 // Give the master thread the newly created valid communicator for
384 communicator = MPI_COMM_WORLD;
385 threadMpiMdrunnerAccessBarrier();
387 GMX_UNUSED_VALUE(numThreadsToLaunch);
388 GMX_UNUSED_VALUE(mdrunner_start_fn);
394 /*! \brief Initialize variables for Verlet scheme simulation */
395 static void prepare_verlet_scheme(FILE* fplog,
399 const gmx_mtop_t* mtop,
401 bool makeGpuPairList,
402 const gmx::CpuInfo& cpuinfo)
404 // We checked the cut-offs in grompp, but double-check here.
405 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
406 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
408 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
409 "With Verlet lists and PME we should have rcoulomb>=rvdw");
413 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
414 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
416 /* For NVE simulations, we will retain the initial list buffer */
417 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
419 /* Update the Verlet buffer size for the current run setup */
421 /* Here we assume SIMD-enabled kernels are being used. But as currently
422 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
423 * and 4x2 gives a larger buffer than 4x4, this is ok.
425 ListSetupType listType =
426 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
427 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
429 const real rlist_new =
430 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
432 if (rlist_new != ir->rlist)
434 if (fplog != nullptr)
437 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
438 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
440 ir->rlist = rlist_new;
444 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
446 gmx_fatal(FARGS, "Can not set nstlist without %s",
447 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
450 if (EI_DYNAMICS(ir->eI))
452 /* Set or try nstlist values */
453 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
457 /*! \brief Override the nslist value in inputrec
459 * with value passed on the command line (if any)
461 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
465 /* override with anything else than the default -2 */
466 if (nsteps_cmdline > -2)
468 char sbuf_steps[STEPSTRSIZE];
469 char sbuf_msg[STRLEN];
471 ir->nsteps = nsteps_cmdline;
472 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
475 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
476 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
480 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
481 gmx_step_str(nsteps_cmdline, sbuf_steps));
484 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
486 else if (nsteps_cmdline < -2)
488 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
490 /* Do nothing if nsteps_cmdline == -2 */
496 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
498 * If not, and if a warning may be issued, logs a warning about
499 * falling back to CPU code. With thread-MPI, only the first
500 * call to this function should have \c issueWarning true. */
501 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
503 bool gpuIsUseful = true;
506 if (ir.opts.ngener - ir.nwall > 1)
508 /* The GPU code does not support more than one energy group.
509 * If the user requested GPUs explicitly, a fatal error is given later.
513 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
514 "For better performance, run on the GPU without energy groups and then do "
515 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
521 warning = "TPI is not implemented for GPUs.";
524 if (!gpuIsUseful && issueWarning)
526 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
532 //! Initializes the logger for mdrun.
533 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
535 gmx::LoggerBuilder builder;
536 if (fplog != nullptr)
538 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
540 if (isSimulationMasterRank)
542 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
544 return builder.build();
547 //! Make a TaskTarget from an mdrun argument string.
548 static TaskTarget findTaskTarget(const char* optionString)
550 TaskTarget returnValue = TaskTarget::Auto;
552 if (strncmp(optionString, "auto", 3) == 0)
554 returnValue = TaskTarget::Auto;
556 else if (strncmp(optionString, "cpu", 3) == 0)
558 returnValue = TaskTarget::Cpu;
560 else if (strncmp(optionString, "gpu", 3) == 0)
562 returnValue = TaskTarget::Gpu;
566 GMX_ASSERT(false, "Option string should have been checked for sanity already");
572 //! Finish run, aggregate data to print performance info.
573 static void finish_run(FILE* fplog,
574 const gmx::MDLogger& mdlog,
576 const t_inputrec* inputrec,
578 gmx_wallcycle_t wcycle,
579 gmx_walltime_accounting_t walltime_accounting,
580 nonbonded_verlet_t* nbv,
581 const gmx_pme_t* pme,
585 double nbfs = 0, mflop = 0;
586 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
587 elapsed_time_over_all_threads_over_all_ranks;
588 /* Control whether it is valid to print a report. Only the
589 simulation master may print, but it should not do so if the run
590 terminated e.g. before a scheduled reset step. This is
591 complicated by the fact that PME ranks are unaware of the
592 reason why they were sent a pmerecvqxFINISH. To avoid
593 communication deadlocks, we always do the communication for the
594 report, even if we've decided not to write the report, because
595 how long it takes to finish the run is not important when we've
596 decided not to report on the simulation performance.
598 Further, we only report performance for dynamical integrators,
599 because those are the only ones for which we plan to
600 consider doing any optimizations. */
601 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
603 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
605 GMX_LOG(mdlog.warning)
607 .appendText("Simulation ended prematurely, no performance report will be written.");
612 std::unique_ptr<t_nrnb> nrnbTotalStorage;
615 nrnbTotalStorage = std::make_unique<t_nrnb>();
616 nrnb_tot = nrnbTotalStorage.get();
618 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
626 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
627 elapsed_time_over_all_threads =
628 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
632 /* reduce elapsed_time over all MPI ranks in the current simulation */
633 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
635 elapsed_time_over_all_ranks /= cr->nnodes;
636 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
637 * current simulation. */
638 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
639 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
644 elapsed_time_over_all_ranks = elapsed_time;
645 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
650 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
653 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
655 print_dd_statistics(cr, inputrec, fplog);
658 /* TODO Move the responsibility for any scaling by thread counts
659 * to the code that handled the thread region, so that there's a
660 * mechanism to keep cycle counting working during the transition
661 * to task parallelism. */
662 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
663 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
664 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
665 nthreads_pp, nthreads_pme);
666 auto cycle_sum(wallcycle_sum(cr, wcycle));
670 auto nbnxn_gpu_timings =
671 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
672 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
674 if (pme_gpu_task_enabled(pme))
676 pme_gpu_get_timings(pme, &pme_gpu_timings);
678 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
679 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
682 if (EI_DYNAMICS(inputrec->eI))
684 delta_t = inputrec->delta_t;
689 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
690 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
691 delta_t, nbfs, mflop);
695 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
696 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
697 delta_t, nbfs, mflop);
702 int Mdrunner::mdrunner()
705 t_forcerec* fr = nullptr;
706 t_fcdata* fcd = nullptr;
707 real ewaldcoeff_q = 0;
708 real ewaldcoeff_lj = 0;
709 int nChargePerturbed = -1, nTypePerturbed = 0;
710 gmx_wallcycle_t wcycle;
711 gmx_walltime_accounting_t walltime_accounting = nullptr;
712 gmx_membed_t* membed = nullptr;
713 gmx_hw_info_t* hwinfo = nullptr;
715 /* CAUTION: threads may be started later on in this function, so
716 cr doesn't reflect the final parallel state right now */
719 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
720 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
721 const bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
722 const bool doRerun = mdrunOptions.rerun;
724 // Handle task-assignment related user options.
725 EmulateGpuNonbonded emulateGpuNonbonded =
726 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
728 std::vector<int> userGpuTaskAssignment;
731 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
733 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
734 auto nonbondedTarget = findTaskTarget(nbpu_opt);
735 auto pmeTarget = findTaskTarget(pme_opt);
736 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
737 auto bondedTarget = findTaskTarget(bonded_opt);
738 auto updateTarget = findTaskTarget(update_opt);
740 FILE* fplog = nullptr;
741 // If we are appending, we don't write log output because we need
742 // to check that the old log file matches what the checkpoint file
743 // expects. Otherwise, we should start to write log output now if
744 // there is a file ready for it.
745 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
747 fplog = gmx_fio_getfp(logFileHandle);
749 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
750 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
751 gmx::MDLogger mdlog(logOwner.logger());
753 // TODO The thread-MPI master rank makes a working
754 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
755 // after the threads have been launched. This works because no use
756 // is made of that communicator until after the execution paths
757 // have rejoined. But it is likely that we can improve the way
758 // this is expressed, e.g. by expressly running detection only the
759 // master rank for thread-MPI, rather than relying on the mutex
760 // and reference count.
761 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
762 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
764 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
766 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
768 // Print citation requests after all software/hardware printing
769 pleaseCiteGromacs(fplog);
771 // TODO Replace this by unique_ptr once t_inputrec is C++
772 t_inputrec inputrecInstance;
773 t_inputrec* inputrec = nullptr;
774 std::unique_ptr<t_state> globalState;
776 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
778 if (isSimulationMasterRank)
780 /* Only the master rank has the global state */
781 globalState = std::make_unique<t_state>();
783 /* Read (nearly) all data required for the simulation
784 * and keep the partly serialized tpr contents to send to other ranks later
786 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
787 &inputrecInstance, globalState.get(), &mtop);
788 inputrec = &inputrecInstance;
791 /* Check and update the hardware options for internal consistency */
792 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
795 if (GMX_THREAD_MPI && isSimulationMasterRank)
797 bool useGpuForNonbonded = false;
798 bool useGpuForPme = false;
801 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
803 // If the user specified the number of ranks, then we must
804 // respect that, but in default mode, we need to allow for
805 // the number of GPUs to choose the number of ranks.
806 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
807 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
808 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
809 canUseGpuForNonbonded,
810 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
811 hw_opt.nthreads_tmpi);
812 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
813 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
814 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
816 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
818 /* Determine how many thread-MPI ranks to start.
820 * TODO Over-writing the user-supplied value here does
821 * prevent any possible subsequent checks from working
823 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
824 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
826 // Now start the threads for thread MPI.
827 spawnThreads(hw_opt.nthreads_tmpi);
828 // The spawned threads enter mdrunner() and execution of
829 // master and spawned threads joins at the end of this block.
830 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
833 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
834 CommrecHandle crHandle = init_commrec(communicator, ms);
835 t_commrec* cr = crHandle.get();
836 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
840 /* now broadcast everything to the non-master nodes/threads: */
841 if (!isSimulationMasterRank)
843 inputrec = &inputrecInstance;
845 init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
847 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
848 partialDeserializedTpr.reset(nullptr);
850 // Now the number of ranks is known to all ranks, and each knows
851 // the inputrec read by the master rank. The ranks can now all run
852 // the task-deciding functions and will agree on the result
853 // without needing to communicate.
854 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
856 // Note that these variables describe only their own node.
858 // Note that when bonded interactions run on a GPU they always run
859 // alongside a nonbonded task, so do not influence task assignment
860 // even though they affect the force calculation workload.
861 bool useGpuForNonbonded = false;
862 bool useGpuForPme = false;
863 bool useGpuForBonded = false;
864 bool useGpuForUpdate = false;
865 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
868 // It's possible that there are different numbers of GPUs on
869 // different nodes, which is the user's responsibility to
870 // handle. If unsuitable, we will notice that during task
872 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
873 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
874 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
875 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
876 useGpuForPme = decideWhetherToUseGpusForPme(
877 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
878 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
879 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
880 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
881 useGpuForBonded = decideWhetherToUseGpusForBonded(
882 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
883 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
884 domdecOptions.numPmeRanks, gpusWereDetected);
886 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
888 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
890 // Initialize development feature flags that enabled by environment variable
891 // and report those features that are enabled.
892 const DevelopmentFeatureFlags devFlags =
893 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
895 const bool useModularSimulator = checkUseModularSimulator(
896 false, inputrec, doRerun, mtop, ms, replExParams, nullptr, doEssentialDynamics, doMembed);
899 // TODO: hide restraint implementation details from Mdrunner.
900 // There is nothing unique about restraints at this point as far as the
901 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
902 // factory functions from the SimulationContext on which to call mdModules_->add().
903 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
904 for (auto&& restraint : restraintManager_->getRestraints())
906 auto module = RestraintMDModule::create(restraint, restraint->sites());
907 mdModules_->add(std::move(module));
910 // TODO: Error handling
911 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
912 // now that the MdModules know their options, they know which callbacks to sign up to
913 mdModules_->subscribeToSimulationSetupNotifications();
914 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
916 if (inputrec->internalParameters != nullptr)
918 mdModulesNotifier.notify(*inputrec->internalParameters);
921 if (fplog != nullptr)
923 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
924 fprintf(fplog, "\n");
929 /* In rerun, set velocities to zero if present */
930 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
932 // rerun does not use velocities
936 "Rerun trajectory contains velocities. Rerun does only evaluate "
937 "potential energy and forces. The velocities will be ignored.");
938 for (int i = 0; i < globalState->natoms; i++)
940 clear_rvec(globalState->v[i]);
942 globalState->flags &= ~(1 << estV);
945 /* now make sure the state is initialized and propagated */
946 set_state_entries(globalState.get(), inputrec, useModularSimulator);
949 /* NM and TPI parallelize over force/energy calculations, not atoms,
950 * so we need to initialize and broadcast the global state.
952 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
956 globalState = std::make_unique<t_state>();
958 broadcastStateWithoutDynamics(cr, globalState.get());
961 /* A parallel command line option consistency check that we can
962 only do after any threads have started. */
964 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
965 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
968 "The -dd or -npme option request a parallel simulation, "
970 "but %s was compiled without threads or MPI enabled",
971 output_env_get_program_display_name(oenv));
973 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
975 "but %s was not started through mpirun/mpiexec or only one rank was requested "
976 "through mpirun/mpiexec",
977 output_env_get_program_display_name(oenv));
981 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
984 "The .mdp file specified an energy mininization or normal mode algorithm, and "
985 "these are not compatible with mdrun -rerun");
988 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
990 if (domdecOptions.numPmeRanks > 0)
992 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
993 "PME-only ranks are requested, but the system does not use PME "
994 "for electrostatics or LJ");
997 domdecOptions.numPmeRanks = 0;
1000 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1002 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1003 * improve performance with many threads per GPU, since our OpenMP
1004 * scaling is bad, but it's difficult to automate the setup.
1006 domdecOptions.numPmeRanks = 0;
1010 if (domdecOptions.numPmeRanks < 0)
1012 domdecOptions.numPmeRanks = 0;
1013 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1017 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1018 "PME GPU decomposition is not supported");
1025 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1029 /* NMR restraints must be initialized before load_checkpoint,
1030 * since with time averaging the history is added to t_state.
1031 * For proper consistency check we therefore need to extend
1033 * So the PME-only nodes (if present) will also initialize
1034 * the distance restraints.
1038 /* This needs to be called before read_checkpoint to extend the state */
1039 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1041 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1043 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1045 ObservablesHistory observablesHistory = {};
1047 if (startingBehavior != StartingBehavior::NewSimulation)
1049 /* Check if checkpoint file exists before doing continuation.
1050 * This way we can use identical input options for the first and subsequent runs...
1052 if (mdrunOptions.numStepsCommandline > -2)
1054 /* Temporarily set the number of steps to unlimited to avoid
1055 * triggering the nsteps check in load_checkpoint().
1056 * This hack will go away soon when the -nsteps option is removed.
1058 inputrec->nsteps = -1;
1061 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1062 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1063 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1065 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1067 // Now we can start normal logging to the truncated log file.
1068 fplog = gmx_fio_getfp(logFileHandle);
1069 prepareLogAppending(fplog);
1070 logOwner = buildLogger(fplog, MASTER(cr));
1071 mdlog = logOwner.logger();
1075 if (mdrunOptions.numStepsCommandline > -2)
1080 "The -nsteps functionality is deprecated, and may be removed in a future "
1082 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1085 /* override nsteps with value set on the commandline */
1086 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1090 copy_mat(globalState->box, box);
1095 gmx_bcast(sizeof(box), box, cr);
1098 if (inputrec->cutoff_scheme != ecutsVERLET)
1101 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1102 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1104 /* Update rlist and nstlist. */
1105 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1106 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1109 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1110 // This builder is necessary while we have multi-part construction
1111 // of DD. Before DD is constructed, we use the existence of
1112 // the builder object to indicate that further construction of DD
1114 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1115 if (useDomainDecomposition)
1117 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1118 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1119 positionsFromStatePointer(globalState.get()));
1123 /* PME, if used, is done on all nodes with 1D decomposition */
1125 cr->duty = (DUTY_PP | DUTY_PME);
1127 if (inputrec->pbcType == PbcType::Screw)
1129 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1133 // Produce the task assignment for this rank.
1134 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1135 GpuTaskAssignments gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1136 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1137 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1138 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1139 // TODO cr->duty & DUTY_PME should imply that a PME
1140 // algorithm is active, but currently does not.
1141 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1143 // Get the device handles for the modules, nullptr when no task is assigned.
1144 // TODO: There should be only one DeviceInformation per rank.
1145 DeviceInformation* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1146 DeviceInformation* pmeDeviceInfo = gpuTaskAssignments.initPmeDevice();
1148 std::unique_ptr<DeviceContext> deviceContext = nullptr;
1151 deviceContext = std::make_unique<DeviceContext>(*pmeDeviceInfo);
1153 else if (nonbondedDeviceInfo)
1155 deviceContext = std::make_unique<DeviceContext>(*nonbondedDeviceInfo);
1158 // TODO Initialize GPU streams here.
1160 // TODO Currently this is always built, yet DD partition code
1161 // checks if it is built before using it. Probably it should
1162 // become an MDModule that is made only when another module
1163 // requires it (e.g. pull, CompEl, density fitting), so that we
1164 // don't update the local atom sets unilaterally every step.
1165 LocalAtomSetManager atomSets;
1168 // TODO Pass the GPU streams to ddBuilder to use in buffer
1169 // transfers (e.g. halo exchange)
1170 cr->dd = ddBuilder->build(&atomSets);
1171 // The builder's job is done, so destruct it
1172 ddBuilder.reset(nullptr);
1173 // Note that local state still does not exist yet.
1176 // The GPU update is decided here because we need to know whether the constraints or
1177 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1178 // defined). This is only known after DD is initialized, hence decision on using GPU
1179 // update is done so late.
1182 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1184 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1185 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1186 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1187 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1188 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1190 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1192 const bool printHostName = (cr->nnodes > 1);
1193 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1195 // If the user chose a task assignment, give them some hints
1196 // where appropriate.
1197 if (!userGpuTaskAssignment.empty())
1199 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1204 /* After possible communicator splitting in make_dd_communicators.
1205 * we can set up the intra/inter node communication.
1207 gmx_setup_nodecomm(fplog, cr);
1213 GMX_LOG(mdlog.warning)
1215 .appendTextFormatted(
1216 "This is simulation %d out of %d running as a composite GROMACS\n"
1217 "multi-simulation job. Setup for this simulation:\n",
1220 GMX_LOG(mdlog.warning)
1221 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1223 cr->nnodes == 1 ? "thread" : "threads"
1225 cr->nnodes == 1 ? "process" : "processes"
1231 // If mdrun -pin auto honors any affinity setting that already
1232 // exists. If so, it is nice to provide feedback about whether
1233 // that existing affinity setting was from OpenMP or something
1234 // else, so we run this code both before and after we initialize
1235 // the OpenMP support.
1236 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1237 /* Check and update the number of OpenMP threads requested */
1238 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1239 pmeRunMode, mtop, *inputrec);
1241 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1242 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1244 // Enable FP exception detection, but not in
1245 // Release mode and not for compilers with known buggy FP
1246 // exception support (clang with any optimization) or suspected
1247 // buggy FP exception support (gcc 7.* with optimization).
1248 #if !defined NDEBUG \
1249 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1250 && defined __OPTIMIZE__)
1251 const bool bEnableFPE = true;
1253 const bool bEnableFPE = false;
1255 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1258 gmx_feenableexcept();
1261 /* Now that we know the setup is consistent, check for efficiency */
1262 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1263 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1265 /* getting number of PP/PME threads on this MPI / tMPI rank.
1266 PME: env variable should be read only on one node to make sure it is
1267 identical everywhere;
1269 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1270 : gmx_omp_nthreads_get(emntPME);
1271 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1272 physicalNodeComm, mdlog);
1274 // Enable Peer access between GPUs where available
1275 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1276 // any of the GPU communication features are active.
1277 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1278 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1280 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1283 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1285 /* Before setting affinity, check whether the affinity has changed
1286 * - which indicates that probably the OpenMP library has changed it
1287 * since we first checked).
1289 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1291 int numThreadsOnThisNode, intraNodeThreadOffset;
1292 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1293 &intraNodeThreadOffset);
1295 /* Set the CPU affinity */
1296 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1297 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1300 if (mdrunOptions.timingOptions.resetStep > -1)
1305 "The -resetstep functionality is deprecated, and may be removed in a "
1308 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1312 /* Master synchronizes its value of reset_counters with all nodes
1313 * including PME only nodes */
1314 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1315 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1316 wcycle_set_reset_counters(wcycle, reset_counters);
1319 // Membrane embedding must be initialized before we call init_forcerec()
1324 fprintf(stderr, "Initializing membed");
1326 /* Note that membed cannot work in parallel because mtop is
1327 * changed here. Fix this if we ever want to make it run with
1328 * multiple ranks. */
1329 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1330 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1333 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1334 std::unique_ptr<MDAtoms> mdAtoms;
1335 std::unique_ptr<gmx_vsite_t> vsite;
1336 std::unique_ptr<GpuBonded> gpuBonded;
1339 if (thisRankHasDuty(cr, DUTY_PP))
1341 mdModulesNotifier.notify(*cr);
1342 mdModulesNotifier.notify(&atomSets);
1343 mdModulesNotifier.notify(inputrec->pbcType);
1344 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1345 /* Initiate forcerecord */
1346 fr = new t_forcerec;
1347 fr->forceProviders = mdModules_->initForceProviders();
1348 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1349 opt2fn("-table", filenames.size(), filenames.data()),
1350 opt2fn("-tablep", filenames.size(), filenames.data()),
1351 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1353 fr->deviceContext = deviceContext.get();
1355 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1358 deviceContext != nullptr,
1359 "Device context can not be nullptr when PME-PP direct communications object.");
1360 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1361 cr->mpi_comm_mysim, cr->dd->pme_nodeid, *deviceContext);
1364 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, nonbondedDeviceInfo,
1365 fr->deviceContext, &mtop, box, wcycle);
1366 if (useGpuForBonded)
1368 auto stream = havePPDomainDecomposition(cr)
1369 ? Nbnxm::gpu_get_command_stream(
1370 fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1371 : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1372 gmx::InteractionLocality::Local);
1374 fr->deviceContext != nullptr,
1375 "Device context can not be nullptr when computing bonded interactions on GPU.");
1376 gpuBonded = std::make_unique<GpuBonded>(mtop.ffparams, *fr->deviceContext, stream, wcycle);
1377 fr->gpuBonded = gpuBonded.get();
1380 /* Initialize the mdAtoms structure.
1381 * mdAtoms is not filled with atom data,
1382 * as this can not be done now with domain decomposition.
1384 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1385 if (globalState && thisRankHasPmeGpuTask)
1387 // The pinning of coordinates in the global state object works, because we only use
1388 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1389 // points to the global state object without DD.
1390 // FIXME: MD and EM separately set up the local state - this should happen in the same
1391 // function, which should also perform the pinning.
1392 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1395 /* Initialize the virtual site communication */
1396 vsite = initVsite(mtop, cr);
1398 calc_shifts(box, fr->shift_vec);
1400 /* With periodic molecules the charge groups should be whole at start up
1401 * and the virtual sites should not be far from their proper positions.
1403 if (!inputrec->bContinuation && MASTER(cr)
1404 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1406 /* Make molecules whole at start of run */
1407 if (fr->pbcType != PbcType::No)
1409 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1413 /* Correct initial vsite positions are required
1414 * for the initial distribution in the domain decomposition
1415 * and for the initial shell prediction.
1417 constructVsitesGlobal(mtop, globalState->x);
1421 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1423 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1424 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1429 /* This is a PME only node */
1431 GMX_ASSERT(globalState == nullptr,
1432 "We don't need the state on a PME only rank and expect it to be unitialized");
1434 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1435 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1438 gmx_pme_t* sepPmeData = nullptr;
1439 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1440 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1441 "Double-checking that only PME-only ranks have no forcerec");
1442 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1444 // TODO should live in ewald module once its testing is improved
1446 // Later, this program could contain kernels that might be later
1447 // re-used as auto-tuning progresses, or subsequent simulations
1449 PmeGpuProgramStorage pmeGpuProgram;
1450 if (thisRankHasPmeGpuTask)
1453 pmeDeviceInfo != nullptr,
1454 "Device information can not be nullptr when building PME GPU program object.");
1456 deviceContext != nullptr,
1457 "Device context can not be nullptr when building PME GPU program object.");
1458 pmeGpuProgram = buildPmeGpuProgram(*pmeDeviceInfo, *deviceContext);
1461 /* Initiate PME if necessary,
1462 * either on all nodes or on dedicated PME nodes only. */
1463 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1465 if (mdAtoms && mdAtoms->mdatoms())
1467 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1468 if (EVDW_PME(inputrec->vdwtype))
1470 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1473 if (cr->npmenodes > 0)
1475 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1476 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1477 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1480 if (thisRankHasDuty(cr, DUTY_PME))
1484 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1485 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1486 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1487 nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1489 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1494 if (EI_DYNAMICS(inputrec->eI))
1496 /* Turn on signal handling on all nodes */
1498 * (A user signal from the PME nodes (if any)
1499 * is communicated to the PP nodes.
1501 signal_handler_install();
1504 pull_t* pull_work = nullptr;
1505 if (thisRankHasDuty(cr, DUTY_PP))
1507 /* Assumes uniform use of the number of OpenMP threads */
1508 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1510 if (inputrec->bPull)
1512 /* Initialize pull code */
1513 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1514 inputrec->fepvals->init_lambda);
1515 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1517 initPullHistory(pull_work, &observablesHistory);
1519 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1521 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1525 std::unique_ptr<EnforcedRotation> enforcedRotation;
1528 /* Initialize enforced rotation code */
1530 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1531 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1534 t_swap* swap = nullptr;
1535 if (inputrec->eSwapCoords != eswapNO)
1537 /* Initialize ion swapping code */
1538 swap = init_swapcoords(fplog, inputrec,
1539 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1540 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1541 oenv, mdrunOptions, startingBehavior);
1544 /* Let makeConstraints know whether we have essential dynamics constraints. */
1545 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1546 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1548 /* Energy terms and groups */
1549 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1550 inputrec->fepvals->n_lambda);
1552 /* Kinetic energy data */
1553 gmx_ekindata_t ekind;
1554 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1556 /* Set up interactive MD (IMD) */
1558 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1559 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1560 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1562 if (DOMAINDECOMP(cr))
1564 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1565 /* This call is not included in init_domain_decomposition mainly
1566 * because fr->cginfo_mb is set later.
1568 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1569 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1572 // TODO This is not the right place to manage the lifetime of
1573 // this data structure, but currently it's the easiest way to
1575 MdrunScheduleWorkload runScheduleWork;
1576 // Also populates the simulation constant workload description.
1577 runScheduleWork.simulationWork =
1578 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1579 useGpuForUpdate, devFlags.enableGpuBufferOps,
1580 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1582 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1583 if (gpusWereDetected
1584 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1585 || runScheduleWork.simulationWork.useGpuBufferOps))
1587 const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1588 const void* localStream =
1589 fr->nbv->gpu_nbv != nullptr
1590 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1592 const void* nonLocalStream =
1593 fr->nbv->gpu_nbv != nullptr
1594 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1596 const int paddingSize = pme_gpu_get_padding_size(fr->pmedata);
1597 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1598 ? GpuApiCallBehavior::Async
1599 : GpuApiCallBehavior::Sync;
1601 deviceContext != nullptr,
1602 "Device context can not be nullptr when building GPU propagator data object.");
1603 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1604 pmeStream, localStream, nonLocalStream, *deviceContext, transferKind,
1605 paddingSize, wcycle);
1606 fr->stateGpu = stateGpu.get();
1609 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1610 SimulatorBuilder simulatorBuilder;
1612 // build and run simulator object based on user-input
1613 auto simulator = simulatorBuilder.build(
1614 useModularSimulator, fplog, cr, ms, mdlog, static_cast<int>(filenames.size()),
1615 filenames.data(), oenv, mdrunOptions, startingBehavior, vsite.get(), constr.get(),
1616 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1617 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1618 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1619 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1620 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1623 if (fr->pmePpCommGpu)
1625 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1626 fr->pmePpCommGpu.reset();
1629 if (inputrec->bPull)
1631 finish_pull(pull_work);
1633 finish_swapcoords(swap);
1637 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1639 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1640 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode,
1641 deviceContext.get());
1644 wallcycle_stop(wcycle, ewcRUN);
1646 /* Finish up, write some stuff
1647 * if rerunMD, don't write last frame again
1649 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1650 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1652 // clean up cycle counter
1653 wallcycle_destroy(wcycle);
1658 gmx_pme_destroy(pmedata);
1662 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1663 // before we destroy the GPU context(s) in free_gpu().
1664 // Pinned buffers are associated with contexts in CUDA.
1665 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1666 mdAtoms.reset(nullptr);
1667 globalState.reset(nullptr);
1668 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1669 gpuBonded.reset(nullptr);
1670 /* Free pinned buffers in *fr */
1674 if (hwinfo->gpu_info.n_dev > 0)
1676 /* stop the GPU profiler (only CUDA) */
1680 /* With tMPI we need to wait for all ranks to finish deallocation before
1681 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1684 * This is not a concern in OpenCL where we use one context per rank which
1685 * is freed in nbnxn_gpu_free().
1687 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1688 * but it is easier and more futureproof to call it on the whole node.
1690 * Note that this function needs to be called even if GPUs are not used
1691 * in this run because the PME ranks have no knowledge of whether GPUs
1692 * are used or not, but all ranks need to enter the barrier below.
1693 * \todo Remove this physical node barrier after making sure
1694 * that it's not needed anymore (with a shared GPU run).
1698 physicalNodeComm.barrier();
1701 free_gpu(nonbondedDeviceInfo);
1702 free_gpu(pmeDeviceInfo);
1703 deviceContext.reset(nullptr);
1708 free_membed(membed);
1711 /* Does what it says */
1712 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1713 walltime_accounting_destroy(walltime_accounting);
1715 // Ensure log file content is written
1718 gmx_fio_flush(logFileHandle);
1721 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1722 * exceptions were enabled before function was called. */
1725 gmx_fedisableexcept();
1728 auto rc = static_cast<int>(gmx_get_stop_condition());
1731 /* we need to join all threads. The sub-threads join when they
1732 exit this function, but the master thread needs to be told to
1734 if (PAR(cr) && MASTER(cr))
1742 Mdrunner::~Mdrunner()
1744 // Clean up of the Manager.
1745 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1746 // but okay as long as threads synchronize some time before adding or accessing
1747 // a new set of restraints.
1748 if (restraintManager_)
1750 restraintManager_->clear();
1751 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1752 "restraints added during runner life time should be cleared at runner "
1757 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1759 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1760 // Not sure if this should be logged through the md logger or something else,
1761 // but it is helpful to have some sort of INFO level message sent somewhere.
1762 // std::cout << "Registering restraint named " << name << std::endl;
1764 // When multiple restraints are used, it may be wasteful to register them separately.
1765 // Maybe instead register an entire Restraint Manager as a force provider.
1766 restraintManager_->addToSpec(std::move(puller), name);
1769 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1771 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1773 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1774 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1776 class Mdrunner::BuilderImplementation
1779 BuilderImplementation() = delete;
1780 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1781 ~BuilderImplementation();
1783 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1784 real forceWarningThreshold,
1785 StartingBehavior startingBehavior);
1787 void addDomdec(const DomdecOptions& options);
1789 void addVerletList(int nstlist);
1791 void addReplicaExchange(const ReplicaExchangeParameters& params);
1793 void addNonBonded(const char* nbpu_opt);
1795 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1797 void addBondedTaskAssignment(const char* bonded_opt);
1799 void addUpdateTaskAssignment(const char* update_opt);
1801 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1803 void addFilenames(ArrayRef<const t_filenm> filenames);
1805 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1807 void addLogFile(t_fileio* logFileHandle);
1809 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1814 // Default parameters copied from runner.h
1815 // \todo Clarify source(s) of default parameters.
1817 const char* nbpu_opt_ = nullptr;
1818 const char* pme_opt_ = nullptr;
1819 const char* pme_fft_opt_ = nullptr;
1820 const char* bonded_opt_ = nullptr;
1821 const char* update_opt_ = nullptr;
1823 MdrunOptions mdrunOptions_;
1825 DomdecOptions domdecOptions_;
1827 ReplicaExchangeParameters replicaExchangeParameters_;
1829 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1832 //! Multisim communicator handle.
1833 gmx_multisim_t* multiSimulation_;
1835 //! mdrun communicator
1836 MPI_Comm communicator_ = MPI_COMM_NULL;
1838 //! Print a warning if any force is larger than this (in kJ/mol nm).
1839 real forceWarningThreshold_ = -1;
1841 //! Whether the simulation will start afresh, or restart with/without appending.
1842 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1844 //! The modules that comprise the functionality of mdrun.
1845 std::unique_ptr<MDModules> mdModules_;
1847 //! \brief Parallelism information.
1848 gmx_hw_opt_t hardwareOptions_;
1850 //! filename options for simulation.
1851 ArrayRef<const t_filenm> filenames_;
1853 /*! \brief Handle to output environment.
1855 * \todo gmx_output_env_t needs lifetime management.
1857 gmx_output_env_t* outputEnvironment_ = nullptr;
1859 /*! \brief Non-owning handle to MD log file.
1861 * \todo Context should own output facilities for client.
1862 * \todo Improve log file handle management.
1864 * Code managing the FILE* relies on the ability to set it to
1865 * nullptr to check whether the filehandle is valid.
1867 t_fileio* logFileHandle_ = nullptr;
1870 * \brief Builder for simulation stop signal handler.
1872 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1875 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1876 compat::not_null<SimulationContext*> context) :
1877 mdModules_(std::move(mdModules))
1879 communicator_ = context->communicator_;
1880 multiSimulation_ = context->multiSimulation_.get();
1883 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1885 Mdrunner::BuilderImplementation&
1886 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1887 const real forceWarningThreshold,
1888 const StartingBehavior startingBehavior)
1890 mdrunOptions_ = options;
1891 forceWarningThreshold_ = forceWarningThreshold;
1892 startingBehavior_ = startingBehavior;
1896 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1898 domdecOptions_ = options;
1901 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1906 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1908 replicaExchangeParameters_ = params;
1911 Mdrunner Mdrunner::BuilderImplementation::build()
1913 auto newRunner = Mdrunner(std::move(mdModules_));
1915 newRunner.mdrunOptions = mdrunOptions_;
1916 newRunner.pforce = forceWarningThreshold_;
1917 newRunner.startingBehavior = startingBehavior_;
1918 newRunner.domdecOptions = domdecOptions_;
1920 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1921 newRunner.hw_opt = hardwareOptions_;
1923 // No invariant to check. This parameter exists to optionally override other behavior.
1924 newRunner.nstlist_cmdline = nstlist_;
1926 newRunner.replExParams = replicaExchangeParameters_;
1928 newRunner.filenames = filenames_;
1930 newRunner.communicator = communicator_;
1932 // nullptr is a valid value for the multisim handle
1933 newRunner.ms = multiSimulation_;
1935 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1936 // \todo Update sanity checking when output environment has clearly specified invariants.
1937 // Initialization and default values for oenv are not well specified in the current version.
1938 if (outputEnvironment_)
1940 newRunner.oenv = outputEnvironment_;
1944 GMX_THROW(gmx::APIError(
1945 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1948 newRunner.logFileHandle = logFileHandle_;
1952 newRunner.nbpu_opt = nbpu_opt_;
1956 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1959 if (pme_opt_ && pme_fft_opt_)
1961 newRunner.pme_opt = pme_opt_;
1962 newRunner.pme_fft_opt = pme_fft_opt_;
1966 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1971 newRunner.bonded_opt = bonded_opt_;
1975 GMX_THROW(gmx::APIError(
1976 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1981 newRunner.update_opt = update_opt_;
1985 GMX_THROW(gmx::APIError(
1986 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1990 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1992 if (stopHandlerBuilder_)
1994 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1998 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2004 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2006 nbpu_opt_ = nbpu_opt;
2009 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2012 pme_fft_opt_ = pme_fft_opt;
2015 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2017 bonded_opt_ = bonded_opt;
2020 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2022 update_opt_ = update_opt;
2025 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2027 hardwareOptions_ = hardwareOptions;
2030 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2032 filenames_ = filenames;
2035 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2037 outputEnvironment_ = outputEnvironment;
2040 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2042 logFileHandle_ = logFileHandle;
2045 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2047 stopHandlerBuilder_ = std::move(builder);
2050 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2051 compat::not_null<SimulationContext*> context) :
2052 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2056 MdrunnerBuilder::~MdrunnerBuilder() = default;
2058 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2059 real forceWarningThreshold,
2060 const StartingBehavior startingBehavior)
2062 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2066 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2068 impl_->addDomdec(options);
2072 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2074 impl_->addVerletList(nstlist);
2078 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2080 impl_->addReplicaExchange(params);
2084 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2086 impl_->addNonBonded(nbpu_opt);
2090 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2092 // The builder method may become more general in the future, but in this version,
2093 // parameters for PME electrostatics are both required and the only parameters
2095 if (pme_opt && pme_fft_opt)
2097 impl_->addPME(pme_opt, pme_fft_opt);
2102 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2107 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2109 impl_->addBondedTaskAssignment(bonded_opt);
2113 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2115 impl_->addUpdateTaskAssignment(update_opt);
2119 Mdrunner MdrunnerBuilder::build()
2121 return impl_->build();
2124 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2126 impl_->addHardwareOptions(hardwareOptions);
2130 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2132 impl_->addFilenames(filenames);
2136 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2138 impl_->addOutputEnvironment(outputEnvironment);
2142 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2144 impl_->addLogFile(logFileHandle);
2148 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2150 impl_->addStopHandlerBuilder(std::move(builder));
2154 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2156 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;