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, 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"
170 /*! \brief Structure that holds boolean flags corresponding to the development
171 * features present enabled through environment variables.
174 struct DevelopmentFeatureFlags
176 //! True if the Buffer ops development feature is enabled
177 // TODO: when the trigger of the buffer ops offload is fully automated this should go away
178 bool enableGpuBufferOps = false;
179 //! If true, forces 'mdrun -update auto' default to 'gpu' when running with DD
180 bool forceGpuUpdateDefaultWithDD = false;
181 //! True if the GPU halo exchange development feature is enabled
182 bool enableGpuHaloExchange = false;
183 //! True if the PME PP direct communication GPU development feature is enabled
184 bool enableGpuPmePPComm = false;
187 /*! \brief Manage any development feature flag variables encountered
189 * The use of dev features indicated by environment variables is
190 * logged in order to ensure that runs with such features enabled can
191 * be identified from their log and standard output. Any cross
192 * dependencies are also checked, and if unsatisfied, a fatal error
195 * Note that some development features overrides are applied already here:
196 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
198 * \param[in] mdlog Logger object.
199 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
200 * \param[in] pmeRunMode The PME run mode for this run
201 * \returns The object populated with development feature flags.
203 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
204 const bool useGpuForNonbonded,
205 const PmeRunMode pmeRunMode)
207 DevelopmentFeatureFlags devFlags;
209 // Some builds of GCC 5 give false positive warnings that these
210 // getenv results are ignored when clearly they are used.
211 #pragma GCC diagnostic push
212 #pragma GCC diagnostic ignored "-Wunused-result"
213 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
214 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
215 devFlags.forceGpuUpdateDefaultWithDD = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
216 devFlags.enableGpuHaloExchange =
217 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
218 devFlags.enableGpuPmePPComm =
219 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
220 #pragma GCC diagnostic pop
222 if (devFlags.enableGpuBufferOps)
224 GMX_LOG(mdlog.warning)
226 .appendTextFormatted(
227 "NOTE: This run uses the 'GPU buffer ops' feature, enabled by the "
228 "GMX_USE_GPU_BUFFER_OPS environment variable.");
231 if (devFlags.forceGpuUpdateDefaultWithDD)
233 GMX_LOG(mdlog.warning)
235 .appendTextFormatted(
236 "NOTE: This run will default to '-update gpu' as requested by the "
237 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable.");
240 if (devFlags.enableGpuHaloExchange)
242 if (useGpuForNonbonded)
244 if (!devFlags.enableGpuBufferOps)
247 "Cannot enable GPU halo exchange without GPU buffer operations, set "
248 "GMX_USE_GPU_BUFFER_OPS=1\n");
250 GMX_LOG(mdlog.warning)
252 .appendTextFormatted(
253 "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the "
254 "GMX_GPU_DD_COMMS environment variable.");
258 GMX_LOG(mdlog.warning)
260 .appendTextFormatted(
261 "NOTE: GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
262 "halo exchange' feature will not be enabled as nonbonded interactions "
263 "are not offloaded.");
264 devFlags.enableGpuHaloExchange = false;
268 if (devFlags.enableGpuPmePPComm)
270 if (pmeRunMode == PmeRunMode::GPU)
272 GMX_LOG(mdlog.warning)
274 .appendTextFormatted(
275 "NOTE: This run uses the 'GPU PME-PP communications' feature, enabled "
276 "by the GMX_GPU_PME_PP_COMMS environment variable.");
280 std::string clarification;
281 if (pmeRunMode == PmeRunMode::Mixed)
284 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
289 clarification = "PME is not offloaded to the GPU.";
291 GMX_LOG(mdlog.warning)
294 "NOTE: GMX_GPU_PME_PP_COMMS environment variable detected, but the "
295 "'GPU PME-PP communications' feature was not enabled as "
297 devFlags.enableGpuPmePPComm = false;
304 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
306 * Used to ensure that the master thread does not modify mdrunner during copy
307 * on the spawned threads. */
308 static void threadMpiMdrunnerAccessBarrier()
311 MPI_Barrier(MPI_COMM_WORLD);
315 Mdrunner Mdrunner::cloneOnSpawnedThread() const
317 auto newRunner = Mdrunner(std::make_unique<MDModules>());
319 // All runners in the same process share a restraint manager resource because it is
320 // part of the interface to the client code, which is associated only with the
321 // original thread. Handles to the same resources can be obtained by copy.
323 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
326 // Copy members of master runner.
327 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
328 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
329 newRunner.hw_opt = hw_opt;
330 newRunner.filenames = filenames;
332 newRunner.oenv = oenv;
333 newRunner.mdrunOptions = mdrunOptions;
334 newRunner.domdecOptions = domdecOptions;
335 newRunner.nbpu_opt = nbpu_opt;
336 newRunner.pme_opt = pme_opt;
337 newRunner.pme_fft_opt = pme_fft_opt;
338 newRunner.bonded_opt = bonded_opt;
339 newRunner.update_opt = update_opt;
340 newRunner.nstlist_cmdline = nstlist_cmdline;
341 newRunner.replExParams = replExParams;
342 newRunner.pforce = pforce;
343 // Give the spawned thread the newly created valid communicator
344 // for the simulation.
345 newRunner.communicator = MPI_COMM_WORLD;
347 newRunner.startingBehavior = startingBehavior;
348 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
350 threadMpiMdrunnerAccessBarrier();
355 /*! \brief The callback used for running on spawned threads.
357 * Obtains the pointer to the master mdrunner object from the one
358 * argument permitted to the thread-launch API call, copies it to make
359 * a new runner for this thread, reinitializes necessary data, and
360 * proceeds to the simulation. */
361 static void mdrunner_start_fn(const void* arg)
365 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
366 /* copy the arg list to make sure that it's thread-local. This
367 doesn't copy pointed-to items, of course; fnm, cr and fplog
368 are reset in the call below, all others should be const. */
369 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
372 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
376 void Mdrunner::spawnThreads(int numThreadsToLaunch)
379 /* now spawn new threads that start mdrunner_start_fn(), while
380 the main thread returns. Thread affinity is handled later. */
381 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
382 static_cast<const void*>(this))
385 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
388 // Give the master thread the newly created valid communicator for
390 communicator = MPI_COMM_WORLD;
391 threadMpiMdrunnerAccessBarrier();
393 GMX_UNUSED_VALUE(numThreadsToLaunch);
394 GMX_UNUSED_VALUE(mdrunner_start_fn);
400 /*! \brief Initialize variables for Verlet scheme simulation */
401 static void prepare_verlet_scheme(FILE* fplog,
405 const gmx_mtop_t* mtop,
407 bool makeGpuPairList,
408 const gmx::CpuInfo& cpuinfo)
410 /* For NVE simulations, we will retain the initial list buffer */
411 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
413 /* Update the Verlet buffer size for the current run setup */
415 /* Here we assume SIMD-enabled kernels are being used. But as currently
416 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
417 * and 4x2 gives a larger buffer than 4x4, this is ok.
419 ListSetupType listType =
420 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
421 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
423 const real rlist_new =
424 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
426 if (rlist_new != ir->rlist)
428 if (fplog != nullptr)
431 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
432 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
434 ir->rlist = rlist_new;
438 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
440 gmx_fatal(FARGS, "Can not set nstlist without %s",
441 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
444 if (EI_DYNAMICS(ir->eI))
446 /* Set or try nstlist values */
447 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
451 /*! \brief Override the nslist value in inputrec
453 * with value passed on the command line (if any)
455 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
459 /* override with anything else than the default -2 */
460 if (nsteps_cmdline > -2)
462 char sbuf_steps[STEPSTRSIZE];
463 char sbuf_msg[STRLEN];
465 ir->nsteps = nsteps_cmdline;
466 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
469 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
470 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
474 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
475 gmx_step_str(nsteps_cmdline, sbuf_steps));
478 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
480 else if (nsteps_cmdline < -2)
482 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
484 /* Do nothing if nsteps_cmdline == -2 */
490 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
492 * If not, and if a warning may be issued, logs a warning about
493 * falling back to CPU code. With thread-MPI, only the first
494 * call to this function should have \c issueWarning true. */
495 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
497 bool gpuIsUseful = true;
500 if (ir.opts.ngener - ir.nwall > 1)
502 /* The GPU code does not support more than one energy group.
503 * If the user requested GPUs explicitly, a fatal error is given later.
507 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
508 "For better performance, run on the GPU without energy groups and then do "
509 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
515 warning = "TPI is not implemented for GPUs.";
518 if (!gpuIsUseful && issueWarning)
520 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
526 //! Initializes the logger for mdrun.
527 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
529 gmx::LoggerBuilder builder;
530 if (fplog != nullptr)
532 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
534 if (isSimulationMasterRank)
536 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
538 return builder.build();
541 //! Make a TaskTarget from an mdrun argument string.
542 static TaskTarget findTaskTarget(const char* optionString)
544 TaskTarget returnValue = TaskTarget::Auto;
546 if (strncmp(optionString, "auto", 3) == 0)
548 returnValue = TaskTarget::Auto;
550 else if (strncmp(optionString, "cpu", 3) == 0)
552 returnValue = TaskTarget::Cpu;
554 else if (strncmp(optionString, "gpu", 3) == 0)
556 returnValue = TaskTarget::Gpu;
560 GMX_ASSERT(false, "Option string should have been checked for sanity already");
566 //! Finish run, aggregate data to print performance info.
567 static void finish_run(FILE* fplog,
568 const gmx::MDLogger& mdlog,
570 const t_inputrec* inputrec,
572 gmx_wallcycle_t wcycle,
573 gmx_walltime_accounting_t walltime_accounting,
574 nonbonded_verlet_t* nbv,
575 const gmx_pme_t* pme,
579 double nbfs = 0, mflop = 0;
580 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
581 elapsed_time_over_all_threads_over_all_ranks;
582 /* Control whether it is valid to print a report. Only the
583 simulation master may print, but it should not do so if the run
584 terminated e.g. before a scheduled reset step. This is
585 complicated by the fact that PME ranks are unaware of the
586 reason why they were sent a pmerecvqxFINISH. To avoid
587 communication deadlocks, we always do the communication for the
588 report, even if we've decided not to write the report, because
589 how long it takes to finish the run is not important when we've
590 decided not to report on the simulation performance.
592 Further, we only report performance for dynamical integrators,
593 because those are the only ones for which we plan to
594 consider doing any optimizations. */
595 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
597 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
599 GMX_LOG(mdlog.warning)
601 .appendText("Simulation ended prematurely, no performance report will be written.");
606 std::unique_ptr<t_nrnb> nrnbTotalStorage;
609 nrnbTotalStorage = std::make_unique<t_nrnb>();
610 nrnb_tot = nrnbTotalStorage.get();
612 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
620 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
621 elapsed_time_over_all_threads =
622 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
626 /* reduce elapsed_time over all MPI ranks in the current simulation */
627 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
629 elapsed_time_over_all_ranks /= cr->nnodes;
630 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
631 * current simulation. */
632 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
633 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
638 elapsed_time_over_all_ranks = elapsed_time;
639 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
644 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
647 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
649 print_dd_statistics(cr, inputrec, fplog);
652 /* TODO Move the responsibility for any scaling by thread counts
653 * to the code that handled the thread region, so that there's a
654 * mechanism to keep cycle counting working during the transition
655 * to task parallelism. */
656 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
657 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
658 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
659 nthreads_pp, nthreads_pme);
660 auto cycle_sum(wallcycle_sum(cr, wcycle));
664 auto nbnxn_gpu_timings =
665 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
666 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
668 if (pme_gpu_task_enabled(pme))
670 pme_gpu_get_timings(pme, &pme_gpu_timings);
672 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
673 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
676 if (EI_DYNAMICS(inputrec->eI))
678 delta_t = inputrec->delta_t;
683 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
684 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
685 delta_t, nbfs, mflop);
689 print_perf(stderr, 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);
696 int Mdrunner::mdrunner()
699 t_forcerec* fr = nullptr;
700 t_fcdata* fcd = nullptr;
701 real ewaldcoeff_q = 0;
702 real ewaldcoeff_lj = 0;
703 int nChargePerturbed = -1, nTypePerturbed = 0;
704 gmx_wallcycle_t wcycle;
705 gmx_walltime_accounting_t walltime_accounting = nullptr;
706 gmx_membed_t* membed = nullptr;
707 gmx_hw_info_t* hwinfo = nullptr;
709 /* CAUTION: threads may be started later on in this function, so
710 cr doesn't reflect the final parallel state right now */
713 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
714 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
715 const bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
716 const bool doRerun = mdrunOptions.rerun;
718 // Handle task-assignment related user options.
719 EmulateGpuNonbonded emulateGpuNonbonded =
720 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
722 std::vector<int> userGpuTaskAssignment;
725 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
727 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
728 auto nonbondedTarget = findTaskTarget(nbpu_opt);
729 auto pmeTarget = findTaskTarget(pme_opt);
730 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
731 auto bondedTarget = findTaskTarget(bonded_opt);
732 auto updateTarget = findTaskTarget(update_opt);
733 PmeRunMode pmeRunMode = PmeRunMode::None;
735 FILE* fplog = nullptr;
736 // If we are appending, we don't write log output because we need
737 // to check that the old log file matches what the checkpoint file
738 // expects. Otherwise, we should start to write log output now if
739 // there is a file ready for it.
740 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
742 fplog = gmx_fio_getfp(logFileHandle);
744 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
745 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
746 gmx::MDLogger mdlog(logOwner.logger());
748 // TODO The thread-MPI master rank makes a working
749 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
750 // after the threads have been launched. This works because no use
751 // is made of that communicator until after the execution paths
752 // have rejoined. But it is likely that we can improve the way
753 // this is expressed, e.g. by expressly running detection only the
754 // master rank for thread-MPI, rather than relying on the mutex
755 // and reference count.
756 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
757 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
759 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
761 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
763 // Print citation requests after all software/hardware printing
764 pleaseCiteGromacs(fplog);
766 // TODO Replace this by unique_ptr once t_inputrec is C++
767 t_inputrec inputrecInstance;
768 t_inputrec* inputrec = nullptr;
769 std::unique_ptr<t_state> globalState;
771 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
773 if (isSimulationMasterRank)
775 /* Only the master rank has the global state */
776 globalState = std::make_unique<t_state>();
778 /* Read (nearly) all data required for the simulation
779 * and keep the partly serialized tpr contents to send to other ranks later
781 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
782 &inputrecInstance, globalState.get(), &mtop);
783 inputrec = &inputrecInstance;
786 /* Check and update the hardware options for internal consistency */
787 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
790 if (GMX_THREAD_MPI && isSimulationMasterRank)
792 bool useGpuForNonbonded = false;
793 bool useGpuForPme = false;
796 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
798 // If the user specified the number of ranks, then we must
799 // respect that, but in default mode, we need to allow for
800 // the number of GPUs to choose the number of ranks.
801 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
802 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
803 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
804 canUseGpuForNonbonded,
805 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
806 hw_opt.nthreads_tmpi);
807 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
808 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
809 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
811 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
813 /* Determine how many thread-MPI ranks to start.
815 * TODO Over-writing the user-supplied value here does
816 * prevent any possible subsequent checks from working
818 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
819 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
821 // Now start the threads for thread MPI.
822 spawnThreads(hw_opt.nthreads_tmpi);
823 // The spawned threads enter mdrunner() and execution of
824 // master and spawned threads joins at the end of this block.
825 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
828 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
829 CommrecHandle crHandle = init_commrec(communicator, ms);
830 t_commrec* cr = crHandle.get();
831 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
835 /* now broadcast everything to the non-master nodes/threads: */
836 if (!isSimulationMasterRank)
838 inputrec = &inputrecInstance;
840 init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
842 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
843 partialDeserializedTpr.reset(nullptr);
845 // Now the number of ranks is known to all ranks, and each knows
846 // the inputrec read by the master rank. The ranks can now all run
847 // the task-deciding functions and will agree on the result
848 // without needing to communicate.
850 // TODO Should we do the communication in debug mode to support
851 // having an assertion?
852 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
854 // Note that these variables describe only their own node.
856 // Note that when bonded interactions run on a GPU they always run
857 // alongside a nonbonded task, so do not influence task assignment
858 // even though they affect the force calculation workload.
859 bool useGpuForNonbonded = false;
860 bool useGpuForPme = false;
861 bool useGpuForBonded = false;
862 bool useGpuForUpdate = false;
863 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
866 // It's possible that there are different numbers of GPUs on
867 // different nodes, which is the user's responsibility to
868 // handle. If unsuitable, we will notice that during task
870 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
871 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
872 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
873 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
874 useGpuForPme = decideWhetherToUseGpusForPme(
875 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
876 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
877 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
878 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
879 useGpuForBonded = decideWhetherToUseGpusForBonded(
880 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
881 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
882 domdecOptions.numPmeRanks, gpusWereDetected);
884 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
885 if (pmeRunMode == PmeRunMode::GPU)
887 if (pmeFftTarget == TaskTarget::Cpu)
889 pmeRunMode = PmeRunMode::Mixed;
892 else if (pmeFftTarget == TaskTarget::Gpu)
895 "Assigning FFTs to GPU requires PME to be assigned to GPU as well. With PME "
896 "on CPU you should not be using -pmefft.");
899 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
901 // Initialize development feature flags that enabled by environment variable
902 // and report those features that are enabled.
903 const DevelopmentFeatureFlags devFlags =
904 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
907 // TODO: hide restraint implementation details from Mdrunner.
908 // There is nothing unique about restraints at this point as far as the
909 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
910 // factory functions from the SimulationContext on which to call mdModules_->add().
911 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
912 for (auto&& restraint : restraintManager_->getRestraints())
914 auto module = RestraintMDModule::create(restraint, restraint->sites());
915 mdModules_->add(std::move(module));
918 // TODO: Error handling
919 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
920 const auto& mdModulesNotifier = mdModules_->notifier().notifier_;
922 if (inputrec->internalParameters != nullptr)
924 mdModulesNotifier.notify(*inputrec->internalParameters);
927 if (fplog != nullptr)
929 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
930 fprintf(fplog, "\n");
935 /* In rerun, set velocities to zero if present */
936 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
938 // rerun does not use velocities
942 "Rerun trajectory contains velocities. Rerun does only evaluate "
943 "potential energy and forces. The velocities will be ignored.");
944 for (int i = 0; i < globalState->natoms; i++)
946 clear_rvec(globalState->v[i]);
948 globalState->flags &= ~(1 << estV);
951 /* now make sure the state is initialized and propagated */
952 set_state_entries(globalState.get(), inputrec);
955 /* NM and TPI parallelize over force/energy calculations, not atoms,
956 * so we need to initialize and broadcast the global state.
958 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
962 globalState = std::make_unique<t_state>();
964 broadcastStateWithoutDynamics(cr, globalState.get());
967 /* A parallel command line option consistency check that we can
968 only do after any threads have started. */
970 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
971 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
974 "The -dd or -npme option request a parallel simulation, "
976 "but %s was compiled without threads or MPI enabled",
977 output_env_get_program_display_name(oenv));
979 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
981 "but %s was not started through mpirun/mpiexec or only one rank was requested "
982 "through mpirun/mpiexec",
983 output_env_get_program_display_name(oenv));
987 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
990 "The .mdp file specified an energy mininization or normal mode algorithm, and "
991 "these are not compatible with mdrun -rerun");
994 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
996 if (domdecOptions.numPmeRanks > 0)
998 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
999 "PME-only ranks are requested, but the system does not use PME "
1000 "for electrostatics or LJ");
1003 domdecOptions.numPmeRanks = 0;
1006 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1008 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1009 * improve performance with many threads per GPU, since our OpenMP
1010 * scaling is bad, but it's difficult to automate the setup.
1012 domdecOptions.numPmeRanks = 0;
1016 if (domdecOptions.numPmeRanks < 0)
1018 domdecOptions.numPmeRanks = 0;
1019 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1023 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1024 "PME GPU decomposition is not supported");
1031 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1035 /* NMR restraints must be initialized before load_checkpoint,
1036 * since with time averaging the history is added to t_state.
1037 * For proper consistency check we therefore need to extend
1039 * So the PME-only nodes (if present) will also initialize
1040 * the distance restraints.
1044 /* This needs to be called before read_checkpoint to extend the state */
1045 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1047 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1049 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1051 ObservablesHistory observablesHistory = {};
1053 if (startingBehavior != StartingBehavior::NewSimulation)
1055 /* Check if checkpoint file exists before doing continuation.
1056 * This way we can use identical input options for the first and subsequent runs...
1058 if (mdrunOptions.numStepsCommandline > -2)
1060 /* Temporarily set the number of steps to unlimited to avoid
1061 * triggering the nsteps check in load_checkpoint().
1062 * This hack will go away soon when the -nsteps option is removed.
1064 inputrec->nsteps = -1;
1067 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1068 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1069 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1071 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1073 // Now we can start normal logging to the truncated log file.
1074 fplog = gmx_fio_getfp(logFileHandle);
1075 prepareLogAppending(fplog);
1076 logOwner = buildLogger(fplog, MASTER(cr));
1077 mdlog = logOwner.logger();
1081 if (mdrunOptions.numStepsCommandline > -2)
1086 "The -nsteps functionality is deprecated, and may be removed in a future "
1088 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1091 /* override nsteps with value set on the commandline */
1092 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1096 copy_mat(globalState->box, box);
1101 gmx_bcast(sizeof(box), box, cr);
1104 if (inputrec->cutoff_scheme != ecutsVERLET)
1107 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1108 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1110 /* Update rlist and nstlist. */
1111 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1112 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1115 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1116 // This builder is necessary while we have multi-part construction
1117 // of DD. Before DD is constructed, we use the existence of
1118 // the builder object to indicate that further construction of DD
1120 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1121 if (useDomainDecomposition)
1123 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1124 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1125 positionsFromStatePointer(globalState.get()));
1129 /* PME, if used, is done on all nodes with 1D decomposition */
1131 cr->duty = (DUTY_PP | DUTY_PME);
1133 if (inputrec->ePBC == epbcSCREW)
1135 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1139 // Produce the task assignment for this rank.
1140 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1141 GpuTaskAssignments gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1142 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, cr, ms, physicalNodeComm, nonbondedTarget,
1143 pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded, useGpuForPme,
1144 thisRankHasDuty(cr, DUTY_PP),
1145 // TODO cr->duty & DUTY_PME should imply that a PME
1146 // algorithm is active, but currently does not.
1147 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1149 const bool printHostName = (cr->nnodes > 1);
1150 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode);
1152 // If the user chose a task assignment, give them some hints
1153 // where appropriate.
1154 if (!userGpuTaskAssignment.empty())
1156 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1159 // Get the device handles for the modules, nullptr when no task is assigned.
1160 gmx_device_info_t* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1161 gmx_device_info_t* pmeDeviceInfo = gpuTaskAssignments.initPmeDevice();
1163 // TODO Initialize GPU streams here.
1165 // TODO Currently this is always built, yet DD partition code
1166 // checks if it is built before using it. Probably it should
1167 // become an MDModule that is made only when another module
1168 // requires it (e.g. pull, CompEl, density fitting), so that we
1169 // don't update the local atom sets unilaterally every step.
1170 LocalAtomSetManager atomSets;
1173 // TODO Pass the GPU streams to ddBuilder to use in buffer
1174 // transfers (e.g. halo exchange)
1175 cr->dd = ddBuilder->build(&atomSets);
1176 // The builder's job is done, so destruct it
1177 ddBuilder.reset(nullptr);
1178 // Note that local state still does not exist yet.
1181 // The GPU update is decided here because we need to know whether the constraints or
1182 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1183 // defined). This is only known after DD is initialized, hence decision on using GPU
1184 // update is done so late.
1187 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1189 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1190 devFlags.forceGpuUpdateDefaultWithDD, useDomainDecomposition, useUpdateGroups,
1191 useGpuForPme, useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1192 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1193 replExParams.exchangeInterval > 0, doRerun);
1195 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1199 /* After possible communicator splitting in make_dd_communicators.
1200 * we can set up the intra/inter node communication.
1202 gmx_setup_nodecomm(fplog, cr);
1208 GMX_LOG(mdlog.warning)
1210 .appendTextFormatted(
1211 "This is simulation %d out of %d running as a composite GROMACS\n"
1212 "multi-simulation job. Setup for this simulation:\n",
1215 GMX_LOG(mdlog.warning)
1216 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1218 cr->nnodes == 1 ? "thread" : "threads"
1220 cr->nnodes == 1 ? "process" : "processes"
1226 // If mdrun -pin auto honors any affinity setting that already
1227 // exists. If so, it is nice to provide feedback about whether
1228 // that existing affinity setting was from OpenMP or something
1229 // else, so we run this code both before and after we initialize
1230 // the OpenMP support.
1231 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1232 /* Check and update the number of OpenMP threads requested */
1233 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1234 pmeRunMode, mtop, *inputrec);
1236 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1237 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1239 // Enable FP exception detection, but not in
1240 // Release mode and not for compilers with known buggy FP
1241 // exception support (clang with any optimization) or suspected
1242 // buggy FP exception support (gcc 7.* with optimization).
1243 #if !defined NDEBUG \
1244 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1245 && defined __OPTIMIZE__)
1246 const bool bEnableFPE = true;
1248 const bool bEnableFPE = false;
1250 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1253 gmx_feenableexcept();
1256 /* Now that we know the setup is consistent, check for efficiency */
1257 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1258 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1260 /* getting number of PP/PME threads on this MPI / tMPI rank.
1261 PME: env variable should be read only on one node to make sure it is
1262 identical everywhere;
1264 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1265 : gmx_omp_nthreads_get(emntPME);
1266 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1267 physicalNodeComm, mdlog);
1269 // Enable Peer access between GPUs where available
1270 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1271 // any of the GPU communication features are active.
1272 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1273 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1275 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1278 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1280 /* Before setting affinity, check whether the affinity has changed
1281 * - which indicates that probably the OpenMP library has changed it
1282 * since we first checked).
1284 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1286 int numThreadsOnThisNode, intraNodeThreadOffset;
1287 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1288 &intraNodeThreadOffset);
1290 /* Set the CPU affinity */
1291 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1292 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1295 if (mdrunOptions.timingOptions.resetStep > -1)
1300 "The -resetstep functionality is deprecated, and may be removed in a "
1303 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1307 /* Master synchronizes its value of reset_counters with all nodes
1308 * including PME only nodes */
1309 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1310 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1311 wcycle_set_reset_counters(wcycle, reset_counters);
1314 // Membrane embedding must be initialized before we call init_forcerec()
1319 fprintf(stderr, "Initializing membed");
1321 /* Note that membed cannot work in parallel because mtop is
1322 * changed here. Fix this if we ever want to make it run with
1323 * multiple ranks. */
1324 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1325 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1328 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1329 std::unique_ptr<MDAtoms> mdAtoms;
1330 std::unique_ptr<gmx_vsite_t> vsite;
1333 if (thisRankHasDuty(cr, DUTY_PP))
1335 mdModulesNotifier.notify(*cr);
1336 mdModulesNotifier.notify(&atomSets);
1337 mdModulesNotifier.notify(PeriodicBoundaryConditionType{ inputrec->ePBC });
1338 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1339 /* Initiate forcerecord */
1340 fr = new t_forcerec;
1341 fr->forceProviders = mdModules_->initForceProviders();
1342 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1343 opt2fn("-table", filenames.size(), filenames.data()),
1344 opt2fn("-tablep", filenames.size(), filenames.data()),
1345 opt2fns("-tableb", filenames.size(), filenames.data()), *hwinfo,
1346 nonbondedDeviceInfo, useGpuForBonded,
1347 pmeRunMode == PmeRunMode::GPU && !thisRankHasDuty(cr, DUTY_PME), pforce, wcycle);
1349 // TODO Move this to happen during domain decomposition setup,
1350 // once stream and event handling works well with that.
1351 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
1352 if (havePPDomainDecomposition(cr) && prefer1DAnd1PulseDD && is1DAnd1PulseDD(*cr->dd))
1354 GMX_RELEASE_ASSERT(devFlags.enableGpuBufferOps,
1355 "Must use GMX_USE_GPU_BUFFER_OPS=1 to use GMX_GPU_DD_COMMS=1");
1357 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local);
1358 void* streamNonLocal =
1359 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal);
1360 GMX_LOG(mdlog.warning)
1362 .appendTextFormatted(
1363 "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the "
1364 "GMX_GPU_DD_COMMS environment variable.");
1365 cr->dd->gpuHaloExchange = std::make_unique<GpuHaloExchange>(
1366 cr->dd, cr->mpi_comm_mysim, streamLocal, streamNonLocal);
1369 /* Initialize the mdAtoms structure.
1370 * mdAtoms is not filled with atom data,
1371 * as this can not be done now with domain decomposition.
1373 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1374 if (globalState && thisRankHasPmeGpuTask)
1376 // The pinning of coordinates in the global state object works, because we only use
1377 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1378 // points to the global state object without DD.
1379 // FIXME: MD and EM separately set up the local state - this should happen in the same
1380 // function, which should also perform the pinning.
1381 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1384 /* Initialize the virtual site communication */
1385 vsite = initVsite(mtop, cr);
1387 calc_shifts(box, fr->shift_vec);
1389 /* With periodic molecules the charge groups should be whole at start up
1390 * and the virtual sites should not be far from their proper positions.
1392 if (!inputrec->bContinuation && MASTER(cr) && !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1394 /* Make molecules whole at start of run */
1395 if (fr->ePBC != epbcNONE)
1397 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1401 /* Correct initial vsite positions are required
1402 * for the initial distribution in the domain decomposition
1403 * and for the initial shell prediction.
1405 constructVsitesGlobal(mtop, globalState->x);
1409 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1411 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1412 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1417 /* This is a PME only node */
1419 GMX_ASSERT(globalState == nullptr,
1420 "We don't need the state on a PME only rank and expect it to be unitialized");
1422 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1423 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1426 gmx_pme_t* sepPmeData = nullptr;
1427 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1428 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1429 "Double-checking that only PME-only ranks have no forcerec");
1430 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1432 // TODO should live in ewald module once its testing is improved
1434 // Later, this program could contain kernels that might be later
1435 // re-used as auto-tuning progresses, or subsequent simulations
1437 PmeGpuProgramStorage pmeGpuProgram;
1438 if (thisRankHasPmeGpuTask)
1440 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1443 /* Initiate PME if necessary,
1444 * either on all nodes or on dedicated PME nodes only. */
1445 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1447 if (mdAtoms && mdAtoms->mdatoms())
1449 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1450 if (EVDW_PME(inputrec->vdwtype))
1452 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1455 if (cr->npmenodes > 0)
1457 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1458 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1459 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1462 if (thisRankHasDuty(cr, DUTY_PME))
1466 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1467 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1468 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1469 nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1471 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1476 if (EI_DYNAMICS(inputrec->eI))
1478 /* Turn on signal handling on all nodes */
1480 * (A user signal from the PME nodes (if any)
1481 * is communicated to the PP nodes.
1483 signal_handler_install();
1486 pull_t* pull_work = nullptr;
1487 if (thisRankHasDuty(cr, DUTY_PP))
1489 /* Assumes uniform use of the number of OpenMP threads */
1490 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1492 if (inputrec->bPull)
1494 /* Initialize pull code */
1495 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1496 inputrec->fepvals->init_lambda);
1497 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1499 initPullHistory(pull_work, &observablesHistory);
1501 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1503 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1507 std::unique_ptr<EnforcedRotation> enforcedRotation;
1510 /* Initialize enforced rotation code */
1512 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1513 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1516 t_swap* swap = nullptr;
1517 if (inputrec->eSwapCoords != eswapNO)
1519 /* Initialize ion swapping code */
1520 swap = init_swapcoords(fplog, inputrec,
1521 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1522 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1523 oenv, mdrunOptions, startingBehavior);
1526 /* Let makeConstraints know whether we have essential dynamics constraints. */
1527 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1528 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1530 /* Energy terms and groups */
1531 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1532 inputrec->fepvals->n_lambda);
1534 /* Kinetic energy data */
1535 gmx_ekindata_t ekind;
1536 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1538 /* Set up interactive MD (IMD) */
1540 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1541 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1542 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1544 if (DOMAINDECOMP(cr))
1546 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1547 /* This call is not included in init_domain_decomposition mainly
1548 * because fr->cginfo_mb is set later.
1550 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1551 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1554 const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
1555 false, inputrec, doRerun, vsite.get(), ms, replExParams, fcd,
1556 static_cast<int>(filenames.size()), filenames.data(), &observablesHistory, membed);
1558 const bool useModularSimulator = inputIsCompatibleWithModularSimulator
1559 && !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
1561 // TODO This is not the right place to manage the lifetime of
1562 // this data structure, but currently it's the easiest way to
1564 MdrunScheduleWorkload runScheduleWork;
1565 // Also populates the simulation constant workload description.
1566 runScheduleWork.simulationWork = createSimulationWorkload(
1567 useGpuForNonbonded, pmeRunMode, useGpuForBonded, useGpuForUpdate,
1568 devFlags.enableGpuBufferOps, devFlags.enableGpuHaloExchange,
1569 devFlags.enableGpuPmePPComm, haveEwaldSurfaceContribution(*inputrec));
1571 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1572 if (gpusWereDetected
1573 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1574 || runScheduleWork.simulationWork.useGpuBufferOps))
1576 const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1577 const void* localStream =
1578 fr->nbv->gpu_nbv != nullptr
1579 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1581 const void* nonLocalStream =
1582 fr->nbv->gpu_nbv != nullptr
1583 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1585 const void* deviceContext = pme_gpu_get_device_context(fr->pmedata);
1586 const int paddingSize = pme_gpu_get_padding_size(fr->pmedata);
1587 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1588 ? GpuApiCallBehavior::Async
1589 : GpuApiCallBehavior::Sync;
1591 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1592 pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
1593 fr->stateGpu = stateGpu.get();
1596 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1597 SimulatorBuilder simulatorBuilder;
1599 // build and run simulator object based on user-input
1600 auto simulator = simulatorBuilder.build(
1601 inputIsCompatibleWithModularSimulator, fplog, cr, ms, mdlog,
1602 static_cast<int>(filenames.size()), filenames.data(), oenv, mdrunOptions,
1603 startingBehavior, vsite.get(), constr.get(),
1604 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1605 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1606 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1607 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1608 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1611 if (fr->pmePpCommGpu)
1613 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1614 fr->pmePpCommGpu.reset();
1617 if (inputrec->bPull)
1619 finish_pull(pull_work);
1621 finish_swapcoords(swap);
1625 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1627 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1628 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1631 wallcycle_stop(wcycle, ewcRUN);
1633 /* Finish up, write some stuff
1634 * if rerunMD, don't write last frame again
1636 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1637 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1639 // clean up cycle counter
1640 wallcycle_destroy(wcycle);
1645 gmx_pme_destroy(pmedata);
1649 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1650 // before we destroy the GPU context(s) in free_gpu_resources().
1651 // Pinned buffers are associated with contexts in CUDA.
1652 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1653 mdAtoms.reset(nullptr);
1654 globalState.reset(nullptr);
1655 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1657 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1658 free_gpu_resources(fr, physicalNodeComm, hwinfo->gpu_info);
1659 free_gpu(nonbondedDeviceInfo);
1660 free_gpu(pmeDeviceInfo);
1661 done_forcerec(fr, mtop.molblock.size());
1666 free_membed(membed);
1669 /* Does what it says */
1670 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1671 walltime_accounting_destroy(walltime_accounting);
1673 // Ensure log file content is written
1676 gmx_fio_flush(logFileHandle);
1679 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1680 * exceptions were enabled before function was called. */
1683 gmx_fedisableexcept();
1686 auto rc = static_cast<int>(gmx_get_stop_condition());
1689 /* we need to join all threads. The sub-threads join when they
1690 exit this function, but the master thread needs to be told to
1692 if (PAR(cr) && MASTER(cr))
1700 Mdrunner::~Mdrunner()
1702 // Clean up of the Manager.
1703 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1704 // but okay as long as threads synchronize some time before adding or accessing
1705 // a new set of restraints.
1706 if (restraintManager_)
1708 restraintManager_->clear();
1709 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1710 "restraints added during runner life time should be cleared at runner "
1715 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1717 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1718 // Not sure if this should be logged through the md logger or something else,
1719 // but it is helpful to have some sort of INFO level message sent somewhere.
1720 // std::cout << "Registering restraint named " << name << std::endl;
1722 // When multiple restraints are used, it may be wasteful to register them separately.
1723 // Maybe instead register an entire Restraint Manager as a force provider.
1724 restraintManager_->addToSpec(std::move(puller), name);
1727 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1729 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1731 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1732 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1734 class Mdrunner::BuilderImplementation
1737 BuilderImplementation() = delete;
1738 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1739 ~BuilderImplementation();
1741 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1742 real forceWarningThreshold,
1743 StartingBehavior startingBehavior);
1745 void addDomdec(const DomdecOptions& options);
1747 void addVerletList(int nstlist);
1749 void addReplicaExchange(const ReplicaExchangeParameters& params);
1751 void addNonBonded(const char* nbpu_opt);
1753 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1755 void addBondedTaskAssignment(const char* bonded_opt);
1757 void addUpdateTaskAssignment(const char* update_opt);
1759 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1761 void addFilenames(ArrayRef<const t_filenm> filenames);
1763 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1765 void addLogFile(t_fileio* logFileHandle);
1767 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1772 // Default parameters copied from runner.h
1773 // \todo Clarify source(s) of default parameters.
1775 const char* nbpu_opt_ = nullptr;
1776 const char* pme_opt_ = nullptr;
1777 const char* pme_fft_opt_ = nullptr;
1778 const char* bonded_opt_ = nullptr;
1779 const char* update_opt_ = nullptr;
1781 MdrunOptions mdrunOptions_;
1783 DomdecOptions domdecOptions_;
1785 ReplicaExchangeParameters replicaExchangeParameters_;
1787 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1790 //! Multisim communicator handle.
1791 gmx_multisim_t* multiSimulation_;
1793 //! mdrun communicator
1794 MPI_Comm communicator_ = MPI_COMM_NULL;
1796 //! Print a warning if any force is larger than this (in kJ/mol nm).
1797 real forceWarningThreshold_ = -1;
1799 //! Whether the simulation will start afresh, or restart with/without appending.
1800 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1802 //! The modules that comprise the functionality of mdrun.
1803 std::unique_ptr<MDModules> mdModules_;
1805 //! \brief Parallelism information.
1806 gmx_hw_opt_t hardwareOptions_;
1808 //! filename options for simulation.
1809 ArrayRef<const t_filenm> filenames_;
1811 /*! \brief Handle to output environment.
1813 * \todo gmx_output_env_t needs lifetime management.
1815 gmx_output_env_t* outputEnvironment_ = nullptr;
1817 /*! \brief Non-owning handle to MD log file.
1819 * \todo Context should own output facilities for client.
1820 * \todo Improve log file handle management.
1822 * Code managing the FILE* relies on the ability to set it to
1823 * nullptr to check whether the filehandle is valid.
1825 t_fileio* logFileHandle_ = nullptr;
1828 * \brief Builder for simulation stop signal handler.
1830 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1833 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1834 compat::not_null<SimulationContext*> context) :
1835 mdModules_(std::move(mdModules))
1837 communicator_ = context->communicator_;
1838 multiSimulation_ = context->multiSimulation_.get();
1841 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1843 Mdrunner::BuilderImplementation&
1844 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1845 const real forceWarningThreshold,
1846 const StartingBehavior startingBehavior)
1848 mdrunOptions_ = options;
1849 forceWarningThreshold_ = forceWarningThreshold;
1850 startingBehavior_ = startingBehavior;
1854 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1856 domdecOptions_ = options;
1859 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1864 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1866 replicaExchangeParameters_ = params;
1869 Mdrunner Mdrunner::BuilderImplementation::build()
1871 auto newRunner = Mdrunner(std::move(mdModules_));
1873 newRunner.mdrunOptions = mdrunOptions_;
1874 newRunner.pforce = forceWarningThreshold_;
1875 newRunner.startingBehavior = startingBehavior_;
1876 newRunner.domdecOptions = domdecOptions_;
1878 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1879 newRunner.hw_opt = hardwareOptions_;
1881 // No invariant to check. This parameter exists to optionally override other behavior.
1882 newRunner.nstlist_cmdline = nstlist_;
1884 newRunner.replExParams = replicaExchangeParameters_;
1886 newRunner.filenames = filenames_;
1888 newRunner.communicator = communicator_;
1890 // nullptr is a valid value for the multisim handle
1891 newRunner.ms = multiSimulation_;
1893 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1894 // \todo Update sanity checking when output environment has clearly specified invariants.
1895 // Initialization and default values for oenv are not well specified in the current version.
1896 if (outputEnvironment_)
1898 newRunner.oenv = outputEnvironment_;
1902 GMX_THROW(gmx::APIError(
1903 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1906 newRunner.logFileHandle = logFileHandle_;
1910 newRunner.nbpu_opt = nbpu_opt_;
1914 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1917 if (pme_opt_ && pme_fft_opt_)
1919 newRunner.pme_opt = pme_opt_;
1920 newRunner.pme_fft_opt = pme_fft_opt_;
1924 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1929 newRunner.bonded_opt = bonded_opt_;
1933 GMX_THROW(gmx::APIError(
1934 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1939 newRunner.update_opt = update_opt_;
1943 GMX_THROW(gmx::APIError(
1944 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1948 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1950 if (stopHandlerBuilder_)
1952 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1956 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1962 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1964 nbpu_opt_ = nbpu_opt;
1967 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
1970 pme_fft_opt_ = pme_fft_opt;
1973 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1975 bonded_opt_ = bonded_opt;
1978 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
1980 update_opt_ = update_opt;
1983 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
1985 hardwareOptions_ = hardwareOptions;
1988 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1990 filenames_ = filenames;
1993 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1995 outputEnvironment_ = outputEnvironment;
1998 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2000 logFileHandle_ = logFileHandle;
2003 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2005 stopHandlerBuilder_ = std::move(builder);
2008 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2009 compat::not_null<SimulationContext*> context) :
2010 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2014 MdrunnerBuilder::~MdrunnerBuilder() = default;
2016 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2017 real forceWarningThreshold,
2018 const StartingBehavior startingBehavior)
2020 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2024 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2026 impl_->addDomdec(options);
2030 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2032 impl_->addVerletList(nstlist);
2036 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2038 impl_->addReplicaExchange(params);
2042 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2044 impl_->addNonBonded(nbpu_opt);
2048 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2050 // The builder method may become more general in the future, but in this version,
2051 // parameters for PME electrostatics are both required and the only parameters
2053 if (pme_opt && pme_fft_opt)
2055 impl_->addPME(pme_opt, pme_fft_opt);
2060 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2065 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2067 impl_->addBondedTaskAssignment(bonded_opt);
2071 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2073 impl_->addUpdateTaskAssignment(update_opt);
2077 Mdrunner MdrunnerBuilder::build()
2079 return impl_->build();
2082 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2084 impl_->addHardwareOptions(hardwareOptions);
2088 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2090 impl_->addFilenames(filenames);
2094 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2096 impl_->addOutputEnvironment(outputEnvironment);
2100 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2102 impl_->addLogFile(logFileHandle);
2106 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2108 impl_->addStopHandlerBuilder(std::move(builder));
2112 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2114 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;