2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/device_context.h"
77 #include "gromacs/gpu_utils/device_stream_manager.h"
78 #include "gromacs/gpu_utils/gpu_utils.h"
79 #include "gromacs/hardware/cpuinfo.h"
80 #include "gromacs/hardware/detecthardware.h"
81 #include "gromacs/hardware/printhardware.h"
82 #include "gromacs/imd/imd.h"
83 #include "gromacs/listed_forces/disre.h"
84 #include "gromacs/listed_forces/gpubonded.h"
85 #include "gromacs/listed_forces/listed_forces.h"
86 #include "gromacs/listed_forces/orires.h"
87 #include "gromacs/math/functions.h"
88 #include "gromacs/math/utilities.h"
89 #include "gromacs/math/vec.h"
90 #include "gromacs/mdlib/boxdeformation.h"
91 #include "gromacs/mdlib/broadcaststructs.h"
92 #include "gromacs/mdlib/calc_verletbuf.h"
93 #include "gromacs/mdlib/dispersioncorrection.h"
94 #include "gromacs/mdlib/enerdata_utils.h"
95 #include "gromacs/mdlib/force.h"
96 #include "gromacs/mdlib/forcerec.h"
97 #include "gromacs/mdlib/gmx_omp_nthreads.h"
98 #include "gromacs/mdlib/makeconstraints.h"
99 #include "gromacs/mdlib/md_support.h"
100 #include "gromacs/mdlib/mdatoms.h"
101 #include "gromacs/mdlib/sighandler.h"
102 #include "gromacs/mdlib/stophandler.h"
103 #include "gromacs/mdlib/tgroup.h"
104 #include "gromacs/mdlib/updategroups.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdrun/mdmodules.h"
107 #include "gromacs/mdrun/simulationcontext.h"
108 #include "gromacs/mdrunutility/handlerestart.h"
109 #include "gromacs/mdrunutility/logging.h"
110 #include "gromacs/mdrunutility/multisim.h"
111 #include "gromacs/mdrunutility/printtime.h"
112 #include "gromacs/mdrunutility/threadaffinity.h"
113 #include "gromacs/mdtypes/commrec.h"
114 #include "gromacs/mdtypes/enerdata.h"
115 #include "gromacs/mdtypes/fcdata.h"
116 #include "gromacs/mdtypes/forcerec.h"
117 #include "gromacs/mdtypes/group.h"
118 #include "gromacs/mdtypes/inputrec.h"
119 #include "gromacs/mdtypes/interaction_const.h"
120 #include "gromacs/mdtypes/md_enums.h"
121 #include "gromacs/mdtypes/mdatom.h"
122 #include "gromacs/mdtypes/mdrunoptions.h"
123 #include "gromacs/mdtypes/observableshistory.h"
124 #include "gromacs/mdtypes/simulation_workload.h"
125 #include "gromacs/mdtypes/state.h"
126 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
127 #include "gromacs/modularsimulator/modularsimulator.h"
128 #include "gromacs/nbnxm/gpu_data_mgmt.h"
129 #include "gromacs/nbnxm/nbnxm.h"
130 #include "gromacs/nbnxm/pairlist_tuning.h"
131 #include "gromacs/pbcutil/pbc.h"
132 #include "gromacs/pulling/output.h"
133 #include "gromacs/pulling/pull.h"
134 #include "gromacs/pulling/pull_rotation.h"
135 #include "gromacs/restraint/manager.h"
136 #include "gromacs/restraint/restraintmdmodule.h"
137 #include "gromacs/restraint/restraintpotential.h"
138 #include "gromacs/swap/swapcoords.h"
139 #include "gromacs/taskassignment/decidegpuusage.h"
140 #include "gromacs/taskassignment/decidesimulationworkload.h"
141 #include "gromacs/taskassignment/resourcedivision.h"
142 #include "gromacs/taskassignment/taskassignment.h"
143 #include "gromacs/taskassignment/usergpuids.h"
144 #include "gromacs/timing/gpu_timing.h"
145 #include "gromacs/timing/wallcycle.h"
146 #include "gromacs/timing/wallcyclereporting.h"
147 #include "gromacs/topology/mtop_util.h"
148 #include "gromacs/trajectory/trajectoryframe.h"
149 #include "gromacs/utility/basenetwork.h"
150 #include "gromacs/utility/cstringutil.h"
151 #include "gromacs/utility/exceptions.h"
152 #include "gromacs/utility/fatalerror.h"
153 #include "gromacs/utility/filestream.h"
154 #include "gromacs/utility/gmxassert.h"
155 #include "gromacs/utility/gmxmpi.h"
156 #include "gromacs/utility/keyvaluetree.h"
157 #include "gromacs/utility/logger.h"
158 #include "gromacs/utility/loggerbuilder.h"
159 #include "gromacs/utility/mdmodulenotification.h"
160 #include "gromacs/utility/physicalnodecommunicator.h"
161 #include "gromacs/utility/pleasecite.h"
162 #include "gromacs/utility/programcontext.h"
163 #include "gromacs/utility/smalloc.h"
164 #include "gromacs/utility/stringutil.h"
166 #include "isimulator.h"
167 #include "membedholder.h"
168 #include "replicaexchange.h"
169 #include "simulatorbuilder.h"
172 # include "corewrap.h"
179 /*! \brief Manage any development feature flag variables encountered
181 * The use of dev features indicated by environment variables is
182 * logged in order to ensure that runs with such features enabled can
183 * be identified from their log and standard output. Any cross
184 * dependencies are also checked, and if unsatisfied, a fatal error
187 * Note that some development features overrides are applied already here:
188 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
190 * \param[in] mdlog Logger object.
191 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
192 * \param[in] pmeRunMode The PME run mode for this run
193 * \returns The object populated with development feature flags.
195 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
196 const bool useGpuForNonbonded,
197 const PmeRunMode pmeRunMode)
199 DevelopmentFeatureFlags devFlags;
201 // Some builds of GCC 5 give false positive warnings that these
202 // getenv results are ignored when clearly they are used.
203 #pragma GCC diagnostic push
204 #pragma GCC diagnostic ignored "-Wunused-result"
205 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
206 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
207 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
208 devFlags.enableGpuHaloExchange =
209 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
210 devFlags.enableGpuPmePPComm =
211 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
212 #pragma GCC diagnostic pop
214 if (devFlags.enableGpuBufferOps)
216 GMX_LOG(mdlog.warning)
218 .appendTextFormatted(
219 "This run uses the 'GPU buffer ops' feature, enabled by the "
220 "GMX_USE_GPU_BUFFER_OPS environment variable.");
223 if (devFlags.forceGpuUpdateDefault)
225 GMX_LOG(mdlog.warning)
227 .appendTextFormatted(
228 "This run will default to '-update gpu' as requested by the "
229 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
230 "decomposition lacks substantial testing and should be used with caution.");
233 if (devFlags.enableGpuHaloExchange)
235 if (useGpuForNonbonded)
237 if (!devFlags.enableGpuBufferOps)
239 GMX_LOG(mdlog.warning)
241 .appendTextFormatted(
242 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
243 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
244 devFlags.enableGpuBufferOps = true;
246 GMX_LOG(mdlog.warning)
248 .appendTextFormatted(
249 "This run has requested the 'GPU halo exchange' feature, enabled by "
251 "GMX_GPU_DD_COMMS environment variable.");
255 GMX_LOG(mdlog.warning)
257 .appendTextFormatted(
258 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
259 "halo exchange' feature will not be enabled as nonbonded interactions "
260 "are not offloaded.");
261 devFlags.enableGpuHaloExchange = false;
265 if (devFlags.enableGpuPmePPComm)
267 if (pmeRunMode == PmeRunMode::GPU)
269 GMX_LOG(mdlog.warning)
271 .appendTextFormatted(
272 "This run uses the 'GPU PME-PP communications' feature, enabled "
273 "by the GMX_GPU_PME_PP_COMMS environment variable.");
277 std::string clarification;
278 if (pmeRunMode == PmeRunMode::Mixed)
281 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
286 clarification = "PME is not offloaded to the GPU.";
288 GMX_LOG(mdlog.warning)
291 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
292 "'GPU PME-PP communications' feature was not enabled as "
294 devFlags.enableGpuPmePPComm = false;
301 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
303 * Used to ensure that the master thread does not modify mdrunner during copy
304 * on the spawned threads. */
305 static void threadMpiMdrunnerAccessBarrier()
308 MPI_Barrier(MPI_COMM_WORLD);
312 Mdrunner Mdrunner::cloneOnSpawnedThread() const
314 auto newRunner = Mdrunner(std::make_unique<MDModules>());
316 // All runners in the same process share a restraint manager resource because it is
317 // part of the interface to the client code, which is associated only with the
318 // original thread. Handles to the same resources can be obtained by copy.
320 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
323 // Copy members of master runner.
324 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
325 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
326 newRunner.hw_opt = hw_opt;
327 newRunner.filenames = filenames;
329 newRunner.oenv = oenv;
330 newRunner.mdrunOptions = mdrunOptions;
331 newRunner.domdecOptions = domdecOptions;
332 newRunner.nbpu_opt = nbpu_opt;
333 newRunner.pme_opt = pme_opt;
334 newRunner.pme_fft_opt = pme_fft_opt;
335 newRunner.bonded_opt = bonded_opt;
336 newRunner.update_opt = update_opt;
337 newRunner.nstlist_cmdline = nstlist_cmdline;
338 newRunner.replExParams = replExParams;
339 newRunner.pforce = pforce;
340 // Give the spawned thread the newly created valid communicator
341 // for the simulation.
342 newRunner.communicator = MPI_COMM_WORLD;
344 newRunner.startingBehavior = startingBehavior;
345 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
347 threadMpiMdrunnerAccessBarrier();
352 /*! \brief The callback used for running on spawned threads.
354 * Obtains the pointer to the master mdrunner object from the one
355 * argument permitted to the thread-launch API call, copies it to make
356 * a new runner for this thread, reinitializes necessary data, and
357 * proceeds to the simulation. */
358 static void mdrunner_start_fn(const void* arg)
362 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
363 /* copy the arg list to make sure that it's thread-local. This
364 doesn't copy pointed-to items, of course; fnm, cr and fplog
365 are reset in the call below, all others should be const. */
366 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
369 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
373 void Mdrunner::spawnThreads(int numThreadsToLaunch)
376 /* now spawn new threads that start mdrunner_start_fn(), while
377 the main thread returns. Thread affinity is handled later. */
378 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
379 static_cast<const void*>(this))
382 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
385 // Give the master thread the newly created valid communicator for
387 communicator = MPI_COMM_WORLD;
388 threadMpiMdrunnerAccessBarrier();
390 GMX_UNUSED_VALUE(numThreadsToLaunch);
391 GMX_UNUSED_VALUE(mdrunner_start_fn);
397 /*! \brief Initialize variables for Verlet scheme simulation */
398 static void prepare_verlet_scheme(FILE* fplog,
402 const gmx_mtop_t* mtop,
404 bool makeGpuPairList,
405 const gmx::CpuInfo& cpuinfo)
407 // We checked the cut-offs in grompp, but double-check here.
408 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
409 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
411 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
412 "With Verlet lists and PME we should have rcoulomb>=rvdw");
416 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
417 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
419 /* For NVE simulations, we will retain the initial list buffer */
420 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
422 /* Update the Verlet buffer size for the current run setup */
424 /* Here we assume SIMD-enabled kernels are being used. But as currently
425 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
426 * and 4x2 gives a larger buffer than 4x4, this is ok.
428 ListSetupType listType =
429 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
430 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
432 const real rlist_new =
433 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
435 if (rlist_new != ir->rlist)
437 if (fplog != nullptr)
440 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
441 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
443 ir->rlist = rlist_new;
447 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
449 gmx_fatal(FARGS, "Can not set nstlist without %s",
450 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
453 if (EI_DYNAMICS(ir->eI))
455 /* Set or try nstlist values */
456 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
460 /*! \brief Override the nslist value in inputrec
462 * with value passed on the command line (if any)
464 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
468 /* override with anything else than the default -2 */
469 if (nsteps_cmdline > -2)
471 char sbuf_steps[STEPSTRSIZE];
472 char sbuf_msg[STRLEN];
474 ir->nsteps = nsteps_cmdline;
475 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
478 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
479 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
483 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
484 gmx_step_str(nsteps_cmdline, sbuf_steps));
487 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
489 else if (nsteps_cmdline < -2)
491 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
493 /* Do nothing if nsteps_cmdline == -2 */
499 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
501 * If not, and if a warning may be issued, logs a warning about
502 * falling back to CPU code. With thread-MPI, only the first
503 * call to this function should have \c issueWarning true. */
504 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
506 bool gpuIsUseful = true;
509 if (ir.opts.ngener - ir.nwall > 1)
511 /* The GPU code does not support more than one energy group.
512 * If the user requested GPUs explicitly, a fatal error is given later.
516 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
517 "For better performance, run on the GPU without energy groups and then do "
518 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
524 warning = "TPI is not implemented for GPUs.";
527 if (!gpuIsUseful && issueWarning)
529 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
535 //! Initializes the logger for mdrun.
536 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
538 gmx::LoggerBuilder builder;
539 if (fplog != nullptr)
541 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
543 if (isSimulationMasterRank)
545 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
547 return builder.build();
550 //! Make a TaskTarget from an mdrun argument string.
551 static TaskTarget findTaskTarget(const char* optionString)
553 TaskTarget returnValue = TaskTarget::Auto;
555 if (strncmp(optionString, "auto", 3) == 0)
557 returnValue = TaskTarget::Auto;
559 else if (strncmp(optionString, "cpu", 3) == 0)
561 returnValue = TaskTarget::Cpu;
563 else if (strncmp(optionString, "gpu", 3) == 0)
565 returnValue = TaskTarget::Gpu;
569 GMX_ASSERT(false, "Option string should have been checked for sanity already");
575 //! Finish run, aggregate data to print performance info.
576 static void finish_run(FILE* fplog,
577 const gmx::MDLogger& mdlog,
579 const t_inputrec* inputrec,
581 gmx_wallcycle_t wcycle,
582 gmx_walltime_accounting_t walltime_accounting,
583 nonbonded_verlet_t* nbv,
584 const gmx_pme_t* pme,
588 double nbfs = 0, mflop = 0;
589 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
590 elapsed_time_over_all_threads_over_all_ranks;
591 /* Control whether it is valid to print a report. Only the
592 simulation master may print, but it should not do so if the run
593 terminated e.g. before a scheduled reset step. This is
594 complicated by the fact that PME ranks are unaware of the
595 reason why they were sent a pmerecvqxFINISH. To avoid
596 communication deadlocks, we always do the communication for the
597 report, even if we've decided not to write the report, because
598 how long it takes to finish the run is not important when we've
599 decided not to report on the simulation performance.
601 Further, we only report performance for dynamical integrators,
602 because those are the only ones for which we plan to
603 consider doing any optimizations. */
604 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
606 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
608 GMX_LOG(mdlog.warning)
610 .appendText("Simulation ended prematurely, no performance report will be written.");
615 std::unique_ptr<t_nrnb> nrnbTotalStorage;
618 nrnbTotalStorage = std::make_unique<t_nrnb>();
619 nrnb_tot = nrnbTotalStorage.get();
621 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
629 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
630 elapsed_time_over_all_threads =
631 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
635 /* reduce elapsed_time over all MPI ranks in the current simulation */
636 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
638 elapsed_time_over_all_ranks /= cr->nnodes;
639 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
640 * current simulation. */
641 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
642 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
647 elapsed_time_over_all_ranks = elapsed_time;
648 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
653 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
656 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
658 print_dd_statistics(cr, inputrec, fplog);
661 /* TODO Move the responsibility for any scaling by thread counts
662 * to the code that handled the thread region, so that there's a
663 * mechanism to keep cycle counting working during the transition
664 * to task parallelism. */
665 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
666 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
667 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
668 nthreads_pp, nthreads_pme);
669 auto cycle_sum(wallcycle_sum(cr, wcycle));
673 auto nbnxn_gpu_timings =
674 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
675 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
677 if (pme_gpu_task_enabled(pme))
679 pme_gpu_get_timings(pme, &pme_gpu_timings);
681 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
682 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
685 if (EI_DYNAMICS(inputrec->eI))
687 delta_t = inputrec->delta_t;
692 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
693 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
694 delta_t, nbfs, mflop);
698 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
699 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
700 delta_t, nbfs, mflop);
705 int Mdrunner::mdrunner()
708 t_forcerec* fr = nullptr;
709 real ewaldcoeff_q = 0;
710 real ewaldcoeff_lj = 0;
711 int nChargePerturbed = -1, nTypePerturbed = 0;
712 gmx_wallcycle_t wcycle;
713 gmx_walltime_accounting_t walltime_accounting = nullptr;
714 MembedHolder membedHolder(filenames.size(), filenames.data());
715 gmx_hw_info_t* hwinfo = nullptr;
717 /* CAUTION: threads may be started later on in this function, so
718 cr doesn't reflect the final parallel state right now */
721 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
722 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
723 const bool doRerun = mdrunOptions.rerun;
725 // Handle task-assignment related user options.
726 EmulateGpuNonbonded emulateGpuNonbonded =
727 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
729 std::vector<int> userGpuTaskAssignment;
732 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
734 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
735 auto nonbondedTarget = findTaskTarget(nbpu_opt);
736 auto pmeTarget = findTaskTarget(pme_opt);
737 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
738 auto bondedTarget = findTaskTarget(bonded_opt);
739 auto updateTarget = findTaskTarget(update_opt);
741 FILE* fplog = nullptr;
742 // If we are appending, we don't write log output because we need
743 // to check that the old log file matches what the checkpoint file
744 // expects. Otherwise, we should start to write log output now if
745 // there is a file ready for it.
746 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
748 fplog = gmx_fio_getfp(logFileHandle);
750 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
751 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
752 gmx::MDLogger mdlog(logOwner.logger());
754 // TODO The thread-MPI master rank makes a working
755 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
756 // after the threads have been launched. This works because no use
757 // is made of that communicator until after the execution paths
758 // have rejoined. But it is likely that we can improve the way
759 // this is expressed, e.g. by expressly running detection only the
760 // master rank for thread-MPI, rather than relying on the mutex
761 // and reference count.
762 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
763 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
765 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
767 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
769 // Print citation requests after all software/hardware printing
770 pleaseCiteGromacs(fplog);
772 // TODO Replace this by unique_ptr once t_inputrec is C++
773 t_inputrec inputrecInstance;
774 t_inputrec* inputrec = nullptr;
775 std::unique_ptr<t_state> globalState;
777 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
779 if (isSimulationMasterRank)
781 /* Only the master rank has the global state */
782 globalState = std::make_unique<t_state>();
784 /* Read (nearly) all data required for the simulation
785 * and keep the partly serialized tpr contents to send to other ranks later
787 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
788 &inputrecInstance, globalState.get(), &mtop);
789 inputrec = &inputrecInstance;
792 /* Check and update the hardware options for internal consistency */
793 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
796 if (GMX_THREAD_MPI && isSimulationMasterRank)
798 bool useGpuForNonbonded = false;
799 bool useGpuForPme = false;
802 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
804 // If the user specified the number of ranks, then we must
805 // respect that, but in default mode, we need to allow for
806 // the number of GPUs to choose the number of ranks.
807 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
808 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
809 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
810 canUseGpuForNonbonded,
811 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
812 hw_opt.nthreads_tmpi);
813 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
814 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
815 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
817 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
819 /* Determine how many thread-MPI ranks to start.
821 * TODO Over-writing the user-supplied value here does
822 * prevent any possible subsequent checks from working
824 hw_opt.nthreads_tmpi =
825 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
826 inputrec, &mtop, mdlog, membedHolder.doMembed());
828 // Now start the threads for thread MPI.
829 spawnThreads(hw_opt.nthreads_tmpi);
830 // The spawned threads enter mdrunner() and execution of
831 // master and spawned threads joins at the end of this block.
832 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
835 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
836 CommrecHandle crHandle = init_commrec(communicator, ms);
837 t_commrec* cr = crHandle.get();
838 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
842 /* now broadcast everything to the non-master nodes/threads: */
843 if (!isSimulationMasterRank)
845 inputrec = &inputrecInstance;
847 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec, &mtop,
848 partialDeserializedTpr.get());
850 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
851 partialDeserializedTpr.reset(nullptr);
853 // Now the number of ranks is known to all ranks, and each knows
854 // the inputrec read by the master rank. The ranks can now all run
855 // the task-deciding functions and will agree on the result
856 // without needing to communicate.
857 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
859 // Note that these variables describe only their own node.
861 // Note that when bonded interactions run on a GPU they always run
862 // alongside a nonbonded task, so do not influence task assignment
863 // even though they affect the force calculation workload.
864 bool useGpuForNonbonded = false;
865 bool useGpuForPme = false;
866 bool useGpuForBonded = false;
867 bool useGpuForUpdate = false;
868 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
871 // It's possible that there are different numbers of GPUs on
872 // different nodes, which is the user's responsibility to
873 // handle. If unsuitable, we will notice that during task
875 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
876 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
877 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
878 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
879 useGpuForPme = decideWhetherToUseGpusForPme(
880 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
881 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
882 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
883 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
884 useGpuForBonded = decideWhetherToUseGpusForBonded(
885 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
886 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
887 domdecOptions.numPmeRanks, gpusWereDetected);
889 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
891 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
893 // Initialize development feature flags that enabled by environment variable
894 // and report those features that are enabled.
895 const DevelopmentFeatureFlags devFlags =
896 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
898 const bool useModularSimulator =
899 checkUseModularSimulator(false, inputrec, doRerun, mtop, ms, replExParams, nullptr,
900 doEssentialDynamics, membedHolder.doMembed());
903 // TODO: hide restraint implementation details from Mdrunner.
904 // There is nothing unique about restraints at this point as far as the
905 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
906 // factory functions from the SimulationContext on which to call mdModules_->add().
907 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
908 for (auto&& restraint : restraintManager_->getRestraints())
910 auto module = RestraintMDModule::create(restraint, restraint->sites());
911 mdModules_->add(std::move(module));
914 // TODO: Error handling
915 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
916 // now that the MdModules know their options, they know which callbacks to sign up to
917 mdModules_->subscribeToSimulationSetupNotifications();
918 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
920 if (inputrec->internalParameters != nullptr)
922 mdModulesNotifier.notify(*inputrec->internalParameters);
925 if (fplog != nullptr)
927 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
928 fprintf(fplog, "\n");
933 /* In rerun, set velocities to zero if present */
934 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
936 // rerun does not use velocities
940 "Rerun trajectory contains velocities. Rerun does only evaluate "
941 "potential energy and forces. The velocities will be ignored.");
942 for (int i = 0; i < globalState->natoms; i++)
944 clear_rvec(globalState->v[i]);
946 globalState->flags &= ~(1 << estV);
949 /* now make sure the state is initialized and propagated */
950 set_state_entries(globalState.get(), inputrec, useModularSimulator);
953 /* NM and TPI parallelize over force/energy calculations, not atoms,
954 * so we need to initialize and broadcast the global state.
956 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
960 globalState = std::make_unique<t_state>();
962 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
966 /* A parallel command line option consistency check that we can
967 only do after any threads have started. */
969 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
970 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
973 "The -dd or -npme option request a parallel simulation, "
975 "but %s was compiled without threads or MPI enabled",
976 output_env_get_program_display_name(oenv));
978 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
980 "but %s was not started through mpirun/mpiexec or only one rank was requested "
981 "through mpirun/mpiexec",
982 output_env_get_program_display_name(oenv));
986 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
989 "The .mdp file specified an energy mininization or normal mode algorithm, and "
990 "these are not compatible with mdrun -rerun");
993 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
995 if (domdecOptions.numPmeRanks > 0)
997 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
998 "PME-only ranks are requested, but the system does not use PME "
999 "for electrostatics or LJ");
1002 domdecOptions.numPmeRanks = 0;
1005 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1007 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1008 * improve performance with many threads per GPU, since our OpenMP
1009 * scaling is bad, but it's difficult to automate the setup.
1011 domdecOptions.numPmeRanks = 0;
1015 if (domdecOptions.numPmeRanks < 0)
1017 domdecOptions.numPmeRanks = 0;
1018 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1022 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1023 "PME GPU decomposition is not supported");
1030 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1034 /* NMR restraints must be initialized before load_checkpoint,
1035 * since with time averaging the history is added to t_state.
1036 * For proper consistency check we therefore need to extend
1038 * So the PME-only nodes (if present) will also initialize
1039 * the distance restraints.
1042 /* This needs to be called before read_checkpoint to extend the state */
1043 t_disresdata* disresdata;
1044 snew(disresdata, 1);
1045 init_disres(fplog, &mtop, inputrec, DisResRunMode::MDRun, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1046 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1047 globalState.get(), replExParams.exchangeInterval > 0);
1049 t_oriresdata* oriresdata;
1050 snew(oriresdata, 1);
1051 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), oriresdata);
1053 auto deform = prepareBoxDeformation(globalState->box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1054 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1055 cr->mpi_comm_mygroup, *inputrec);
1057 ObservablesHistory observablesHistory = {};
1059 if (startingBehavior != StartingBehavior::NewSimulation)
1061 /* Check if checkpoint file exists before doing continuation.
1062 * This way we can use identical input options for the first and subsequent runs...
1064 if (mdrunOptions.numStepsCommandline > -2)
1066 /* Temporarily set the number of steps to unlimited to avoid
1067 * triggering the nsteps check in load_checkpoint().
1068 * This hack will go away soon when the -nsteps option is removed.
1070 inputrec->nsteps = -1;
1073 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1074 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1075 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1077 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1079 // Now we can start normal logging to the truncated log file.
1080 fplog = gmx_fio_getfp(logFileHandle);
1081 prepareLogAppending(fplog);
1082 logOwner = buildLogger(fplog, MASTER(cr));
1083 mdlog = logOwner.logger();
1087 if (mdrunOptions.numStepsCommandline > -2)
1092 "The -nsteps functionality is deprecated, and may be removed in a future "
1094 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1097 /* override nsteps with value set on the commandline */
1098 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1102 copy_mat(globalState->box, box);
1107 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1110 if (inputrec->cutoff_scheme != ecutsVERLET)
1113 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1114 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1116 /* Update rlist and nstlist. */
1117 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1118 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1119 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1121 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1122 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1125 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1126 // This builder is necessary while we have multi-part construction
1127 // of DD. Before DD is constructed, we use the existence of
1128 // the builder object to indicate that further construction of DD
1130 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1131 if (useDomainDecomposition)
1133 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1134 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1135 positionsFromStatePointer(globalState.get()));
1139 /* PME, if used, is done on all nodes with 1D decomposition */
1140 cr->nnodes = cr->sizeOfDefaultCommunicator;
1141 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1142 cr->nodeid = cr->rankInDefaultCommunicator;
1144 cr->duty = (DUTY_PP | DUTY_PME);
1146 if (inputrec->pbcType == PbcType::Screw)
1148 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1152 // Produce the task assignment for this rank - done after DD is constructed
1153 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1154 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1155 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1156 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1157 // TODO cr->duty & DUTY_PME should imply that a PME
1158 // algorithm is active, but currently does not.
1159 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1161 // Get the device handles for the modules, nullptr when no task is assigned.
1163 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1165 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1166 bool useTiming = true;
1167 if (GMX_GPU == GMX_GPU_CUDA)
1169 /* WARNING: CUDA timings are incorrect with multiple streams.
1170 * This is the main reason why they are disabled by default.
1172 // TODO: Consider turning on by default when we can detect nr of streams.
1173 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1175 else if (GMX_GPU == GMX_GPU_OPENCL)
1177 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1180 // TODO Currently this is always built, yet DD partition code
1181 // checks if it is built before using it. Probably it should
1182 // become an MDModule that is made only when another module
1183 // requires it (e.g. pull, CompEl, density fitting), so that we
1184 // don't update the local atom sets unilaterally every step.
1185 LocalAtomSetManager atomSets;
1188 // TODO Pass the GPU streams to ddBuilder to use in buffer
1189 // transfers (e.g. halo exchange)
1190 cr->dd = ddBuilder->build(&atomSets);
1191 // The builder's job is done, so destruct it
1192 ddBuilder.reset(nullptr);
1193 // Note that local state still does not exist yet.
1196 // The GPU update is decided here because we need to know whether the constraints or
1197 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1198 // defined). This is only known after DD is initialized, hence decision on using GPU
1199 // update is done so late.
1202 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1204 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1205 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1206 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1207 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1208 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1210 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1212 const bool printHostName = (cr->nnodes > 1);
1213 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1215 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1217 if (deviceInfo != nullptr)
1219 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1221 dd_setup_dlb_resource_sharing(cr, deviceId);
1223 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1224 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1225 useGpuForUpdate, useTiming);
1228 // If the user chose a task assignment, give them some hints
1229 // where appropriate.
1230 if (!userGpuTaskAssignment.empty())
1232 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1237 /* After possible communicator splitting in make_dd_communicators.
1238 * we can set up the intra/inter node communication.
1240 gmx_setup_nodecomm(fplog, cr);
1246 GMX_LOG(mdlog.warning)
1248 .appendTextFormatted(
1249 "This is simulation %d out of %d running as a composite GROMACS\n"
1250 "multi-simulation job. Setup for this simulation:\n",
1253 GMX_LOG(mdlog.warning)
1254 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1256 cr->nnodes == 1 ? "thread" : "threads"
1258 cr->nnodes == 1 ? "process" : "processes"
1264 // If mdrun -pin auto honors any affinity setting that already
1265 // exists. If so, it is nice to provide feedback about whether
1266 // that existing affinity setting was from OpenMP or something
1267 // else, so we run this code both before and after we initialize
1268 // the OpenMP support.
1269 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1270 /* Check and update the number of OpenMP threads requested */
1271 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1272 pmeRunMode, mtop, *inputrec);
1274 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1275 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1277 // Enable FP exception detection, but not in
1278 // Release mode and not for compilers with known buggy FP
1279 // exception support (clang with any optimization) or suspected
1280 // buggy FP exception support (gcc 7.* with optimization).
1281 #if !defined NDEBUG \
1282 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1283 && defined __OPTIMIZE__)
1284 const bool bEnableFPE = true;
1286 const bool bEnableFPE = false;
1288 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1291 gmx_feenableexcept();
1294 /* Now that we know the setup is consistent, check for efficiency */
1295 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1296 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1298 /* getting number of PP/PME threads on this MPI / tMPI rank.
1299 PME: env variable should be read only on one node to make sure it is
1300 identical everywhere;
1302 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1303 : gmx_omp_nthreads_get(emntPME);
1304 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1305 physicalNodeComm, mdlog);
1307 // Enable Peer access between GPUs where available
1308 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1309 // any of the GPU communication features are active.
1310 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1311 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1313 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1316 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1318 /* Before setting affinity, check whether the affinity has changed
1319 * - which indicates that probably the OpenMP library has changed it
1320 * since we first checked).
1322 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1324 int numThreadsOnThisNode, intraNodeThreadOffset;
1325 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1326 &intraNodeThreadOffset);
1328 /* Set the CPU affinity */
1329 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1330 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1333 if (mdrunOptions.timingOptions.resetStep > -1)
1338 "The -resetstep functionality is deprecated, and may be removed in a "
1341 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1345 /* Master synchronizes its value of reset_counters with all nodes
1346 * including PME only nodes */
1347 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1348 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1349 wcycle_set_reset_counters(wcycle, reset_counters);
1352 // Membrane embedding must be initialized before we call init_forcerec()
1353 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1354 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1356 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1357 std::unique_ptr<MDAtoms> mdAtoms;
1358 std::unique_ptr<VirtualSitesHandler> vsite;
1359 std::unique_ptr<GpuBonded> gpuBonded;
1362 if (thisRankHasDuty(cr, DUTY_PP))
1364 mdModulesNotifier.notify(*cr);
1365 mdModulesNotifier.notify(&atomSets);
1366 mdModulesNotifier.notify(inputrec->pbcType);
1367 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1368 /* Initiate forcerecord */
1369 fr = new t_forcerec;
1370 fr->forceProviders = mdModules_->initForceProviders();
1371 init_forcerec(fplog, mdlog, fr, inputrec, &mtop, cr, box,
1372 opt2fn("-table", filenames.size(), filenames.data()),
1373 opt2fn("-tablep", filenames.size(), filenames.data()),
1374 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1375 // Dirty hack, for fixing disres and orires should be made mdmodules
1376 fr->listedForces->fcdata().disres = disresdata;
1377 fr->listedForces->fcdata().orires = oriresdata;
1379 // Save a handle to device stream manager to use elsewhere in the code
1380 // TODO: Forcerec is not a correct place to store it.
1381 fr->deviceStreamManager = deviceStreamManager.get();
1383 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1386 deviceStreamManager != nullptr,
1387 "GPU device stream manager should be valid in order to use PME-PP direct "
1390 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1391 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1393 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1394 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1395 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1398 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, useGpuForNonbonded,
1399 deviceStreamManager.get(), &mtop, box, wcycle);
1400 // TODO: Move the logic below to a GPU bonded builder
1401 if (useGpuForBonded)
1403 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1404 "GPU device stream manager should be valid in order to use GPU "
1405 "version of bonded forces.");
1406 gpuBonded = std::make_unique<GpuBonded>(
1407 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1408 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1409 fr->gpuBonded = gpuBonded.get();
1412 /* Initialize the mdAtoms structure.
1413 * mdAtoms is not filled with atom data,
1414 * as this can not be done now with domain decomposition.
1416 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1417 if (globalState && thisRankHasPmeGpuTask)
1419 // The pinning of coordinates in the global state object works, because we only use
1420 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1421 // points to the global state object without DD.
1422 // FIXME: MD and EM separately set up the local state - this should happen in the same
1423 // function, which should also perform the pinning.
1424 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1427 /* Initialize the virtual site communication */
1428 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1430 calc_shifts(box, fr->shift_vec);
1432 /* With periodic molecules the charge groups should be whole at start up
1433 * and the virtual sites should not be far from their proper positions.
1435 if (!inputrec->bContinuation && MASTER(cr)
1436 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1438 /* Make molecules whole at start of run */
1439 if (fr->pbcType != PbcType::No)
1441 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1445 /* Correct initial vsite positions are required
1446 * for the initial distribution in the domain decomposition
1447 * and for the initial shell prediction.
1449 constructVirtualSitesGlobal(mtop, globalState->x);
1453 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1455 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1456 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1461 /* This is a PME only node */
1463 GMX_ASSERT(globalState == nullptr,
1464 "We don't need the state on a PME only rank and expect it to be unitialized");
1466 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1467 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1470 gmx_pme_t* sepPmeData = nullptr;
1471 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1472 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1473 "Double-checking that only PME-only ranks have no forcerec");
1474 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1476 // TODO should live in ewald module once its testing is improved
1478 // Later, this program could contain kernels that might be later
1479 // re-used as auto-tuning progresses, or subsequent simulations
1481 PmeGpuProgramStorage pmeGpuProgram;
1482 if (thisRankHasPmeGpuTask)
1485 (deviceStreamManager != nullptr),
1486 "GPU device stream manager should be initialized in order to use GPU for PME.");
1487 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1488 "GPU device should be initialized in order to use GPU for PME.");
1489 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1492 /* Initiate PME if necessary,
1493 * either on all nodes or on dedicated PME nodes only. */
1494 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1496 if (mdAtoms && mdAtoms->mdatoms())
1498 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1499 if (EVDW_PME(inputrec->vdwtype))
1501 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1504 if (cr->npmenodes > 0)
1506 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1507 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1508 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1511 if (thisRankHasDuty(cr, DUTY_PME))
1515 // TODO: This should be in the builder.
1516 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1517 "Device stream manager should be valid in order to use GPU "
1520 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1521 "GPU PME stream should be valid in order to use GPU version of PME.");
1523 const DeviceContext* deviceContext =
1524 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1525 const DeviceStream* pmeStream =
1526 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1528 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1529 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1530 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1531 nullptr, deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1533 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1538 if (EI_DYNAMICS(inputrec->eI))
1540 /* Turn on signal handling on all nodes */
1542 * (A user signal from the PME nodes (if any)
1543 * is communicated to the PP nodes.
1545 signal_handler_install();
1548 pull_t* pull_work = nullptr;
1549 if (thisRankHasDuty(cr, DUTY_PP))
1551 /* Assumes uniform use of the number of OpenMP threads */
1552 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1554 if (inputrec->bPull)
1556 /* Initialize pull code */
1557 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1558 inputrec->fepvals->init_lambda);
1559 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1561 initPullHistory(pull_work, &observablesHistory);
1563 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1565 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1569 std::unique_ptr<EnforcedRotation> enforcedRotation;
1572 /* Initialize enforced rotation code */
1574 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1575 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1578 t_swap* swap = nullptr;
1579 if (inputrec->eSwapCoords != eswapNO)
1581 /* Initialize ion swapping code */
1582 swap = init_swapcoords(fplog, inputrec,
1583 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1584 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1585 oenv, mdrunOptions, startingBehavior);
1588 /* Let makeConstraints know whether we have essential dynamics constraints. */
1589 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1590 ms, &nrnb, wcycle, fr->bMolPBC);
1592 /* Energy terms and groups */
1593 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1594 inputrec->fepvals->n_lambda);
1596 /* Kinetic energy data */
1597 gmx_ekindata_t ekind;
1598 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1600 /* Set up interactive MD (IMD) */
1602 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1603 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1604 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1606 if (DOMAINDECOMP(cr))
1608 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1609 /* This call is not included in init_domain_decomposition mainly
1610 * because fr->cginfo_mb is set later.
1612 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1613 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1616 // TODO This is not the right place to manage the lifetime of
1617 // this data structure, but currently it's the easiest way to
1619 MdrunScheduleWorkload runScheduleWork;
1620 // Also populates the simulation constant workload description.
1621 runScheduleWork.simulationWork =
1622 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1623 useGpuForUpdate, devFlags.enableGpuBufferOps,
1624 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1626 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1627 if (gpusWereDetected
1628 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1629 || runScheduleWork.simulationWork.useGpuBufferOps))
1631 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1632 ? GpuApiCallBehavior::Async
1633 : GpuApiCallBehavior::Sync;
1634 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1635 "GPU device stream manager should be initialized to use GPU.");
1636 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1637 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1638 fr->stateGpu = stateGpu.get();
1641 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1642 SimulatorBuilder simulatorBuilder;
1644 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1645 simulatorBuilder.add(std::move(membedHolder));
1646 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1647 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1650 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1651 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1652 simulatorBuilder.add(ConstraintsParam(
1653 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1655 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1656 simulatorBuilder.add(
1657 LegacyInput(static_cast<int>(filenames.size()), filenames.data(), inputrec, fr));
1658 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1659 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1660 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1661 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1662 // Todo move to an MDModule
1663 simulatorBuilder.add(IonSwapping(swap));
1664 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1665 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1667 // build and run simulator object based on user-input
1668 auto simulator = simulatorBuilder.build(useModularSimulator);
1671 if (fr->pmePpCommGpu)
1673 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1674 fr->pmePpCommGpu.reset();
1677 if (inputrec->bPull)
1679 finish_pull(pull_work);
1681 finish_swapcoords(swap);
1685 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1687 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1688 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode,
1689 deviceStreamManager.get());
1692 wallcycle_stop(wcycle, ewcRUN);
1694 /* Finish up, write some stuff
1695 * if rerunMD, don't write last frame again
1697 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1698 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1700 // clean up cycle counter
1701 wallcycle_destroy(wcycle);
1703 deviceStreamManager.reset(nullptr);
1707 gmx_pme_destroy(pmedata);
1711 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1712 // before we destroy the GPU context(s) in free_gpu().
1713 // Pinned buffers are associated with contexts in CUDA.
1714 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1715 mdAtoms.reset(nullptr);
1716 globalState.reset(nullptr);
1717 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1718 gpuBonded.reset(nullptr);
1719 /* Free pinned buffers in *fr */
1722 // TODO convert to C++ so we can get rid of these frees
1726 if (hwinfo->gpu_info.n_dev > 0)
1728 /* stop the GPU profiler (only CUDA) */
1732 /* With tMPI we need to wait for all ranks to finish deallocation before
1733 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1736 * This is not a concern in OpenCL where we use one context per rank which
1737 * is freed in nbnxn_gpu_free().
1739 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1740 * but it is easier and more futureproof to call it on the whole node.
1742 * Note that this function needs to be called even if GPUs are not used
1743 * in this run because the PME ranks have no knowledge of whether GPUs
1744 * are used or not, but all ranks need to enter the barrier below.
1745 * \todo Remove this physical node barrier after making sure
1746 * that it's not needed anymore (with a shared GPU run).
1750 physicalNodeComm.barrier();
1753 free_gpu(deviceInfo);
1755 /* Does what it says */
1756 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1757 walltime_accounting_destroy(walltime_accounting);
1759 // Ensure log file content is written
1762 gmx_fio_flush(logFileHandle);
1765 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1766 * exceptions were enabled before function was called. */
1769 gmx_fedisableexcept();
1772 auto rc = static_cast<int>(gmx_get_stop_condition());
1775 /* we need to join all threads. The sub-threads join when they
1776 exit this function, but the master thread needs to be told to
1778 if (PAR(cr) && MASTER(cr))
1786 Mdrunner::~Mdrunner()
1788 // Clean up of the Manager.
1789 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1790 // but okay as long as threads synchronize some time before adding or accessing
1791 // a new set of restraints.
1792 if (restraintManager_)
1794 restraintManager_->clear();
1795 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1796 "restraints added during runner life time should be cleared at runner "
1801 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1803 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1804 // Not sure if this should be logged through the md logger or something else,
1805 // but it is helpful to have some sort of INFO level message sent somewhere.
1806 // std::cout << "Registering restraint named " << name << std::endl;
1808 // When multiple restraints are used, it may be wasteful to register them separately.
1809 // Maybe instead register an entire Restraint Manager as a force provider.
1810 restraintManager_->addToSpec(std::move(puller), name);
1813 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1815 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1817 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1819 class Mdrunner::BuilderImplementation
1822 BuilderImplementation() = delete;
1823 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1824 ~BuilderImplementation();
1826 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1827 real forceWarningThreshold,
1828 StartingBehavior startingBehavior);
1830 void addDomdec(const DomdecOptions& options);
1832 void addVerletList(int nstlist);
1834 void addReplicaExchange(const ReplicaExchangeParameters& params);
1836 void addNonBonded(const char* nbpu_opt);
1838 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1840 void addBondedTaskAssignment(const char* bonded_opt);
1842 void addUpdateTaskAssignment(const char* update_opt);
1844 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1846 void addFilenames(ArrayRef<const t_filenm> filenames);
1848 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1850 void addLogFile(t_fileio* logFileHandle);
1852 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1857 // Default parameters copied from runner.h
1858 // \todo Clarify source(s) of default parameters.
1860 const char* nbpu_opt_ = nullptr;
1861 const char* pme_opt_ = nullptr;
1862 const char* pme_fft_opt_ = nullptr;
1863 const char* bonded_opt_ = nullptr;
1864 const char* update_opt_ = nullptr;
1866 MdrunOptions mdrunOptions_;
1868 DomdecOptions domdecOptions_;
1870 ReplicaExchangeParameters replicaExchangeParameters_;
1872 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1875 //! Multisim communicator handle.
1876 gmx_multisim_t* multiSimulation_;
1878 //! mdrun communicator
1879 MPI_Comm communicator_ = MPI_COMM_NULL;
1881 //! Print a warning if any force is larger than this (in kJ/mol nm).
1882 real forceWarningThreshold_ = -1;
1884 //! Whether the simulation will start afresh, or restart with/without appending.
1885 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1887 //! The modules that comprise the functionality of mdrun.
1888 std::unique_ptr<MDModules> mdModules_;
1890 //! \brief Parallelism information.
1891 gmx_hw_opt_t hardwareOptions_;
1893 //! filename options for simulation.
1894 ArrayRef<const t_filenm> filenames_;
1896 /*! \brief Handle to output environment.
1898 * \todo gmx_output_env_t needs lifetime management.
1900 gmx_output_env_t* outputEnvironment_ = nullptr;
1902 /*! \brief Non-owning handle to MD log file.
1904 * \todo Context should own output facilities for client.
1905 * \todo Improve log file handle management.
1907 * Code managing the FILE* relies on the ability to set it to
1908 * nullptr to check whether the filehandle is valid.
1910 t_fileio* logFileHandle_ = nullptr;
1913 * \brief Builder for simulation stop signal handler.
1915 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1918 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1919 compat::not_null<SimulationContext*> context) :
1920 mdModules_(std::move(mdModules))
1922 communicator_ = context->communicator_;
1923 multiSimulation_ = context->multiSimulation_.get();
1926 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1928 Mdrunner::BuilderImplementation&
1929 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1930 const real forceWarningThreshold,
1931 const StartingBehavior startingBehavior)
1933 mdrunOptions_ = options;
1934 forceWarningThreshold_ = forceWarningThreshold;
1935 startingBehavior_ = startingBehavior;
1939 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1941 domdecOptions_ = options;
1944 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1949 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1951 replicaExchangeParameters_ = params;
1954 Mdrunner Mdrunner::BuilderImplementation::build()
1956 auto newRunner = Mdrunner(std::move(mdModules_));
1958 newRunner.mdrunOptions = mdrunOptions_;
1959 newRunner.pforce = forceWarningThreshold_;
1960 newRunner.startingBehavior = startingBehavior_;
1961 newRunner.domdecOptions = domdecOptions_;
1963 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1964 newRunner.hw_opt = hardwareOptions_;
1966 // No invariant to check. This parameter exists to optionally override other behavior.
1967 newRunner.nstlist_cmdline = nstlist_;
1969 newRunner.replExParams = replicaExchangeParameters_;
1971 newRunner.filenames = filenames_;
1973 newRunner.communicator = communicator_;
1975 // nullptr is a valid value for the multisim handle
1976 newRunner.ms = multiSimulation_;
1978 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1979 // \todo Update sanity checking when output environment has clearly specified invariants.
1980 // Initialization and default values for oenv are not well specified in the current version.
1981 if (outputEnvironment_)
1983 newRunner.oenv = outputEnvironment_;
1987 GMX_THROW(gmx::APIError(
1988 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1991 newRunner.logFileHandle = logFileHandle_;
1995 newRunner.nbpu_opt = nbpu_opt_;
1999 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2002 if (pme_opt_ && pme_fft_opt_)
2004 newRunner.pme_opt = pme_opt_;
2005 newRunner.pme_fft_opt = pme_fft_opt_;
2009 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2014 newRunner.bonded_opt = bonded_opt_;
2018 GMX_THROW(gmx::APIError(
2019 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2024 newRunner.update_opt = update_opt_;
2028 GMX_THROW(gmx::APIError(
2029 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2033 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2035 if (stopHandlerBuilder_)
2037 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2041 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2047 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2049 nbpu_opt_ = nbpu_opt;
2052 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2055 pme_fft_opt_ = pme_fft_opt;
2058 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2060 bonded_opt_ = bonded_opt;
2063 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2065 update_opt_ = update_opt;
2068 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2070 hardwareOptions_ = hardwareOptions;
2073 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2075 filenames_ = filenames;
2078 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2080 outputEnvironment_ = outputEnvironment;
2083 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2085 logFileHandle_ = logFileHandle;
2088 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2090 stopHandlerBuilder_ = std::move(builder);
2093 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2094 compat::not_null<SimulationContext*> context) :
2095 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2099 MdrunnerBuilder::~MdrunnerBuilder() = default;
2101 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2102 real forceWarningThreshold,
2103 const StartingBehavior startingBehavior)
2105 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2109 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2111 impl_->addDomdec(options);
2115 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2117 impl_->addVerletList(nstlist);
2121 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2123 impl_->addReplicaExchange(params);
2127 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2129 impl_->addNonBonded(nbpu_opt);
2133 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2135 // The builder method may become more general in the future, but in this version,
2136 // parameters for PME electrostatics are both required and the only parameters
2138 if (pme_opt && pme_fft_opt)
2140 impl_->addPME(pme_opt, pme_fft_opt);
2145 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2150 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2152 impl_->addBondedTaskAssignment(bonded_opt);
2156 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2158 impl_->addUpdateTaskAssignment(update_opt);
2162 Mdrunner MdrunnerBuilder::build()
2164 return impl_->build();
2167 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2169 impl_->addHardwareOptions(hardwareOptions);
2173 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2175 impl_->addFilenames(filenames);
2179 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2181 impl_->addOutputEnvironment(outputEnvironment);
2185 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2187 impl_->addLogFile(logFileHandle);
2191 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2193 impl_->addStopHandlerBuilder(std::move(builder));
2197 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2199 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;