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.h"
68 #include "gromacs/ewald/pme_gpu_program.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/mdrun/mdmodules.h"
104 #include "gromacs/mdrun/simulationcontext.h"
105 #include "gromacs/mdrunutility/handlerestart.h"
106 #include "gromacs/mdrunutility/logging.h"
107 #include "gromacs/mdrunutility/multisim.h"
108 #include "gromacs/mdrunutility/printtime.h"
109 #include "gromacs/mdrunutility/threadaffinity.h"
110 #include "gromacs/mdtypes/commrec.h"
111 #include "gromacs/mdtypes/enerdata.h"
112 #include "gromacs/mdtypes/fcdata.h"
113 #include "gromacs/mdtypes/group.h"
114 #include "gromacs/mdtypes/inputrec.h"
115 #include "gromacs/mdtypes/md_enums.h"
116 #include "gromacs/mdtypes/mdrunoptions.h"
117 #include "gromacs/mdtypes/observableshistory.h"
118 #include "gromacs/mdtypes/simulation_workload.h"
119 #include "gromacs/mdtypes/state.h"
120 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
121 #include "gromacs/nbnxm/gpu_data_mgmt.h"
122 #include "gromacs/nbnxm/nbnxm.h"
123 #include "gromacs/nbnxm/pairlist_tuning.h"
124 #include "gromacs/pbcutil/pbc.h"
125 #include "gromacs/pulling/output.h"
126 #include "gromacs/pulling/pull.h"
127 #include "gromacs/pulling/pull_rotation.h"
128 #include "gromacs/restraint/manager.h"
129 #include "gromacs/restraint/restraintmdmodule.h"
130 #include "gromacs/restraint/restraintpotential.h"
131 #include "gromacs/swap/swapcoords.h"
132 #include "gromacs/taskassignment/decidegpuusage.h"
133 #include "gromacs/taskassignment/decidesimulationworkload.h"
134 #include "gromacs/taskassignment/resourcedivision.h"
135 #include "gromacs/taskassignment/taskassignment.h"
136 #include "gromacs/taskassignment/usergpuids.h"
137 #include "gromacs/timing/gpu_timing.h"
138 #include "gromacs/timing/wallcycle.h"
139 #include "gromacs/timing/wallcyclereporting.h"
140 #include "gromacs/topology/mtop_util.h"
141 #include "gromacs/trajectory/trajectoryframe.h"
142 #include "gromacs/utility/basenetwork.h"
143 #include "gromacs/utility/cstringutil.h"
144 #include "gromacs/utility/exceptions.h"
145 #include "gromacs/utility/fatalerror.h"
146 #include "gromacs/utility/filestream.h"
147 #include "gromacs/utility/gmxassert.h"
148 #include "gromacs/utility/gmxmpi.h"
149 #include "gromacs/utility/keyvaluetree.h"
150 #include "gromacs/utility/logger.h"
151 #include "gromacs/utility/loggerbuilder.h"
152 #include "gromacs/utility/mdmodulenotification.h"
153 #include "gromacs/utility/physicalnodecommunicator.h"
154 #include "gromacs/utility/pleasecite.h"
155 #include "gromacs/utility/programcontext.h"
156 #include "gromacs/utility/smalloc.h"
157 #include "gromacs/utility/stringutil.h"
159 #include "isimulator.h"
160 #include "replicaexchange.h"
161 #include "simulatorbuilder.h"
164 # include "corewrap.h"
171 /*! \brief Manage any development feature flag variables encountered
173 * The use of dev features indicated by environment variables is
174 * logged in order to ensure that runs with such features enabled can
175 * be identified from their log and standard output. Any cross
176 * dependencies are also checked, and if unsatisfied, a fatal error
179 * Note that some development features overrides are applied already here:
180 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
182 * \param[in] mdlog Logger object.
183 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
184 * \param[in] pmeRunMode The PME run mode for this run
185 * \returns The object populated with development feature flags.
187 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
188 const bool useGpuForNonbonded,
189 const PmeRunMode pmeRunMode)
191 DevelopmentFeatureFlags devFlags;
193 // Some builds of GCC 5 give false positive warnings that these
194 // getenv results are ignored when clearly they are used.
195 #pragma GCC diagnostic push
196 #pragma GCC diagnostic ignored "-Wunused-result"
197 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
198 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
199 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
200 devFlags.enableGpuHaloExchange =
201 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
202 devFlags.enableGpuPmePPComm =
203 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
204 #pragma GCC diagnostic pop
206 if (devFlags.enableGpuBufferOps)
208 GMX_LOG(mdlog.warning)
210 .appendTextFormatted(
211 "This run uses the 'GPU buffer ops' feature, enabled by the "
212 "GMX_USE_GPU_BUFFER_OPS environment variable.");
215 if (devFlags.forceGpuUpdateDefault)
217 GMX_LOG(mdlog.warning)
219 .appendTextFormatted(
220 "This run will default to '-update gpu' as requested by the "
221 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
222 "decomposition lacks substantial testing and should be used with caution.");
225 if (devFlags.enableGpuHaloExchange)
227 if (useGpuForNonbonded)
229 if (!devFlags.enableGpuBufferOps)
231 GMX_LOG(mdlog.warning)
233 .appendTextFormatted(
234 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
235 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
236 devFlags.enableGpuBufferOps = true;
238 GMX_LOG(mdlog.warning)
240 .appendTextFormatted(
241 "This run uses the 'GPU halo exchange' feature, enabled by the "
242 "GMX_GPU_DD_COMMS environment variable.");
246 GMX_LOG(mdlog.warning)
248 .appendTextFormatted(
249 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
250 "halo exchange' feature will not be enabled as nonbonded interactions "
251 "are not offloaded.");
252 devFlags.enableGpuHaloExchange = false;
256 if (devFlags.enableGpuPmePPComm)
258 if (pmeRunMode == PmeRunMode::GPU)
260 GMX_LOG(mdlog.warning)
262 .appendTextFormatted(
263 "This run uses the 'GPU PME-PP communications' feature, enabled "
264 "by the GMX_GPU_PME_PP_COMMS environment variable.");
268 std::string clarification;
269 if (pmeRunMode == PmeRunMode::Mixed)
272 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
277 clarification = "PME is not offloaded to the GPU.";
279 GMX_LOG(mdlog.warning)
282 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
283 "'GPU PME-PP communications' feature was not enabled as "
285 devFlags.enableGpuPmePPComm = false;
292 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
294 * Used to ensure that the master thread does not modify mdrunner during copy
295 * on the spawned threads. */
296 static void threadMpiMdrunnerAccessBarrier()
299 MPI_Barrier(MPI_COMM_WORLD);
303 Mdrunner Mdrunner::cloneOnSpawnedThread() const
305 auto newRunner = Mdrunner(std::make_unique<MDModules>());
307 // All runners in the same process share a restraint manager resource because it is
308 // part of the interface to the client code, which is associated only with the
309 // original thread. Handles to the same resources can be obtained by copy.
311 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
314 // Copy members of master runner.
315 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
316 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
317 newRunner.hw_opt = hw_opt;
318 newRunner.filenames = filenames;
320 newRunner.oenv = oenv;
321 newRunner.mdrunOptions = mdrunOptions;
322 newRunner.domdecOptions = domdecOptions;
323 newRunner.nbpu_opt = nbpu_opt;
324 newRunner.pme_opt = pme_opt;
325 newRunner.pme_fft_opt = pme_fft_opt;
326 newRunner.bonded_opt = bonded_opt;
327 newRunner.update_opt = update_opt;
328 newRunner.nstlist_cmdline = nstlist_cmdline;
329 newRunner.replExParams = replExParams;
330 newRunner.pforce = pforce;
331 // Give the spawned thread the newly created valid communicator
332 // for the simulation.
333 newRunner.communicator = MPI_COMM_WORLD;
335 newRunner.startingBehavior = startingBehavior;
336 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
338 threadMpiMdrunnerAccessBarrier();
343 /*! \brief The callback used for running on spawned threads.
345 * Obtains the pointer to the master mdrunner object from the one
346 * argument permitted to the thread-launch API call, copies it to make
347 * a new runner for this thread, reinitializes necessary data, and
348 * proceeds to the simulation. */
349 static void mdrunner_start_fn(const void* arg)
353 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
354 /* copy the arg list to make sure that it's thread-local. This
355 doesn't copy pointed-to items, of course; fnm, cr and fplog
356 are reset in the call below, all others should be const. */
357 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
360 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
364 void Mdrunner::spawnThreads(int numThreadsToLaunch)
367 /* now spawn new threads that start mdrunner_start_fn(), while
368 the main thread returns. Thread affinity is handled later. */
369 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
370 static_cast<const void*>(this))
373 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
376 // Give the master thread the newly created valid communicator for
378 communicator = MPI_COMM_WORLD;
379 threadMpiMdrunnerAccessBarrier();
381 GMX_UNUSED_VALUE(numThreadsToLaunch);
382 GMX_UNUSED_VALUE(mdrunner_start_fn);
388 /*! \brief Initialize variables for Verlet scheme simulation */
389 static void prepare_verlet_scheme(FILE* fplog,
393 const gmx_mtop_t* mtop,
395 bool makeGpuPairList,
396 const gmx::CpuInfo& cpuinfo)
398 /* For NVE simulations, we will retain the initial list buffer */
399 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
401 /* Update the Verlet buffer size for the current run setup */
403 /* Here we assume SIMD-enabled kernels are being used. But as currently
404 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
405 * and 4x2 gives a larger buffer than 4x4, this is ok.
407 ListSetupType listType =
408 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
409 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
411 const real rlist_new =
412 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
414 if (rlist_new != ir->rlist)
416 if (fplog != nullptr)
419 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
420 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
422 ir->rlist = rlist_new;
426 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
428 gmx_fatal(FARGS, "Can not set nstlist without %s",
429 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
432 if (EI_DYNAMICS(ir->eI))
434 /* Set or try nstlist values */
435 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
439 /*! \brief Override the nslist value in inputrec
441 * with value passed on the command line (if any)
443 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
447 /* override with anything else than the default -2 */
448 if (nsteps_cmdline > -2)
450 char sbuf_steps[STEPSTRSIZE];
451 char sbuf_msg[STRLEN];
453 ir->nsteps = nsteps_cmdline;
454 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
457 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
458 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
462 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
463 gmx_step_str(nsteps_cmdline, sbuf_steps));
466 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
468 else if (nsteps_cmdline < -2)
470 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
472 /* Do nothing if nsteps_cmdline == -2 */
478 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
480 * If not, and if a warning may be issued, logs a warning about
481 * falling back to CPU code. With thread-MPI, only the first
482 * call to this function should have \c issueWarning true. */
483 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
485 bool gpuIsUseful = true;
488 if (ir.opts.ngener - ir.nwall > 1)
490 /* The GPU code does not support more than one energy group.
491 * If the user requested GPUs explicitly, a fatal error is given later.
495 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
496 "For better performance, run on the GPU without energy groups and then do "
497 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
503 warning = "TPI is not implemented for GPUs.";
506 if (!gpuIsUseful && issueWarning)
508 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
514 //! Initializes the logger for mdrun.
515 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
517 gmx::LoggerBuilder builder;
518 if (fplog != nullptr)
520 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
522 if (isSimulationMasterRank)
524 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
526 return builder.build();
529 //! Make a TaskTarget from an mdrun argument string.
530 static TaskTarget findTaskTarget(const char* optionString)
532 TaskTarget returnValue = TaskTarget::Auto;
534 if (strncmp(optionString, "auto", 3) == 0)
536 returnValue = TaskTarget::Auto;
538 else if (strncmp(optionString, "cpu", 3) == 0)
540 returnValue = TaskTarget::Cpu;
542 else if (strncmp(optionString, "gpu", 3) == 0)
544 returnValue = TaskTarget::Gpu;
548 GMX_ASSERT(false, "Option string should have been checked for sanity already");
554 //! Finish run, aggregate data to print performance info.
555 static void finish_run(FILE* fplog,
556 const gmx::MDLogger& mdlog,
558 const t_inputrec* inputrec,
560 gmx_wallcycle_t wcycle,
561 gmx_walltime_accounting_t walltime_accounting,
562 nonbonded_verlet_t* nbv,
563 const gmx_pme_t* pme,
567 double nbfs = 0, mflop = 0;
568 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
569 elapsed_time_over_all_threads_over_all_ranks;
570 /* Control whether it is valid to print a report. Only the
571 simulation master may print, but it should not do so if the run
572 terminated e.g. before a scheduled reset step. This is
573 complicated by the fact that PME ranks are unaware of the
574 reason why they were sent a pmerecvqxFINISH. To avoid
575 communication deadlocks, we always do the communication for the
576 report, even if we've decided not to write the report, because
577 how long it takes to finish the run is not important when we've
578 decided not to report on the simulation performance.
580 Further, we only report performance for dynamical integrators,
581 because those are the only ones for which we plan to
582 consider doing any optimizations. */
583 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
585 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
587 GMX_LOG(mdlog.warning)
589 .appendText("Simulation ended prematurely, no performance report will be written.");
594 std::unique_ptr<t_nrnb> nrnbTotalStorage;
597 nrnbTotalStorage = std::make_unique<t_nrnb>();
598 nrnb_tot = nrnbTotalStorage.get();
600 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
608 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
609 elapsed_time_over_all_threads =
610 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
614 /* reduce elapsed_time over all MPI ranks in the current simulation */
615 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
617 elapsed_time_over_all_ranks /= cr->nnodes;
618 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
619 * current simulation. */
620 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
621 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
626 elapsed_time_over_all_ranks = elapsed_time;
627 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
632 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
635 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
637 print_dd_statistics(cr, inputrec, fplog);
640 /* TODO Move the responsibility for any scaling by thread counts
641 * to the code that handled the thread region, so that there's a
642 * mechanism to keep cycle counting working during the transition
643 * to task parallelism. */
644 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
645 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
646 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
647 nthreads_pp, nthreads_pme);
648 auto cycle_sum(wallcycle_sum(cr, wcycle));
652 auto nbnxn_gpu_timings =
653 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
654 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
656 if (pme_gpu_task_enabled(pme))
658 pme_gpu_get_timings(pme, &pme_gpu_timings);
660 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
661 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
664 if (EI_DYNAMICS(inputrec->eI))
666 delta_t = inputrec->delta_t;
671 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
672 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
673 delta_t, nbfs, mflop);
677 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
678 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
679 delta_t, nbfs, mflop);
684 int Mdrunner::mdrunner()
687 t_forcerec* fr = nullptr;
688 t_fcdata* fcd = nullptr;
689 real ewaldcoeff_q = 0;
690 real ewaldcoeff_lj = 0;
691 int nChargePerturbed = -1, nTypePerturbed = 0;
692 gmx_wallcycle_t wcycle;
693 gmx_walltime_accounting_t walltime_accounting = nullptr;
694 gmx_membed_t* membed = nullptr;
695 gmx_hw_info_t* hwinfo = nullptr;
697 /* CAUTION: threads may be started later on in this function, so
698 cr doesn't reflect the final parallel state right now */
701 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
702 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
703 const bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
704 const bool doRerun = mdrunOptions.rerun;
706 // Handle task-assignment related user options.
707 EmulateGpuNonbonded emulateGpuNonbonded =
708 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
710 std::vector<int> userGpuTaskAssignment;
713 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
715 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
716 auto nonbondedTarget = findTaskTarget(nbpu_opt);
717 auto pmeTarget = findTaskTarget(pme_opt);
718 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
719 auto bondedTarget = findTaskTarget(bonded_opt);
720 auto updateTarget = findTaskTarget(update_opt);
722 FILE* fplog = nullptr;
723 // If we are appending, we don't write log output because we need
724 // to check that the old log file matches what the checkpoint file
725 // expects. Otherwise, we should start to write log output now if
726 // there is a file ready for it.
727 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
729 fplog = gmx_fio_getfp(logFileHandle);
731 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
732 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
733 gmx::MDLogger mdlog(logOwner.logger());
735 // TODO The thread-MPI master rank makes a working
736 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
737 // after the threads have been launched. This works because no use
738 // is made of that communicator until after the execution paths
739 // have rejoined. But it is likely that we can improve the way
740 // this is expressed, e.g. by expressly running detection only the
741 // master rank for thread-MPI, rather than relying on the mutex
742 // and reference count.
743 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
744 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
746 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
748 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
750 // Print citation requests after all software/hardware printing
751 pleaseCiteGromacs(fplog);
753 // TODO Replace this by unique_ptr once t_inputrec is C++
754 t_inputrec inputrecInstance;
755 t_inputrec* inputrec = nullptr;
756 std::unique_ptr<t_state> globalState;
758 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
760 if (isSimulationMasterRank)
762 /* Only the master rank has the global state */
763 globalState = std::make_unique<t_state>();
765 /* Read (nearly) all data required for the simulation
766 * and keep the partly serialized tpr contents to send to other ranks later
768 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
769 &inputrecInstance, globalState.get(), &mtop);
770 inputrec = &inputrecInstance;
773 /* Check and update the hardware options for internal consistency */
774 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
777 if (GMX_THREAD_MPI && isSimulationMasterRank)
779 bool useGpuForNonbonded = false;
780 bool useGpuForPme = false;
783 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
785 // If the user specified the number of ranks, then we must
786 // respect that, but in default mode, we need to allow for
787 // the number of GPUs to choose the number of ranks.
788 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
789 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
790 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
791 canUseGpuForNonbonded,
792 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
793 hw_opt.nthreads_tmpi);
794 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
795 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
796 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
798 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
800 /* Determine how many thread-MPI ranks to start.
802 * TODO Over-writing the user-supplied value here does
803 * prevent any possible subsequent checks from working
805 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
806 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
808 // Now start the threads for thread MPI.
809 spawnThreads(hw_opt.nthreads_tmpi);
810 // The spawned threads enter mdrunner() and execution of
811 // master and spawned threads joins at the end of this block.
812 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
815 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
816 CommrecHandle crHandle = init_commrec(communicator, ms);
817 t_commrec* cr = crHandle.get();
818 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
822 /* now broadcast everything to the non-master nodes/threads: */
823 if (!isSimulationMasterRank)
825 inputrec = &inputrecInstance;
827 init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
829 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
830 partialDeserializedTpr.reset(nullptr);
832 // Now the number of ranks is known to all ranks, and each knows
833 // the inputrec read by the master rank. The ranks can now all run
834 // the task-deciding functions and will agree on the result
835 // without needing to communicate.
837 // TODO Should we do the communication in debug mode to support
838 // having an assertion?
839 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
841 // Note that these variables describe only their own node.
843 // Note that when bonded interactions run on a GPU they always run
844 // alongside a nonbonded task, so do not influence task assignment
845 // even though they affect the force calculation workload.
846 bool useGpuForNonbonded = false;
847 bool useGpuForPme = false;
848 bool useGpuForBonded = false;
849 bool useGpuForUpdate = false;
850 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
853 // It's possible that there are different numbers of GPUs on
854 // different nodes, which is the user's responsibility to
855 // handle. If unsuitable, we will notice that during task
857 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
858 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
859 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
860 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
861 useGpuForPme = decideWhetherToUseGpusForPme(
862 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
863 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
864 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
865 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
866 useGpuForBonded = decideWhetherToUseGpusForBonded(
867 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
868 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
869 domdecOptions.numPmeRanks, gpusWereDetected);
871 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
873 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
875 // Initialize development feature flags that enabled by environment variable
876 // and report those features that are enabled.
877 const DevelopmentFeatureFlags devFlags =
878 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
880 const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
881 false, inputrec, doRerun, mtop, ms, replExParams, nullptr, doEssentialDynamics, doMembed);
882 const bool useModularSimulator = inputIsCompatibleWithModularSimulator
883 && !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
886 // TODO: hide restraint implementation details from Mdrunner.
887 // There is nothing unique about restraints at this point as far as the
888 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
889 // factory functions from the SimulationContext on which to call mdModules_->add().
890 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
891 for (auto&& restraint : restraintManager_->getRestraints())
893 auto module = RestraintMDModule::create(restraint, restraint->sites());
894 mdModules_->add(std::move(module));
897 // TODO: Error handling
898 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
899 const auto& mdModulesNotifier = mdModules_->notifier().notifier_;
901 if (inputrec->internalParameters != nullptr)
903 mdModulesNotifier.notify(*inputrec->internalParameters);
906 if (fplog != nullptr)
908 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
909 fprintf(fplog, "\n");
914 /* In rerun, set velocities to zero if present */
915 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
917 // rerun does not use velocities
921 "Rerun trajectory contains velocities. Rerun does only evaluate "
922 "potential energy and forces. The velocities will be ignored.");
923 for (int i = 0; i < globalState->natoms; i++)
925 clear_rvec(globalState->v[i]);
927 globalState->flags &= ~(1 << estV);
930 /* now make sure the state is initialized and propagated */
931 set_state_entries(globalState.get(), inputrec, useModularSimulator);
934 /* NM and TPI parallelize over force/energy calculations, not atoms,
935 * so we need to initialize and broadcast the global state.
937 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
941 globalState = std::make_unique<t_state>();
943 broadcastStateWithoutDynamics(cr, globalState.get());
946 /* A parallel command line option consistency check that we can
947 only do after any threads have started. */
949 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
950 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
953 "The -dd or -npme option request a parallel simulation, "
955 "but %s was compiled without threads or MPI enabled",
956 output_env_get_program_display_name(oenv));
958 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
960 "but %s was not started through mpirun/mpiexec or only one rank was requested "
961 "through mpirun/mpiexec",
962 output_env_get_program_display_name(oenv));
966 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
969 "The .mdp file specified an energy mininization or normal mode algorithm, and "
970 "these are not compatible with mdrun -rerun");
973 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
975 if (domdecOptions.numPmeRanks > 0)
977 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
978 "PME-only ranks are requested, but the system does not use PME "
979 "for electrostatics or LJ");
982 domdecOptions.numPmeRanks = 0;
985 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
987 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
988 * improve performance with many threads per GPU, since our OpenMP
989 * scaling is bad, but it's difficult to automate the setup.
991 domdecOptions.numPmeRanks = 0;
995 if (domdecOptions.numPmeRanks < 0)
997 domdecOptions.numPmeRanks = 0;
998 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1002 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1003 "PME GPU decomposition is not supported");
1010 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1014 /* NMR restraints must be initialized before load_checkpoint,
1015 * since with time averaging the history is added to t_state.
1016 * For proper consistency check we therefore need to extend
1018 * So the PME-only nodes (if present) will also initialize
1019 * the distance restraints.
1023 /* This needs to be called before read_checkpoint to extend the state */
1024 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1026 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1028 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1030 ObservablesHistory observablesHistory = {};
1032 if (startingBehavior != StartingBehavior::NewSimulation)
1034 /* Check if checkpoint file exists before doing continuation.
1035 * This way we can use identical input options for the first and subsequent runs...
1037 if (mdrunOptions.numStepsCommandline > -2)
1039 /* Temporarily set the number of steps to unlimited to avoid
1040 * triggering the nsteps check in load_checkpoint().
1041 * This hack will go away soon when the -nsteps option is removed.
1043 inputrec->nsteps = -1;
1046 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1047 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1048 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1050 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1052 // Now we can start normal logging to the truncated log file.
1053 fplog = gmx_fio_getfp(logFileHandle);
1054 prepareLogAppending(fplog);
1055 logOwner = buildLogger(fplog, MASTER(cr));
1056 mdlog = logOwner.logger();
1060 if (mdrunOptions.numStepsCommandline > -2)
1065 "The -nsteps functionality is deprecated, and may be removed in a future "
1067 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1070 /* override nsteps with value set on the commandline */
1071 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1075 copy_mat(globalState->box, box);
1080 gmx_bcast(sizeof(box), box, cr);
1083 if (inputrec->cutoff_scheme != ecutsVERLET)
1086 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1087 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1089 /* Update rlist and nstlist. */
1090 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1091 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1094 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1095 // This builder is necessary while we have multi-part construction
1096 // of DD. Before DD is constructed, we use the existence of
1097 // the builder object to indicate that further construction of DD
1099 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1100 if (useDomainDecomposition)
1102 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1103 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1104 positionsFromStatePointer(globalState.get()));
1108 /* PME, if used, is done on all nodes with 1D decomposition */
1110 cr->duty = (DUTY_PP | DUTY_PME);
1112 if (inputrec->ePBC == epbcSCREW)
1114 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1118 // Produce the task assignment for this rank.
1119 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1120 GpuTaskAssignments gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1121 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1122 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1123 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1124 // TODO cr->duty & DUTY_PME should imply that a PME
1125 // algorithm is active, but currently does not.
1126 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1128 // Get the device handles for the modules, nullptr when no task is assigned.
1129 gmx_device_info_t* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1130 gmx_device_info_t* pmeDeviceInfo = gpuTaskAssignments.initPmeDevice();
1132 // TODO Initialize GPU streams here.
1134 // TODO Currently this is always built, yet DD partition code
1135 // checks if it is built before using it. Probably it should
1136 // become an MDModule that is made only when another module
1137 // requires it (e.g. pull, CompEl, density fitting), so that we
1138 // don't update the local atom sets unilaterally every step.
1139 LocalAtomSetManager atomSets;
1142 // TODO Pass the GPU streams to ddBuilder to use in buffer
1143 // transfers (e.g. halo exchange)
1144 cr->dd = ddBuilder->build(&atomSets);
1145 // The builder's job is done, so destruct it
1146 ddBuilder.reset(nullptr);
1147 // Note that local state still does not exist yet.
1150 // The GPU update is decided here because we need to know whether the constraints or
1151 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1152 // defined). This is only known after DD is initialized, hence decision on using GPU
1153 // update is done so late.
1156 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1158 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1159 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1160 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1161 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1162 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1164 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1166 const bool printHostName = (cr->nnodes > 1);
1167 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1169 // If the user chose a task assignment, give them some hints
1170 // where appropriate.
1171 if (!userGpuTaskAssignment.empty())
1173 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1178 /* After possible communicator splitting in make_dd_communicators.
1179 * we can set up the intra/inter node communication.
1181 gmx_setup_nodecomm(fplog, cr);
1187 GMX_LOG(mdlog.warning)
1189 .appendTextFormatted(
1190 "This is simulation %d out of %d running as a composite GROMACS\n"
1191 "multi-simulation job. Setup for this simulation:\n",
1194 GMX_LOG(mdlog.warning)
1195 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1197 cr->nnodes == 1 ? "thread" : "threads"
1199 cr->nnodes == 1 ? "process" : "processes"
1205 // If mdrun -pin auto honors any affinity setting that already
1206 // exists. If so, it is nice to provide feedback about whether
1207 // that existing affinity setting was from OpenMP or something
1208 // else, so we run this code both before and after we initialize
1209 // the OpenMP support.
1210 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1211 /* Check and update the number of OpenMP threads requested */
1212 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1213 pmeRunMode, mtop, *inputrec);
1215 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1216 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1218 // Enable FP exception detection, but not in
1219 // Release mode and not for compilers with known buggy FP
1220 // exception support (clang with any optimization) or suspected
1221 // buggy FP exception support (gcc 7.* with optimization).
1222 #if !defined NDEBUG \
1223 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1224 && defined __OPTIMIZE__)
1225 const bool bEnableFPE = true;
1227 const bool bEnableFPE = false;
1229 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1232 gmx_feenableexcept();
1235 /* Now that we know the setup is consistent, check for efficiency */
1236 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1237 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1239 /* getting number of PP/PME threads on this MPI / tMPI rank.
1240 PME: env variable should be read only on one node to make sure it is
1241 identical everywhere;
1243 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1244 : gmx_omp_nthreads_get(emntPME);
1245 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1246 physicalNodeComm, mdlog);
1248 // Enable Peer access between GPUs where available
1249 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1250 // any of the GPU communication features are active.
1251 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1252 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1254 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1257 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1259 /* Before setting affinity, check whether the affinity has changed
1260 * - which indicates that probably the OpenMP library has changed it
1261 * since we first checked).
1263 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1265 int numThreadsOnThisNode, intraNodeThreadOffset;
1266 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1267 &intraNodeThreadOffset);
1269 /* Set the CPU affinity */
1270 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1271 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1274 if (mdrunOptions.timingOptions.resetStep > -1)
1279 "The -resetstep functionality is deprecated, and may be removed in a "
1282 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1286 /* Master synchronizes its value of reset_counters with all nodes
1287 * including PME only nodes */
1288 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1289 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1290 wcycle_set_reset_counters(wcycle, reset_counters);
1293 // Membrane embedding must be initialized before we call init_forcerec()
1298 fprintf(stderr, "Initializing membed");
1300 /* Note that membed cannot work in parallel because mtop is
1301 * changed here. Fix this if we ever want to make it run with
1302 * multiple ranks. */
1303 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1304 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1307 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1308 std::unique_ptr<MDAtoms> mdAtoms;
1309 std::unique_ptr<gmx_vsite_t> vsite;
1312 if (thisRankHasDuty(cr, DUTY_PP))
1314 mdModulesNotifier.notify(*cr);
1315 mdModulesNotifier.notify(&atomSets);
1316 mdModulesNotifier.notify(PeriodicBoundaryConditionType{ inputrec->ePBC });
1317 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1318 /* Initiate forcerecord */
1319 fr = new t_forcerec;
1320 fr->forceProviders = mdModules_->initForceProviders();
1321 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1322 opt2fn("-table", filenames.size(), filenames.data()),
1323 opt2fn("-tablep", filenames.size(), filenames.data()),
1324 opt2fns("-tableb", filenames.size(), filenames.data()), *hwinfo,
1325 nonbondedDeviceInfo, useGpuForBonded,
1326 pmeRunMode == PmeRunMode::GPU && !thisRankHasDuty(cr, DUTY_PME), pforce, wcycle);
1328 // TODO Move this to happen during domain decomposition setup,
1329 // once stream and event handling works well with that.
1330 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
1331 if (havePPDomainDecomposition(cr) && prefer1DAnd1PulseDD && is1DAnd1PulseDD(*cr->dd))
1333 GMX_RELEASE_ASSERT(devFlags.enableGpuBufferOps,
1334 "Must use GMX_USE_GPU_BUFFER_OPS=1 to use GMX_GPU_DD_COMMS=1");
1336 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local);
1337 void* streamNonLocal =
1338 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal);
1339 GMX_LOG(mdlog.warning)
1341 .appendTextFormatted(
1342 "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the "
1343 "GMX_GPU_DD_COMMS environment variable.");
1344 cr->dd->gpuHaloExchange = std::make_unique<GpuHaloExchange>(
1345 cr->dd, cr->mpi_comm_mysim, streamLocal, streamNonLocal);
1348 /* Initialize the mdAtoms structure.
1349 * mdAtoms is not filled with atom data,
1350 * as this can not be done now with domain decomposition.
1352 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1353 if (globalState && thisRankHasPmeGpuTask)
1355 // The pinning of coordinates in the global state object works, because we only use
1356 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1357 // points to the global state object without DD.
1358 // FIXME: MD and EM separately set up the local state - this should happen in the same
1359 // function, which should also perform the pinning.
1360 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1363 /* Initialize the virtual site communication */
1364 vsite = initVsite(mtop, cr);
1366 calc_shifts(box, fr->shift_vec);
1368 /* With periodic molecules the charge groups should be whole at start up
1369 * and the virtual sites should not be far from their proper positions.
1371 if (!inputrec->bContinuation && MASTER(cr) && !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1373 /* Make molecules whole at start of run */
1374 if (fr->ePBC != epbcNONE)
1376 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1380 /* Correct initial vsite positions are required
1381 * for the initial distribution in the domain decomposition
1382 * and for the initial shell prediction.
1384 constructVsitesGlobal(mtop, globalState->x);
1388 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1390 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1391 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1396 /* This is a PME only node */
1398 GMX_ASSERT(globalState == nullptr,
1399 "We don't need the state on a PME only rank and expect it to be unitialized");
1401 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1402 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1405 gmx_pme_t* sepPmeData = nullptr;
1406 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1407 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1408 "Double-checking that only PME-only ranks have no forcerec");
1409 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1411 // TODO should live in ewald module once its testing is improved
1413 // Later, this program could contain kernels that might be later
1414 // re-used as auto-tuning progresses, or subsequent simulations
1416 PmeGpuProgramStorage pmeGpuProgram;
1417 if (thisRankHasPmeGpuTask)
1419 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1422 /* Initiate PME if necessary,
1423 * either on all nodes or on dedicated PME nodes only. */
1424 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1426 if (mdAtoms && mdAtoms->mdatoms())
1428 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1429 if (EVDW_PME(inputrec->vdwtype))
1431 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1434 if (cr->npmenodes > 0)
1436 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1437 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1438 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1441 if (thisRankHasDuty(cr, DUTY_PME))
1445 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1446 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1447 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1448 nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1450 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1455 if (EI_DYNAMICS(inputrec->eI))
1457 /* Turn on signal handling on all nodes */
1459 * (A user signal from the PME nodes (if any)
1460 * is communicated to the PP nodes.
1462 signal_handler_install();
1465 pull_t* pull_work = nullptr;
1466 if (thisRankHasDuty(cr, DUTY_PP))
1468 /* Assumes uniform use of the number of OpenMP threads */
1469 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1471 if (inputrec->bPull)
1473 /* Initialize pull code */
1474 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1475 inputrec->fepvals->init_lambda);
1476 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1478 initPullHistory(pull_work, &observablesHistory);
1480 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1482 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1486 std::unique_ptr<EnforcedRotation> enforcedRotation;
1489 /* Initialize enforced rotation code */
1491 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1492 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1495 t_swap* swap = nullptr;
1496 if (inputrec->eSwapCoords != eswapNO)
1498 /* Initialize ion swapping code */
1499 swap = init_swapcoords(fplog, inputrec,
1500 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1501 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1502 oenv, mdrunOptions, startingBehavior);
1505 /* Let makeConstraints know whether we have essential dynamics constraints. */
1506 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1507 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1509 /* Energy terms and groups */
1510 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1511 inputrec->fepvals->n_lambda);
1513 /* Kinetic energy data */
1514 gmx_ekindata_t ekind;
1515 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1517 /* Set up interactive MD (IMD) */
1519 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1520 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1521 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1523 if (DOMAINDECOMP(cr))
1525 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1526 /* This call is not included in init_domain_decomposition mainly
1527 * because fr->cginfo_mb is set later.
1529 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1530 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1533 // TODO This is not the right place to manage the lifetime of
1534 // this data structure, but currently it's the easiest way to
1536 MdrunScheduleWorkload runScheduleWork;
1537 // Also populates the simulation constant workload description.
1538 runScheduleWork.simulationWork = createSimulationWorkload(
1539 useGpuForNonbonded, pmeRunMode, useGpuForBonded, useGpuForUpdate,
1540 devFlags.enableGpuBufferOps, devFlags.enableGpuHaloExchange,
1541 devFlags.enableGpuPmePPComm, haveEwaldSurfaceContribution(*inputrec));
1543 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1544 if (gpusWereDetected
1545 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1546 || runScheduleWork.simulationWork.useGpuBufferOps))
1548 const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1549 const void* localStream =
1550 fr->nbv->gpu_nbv != nullptr
1551 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1553 const void* nonLocalStream =
1554 fr->nbv->gpu_nbv != nullptr
1555 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1557 const void* deviceContext = pme_gpu_get_device_context(fr->pmedata);
1558 const int paddingSize = pme_gpu_get_padding_size(fr->pmedata);
1559 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1560 ? GpuApiCallBehavior::Async
1561 : GpuApiCallBehavior::Sync;
1563 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1564 pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
1565 fr->stateGpu = stateGpu.get();
1568 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1569 SimulatorBuilder simulatorBuilder;
1571 // build and run simulator object based on user-input
1572 auto simulator = simulatorBuilder.build(
1573 inputIsCompatibleWithModularSimulator, fplog, cr, ms, mdlog,
1574 static_cast<int>(filenames.size()), filenames.data(), oenv, mdrunOptions,
1575 startingBehavior, vsite.get(), constr.get(),
1576 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1577 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1578 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1579 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1580 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1583 if (fr->pmePpCommGpu)
1585 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1586 fr->pmePpCommGpu.reset();
1589 if (inputrec->bPull)
1591 finish_pull(pull_work);
1593 finish_swapcoords(swap);
1597 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1599 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1600 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1603 wallcycle_stop(wcycle, ewcRUN);
1605 /* Finish up, write some stuff
1606 * if rerunMD, don't write last frame again
1608 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1609 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1611 // clean up cycle counter
1612 wallcycle_destroy(wcycle);
1617 gmx_pme_destroy(pmedata);
1621 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1622 // before we destroy the GPU context(s) in free_gpu_resources().
1623 // Pinned buffers are associated with contexts in CUDA.
1624 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1625 mdAtoms.reset(nullptr);
1626 globalState.reset(nullptr);
1627 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1629 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1630 free_gpu_resources(fr, physicalNodeComm, hwinfo->gpu_info);
1631 free_gpu(nonbondedDeviceInfo);
1632 free_gpu(pmeDeviceInfo);
1633 done_forcerec(fr, mtop.molblock.size());
1638 free_membed(membed);
1641 /* Does what it says */
1642 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1643 walltime_accounting_destroy(walltime_accounting);
1645 // Ensure log file content is written
1648 gmx_fio_flush(logFileHandle);
1651 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1652 * exceptions were enabled before function was called. */
1655 gmx_fedisableexcept();
1658 auto rc = static_cast<int>(gmx_get_stop_condition());
1661 /* we need to join all threads. The sub-threads join when they
1662 exit this function, but the master thread needs to be told to
1664 if (PAR(cr) && MASTER(cr))
1672 Mdrunner::~Mdrunner()
1674 // Clean up of the Manager.
1675 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1676 // but okay as long as threads synchronize some time before adding or accessing
1677 // a new set of restraints.
1678 if (restraintManager_)
1680 restraintManager_->clear();
1681 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1682 "restraints added during runner life time should be cleared at runner "
1687 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1689 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1690 // Not sure if this should be logged through the md logger or something else,
1691 // but it is helpful to have some sort of INFO level message sent somewhere.
1692 // std::cout << "Registering restraint named " << name << std::endl;
1694 // When multiple restraints are used, it may be wasteful to register them separately.
1695 // Maybe instead register an entire Restraint Manager as a force provider.
1696 restraintManager_->addToSpec(std::move(puller), name);
1699 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1701 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1703 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1704 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1706 class Mdrunner::BuilderImplementation
1709 BuilderImplementation() = delete;
1710 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1711 ~BuilderImplementation();
1713 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1714 real forceWarningThreshold,
1715 StartingBehavior startingBehavior);
1717 void addDomdec(const DomdecOptions& options);
1719 void addVerletList(int nstlist);
1721 void addReplicaExchange(const ReplicaExchangeParameters& params);
1723 void addNonBonded(const char* nbpu_opt);
1725 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1727 void addBondedTaskAssignment(const char* bonded_opt);
1729 void addUpdateTaskAssignment(const char* update_opt);
1731 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1733 void addFilenames(ArrayRef<const t_filenm> filenames);
1735 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1737 void addLogFile(t_fileio* logFileHandle);
1739 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1744 // Default parameters copied from runner.h
1745 // \todo Clarify source(s) of default parameters.
1747 const char* nbpu_opt_ = nullptr;
1748 const char* pme_opt_ = nullptr;
1749 const char* pme_fft_opt_ = nullptr;
1750 const char* bonded_opt_ = nullptr;
1751 const char* update_opt_ = nullptr;
1753 MdrunOptions mdrunOptions_;
1755 DomdecOptions domdecOptions_;
1757 ReplicaExchangeParameters replicaExchangeParameters_;
1759 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1762 //! Multisim communicator handle.
1763 gmx_multisim_t* multiSimulation_;
1765 //! mdrun communicator
1766 MPI_Comm communicator_ = MPI_COMM_NULL;
1768 //! Print a warning if any force is larger than this (in kJ/mol nm).
1769 real forceWarningThreshold_ = -1;
1771 //! Whether the simulation will start afresh, or restart with/without appending.
1772 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1774 //! The modules that comprise the functionality of mdrun.
1775 std::unique_ptr<MDModules> mdModules_;
1777 //! \brief Parallelism information.
1778 gmx_hw_opt_t hardwareOptions_;
1780 //! filename options for simulation.
1781 ArrayRef<const t_filenm> filenames_;
1783 /*! \brief Handle to output environment.
1785 * \todo gmx_output_env_t needs lifetime management.
1787 gmx_output_env_t* outputEnvironment_ = nullptr;
1789 /*! \brief Non-owning handle to MD log file.
1791 * \todo Context should own output facilities for client.
1792 * \todo Improve log file handle management.
1794 * Code managing the FILE* relies on the ability to set it to
1795 * nullptr to check whether the filehandle is valid.
1797 t_fileio* logFileHandle_ = nullptr;
1800 * \brief Builder for simulation stop signal handler.
1802 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1805 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1806 compat::not_null<SimulationContext*> context) :
1807 mdModules_(std::move(mdModules))
1809 communicator_ = context->communicator_;
1810 multiSimulation_ = context->multiSimulation_.get();
1813 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1815 Mdrunner::BuilderImplementation&
1816 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1817 const real forceWarningThreshold,
1818 const StartingBehavior startingBehavior)
1820 mdrunOptions_ = options;
1821 forceWarningThreshold_ = forceWarningThreshold;
1822 startingBehavior_ = startingBehavior;
1826 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1828 domdecOptions_ = options;
1831 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1836 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1838 replicaExchangeParameters_ = params;
1841 Mdrunner Mdrunner::BuilderImplementation::build()
1843 auto newRunner = Mdrunner(std::move(mdModules_));
1845 newRunner.mdrunOptions = mdrunOptions_;
1846 newRunner.pforce = forceWarningThreshold_;
1847 newRunner.startingBehavior = startingBehavior_;
1848 newRunner.domdecOptions = domdecOptions_;
1850 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1851 newRunner.hw_opt = hardwareOptions_;
1853 // No invariant to check. This parameter exists to optionally override other behavior.
1854 newRunner.nstlist_cmdline = nstlist_;
1856 newRunner.replExParams = replicaExchangeParameters_;
1858 newRunner.filenames = filenames_;
1860 newRunner.communicator = communicator_;
1862 // nullptr is a valid value for the multisim handle
1863 newRunner.ms = multiSimulation_;
1865 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1866 // \todo Update sanity checking when output environment has clearly specified invariants.
1867 // Initialization and default values for oenv are not well specified in the current version.
1868 if (outputEnvironment_)
1870 newRunner.oenv = outputEnvironment_;
1874 GMX_THROW(gmx::APIError(
1875 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1878 newRunner.logFileHandle = logFileHandle_;
1882 newRunner.nbpu_opt = nbpu_opt_;
1886 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1889 if (pme_opt_ && pme_fft_opt_)
1891 newRunner.pme_opt = pme_opt_;
1892 newRunner.pme_fft_opt = pme_fft_opt_;
1896 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1901 newRunner.bonded_opt = bonded_opt_;
1905 GMX_THROW(gmx::APIError(
1906 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1911 newRunner.update_opt = update_opt_;
1915 GMX_THROW(gmx::APIError(
1916 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1920 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1922 if (stopHandlerBuilder_)
1924 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1928 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1934 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1936 nbpu_opt_ = nbpu_opt;
1939 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
1942 pme_fft_opt_ = pme_fft_opt;
1945 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1947 bonded_opt_ = bonded_opt;
1950 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
1952 update_opt_ = update_opt;
1955 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
1957 hardwareOptions_ = hardwareOptions;
1960 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1962 filenames_ = filenames;
1965 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1967 outputEnvironment_ = outputEnvironment;
1970 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
1972 logFileHandle_ = logFileHandle;
1975 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1977 stopHandlerBuilder_ = std::move(builder);
1980 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
1981 compat::not_null<SimulationContext*> context) :
1982 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
1986 MdrunnerBuilder::~MdrunnerBuilder() = default;
1988 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
1989 real forceWarningThreshold,
1990 const StartingBehavior startingBehavior)
1992 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
1996 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
1998 impl_->addDomdec(options);
2002 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2004 impl_->addVerletList(nstlist);
2008 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2010 impl_->addReplicaExchange(params);
2014 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2016 impl_->addNonBonded(nbpu_opt);
2020 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2022 // The builder method may become more general in the future, but in this version,
2023 // parameters for PME electrostatics are both required and the only parameters
2025 if (pme_opt && pme_fft_opt)
2027 impl_->addPME(pme_opt, pme_fft_opt);
2032 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2037 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2039 impl_->addBondedTaskAssignment(bonded_opt);
2043 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2045 impl_->addUpdateTaskAssignment(update_opt);
2049 Mdrunner MdrunnerBuilder::build()
2051 return impl_->build();
2054 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2056 impl_->addHardwareOptions(hardwareOptions);
2060 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2062 impl_->addFilenames(filenames);
2066 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2068 impl_->addOutputEnvironment(outputEnvironment);
2072 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2074 impl_->addLogFile(logFileHandle);
2078 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2080 impl_->addStopHandlerBuilder(std::move(builder));
2084 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2086 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;