Two sets of coefficients for Coulomb FEP PME on GPU
[gromacs.git] / src / gromacs / mdrun / runner.cpp
blobb90bc9aea50cf0cafa4a174e934731c015eaf20a
1 /*
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.
37 /*! \internal \file
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
44 #include "gmxpre.h"
46 #include "runner.h"
48 #include "config.h"
50 #include <cassert>
51 #include <cinttypes>
52 #include <csignal>
53 #include <cstdlib>
54 #include <cstring>
56 #include <algorithm>
57 #include <memory>
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/device_context.h"
77 #include "gromacs/gpu_utils/device_stream_manager.h"
78 #include "gromacs/gpu_utils/gpu_utils.h"
79 #include "gromacs/hardware/cpuinfo.h"
80 #include "gromacs/hardware/detecthardware.h"
81 #include "gromacs/hardware/printhardware.h"
82 #include "gromacs/imd/imd.h"
83 #include "gromacs/listed_forces/disre.h"
84 #include "gromacs/listed_forces/gpubonded.h"
85 #include "gromacs/listed_forces/listed_forces.h"
86 #include "gromacs/listed_forces/orires.h"
87 #include "gromacs/math/functions.h"
88 #include "gromacs/math/utilities.h"
89 #include "gromacs/math/vec.h"
90 #include "gromacs/mdlib/boxdeformation.h"
91 #include "gromacs/mdlib/broadcaststructs.h"
92 #include "gromacs/mdlib/calc_verletbuf.h"
93 #include "gromacs/mdlib/dispersioncorrection.h"
94 #include "gromacs/mdlib/enerdata_utils.h"
95 #include "gromacs/mdlib/force.h"
96 #include "gromacs/mdlib/forcerec.h"
97 #include "gromacs/mdlib/gmx_omp_nthreads.h"
98 #include "gromacs/mdlib/makeconstraints.h"
99 #include "gromacs/mdlib/md_support.h"
100 #include "gromacs/mdlib/mdatoms.h"
101 #include "gromacs/mdlib/sighandler.h"
102 #include "gromacs/mdlib/stophandler.h"
103 #include "gromacs/mdlib/tgroup.h"
104 #include "gromacs/mdlib/updategroups.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdrun/mdmodules.h"
107 #include "gromacs/mdrun/simulationcontext.h"
108 #include "gromacs/mdrunutility/handlerestart.h"
109 #include "gromacs/mdrunutility/logging.h"
110 #include "gromacs/mdrunutility/multisim.h"
111 #include "gromacs/mdrunutility/printtime.h"
112 #include "gromacs/mdrunutility/threadaffinity.h"
113 #include "gromacs/mdtypes/commrec.h"
114 #include "gromacs/mdtypes/enerdata.h"
115 #include "gromacs/mdtypes/fcdata.h"
116 #include "gromacs/mdtypes/forcerec.h"
117 #include "gromacs/mdtypes/group.h"
118 #include "gromacs/mdtypes/inputrec.h"
119 #include "gromacs/mdtypes/interaction_const.h"
120 #include "gromacs/mdtypes/md_enums.h"
121 #include "gromacs/mdtypes/mdatom.h"
122 #include "gromacs/mdtypes/mdrunoptions.h"
123 #include "gromacs/mdtypes/observableshistory.h"
124 #include "gromacs/mdtypes/simulation_workload.h"
125 #include "gromacs/mdtypes/state.h"
126 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
127 #include "gromacs/modularsimulator/modularsimulator.h"
128 #include "gromacs/nbnxm/gpu_data_mgmt.h"
129 #include "gromacs/nbnxm/nbnxm.h"
130 #include "gromacs/nbnxm/pairlist_tuning.h"
131 #include "gromacs/pbcutil/pbc.h"
132 #include "gromacs/pulling/output.h"
133 #include "gromacs/pulling/pull.h"
134 #include "gromacs/pulling/pull_rotation.h"
135 #include "gromacs/restraint/manager.h"
136 #include "gromacs/restraint/restraintmdmodule.h"
137 #include "gromacs/restraint/restraintpotential.h"
138 #include "gromacs/swap/swapcoords.h"
139 #include "gromacs/taskassignment/decidegpuusage.h"
140 #include "gromacs/taskassignment/decidesimulationworkload.h"
141 #include "gromacs/taskassignment/resourcedivision.h"
142 #include "gromacs/taskassignment/taskassignment.h"
143 #include "gromacs/taskassignment/usergpuids.h"
144 #include "gromacs/timing/gpu_timing.h"
145 #include "gromacs/timing/wallcycle.h"
146 #include "gromacs/timing/wallcyclereporting.h"
147 #include "gromacs/topology/mtop_util.h"
148 #include "gromacs/trajectory/trajectoryframe.h"
149 #include "gromacs/utility/basenetwork.h"
150 #include "gromacs/utility/cstringutil.h"
151 #include "gromacs/utility/exceptions.h"
152 #include "gromacs/utility/fatalerror.h"
153 #include "gromacs/utility/filestream.h"
154 #include "gromacs/utility/gmxassert.h"
155 #include "gromacs/utility/gmxmpi.h"
156 #include "gromacs/utility/keyvaluetree.h"
157 #include "gromacs/utility/logger.h"
158 #include "gromacs/utility/loggerbuilder.h"
159 #include "gromacs/utility/mdmodulenotification.h"
160 #include "gromacs/utility/physicalnodecommunicator.h"
161 #include "gromacs/utility/pleasecite.h"
162 #include "gromacs/utility/programcontext.h"
163 #include "gromacs/utility/smalloc.h"
164 #include "gromacs/utility/stringutil.h"
166 #include "isimulator.h"
167 #include "membedholder.h"
168 #include "replicaexchange.h"
169 #include "simulatorbuilder.h"
171 #if GMX_FAHCORE
172 # include "corewrap.h"
173 #endif
175 namespace gmx
179 /*! \brief Manage any development feature flag variables encountered
181 * The use of dev features indicated by environment variables is
182 * logged in order to ensure that runs with such features enabled can
183 * be identified from their log and standard output. Any cross
184 * dependencies are also checked, and if unsatisfied, a fatal error
185 * issued.
187 * Note that some development features overrides are applied already here:
188 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
190 * \param[in] mdlog Logger object.
191 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
192 * \param[in] pmeRunMode The PME run mode for this run
193 * \returns The object populated with development feature flags.
195 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
196 const bool useGpuForNonbonded,
197 const PmeRunMode pmeRunMode)
199 DevelopmentFeatureFlags devFlags;
201 // Some builds of GCC 5 give false positive warnings that these
202 // getenv results are ignored when clearly they are used.
203 #pragma GCC diagnostic push
204 #pragma GCC diagnostic ignored "-Wunused-result"
205 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
206 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
207 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
208 devFlags.enableGpuHaloExchange =
209 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
210 devFlags.enableGpuPmePPComm =
211 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
212 #pragma GCC diagnostic pop
214 if (devFlags.enableGpuBufferOps)
216 GMX_LOG(mdlog.warning)
217 .asParagraph()
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)
226 .asParagraph()
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)
240 .asParagraph()
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)
247 .asParagraph()
248 .appendTextFormatted(
249 "This run has requested the 'GPU halo exchange' feature, enabled by "
250 "the "
251 "GMX_GPU_DD_COMMS environment variable.");
253 else
255 GMX_LOG(mdlog.warning)
256 .asParagraph()
257 .appendTextFormatted(
258 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
259 "halo exchange' feature will not be enabled as nonbonded interactions "
260 "are not offloaded.");
261 devFlags.enableGpuHaloExchange = false;
265 if (devFlags.enableGpuPmePPComm)
267 if (pmeRunMode == PmeRunMode::GPU)
269 GMX_LOG(mdlog.warning)
270 .asParagraph()
271 .appendTextFormatted(
272 "This run uses the 'GPU PME-PP communications' feature, enabled "
273 "by the GMX_GPU_PME_PP_COMMS environment variable.");
275 else
277 std::string clarification;
278 if (pmeRunMode == PmeRunMode::Mixed)
280 clarification =
281 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
282 "mode).";
284 else
286 clarification = "PME is not offloaded to the GPU.";
288 GMX_LOG(mdlog.warning)
289 .asParagraph()
290 .appendText(
291 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
292 "'GPU PME-PP communications' feature was not enabled as "
293 + clarification);
294 devFlags.enableGpuPmePPComm = false;
298 return devFlags;
301 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
303 * Used to ensure that the master thread does not modify mdrunner during copy
304 * on the spawned threads. */
305 static void threadMpiMdrunnerAccessBarrier()
307 #if GMX_THREAD_MPI
308 MPI_Barrier(MPI_COMM_WORLD);
309 #endif
312 Mdrunner Mdrunner::cloneOnSpawnedThread() const
314 auto newRunner = Mdrunner(std::make_unique<MDModules>());
316 // All runners in the same process share a restraint manager resource because it is
317 // part of the interface to the client code, which is associated only with the
318 // original thread. Handles to the same resources can be obtained by copy.
320 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
323 // Copy members of master runner.
324 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
325 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
326 newRunner.hw_opt = hw_opt;
327 newRunner.filenames = filenames;
329 newRunner.oenv = oenv;
330 newRunner.mdrunOptions = mdrunOptions;
331 newRunner.domdecOptions = domdecOptions;
332 newRunner.nbpu_opt = nbpu_opt;
333 newRunner.pme_opt = pme_opt;
334 newRunner.pme_fft_opt = pme_fft_opt;
335 newRunner.bonded_opt = bonded_opt;
336 newRunner.update_opt = update_opt;
337 newRunner.nstlist_cmdline = nstlist_cmdline;
338 newRunner.replExParams = replExParams;
339 newRunner.pforce = pforce;
340 // Give the spawned thread the newly created valid communicator
341 // for the simulation.
342 newRunner.communicator = MPI_COMM_WORLD;
343 newRunner.ms = ms;
344 newRunner.startingBehavior = startingBehavior;
345 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
347 threadMpiMdrunnerAccessBarrier();
349 return newRunner;
352 /*! \brief The callback used for running on spawned threads.
354 * Obtains the pointer to the master mdrunner object from the one
355 * argument permitted to the thread-launch API call, copies it to make
356 * a new runner for this thread, reinitializes necessary data, and
357 * proceeds to the simulation. */
358 static void mdrunner_start_fn(const void* arg)
362 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
363 /* copy the arg list to make sure that it's thread-local. This
364 doesn't copy pointed-to items, of course; fnm, cr and fplog
365 are reset in the call below, all others should be const. */
366 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
367 mdrunner.mdrunner();
369 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
373 void Mdrunner::spawnThreads(int numThreadsToLaunch)
375 #if GMX_THREAD_MPI
376 /* now spawn new threads that start mdrunner_start_fn(), while
377 the main thread returns. Thread affinity is handled later. */
378 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
379 static_cast<const void*>(this))
380 != TMPI_SUCCESS)
382 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
385 // Give the master thread the newly created valid communicator for
386 // the simulation.
387 communicator = MPI_COMM_WORLD;
388 threadMpiMdrunnerAccessBarrier();
389 #else
390 GMX_UNUSED_VALUE(numThreadsToLaunch);
391 GMX_UNUSED_VALUE(mdrunner_start_fn);
392 #endif
395 } // namespace gmx
397 /*! \brief Initialize variables for Verlet scheme simulation */
398 static void prepare_verlet_scheme(FILE* fplog,
399 t_commrec* cr,
400 t_inputrec* ir,
401 int nstlist_cmdline,
402 const gmx_mtop_t* mtop,
403 const matrix box,
404 bool makeGpuPairList,
405 const gmx::CpuInfo& cpuinfo)
407 // We checked the cut-offs in grompp, but double-check here.
408 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
409 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
411 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
412 "With Verlet lists and PME we should have rcoulomb>=rvdw");
414 else
416 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
417 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
419 /* For NVE simulations, we will retain the initial list buffer */
420 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
422 /* Update the Verlet buffer size for the current run setup */
424 /* Here we assume SIMD-enabled kernels are being used. But as currently
425 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
426 * and 4x2 gives a larger buffer than 4x4, this is ok.
428 ListSetupType listType =
429 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
430 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
432 const real rlist_new =
433 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
435 if (rlist_new != ir->rlist)
437 if (fplog != nullptr)
439 fprintf(fplog,
440 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
441 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
443 ir->rlist = rlist_new;
447 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
449 gmx_fatal(FARGS, "Can not set nstlist without %s",
450 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
453 if (EI_DYNAMICS(ir->eI))
455 /* Set or try nstlist values */
456 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
460 /*! \brief Override the nslist value in inputrec
462 * with value passed on the command line (if any)
464 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
466 assert(ir);
468 /* override with anything else than the default -2 */
469 if (nsteps_cmdline > -2)
471 char sbuf_steps[STEPSTRSIZE];
472 char sbuf_msg[STRLEN];
474 ir->nsteps = nsteps_cmdline;
475 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
477 sprintf(sbuf_msg,
478 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
479 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
481 else
483 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
484 gmx_step_str(nsteps_cmdline, sbuf_steps));
487 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
489 else if (nsteps_cmdline < -2)
491 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
493 /* Do nothing if nsteps_cmdline == -2 */
496 namespace gmx
499 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
501 * If not, and if a warning may be issued, logs a warning about
502 * falling back to CPU code. With thread-MPI, only the first
503 * call to this function should have \c issueWarning true. */
504 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
506 bool gpuIsUseful = true;
507 std::string warning;
509 if (ir.opts.ngener - ir.nwall > 1)
511 /* The GPU code does not support more than one energy group.
512 * If the user requested GPUs explicitly, a fatal error is given later.
514 gpuIsUseful = false;
515 warning =
516 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
517 "For better performance, run on the GPU without energy groups and then do "
518 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
521 if (EI_TPI(ir.eI))
523 gpuIsUseful = false;
524 warning = "TPI is not implemented for GPUs.";
527 if (!gpuIsUseful && issueWarning)
529 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
532 return gpuIsUseful;
535 //! Initializes the logger for mdrun.
536 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
538 gmx::LoggerBuilder builder;
539 if (fplog != nullptr)
541 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
543 if (isSimulationMasterRank)
545 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
547 return builder.build();
550 //! Make a TaskTarget from an mdrun argument string.
551 static TaskTarget findTaskTarget(const char* optionString)
553 TaskTarget returnValue = TaskTarget::Auto;
555 if (strncmp(optionString, "auto", 3) == 0)
557 returnValue = TaskTarget::Auto;
559 else if (strncmp(optionString, "cpu", 3) == 0)
561 returnValue = TaskTarget::Cpu;
563 else if (strncmp(optionString, "gpu", 3) == 0)
565 returnValue = TaskTarget::Gpu;
567 else
569 GMX_ASSERT(false, "Option string should have been checked for sanity already");
572 return returnValue;
575 //! Finish run, aggregate data to print performance info.
576 static void finish_run(FILE* fplog,
577 const gmx::MDLogger& mdlog,
578 const t_commrec* cr,
579 const t_inputrec* inputrec,
580 t_nrnb nrnb[],
581 gmx_wallcycle_t wcycle,
582 gmx_walltime_accounting_t walltime_accounting,
583 nonbonded_verlet_t* nbv,
584 const gmx_pme_t* pme,
585 gmx_bool bWriteStat)
587 double delta_t = 0;
588 double nbfs = 0, mflop = 0;
589 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
590 elapsed_time_over_all_threads_over_all_ranks;
591 /* Control whether it is valid to print a report. Only the
592 simulation master may print, but it should not do so if the run
593 terminated e.g. before a scheduled reset step. This is
594 complicated by the fact that PME ranks are unaware of the
595 reason why they were sent a pmerecvqxFINISH. To avoid
596 communication deadlocks, we always do the communication for the
597 report, even if we've decided not to write the report, because
598 how long it takes to finish the run is not important when we've
599 decided not to report on the simulation performance.
601 Further, we only report performance for dynamical integrators,
602 because those are the only ones for which we plan to
603 consider doing any optimizations. */
604 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
606 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
608 GMX_LOG(mdlog.warning)
609 .asParagraph()
610 .appendText("Simulation ended prematurely, no performance report will be written.");
611 printReport = false;
614 t_nrnb* nrnb_tot;
615 std::unique_ptr<t_nrnb> nrnbTotalStorage;
616 if (cr->nnodes > 1)
618 nrnbTotalStorage = std::make_unique<t_nrnb>();
619 nrnb_tot = nrnbTotalStorage.get();
620 #if GMX_MPI
621 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
622 #endif
624 else
626 nrnb_tot = nrnb;
629 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
630 elapsed_time_over_all_threads =
631 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
632 if (cr->nnodes > 1)
634 #if GMX_MPI
635 /* reduce elapsed_time over all MPI ranks in the current simulation */
636 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
637 cr->mpi_comm_mysim);
638 elapsed_time_over_all_ranks /= cr->nnodes;
639 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
640 * current simulation. */
641 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
642 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
643 #endif
645 else
647 elapsed_time_over_all_ranks = elapsed_time;
648 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
651 if (printReport)
653 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
656 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
658 print_dd_statistics(cr, inputrec, fplog);
661 /* TODO Move the responsibility for any scaling by thread counts
662 * to the code that handled the thread region, so that there's a
663 * mechanism to keep cycle counting working during the transition
664 * to task parallelism. */
665 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
666 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
667 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
668 nthreads_pp, nthreads_pme);
669 auto cycle_sum(wallcycle_sum(cr, wcycle));
671 if (printReport)
673 auto nbnxn_gpu_timings =
674 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
675 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
677 if (pme_gpu_task_enabled(pme))
679 pme_gpu_get_timings(pme, &pme_gpu_timings);
681 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
682 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
683 &pme_gpu_timings);
685 if (EI_DYNAMICS(inputrec->eI))
687 delta_t = inputrec->delta_t;
690 if (fplog)
692 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
693 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
694 delta_t, nbfs, mflop);
696 if (bWriteStat)
698 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
699 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
700 delta_t, nbfs, mflop);
705 int Mdrunner::mdrunner()
707 matrix box;
708 t_forcerec* fr = nullptr;
709 real ewaldcoeff_q = 0;
710 real ewaldcoeff_lj = 0;
711 int nChargePerturbed = -1, nTypePerturbed = 0;
712 gmx_wallcycle_t wcycle;
713 gmx_walltime_accounting_t walltime_accounting = nullptr;
714 MembedHolder membedHolder(filenames.size(), filenames.data());
715 gmx_hw_info_t* hwinfo = nullptr;
717 /* CAUTION: threads may be started later on in this function, so
718 cr doesn't reflect the final parallel state right now */
719 gmx_mtop_t mtop;
721 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
722 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
723 const bool doRerun = mdrunOptions.rerun;
725 // Handle task-assignment related user options.
726 EmulateGpuNonbonded emulateGpuNonbonded =
727 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
729 std::vector<int> userGpuTaskAssignment;
732 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
734 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
735 auto nonbondedTarget = findTaskTarget(nbpu_opt);
736 auto pmeTarget = findTaskTarget(pme_opt);
737 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
738 auto bondedTarget = findTaskTarget(bonded_opt);
739 auto updateTarget = findTaskTarget(update_opt);
741 FILE* fplog = nullptr;
742 // If we are appending, we don't write log output because we need
743 // to check that the old log file matches what the checkpoint file
744 // expects. Otherwise, we should start to write log output now if
745 // there is a file ready for it.
746 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
748 fplog = gmx_fio_getfp(logFileHandle);
750 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
751 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
752 gmx::MDLogger mdlog(logOwner.logger());
754 // TODO The thread-MPI master rank makes a working
755 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
756 // after the threads have been launched. This works because no use
757 // is made of that communicator until after the execution paths
758 // have rejoined. But it is likely that we can improve the way
759 // this is expressed, e.g. by expressly running detection only the
760 // master rank for thread-MPI, rather than relying on the mutex
761 // and reference count.
762 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
763 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
765 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
767 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
769 // Print citation requests after all software/hardware printing
770 pleaseCiteGromacs(fplog);
772 // TODO Replace this by unique_ptr once t_inputrec is C++
773 t_inputrec inputrecInstance;
774 t_inputrec* inputrec = nullptr;
775 std::unique_ptr<t_state> globalState;
777 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
779 if (isSimulationMasterRank)
781 /* Only the master rank has the global state */
782 globalState = std::make_unique<t_state>();
784 /* Read (nearly) all data required for the simulation
785 * and keep the partly serialized tpr contents to send to other ranks later
787 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
788 &inputrecInstance, globalState.get(), &mtop);
789 inputrec = &inputrecInstance;
792 /* Check and update the hardware options for internal consistency */
793 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
794 inputrec);
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, 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
823 * correctly. */
824 hw_opt.nthreads_tmpi =
825 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
826 inputrec, &mtop, mdlog, membedHolder.doMembed());
828 // Now start the threads for thread MPI.
829 spawnThreads(hw_opt.nthreads_tmpi);
830 // The spawned threads enter mdrunner() and execution of
831 // master and spawned threads joins at the end of this block.
832 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
835 GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
836 "Must have valid world communicator unless running a multi-simulation");
837 CommrecHandle crHandle = init_commrec(communicator);
838 t_commrec* cr = crHandle.get();
839 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
841 if (PAR(cr))
843 /* now broadcast everything to the non-master nodes/threads: */
844 if (!isSimulationMasterRank)
846 inputrec = &inputrecInstance;
848 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec, &mtop,
849 partialDeserializedTpr.get());
851 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
852 partialDeserializedTpr.reset(nullptr);
854 // Now the number of ranks is known to all ranks, and each knows
855 // the inputrec read by the master rank. The ranks can now all run
856 // the task-deciding functions and will agree on the result
857 // without needing to communicate.
858 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
860 // Note that these variables describe only their own node.
862 // Note that when bonded interactions run on a GPU they always run
863 // alongside a nonbonded task, so do not influence task assignment
864 // even though they affect the force calculation workload.
865 bool useGpuForNonbonded = false;
866 bool useGpuForPme = false;
867 bool useGpuForBonded = false;
868 bool useGpuForUpdate = false;
869 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
872 // It's possible that there are different numbers of GPUs on
873 // different nodes, which is the user's responsibility to
874 // handle. If unsuitable, we will notice that during task
875 // assignment.
876 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
877 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
878 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
879 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
880 useGpuForPme = decideWhetherToUseGpusForPme(
881 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec,
882 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
883 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
884 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
885 useGpuForBonded = decideWhetherToUseGpusForBonded(
886 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
887 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
888 domdecOptions.numPmeRanks, gpusWereDetected);
890 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
892 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
894 // Initialize development feature flags that enabled by environment variable
895 // and report those features that are enabled.
896 const DevelopmentFeatureFlags devFlags =
897 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
899 const bool useModularSimulator =
900 checkUseModularSimulator(false, inputrec, doRerun, mtop, ms, replExParams, nullptr,
901 doEssentialDynamics, membedHolder.doMembed());
903 // Build restraints.
904 // TODO: hide restraint implementation details from Mdrunner.
905 // There is nothing unique about restraints at this point as far as the
906 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
907 // factory functions from the SimulationContext on which to call mdModules_->add().
908 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
909 for (auto&& restraint : restraintManager_->getRestraints())
911 auto module = RestraintMDModule::create(restraint, restraint->sites());
912 mdModules_->add(std::move(module));
915 // TODO: Error handling
916 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
917 // now that the MdModules know their options, they know which callbacks to sign up to
918 mdModules_->subscribeToSimulationSetupNotifications();
919 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
921 if (inputrec->internalParameters != nullptr)
923 mdModulesNotifier.notify(*inputrec->internalParameters);
926 if (fplog != nullptr)
928 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
929 fprintf(fplog, "\n");
932 if (SIMMASTER(cr))
934 /* In rerun, set velocities to zero if present */
935 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
937 // rerun does not use velocities
938 GMX_LOG(mdlog.info)
939 .asParagraph()
940 .appendText(
941 "Rerun trajectory contains velocities. Rerun does only evaluate "
942 "potential energy and forces. The velocities will be ignored.");
943 for (int i = 0; i < globalState->natoms; i++)
945 clear_rvec(globalState->v[i]);
947 globalState->flags &= ~(1 << estV);
950 /* now make sure the state is initialized and propagated */
951 set_state_entries(globalState.get(), inputrec, useModularSimulator);
954 /* NM and TPI parallelize over force/energy calculations, not atoms,
955 * so we need to initialize and broadcast the global state.
957 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
959 if (!MASTER(cr))
961 globalState = std::make_unique<t_state>();
963 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
964 globalState.get());
967 /* A parallel command line option consistency check that we can
968 only do after any threads have started. */
969 if (!PAR(cr)
970 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
971 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
973 gmx_fatal(FARGS,
974 "The -dd or -npme option request a parallel simulation, "
975 #if !GMX_MPI
976 "but %s was compiled without threads or MPI enabled",
977 output_env_get_program_display_name(oenv));
978 #elif GMX_THREAD_MPI
979 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
980 #else
981 "but %s was not started through mpirun/mpiexec or only one rank was requested "
982 "through mpirun/mpiexec",
983 output_env_get_program_display_name(oenv));
984 #endif
987 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
989 gmx_fatal(FARGS,
990 "The .mdp file specified an energy mininization or normal mode algorithm, and "
991 "these are not compatible with mdrun -rerun");
994 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
996 if (domdecOptions.numPmeRanks > 0)
998 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
999 "PME-only ranks are requested, but the system does not use PME "
1000 "for electrostatics or LJ");
1003 domdecOptions.numPmeRanks = 0;
1006 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1008 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1009 * improve performance with many threads per GPU, since our OpenMP
1010 * scaling is bad, but it's difficult to automate the setup.
1012 domdecOptions.numPmeRanks = 0;
1014 if (useGpuForPme)
1016 if (domdecOptions.numPmeRanks < 0)
1018 domdecOptions.numPmeRanks = 0;
1019 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1021 else
1023 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1024 "PME GPU decomposition is not supported");
1028 #if GMX_FAHCORE
1029 if (MASTER(cr))
1031 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1033 #endif
1035 /* NMR restraints must be initialized before load_checkpoint,
1036 * since with time averaging the history is added to t_state.
1037 * For proper consistency check we therefore need to extend
1038 * t_state here.
1039 * So the PME-only nodes (if present) will also initialize
1040 * the distance restraints.
1043 /* This needs to be called before read_checkpoint to extend the state */
1044 t_disresdata* disresdata;
1045 snew(disresdata, 1);
1046 init_disres(fplog, &mtop, inputrec, DisResRunMode::MDRun, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1047 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1048 globalState.get(), replExParams.exchangeInterval > 0);
1050 t_oriresdata* oriresdata;
1051 snew(oriresdata, 1);
1052 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), oriresdata);
1054 auto deform = prepareBoxDeformation(globalState->box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1055 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1056 cr->mpi_comm_mygroup, *inputrec);
1058 ObservablesHistory observablesHistory = {};
1060 if (startingBehavior != StartingBehavior::NewSimulation)
1062 /* Check if checkpoint file exists before doing continuation.
1063 * This way we can use identical input options for the first and subsequent runs...
1065 if (mdrunOptions.numStepsCommandline > -2)
1067 /* Temporarily set the number of steps to unlimited to avoid
1068 * triggering the nsteps check in load_checkpoint().
1069 * This hack will go away soon when the -nsteps option is removed.
1071 inputrec->nsteps = -1;
1074 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1075 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1076 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1078 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1080 // Now we can start normal logging to the truncated log file.
1081 fplog = gmx_fio_getfp(logFileHandle);
1082 prepareLogAppending(fplog);
1083 logOwner = buildLogger(fplog, MASTER(cr));
1084 mdlog = logOwner.logger();
1088 if (mdrunOptions.numStepsCommandline > -2)
1090 GMX_LOG(mdlog.info)
1091 .asParagraph()
1092 .appendText(
1093 "The -nsteps functionality is deprecated, and may be removed in a future "
1094 "version. "
1095 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1096 "file field.");
1098 /* override nsteps with value set on the commandline */
1099 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1101 if (SIMMASTER(cr))
1103 copy_mat(globalState->box, box);
1106 if (PAR(cr))
1108 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1111 if (inputrec->cutoff_scheme != ecutsVERLET)
1113 gmx_fatal(FARGS,
1114 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1115 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1117 /* Update rlist and nstlist. */
1118 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1119 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1120 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1122 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1123 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1124 *hwinfo->cpuInfo);
1126 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1127 // This builder is necessary while we have multi-part construction
1128 // of DD. Before DD is constructed, we use the existence of
1129 // the builder object to indicate that further construction of DD
1130 // is needed.
1131 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1132 if (useDomainDecomposition)
1134 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1135 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1136 positionsFromStatePointer(globalState.get()));
1138 else
1140 /* PME, if used, is done on all nodes with 1D decomposition */
1141 cr->nnodes = cr->sizeOfDefaultCommunicator;
1142 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1143 cr->nodeid = cr->rankInDefaultCommunicator;
1144 cr->npmenodes = 0;
1145 cr->duty = (DUTY_PP | DUTY_PME);
1147 if (inputrec->pbcType == PbcType::Screw)
1149 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1153 // Produce the task assignment for this rank - done after DD is constructed
1154 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1155 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1156 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1157 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1158 // TODO cr->duty & DUTY_PME should imply that a PME
1159 // algorithm is active, but currently does not.
1160 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1162 // Get the device handles for the modules, nullptr when no task is assigned.
1163 int deviceId = -1;
1164 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1166 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1167 bool useTiming = true;
1168 if (GMX_GPU == GMX_GPU_CUDA)
1170 /* WARNING: CUDA timings are incorrect with multiple streams.
1171 * This is the main reason why they are disabled by default.
1173 // TODO: Consider turning on by default when we can detect nr of streams.
1174 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1176 else if (GMX_GPU == GMX_GPU_OPENCL)
1178 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1181 // TODO Currently this is always built, yet DD partition code
1182 // checks if it is built before using it. Probably it should
1183 // become an MDModule that is made only when another module
1184 // requires it (e.g. pull, CompEl, density fitting), so that we
1185 // don't update the local atom sets unilaterally every step.
1186 LocalAtomSetManager atomSets;
1187 if (ddBuilder)
1189 // TODO Pass the GPU streams to ddBuilder to use in buffer
1190 // transfers (e.g. halo exchange)
1191 cr->dd = ddBuilder->build(&atomSets);
1192 // The builder's job is done, so destruct it
1193 ddBuilder.reset(nullptr);
1194 // Note that local state still does not exist yet.
1197 // The GPU update is decided here because we need to know whether the constraints or
1198 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1199 // defined). This is only known after DD is initialized, hence decision on using GPU
1200 // update is done so late.
1203 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1205 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1206 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1207 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1208 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1209 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1211 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1213 const bool printHostName = (cr->nnodes > 1);
1214 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1216 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1218 if (deviceInfo != nullptr)
1220 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1222 dd_setup_dlb_resource_sharing(cr, deviceId);
1224 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1225 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1226 useGpuForUpdate, useTiming);
1229 // If the user chose a task assignment, give them some hints
1230 // where appropriate.
1231 if (!userGpuTaskAssignment.empty())
1233 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1236 if (PAR(cr))
1238 /* After possible communicator splitting in make_dd_communicators.
1239 * we can set up the intra/inter node communication.
1241 gmx_setup_nodecomm(fplog, cr);
1244 #if GMX_MPI
1245 if (isMultiSim(ms))
1247 GMX_LOG(mdlog.warning)
1248 .asParagraph()
1249 .appendTextFormatted(
1250 "This is simulation %d out of %d running as a composite GROMACS\n"
1251 "multi-simulation job. Setup for this simulation:\n",
1252 ms->simulationIndex_, ms->numSimulations_);
1254 GMX_LOG(mdlog.warning)
1255 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1256 # if GMX_THREAD_MPI
1257 cr->nnodes == 1 ? "thread" : "threads"
1258 # else
1259 cr->nnodes == 1 ? "process" : "processes"
1260 # endif
1262 fflush(stderr);
1263 #endif
1265 // If mdrun -pin auto honors any affinity setting that already
1266 // exists. If so, it is nice to provide feedback about whether
1267 // that existing affinity setting was from OpenMP or something
1268 // else, so we run this code both before and after we initialize
1269 // the OpenMP support.
1270 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1271 /* Check and update the number of OpenMP threads requested */
1272 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1273 pmeRunMode, mtop, *inputrec);
1275 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1276 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1278 // Enable FP exception detection, but not in
1279 // Release mode and not for compilers with known buggy FP
1280 // exception support (clang with any optimization) or suspected
1281 // buggy FP exception support (gcc 7.* with optimization).
1282 #if !defined NDEBUG \
1283 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1284 && defined __OPTIMIZE__)
1285 const bool bEnableFPE = true;
1286 #else
1287 const bool bEnableFPE = false;
1288 #endif
1289 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1290 if (bEnableFPE)
1292 gmx_feenableexcept();
1295 /* Now that we know the setup is consistent, check for efficiency */
1296 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1297 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1299 /* getting number of PP/PME threads on this MPI / tMPI rank.
1300 PME: env variable should be read only on one node to make sure it is
1301 identical everywhere;
1303 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1304 : gmx_omp_nthreads_get(emntPME);
1305 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1306 physicalNodeComm, mdlog);
1308 // Enable Peer access between GPUs where available
1309 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1310 // any of the GPU communication features are active.
1311 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1312 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1314 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1317 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1319 /* Before setting affinity, check whether the affinity has changed
1320 * - which indicates that probably the OpenMP library has changed it
1321 * since we first checked).
1323 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1325 int numThreadsOnThisNode, intraNodeThreadOffset;
1326 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1327 &intraNodeThreadOffset);
1329 /* Set the CPU affinity */
1330 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1331 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1334 if (mdrunOptions.timingOptions.resetStep > -1)
1336 GMX_LOG(mdlog.info)
1337 .asParagraph()
1338 .appendText(
1339 "The -resetstep functionality is deprecated, and may be removed in a "
1340 "future version.");
1342 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1344 if (PAR(cr))
1346 /* Master synchronizes its value of reset_counters with all nodes
1347 * including PME only nodes */
1348 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1349 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1350 wcycle_set_reset_counters(wcycle, reset_counters);
1353 // Membrane embedding must be initialized before we call init_forcerec()
1354 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1355 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1357 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1358 std::unique_ptr<MDAtoms> mdAtoms;
1359 std::unique_ptr<VirtualSitesHandler> vsite;
1360 std::unique_ptr<GpuBonded> gpuBonded;
1362 t_nrnb nrnb;
1363 if (thisRankHasDuty(cr, DUTY_PP))
1365 mdModulesNotifier.notify(*cr);
1366 mdModulesNotifier.notify(&atomSets);
1367 mdModulesNotifier.notify(inputrec->pbcType);
1368 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1369 /* Initiate forcerecord */
1370 fr = new t_forcerec;
1371 fr->forceProviders = mdModules_->initForceProviders();
1372 init_forcerec(fplog, mdlog, fr, inputrec, &mtop, cr, box,
1373 opt2fn("-table", filenames.size(), filenames.data()),
1374 opt2fn("-tablep", filenames.size(), filenames.data()),
1375 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1376 // Dirty hack, for fixing disres and orires should be made mdmodules
1377 fr->listedForces->fcdata().disres = disresdata;
1378 fr->listedForces->fcdata().orires = oriresdata;
1380 // Save a handle to device stream manager to use elsewhere in the code
1381 // TODO: Forcerec is not a correct place to store it.
1382 fr->deviceStreamManager = deviceStreamManager.get();
1384 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1386 GMX_RELEASE_ASSERT(
1387 deviceStreamManager != nullptr,
1388 "GPU device stream manager should be valid in order to use PME-PP direct "
1389 "communications.");
1390 GMX_RELEASE_ASSERT(
1391 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1392 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1393 "communications.");
1394 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1395 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1396 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1399 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, useGpuForNonbonded,
1400 deviceStreamManager.get(), &mtop, box, wcycle);
1401 // TODO: Move the logic below to a GPU bonded builder
1402 if (useGpuForBonded)
1404 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1405 "GPU device stream manager should be valid in order to use GPU "
1406 "version of bonded forces.");
1407 gpuBonded = std::make_unique<GpuBonded>(
1408 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1409 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1410 fr->gpuBonded = gpuBonded.get();
1413 /* Initialize the mdAtoms structure.
1414 * mdAtoms is not filled with atom data,
1415 * as this can not be done now with domain decomposition.
1417 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1418 if (globalState && thisRankHasPmeGpuTask)
1420 // The pinning of coordinates in the global state object works, because we only use
1421 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1422 // points to the global state object without DD.
1423 // FIXME: MD and EM separately set up the local state - this should happen in the same
1424 // function, which should also perform the pinning.
1425 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1428 /* Initialize the virtual site communication */
1429 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1431 calc_shifts(box, fr->shift_vec);
1433 /* With periodic molecules the charge groups should be whole at start up
1434 * and the virtual sites should not be far from their proper positions.
1436 if (!inputrec->bContinuation && MASTER(cr)
1437 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1439 /* Make molecules whole at start of run */
1440 if (fr->pbcType != PbcType::No)
1442 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1444 if (vsite)
1446 /* Correct initial vsite positions are required
1447 * for the initial distribution in the domain decomposition
1448 * and for the initial shell prediction.
1450 constructVirtualSitesGlobal(mtop, globalState->x);
1454 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1456 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1457 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1460 else
1462 /* This is a PME only node */
1464 GMX_ASSERT(globalState == nullptr,
1465 "We don't need the state on a PME only rank and expect it to be unitialized");
1467 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1468 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1471 gmx_pme_t* sepPmeData = nullptr;
1472 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1473 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1474 "Double-checking that only PME-only ranks have no forcerec");
1475 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1477 // TODO should live in ewald module once its testing is improved
1479 // Later, this program could contain kernels that might be later
1480 // re-used as auto-tuning progresses, or subsequent simulations
1481 // are invoked.
1482 PmeGpuProgramStorage pmeGpuProgram;
1483 if (thisRankHasPmeGpuTask)
1485 GMX_RELEASE_ASSERT(
1486 (deviceStreamManager != nullptr),
1487 "GPU device stream manager should be initialized in order to use GPU for PME.");
1488 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1489 "GPU device should be initialized in order to use GPU for PME.");
1490 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1493 /* Initiate PME if necessary,
1494 * either on all nodes or on dedicated PME nodes only. */
1495 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1497 if (mdAtoms && mdAtoms->mdatoms())
1499 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1500 if (EVDW_PME(inputrec->vdwtype))
1502 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1505 if (cr->npmenodes > 0)
1507 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1508 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1509 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1512 if (thisRankHasDuty(cr, DUTY_PME))
1516 // TODO: This should be in the builder.
1517 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1518 "Device stream manager should be valid in order to use GPU "
1519 "version of PME.");
1520 GMX_RELEASE_ASSERT(
1521 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1522 "GPU PME stream should be valid in order to use GPU version of PME.");
1524 const DeviceContext* deviceContext =
1525 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1526 const DeviceStream* pmeStream =
1527 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1529 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1530 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1531 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1532 nullptr, deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1534 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1539 if (EI_DYNAMICS(inputrec->eI))
1541 /* Turn on signal handling on all nodes */
1543 * (A user signal from the PME nodes (if any)
1544 * is communicated to the PP nodes.
1546 signal_handler_install();
1549 pull_t* pull_work = nullptr;
1550 if (thisRankHasDuty(cr, DUTY_PP))
1552 /* Assumes uniform use of the number of OpenMP threads */
1553 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1555 if (inputrec->bPull)
1557 /* Initialize pull code */
1558 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1559 inputrec->fepvals->init_lambda);
1560 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1562 initPullHistory(pull_work, &observablesHistory);
1564 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1566 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1570 std::unique_ptr<EnforcedRotation> enforcedRotation;
1571 if (inputrec->bRot)
1573 /* Initialize enforced rotation code */
1574 enforcedRotation =
1575 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1576 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1579 t_swap* swap = nullptr;
1580 if (inputrec->eSwapCoords != eswapNO)
1582 /* Initialize ion swapping code */
1583 swap = init_swapcoords(fplog, inputrec,
1584 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1585 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1586 oenv, mdrunOptions, startingBehavior);
1589 /* Let makeConstraints know whether we have essential dynamics constraints. */
1590 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1591 ms, &nrnb, wcycle, fr->bMolPBC);
1593 /* Energy terms and groups */
1594 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1595 inputrec->fepvals->n_lambda);
1597 /* Kinetic energy data */
1598 gmx_ekindata_t ekind;
1599 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1601 /* Set up interactive MD (IMD) */
1602 auto imdSession =
1603 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1604 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1605 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1607 if (DOMAINDECOMP(cr))
1609 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1610 /* This call is not included in init_domain_decomposition mainly
1611 * because fr->cginfo_mb is set later.
1613 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1614 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1617 // TODO This is not the right place to manage the lifetime of
1618 // this data structure, but currently it's the easiest way to
1619 // make it work.
1620 MdrunScheduleWorkload runScheduleWork;
1621 // Also populates the simulation constant workload description.
1622 runScheduleWork.simulationWork =
1623 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1624 useGpuForUpdate, devFlags.enableGpuBufferOps,
1625 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1627 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1628 if (gpusWereDetected
1629 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1630 || runScheduleWork.simulationWork.useGpuBufferOps))
1632 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1633 ? GpuApiCallBehavior::Async
1634 : GpuApiCallBehavior::Sync;
1635 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1636 "GPU device stream manager should be initialized to use GPU.");
1637 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1638 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1639 fr->stateGpu = stateGpu.get();
1642 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1643 SimulatorBuilder simulatorBuilder;
1645 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1646 simulatorBuilder.add(std::move(membedHolder));
1647 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1648 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1651 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1652 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1653 simulatorBuilder.add(ConstraintsParam(
1654 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1655 vsite.get()));
1656 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1657 simulatorBuilder.add(
1658 LegacyInput(static_cast<int>(filenames.size()), filenames.data(), inputrec, fr));
1659 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1660 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1661 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1662 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1663 // Todo move to an MDModule
1664 simulatorBuilder.add(IonSwapping(swap));
1665 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1666 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1668 // build and run simulator object based on user-input
1669 auto simulator = simulatorBuilder.build(useModularSimulator);
1670 simulator->run();
1672 if (fr->pmePpCommGpu)
1674 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1675 fr->pmePpCommGpu.reset();
1678 if (inputrec->bPull)
1680 finish_pull(pull_work);
1682 finish_swapcoords(swap);
1684 else
1686 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1687 /* do PME only */
1688 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1689 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode,
1690 deviceStreamManager.get());
1693 wallcycle_stop(wcycle, ewcRUN);
1695 /* Finish up, write some stuff
1696 * if rerunMD, don't write last frame again
1698 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1699 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1701 // clean up cycle counter
1702 wallcycle_destroy(wcycle);
1704 deviceStreamManager.reset(nullptr);
1705 // Free PME data
1706 if (pmedata)
1708 gmx_pme_destroy(pmedata);
1709 pmedata = nullptr;
1712 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1713 // before we destroy the GPU context(s) in free_gpu().
1714 // Pinned buffers are associated with contexts in CUDA.
1715 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1716 mdAtoms.reset(nullptr);
1717 globalState.reset(nullptr);
1718 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1719 gpuBonded.reset(nullptr);
1720 /* Free pinned buffers in *fr */
1721 delete fr;
1722 fr = nullptr;
1723 // TODO convert to C++ so we can get rid of these frees
1724 sfree(disresdata);
1725 sfree(oriresdata);
1727 if (hwinfo->gpu_info.n_dev > 0)
1729 /* stop the GPU profiler (only CUDA) */
1730 stopGpuProfiler();
1733 /* With tMPI we need to wait for all ranks to finish deallocation before
1734 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1735 * GPU and context.
1737 * This is not a concern in OpenCL where we use one context per rank which
1738 * is freed in nbnxn_gpu_free().
1740 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1741 * but it is easier and more futureproof to call it on the whole node.
1743 * Note that this function needs to be called even if GPUs are not used
1744 * in this run because the PME ranks have no knowledge of whether GPUs
1745 * are used or not, but all ranks need to enter the barrier below.
1746 * \todo Remove this physical node barrier after making sure
1747 * that it's not needed anymore (with a shared GPU run).
1749 if (GMX_THREAD_MPI)
1751 physicalNodeComm.barrier();
1754 free_gpu(deviceInfo);
1756 /* Does what it says */
1757 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1758 walltime_accounting_destroy(walltime_accounting);
1760 // Ensure log file content is written
1761 if (logFileHandle)
1763 gmx_fio_flush(logFileHandle);
1766 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1767 * exceptions were enabled before function was called. */
1768 if (bEnableFPE)
1770 gmx_fedisableexcept();
1773 auto rc = static_cast<int>(gmx_get_stop_condition());
1775 #if GMX_THREAD_MPI
1776 /* we need to join all threads. The sub-threads join when they
1777 exit this function, but the master thread needs to be told to
1778 wait for that. */
1779 if (PAR(cr) && MASTER(cr))
1781 tMPI_Finalize();
1783 #endif
1784 return rc;
1785 } // namespace gmx
1787 Mdrunner::~Mdrunner()
1789 // Clean up of the Manager.
1790 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1791 // but okay as long as threads synchronize some time before adding or accessing
1792 // a new set of restraints.
1793 if (restraintManager_)
1795 restraintManager_->clear();
1796 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1797 "restraints added during runner life time should be cleared at runner "
1798 "destruction.");
1802 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1804 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1805 // Not sure if this should be logged through the md logger or something else,
1806 // but it is helpful to have some sort of INFO level message sent somewhere.
1807 // std::cout << "Registering restraint named " << name << std::endl;
1809 // When multiple restraints are used, it may be wasteful to register them separately.
1810 // Maybe instead register an entire Restraint Manager as a force provider.
1811 restraintManager_->addToSpec(std::move(puller), name);
1814 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1816 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1818 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1820 class Mdrunner::BuilderImplementation
1822 public:
1823 BuilderImplementation() = delete;
1824 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1825 ~BuilderImplementation();
1827 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1828 real forceWarningThreshold,
1829 StartingBehavior startingBehavior);
1831 void addDomdec(const DomdecOptions& options);
1833 void addVerletList(int nstlist);
1835 void addReplicaExchange(const ReplicaExchangeParameters& params);
1837 void addNonBonded(const char* nbpu_opt);
1839 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1841 void addBondedTaskAssignment(const char* bonded_opt);
1843 void addUpdateTaskAssignment(const char* update_opt);
1845 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1847 void addFilenames(ArrayRef<const t_filenm> filenames);
1849 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1851 void addLogFile(t_fileio* logFileHandle);
1853 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1855 Mdrunner build();
1857 private:
1858 // Default parameters copied from runner.h
1859 // \todo Clarify source(s) of default parameters.
1861 const char* nbpu_opt_ = nullptr;
1862 const char* pme_opt_ = nullptr;
1863 const char* pme_fft_opt_ = nullptr;
1864 const char* bonded_opt_ = nullptr;
1865 const char* update_opt_ = nullptr;
1867 MdrunOptions mdrunOptions_;
1869 DomdecOptions domdecOptions_;
1871 ReplicaExchangeParameters replicaExchangeParameters_;
1873 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1874 int nstlist_ = 0;
1876 //! Multisim communicator handle.
1877 gmx_multisim_t* multiSimulation_;
1879 //! mdrun communicator
1880 MPI_Comm communicator_ = MPI_COMM_NULL;
1882 //! Print a warning if any force is larger than this (in kJ/mol nm).
1883 real forceWarningThreshold_ = -1;
1885 //! Whether the simulation will start afresh, or restart with/without appending.
1886 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1888 //! The modules that comprise the functionality of mdrun.
1889 std::unique_ptr<MDModules> mdModules_;
1891 //! \brief Parallelism information.
1892 gmx_hw_opt_t hardwareOptions_;
1894 //! filename options for simulation.
1895 ArrayRef<const t_filenm> filenames_;
1897 /*! \brief Handle to output environment.
1899 * \todo gmx_output_env_t needs lifetime management.
1901 gmx_output_env_t* outputEnvironment_ = nullptr;
1903 /*! \brief Non-owning handle to MD log file.
1905 * \todo Context should own output facilities for client.
1906 * \todo Improve log file handle management.
1907 * \internal
1908 * Code managing the FILE* relies on the ability to set it to
1909 * nullptr to check whether the filehandle is valid.
1911 t_fileio* logFileHandle_ = nullptr;
1914 * \brief Builder for simulation stop signal handler.
1916 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1919 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1920 compat::not_null<SimulationContext*> context) :
1921 mdModules_(std::move(mdModules))
1923 communicator_ = context->communicator_;
1924 multiSimulation_ = context->multiSimulation_.get();
1927 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1929 Mdrunner::BuilderImplementation&
1930 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1931 const real forceWarningThreshold,
1932 const StartingBehavior startingBehavior)
1934 mdrunOptions_ = options;
1935 forceWarningThreshold_ = forceWarningThreshold;
1936 startingBehavior_ = startingBehavior;
1937 return *this;
1940 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1942 domdecOptions_ = options;
1945 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1947 nstlist_ = nstlist;
1950 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1952 replicaExchangeParameters_ = params;
1955 Mdrunner Mdrunner::BuilderImplementation::build()
1957 auto newRunner = Mdrunner(std::move(mdModules_));
1959 newRunner.mdrunOptions = mdrunOptions_;
1960 newRunner.pforce = forceWarningThreshold_;
1961 newRunner.startingBehavior = startingBehavior_;
1962 newRunner.domdecOptions = domdecOptions_;
1964 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1965 newRunner.hw_opt = hardwareOptions_;
1967 // No invariant to check. This parameter exists to optionally override other behavior.
1968 newRunner.nstlist_cmdline = nstlist_;
1970 newRunner.replExParams = replicaExchangeParameters_;
1972 newRunner.filenames = filenames_;
1974 newRunner.communicator = communicator_;
1976 // nullptr is a valid value for the multisim handle
1977 newRunner.ms = multiSimulation_;
1979 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1980 // \todo Update sanity checking when output environment has clearly specified invariants.
1981 // Initialization and default values for oenv are not well specified in the current version.
1982 if (outputEnvironment_)
1984 newRunner.oenv = outputEnvironment_;
1986 else
1988 GMX_THROW(gmx::APIError(
1989 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1992 newRunner.logFileHandle = logFileHandle_;
1994 if (nbpu_opt_)
1996 newRunner.nbpu_opt = nbpu_opt_;
1998 else
2000 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2003 if (pme_opt_ && pme_fft_opt_)
2005 newRunner.pme_opt = pme_opt_;
2006 newRunner.pme_fft_opt = pme_fft_opt_;
2008 else
2010 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2013 if (bonded_opt_)
2015 newRunner.bonded_opt = bonded_opt_;
2017 else
2019 GMX_THROW(gmx::APIError(
2020 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2023 if (update_opt_)
2025 newRunner.update_opt = update_opt_;
2027 else
2029 GMX_THROW(gmx::APIError(
2030 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2034 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2036 if (stopHandlerBuilder_)
2038 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2040 else
2042 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2045 return newRunner;
2048 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2050 nbpu_opt_ = nbpu_opt;
2053 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2055 pme_opt_ = pme_opt;
2056 pme_fft_opt_ = pme_fft_opt;
2059 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2061 bonded_opt_ = bonded_opt;
2064 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2066 update_opt_ = update_opt;
2069 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2071 hardwareOptions_ = hardwareOptions;
2074 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2076 filenames_ = filenames;
2079 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2081 outputEnvironment_ = outputEnvironment;
2084 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2086 logFileHandle_ = logFileHandle;
2089 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2091 stopHandlerBuilder_ = std::move(builder);
2094 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2095 compat::not_null<SimulationContext*> context) :
2096 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2100 MdrunnerBuilder::~MdrunnerBuilder() = default;
2102 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2103 real forceWarningThreshold,
2104 const StartingBehavior startingBehavior)
2106 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2107 return *this;
2110 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2112 impl_->addDomdec(options);
2113 return *this;
2116 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2118 impl_->addVerletList(nstlist);
2119 return *this;
2122 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2124 impl_->addReplicaExchange(params);
2125 return *this;
2128 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2130 impl_->addNonBonded(nbpu_opt);
2131 return *this;
2134 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2136 // The builder method may become more general in the future, but in this version,
2137 // parameters for PME electrostatics are both required and the only parameters
2138 // available.
2139 if (pme_opt && pme_fft_opt)
2141 impl_->addPME(pme_opt, pme_fft_opt);
2143 else
2145 GMX_THROW(
2146 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2148 return *this;
2151 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2153 impl_->addBondedTaskAssignment(bonded_opt);
2154 return *this;
2157 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2159 impl_->addUpdateTaskAssignment(update_opt);
2160 return *this;
2163 Mdrunner MdrunnerBuilder::build()
2165 return impl_->build();
2168 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2170 impl_->addHardwareOptions(hardwareOptions);
2171 return *this;
2174 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2176 impl_->addFilenames(filenames);
2177 return *this;
2180 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2182 impl_->addOutputEnvironment(outputEnvironment);
2183 return *this;
2186 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2188 impl_->addLogFile(logFileHandle);
2189 return *this;
2192 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2194 impl_->addStopHandlerBuilder(std::move(builder));
2195 return *this;
2198 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2200 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;
2202 } // namespace gmx