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/orires.h"
86 #include "gromacs/math/functions.h"
87 #include "gromacs/math/utilities.h"
88 #include "gromacs/math/vec.h"
89 #include "gromacs/mdlib/boxdeformation.h"
90 #include "gromacs/mdlib/broadcaststructs.h"
91 #include "gromacs/mdlib/calc_verletbuf.h"
92 #include "gromacs/mdlib/dispersioncorrection.h"
93 #include "gromacs/mdlib/enerdata_utils.h"
94 #include "gromacs/mdlib/force.h"
95 #include "gromacs/mdlib/forcerec.h"
96 #include "gromacs/mdlib/gmx_omp_nthreads.h"
97 #include "gromacs/mdlib/makeconstraints.h"
98 #include "gromacs/mdlib/md_support.h"
99 #include "gromacs/mdlib/mdatoms.h"
100 #include "gromacs/mdlib/membed.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/nbnxm/gpu_data_mgmt.h"
128 #include "gromacs/nbnxm/nbnxm.h"
129 #include "gromacs/nbnxm/pairlist_tuning.h"
130 #include "gromacs/pbcutil/pbc.h"
131 #include "gromacs/pulling/output.h"
132 #include "gromacs/pulling/pull.h"
133 #include "gromacs/pulling/pull_rotation.h"
134 #include "gromacs/restraint/manager.h"
135 #include "gromacs/restraint/restraintmdmodule.h"
136 #include "gromacs/restraint/restraintpotential.h"
137 #include "gromacs/swap/swapcoords.h"
138 #include "gromacs/taskassignment/decidegpuusage.h"
139 #include "gromacs/taskassignment/decidesimulationworkload.h"
140 #include "gromacs/taskassignment/resourcedivision.h"
141 #include "gromacs/taskassignment/taskassignment.h"
142 #include "gromacs/taskassignment/usergpuids.h"
143 #include "gromacs/timing/gpu_timing.h"
144 #include "gromacs/timing/wallcycle.h"
145 #include "gromacs/timing/wallcyclereporting.h"
146 #include "gromacs/topology/mtop_util.h"
147 #include "gromacs/trajectory/trajectoryframe.h"
148 #include "gromacs/utility/basenetwork.h"
149 #include "gromacs/utility/cstringutil.h"
150 #include "gromacs/utility/exceptions.h"
151 #include "gromacs/utility/fatalerror.h"
152 #include "gromacs/utility/filestream.h"
153 #include "gromacs/utility/gmxassert.h"
154 #include "gromacs/utility/gmxmpi.h"
155 #include "gromacs/utility/keyvaluetree.h"
156 #include "gromacs/utility/logger.h"
157 #include "gromacs/utility/loggerbuilder.h"
158 #include "gromacs/utility/mdmodulenotification.h"
159 #include "gromacs/utility/physicalnodecommunicator.h"
160 #include "gromacs/utility/pleasecite.h"
161 #include "gromacs/utility/programcontext.h"
162 #include "gromacs/utility/smalloc.h"
163 #include "gromacs/utility/stringutil.h"
165 #include "isimulator.h"
166 #include "replicaexchange.h"
167 #include "simulatorbuilder.h"
170 # include "corewrap.h"
177 /*! \brief Manage any development feature flag variables encountered
179 * The use of dev features indicated by environment variables is
180 * logged in order to ensure that runs with such features enabled can
181 * be identified from their log and standard output. Any cross
182 * dependencies are also checked, and if unsatisfied, a fatal error
185 * Note that some development features overrides are applied already here:
186 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
188 * \param[in] mdlog Logger object.
189 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
190 * \param[in] pmeRunMode The PME run mode for this run
191 * \returns The object populated with development feature flags.
193 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
194 const bool useGpuForNonbonded,
195 const PmeRunMode pmeRunMode)
197 DevelopmentFeatureFlags devFlags;
199 // Some builds of GCC 5 give false positive warnings that these
200 // getenv results are ignored when clearly they are used.
201 #pragma GCC diagnostic push
202 #pragma GCC diagnostic ignored "-Wunused-result"
203 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
204 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
205 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
206 devFlags.enableGpuHaloExchange =
207 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
208 devFlags.enableGpuPmePPComm =
209 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
210 #pragma GCC diagnostic pop
212 if (devFlags.enableGpuBufferOps)
214 GMX_LOG(mdlog.warning)
216 .appendTextFormatted(
217 "This run uses the 'GPU buffer ops' feature, enabled by the "
218 "GMX_USE_GPU_BUFFER_OPS environment variable.");
221 if (devFlags.forceGpuUpdateDefault)
223 GMX_LOG(mdlog.warning)
225 .appendTextFormatted(
226 "This run will default to '-update gpu' as requested by the "
227 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
228 "decomposition lacks substantial testing and should be used with caution.");
231 if (devFlags.enableGpuHaloExchange)
233 if (useGpuForNonbonded)
235 if (!devFlags.enableGpuBufferOps)
237 GMX_LOG(mdlog.warning)
239 .appendTextFormatted(
240 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
241 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
242 devFlags.enableGpuBufferOps = true;
244 GMX_LOG(mdlog.warning)
246 .appendTextFormatted(
247 "This run has requested the 'GPU halo exchange' feature, enabled by "
249 "GMX_GPU_DD_COMMS environment variable.");
253 GMX_LOG(mdlog.warning)
255 .appendTextFormatted(
256 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
257 "halo exchange' feature will not be enabled as nonbonded interactions "
258 "are not offloaded.");
259 devFlags.enableGpuHaloExchange = false;
263 if (devFlags.enableGpuPmePPComm)
265 if (pmeRunMode == PmeRunMode::GPU)
267 GMX_LOG(mdlog.warning)
269 .appendTextFormatted(
270 "This run uses the 'GPU PME-PP communications' feature, enabled "
271 "by the GMX_GPU_PME_PP_COMMS environment variable.");
275 std::string clarification;
276 if (pmeRunMode == PmeRunMode::Mixed)
279 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
284 clarification = "PME is not offloaded to the GPU.";
286 GMX_LOG(mdlog.warning)
289 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
290 "'GPU PME-PP communications' feature was not enabled as "
292 devFlags.enableGpuPmePPComm = false;
299 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
301 * Used to ensure that the master thread does not modify mdrunner during copy
302 * on the spawned threads. */
303 static void threadMpiMdrunnerAccessBarrier()
306 MPI_Barrier(MPI_COMM_WORLD);
310 Mdrunner Mdrunner::cloneOnSpawnedThread() const
312 auto newRunner = Mdrunner(std::make_unique<MDModules>());
314 // All runners in the same process share a restraint manager resource because it is
315 // part of the interface to the client code, which is associated only with the
316 // original thread. Handles to the same resources can be obtained by copy.
318 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
321 // Copy members of master runner.
322 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
323 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
324 newRunner.hw_opt = hw_opt;
325 newRunner.filenames = filenames;
327 newRunner.oenv = oenv;
328 newRunner.mdrunOptions = mdrunOptions;
329 newRunner.domdecOptions = domdecOptions;
330 newRunner.nbpu_opt = nbpu_opt;
331 newRunner.pme_opt = pme_opt;
332 newRunner.pme_fft_opt = pme_fft_opt;
333 newRunner.bonded_opt = bonded_opt;
334 newRunner.update_opt = update_opt;
335 newRunner.nstlist_cmdline = nstlist_cmdline;
336 newRunner.replExParams = replExParams;
337 newRunner.pforce = pforce;
338 // Give the spawned thread the newly created valid communicator
339 // for the simulation.
340 newRunner.communicator = MPI_COMM_WORLD;
342 newRunner.startingBehavior = startingBehavior;
343 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
345 threadMpiMdrunnerAccessBarrier();
350 /*! \brief The callback used for running on spawned threads.
352 * Obtains the pointer to the master mdrunner object from the one
353 * argument permitted to the thread-launch API call, copies it to make
354 * a new runner for this thread, reinitializes necessary data, and
355 * proceeds to the simulation. */
356 static void mdrunner_start_fn(const void* arg)
360 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
361 /* copy the arg list to make sure that it's thread-local. This
362 doesn't copy pointed-to items, of course; fnm, cr and fplog
363 are reset in the call below, all others should be const. */
364 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
367 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
371 void Mdrunner::spawnThreads(int numThreadsToLaunch)
374 /* now spawn new threads that start mdrunner_start_fn(), while
375 the main thread returns. Thread affinity is handled later. */
376 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
377 static_cast<const void*>(this))
380 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
383 // Give the master thread the newly created valid communicator for
385 communicator = MPI_COMM_WORLD;
386 threadMpiMdrunnerAccessBarrier();
388 GMX_UNUSED_VALUE(numThreadsToLaunch);
389 GMX_UNUSED_VALUE(mdrunner_start_fn);
395 /*! \brief Initialize variables for Verlet scheme simulation */
396 static void prepare_verlet_scheme(FILE* fplog,
400 const gmx_mtop_t* mtop,
402 bool makeGpuPairList,
403 const gmx::CpuInfo& cpuinfo)
405 // We checked the cut-offs in grompp, but double-check here.
406 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
407 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
409 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
410 "With Verlet lists and PME we should have rcoulomb>=rvdw");
414 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
415 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
417 /* For NVE simulations, we will retain the initial list buffer */
418 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
420 /* Update the Verlet buffer size for the current run setup */
422 /* Here we assume SIMD-enabled kernels are being used. But as currently
423 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
424 * and 4x2 gives a larger buffer than 4x4, this is ok.
426 ListSetupType listType =
427 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
428 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
430 const real rlist_new =
431 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
433 if (rlist_new != ir->rlist)
435 if (fplog != nullptr)
438 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
439 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
441 ir->rlist = rlist_new;
445 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
447 gmx_fatal(FARGS, "Can not set nstlist without %s",
448 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
451 if (EI_DYNAMICS(ir->eI))
453 /* Set or try nstlist values */
454 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
458 /*! \brief Override the nslist value in inputrec
460 * with value passed on the command line (if any)
462 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
466 /* override with anything else than the default -2 */
467 if (nsteps_cmdline > -2)
469 char sbuf_steps[STEPSTRSIZE];
470 char sbuf_msg[STRLEN];
472 ir->nsteps = nsteps_cmdline;
473 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
476 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
477 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
481 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
482 gmx_step_str(nsteps_cmdline, sbuf_steps));
485 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
487 else if (nsteps_cmdline < -2)
489 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
491 /* Do nothing if nsteps_cmdline == -2 */
497 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
499 * If not, and if a warning may be issued, logs a warning about
500 * falling back to CPU code. With thread-MPI, only the first
501 * call to this function should have \c issueWarning true. */
502 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
504 bool gpuIsUseful = true;
507 if (ir.opts.ngener - ir.nwall > 1)
509 /* The GPU code does not support more than one energy group.
510 * If the user requested GPUs explicitly, a fatal error is given later.
514 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
515 "For better performance, run on the GPU without energy groups and then do "
516 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
522 warning = "TPI is not implemented for GPUs.";
525 if (!gpuIsUseful && issueWarning)
527 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
533 //! Initializes the logger for mdrun.
534 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
536 gmx::LoggerBuilder builder;
537 if (fplog != nullptr)
539 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
541 if (isSimulationMasterRank)
543 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
545 return builder.build();
548 //! Make a TaskTarget from an mdrun argument string.
549 static TaskTarget findTaskTarget(const char* optionString)
551 TaskTarget returnValue = TaskTarget::Auto;
553 if (strncmp(optionString, "auto", 3) == 0)
555 returnValue = TaskTarget::Auto;
557 else if (strncmp(optionString, "cpu", 3) == 0)
559 returnValue = TaskTarget::Cpu;
561 else if (strncmp(optionString, "gpu", 3) == 0)
563 returnValue = TaskTarget::Gpu;
567 GMX_ASSERT(false, "Option string should have been checked for sanity already");
573 //! Finish run, aggregate data to print performance info.
574 static void finish_run(FILE* fplog,
575 const gmx::MDLogger& mdlog,
577 const t_inputrec* inputrec,
579 gmx_wallcycle_t wcycle,
580 gmx_walltime_accounting_t walltime_accounting,
581 nonbonded_verlet_t* nbv,
582 const gmx_pme_t* pme,
586 double nbfs = 0, mflop = 0;
587 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
588 elapsed_time_over_all_threads_over_all_ranks;
589 /* Control whether it is valid to print a report. Only the
590 simulation master may print, but it should not do so if the run
591 terminated e.g. before a scheduled reset step. This is
592 complicated by the fact that PME ranks are unaware of the
593 reason why they were sent a pmerecvqxFINISH. To avoid
594 communication deadlocks, we always do the communication for the
595 report, even if we've decided not to write the report, because
596 how long it takes to finish the run is not important when we've
597 decided not to report on the simulation performance.
599 Further, we only report performance for dynamical integrators,
600 because those are the only ones for which we plan to
601 consider doing any optimizations. */
602 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
604 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
606 GMX_LOG(mdlog.warning)
608 .appendText("Simulation ended prematurely, no performance report will be written.");
613 std::unique_ptr<t_nrnb> nrnbTotalStorage;
616 nrnbTotalStorage = std::make_unique<t_nrnb>();
617 nrnb_tot = nrnbTotalStorage.get();
619 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
627 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
628 elapsed_time_over_all_threads =
629 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
633 /* reduce elapsed_time over all MPI ranks in the current simulation */
634 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
636 elapsed_time_over_all_ranks /= cr->nnodes;
637 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
638 * current simulation. */
639 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
640 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
645 elapsed_time_over_all_ranks = elapsed_time;
646 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
651 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
654 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
656 print_dd_statistics(cr, inputrec, fplog);
659 /* TODO Move the responsibility for any scaling by thread counts
660 * to the code that handled the thread region, so that there's a
661 * mechanism to keep cycle counting working during the transition
662 * to task parallelism. */
663 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
664 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
665 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
666 nthreads_pp, nthreads_pme);
667 auto cycle_sum(wallcycle_sum(cr, wcycle));
671 auto nbnxn_gpu_timings =
672 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
673 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
675 if (pme_gpu_task_enabled(pme))
677 pme_gpu_get_timings(pme, &pme_gpu_timings);
679 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
680 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
683 if (EI_DYNAMICS(inputrec->eI))
685 delta_t = inputrec->delta_t;
690 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
691 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
692 delta_t, nbfs, mflop);
696 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
697 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
698 delta_t, nbfs, mflop);
703 int Mdrunner::mdrunner()
706 t_forcerec* fr = nullptr;
707 t_fcdata* fcd = nullptr;
708 real ewaldcoeff_q = 0;
709 real ewaldcoeff_lj = 0;
710 int nChargePerturbed = -1, nTypePerturbed = 0;
711 gmx_wallcycle_t wcycle;
712 gmx_walltime_accounting_t walltime_accounting = nullptr;
713 gmx_membed_t* membed = nullptr;
714 gmx_hw_info_t* hwinfo = nullptr;
716 /* CAUTION: threads may be started later on in this function, so
717 cr doesn't reflect the final parallel state right now */
720 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
721 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
722 const bool doMembed = opt2bSet("-membed", 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 = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
825 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
827 // Now start the threads for thread MPI.
828 spawnThreads(hw_opt.nthreads_tmpi);
829 // The spawned threads enter mdrunner() and execution of
830 // master and spawned threads joins at the end of this block.
831 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
834 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
835 CommrecHandle crHandle = init_commrec(communicator, ms);
836 t_commrec* cr = crHandle.get();
837 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
841 /* now broadcast everything to the non-master nodes/threads: */
842 if (!isSimulationMasterRank)
844 inputrec = &inputrecInstance;
846 init_parallel(cr->mpi_comm_mygroup, MASTER(cr), inputrec, &mtop, partialDeserializedTpr.get());
848 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
849 partialDeserializedTpr.reset(nullptr);
851 // Now the number of ranks is known to all ranks, and each knows
852 // the inputrec read by the master rank. The ranks can now all run
853 // the task-deciding functions and will agree on the result
854 // without needing to communicate.
855 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
857 // Note that these variables describe only their own node.
859 // Note that when bonded interactions run on a GPU they always run
860 // alongside a nonbonded task, so do not influence task assignment
861 // even though they affect the force calculation workload.
862 bool useGpuForNonbonded = false;
863 bool useGpuForPme = false;
864 bool useGpuForBonded = false;
865 bool useGpuForUpdate = false;
866 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
869 // It's possible that there are different numbers of GPUs on
870 // different nodes, which is the user's responsibility to
871 // handle. If unsuitable, we will notice that during task
873 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
874 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
875 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
876 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
877 useGpuForPme = decideWhetherToUseGpusForPme(
878 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
879 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
880 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
881 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
882 useGpuForBonded = decideWhetherToUseGpusForBonded(
883 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
884 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
885 domdecOptions.numPmeRanks, gpusWereDetected);
887 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
889 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
891 // Initialize development feature flags that enabled by environment variable
892 // and report those features that are enabled.
893 const DevelopmentFeatureFlags devFlags =
894 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
896 const bool useModularSimulator = checkUseModularSimulator(
897 false, inputrec, doRerun, mtop, ms, replExParams, nullptr, doEssentialDynamics, doMembed);
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->mpi_comm_mygroup, DOMAINDECOMP(cr), PAR(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, DisResRunMode::MDRun, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1041 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, fcd,
1042 globalState.get(), replExParams.exchangeInterval > 0);
1044 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1046 auto deform = prepareBoxDeformation(globalState->box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1047 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1048 cr->mpi_comm_mygroup, *inputrec);
1050 ObservablesHistory observablesHistory = {};
1052 if (startingBehavior != StartingBehavior::NewSimulation)
1054 /* Check if checkpoint file exists before doing continuation.
1055 * This way we can use identical input options for the first and subsequent runs...
1057 if (mdrunOptions.numStepsCommandline > -2)
1059 /* Temporarily set the number of steps to unlimited to avoid
1060 * triggering the nsteps check in load_checkpoint().
1061 * This hack will go away soon when the -nsteps option is removed.
1063 inputrec->nsteps = -1;
1066 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1067 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1068 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1070 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1072 // Now we can start normal logging to the truncated log file.
1073 fplog = gmx_fio_getfp(logFileHandle);
1074 prepareLogAppending(fplog);
1075 logOwner = buildLogger(fplog, MASTER(cr));
1076 mdlog = logOwner.logger();
1080 if (mdrunOptions.numStepsCommandline > -2)
1085 "The -nsteps functionality is deprecated, and may be removed in a future "
1087 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1090 /* override nsteps with value set on the commandline */
1091 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1095 copy_mat(globalState->box, box);
1100 gmx_bcast(sizeof(box), box, cr->mpi_comm_mygroup);
1103 if (inputrec->cutoff_scheme != ecutsVERLET)
1106 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1107 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1109 /* Update rlist and nstlist. */
1110 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1111 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1114 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1115 // This builder is necessary while we have multi-part construction
1116 // of DD. Before DD is constructed, we use the existence of
1117 // the builder object to indicate that further construction of DD
1119 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1120 if (useDomainDecomposition)
1122 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1123 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1124 positionsFromStatePointer(globalState.get()));
1128 /* PME, if used, is done on all nodes with 1D decomposition */
1130 cr->duty = (DUTY_PP | DUTY_PME);
1132 if (inputrec->pbcType == PbcType::Screw)
1134 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1138 // Produce the task assignment for this rank - done after DD is constructed
1139 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1140 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1141 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1142 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1143 // TODO cr->duty & DUTY_PME should imply that a PME
1144 // algorithm is active, but currently does not.
1145 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1147 // Get the device handles for the modules, nullptr when no task is assigned.
1149 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1151 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1152 bool useTiming = true;
1153 if (GMX_GPU == GMX_GPU_CUDA)
1155 /* WARNING: CUDA timings are incorrect with multiple streams.
1156 * This is the main reason why they are disabled by default.
1158 // TODO: Consider turning on by default when we can detect nr of streams.
1159 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1161 else if (GMX_GPU == GMX_GPU_OPENCL)
1163 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1166 // TODO Currently this is always built, yet DD partition code
1167 // checks if it is built before using it. Probably it should
1168 // become an MDModule that is made only when another module
1169 // requires it (e.g. pull, CompEl, density fitting), so that we
1170 // don't update the local atom sets unilaterally every step.
1171 LocalAtomSetManager atomSets;
1174 // TODO Pass the GPU streams to ddBuilder to use in buffer
1175 // transfers (e.g. halo exchange)
1176 cr->dd = ddBuilder->build(&atomSets);
1177 // The builder's job is done, so destruct it
1178 ddBuilder.reset(nullptr);
1179 // Note that local state still does not exist yet.
1182 // The GPU update is decided here because we need to know whether the constraints or
1183 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1184 // defined). This is only known after DD is initialized, hence decision on using GPU
1185 // update is done so late.
1188 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1190 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1191 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1192 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1193 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1194 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1196 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1198 const bool printHostName = (cr->nnodes > 1);
1199 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1201 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1203 if (deviceInfo != nullptr)
1205 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1207 dd_setup_dlb_resource_sharing(cr, deviceId);
1209 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1210 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1211 useGpuForUpdate, useTiming);
1214 // If the user chose a task assignment, give them some hints
1215 // where appropriate.
1216 if (!userGpuTaskAssignment.empty())
1218 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1223 /* After possible communicator splitting in make_dd_communicators.
1224 * we can set up the intra/inter node communication.
1226 gmx_setup_nodecomm(fplog, cr);
1232 GMX_LOG(mdlog.warning)
1234 .appendTextFormatted(
1235 "This is simulation %d out of %d running as a composite GROMACS\n"
1236 "multi-simulation job. Setup for this simulation:\n",
1239 GMX_LOG(mdlog.warning)
1240 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1242 cr->nnodes == 1 ? "thread" : "threads"
1244 cr->nnodes == 1 ? "process" : "processes"
1250 // If mdrun -pin auto honors any affinity setting that already
1251 // exists. If so, it is nice to provide feedback about whether
1252 // that existing affinity setting was from OpenMP or something
1253 // else, so we run this code both before and after we initialize
1254 // the OpenMP support.
1255 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1256 /* Check and update the number of OpenMP threads requested */
1257 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1258 pmeRunMode, mtop, *inputrec);
1260 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1261 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1263 // Enable FP exception detection, but not in
1264 // Release mode and not for compilers with known buggy FP
1265 // exception support (clang with any optimization) or suspected
1266 // buggy FP exception support (gcc 7.* with optimization).
1267 #if !defined NDEBUG \
1268 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1269 && defined __OPTIMIZE__)
1270 const bool bEnableFPE = true;
1272 const bool bEnableFPE = false;
1274 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1277 gmx_feenableexcept();
1280 /* Now that we know the setup is consistent, check for efficiency */
1281 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1282 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1284 /* getting number of PP/PME threads on this MPI / tMPI rank.
1285 PME: env variable should be read only on one node to make sure it is
1286 identical everywhere;
1288 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1289 : gmx_omp_nthreads_get(emntPME);
1290 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1291 physicalNodeComm, mdlog);
1293 // Enable Peer access between GPUs where available
1294 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1295 // any of the GPU communication features are active.
1296 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1297 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1299 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1302 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1304 /* Before setting affinity, check whether the affinity has changed
1305 * - which indicates that probably the OpenMP library has changed it
1306 * since we first checked).
1308 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1310 int numThreadsOnThisNode, intraNodeThreadOffset;
1311 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1312 &intraNodeThreadOffset);
1314 /* Set the CPU affinity */
1315 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1316 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1319 if (mdrunOptions.timingOptions.resetStep > -1)
1324 "The -resetstep functionality is deprecated, and may be removed in a "
1327 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1331 /* Master synchronizes its value of reset_counters with all nodes
1332 * including PME only nodes */
1333 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1334 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1335 wcycle_set_reset_counters(wcycle, reset_counters);
1338 // Membrane embedding must be initialized before we call init_forcerec()
1343 fprintf(stderr, "Initializing membed");
1345 /* Note that membed cannot work in parallel because mtop is
1346 * changed here. Fix this if we ever want to make it run with
1347 * multiple ranks. */
1348 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1349 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1352 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1353 std::unique_ptr<MDAtoms> mdAtoms;
1354 std::unique_ptr<gmx_vsite_t> vsite;
1355 std::unique_ptr<GpuBonded> gpuBonded;
1358 if (thisRankHasDuty(cr, DUTY_PP))
1360 mdModulesNotifier.notify(*cr);
1361 mdModulesNotifier.notify(&atomSets);
1362 mdModulesNotifier.notify(inputrec->pbcType);
1363 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1364 /* Initiate forcerecord */
1365 fr = new t_forcerec;
1366 fr->forceProviders = mdModules_->initForceProviders();
1367 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1368 opt2fn("-table", filenames.size(), filenames.data()),
1369 opt2fn("-tablep", filenames.size(), filenames.data()),
1370 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1372 // Save a handle to device stream manager to use elsewhere in the code
1373 // TODO: Forcerec is not a correct place to store it.
1374 fr->deviceStreamManager = deviceStreamManager.get();
1376 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1379 deviceStreamManager != nullptr,
1380 "GPU device stream manager should be valid in order to use PME-PP direct "
1383 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1384 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1386 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1387 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1388 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1391 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, useGpuForNonbonded,
1392 deviceStreamManager.get(), &mtop, box, wcycle);
1393 // TODO: Move the logic below to a GPU bonded builder
1394 if (useGpuForBonded)
1396 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1397 "GPU device stream manager should be valid in order to use GPU "
1398 "version of bonded forces.");
1399 gpuBonded = std::make_unique<GpuBonded>(
1400 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1401 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1402 fr->gpuBonded = gpuBonded.get();
1405 /* Initialize the mdAtoms structure.
1406 * mdAtoms is not filled with atom data,
1407 * as this can not be done now with domain decomposition.
1409 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1410 if (globalState && thisRankHasPmeGpuTask)
1412 // The pinning of coordinates in the global state object works, because we only use
1413 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1414 // points to the global state object without DD.
1415 // FIXME: MD and EM separately set up the local state - this should happen in the same
1416 // function, which should also perform the pinning.
1417 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1420 /* Initialize the virtual site communication */
1421 vsite = initVsite(mtop, cr);
1423 calc_shifts(box, fr->shift_vec);
1425 /* With periodic molecules the charge groups should be whole at start up
1426 * and the virtual sites should not be far from their proper positions.
1428 if (!inputrec->bContinuation && MASTER(cr)
1429 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1431 /* Make molecules whole at start of run */
1432 if (fr->pbcType != PbcType::No)
1434 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1438 /* Correct initial vsite positions are required
1439 * for the initial distribution in the domain decomposition
1440 * and for the initial shell prediction.
1442 constructVsitesGlobal(mtop, globalState->x);
1446 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1448 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1449 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1454 /* This is a PME only node */
1456 GMX_ASSERT(globalState == nullptr,
1457 "We don't need the state on a PME only rank and expect it to be unitialized");
1459 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1460 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1463 gmx_pme_t* sepPmeData = nullptr;
1464 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1465 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1466 "Double-checking that only PME-only ranks have no forcerec");
1467 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1469 // TODO should live in ewald module once its testing is improved
1471 // Later, this program could contain kernels that might be later
1472 // re-used as auto-tuning progresses, or subsequent simulations
1474 PmeGpuProgramStorage pmeGpuProgram;
1475 if (thisRankHasPmeGpuTask)
1478 (deviceStreamManager != nullptr),
1479 "GPU device stream manager should be initialized in order to use GPU for PME.");
1480 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1481 "GPU device should be initialized in order to use GPU for PME.");
1482 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1485 /* Initiate PME if necessary,
1486 * either on all nodes or on dedicated PME nodes only. */
1487 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1489 if (mdAtoms && mdAtoms->mdatoms())
1491 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1492 if (EVDW_PME(inputrec->vdwtype))
1494 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1497 if (cr->npmenodes > 0)
1499 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1500 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1501 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1504 if (thisRankHasDuty(cr, DUTY_PME))
1508 // TODO: This should be in the builder.
1509 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1510 "Device stream manager should be valid in order to use GPU "
1513 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1514 "GPU PME stream should be valid in order to use GPU version of PME.");
1516 const DeviceContext* deviceContext =
1517 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1518 const DeviceStream* pmeStream =
1519 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1521 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1522 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1523 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1524 nullptr, deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1526 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1531 if (EI_DYNAMICS(inputrec->eI))
1533 /* Turn on signal handling on all nodes */
1535 * (A user signal from the PME nodes (if any)
1536 * is communicated to the PP nodes.
1538 signal_handler_install();
1541 pull_t* pull_work = nullptr;
1542 if (thisRankHasDuty(cr, DUTY_PP))
1544 /* Assumes uniform use of the number of OpenMP threads */
1545 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1547 if (inputrec->bPull)
1549 /* Initialize pull code */
1550 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1551 inputrec->fepvals->init_lambda);
1552 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1554 initPullHistory(pull_work, &observablesHistory);
1556 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1558 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1562 std::unique_ptr<EnforcedRotation> enforcedRotation;
1565 /* Initialize enforced rotation code */
1567 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1568 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1571 t_swap* swap = nullptr;
1572 if (inputrec->eSwapCoords != eswapNO)
1574 /* Initialize ion swapping code */
1575 swap = init_swapcoords(fplog, inputrec,
1576 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1577 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1578 oenv, mdrunOptions, startingBehavior);
1581 /* Let makeConstraints know whether we have essential dynamics constraints. */
1582 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1583 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1585 /* Energy terms and groups */
1586 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1587 inputrec->fepvals->n_lambda);
1589 /* Kinetic energy data */
1590 gmx_ekindata_t ekind;
1591 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1593 /* Set up interactive MD (IMD) */
1595 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1596 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1597 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1599 if (DOMAINDECOMP(cr))
1601 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1602 /* This call is not included in init_domain_decomposition mainly
1603 * because fr->cginfo_mb is set later.
1605 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1606 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1609 // TODO This is not the right place to manage the lifetime of
1610 // this data structure, but currently it's the easiest way to
1612 MdrunScheduleWorkload runScheduleWork;
1613 // Also populates the simulation constant workload description.
1614 runScheduleWork.simulationWork =
1615 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1616 useGpuForUpdate, devFlags.enableGpuBufferOps,
1617 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1619 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1620 if (gpusWereDetected
1621 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1622 || runScheduleWork.simulationWork.useGpuBufferOps))
1624 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1625 ? GpuApiCallBehavior::Async
1626 : GpuApiCallBehavior::Sync;
1627 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1628 "GPU device stream manager should be initialized to use GPU.");
1629 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1630 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1631 fr->stateGpu = stateGpu.get();
1634 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1635 SimulatorBuilder simulatorBuilder;
1637 // build and run simulator object based on user-input
1638 auto simulator = simulatorBuilder.build(
1639 useModularSimulator, fplog, cr, ms, mdlog, static_cast<int>(filenames.size()),
1640 filenames.data(), oenv, mdrunOptions, startingBehavior, vsite.get(), constr.get(),
1641 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1642 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1643 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1644 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1645 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1648 if (fr->pmePpCommGpu)
1650 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1651 fr->pmePpCommGpu.reset();
1654 if (inputrec->bPull)
1656 finish_pull(pull_work);
1658 finish_swapcoords(swap);
1662 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1664 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1665 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode,
1666 deviceStreamManager.get());
1669 wallcycle_stop(wcycle, ewcRUN);
1671 /* Finish up, write some stuff
1672 * if rerunMD, don't write last frame again
1674 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1675 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1677 // clean up cycle counter
1678 wallcycle_destroy(wcycle);
1680 deviceStreamManager.reset(nullptr);
1684 gmx_pme_destroy(pmedata);
1688 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1689 // before we destroy the GPU context(s) in free_gpu().
1690 // Pinned buffers are associated with contexts in CUDA.
1691 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1692 mdAtoms.reset(nullptr);
1693 globalState.reset(nullptr);
1694 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1695 gpuBonded.reset(nullptr);
1696 /* Free pinned buffers in *fr */
1700 if (hwinfo->gpu_info.n_dev > 0)
1702 /* stop the GPU profiler (only CUDA) */
1706 /* With tMPI we need to wait for all ranks to finish deallocation before
1707 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1710 * This is not a concern in OpenCL where we use one context per rank which
1711 * is freed in nbnxn_gpu_free().
1713 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1714 * but it is easier and more futureproof to call it on the whole node.
1716 * Note that this function needs to be called even if GPUs are not used
1717 * in this run because the PME ranks have no knowledge of whether GPUs
1718 * are used or not, but all ranks need to enter the barrier below.
1719 * \todo Remove this physical node barrier after making sure
1720 * that it's not needed anymore (with a shared GPU run).
1724 physicalNodeComm.barrier();
1727 free_gpu(deviceInfo);
1732 free_membed(membed);
1735 /* Does what it says */
1736 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1737 walltime_accounting_destroy(walltime_accounting);
1739 // Ensure log file content is written
1742 gmx_fio_flush(logFileHandle);
1745 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1746 * exceptions were enabled before function was called. */
1749 gmx_fedisableexcept();
1752 auto rc = static_cast<int>(gmx_get_stop_condition());
1755 /* we need to join all threads. The sub-threads join when they
1756 exit this function, but the master thread needs to be told to
1758 if (PAR(cr) && MASTER(cr))
1766 Mdrunner::~Mdrunner()
1768 // Clean up of the Manager.
1769 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1770 // but okay as long as threads synchronize some time before adding or accessing
1771 // a new set of restraints.
1772 if (restraintManager_)
1774 restraintManager_->clear();
1775 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1776 "restraints added during runner life time should be cleared at runner "
1781 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1783 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1784 // Not sure if this should be logged through the md logger or something else,
1785 // but it is helpful to have some sort of INFO level message sent somewhere.
1786 // std::cout << "Registering restraint named " << name << std::endl;
1788 // When multiple restraints are used, it may be wasteful to register them separately.
1789 // Maybe instead register an entire Restraint Manager as a force provider.
1790 restraintManager_->addToSpec(std::move(puller), name);
1793 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1795 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1797 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1798 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1800 class Mdrunner::BuilderImplementation
1803 BuilderImplementation() = delete;
1804 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1805 ~BuilderImplementation();
1807 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1808 real forceWarningThreshold,
1809 StartingBehavior startingBehavior);
1811 void addDomdec(const DomdecOptions& options);
1813 void addVerletList(int nstlist);
1815 void addReplicaExchange(const ReplicaExchangeParameters& params);
1817 void addNonBonded(const char* nbpu_opt);
1819 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1821 void addBondedTaskAssignment(const char* bonded_opt);
1823 void addUpdateTaskAssignment(const char* update_opt);
1825 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1827 void addFilenames(ArrayRef<const t_filenm> filenames);
1829 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1831 void addLogFile(t_fileio* logFileHandle);
1833 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1838 // Default parameters copied from runner.h
1839 // \todo Clarify source(s) of default parameters.
1841 const char* nbpu_opt_ = nullptr;
1842 const char* pme_opt_ = nullptr;
1843 const char* pme_fft_opt_ = nullptr;
1844 const char* bonded_opt_ = nullptr;
1845 const char* update_opt_ = nullptr;
1847 MdrunOptions mdrunOptions_;
1849 DomdecOptions domdecOptions_;
1851 ReplicaExchangeParameters replicaExchangeParameters_;
1853 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1856 //! Multisim communicator handle.
1857 gmx_multisim_t* multiSimulation_;
1859 //! mdrun communicator
1860 MPI_Comm communicator_ = MPI_COMM_NULL;
1862 //! Print a warning if any force is larger than this (in kJ/mol nm).
1863 real forceWarningThreshold_ = -1;
1865 //! Whether the simulation will start afresh, or restart with/without appending.
1866 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1868 //! The modules that comprise the functionality of mdrun.
1869 std::unique_ptr<MDModules> mdModules_;
1871 //! \brief Parallelism information.
1872 gmx_hw_opt_t hardwareOptions_;
1874 //! filename options for simulation.
1875 ArrayRef<const t_filenm> filenames_;
1877 /*! \brief Handle to output environment.
1879 * \todo gmx_output_env_t needs lifetime management.
1881 gmx_output_env_t* outputEnvironment_ = nullptr;
1883 /*! \brief Non-owning handle to MD log file.
1885 * \todo Context should own output facilities for client.
1886 * \todo Improve log file handle management.
1888 * Code managing the FILE* relies on the ability to set it to
1889 * nullptr to check whether the filehandle is valid.
1891 t_fileio* logFileHandle_ = nullptr;
1894 * \brief Builder for simulation stop signal handler.
1896 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1899 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1900 compat::not_null<SimulationContext*> context) :
1901 mdModules_(std::move(mdModules))
1903 communicator_ = context->communicator_;
1904 multiSimulation_ = context->multiSimulation_.get();
1907 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1909 Mdrunner::BuilderImplementation&
1910 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1911 const real forceWarningThreshold,
1912 const StartingBehavior startingBehavior)
1914 mdrunOptions_ = options;
1915 forceWarningThreshold_ = forceWarningThreshold;
1916 startingBehavior_ = startingBehavior;
1920 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1922 domdecOptions_ = options;
1925 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1930 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1932 replicaExchangeParameters_ = params;
1935 Mdrunner Mdrunner::BuilderImplementation::build()
1937 auto newRunner = Mdrunner(std::move(mdModules_));
1939 newRunner.mdrunOptions = mdrunOptions_;
1940 newRunner.pforce = forceWarningThreshold_;
1941 newRunner.startingBehavior = startingBehavior_;
1942 newRunner.domdecOptions = domdecOptions_;
1944 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1945 newRunner.hw_opt = hardwareOptions_;
1947 // No invariant to check. This parameter exists to optionally override other behavior.
1948 newRunner.nstlist_cmdline = nstlist_;
1950 newRunner.replExParams = replicaExchangeParameters_;
1952 newRunner.filenames = filenames_;
1954 newRunner.communicator = communicator_;
1956 // nullptr is a valid value for the multisim handle
1957 newRunner.ms = multiSimulation_;
1959 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1960 // \todo Update sanity checking when output environment has clearly specified invariants.
1961 // Initialization and default values for oenv are not well specified in the current version.
1962 if (outputEnvironment_)
1964 newRunner.oenv = outputEnvironment_;
1968 GMX_THROW(gmx::APIError(
1969 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1972 newRunner.logFileHandle = logFileHandle_;
1976 newRunner.nbpu_opt = nbpu_opt_;
1980 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1983 if (pme_opt_ && pme_fft_opt_)
1985 newRunner.pme_opt = pme_opt_;
1986 newRunner.pme_fft_opt = pme_fft_opt_;
1990 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1995 newRunner.bonded_opt = bonded_opt_;
1999 GMX_THROW(gmx::APIError(
2000 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2005 newRunner.update_opt = update_opt_;
2009 GMX_THROW(gmx::APIError(
2010 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2014 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2016 if (stopHandlerBuilder_)
2018 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2022 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2028 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2030 nbpu_opt_ = nbpu_opt;
2033 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2036 pme_fft_opt_ = pme_fft_opt;
2039 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2041 bonded_opt_ = bonded_opt;
2044 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2046 update_opt_ = update_opt;
2049 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2051 hardwareOptions_ = hardwareOptions;
2054 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2056 filenames_ = filenames;
2059 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2061 outputEnvironment_ = outputEnvironment;
2064 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2066 logFileHandle_ = logFileHandle;
2069 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2071 stopHandlerBuilder_ = std::move(builder);
2074 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2075 compat::not_null<SimulationContext*> context) :
2076 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2080 MdrunnerBuilder::~MdrunnerBuilder() = default;
2082 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2083 real forceWarningThreshold,
2084 const StartingBehavior startingBehavior)
2086 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2090 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2092 impl_->addDomdec(options);
2096 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2098 impl_->addVerletList(nstlist);
2102 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2104 impl_->addReplicaExchange(params);
2108 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2110 impl_->addNonBonded(nbpu_opt);
2114 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2116 // The builder method may become more general in the future, but in this version,
2117 // parameters for PME electrostatics are both required and the only parameters
2119 if (pme_opt && pme_fft_opt)
2121 impl_->addPME(pme_opt, pme_fft_opt);
2126 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2131 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2133 impl_->addBondedTaskAssignment(bonded_opt);
2137 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2139 impl_->addUpdateTaskAssignment(update_opt);
2143 Mdrunner MdrunnerBuilder::build()
2145 return impl_->build();
2148 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2150 impl_->addHardwareOptions(hardwareOptions);
2154 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2156 impl_->addFilenames(filenames);
2160 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2162 impl_->addOutputEnvironment(outputEnvironment);
2166 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2168 impl_->addLogFile(logFileHandle);
2172 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2174 impl_->addStopHandlerBuilder(std::move(builder));
2178 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2180 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;