2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/gpu_utils.h"
77 #include "gromacs/hardware/cpuinfo.h"
78 #include "gromacs/hardware/detecthardware.h"
79 #include "gromacs/hardware/printhardware.h"
80 #include "gromacs/imd/imd.h"
81 #include "gromacs/listed_forces/disre.h"
82 #include "gromacs/listed_forces/gpubonded.h"
83 #include "gromacs/listed_forces/orires.h"
84 #include "gromacs/math/functions.h"
85 #include "gromacs/math/utilities.h"
86 #include "gromacs/math/vec.h"
87 #include "gromacs/mdlib/boxdeformation.h"
88 #include "gromacs/mdlib/broadcaststructs.h"
89 #include "gromacs/mdlib/calc_verletbuf.h"
90 #include "gromacs/mdlib/dispersioncorrection.h"
91 #include "gromacs/mdlib/enerdata_utils.h"
92 #include "gromacs/mdlib/force.h"
93 #include "gromacs/mdlib/forcerec.h"
94 #include "gromacs/mdlib/gmx_omp_nthreads.h"
95 #include "gromacs/mdlib/makeconstraints.h"
96 #include "gromacs/mdlib/md_support.h"
97 #include "gromacs/mdlib/mdatoms.h"
98 #include "gromacs/mdlib/membed.h"
99 #include "gromacs/mdlib/qmmm.h"
100 #include "gromacs/mdlib/sighandler.h"
101 #include "gromacs/mdlib/stophandler.h"
102 #include "gromacs/mdlib/updategroups.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdrun/mdmodules.h"
105 #include "gromacs/mdrun/simulationcontext.h"
106 #include "gromacs/mdrunutility/handlerestart.h"
107 #include "gromacs/mdrunutility/logging.h"
108 #include "gromacs/mdrunutility/multisim.h"
109 #include "gromacs/mdrunutility/printtime.h"
110 #include "gromacs/mdrunutility/threadaffinity.h"
111 #include "gromacs/mdtypes/commrec.h"
112 #include "gromacs/mdtypes/enerdata.h"
113 #include "gromacs/mdtypes/fcdata.h"
114 #include "gromacs/mdtypes/forcerec.h"
115 #include "gromacs/mdtypes/group.h"
116 #include "gromacs/mdtypes/inputrec.h"
117 #include "gromacs/mdtypes/interaction_const.h"
118 #include "gromacs/mdtypes/md_enums.h"
119 #include "gromacs/mdtypes/mdatom.h"
120 #include "gromacs/mdtypes/mdrunoptions.h"
121 #include "gromacs/mdtypes/observableshistory.h"
122 #include "gromacs/mdtypes/simulation_workload.h"
123 #include "gromacs/mdtypes/state.h"
124 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
125 #include "gromacs/nbnxm/gpu_data_mgmt.h"
126 #include "gromacs/nbnxm/nbnxm.h"
127 #include "gromacs/nbnxm/pairlist_tuning.h"
128 #include "gromacs/pbcutil/pbc.h"
129 #include "gromacs/pulling/output.h"
130 #include "gromacs/pulling/pull.h"
131 #include "gromacs/pulling/pull_rotation.h"
132 #include "gromacs/restraint/manager.h"
133 #include "gromacs/restraint/restraintmdmodule.h"
134 #include "gromacs/restraint/restraintpotential.h"
135 #include "gromacs/swap/swapcoords.h"
136 #include "gromacs/taskassignment/decidegpuusage.h"
137 #include "gromacs/taskassignment/decidesimulationworkload.h"
138 #include "gromacs/taskassignment/resourcedivision.h"
139 #include "gromacs/taskassignment/taskassignment.h"
140 #include "gromacs/taskassignment/usergpuids.h"
141 #include "gromacs/timing/gpu_timing.h"
142 #include "gromacs/timing/wallcycle.h"
143 #include "gromacs/timing/wallcyclereporting.h"
144 #include "gromacs/topology/mtop_util.h"
145 #include "gromacs/trajectory/trajectoryframe.h"
146 #include "gromacs/utility/basenetwork.h"
147 #include "gromacs/utility/cstringutil.h"
148 #include "gromacs/utility/exceptions.h"
149 #include "gromacs/utility/fatalerror.h"
150 #include "gromacs/utility/filestream.h"
151 #include "gromacs/utility/gmxassert.h"
152 #include "gromacs/utility/gmxmpi.h"
153 #include "gromacs/utility/keyvaluetree.h"
154 #include "gromacs/utility/logger.h"
155 #include "gromacs/utility/loggerbuilder.h"
156 #include "gromacs/utility/mdmodulenotification.h"
157 #include "gromacs/utility/physicalnodecommunicator.h"
158 #include "gromacs/utility/pleasecite.h"
159 #include "gromacs/utility/programcontext.h"
160 #include "gromacs/utility/smalloc.h"
161 #include "gromacs/utility/stringutil.h"
163 #include "isimulator.h"
164 #include "replicaexchange.h"
165 #include "simulatorbuilder.h"
168 # include "corewrap.h"
175 /*! \brief Manage any development feature flag variables encountered
177 * The use of dev features indicated by environment variables is
178 * logged in order to ensure that runs with such features enabled can
179 * be identified from their log and standard output. Any cross
180 * dependencies are also checked, and if unsatisfied, a fatal error
183 * Note that some development features overrides are applied already here:
184 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
186 * \param[in] mdlog Logger object.
187 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
188 * \param[in] pmeRunMode The PME run mode for this run
189 * \returns The object populated with development feature flags.
191 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
192 const bool useGpuForNonbonded,
193 const PmeRunMode pmeRunMode)
195 DevelopmentFeatureFlags devFlags;
197 // Some builds of GCC 5 give false positive warnings that these
198 // getenv results are ignored when clearly they are used.
199 #pragma GCC diagnostic push
200 #pragma GCC diagnostic ignored "-Wunused-result"
201 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
202 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
203 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
204 devFlags.enableGpuHaloExchange =
205 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
206 devFlags.enableGpuPmePPComm =
207 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
208 #pragma GCC diagnostic pop
210 if (devFlags.enableGpuBufferOps)
212 GMX_LOG(mdlog.warning)
214 .appendTextFormatted(
215 "This run uses the 'GPU buffer ops' feature, enabled by the "
216 "GMX_USE_GPU_BUFFER_OPS environment variable.");
219 if (devFlags.forceGpuUpdateDefault)
221 GMX_LOG(mdlog.warning)
223 .appendTextFormatted(
224 "This run will default to '-update gpu' as requested by the "
225 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
226 "decomposition lacks substantial testing and should be used with caution.");
229 if (devFlags.enableGpuHaloExchange)
231 if (useGpuForNonbonded)
233 if (!devFlags.enableGpuBufferOps)
235 GMX_LOG(mdlog.warning)
237 .appendTextFormatted(
238 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
239 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
240 devFlags.enableGpuBufferOps = true;
242 GMX_LOG(mdlog.warning)
244 .appendTextFormatted(
245 "This run has requested the 'GPU halo exchange' feature, enabled by "
247 "GMX_GPU_DD_COMMS environment variable.");
251 GMX_LOG(mdlog.warning)
253 .appendTextFormatted(
254 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
255 "halo exchange' feature will not be enabled as nonbonded interactions "
256 "are not offloaded.");
257 devFlags.enableGpuHaloExchange = false;
261 if (devFlags.enableGpuPmePPComm)
263 if (pmeRunMode == PmeRunMode::GPU)
265 GMX_LOG(mdlog.warning)
267 .appendTextFormatted(
268 "This run uses the 'GPU PME-PP communications' feature, enabled "
269 "by the GMX_GPU_PME_PP_COMMS environment variable.");
273 std::string clarification;
274 if (pmeRunMode == PmeRunMode::Mixed)
277 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
282 clarification = "PME is not offloaded to the GPU.";
284 GMX_LOG(mdlog.warning)
287 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
288 "'GPU PME-PP communications' feature was not enabled as "
290 devFlags.enableGpuPmePPComm = false;
297 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
299 * Used to ensure that the master thread does not modify mdrunner during copy
300 * on the spawned threads. */
301 static void threadMpiMdrunnerAccessBarrier()
304 MPI_Barrier(MPI_COMM_WORLD);
308 Mdrunner Mdrunner::cloneOnSpawnedThread() const
310 auto newRunner = Mdrunner(std::make_unique<MDModules>());
312 // All runners in the same process share a restraint manager resource because it is
313 // part of the interface to the client code, which is associated only with the
314 // original thread. Handles to the same resources can be obtained by copy.
316 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
319 // Copy members of master runner.
320 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
321 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
322 newRunner.hw_opt = hw_opt;
323 newRunner.filenames = filenames;
325 newRunner.oenv = oenv;
326 newRunner.mdrunOptions = mdrunOptions;
327 newRunner.domdecOptions = domdecOptions;
328 newRunner.nbpu_opt = nbpu_opt;
329 newRunner.pme_opt = pme_opt;
330 newRunner.pme_fft_opt = pme_fft_opt;
331 newRunner.bonded_opt = bonded_opt;
332 newRunner.update_opt = update_opt;
333 newRunner.nstlist_cmdline = nstlist_cmdline;
334 newRunner.replExParams = replExParams;
335 newRunner.pforce = pforce;
336 // Give the spawned thread the newly created valid communicator
337 // for the simulation.
338 newRunner.communicator = MPI_COMM_WORLD;
340 newRunner.startingBehavior = startingBehavior;
341 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
343 threadMpiMdrunnerAccessBarrier();
348 /*! \brief The callback used for running on spawned threads.
350 * Obtains the pointer to the master mdrunner object from the one
351 * argument permitted to the thread-launch API call, copies it to make
352 * a new runner for this thread, reinitializes necessary data, and
353 * proceeds to the simulation. */
354 static void mdrunner_start_fn(const void* arg)
358 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
359 /* copy the arg list to make sure that it's thread-local. This
360 doesn't copy pointed-to items, of course; fnm, cr and fplog
361 are reset in the call below, all others should be const. */
362 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
365 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
369 void Mdrunner::spawnThreads(int numThreadsToLaunch)
372 /* now spawn new threads that start mdrunner_start_fn(), while
373 the main thread returns. Thread affinity is handled later. */
374 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
375 static_cast<const void*>(this))
378 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
381 // Give the master thread the newly created valid communicator for
383 communicator = MPI_COMM_WORLD;
384 threadMpiMdrunnerAccessBarrier();
386 GMX_UNUSED_VALUE(numThreadsToLaunch);
387 GMX_UNUSED_VALUE(mdrunner_start_fn);
393 /*! \brief Initialize variables for Verlet scheme simulation */
394 static void prepare_verlet_scheme(FILE* fplog,
398 const gmx_mtop_t* mtop,
400 bool makeGpuPairList,
401 const gmx::CpuInfo& cpuinfo)
403 // We checked the cut-offs in grompp, but double-check here.
404 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
405 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
407 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
408 "With Verlet lists and PME we should have rcoulomb>=rvdw");
412 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
413 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
415 /* For NVE simulations, we will retain the initial list buffer */
416 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
418 /* Update the Verlet buffer size for the current run setup */
420 /* Here we assume SIMD-enabled kernels are being used. But as currently
421 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
422 * and 4x2 gives a larger buffer than 4x4, this is ok.
424 ListSetupType listType =
425 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
426 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
428 const real rlist_new =
429 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
431 if (rlist_new != ir->rlist)
433 if (fplog != nullptr)
436 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
437 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
439 ir->rlist = rlist_new;
443 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
445 gmx_fatal(FARGS, "Can not set nstlist without %s",
446 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
449 if (EI_DYNAMICS(ir->eI))
451 /* Set or try nstlist values */
452 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
456 /*! \brief Override the nslist value in inputrec
458 * with value passed on the command line (if any)
460 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
464 /* override with anything else than the default -2 */
465 if (nsteps_cmdline > -2)
467 char sbuf_steps[STEPSTRSIZE];
468 char sbuf_msg[STRLEN];
470 ir->nsteps = nsteps_cmdline;
471 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
474 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
475 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
479 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
480 gmx_step_str(nsteps_cmdline, sbuf_steps));
483 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
485 else if (nsteps_cmdline < -2)
487 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
489 /* Do nothing if nsteps_cmdline == -2 */
495 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
497 * If not, and if a warning may be issued, logs a warning about
498 * falling back to CPU code. With thread-MPI, only the first
499 * call to this function should have \c issueWarning true. */
500 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
502 bool gpuIsUseful = true;
505 if (ir.opts.ngener - ir.nwall > 1)
507 /* The GPU code does not support more than one energy group.
508 * If the user requested GPUs explicitly, a fatal error is given later.
512 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
513 "For better performance, run on the GPU without energy groups and then do "
514 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
520 warning = "TPI is not implemented for GPUs.";
523 if (!gpuIsUseful && issueWarning)
525 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
531 //! Initializes the logger for mdrun.
532 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
534 gmx::LoggerBuilder builder;
535 if (fplog != nullptr)
537 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
539 if (isSimulationMasterRank)
541 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
543 return builder.build();
546 //! Make a TaskTarget from an mdrun argument string.
547 static TaskTarget findTaskTarget(const char* optionString)
549 TaskTarget returnValue = TaskTarget::Auto;
551 if (strncmp(optionString, "auto", 3) == 0)
553 returnValue = TaskTarget::Auto;
555 else if (strncmp(optionString, "cpu", 3) == 0)
557 returnValue = TaskTarget::Cpu;
559 else if (strncmp(optionString, "gpu", 3) == 0)
561 returnValue = TaskTarget::Gpu;
565 GMX_ASSERT(false, "Option string should have been checked for sanity already");
571 //! Finish run, aggregate data to print performance info.
572 static void finish_run(FILE* fplog,
573 const gmx::MDLogger& mdlog,
575 const t_inputrec* inputrec,
577 gmx_wallcycle_t wcycle,
578 gmx_walltime_accounting_t walltime_accounting,
579 nonbonded_verlet_t* nbv,
580 const gmx_pme_t* pme,
584 double nbfs = 0, mflop = 0;
585 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
586 elapsed_time_over_all_threads_over_all_ranks;
587 /* Control whether it is valid to print a report. Only the
588 simulation master may print, but it should not do so if the run
589 terminated e.g. before a scheduled reset step. This is
590 complicated by the fact that PME ranks are unaware of the
591 reason why they were sent a pmerecvqxFINISH. To avoid
592 communication deadlocks, we always do the communication for the
593 report, even if we've decided not to write the report, because
594 how long it takes to finish the run is not important when we've
595 decided not to report on the simulation performance.
597 Further, we only report performance for dynamical integrators,
598 because those are the only ones for which we plan to
599 consider doing any optimizations. */
600 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
602 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
604 GMX_LOG(mdlog.warning)
606 .appendText("Simulation ended prematurely, no performance report will be written.");
611 std::unique_ptr<t_nrnb> nrnbTotalStorage;
614 nrnbTotalStorage = std::make_unique<t_nrnb>();
615 nrnb_tot = nrnbTotalStorage.get();
617 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
625 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
626 elapsed_time_over_all_threads =
627 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
631 /* reduce elapsed_time over all MPI ranks in the current simulation */
632 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
634 elapsed_time_over_all_ranks /= cr->nnodes;
635 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
636 * current simulation. */
637 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
638 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
643 elapsed_time_over_all_ranks = elapsed_time;
644 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
649 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
652 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
654 print_dd_statistics(cr, inputrec, fplog);
657 /* TODO Move the responsibility for any scaling by thread counts
658 * to the code that handled the thread region, so that there's a
659 * mechanism to keep cycle counting working during the transition
660 * to task parallelism. */
661 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
662 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
663 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
664 nthreads_pp, nthreads_pme);
665 auto cycle_sum(wallcycle_sum(cr, wcycle));
669 auto nbnxn_gpu_timings =
670 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
671 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
673 if (pme_gpu_task_enabled(pme))
675 pme_gpu_get_timings(pme, &pme_gpu_timings);
677 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
678 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
681 if (EI_DYNAMICS(inputrec->eI))
683 delta_t = inputrec->delta_t;
688 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
689 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
690 delta_t, nbfs, mflop);
694 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
695 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
696 delta_t, nbfs, mflop);
701 int Mdrunner::mdrunner()
704 t_forcerec* fr = nullptr;
705 t_fcdata* fcd = nullptr;
706 real ewaldcoeff_q = 0;
707 real ewaldcoeff_lj = 0;
708 int nChargePerturbed = -1, nTypePerturbed = 0;
709 gmx_wallcycle_t wcycle;
710 gmx_walltime_accounting_t walltime_accounting = nullptr;
711 gmx_membed_t* membed = nullptr;
712 gmx_hw_info_t* hwinfo = nullptr;
714 /* CAUTION: threads may be started later on in this function, so
715 cr doesn't reflect the final parallel state right now */
718 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
719 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
720 const bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
721 const bool doRerun = mdrunOptions.rerun;
723 // Handle task-assignment related user options.
724 EmulateGpuNonbonded emulateGpuNonbonded =
725 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
727 std::vector<int> userGpuTaskAssignment;
730 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
732 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
733 auto nonbondedTarget = findTaskTarget(nbpu_opt);
734 auto pmeTarget = findTaskTarget(pme_opt);
735 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
736 auto bondedTarget = findTaskTarget(bonded_opt);
737 auto updateTarget = findTaskTarget(update_opt);
739 FILE* fplog = nullptr;
740 // If we are appending, we don't write log output because we need
741 // to check that the old log file matches what the checkpoint file
742 // expects. Otherwise, we should start to write log output now if
743 // there is a file ready for it.
744 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
746 fplog = gmx_fio_getfp(logFileHandle);
748 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
749 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
750 gmx::MDLogger mdlog(logOwner.logger());
752 // TODO The thread-MPI master rank makes a working
753 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
754 // after the threads have been launched. This works because no use
755 // is made of that communicator until after the execution paths
756 // have rejoined. But it is likely that we can improve the way
757 // this is expressed, e.g. by expressly running detection only the
758 // master rank for thread-MPI, rather than relying on the mutex
759 // and reference count.
760 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
761 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
763 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
765 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
767 // Print citation requests after all software/hardware printing
768 pleaseCiteGromacs(fplog);
770 // TODO Replace this by unique_ptr once t_inputrec is C++
771 t_inputrec inputrecInstance;
772 t_inputrec* inputrec = nullptr;
773 std::unique_ptr<t_state> globalState;
775 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
777 if (isSimulationMasterRank)
779 /* Only the master rank has the global state */
780 globalState = std::make_unique<t_state>();
782 /* Read (nearly) all data required for the simulation
783 * and keep the partly serialized tpr contents to send to other ranks later
785 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
786 &inputrecInstance, globalState.get(), &mtop);
787 inputrec = &inputrecInstance;
790 /* Check and update the hardware options for internal consistency */
791 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
794 if (GMX_THREAD_MPI && isSimulationMasterRank)
796 bool useGpuForNonbonded = false;
797 bool useGpuForPme = false;
800 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
802 // If the user specified the number of ranks, then we must
803 // respect that, but in default mode, we need to allow for
804 // the number of GPUs to choose the number of ranks.
805 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
806 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
807 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
808 canUseGpuForNonbonded,
809 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
810 hw_opt.nthreads_tmpi);
811 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
812 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
813 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
815 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
817 /* Determine how many thread-MPI ranks to start.
819 * TODO Over-writing the user-supplied value here does
820 * prevent any possible subsequent checks from working
822 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
823 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
825 // Now start the threads for thread MPI.
826 spawnThreads(hw_opt.nthreads_tmpi);
827 // The spawned threads enter mdrunner() and execution of
828 // master and spawned threads joins at the end of this block.
829 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
832 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
833 CommrecHandle crHandle = init_commrec(communicator, ms);
834 t_commrec* cr = crHandle.get();
835 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
839 /* now broadcast everything to the non-master nodes/threads: */
840 if (!isSimulationMasterRank)
842 inputrec = &inputrecInstance;
844 init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
846 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
847 partialDeserializedTpr.reset(nullptr);
849 // Now the number of ranks is known to all ranks, and each knows
850 // the inputrec read by the master rank. The ranks can now all run
851 // the task-deciding functions and will agree on the result
852 // without needing to communicate.
853 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
855 // Note that these variables describe only their own node.
857 // Note that when bonded interactions run on a GPU they always run
858 // alongside a nonbonded task, so do not influence task assignment
859 // even though they affect the force calculation workload.
860 bool useGpuForNonbonded = false;
861 bool useGpuForPme = false;
862 bool useGpuForBonded = false;
863 bool useGpuForUpdate = false;
864 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
867 // It's possible that there are different numbers of GPUs on
868 // different nodes, which is the user's responsibility to
869 // handle. If unsuitable, we will notice that during task
871 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
872 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
873 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
874 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
875 useGpuForPme = decideWhetherToUseGpusForPme(
876 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
877 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
878 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
879 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
880 useGpuForBonded = decideWhetherToUseGpusForBonded(
881 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
882 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
883 domdecOptions.numPmeRanks, gpusWereDetected);
885 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
887 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
889 // Initialize development feature flags that enabled by environment variable
890 // and report those features that are enabled.
891 const DevelopmentFeatureFlags devFlags =
892 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
894 const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
895 false, inputrec, doRerun, mtop, ms, replExParams, nullptr, doEssentialDynamics, doMembed);
896 const bool useModularSimulator = inputIsCompatibleWithModularSimulator
897 && !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
900 // TODO: hide restraint implementation details from Mdrunner.
901 // There is nothing unique about restraints at this point as far as the
902 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
903 // factory functions from the SimulationContext on which to call mdModules_->add().
904 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
905 for (auto&& restraint : restraintManager_->getRestraints())
907 auto module = RestraintMDModule::create(restraint, restraint->sites());
908 mdModules_->add(std::move(module));
911 // TODO: Error handling
912 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
913 // now that the MdModules know their options, they know which callbacks to sign up to
914 mdModules_->subscribeToSimulationSetupNotifications();
915 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
917 if (inputrec->internalParameters != nullptr)
919 mdModulesNotifier.notify(*inputrec->internalParameters);
922 if (fplog != nullptr)
924 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
925 fprintf(fplog, "\n");
930 /* In rerun, set velocities to zero if present */
931 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
933 // rerun does not use velocities
937 "Rerun trajectory contains velocities. Rerun does only evaluate "
938 "potential energy and forces. The velocities will be ignored.");
939 for (int i = 0; i < globalState->natoms; i++)
941 clear_rvec(globalState->v[i]);
943 globalState->flags &= ~(1 << estV);
946 /* now make sure the state is initialized and propagated */
947 set_state_entries(globalState.get(), inputrec, useModularSimulator);
950 /* NM and TPI parallelize over force/energy calculations, not atoms,
951 * so we need to initialize and broadcast the global state.
953 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
957 globalState = std::make_unique<t_state>();
959 broadcastStateWithoutDynamics(cr, globalState.get());
962 /* A parallel command line option consistency check that we can
963 only do after any threads have started. */
965 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
966 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
969 "The -dd or -npme option request a parallel simulation, "
971 "but %s was compiled without threads or MPI enabled",
972 output_env_get_program_display_name(oenv));
974 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
976 "but %s was not started through mpirun/mpiexec or only one rank was requested "
977 "through mpirun/mpiexec",
978 output_env_get_program_display_name(oenv));
982 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
985 "The .mdp file specified an energy mininization or normal mode algorithm, and "
986 "these are not compatible with mdrun -rerun");
989 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
991 if (domdecOptions.numPmeRanks > 0)
993 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
994 "PME-only ranks are requested, but the system does not use PME "
995 "for electrostatics or LJ");
998 domdecOptions.numPmeRanks = 0;
1001 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1003 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1004 * improve performance with many threads per GPU, since our OpenMP
1005 * scaling is bad, but it's difficult to automate the setup.
1007 domdecOptions.numPmeRanks = 0;
1011 if (domdecOptions.numPmeRanks < 0)
1013 domdecOptions.numPmeRanks = 0;
1014 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1018 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1019 "PME GPU decomposition is not supported");
1026 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1030 /* NMR restraints must be initialized before load_checkpoint,
1031 * since with time averaging the history is added to t_state.
1032 * For proper consistency check we therefore need to extend
1034 * So the PME-only nodes (if present) will also initialize
1035 * the distance restraints.
1039 /* This needs to be called before read_checkpoint to extend the state */
1040 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1042 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1044 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1046 ObservablesHistory observablesHistory = {};
1048 if (startingBehavior != StartingBehavior::NewSimulation)
1050 /* Check if checkpoint file exists before doing continuation.
1051 * This way we can use identical input options for the first and subsequent runs...
1053 if (mdrunOptions.numStepsCommandline > -2)
1055 /* Temporarily set the number of steps to unlimited to avoid
1056 * triggering the nsteps check in load_checkpoint().
1057 * This hack will go away soon when the -nsteps option is removed.
1059 inputrec->nsteps = -1;
1062 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1063 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1064 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1066 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1068 // Now we can start normal logging to the truncated log file.
1069 fplog = gmx_fio_getfp(logFileHandle);
1070 prepareLogAppending(fplog);
1071 logOwner = buildLogger(fplog, MASTER(cr));
1072 mdlog = logOwner.logger();
1076 if (mdrunOptions.numStepsCommandline > -2)
1081 "The -nsteps functionality is deprecated, and may be removed in a future "
1083 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1086 /* override nsteps with value set on the commandline */
1087 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1091 copy_mat(globalState->box, box);
1096 gmx_bcast(sizeof(box), box, cr);
1099 if (inputrec->cutoff_scheme != ecutsVERLET)
1102 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1103 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1105 /* Update rlist and nstlist. */
1106 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1107 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1110 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1111 // This builder is necessary while we have multi-part construction
1112 // of DD. Before DD is constructed, we use the existence of
1113 // the builder object to indicate that further construction of DD
1115 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1116 if (useDomainDecomposition)
1118 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1119 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1120 positionsFromStatePointer(globalState.get()));
1124 /* PME, if used, is done on all nodes with 1D decomposition */
1126 cr->duty = (DUTY_PP | DUTY_PME);
1128 if (inputrec->pbcType == PbcType::Screw)
1130 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1134 // Produce the task assignment for this rank.
1135 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1136 GpuTaskAssignments gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1137 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1138 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1139 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1140 // TODO cr->duty & DUTY_PME should imply that a PME
1141 // algorithm is active, but currently does not.
1142 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1144 // Get the device handles for the modules, nullptr when no task is assigned.
1145 DeviceInformation* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1146 DeviceInformation* pmeDeviceInfo = gpuTaskAssignments.initPmeDevice();
1148 // TODO Initialize GPU streams here.
1150 // TODO Currently this is always built, yet DD partition code
1151 // checks if it is built before using it. Probably it should
1152 // become an MDModule that is made only when another module
1153 // requires it (e.g. pull, CompEl, density fitting), so that we
1154 // don't update the local atom sets unilaterally every step.
1155 LocalAtomSetManager atomSets;
1158 // TODO Pass the GPU streams to ddBuilder to use in buffer
1159 // transfers (e.g. halo exchange)
1160 cr->dd = ddBuilder->build(&atomSets);
1161 // The builder's job is done, so destruct it
1162 ddBuilder.reset(nullptr);
1163 // Note that local state still does not exist yet.
1166 // The GPU update is decided here because we need to know whether the constraints or
1167 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1168 // defined). This is only known after DD is initialized, hence decision on using GPU
1169 // update is done so late.
1172 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1174 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1175 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1176 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1177 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1178 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1180 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1182 const bool printHostName = (cr->nnodes > 1);
1183 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1185 // If the user chose a task assignment, give them some hints
1186 // where appropriate.
1187 if (!userGpuTaskAssignment.empty())
1189 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1194 /* After possible communicator splitting in make_dd_communicators.
1195 * we can set up the intra/inter node communication.
1197 gmx_setup_nodecomm(fplog, cr);
1203 GMX_LOG(mdlog.warning)
1205 .appendTextFormatted(
1206 "This is simulation %d out of %d running as a composite GROMACS\n"
1207 "multi-simulation job. Setup for this simulation:\n",
1210 GMX_LOG(mdlog.warning)
1211 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1213 cr->nnodes == 1 ? "thread" : "threads"
1215 cr->nnodes == 1 ? "process" : "processes"
1221 // If mdrun -pin auto honors any affinity setting that already
1222 // exists. If so, it is nice to provide feedback about whether
1223 // that existing affinity setting was from OpenMP or something
1224 // else, so we run this code both before and after we initialize
1225 // the OpenMP support.
1226 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1227 /* Check and update the number of OpenMP threads requested */
1228 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1229 pmeRunMode, mtop, *inputrec);
1231 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1232 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1234 // Enable FP exception detection, but not in
1235 // Release mode and not for compilers with known buggy FP
1236 // exception support (clang with any optimization) or suspected
1237 // buggy FP exception support (gcc 7.* with optimization).
1238 #if !defined NDEBUG \
1239 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1240 && defined __OPTIMIZE__)
1241 const bool bEnableFPE = true;
1243 const bool bEnableFPE = false;
1245 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1248 gmx_feenableexcept();
1251 /* Now that we know the setup is consistent, check for efficiency */
1252 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1253 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1255 /* getting number of PP/PME threads on this MPI / tMPI rank.
1256 PME: env variable should be read only on one node to make sure it is
1257 identical everywhere;
1259 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1260 : gmx_omp_nthreads_get(emntPME);
1261 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1262 physicalNodeComm, mdlog);
1264 // Enable Peer access between GPUs where available
1265 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1266 // any of the GPU communication features are active.
1267 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1268 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1270 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1273 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1275 /* Before setting affinity, check whether the affinity has changed
1276 * - which indicates that probably the OpenMP library has changed it
1277 * since we first checked).
1279 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1281 int numThreadsOnThisNode, intraNodeThreadOffset;
1282 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1283 &intraNodeThreadOffset);
1285 /* Set the CPU affinity */
1286 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1287 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1290 if (mdrunOptions.timingOptions.resetStep > -1)
1295 "The -resetstep functionality is deprecated, and may be removed in a "
1298 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1302 /* Master synchronizes its value of reset_counters with all nodes
1303 * including PME only nodes */
1304 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1305 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1306 wcycle_set_reset_counters(wcycle, reset_counters);
1309 // Membrane embedding must be initialized before we call init_forcerec()
1314 fprintf(stderr, "Initializing membed");
1316 /* Note that membed cannot work in parallel because mtop is
1317 * changed here. Fix this if we ever want to make it run with
1318 * multiple ranks. */
1319 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1320 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1323 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1324 std::unique_ptr<MDAtoms> mdAtoms;
1325 std::unique_ptr<gmx_vsite_t> vsite;
1326 std::unique_ptr<GpuBonded> gpuBonded;
1329 if (thisRankHasDuty(cr, DUTY_PP))
1331 mdModulesNotifier.notify(*cr);
1332 mdModulesNotifier.notify(&atomSets);
1333 mdModulesNotifier.notify(inputrec->pbcType);
1334 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1335 /* Initiate forcerecord */
1336 fr = new t_forcerec;
1337 fr->forceProviders = mdModules_->initForceProviders();
1338 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1339 opt2fn("-table", filenames.size(), filenames.data()),
1340 opt2fn("-tablep", filenames.size(), filenames.data()),
1341 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1343 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1345 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(cr->mpi_comm_mysim, cr->dd->pme_nodeid);
1348 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, nonbondedDeviceInfo,
1349 &mtop, box, wcycle);
1350 if (useGpuForBonded)
1352 auto stream = havePPDomainDecomposition(cr)
1353 ? Nbnxm::gpu_get_command_stream(
1354 fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1355 : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1356 gmx::InteractionLocality::Local);
1357 gpuBonded = std::make_unique<GpuBonded>(mtop.ffparams, stream, wcycle);
1358 fr->gpuBonded = gpuBonded.get();
1361 /* Initialize the mdAtoms structure.
1362 * mdAtoms is not filled with atom data,
1363 * as this can not be done now with domain decomposition.
1365 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1366 if (globalState && thisRankHasPmeGpuTask)
1368 // The pinning of coordinates in the global state object works, because we only use
1369 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1370 // points to the global state object without DD.
1371 // FIXME: MD and EM separately set up the local state - this should happen in the same
1372 // function, which should also perform the pinning.
1373 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1376 /* Initialize the virtual site communication */
1377 vsite = initVsite(mtop, cr);
1379 calc_shifts(box, fr->shift_vec);
1381 /* With periodic molecules the charge groups should be whole at start up
1382 * and the virtual sites should not be far from their proper positions.
1384 if (!inputrec->bContinuation && MASTER(cr)
1385 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1387 /* Make molecules whole at start of run */
1388 if (fr->pbcType != PbcType::No)
1390 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1394 /* Correct initial vsite positions are required
1395 * for the initial distribution in the domain decomposition
1396 * and for the initial shell prediction.
1398 constructVsitesGlobal(mtop, globalState->x);
1402 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1404 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1405 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1410 /* This is a PME only node */
1412 GMX_ASSERT(globalState == nullptr,
1413 "We don't need the state on a PME only rank and expect it to be unitialized");
1415 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1416 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1419 gmx_pme_t* sepPmeData = nullptr;
1420 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1421 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1422 "Double-checking that only PME-only ranks have no forcerec");
1423 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1425 // TODO should live in ewald module once its testing is improved
1427 // Later, this program could contain kernels that might be later
1428 // re-used as auto-tuning progresses, or subsequent simulations
1430 PmeGpuProgramStorage pmeGpuProgram;
1431 if (thisRankHasPmeGpuTask)
1433 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1436 /* Initiate PME if necessary,
1437 * either on all nodes or on dedicated PME nodes only. */
1438 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1440 if (mdAtoms && mdAtoms->mdatoms())
1442 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1443 if (EVDW_PME(inputrec->vdwtype))
1445 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1448 if (cr->npmenodes > 0)
1450 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1451 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1452 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1455 if (thisRankHasDuty(cr, DUTY_PME))
1459 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1460 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1461 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1462 nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1464 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1469 if (EI_DYNAMICS(inputrec->eI))
1471 /* Turn on signal handling on all nodes */
1473 * (A user signal from the PME nodes (if any)
1474 * is communicated to the PP nodes.
1476 signal_handler_install();
1479 pull_t* pull_work = nullptr;
1480 if (thisRankHasDuty(cr, DUTY_PP))
1482 /* Assumes uniform use of the number of OpenMP threads */
1483 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1485 if (inputrec->bPull)
1487 /* Initialize pull code */
1488 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1489 inputrec->fepvals->init_lambda);
1490 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1492 initPullHistory(pull_work, &observablesHistory);
1494 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1496 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1500 std::unique_ptr<EnforcedRotation> enforcedRotation;
1503 /* Initialize enforced rotation code */
1505 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1506 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1509 t_swap* swap = nullptr;
1510 if (inputrec->eSwapCoords != eswapNO)
1512 /* Initialize ion swapping code */
1513 swap = init_swapcoords(fplog, inputrec,
1514 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1515 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1516 oenv, mdrunOptions, startingBehavior);
1519 /* Let makeConstraints know whether we have essential dynamics constraints. */
1520 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1521 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1523 /* Energy terms and groups */
1524 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1525 inputrec->fepvals->n_lambda);
1527 /* Kinetic energy data */
1528 gmx_ekindata_t ekind;
1529 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1531 /* Set up interactive MD (IMD) */
1533 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1534 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1535 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1537 if (DOMAINDECOMP(cr))
1539 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1540 /* This call is not included in init_domain_decomposition mainly
1541 * because fr->cginfo_mb is set later.
1543 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1544 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1547 // TODO This is not the right place to manage the lifetime of
1548 // this data structure, but currently it's the easiest way to
1550 MdrunScheduleWorkload runScheduleWork;
1551 // Also populates the simulation constant workload description.
1552 runScheduleWork.simulationWork =
1553 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1554 useGpuForUpdate, devFlags.enableGpuBufferOps,
1555 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1557 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1558 if (gpusWereDetected
1559 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1560 || runScheduleWork.simulationWork.useGpuBufferOps))
1562 const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1563 const void* localStream =
1564 fr->nbv->gpu_nbv != nullptr
1565 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1567 const void* nonLocalStream =
1568 fr->nbv->gpu_nbv != nullptr
1569 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1571 const void* deviceContext = pme_gpu_get_device_context(fr->pmedata);
1572 const int paddingSize = pme_gpu_get_padding_size(fr->pmedata);
1573 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1574 ? GpuApiCallBehavior::Async
1575 : GpuApiCallBehavior::Sync;
1577 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1578 pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
1579 fr->stateGpu = stateGpu.get();
1582 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1583 SimulatorBuilder simulatorBuilder;
1585 // build and run simulator object based on user-input
1586 auto simulator = simulatorBuilder.build(
1587 inputIsCompatibleWithModularSimulator, fplog, cr, ms, mdlog,
1588 static_cast<int>(filenames.size()), filenames.data(), oenv, mdrunOptions,
1589 startingBehavior, vsite.get(), constr.get(),
1590 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1591 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1592 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1593 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1594 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1597 if (fr->pmePpCommGpu)
1599 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1600 fr->pmePpCommGpu.reset();
1603 if (inputrec->bPull)
1605 finish_pull(pull_work);
1607 finish_swapcoords(swap);
1611 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1613 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1614 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1617 wallcycle_stop(wcycle, ewcRUN);
1619 /* Finish up, write some stuff
1620 * if rerunMD, don't write last frame again
1622 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1623 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1625 // clean up cycle counter
1626 wallcycle_destroy(wcycle);
1631 gmx_pme_destroy(pmedata);
1635 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1636 // before we destroy the GPU context(s) in free_gpu().
1637 // Pinned buffers are associated with contexts in CUDA.
1638 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1639 mdAtoms.reset(nullptr);
1640 globalState.reset(nullptr);
1641 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1642 gpuBonded.reset(nullptr);
1643 /* Free pinned buffers in *fr */
1647 if (hwinfo->gpu_info.n_dev > 0)
1649 /* stop the GPU profiler (only CUDA) */
1653 /* With tMPI we need to wait for all ranks to finish deallocation before
1654 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1657 * This is not a concern in OpenCL where we use one context per rank which
1658 * is freed in nbnxn_gpu_free().
1660 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1661 * but it is easier and more futureproof to call it on the whole node.
1663 * Note that this function needs to be called even if GPUs are not used
1664 * in this run because the PME ranks have no knowledge of whether GPUs
1665 * are used or not, but all ranks need to enter the barrier below.
1666 * \todo Remove this physical node barrier after making sure
1667 * that it's not needed anymore (with a shared GPU run).
1671 physicalNodeComm.barrier();
1674 free_gpu(nonbondedDeviceInfo);
1675 free_gpu(pmeDeviceInfo);
1680 free_membed(membed);
1683 /* Does what it says */
1684 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1685 walltime_accounting_destroy(walltime_accounting);
1687 // Ensure log file content is written
1690 gmx_fio_flush(logFileHandle);
1693 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1694 * exceptions were enabled before function was called. */
1697 gmx_fedisableexcept();
1700 auto rc = static_cast<int>(gmx_get_stop_condition());
1703 /* we need to join all threads. The sub-threads join when they
1704 exit this function, but the master thread needs to be told to
1706 if (PAR(cr) && MASTER(cr))
1714 Mdrunner::~Mdrunner()
1716 // Clean up of the Manager.
1717 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1718 // but okay as long as threads synchronize some time before adding or accessing
1719 // a new set of restraints.
1720 if (restraintManager_)
1722 restraintManager_->clear();
1723 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1724 "restraints added during runner life time should be cleared at runner "
1729 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1731 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1732 // Not sure if this should be logged through the md logger or something else,
1733 // but it is helpful to have some sort of INFO level message sent somewhere.
1734 // std::cout << "Registering restraint named " << name << std::endl;
1736 // When multiple restraints are used, it may be wasteful to register them separately.
1737 // Maybe instead register an entire Restraint Manager as a force provider.
1738 restraintManager_->addToSpec(std::move(puller), name);
1741 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1743 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1745 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1746 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1748 class Mdrunner::BuilderImplementation
1751 BuilderImplementation() = delete;
1752 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1753 ~BuilderImplementation();
1755 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1756 real forceWarningThreshold,
1757 StartingBehavior startingBehavior);
1759 void addDomdec(const DomdecOptions& options);
1761 void addVerletList(int nstlist);
1763 void addReplicaExchange(const ReplicaExchangeParameters& params);
1765 void addNonBonded(const char* nbpu_opt);
1767 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1769 void addBondedTaskAssignment(const char* bonded_opt);
1771 void addUpdateTaskAssignment(const char* update_opt);
1773 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1775 void addFilenames(ArrayRef<const t_filenm> filenames);
1777 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1779 void addLogFile(t_fileio* logFileHandle);
1781 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1786 // Default parameters copied from runner.h
1787 // \todo Clarify source(s) of default parameters.
1789 const char* nbpu_opt_ = nullptr;
1790 const char* pme_opt_ = nullptr;
1791 const char* pme_fft_opt_ = nullptr;
1792 const char* bonded_opt_ = nullptr;
1793 const char* update_opt_ = nullptr;
1795 MdrunOptions mdrunOptions_;
1797 DomdecOptions domdecOptions_;
1799 ReplicaExchangeParameters replicaExchangeParameters_;
1801 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1804 //! Multisim communicator handle.
1805 gmx_multisim_t* multiSimulation_;
1807 //! mdrun communicator
1808 MPI_Comm communicator_ = MPI_COMM_NULL;
1810 //! Print a warning if any force is larger than this (in kJ/mol nm).
1811 real forceWarningThreshold_ = -1;
1813 //! Whether the simulation will start afresh, or restart with/without appending.
1814 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1816 //! The modules that comprise the functionality of mdrun.
1817 std::unique_ptr<MDModules> mdModules_;
1819 //! \brief Parallelism information.
1820 gmx_hw_opt_t hardwareOptions_;
1822 //! filename options for simulation.
1823 ArrayRef<const t_filenm> filenames_;
1825 /*! \brief Handle to output environment.
1827 * \todo gmx_output_env_t needs lifetime management.
1829 gmx_output_env_t* outputEnvironment_ = nullptr;
1831 /*! \brief Non-owning handle to MD log file.
1833 * \todo Context should own output facilities for client.
1834 * \todo Improve log file handle management.
1836 * Code managing the FILE* relies on the ability to set it to
1837 * nullptr to check whether the filehandle is valid.
1839 t_fileio* logFileHandle_ = nullptr;
1842 * \brief Builder for simulation stop signal handler.
1844 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1847 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1848 compat::not_null<SimulationContext*> context) :
1849 mdModules_(std::move(mdModules))
1851 communicator_ = context->communicator_;
1852 multiSimulation_ = context->multiSimulation_.get();
1855 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1857 Mdrunner::BuilderImplementation&
1858 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1859 const real forceWarningThreshold,
1860 const StartingBehavior startingBehavior)
1862 mdrunOptions_ = options;
1863 forceWarningThreshold_ = forceWarningThreshold;
1864 startingBehavior_ = startingBehavior;
1868 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1870 domdecOptions_ = options;
1873 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1878 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1880 replicaExchangeParameters_ = params;
1883 Mdrunner Mdrunner::BuilderImplementation::build()
1885 auto newRunner = Mdrunner(std::move(mdModules_));
1887 newRunner.mdrunOptions = mdrunOptions_;
1888 newRunner.pforce = forceWarningThreshold_;
1889 newRunner.startingBehavior = startingBehavior_;
1890 newRunner.domdecOptions = domdecOptions_;
1892 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1893 newRunner.hw_opt = hardwareOptions_;
1895 // No invariant to check. This parameter exists to optionally override other behavior.
1896 newRunner.nstlist_cmdline = nstlist_;
1898 newRunner.replExParams = replicaExchangeParameters_;
1900 newRunner.filenames = filenames_;
1902 newRunner.communicator = communicator_;
1904 // nullptr is a valid value for the multisim handle
1905 newRunner.ms = multiSimulation_;
1907 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1908 // \todo Update sanity checking when output environment has clearly specified invariants.
1909 // Initialization and default values for oenv are not well specified in the current version.
1910 if (outputEnvironment_)
1912 newRunner.oenv = outputEnvironment_;
1916 GMX_THROW(gmx::APIError(
1917 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1920 newRunner.logFileHandle = logFileHandle_;
1924 newRunner.nbpu_opt = nbpu_opt_;
1928 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1931 if (pme_opt_ && pme_fft_opt_)
1933 newRunner.pme_opt = pme_opt_;
1934 newRunner.pme_fft_opt = pme_fft_opt_;
1938 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1943 newRunner.bonded_opt = bonded_opt_;
1947 GMX_THROW(gmx::APIError(
1948 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1953 newRunner.update_opt = update_opt_;
1957 GMX_THROW(gmx::APIError(
1958 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1962 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1964 if (stopHandlerBuilder_)
1966 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1970 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1976 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1978 nbpu_opt_ = nbpu_opt;
1981 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
1984 pme_fft_opt_ = pme_fft_opt;
1987 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1989 bonded_opt_ = bonded_opt;
1992 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
1994 update_opt_ = update_opt;
1997 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
1999 hardwareOptions_ = hardwareOptions;
2002 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2004 filenames_ = filenames;
2007 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2009 outputEnvironment_ = outputEnvironment;
2012 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2014 logFileHandle_ = logFileHandle;
2017 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2019 stopHandlerBuilder_ = std::move(builder);
2022 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2023 compat::not_null<SimulationContext*> context) :
2024 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2028 MdrunnerBuilder::~MdrunnerBuilder() = default;
2030 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2031 real forceWarningThreshold,
2032 const StartingBehavior startingBehavior)
2034 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2038 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2040 impl_->addDomdec(options);
2044 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2046 impl_->addVerletList(nstlist);
2050 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2052 impl_->addReplicaExchange(params);
2056 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2058 impl_->addNonBonded(nbpu_opt);
2062 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2064 // The builder method may become more general in the future, but in this version,
2065 // parameters for PME electrostatics are both required and the only parameters
2067 if (pme_opt && pme_fft_opt)
2069 impl_->addPME(pme_opt, pme_fft_opt);
2074 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2079 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2081 impl_->addBondedTaskAssignment(bonded_opt);
2085 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2087 impl_->addUpdateTaskAssignment(update_opt);
2091 Mdrunner MdrunnerBuilder::build()
2093 return impl_->build();
2096 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2098 impl_->addHardwareOptions(hardwareOptions);
2102 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2104 impl_->addFilenames(filenames);
2108 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2110 impl_->addOutputEnvironment(outputEnvironment);
2114 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2116 impl_->addLogFile(logFileHandle);
2120 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2122 impl_->addStopHandlerBuilder(std::move(builder));
2126 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2128 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;