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_stream_manager.h"
77 #include "gromacs/hardware/cpuinfo.h"
78 #include "gromacs/hardware/detecthardware.h"
79 #include "gromacs/hardware/device_management.h"
80 #include "gromacs/hardware/printhardware.h"
81 #include "gromacs/imd/imd.h"
82 #include "gromacs/listed_forces/disre.h"
83 #include "gromacs/listed_forces/gpubonded.h"
84 #include "gromacs/listed_forces/listed_forces.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/gpuforcereduction.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/mdrun/simulationinput.h"
109 #include "gromacs/mdrun/simulationinputhandle.h"
110 #include "gromacs/mdrunutility/handlerestart.h"
111 #include "gromacs/mdrunutility/logging.h"
112 #include "gromacs/mdrunutility/multisim.h"
113 #include "gromacs/mdrunutility/printtime.h"
114 #include "gromacs/mdrunutility/threadaffinity.h"
115 #include "gromacs/mdtypes/checkpointdata.h"
116 #include "gromacs/mdtypes/commrec.h"
117 #include "gromacs/mdtypes/enerdata.h"
118 #include "gromacs/mdtypes/fcdata.h"
119 #include "gromacs/mdtypes/forcerec.h"
120 #include "gromacs/mdtypes/group.h"
121 #include "gromacs/mdtypes/inputrec.h"
122 #include "gromacs/mdtypes/interaction_const.h"
123 #include "gromacs/mdtypes/md_enums.h"
124 #include "gromacs/mdtypes/mdatom.h"
125 #include "gromacs/mdtypes/mdrunoptions.h"
126 #include "gromacs/mdtypes/observableshistory.h"
127 #include "gromacs/mdtypes/simulation_workload.h"
128 #include "gromacs/mdtypes/state.h"
129 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
130 #include "gromacs/modularsimulator/modularsimulator.h"
131 #include "gromacs/nbnxm/gpu_data_mgmt.h"
132 #include "gromacs/nbnxm/nbnxm.h"
133 #include "gromacs/nbnxm/pairlist_tuning.h"
134 #include "gromacs/pbcutil/pbc.h"
135 #include "gromacs/pulling/output.h"
136 #include "gromacs/pulling/pull.h"
137 #include "gromacs/pulling/pull_rotation.h"
138 #include "gromacs/restraint/manager.h"
139 #include "gromacs/restraint/restraintmdmodule.h"
140 #include "gromacs/restraint/restraintpotential.h"
141 #include "gromacs/swap/swapcoords.h"
142 #include "gromacs/taskassignment/decidegpuusage.h"
143 #include "gromacs/taskassignment/decidesimulationworkload.h"
144 #include "gromacs/taskassignment/resourcedivision.h"
145 #include "gromacs/taskassignment/taskassignment.h"
146 #include "gromacs/taskassignment/usergpuids.h"
147 #include "gromacs/timing/gpu_timing.h"
148 #include "gromacs/timing/wallcycle.h"
149 #include "gromacs/timing/wallcyclereporting.h"
150 #include "gromacs/topology/mtop_util.h"
151 #include "gromacs/trajectory/trajectoryframe.h"
152 #include "gromacs/utility/basenetwork.h"
153 #include "gromacs/utility/cstringutil.h"
154 #include "gromacs/utility/exceptions.h"
155 #include "gromacs/utility/fatalerror.h"
156 #include "gromacs/utility/filestream.h"
157 #include "gromacs/utility/gmxassert.h"
158 #include "gromacs/utility/gmxmpi.h"
159 #include "gromacs/utility/keyvaluetree.h"
160 #include "gromacs/utility/logger.h"
161 #include "gromacs/utility/loggerbuilder.h"
162 #include "gromacs/utility/mdmodulenotification.h"
163 #include "gromacs/utility/physicalnodecommunicator.h"
164 #include "gromacs/utility/pleasecite.h"
165 #include "gromacs/utility/programcontext.h"
166 #include "gromacs/utility/smalloc.h"
167 #include "gromacs/utility/stringutil.h"
169 #include "isimulator.h"
170 #include "membedholder.h"
171 #include "replicaexchange.h"
172 #include "simulatorbuilder.h"
178 /*! \brief Manage any development feature flag variables encountered
180 * The use of dev features indicated by environment variables is
181 * logged in order to ensure that runs with such features enabled can
182 * be identified from their log and standard output. Any cross
183 * dependencies are also checked, and if unsatisfied, a fatal error
186 * Note that some development features overrides are applied already here:
187 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
189 * \param[in] mdlog Logger object.
190 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
191 * \param[in] pmeRunMode The PME run mode for this run
192 * \returns The object populated with development feature flags.
194 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
195 const bool useGpuForNonbonded,
196 const PmeRunMode pmeRunMode)
198 DevelopmentFeatureFlags devFlags;
200 // Some builds of GCC 5 give false positive warnings that these
201 // getenv results are ignored when clearly they are used.
202 #pragma GCC diagnostic push
203 #pragma GCC diagnostic ignored "-Wunused-result"
205 devFlags.enableGpuBufferOps =
206 GMX_GPU_CUDA && useGpuForNonbonded && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
207 devFlags.enableGpuHaloExchange = GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_DD_COMMS") != nullptr;
208 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr) || GMX_FAHCORE;
209 devFlags.enableGpuPmePPComm =
210 GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
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 if (!devFlags.enableGpuBufferOps)
271 GMX_LOG(mdlog.warning)
273 .appendTextFormatted(
274 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
275 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
276 devFlags.enableGpuBufferOps = true;
278 GMX_LOG(mdlog.warning)
280 .appendTextFormatted(
281 "This run uses the 'GPU PME-PP communications' feature, enabled "
282 "by the GMX_GPU_PME_PP_COMMS environment variable.");
286 std::string clarification;
287 if (pmeRunMode == PmeRunMode::Mixed)
290 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
295 clarification = "PME is not offloaded to the GPU.";
297 GMX_LOG(mdlog.warning)
300 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
301 "'GPU PME-PP communications' feature was not enabled as "
303 devFlags.enableGpuPmePPComm = false;
310 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
312 * Used to ensure that the master thread does not modify mdrunner during copy
313 * on the spawned threads. */
314 static void threadMpiMdrunnerAccessBarrier()
317 MPI_Barrier(MPI_COMM_WORLD);
321 Mdrunner Mdrunner::cloneOnSpawnedThread() const
323 auto newRunner = Mdrunner(std::make_unique<MDModules>());
325 // All runners in the same process share a restraint manager resource because it is
326 // part of the interface to the client code, which is associated only with the
327 // original thread. Handles to the same resources can be obtained by copy.
329 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
332 // Copy members of master runner.
333 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
334 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
335 newRunner.hw_opt = hw_opt;
336 newRunner.filenames = filenames;
338 newRunner.oenv = oenv;
339 newRunner.mdrunOptions = mdrunOptions;
340 newRunner.domdecOptions = domdecOptions;
341 newRunner.nbpu_opt = nbpu_opt;
342 newRunner.pme_opt = pme_opt;
343 newRunner.pme_fft_opt = pme_fft_opt;
344 newRunner.bonded_opt = bonded_opt;
345 newRunner.update_opt = update_opt;
346 newRunner.nstlist_cmdline = nstlist_cmdline;
347 newRunner.replExParams = replExParams;
348 newRunner.pforce = pforce;
349 // Give the spawned thread the newly created valid communicator
350 // for the simulation.
351 newRunner.libraryWorldCommunicator = MPI_COMM_WORLD;
352 newRunner.simulationCommunicator = MPI_COMM_WORLD;
354 newRunner.startingBehavior = startingBehavior;
355 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
356 newRunner.inputHolder_ = inputHolder_;
358 threadMpiMdrunnerAccessBarrier();
363 /*! \brief The callback used for running on spawned threads.
365 * Obtains the pointer to the master mdrunner object from the one
366 * argument permitted to the thread-launch API call, copies it to make
367 * a new runner for this thread, reinitializes necessary data, and
368 * proceeds to the simulation. */
369 static void mdrunner_start_fn(const void* arg)
373 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
374 /* copy the arg list to make sure that it's thread-local. This
375 doesn't copy pointed-to items, of course; fnm, cr and fplog
376 are reset in the call below, all others should be const. */
377 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
380 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
384 void Mdrunner::spawnThreads(int numThreadsToLaunch)
387 /* now spawn new threads that start mdrunner_start_fn(), while
388 the main thread returns. Thread affinity is handled later. */
389 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
390 static_cast<const void*>(this))
393 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
396 // Give the master thread the newly created valid communicator for
398 libraryWorldCommunicator = MPI_COMM_WORLD;
399 simulationCommunicator = MPI_COMM_WORLD;
400 threadMpiMdrunnerAccessBarrier();
402 GMX_UNUSED_VALUE(numThreadsToLaunch);
403 GMX_UNUSED_VALUE(mdrunner_start_fn);
409 /*! \brief Initialize variables for Verlet scheme simulation */
410 static void prepare_verlet_scheme(FILE* fplog,
414 const gmx_mtop_t* mtop,
416 bool makeGpuPairList,
417 const gmx::CpuInfo& cpuinfo)
419 // We checked the cut-offs in grompp, but double-check here.
420 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
421 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
423 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
424 "With Verlet lists and PME we should have rcoulomb>=rvdw");
428 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
429 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
431 /* For NVE simulations, we will retain the initial list buffer */
432 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
434 /* Update the Verlet buffer size for the current run setup */
436 /* Here we assume SIMD-enabled kernels are being used. But as currently
437 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
438 * and 4x2 gives a larger buffer than 4x4, this is ok.
440 ListSetupType listType =
441 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
442 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
444 const real rlist_new =
445 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
447 if (rlist_new != ir->rlist)
449 if (fplog != nullptr)
452 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
453 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
455 ir->rlist = rlist_new;
459 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
461 gmx_fatal(FARGS, "Can not set nstlist without %s",
462 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
465 if (EI_DYNAMICS(ir->eI))
467 /* Set or try nstlist values */
468 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
472 /*! \brief Override the nslist value in inputrec
474 * with value passed on the command line (if any)
476 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
480 /* override with anything else than the default -2 */
481 if (nsteps_cmdline > -2)
483 char sbuf_steps[STEPSTRSIZE];
484 char sbuf_msg[STRLEN];
486 ir->nsteps = nsteps_cmdline;
487 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
490 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
491 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
495 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
496 gmx_step_str(nsteps_cmdline, sbuf_steps));
499 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
501 else if (nsteps_cmdline < -2)
503 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
505 /* Do nothing if nsteps_cmdline == -2 */
511 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
513 * If not, and if a warning may be issued, logs a warning about
514 * falling back to CPU code. With thread-MPI, only the first
515 * call to this function should have \c issueWarning true. */
516 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
518 bool gpuIsUseful = true;
521 if (ir.opts.ngener - ir.nwall > 1)
523 /* The GPU code does not support more than one energy group.
524 * If the user requested GPUs explicitly, a fatal error is given later.
528 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
529 "For better performance, run on the GPU without energy groups and then do "
530 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
536 warning = "TPI is not implemented for GPUs.";
539 if (!gpuIsUseful && issueWarning)
541 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
547 //! Initializes the logger for mdrun.
548 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
550 gmx::LoggerBuilder builder;
551 if (fplog != nullptr)
553 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
555 if (isSimulationMasterRank)
557 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
559 return builder.build();
562 //! Make a TaskTarget from an mdrun argument string.
563 static TaskTarget findTaskTarget(const char* optionString)
565 TaskTarget returnValue = TaskTarget::Auto;
567 if (strncmp(optionString, "auto", 3) == 0)
569 returnValue = TaskTarget::Auto;
571 else if (strncmp(optionString, "cpu", 3) == 0)
573 returnValue = TaskTarget::Cpu;
575 else if (strncmp(optionString, "gpu", 3) == 0)
577 returnValue = TaskTarget::Gpu;
581 GMX_ASSERT(false, "Option string should have been checked for sanity already");
587 //! Finish run, aggregate data to print performance info.
588 static void finish_run(FILE* fplog,
589 const gmx::MDLogger& mdlog,
591 const t_inputrec* inputrec,
593 gmx_wallcycle_t wcycle,
594 gmx_walltime_accounting_t walltime_accounting,
595 nonbonded_verlet_t* nbv,
596 const gmx_pme_t* pme,
600 double nbfs = 0, mflop = 0;
601 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
602 elapsed_time_over_all_threads_over_all_ranks;
603 /* Control whether it is valid to print a report. Only the
604 simulation master may print, but it should not do so if the run
605 terminated e.g. before a scheduled reset step. This is
606 complicated by the fact that PME ranks are unaware of the
607 reason why they were sent a pmerecvqxFINISH. To avoid
608 communication deadlocks, we always do the communication for the
609 report, even if we've decided not to write the report, because
610 how long it takes to finish the run is not important when we've
611 decided not to report on the simulation performance.
613 Further, we only report performance for dynamical integrators,
614 because those are the only ones for which we plan to
615 consider doing any optimizations. */
616 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
618 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
620 GMX_LOG(mdlog.warning)
622 .appendText("Simulation ended prematurely, no performance report will be written.");
627 std::unique_ptr<t_nrnb> nrnbTotalStorage;
630 nrnbTotalStorage = std::make_unique<t_nrnb>();
631 nrnb_tot = nrnbTotalStorage.get();
633 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
641 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
642 elapsed_time_over_all_threads =
643 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
647 /* reduce elapsed_time over all MPI ranks in the current simulation */
648 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
650 elapsed_time_over_all_ranks /= cr->nnodes;
651 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
652 * current simulation. */
653 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
654 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
659 elapsed_time_over_all_ranks = elapsed_time;
660 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
665 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
668 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
670 print_dd_statistics(cr, inputrec, fplog);
673 /* TODO Move the responsibility for any scaling by thread counts
674 * to the code that handled the thread region, so that there's a
675 * mechanism to keep cycle counting working during the transition
676 * to task parallelism. */
677 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
678 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
679 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
680 nthreads_pp, nthreads_pme);
681 auto cycle_sum(wallcycle_sum(cr, wcycle));
685 auto nbnxn_gpu_timings =
686 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
687 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
689 if (pme_gpu_task_enabled(pme))
691 pme_gpu_get_timings(pme, &pme_gpu_timings);
693 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
694 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
697 if (EI_DYNAMICS(inputrec->eI))
699 delta_t = inputrec->delta_t;
704 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
705 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
706 delta_t, nbfs, mflop);
710 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
711 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
712 delta_t, nbfs, mflop);
717 int Mdrunner::mdrunner()
720 t_forcerec* fr = nullptr;
721 real ewaldcoeff_q = 0;
722 real ewaldcoeff_lj = 0;
723 int nChargePerturbed = -1, nTypePerturbed = 0;
724 gmx_wallcycle_t wcycle;
725 gmx_walltime_accounting_t walltime_accounting = nullptr;
726 MembedHolder membedHolder(filenames.size(), filenames.data());
727 gmx_hw_info_t* hwinfo = nullptr;
729 /* CAUTION: threads may be started later on in this function, so
730 cr doesn't reflect the final parallel state right now */
733 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
734 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
735 const bool doRerun = mdrunOptions.rerun;
737 // Handle task-assignment related user options.
738 EmulateGpuNonbonded emulateGpuNonbonded =
739 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
741 std::vector<int> userGpuTaskAssignment;
744 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
746 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
747 auto nonbondedTarget = findTaskTarget(nbpu_opt);
748 auto pmeTarget = findTaskTarget(pme_opt);
749 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
750 auto bondedTarget = findTaskTarget(bonded_opt);
751 auto updateTarget = findTaskTarget(update_opt);
753 FILE* fplog = nullptr;
754 // If we are appending, we don't write log output because we need
755 // to check that the old log file matches what the checkpoint file
756 // expects. Otherwise, we should start to write log output now if
757 // there is a file ready for it.
758 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
760 fplog = gmx_fio_getfp(logFileHandle);
762 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, simulationCommunicator);
763 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
764 gmx::MDLogger mdlog(logOwner.logger());
766 // TODO The thread-MPI master rank makes a working
767 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
768 // after the threads have been launched. This works because no use
769 // is made of that communicator until after the execution paths
770 // have rejoined. But it is likely that we can improve the way
771 // this is expressed, e.g. by expressly running detection only the
772 // master rank for thread-MPI, rather than relying on the mutex
773 // and reference count.
774 PhysicalNodeCommunicator physicalNodeComm(libraryWorldCommunicator, gmx_physicalnode_id_hash());
775 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
777 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
779 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->deviceInfoList, hw_opt.gpuIdsAvailable);
780 const int numDevicesToUse = gmx::ssize(gpuIdsToUse);
782 // Print citation requests after all software/hardware printing
783 pleaseCiteGromacs(fplog);
785 // Note: legacy program logic relies on checking whether these pointers are assigned.
786 // Objects may or may not be allocated later.
787 std::unique_ptr<t_inputrec> inputrec;
788 std::unique_ptr<t_state> globalState;
790 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
792 if (isSimulationMasterRank)
794 // Allocate objects to be initialized by later function calls.
795 /* Only the master rank has the global state */
796 globalState = std::make_unique<t_state>();
797 inputrec = std::make_unique<t_inputrec>();
799 /* Read (nearly) all data required for the simulation
800 * and keep the partly serialized tpr contents to send to other ranks later
802 applyGlobalSimulationState(*inputHolder_.get(), partialDeserializedTpr.get(),
803 globalState.get(), inputrec.get(), &mtop);
806 /* Check and update the hardware options for internal consistency */
807 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
810 if (GMX_THREAD_MPI && isSimulationMasterRank)
812 bool useGpuForNonbonded = false;
813 bool useGpuForPme = false;
816 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
818 // If the user specified the number of ranks, then we must
819 // respect that, but in default mode, we need to allow for
820 // the number of GPUs to choose the number of ranks.
821 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
822 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
823 nonbondedTarget, numDevicesToUse, userGpuTaskAssignment, emulateGpuNonbonded,
824 canUseGpuForNonbonded,
825 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
826 hw_opt.nthreads_tmpi);
827 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
828 useGpuForNonbonded, pmeTarget, numDevicesToUse, userGpuTaskAssignment, *hwinfo,
829 *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
831 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
833 /* Determine how many thread-MPI ranks to start.
835 * TODO Over-writing the user-supplied value here does
836 * prevent any possible subsequent checks from working
838 hw_opt.nthreads_tmpi =
839 get_nthreads_mpi(hwinfo, &hw_opt, numDevicesToUse, useGpuForNonbonded, useGpuForPme,
840 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
842 // Now start the threads for thread MPI.
843 spawnThreads(hw_opt.nthreads_tmpi);
844 // The spawned threads enter mdrunner() and execution of
845 // master and spawned threads joins at the end of this block.
847 PhysicalNodeCommunicator(libraryWorldCommunicator, gmx_physicalnode_id_hash());
850 GMX_RELEASE_ASSERT(ms || simulationCommunicator != MPI_COMM_NULL,
851 "Must have valid communicator unless running a multi-simulation");
852 CommrecHandle crHandle = init_commrec(simulationCommunicator);
853 t_commrec* cr = crHandle.get();
854 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
858 /* now broadcast everything to the non-master nodes/threads: */
859 if (!isSimulationMasterRank)
861 // Until now, only the master rank has a non-null pointer.
862 // On non-master ranks, allocate the object that will receive data in the following call.
863 inputrec = std::make_unique<t_inputrec>();
865 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec.get(), &mtop,
866 partialDeserializedTpr.get());
868 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
869 partialDeserializedTpr.reset(nullptr);
871 // Now the number of ranks is known to all ranks, and each knows
872 // the inputrec read by the master rank. The ranks can now all run
873 // the task-deciding functions and will agree on the result
874 // without needing to communicate.
875 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
877 // Note that these variables describe only their own node.
879 // Note that when bonded interactions run on a GPU they always run
880 // alongside a nonbonded task, so do not influence task assignment
881 // even though they affect the force calculation workload.
882 bool useGpuForNonbonded = false;
883 bool useGpuForPme = false;
884 bool useGpuForBonded = false;
885 bool useGpuForUpdate = false;
886 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
889 // It's possible that there are different numbers of GPUs on
890 // different nodes, which is the user's responsibility to
891 // handle. If unsuitable, we will notice that during task
893 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
894 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
895 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
896 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
897 useGpuForPme = decideWhetherToUseGpusForPme(
898 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec,
899 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
900 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
901 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
902 useGpuForBonded = decideWhetherToUseGpusForBonded(
903 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
904 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
905 domdecOptions.numPmeRanks, gpusWereDetected);
907 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
909 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
911 // Initialize development feature flags that enabled by environment variable
912 // and report those features that are enabled.
913 const DevelopmentFeatureFlags devFlags =
914 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
916 const bool useModularSimulator =
917 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
918 nullptr, doEssentialDynamics, membedHolder.doMembed());
921 // TODO: hide restraint implementation details from Mdrunner.
922 // There is nothing unique about restraints at this point as far as the
923 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
924 // factory functions from the SimulationContext on which to call mdModules_->add().
925 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
926 for (auto&& restraint : restraintManager_->getRestraints())
928 auto module = RestraintMDModule::create(restraint, restraint->sites());
929 mdModules_->add(std::move(module));
932 // TODO: Error handling
933 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
934 // now that the MdModules know their options, they know which callbacks to sign up to
935 mdModules_->subscribeToSimulationSetupNotifications();
936 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
938 if (inputrec->internalParameters != nullptr)
940 mdModulesNotifier.notify(*inputrec->internalParameters);
943 if (fplog != nullptr)
945 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
946 fprintf(fplog, "\n");
951 /* In rerun, set velocities to zero if present */
952 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
954 // rerun does not use velocities
958 "Rerun trajectory contains velocities. Rerun does only evaluate "
959 "potential energy and forces. The velocities will be ignored.");
960 for (int i = 0; i < globalState->natoms; i++)
962 clear_rvec(globalState->v[i]);
964 globalState->flags &= ~(1 << estV);
967 /* now make sure the state is initialized and propagated */
968 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
971 /* NM and TPI parallelize over force/energy calculations, not atoms,
972 * so we need to initialize and broadcast the global state.
974 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
978 globalState = std::make_unique<t_state>();
980 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
984 /* A parallel command line option consistency check that we can
985 only do after any threads have started. */
987 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
988 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
991 "The -dd or -npme option request a parallel simulation, "
993 "but %s was compiled without threads or MPI enabled",
994 output_env_get_program_display_name(oenv));
996 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
998 "but %s was not started through mpirun/mpiexec or only one rank was requested "
999 "through mpirun/mpiexec",
1000 output_env_get_program_display_name(oenv));
1004 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1007 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1008 "these are not compatible with mdrun -rerun");
1011 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1013 if (domdecOptions.numPmeRanks > 0)
1015 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1016 "PME-only ranks are requested, but the system does not use PME "
1017 "for electrostatics or LJ");
1020 domdecOptions.numPmeRanks = 0;
1023 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1025 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1026 * improve performance with many threads per GPU, since our OpenMP
1027 * scaling is bad, but it's difficult to automate the setup.
1029 domdecOptions.numPmeRanks = 0;
1033 if (domdecOptions.numPmeRanks < 0)
1035 domdecOptions.numPmeRanks = 0;
1036 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1040 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1041 "PME GPU decomposition is not supported");
1045 /* NMR restraints must be initialized before load_checkpoint,
1046 * since with time averaging the history is added to t_state.
1047 * For proper consistency check we therefore need to extend
1049 * So the PME-only nodes (if present) will also initialize
1050 * the distance restraints.
1053 /* This needs to be called before read_checkpoint to extend the state */
1054 t_disresdata* disresdata;
1055 snew(disresdata, 1);
1056 init_disres(fplog, &mtop, inputrec.get(), DisResRunMode::MDRun,
1057 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1058 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1059 globalState.get(), replExParams.exchangeInterval > 0);
1061 t_oriresdata* oriresdata;
1062 snew(oriresdata, 1);
1063 init_orires(fplog, &mtop, inputrec.get(), cr, ms, globalState.get(), oriresdata);
1065 auto deform = prepareBoxDeformation(
1066 globalState != nullptr ? globalState->box : box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1067 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mygroup, *inputrec);
1070 /* We have to remember the generation's first step before reading checkpoint.
1071 This way, we can report to the F@H core both the generation's first step
1072 and the restored first step, thus making it able to distinguish between
1073 an interruption/resume and start of the n-th generation simulation.
1074 Having this information, the F@H core can correctly calculate and report
1077 int gen_first_step = 0;
1080 gen_first_step = inputrec->init_step;
1084 ObservablesHistory observablesHistory = {};
1086 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1087 if (startingBehavior != StartingBehavior::NewSimulation)
1089 /* Check if checkpoint file exists before doing continuation.
1090 * This way we can use identical input options for the first and subsequent runs...
1092 if (mdrunOptions.numStepsCommandline > -2)
1094 /* Temporarily set the number of steps to unlimited to avoid
1095 * triggering the nsteps check in load_checkpoint().
1096 * This hack will go away soon when the -nsteps option is removed.
1098 inputrec->nsteps = -1;
1101 // Finish applying initial simulation state information from external sources on all ranks.
1102 // Reconcile checkpoint file data with Mdrunner state established up to this point.
1103 applyLocalState(*inputHolder_.get(), logFileHandle, cr, domdecOptions.numCells,
1104 inputrec.get(), globalState.get(), &observablesHistory,
1105 mdrunOptions.reproducible, mdModules_->notifier(),
1106 modularSimulatorCheckpointData.get(), useModularSimulator);
1107 // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1109 // on all code paths.
1110 // Write checkpoint or provide hook to update SimulationInput.
1111 // If there was a checkpoint file, SimulationInput contains more information
1112 // than if there wasn't. At this point, we have synchronized the in-memory
1113 // state with the filesystem state only for restarted simulations. We should
1114 // be calling applyLocalState unconditionally and expect that the completeness
1115 // of SimulationInput is not dependent on its creation method.
1117 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1119 // Now we can start normal logging to the truncated log file.
1120 fplog = gmx_fio_getfp(logFileHandle);
1121 prepareLogAppending(fplog);
1122 logOwner = buildLogger(fplog, MASTER(cr));
1123 mdlog = logOwner.logger();
1130 fcRegisterSteps(inputrec->nsteps + inputrec->init_step, gen_first_step);
1134 if (mdrunOptions.numStepsCommandline > -2)
1139 "The -nsteps functionality is deprecated, and may be removed in a future "
1141 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1144 /* override nsteps with value set on the commandline */
1145 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1147 if (isSimulationMasterRank)
1149 copy_mat(globalState->box, box);
1154 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1157 if (inputrec->cutoff_scheme != ecutsVERLET)
1160 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1161 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1163 /* Update rlist and nstlist. */
1164 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1165 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1166 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1168 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1169 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1172 // This builder is necessary while we have multi-part construction
1173 // of DD. Before DD is constructed, we use the existence of
1174 // the builder object to indicate that further construction of DD
1176 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1177 if (useDomainDecomposition)
1179 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1180 mdlog, cr, domdecOptions, mdrunOptions, mtop, *inputrec, box,
1181 positionsFromStatePointer(globalState.get()));
1185 /* PME, if used, is done on all nodes with 1D decomposition */
1186 cr->nnodes = cr->sizeOfDefaultCommunicator;
1187 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1188 cr->nodeid = cr->rankInDefaultCommunicator;
1190 cr->duty = (DUTY_PP | DUTY_PME);
1192 if (inputrec->pbcType == PbcType::Screw)
1194 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1198 // Produce the task assignment for this rank - done after DD is constructed
1199 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1200 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, simulationCommunicator, physicalNodeComm,
1201 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1202 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1203 // TODO cr->duty & DUTY_PME should imply that a PME
1204 // algorithm is active, but currently does not.
1205 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1207 // Get the device handles for the modules, nullptr when no task is assigned.
1209 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1211 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1212 bool useTiming = true;
1216 /* WARNING: CUDA timings are incorrect with multiple streams.
1217 * This is the main reason why they are disabled by default.
1219 // TODO: Consider turning on by default when we can detect nr of streams.
1220 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1222 else if (GMX_GPU_OPENCL)
1224 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1227 // TODO Currently this is always built, yet DD partition code
1228 // checks if it is built before using it. Probably it should
1229 // become an MDModule that is made only when another module
1230 // requires it (e.g. pull, CompEl, density fitting), so that we
1231 // don't update the local atom sets unilaterally every step.
1232 LocalAtomSetManager atomSets;
1235 // TODO Pass the GPU streams to ddBuilder to use in buffer
1236 // transfers (e.g. halo exchange)
1237 cr->dd = ddBuilder->build(&atomSets);
1238 // The builder's job is done, so destruct it
1239 ddBuilder.reset(nullptr);
1240 // Note that local state still does not exist yet.
1243 // The GPU update is decided here because we need to know whether the constraints or
1244 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1245 // defined). This is only known after DD is initialized, hence decision on using GPU
1246 // update is done so late.
1249 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1251 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1252 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1253 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1254 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1255 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1257 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1259 const bool printHostName = (cr->nnodes > 1);
1260 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1262 const bool disableNonbondedCalculation = (getenv("GMX_NO_NONBONDED") != nullptr);
1263 if (disableNonbondedCalculation)
1265 /* turn off non-bonded calculations */
1266 GMX_LOG(mdlog.warning)
1269 "Found environment variable GMX_NO_NONBONDED.\n"
1270 "Disabling nonbonded calculations.");
1273 MdrunScheduleWorkload runScheduleWork;
1275 bool useGpuDirectHalo = decideWhetherToUseGpuForHalo(
1276 devFlags, havePPDomainDecomposition(cr), useGpuForNonbonded, useModularSimulator,
1277 doRerun, EI_ENERGY_MINIMIZATION(inputrec->eI));
1279 // Also populates the simulation constant workload description.
1280 runScheduleWork.simulationWork = createSimulationWorkload(
1281 *inputrec, disableNonbondedCalculation, devFlags, useGpuForNonbonded, pmeRunMode,
1282 useGpuForBonded, useGpuForUpdate, useGpuDirectHalo);
1284 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1286 if (deviceInfo != nullptr)
1288 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1290 dd_setup_dlb_resource_sharing(cr, deviceId);
1292 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1293 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, useTiming);
1296 // If the user chose a task assignment, give them some hints
1297 // where appropriate.
1298 if (!userGpuTaskAssignment.empty())
1300 gpuTaskAssignments.logPerformanceHints(mdlog, numDevicesToUse);
1305 /* After possible communicator splitting in make_dd_communicators.
1306 * we can set up the intra/inter node communication.
1308 gmx_setup_nodecomm(fplog, cr);
1314 GMX_LOG(mdlog.warning)
1316 .appendTextFormatted(
1317 "This is simulation %d out of %d running as a composite GROMACS\n"
1318 "multi-simulation job. Setup for this simulation:\n",
1319 ms->simulationIndex_, ms->numSimulations_);
1321 GMX_LOG(mdlog.warning)
1322 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1324 cr->nnodes == 1 ? "thread" : "threads"
1326 cr->nnodes == 1 ? "process" : "processes"
1332 // If mdrun -pin auto honors any affinity setting that already
1333 // exists. If so, it is nice to provide feedback about whether
1334 // that existing affinity setting was from OpenMP or something
1335 // else, so we run this code both before and after we initialize
1336 // the OpenMP support.
1337 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1338 /* Check and update the number of OpenMP threads requested */
1339 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1340 pmeRunMode, mtop, *inputrec);
1342 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1343 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1345 // Enable FP exception detection, but not in
1346 // Release mode and not for compilers with known buggy FP
1347 // exception support (clang with any optimization) or suspected
1348 // buggy FP exception support (gcc 7.* with optimization).
1349 #if !defined NDEBUG \
1350 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1351 && defined __OPTIMIZE__)
1352 const bool bEnableFPE = true;
1354 const bool bEnableFPE = false;
1356 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1359 gmx_feenableexcept();
1362 /* Now that we know the setup is consistent, check for efficiency */
1363 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1364 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1366 /* getting number of PP/PME threads on this MPI / tMPI rank.
1367 PME: env variable should be read only on one node to make sure it is
1368 identical everywhere;
1370 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1371 : gmx_omp_nthreads_get(emntPME);
1372 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1373 physicalNodeComm, mdlog);
1375 // Enable Peer access between GPUs where available
1376 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1377 // any of the GPU communication features are active.
1378 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1379 && (runScheduleWork.simulationWork.useGpuHaloExchange
1380 || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1382 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1385 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1387 /* Before setting affinity, check whether the affinity has changed
1388 * - which indicates that probably the OpenMP library has changed it
1389 * since we first checked).
1391 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1393 int numThreadsOnThisNode, intraNodeThreadOffset;
1394 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1395 &intraNodeThreadOffset);
1397 /* Set the CPU affinity */
1398 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1399 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1402 if (mdrunOptions.timingOptions.resetStep > -1)
1407 "The -resetstep functionality is deprecated, and may be removed in a "
1410 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1414 /* Master synchronizes its value of reset_counters with all nodes
1415 * including PME only nodes */
1416 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1417 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1418 wcycle_set_reset_counters(wcycle, reset_counters);
1421 // Membrane embedding must be initialized before we call init_forcerec()
1422 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1423 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1425 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1426 std::unique_ptr<MDAtoms> mdAtoms;
1427 std::unique_ptr<VirtualSitesHandler> vsite;
1428 std::unique_ptr<GpuBonded> gpuBonded;
1431 if (thisRankHasDuty(cr, DUTY_PP))
1433 mdModulesNotifier.notify(*cr);
1434 mdModulesNotifier.notify(&atomSets);
1435 mdModulesNotifier.notify(inputrec->pbcType);
1436 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1437 /* Initiate forcerecord */
1438 fr = new t_forcerec;
1439 fr->forceProviders = mdModules_->initForceProviders();
1440 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1441 opt2fn("-table", filenames.size(), filenames.data()),
1442 opt2fn("-tablep", filenames.size(), filenames.data()),
1443 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1444 // Dirty hack, for fixing disres and orires should be made mdmodules
1445 fr->fcdata->disres = disresdata;
1446 fr->fcdata->orires = oriresdata;
1448 // Save a handle to device stream manager to use elsewhere in the code
1449 // TODO: Forcerec is not a correct place to store it.
1450 fr->deviceStreamManager = deviceStreamManager.get();
1452 if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1455 deviceStreamManager != nullptr,
1456 "GPU device stream manager should be valid in order to use PME-PP direct "
1459 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1460 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1462 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1463 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1464 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1467 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo,
1468 runScheduleWork.simulationWork.useGpuNonbonded,
1469 deviceStreamManager.get(), &mtop, box, wcycle);
1470 // TODO: Move the logic below to a GPU bonded builder
1471 if (runScheduleWork.simulationWork.useGpuBonded)
1473 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1474 "GPU device stream manager should be valid in order to use GPU "
1475 "version of bonded forces.");
1476 gpuBonded = std::make_unique<GpuBonded>(
1477 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1478 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1479 fr->gpuBonded = gpuBonded.get();
1482 /* Initialize the mdAtoms structure.
1483 * mdAtoms is not filled with atom data,
1484 * as this can not be done now with domain decomposition.
1486 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1487 if (globalState && thisRankHasPmeGpuTask)
1489 // The pinning of coordinates in the global state object works, because we only use
1490 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1491 // points to the global state object without DD.
1492 // FIXME: MD and EM separately set up the local state - this should happen in the same
1493 // function, which should also perform the pinning.
1494 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1497 /* Initialize the virtual site communication */
1498 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1500 calc_shifts(box, fr->shift_vec);
1502 /* With periodic molecules the charge groups should be whole at start up
1503 * and the virtual sites should not be far from their proper positions.
1505 if (!inputrec->bContinuation && MASTER(cr)
1506 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1508 /* Make molecules whole at start of run */
1509 if (fr->pbcType != PbcType::No)
1511 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1515 /* Correct initial vsite positions are required
1516 * for the initial distribution in the domain decomposition
1517 * and for the initial shell prediction.
1519 constructVirtualSitesGlobal(mtop, globalState->x);
1523 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1525 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1526 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1531 /* This is a PME only node */
1533 GMX_ASSERT(globalState == nullptr,
1534 "We don't need the state on a PME only rank and expect it to be unitialized");
1536 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1537 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1540 gmx_pme_t* sepPmeData = nullptr;
1541 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1542 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1543 "Double-checking that only PME-only ranks have no forcerec");
1544 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1546 // TODO should live in ewald module once its testing is improved
1548 // Later, this program could contain kernels that might be later
1549 // re-used as auto-tuning progresses, or subsequent simulations
1551 PmeGpuProgramStorage pmeGpuProgram;
1552 if (thisRankHasPmeGpuTask)
1555 (deviceStreamManager != nullptr),
1556 "GPU device stream manager should be initialized in order to use GPU for PME.");
1557 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1558 "GPU device should be initialized in order to use GPU for PME.");
1559 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1562 /* Initiate PME if necessary,
1563 * either on all nodes or on dedicated PME nodes only. */
1564 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1566 if (mdAtoms && mdAtoms->mdatoms())
1568 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1569 if (EVDW_PME(inputrec->vdwtype))
1571 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1574 if (cr->npmenodes > 0)
1576 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1577 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1578 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1581 if (thisRankHasDuty(cr, DUTY_PME))
1585 // TODO: This should be in the builder.
1586 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1587 || (deviceStreamManager != nullptr),
1588 "Device stream manager should be valid in order to use GPU "
1591 !runScheduleWork.simulationWork.useGpuPme
1592 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1593 "GPU PME stream should be valid in order to use GPU version of PME.");
1595 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1596 ? &deviceStreamManager->context()
1598 const DeviceStream* pmeStream =
1599 runScheduleWork.simulationWork.useGpuPme
1600 ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1603 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1604 nChargePerturbed != 0, nTypePerturbed != 0,
1605 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1606 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1607 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1609 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1614 if (EI_DYNAMICS(inputrec->eI))
1616 /* Turn on signal handling on all nodes */
1618 * (A user signal from the PME nodes (if any)
1619 * is communicated to the PP nodes.
1621 signal_handler_install();
1624 pull_t* pull_work = nullptr;
1625 if (thisRankHasDuty(cr, DUTY_PP))
1627 /* Assumes uniform use of the number of OpenMP threads */
1628 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1630 if (inputrec->bPull)
1632 /* Initialize pull code */
1633 pull_work = init_pull(fplog, inputrec->pull.get(), inputrec.get(), &mtop, cr, &atomSets,
1634 inputrec->fepvals->init_lambda);
1635 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1637 initPullHistory(pull_work, &observablesHistory);
1639 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1641 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1645 std::unique_ptr<EnforcedRotation> enforcedRotation;
1648 /* Initialize enforced rotation code */
1649 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1650 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1654 t_swap* swap = nullptr;
1655 if (inputrec->eSwapCoords != eswapNO)
1657 /* Initialize ion swapping code */
1658 swap = init_swapcoords(fplog, inputrec.get(),
1659 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1660 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1661 oenv, mdrunOptions, startingBehavior);
1664 /* Let makeConstraints know whether we have essential dynamics constraints. */
1665 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1666 ms, &nrnb, wcycle, fr->bMolPBC);
1668 /* Energy terms and groups */
1669 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1670 inputrec->fepvals->n_lambda);
1672 /* Kinetic energy data */
1673 gmx_ekindata_t ekind;
1674 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1676 /* Set up interactive MD (IMD) */
1678 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1679 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1680 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1682 if (DOMAINDECOMP(cr))
1684 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1685 /* This call is not included in init_domain_decomposition mainly
1686 * because fr->cginfo_mb is set later.
1688 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1689 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1692 if (runScheduleWork.simulationWork.useGpuBufferOps)
1694 fr->gpuForceReduction[gmx::AtomLocality::Local] = std::make_unique<gmx::GpuForceReduction>(
1695 deviceStreamManager->context(),
1696 deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedLocal), wcycle);
1697 fr->gpuForceReduction[gmx::AtomLocality::NonLocal] = std::make_unique<gmx::GpuForceReduction>(
1698 deviceStreamManager->context(),
1699 deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedNonLocal), wcycle);
1702 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1703 if (gpusWereDetected
1704 && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
1705 || runScheduleWork.simulationWork.useGpuBufferOps))
1707 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1708 ? GpuApiCallBehavior::Async
1709 : GpuApiCallBehavior::Sync;
1710 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1711 "GPU device stream manager should be initialized to use GPU.");
1712 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1713 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1714 fr->stateGpu = stateGpu.get();
1717 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1718 SimulatorBuilder simulatorBuilder;
1720 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1721 simulatorBuilder.add(std::move(membedHolder));
1722 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1723 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1726 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1727 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1728 simulatorBuilder.add(ConstraintsParam(
1729 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1731 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1732 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1733 inputrec.get(), fr));
1734 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1735 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1736 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1737 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1738 // Todo move to an MDModule
1739 simulatorBuilder.add(IonSwapping(swap));
1740 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1741 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1742 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
1744 // build and run simulator object based on user-input
1745 auto simulator = simulatorBuilder.build(useModularSimulator);
1748 if (fr->pmePpCommGpu)
1750 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1751 fr->pmePpCommGpu.reset();
1754 if (inputrec->bPull)
1756 finish_pull(pull_work);
1758 finish_swapcoords(swap);
1762 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1764 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1765 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1766 deviceStreamManager.get());
1769 wallcycle_stop(wcycle, ewcRUN);
1771 /* Finish up, write some stuff
1772 * if rerunMD, don't write last frame again
1774 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1775 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1777 // clean up cycle counter
1778 wallcycle_destroy(wcycle);
1780 deviceStreamManager.reset(nullptr);
1784 gmx_pme_destroy(pmedata);
1788 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1789 // before we destroy the GPU context(s)
1790 // Pinned buffers are associated with contexts in CUDA.
1791 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1792 mdAtoms.reset(nullptr);
1793 globalState.reset(nullptr);
1794 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1795 gpuBonded.reset(nullptr);
1796 /* Free pinned buffers in *fr */
1799 // TODO convert to C++ so we can get rid of these frees
1803 if (!hwinfo->deviceInfoList.empty())
1805 /* stop the GPU profiler (only CUDA) */
1809 /* With tMPI we need to wait for all ranks to finish deallocation before
1810 * destroying the CUDA context as some tMPI ranks may be sharing
1813 * This is not a concern in OpenCL where we use one context per rank.
1815 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1816 * but it is easier and more futureproof to call it on the whole node.
1818 * Note that this function needs to be called even if GPUs are not used
1819 * in this run because the PME ranks have no knowledge of whether GPUs
1820 * are used or not, but all ranks need to enter the barrier below.
1821 * \todo Remove this physical node barrier after making sure
1822 * that it's not needed anymore (with a shared GPU run).
1826 physicalNodeComm.barrier();
1828 releaseDevice(deviceInfo);
1830 /* Does what it says */
1831 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1832 walltime_accounting_destroy(walltime_accounting);
1834 // Ensure log file content is written
1837 gmx_fio_flush(logFileHandle);
1840 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1841 * exceptions were enabled before function was called. */
1844 gmx_fedisableexcept();
1847 auto rc = static_cast<int>(gmx_get_stop_condition());
1850 /* we need to join all threads. The sub-threads join when they
1851 exit this function, but the master thread needs to be told to
1861 Mdrunner::~Mdrunner()
1863 // Clean up of the Manager.
1864 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1865 // but okay as long as threads synchronize some time before adding or accessing
1866 // a new set of restraints.
1867 if (restraintManager_)
1869 restraintManager_->clear();
1870 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1871 "restraints added during runner life time should be cleared at runner "
1876 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1878 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1879 // Not sure if this should be logged through the md logger or something else,
1880 // but it is helpful to have some sort of INFO level message sent somewhere.
1881 // std::cout << "Registering restraint named " << name << std::endl;
1883 // When multiple restraints are used, it may be wasteful to register them separately.
1884 // Maybe instead register an entire Restraint Manager as a force provider.
1885 restraintManager_->addToSpec(std::move(puller), name);
1888 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1890 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1892 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1894 class Mdrunner::BuilderImplementation
1897 BuilderImplementation() = delete;
1898 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1899 ~BuilderImplementation();
1901 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1902 real forceWarningThreshold,
1903 StartingBehavior startingBehavior);
1905 void addDomdec(const DomdecOptions& options);
1907 void addInput(SimulationInputHandle inputHolder);
1909 void addVerletList(int nstlist);
1911 void addReplicaExchange(const ReplicaExchangeParameters& params);
1913 void addNonBonded(const char* nbpu_opt);
1915 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1917 void addBondedTaskAssignment(const char* bonded_opt);
1919 void addUpdateTaskAssignment(const char* update_opt);
1921 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1923 void addFilenames(ArrayRef<const t_filenm> filenames);
1925 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1927 void addLogFile(t_fileio* logFileHandle);
1929 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1934 // Default parameters copied from runner.h
1935 // \todo Clarify source(s) of default parameters.
1937 const char* nbpu_opt_ = nullptr;
1938 const char* pme_opt_ = nullptr;
1939 const char* pme_fft_opt_ = nullptr;
1940 const char* bonded_opt_ = nullptr;
1941 const char* update_opt_ = nullptr;
1943 MdrunOptions mdrunOptions_;
1945 DomdecOptions domdecOptions_;
1947 ReplicaExchangeParameters replicaExchangeParameters_;
1949 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1952 //! World communicator, used for hardware detection and task assignment
1953 MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
1955 //! Multisim communicator handle.
1956 gmx_multisim_t* multiSimulation_;
1958 //! mdrun communicator
1959 MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
1961 //! Print a warning if any force is larger than this (in kJ/mol nm).
1962 real forceWarningThreshold_ = -1;
1964 //! Whether the simulation will start afresh, or restart with/without appending.
1965 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1967 //! The modules that comprise the functionality of mdrun.
1968 std::unique_ptr<MDModules> mdModules_;
1970 //! \brief Parallelism information.
1971 gmx_hw_opt_t hardwareOptions_;
1973 //! filename options for simulation.
1974 ArrayRef<const t_filenm> filenames_;
1976 /*! \brief Handle to output environment.
1978 * \todo gmx_output_env_t needs lifetime management.
1980 gmx_output_env_t* outputEnvironment_ = nullptr;
1982 /*! \brief Non-owning handle to MD log file.
1984 * \todo Context should own output facilities for client.
1985 * \todo Improve log file handle management.
1987 * Code managing the FILE* relies on the ability to set it to
1988 * nullptr to check whether the filehandle is valid.
1990 t_fileio* logFileHandle_ = nullptr;
1993 * \brief Builder for simulation stop signal handler.
1995 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1998 * \brief Sources for initial simulation state.
2000 * See issue #3652 for near-term refinements to the SimulationInput interface.
2002 * See issue #3379 for broader discussion on API aspects of simulation inputs and outputs.
2004 SimulationInputHandle inputHolder_;
2007 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
2008 compat::not_null<SimulationContext*> context) :
2009 mdModules_(std::move(mdModules))
2011 libraryWorldCommunicator_ = context->libraryWorldCommunicator_;
2012 simulationCommunicator_ = context->simulationCommunicator_;
2013 multiSimulation_ = context->multiSimulation_.get();
2016 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
2018 Mdrunner::BuilderImplementation&
2019 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
2020 const real forceWarningThreshold,
2021 const StartingBehavior startingBehavior)
2023 mdrunOptions_ = options;
2024 forceWarningThreshold_ = forceWarningThreshold;
2025 startingBehavior_ = startingBehavior;
2029 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
2031 domdecOptions_ = options;
2034 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
2039 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
2041 replicaExchangeParameters_ = params;
2044 Mdrunner Mdrunner::BuilderImplementation::build()
2046 auto newRunner = Mdrunner(std::move(mdModules_));
2048 newRunner.mdrunOptions = mdrunOptions_;
2049 newRunner.pforce = forceWarningThreshold_;
2050 newRunner.startingBehavior = startingBehavior_;
2051 newRunner.domdecOptions = domdecOptions_;
2053 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2054 newRunner.hw_opt = hardwareOptions_;
2056 // No invariant to check. This parameter exists to optionally override other behavior.
2057 newRunner.nstlist_cmdline = nstlist_;
2059 newRunner.replExParams = replicaExchangeParameters_;
2061 newRunner.filenames = filenames_;
2063 newRunner.libraryWorldCommunicator = libraryWorldCommunicator_;
2065 newRunner.simulationCommunicator = simulationCommunicator_;
2067 // nullptr is a valid value for the multisim handle
2068 newRunner.ms = multiSimulation_;
2072 newRunner.inputHolder_ = std::move(inputHolder_);
2076 GMX_THROW(gmx::APIError("MdrunnerBuilder::addInput() is required before build()."));
2079 // \todo Clarify ownership and lifetime management for gmx_output_env_t
2080 // \todo Update sanity checking when output environment has clearly specified invariants.
2081 // Initialization and default values for oenv are not well specified in the current version.
2082 if (outputEnvironment_)
2084 newRunner.oenv = outputEnvironment_;
2088 GMX_THROW(gmx::APIError(
2089 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2092 newRunner.logFileHandle = logFileHandle_;
2096 newRunner.nbpu_opt = nbpu_opt_;
2100 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2103 if (pme_opt_ && pme_fft_opt_)
2105 newRunner.pme_opt = pme_opt_;
2106 newRunner.pme_fft_opt = pme_fft_opt_;
2110 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2115 newRunner.bonded_opt = bonded_opt_;
2119 GMX_THROW(gmx::APIError(
2120 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2125 newRunner.update_opt = update_opt_;
2129 GMX_THROW(gmx::APIError(
2130 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2134 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2136 if (stopHandlerBuilder_)
2138 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2142 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2148 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2150 nbpu_opt_ = nbpu_opt;
2153 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2156 pme_fft_opt_ = pme_fft_opt;
2159 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2161 bonded_opt_ = bonded_opt;
2164 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2166 update_opt_ = update_opt;
2169 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2171 hardwareOptions_ = hardwareOptions;
2174 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2176 filenames_ = filenames;
2179 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2181 outputEnvironment_ = outputEnvironment;
2184 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2186 logFileHandle_ = logFileHandle;
2189 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2191 stopHandlerBuilder_ = std::move(builder);
2194 void Mdrunner::BuilderImplementation::addInput(SimulationInputHandle inputHolder)
2196 inputHolder_ = std::move(inputHolder);
2199 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2200 compat::not_null<SimulationContext*> context) :
2201 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2205 MdrunnerBuilder::~MdrunnerBuilder() = default;
2207 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2208 real forceWarningThreshold,
2209 const StartingBehavior startingBehavior)
2211 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2215 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2217 impl_->addDomdec(options);
2221 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2223 impl_->addVerletList(nstlist);
2227 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2229 impl_->addReplicaExchange(params);
2233 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2235 impl_->addNonBonded(nbpu_opt);
2239 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2241 // The builder method may become more general in the future, but in this version,
2242 // parameters for PME electrostatics are both required and the only parameters
2244 if (pme_opt && pme_fft_opt)
2246 impl_->addPME(pme_opt, pme_fft_opt);
2251 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2256 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2258 impl_->addBondedTaskAssignment(bonded_opt);
2262 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2264 impl_->addUpdateTaskAssignment(update_opt);
2268 Mdrunner MdrunnerBuilder::build()
2270 return impl_->build();
2273 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2275 impl_->addHardwareOptions(hardwareOptions);
2279 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2281 impl_->addFilenames(filenames);
2285 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2287 impl_->addOutputEnvironment(outputEnvironment);
2291 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2293 impl_->addLogFile(logFileHandle);
2297 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2299 impl_->addStopHandlerBuilder(std::move(builder));
2303 MdrunnerBuilder& MdrunnerBuilder::addInput(SimulationInputHandle input)
2305 impl_->addInput(std::move(input));
2309 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2311 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;