Fix clang-format
[gromacs.git] / src / gromacs / mdrun / runner.cpp
blobcc76ce80d7e2adb580ab19b2f3dbbce0e8019d25
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/hardware/cpuinfo.h"
79 #include "gromacs/hardware/detecthardware.h"
80 #include "gromacs/hardware/device_management.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/checkpointdata.h"
114 #include "gromacs/mdtypes/commrec.h"
115 #include "gromacs/mdtypes/enerdata.h"
116 #include "gromacs/mdtypes/fcdata.h"
117 #include "gromacs/mdtypes/forcerec.h"
118 #include "gromacs/mdtypes/group.h"
119 #include "gromacs/mdtypes/inputrec.h"
120 #include "gromacs/mdtypes/interaction_const.h"
121 #include "gromacs/mdtypes/md_enums.h"
122 #include "gromacs/mdtypes/mdatom.h"
123 #include "gromacs/mdtypes/mdrunoptions.h"
124 #include "gromacs/mdtypes/observableshistory.h"
125 #include "gromacs/mdtypes/simulation_workload.h"
126 #include "gromacs/mdtypes/state.h"
127 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
128 #include "gromacs/modularsimulator/modularsimulator.h"
129 #include "gromacs/nbnxm/gpu_data_mgmt.h"
130 #include "gromacs/nbnxm/nbnxm.h"
131 #include "gromacs/nbnxm/pairlist_tuning.h"
132 #include "gromacs/pbcutil/pbc.h"
133 #include "gromacs/pulling/output.h"
134 #include "gromacs/pulling/pull.h"
135 #include "gromacs/pulling/pull_rotation.h"
136 #include "gromacs/restraint/manager.h"
137 #include "gromacs/restraint/restraintmdmodule.h"
138 #include "gromacs/restraint/restraintpotential.h"
139 #include "gromacs/swap/swapcoords.h"
140 #include "gromacs/taskassignment/decidegpuusage.h"
141 #include "gromacs/taskassignment/decidesimulationworkload.h"
142 #include "gromacs/taskassignment/resourcedivision.h"
143 #include "gromacs/taskassignment/taskassignment.h"
144 #include "gromacs/taskassignment/usergpuids.h"
145 #include "gromacs/timing/gpu_timing.h"
146 #include "gromacs/timing/wallcycle.h"
147 #include "gromacs/timing/wallcyclereporting.h"
148 #include "gromacs/topology/mtop_util.h"
149 #include "gromacs/trajectory/trajectoryframe.h"
150 #include "gromacs/utility/basenetwork.h"
151 #include "gromacs/utility/cstringutil.h"
152 #include "gromacs/utility/exceptions.h"
153 #include "gromacs/utility/fatalerror.h"
154 #include "gromacs/utility/filestream.h"
155 #include "gromacs/utility/gmxassert.h"
156 #include "gromacs/utility/gmxmpi.h"
157 #include "gromacs/utility/keyvaluetree.h"
158 #include "gromacs/utility/logger.h"
159 #include "gromacs/utility/loggerbuilder.h"
160 #include "gromacs/utility/mdmodulenotification.h"
161 #include "gromacs/utility/physicalnodecommunicator.h"
162 #include "gromacs/utility/pleasecite.h"
163 #include "gromacs/utility/programcontext.h"
164 #include "gromacs/utility/smalloc.h"
165 #include "gromacs/utility/stringutil.h"
167 #include "isimulator.h"
168 #include "membedholder.h"
169 #include "replicaexchange.h"
170 #include "simulationinput.h"
171 #include "simulationinpututility.h"
172 #include "simulatorbuilder.h"
174 #if GMX_FAHCORE
175 # include "corewrap.h"
176 #endif
178 namespace gmx
182 /*! \brief Manage any development feature flag variables encountered
184 * The use of dev features indicated by environment variables is
185 * logged in order to ensure that runs with such features enabled can
186 * be identified from their log and standard output. Any cross
187 * dependencies are also checked, and if unsatisfied, a fatal error
188 * issued.
190 * Note that some development features overrides are applied already here:
191 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
193 * \param[in] mdlog Logger object.
194 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
195 * \param[in] pmeRunMode The PME run mode for this run
196 * \returns The object populated with development feature flags.
198 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
199 const bool useGpuForNonbonded,
200 const PmeRunMode pmeRunMode)
202 DevelopmentFeatureFlags devFlags;
204 // Some builds of GCC 5 give false positive warnings that these
205 // getenv results are ignored when clearly they are used.
206 #pragma GCC diagnostic push
207 #pragma GCC diagnostic ignored "-Wunused-result"
209 devFlags.enableGpuBufferOps =
210 GMX_GPU_CUDA && useGpuForNonbonded && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
211 devFlags.enableGpuHaloExchange = GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_DD_COMMS") != nullptr;
212 devFlags.enableGpuPmePPComm =
213 GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
215 #pragma GCC diagnostic pop
217 if (devFlags.enableGpuBufferOps)
219 GMX_LOG(mdlog.warning)
220 .asParagraph()
221 .appendTextFormatted(
222 "This run uses the 'GPU buffer ops' feature, enabled by the "
223 "GMX_USE_GPU_BUFFER_OPS environment variable.");
226 if (devFlags.forceGpuUpdateDefault)
228 GMX_LOG(mdlog.warning)
229 .asParagraph()
230 .appendTextFormatted(
231 "This run will default to '-update gpu' as requested by the "
232 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
233 "decomposition lacks substantial testing and should be used with caution.");
236 if (devFlags.enableGpuHaloExchange)
238 if (useGpuForNonbonded)
240 if (!devFlags.enableGpuBufferOps)
242 GMX_LOG(mdlog.warning)
243 .asParagraph()
244 .appendTextFormatted(
245 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
246 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
247 devFlags.enableGpuBufferOps = true;
249 GMX_LOG(mdlog.warning)
250 .asParagraph()
251 .appendTextFormatted(
252 "This run has requested the 'GPU halo exchange' feature, enabled by "
253 "the "
254 "GMX_GPU_DD_COMMS environment variable.");
256 else
258 GMX_LOG(mdlog.warning)
259 .asParagraph()
260 .appendTextFormatted(
261 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
262 "halo exchange' feature will not be enabled as nonbonded interactions "
263 "are not offloaded.");
264 devFlags.enableGpuHaloExchange = false;
268 if (devFlags.enableGpuPmePPComm)
270 if (pmeRunMode == PmeRunMode::GPU)
272 if (!devFlags.enableGpuBufferOps)
274 GMX_LOG(mdlog.warning)
275 .asParagraph()
276 .appendTextFormatted(
277 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
278 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
279 devFlags.enableGpuBufferOps = true;
281 GMX_LOG(mdlog.warning)
282 .asParagraph()
283 .appendTextFormatted(
284 "This run uses the 'GPU PME-PP communications' feature, enabled "
285 "by the GMX_GPU_PME_PP_COMMS environment variable.");
287 else
289 std::string clarification;
290 if (pmeRunMode == PmeRunMode::Mixed)
292 clarification =
293 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
294 "mode).";
296 else
298 clarification = "PME is not offloaded to the GPU.";
300 GMX_LOG(mdlog.warning)
301 .asParagraph()
302 .appendText(
303 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
304 "'GPU PME-PP communications' feature was not enabled as "
305 + clarification);
306 devFlags.enableGpuPmePPComm = false;
310 return devFlags;
313 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
315 * Used to ensure that the master thread does not modify mdrunner during copy
316 * on the spawned threads. */
317 static void threadMpiMdrunnerAccessBarrier()
319 #if GMX_THREAD_MPI
320 MPI_Barrier(MPI_COMM_WORLD);
321 #endif
324 Mdrunner Mdrunner::cloneOnSpawnedThread() const
326 auto newRunner = Mdrunner(std::make_unique<MDModules>());
328 // All runners in the same process share a restraint manager resource because it is
329 // part of the interface to the client code, which is associated only with the
330 // original thread. Handles to the same resources can be obtained by copy.
332 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
335 // Copy members of master runner.
336 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
337 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
338 newRunner.hw_opt = hw_opt;
339 newRunner.filenames = filenames;
341 newRunner.oenv = oenv;
342 newRunner.mdrunOptions = mdrunOptions;
343 newRunner.domdecOptions = domdecOptions;
344 newRunner.nbpu_opt = nbpu_opt;
345 newRunner.pme_opt = pme_opt;
346 newRunner.pme_fft_opt = pme_fft_opt;
347 newRunner.bonded_opt = bonded_opt;
348 newRunner.update_opt = update_opt;
349 newRunner.nstlist_cmdline = nstlist_cmdline;
350 newRunner.replExParams = replExParams;
351 newRunner.pforce = pforce;
352 // Give the spawned thread the newly created valid communicator
353 // for the simulation.
354 newRunner.communicator = MPI_COMM_WORLD;
355 newRunner.ms = ms;
356 newRunner.startingBehavior = startingBehavior;
357 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
359 threadMpiMdrunnerAccessBarrier();
361 return newRunner;
364 /*! \brief The callback used for running on spawned threads.
366 * Obtains the pointer to the master mdrunner object from the one
367 * argument permitted to the thread-launch API call, copies it to make
368 * a new runner for this thread, reinitializes necessary data, and
369 * proceeds to the simulation. */
370 static void mdrunner_start_fn(const void* arg)
374 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
375 /* copy the arg list to make sure that it's thread-local. This
376 doesn't copy pointed-to items, of course; fnm, cr and fplog
377 are reset in the call below, all others should be const. */
378 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
379 mdrunner.mdrunner();
381 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
385 void Mdrunner::spawnThreads(int numThreadsToLaunch)
387 #if GMX_THREAD_MPI
388 /* now spawn new threads that start mdrunner_start_fn(), while
389 the main thread returns. Thread affinity is handled later. */
390 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
391 static_cast<const void*>(this))
392 != TMPI_SUCCESS)
394 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
397 // Give the master thread the newly created valid communicator for
398 // the simulation.
399 communicator = MPI_COMM_WORLD;
400 threadMpiMdrunnerAccessBarrier();
401 #else
402 GMX_UNUSED_VALUE(numThreadsToLaunch);
403 GMX_UNUSED_VALUE(mdrunner_start_fn);
404 #endif
407 } // namespace gmx
409 /*! \brief Initialize variables for Verlet scheme simulation */
410 static void prepare_verlet_scheme(FILE* fplog,
411 t_commrec* cr,
412 t_inputrec* ir,
413 int nstlist_cmdline,
414 const gmx_mtop_t* mtop,
415 const matrix box,
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");
426 else
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)
451 fprintf(fplog,
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)
478 assert(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)
489 sprintf(sbuf_msg,
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));
493 else
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 */
508 namespace gmx
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;
519 std::string warning;
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.
526 gpuIsUseful = false;
527 warning =
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.";
533 if (EI_TPI(ir.eI))
535 gpuIsUseful = false;
536 warning = "TPI is not implemented for GPUs.";
539 if (!gpuIsUseful && issueWarning)
541 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
544 return gpuIsUseful;
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;
579 else
581 GMX_ASSERT(false, "Option string should have been checked for sanity already");
584 return returnValue;
587 //! Finish run, aggregate data to print performance info.
588 static void finish_run(FILE* fplog,
589 const gmx::MDLogger& mdlog,
590 const t_commrec* cr,
591 const t_inputrec* inputrec,
592 t_nrnb nrnb[],
593 gmx_wallcycle_t wcycle,
594 gmx_walltime_accounting_t walltime_accounting,
595 nonbonded_verlet_t* nbv,
596 const gmx_pme_t* pme,
597 gmx_bool bWriteStat)
599 double delta_t = 0;
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)
621 .asParagraph()
622 .appendText("Simulation ended prematurely, no performance report will be written.");
623 printReport = false;
626 t_nrnb* nrnb_tot;
627 std::unique_ptr<t_nrnb> nrnbTotalStorage;
628 if (cr->nnodes > 1)
630 nrnbTotalStorage = std::make_unique<t_nrnb>();
631 nrnb_tot = nrnbTotalStorage.get();
632 #if GMX_MPI
633 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
634 #endif
636 else
638 nrnb_tot = nrnb;
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);
644 if (cr->nnodes > 1)
646 #if GMX_MPI
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,
649 cr->mpi_comm_mysim);
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);
655 #endif
657 else
659 elapsed_time_over_all_ranks = elapsed_time;
660 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
663 if (printReport)
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));
683 if (printReport)
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,
695 &pme_gpu_timings);
697 if (EI_DYNAMICS(inputrec->eI))
699 delta_t = inputrec->delta_t;
702 if (fplog)
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);
708 if (bWriteStat)
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()
719 matrix box;
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 */
731 gmx_mtop_t mtop;
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, communicator);
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(communicator, 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);
781 // Print citation requests after all software/hardware printing
782 pleaseCiteGromacs(fplog);
784 // TODO: Use SimulationInputHolder member to access SimulationInput. Issue #3374.
785 const auto* const tprFilename = ftp2fn(efTPR, filenames.size(), filenames.data());
786 const auto* const cpiFilename = opt2fn("-cpi", filenames.size(), filenames.data());
787 // Note that, as of this change, there is no public API for simulationInput or its creation.
788 // TODO: (#3374) Public API for providing simulationInput from client.
789 auto simulationInput = detail::makeSimulationInput(tprFilename, cpiFilename);
791 // Note: legacy program logic relies on checking whether these pointers are assigned.
792 // Objects may or may not be allocated later.
793 std::unique_ptr<t_inputrec> inputrec;
794 std::unique_ptr<t_state> globalState;
796 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
798 if (isSimulationMasterRank)
800 // Allocate objects to be initialized by later function calls.
801 /* Only the master rank has the global state */
802 globalState = std::make_unique<t_state>();
803 inputrec = std::make_unique<t_inputrec>();
805 /* Read (nearly) all data required for the simulation
806 * and keep the partly serialized tpr contents to send to other ranks later
808 applyGlobalSimulationState(*simulationInput.object_, partialDeserializedTpr.get(),
809 globalState.get(), inputrec.get(), &mtop);
812 /* Check and update the hardware options for internal consistency */
813 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
814 inputrec.get());
816 if (GMX_THREAD_MPI && isSimulationMasterRank)
818 bool useGpuForNonbonded = false;
819 bool useGpuForPme = false;
822 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
824 // If the user specified the number of ranks, then we must
825 // respect that, but in default mode, we need to allow for
826 // the number of GPUs to choose the number of ranks.
827 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
828 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
829 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
830 canUseGpuForNonbonded,
831 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
832 hw_opt.nthreads_tmpi);
833 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
834 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
835 *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
837 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
839 /* Determine how many thread-MPI ranks to start.
841 * TODO Over-writing the user-supplied value here does
842 * prevent any possible subsequent checks from working
843 * correctly. */
844 hw_opt.nthreads_tmpi =
845 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
846 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
848 // Now start the threads for thread MPI.
849 spawnThreads(hw_opt.nthreads_tmpi);
850 // The spawned threads enter mdrunner() and execution of
851 // master and spawned threads joins at the end of this block.
852 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
855 GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
856 "Must have valid world communicator unless running a multi-simulation");
857 CommrecHandle crHandle = init_commrec(communicator);
858 t_commrec* cr = crHandle.get();
859 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
861 if (PAR(cr))
863 /* now broadcast everything to the non-master nodes/threads: */
864 if (!isSimulationMasterRank)
866 // Until now, only the master rank has a non-null pointer.
867 // On non-master ranks, allocate the object that will receive data in the following call.
868 inputrec = std::make_unique<t_inputrec>();
870 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec.get(), &mtop,
871 partialDeserializedTpr.get());
873 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
874 partialDeserializedTpr.reset(nullptr);
876 // Now the number of ranks is known to all ranks, and each knows
877 // the inputrec read by the master rank. The ranks can now all run
878 // the task-deciding functions and will agree on the result
879 // without needing to communicate.
880 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
882 // Note that these variables describe only their own node.
884 // Note that when bonded interactions run on a GPU they always run
885 // alongside a nonbonded task, so do not influence task assignment
886 // even though they affect the force calculation workload.
887 bool useGpuForNonbonded = false;
888 bool useGpuForPme = false;
889 bool useGpuForBonded = false;
890 bool useGpuForUpdate = false;
891 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
894 // It's possible that there are different numbers of GPUs on
895 // different nodes, which is the user's responsibility to
896 // handle. If unsuitable, we will notice that during task
897 // assignment.
898 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
899 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
900 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
901 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
902 useGpuForPme = decideWhetherToUseGpusForPme(
903 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec,
904 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
905 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
906 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
907 useGpuForBonded = decideWhetherToUseGpusForBonded(
908 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
909 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
910 domdecOptions.numPmeRanks, gpusWereDetected);
912 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
914 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
916 // Initialize development feature flags that enabled by environment variable
917 // and report those features that are enabled.
918 const DevelopmentFeatureFlags devFlags =
919 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
921 const bool useModularSimulator =
922 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
923 nullptr, doEssentialDynamics, membedHolder.doMembed());
925 // Build restraints.
926 // TODO: hide restraint implementation details from Mdrunner.
927 // There is nothing unique about restraints at this point as far as the
928 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
929 // factory functions from the SimulationContext on which to call mdModules_->add().
930 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
931 for (auto&& restraint : restraintManager_->getRestraints())
933 auto module = RestraintMDModule::create(restraint, restraint->sites());
934 mdModules_->add(std::move(module));
937 // TODO: Error handling
938 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
939 // now that the MdModules know their options, they know which callbacks to sign up to
940 mdModules_->subscribeToSimulationSetupNotifications();
941 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
943 if (inputrec->internalParameters != nullptr)
945 mdModulesNotifier.notify(*inputrec->internalParameters);
948 if (fplog != nullptr)
950 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
951 fprintf(fplog, "\n");
954 if (SIMMASTER(cr))
956 /* In rerun, set velocities to zero if present */
957 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
959 // rerun does not use velocities
960 GMX_LOG(mdlog.info)
961 .asParagraph()
962 .appendText(
963 "Rerun trajectory contains velocities. Rerun does only evaluate "
964 "potential energy and forces. The velocities will be ignored.");
965 for (int i = 0; i < globalState->natoms; i++)
967 clear_rvec(globalState->v[i]);
969 globalState->flags &= ~(1 << estV);
972 /* now make sure the state is initialized and propagated */
973 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
976 /* NM and TPI parallelize over force/energy calculations, not atoms,
977 * so we need to initialize and broadcast the global state.
979 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
981 if (!MASTER(cr))
983 globalState = std::make_unique<t_state>();
985 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
986 globalState.get());
989 /* A parallel command line option consistency check that we can
990 only do after any threads have started. */
991 if (!PAR(cr)
992 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
993 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
995 gmx_fatal(FARGS,
996 "The -dd or -npme option request a parallel simulation, "
997 #if !GMX_MPI
998 "but %s was compiled without threads or MPI enabled",
999 output_env_get_program_display_name(oenv));
1000 #elif GMX_THREAD_MPI
1001 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
1002 #else
1003 "but %s was not started through mpirun/mpiexec or only one rank was requested "
1004 "through mpirun/mpiexec",
1005 output_env_get_program_display_name(oenv));
1006 #endif
1009 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1011 gmx_fatal(FARGS,
1012 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1013 "these are not compatible with mdrun -rerun");
1016 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1018 if (domdecOptions.numPmeRanks > 0)
1020 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1021 "PME-only ranks are requested, but the system does not use PME "
1022 "for electrostatics or LJ");
1025 domdecOptions.numPmeRanks = 0;
1028 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1030 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1031 * improve performance with many threads per GPU, since our OpenMP
1032 * scaling is bad, but it's difficult to automate the setup.
1034 domdecOptions.numPmeRanks = 0;
1036 if (useGpuForPme)
1038 if (domdecOptions.numPmeRanks < 0)
1040 domdecOptions.numPmeRanks = 0;
1041 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1043 else
1045 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1046 "PME GPU decomposition is not supported");
1050 #if GMX_FAHCORE
1051 if (MASTER(cr))
1053 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1055 #endif
1057 /* NMR restraints must be initialized before load_checkpoint,
1058 * since with time averaging the history is added to t_state.
1059 * For proper consistency check we therefore need to extend
1060 * t_state here.
1061 * So the PME-only nodes (if present) will also initialize
1062 * the distance restraints.
1065 /* This needs to be called before read_checkpoint to extend the state */
1066 t_disresdata* disresdata;
1067 snew(disresdata, 1);
1068 init_disres(fplog, &mtop, inputrec.get(), DisResRunMode::MDRun,
1069 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1070 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1071 globalState.get(), replExParams.exchangeInterval > 0);
1073 t_oriresdata* oriresdata;
1074 snew(oriresdata, 1);
1075 init_orires(fplog, &mtop, inputrec.get(), cr, ms, globalState.get(), oriresdata);
1077 auto deform = prepareBoxDeformation(
1078 globalState != nullptr ? globalState->box : box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1079 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mygroup, *inputrec);
1081 ObservablesHistory observablesHistory = {};
1083 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1084 if (startingBehavior != StartingBehavior::NewSimulation)
1086 /* Check if checkpoint file exists before doing continuation.
1087 * This way we can use identical input options for the first and subsequent runs...
1089 if (mdrunOptions.numStepsCommandline > -2)
1091 /* Temporarily set the number of steps to unlimited to avoid
1092 * triggering the nsteps check in load_checkpoint().
1093 * This hack will go away soon when the -nsteps option is removed.
1095 inputrec->nsteps = -1;
1098 // Finish applying initial simulation state information from external sources on all ranks.
1099 // Reconcile checkpoint file data with Mdrunner state established up to this point.
1100 applyLocalState(*simulationInput.object_, logFileHandle, cr, domdecOptions.numCells,
1101 inputrec.get(), globalState.get(), &observablesHistory,
1102 mdrunOptions.reproducible, mdModules_->notifier(),
1103 modularSimulatorCheckpointData.get(), useModularSimulator);
1104 // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1105 // invariants
1106 // on all code paths.
1107 // Write checkpoint or provide hook to update SimulationInput.
1108 // If there was a checkpoint file, SimulationInput contains more information
1109 // than if there wasn't. At this point, we have synchronized the in-memory
1110 // state with the filesystem state only for restarted simulations. We should
1111 // be calling applyLocalState unconditionally and expect that the completeness
1112 // of SimulationInput is not dependent on its creation method.
1114 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1116 // Now we can start normal logging to the truncated log file.
1117 fplog = gmx_fio_getfp(logFileHandle);
1118 prepareLogAppending(fplog);
1119 logOwner = buildLogger(fplog, MASTER(cr));
1120 mdlog = logOwner.logger();
1124 if (mdrunOptions.numStepsCommandline > -2)
1126 GMX_LOG(mdlog.info)
1127 .asParagraph()
1128 .appendText(
1129 "The -nsteps functionality is deprecated, and may be removed in a future "
1130 "version. "
1131 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1132 "file field.");
1134 /* override nsteps with value set on the commandline */
1135 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1137 if (isSimulationMasterRank)
1139 copy_mat(globalState->box, box);
1142 if (PAR(cr))
1144 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1147 if (inputrec->cutoff_scheme != ecutsVERLET)
1149 gmx_fatal(FARGS,
1150 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1151 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1153 /* Update rlist and nstlist. */
1154 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1155 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1156 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1158 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1159 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1160 *hwinfo->cpuInfo);
1162 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1163 // This builder is necessary while we have multi-part construction
1164 // of DD. Before DD is constructed, we use the existence of
1165 // the builder object to indicate that further construction of DD
1166 // is needed.
1167 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1168 if (useDomainDecomposition)
1170 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1171 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1172 positionsFromStatePointer(globalState.get()));
1174 else
1176 /* PME, if used, is done on all nodes with 1D decomposition */
1177 cr->nnodes = cr->sizeOfDefaultCommunicator;
1178 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1179 cr->nodeid = cr->rankInDefaultCommunicator;
1180 cr->npmenodes = 0;
1181 cr->duty = (DUTY_PP | DUTY_PME);
1183 if (inputrec->pbcType == PbcType::Screw)
1185 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1189 // Produce the task assignment for this rank - done after DD is constructed
1190 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1191 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1192 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1193 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1194 // TODO cr->duty & DUTY_PME should imply that a PME
1195 // algorithm is active, but currently does not.
1196 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1198 // Get the device handles for the modules, nullptr when no task is assigned.
1199 int deviceId = -1;
1200 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1202 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1203 bool useTiming = true;
1205 if (GMX_GPU_CUDA)
1207 /* WARNING: CUDA timings are incorrect with multiple streams.
1208 * This is the main reason why they are disabled by default.
1210 // TODO: Consider turning on by default when we can detect nr of streams.
1211 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1213 else if (GMX_GPU_OPENCL)
1215 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1218 // TODO Currently this is always built, yet DD partition code
1219 // checks if it is built before using it. Probably it should
1220 // become an MDModule that is made only when another module
1221 // requires it (e.g. pull, CompEl, density fitting), so that we
1222 // don't update the local atom sets unilaterally every step.
1223 LocalAtomSetManager atomSets;
1224 if (ddBuilder)
1226 // TODO Pass the GPU streams to ddBuilder to use in buffer
1227 // transfers (e.g. halo exchange)
1228 cr->dd = ddBuilder->build(&atomSets);
1229 // The builder's job is done, so destruct it
1230 ddBuilder.reset(nullptr);
1231 // Note that local state still does not exist yet.
1234 // The GPU update is decided here because we need to know whether the constraints or
1235 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1236 // defined). This is only known after DD is initialized, hence decision on using GPU
1237 // update is done so late.
1240 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1242 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1243 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1244 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1245 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1246 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1248 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1250 const bool printHostName = (cr->nnodes > 1);
1251 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1253 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1255 if (deviceInfo != nullptr)
1257 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1259 dd_setup_dlb_resource_sharing(cr, deviceId);
1261 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1262 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1263 useGpuForUpdate, useTiming);
1266 // If the user chose a task assignment, give them some hints
1267 // where appropriate.
1268 if (!userGpuTaskAssignment.empty())
1270 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1273 if (PAR(cr))
1275 /* After possible communicator splitting in make_dd_communicators.
1276 * we can set up the intra/inter node communication.
1278 gmx_setup_nodecomm(fplog, cr);
1281 #if GMX_MPI
1282 if (isMultiSim(ms))
1284 GMX_LOG(mdlog.warning)
1285 .asParagraph()
1286 .appendTextFormatted(
1287 "This is simulation %d out of %d running as a composite GROMACS\n"
1288 "multi-simulation job. Setup for this simulation:\n",
1289 ms->simulationIndex_, ms->numSimulations_);
1291 GMX_LOG(mdlog.warning)
1292 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1293 # if GMX_THREAD_MPI
1294 cr->nnodes == 1 ? "thread" : "threads"
1295 # else
1296 cr->nnodes == 1 ? "process" : "processes"
1297 # endif
1299 fflush(stderr);
1300 #endif
1302 // If mdrun -pin auto honors any affinity setting that already
1303 // exists. If so, it is nice to provide feedback about whether
1304 // that existing affinity setting was from OpenMP or something
1305 // else, so we run this code both before and after we initialize
1306 // the OpenMP support.
1307 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1308 /* Check and update the number of OpenMP threads requested */
1309 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1310 pmeRunMode, mtop, *inputrec);
1312 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1313 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1315 // Enable FP exception detection, but not in
1316 // Release mode and not for compilers with known buggy FP
1317 // exception support (clang with any optimization) or suspected
1318 // buggy FP exception support (gcc 7.* with optimization).
1319 #if !defined NDEBUG \
1320 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1321 && defined __OPTIMIZE__)
1322 const bool bEnableFPE = true;
1323 #else
1324 const bool bEnableFPE = false;
1325 #endif
1326 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1327 if (bEnableFPE)
1329 gmx_feenableexcept();
1332 /* Now that we know the setup is consistent, check for efficiency */
1333 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1334 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1336 /* getting number of PP/PME threads on this MPI / tMPI rank.
1337 PME: env variable should be read only on one node to make sure it is
1338 identical everywhere;
1340 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1341 : gmx_omp_nthreads_get(emntPME);
1342 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1343 physicalNodeComm, mdlog);
1345 // Enable Peer access between GPUs where available
1346 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1347 // any of the GPU communication features are active.
1348 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1349 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1351 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1354 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1356 /* Before setting affinity, check whether the affinity has changed
1357 * - which indicates that probably the OpenMP library has changed it
1358 * since we first checked).
1360 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1362 int numThreadsOnThisNode, intraNodeThreadOffset;
1363 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1364 &intraNodeThreadOffset);
1366 /* Set the CPU affinity */
1367 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1368 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1371 if (mdrunOptions.timingOptions.resetStep > -1)
1373 GMX_LOG(mdlog.info)
1374 .asParagraph()
1375 .appendText(
1376 "The -resetstep functionality is deprecated, and may be removed in a "
1377 "future version.");
1379 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1381 if (PAR(cr))
1383 /* Master synchronizes its value of reset_counters with all nodes
1384 * including PME only nodes */
1385 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1386 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1387 wcycle_set_reset_counters(wcycle, reset_counters);
1390 // Membrane embedding must be initialized before we call init_forcerec()
1391 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1392 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1394 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1395 std::unique_ptr<MDAtoms> mdAtoms;
1396 std::unique_ptr<VirtualSitesHandler> vsite;
1397 std::unique_ptr<GpuBonded> gpuBonded;
1399 t_nrnb nrnb;
1400 if (thisRankHasDuty(cr, DUTY_PP))
1402 mdModulesNotifier.notify(*cr);
1403 mdModulesNotifier.notify(&atomSets);
1404 mdModulesNotifier.notify(inputrec->pbcType);
1405 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1406 /* Initiate forcerecord */
1407 fr = new t_forcerec;
1408 fr->forceProviders = mdModules_->initForceProviders();
1409 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1410 opt2fn("-table", filenames.size(), filenames.data()),
1411 opt2fn("-tablep", filenames.size(), filenames.data()),
1412 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1413 // Dirty hack, for fixing disres and orires should be made mdmodules
1414 fr->fcdata->disres = disresdata;
1415 fr->fcdata->orires = oriresdata;
1417 // Save a handle to device stream manager to use elsewhere in the code
1418 // TODO: Forcerec is not a correct place to store it.
1419 fr->deviceStreamManager = deviceStreamManager.get();
1421 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1423 GMX_RELEASE_ASSERT(
1424 deviceStreamManager != nullptr,
1425 "GPU device stream manager should be valid in order to use PME-PP direct "
1426 "communications.");
1427 GMX_RELEASE_ASSERT(
1428 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1429 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1430 "communications.");
1431 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1432 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1433 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1436 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo, useGpuForNonbonded,
1437 deviceStreamManager.get(), &mtop, box, wcycle);
1438 // TODO: Move the logic below to a GPU bonded builder
1439 if (useGpuForBonded)
1441 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1442 "GPU device stream manager should be valid in order to use GPU "
1443 "version of bonded forces.");
1444 gpuBonded = std::make_unique<GpuBonded>(
1445 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1446 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1447 fr->gpuBonded = gpuBonded.get();
1450 /* Initialize the mdAtoms structure.
1451 * mdAtoms is not filled with atom data,
1452 * as this can not be done now with domain decomposition.
1454 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1455 if (globalState && thisRankHasPmeGpuTask)
1457 // The pinning of coordinates in the global state object works, because we only use
1458 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1459 // points to the global state object without DD.
1460 // FIXME: MD and EM separately set up the local state - this should happen in the same
1461 // function, which should also perform the pinning.
1462 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1465 /* Initialize the virtual site communication */
1466 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1468 calc_shifts(box, fr->shift_vec);
1470 /* With periodic molecules the charge groups should be whole at start up
1471 * and the virtual sites should not be far from their proper positions.
1473 if (!inputrec->bContinuation && MASTER(cr)
1474 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1476 /* Make molecules whole at start of run */
1477 if (fr->pbcType != PbcType::No)
1479 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1481 if (vsite)
1483 /* Correct initial vsite positions are required
1484 * for the initial distribution in the domain decomposition
1485 * and for the initial shell prediction.
1487 constructVirtualSitesGlobal(mtop, globalState->x);
1491 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1493 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1494 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1497 else
1499 /* This is a PME only node */
1501 GMX_ASSERT(globalState == nullptr,
1502 "We don't need the state on a PME only rank and expect it to be unitialized");
1504 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1505 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1508 gmx_pme_t* sepPmeData = nullptr;
1509 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1510 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1511 "Double-checking that only PME-only ranks have no forcerec");
1512 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1514 // TODO should live in ewald module once its testing is improved
1516 // Later, this program could contain kernels that might be later
1517 // re-used as auto-tuning progresses, or subsequent simulations
1518 // are invoked.
1519 PmeGpuProgramStorage pmeGpuProgram;
1520 if (thisRankHasPmeGpuTask)
1522 GMX_RELEASE_ASSERT(
1523 (deviceStreamManager != nullptr),
1524 "GPU device stream manager should be initialized in order to use GPU for PME.");
1525 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1526 "GPU device should be initialized in order to use GPU for PME.");
1527 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1530 /* Initiate PME if necessary,
1531 * either on all nodes or on dedicated PME nodes only. */
1532 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1534 if (mdAtoms && mdAtoms->mdatoms())
1536 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1537 if (EVDW_PME(inputrec->vdwtype))
1539 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1542 if (cr->npmenodes > 0)
1544 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1545 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1546 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1549 if (thisRankHasDuty(cr, DUTY_PME))
1553 // TODO: This should be in the builder.
1554 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1555 "Device stream manager should be valid in order to use GPU "
1556 "version of PME.");
1557 GMX_RELEASE_ASSERT(
1558 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1559 "GPU PME stream should be valid in order to use GPU version of PME.");
1561 const DeviceContext* deviceContext =
1562 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1563 const DeviceStream* pmeStream =
1564 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1566 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1567 nChargePerturbed != 0, nTypePerturbed != 0,
1568 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1569 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1570 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1572 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1577 if (EI_DYNAMICS(inputrec->eI))
1579 /* Turn on signal handling on all nodes */
1581 * (A user signal from the PME nodes (if any)
1582 * is communicated to the PP nodes.
1584 signal_handler_install();
1587 pull_t* pull_work = nullptr;
1588 if (thisRankHasDuty(cr, DUTY_PP))
1590 /* Assumes uniform use of the number of OpenMP threads */
1591 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1593 if (inputrec->bPull)
1595 /* Initialize pull code */
1596 pull_work = init_pull(fplog, inputrec->pull, inputrec.get(), &mtop, cr, &atomSets,
1597 inputrec->fepvals->init_lambda);
1598 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1600 initPullHistory(pull_work, &observablesHistory);
1602 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1604 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1608 std::unique_ptr<EnforcedRotation> enforcedRotation;
1609 if (inputrec->bRot)
1611 /* Initialize enforced rotation code */
1612 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1613 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1614 startingBehavior);
1617 t_swap* swap = nullptr;
1618 if (inputrec->eSwapCoords != eswapNO)
1620 /* Initialize ion swapping code */
1621 swap = init_swapcoords(fplog, inputrec.get(),
1622 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1623 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1624 oenv, mdrunOptions, startingBehavior);
1627 /* Let makeConstraints know whether we have essential dynamics constraints. */
1628 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1629 ms, &nrnb, wcycle, fr->bMolPBC);
1631 /* Energy terms and groups */
1632 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1633 inputrec->fepvals->n_lambda);
1635 /* Kinetic energy data */
1636 gmx_ekindata_t ekind;
1637 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1639 /* Set up interactive MD (IMD) */
1640 auto imdSession =
1641 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1642 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1643 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1645 if (DOMAINDECOMP(cr))
1647 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1648 /* This call is not included in init_domain_decomposition mainly
1649 * because fr->cginfo_mb is set later.
1651 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1652 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1655 // TODO This is not the right place to manage the lifetime of
1656 // this data structure, but currently it's the easiest way to
1657 // make it work.
1658 MdrunScheduleWorkload runScheduleWork;
1659 // Also populates the simulation constant workload description.
1660 runScheduleWork.simulationWork =
1661 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1662 useGpuForUpdate, devFlags.enableGpuBufferOps,
1663 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1665 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1666 if (gpusWereDetected
1667 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1668 || runScheduleWork.simulationWork.useGpuBufferOps))
1670 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1671 ? GpuApiCallBehavior::Async
1672 : GpuApiCallBehavior::Sync;
1673 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1674 "GPU device stream manager should be initialized to use GPU.");
1675 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1676 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1677 fr->stateGpu = stateGpu.get();
1680 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1681 SimulatorBuilder simulatorBuilder;
1683 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1684 simulatorBuilder.add(std::move(membedHolder));
1685 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1686 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1689 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1690 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1691 simulatorBuilder.add(ConstraintsParam(
1692 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1693 vsite.get()));
1694 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1695 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1696 inputrec.get(), fr));
1697 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1698 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1699 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1700 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1701 // Todo move to an MDModule
1702 simulatorBuilder.add(IonSwapping(swap));
1703 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1704 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1705 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
1707 // build and run simulator object based on user-input
1708 auto simulator = simulatorBuilder.build(useModularSimulator);
1709 simulator->run();
1711 if (fr->pmePpCommGpu)
1713 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1714 fr->pmePpCommGpu.reset();
1717 if (inputrec->bPull)
1719 finish_pull(pull_work);
1721 finish_swapcoords(swap);
1723 else
1725 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1726 /* do PME only */
1727 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1728 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1729 deviceStreamManager.get());
1732 wallcycle_stop(wcycle, ewcRUN);
1734 /* Finish up, write some stuff
1735 * if rerunMD, don't write last frame again
1737 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1738 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1740 // clean up cycle counter
1741 wallcycle_destroy(wcycle);
1743 deviceStreamManager.reset(nullptr);
1744 // Free PME data
1745 if (pmedata)
1747 gmx_pme_destroy(pmedata);
1748 pmedata = nullptr;
1751 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1752 // before we destroy the GPU context(s)
1753 // Pinned buffers are associated with contexts in CUDA.
1754 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1755 mdAtoms.reset(nullptr);
1756 globalState.reset(nullptr);
1757 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1758 gpuBonded.reset(nullptr);
1759 /* Free pinned buffers in *fr */
1760 delete fr;
1761 fr = nullptr;
1762 // TODO convert to C++ so we can get rid of these frees
1763 sfree(disresdata);
1764 sfree(oriresdata);
1766 if (!hwinfo->deviceInfoList.empty())
1768 /* stop the GPU profiler (only CUDA) */
1769 stopGpuProfiler();
1772 /* With tMPI we need to wait for all ranks to finish deallocation before
1773 * destroying the CUDA context as some tMPI ranks may be sharing
1774 * GPU and context.
1776 * This is not a concern in OpenCL where we use one context per rank.
1778 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1779 * but it is easier and more futureproof to call it on the whole node.
1781 * Note that this function needs to be called even if GPUs are not used
1782 * in this run because the PME ranks have no knowledge of whether GPUs
1783 * are used or not, but all ranks need to enter the barrier below.
1784 * \todo Remove this physical node barrier after making sure
1785 * that it's not needed anymore (with a shared GPU run).
1787 if (GMX_THREAD_MPI)
1789 physicalNodeComm.barrier();
1791 releaseDevice(deviceInfo);
1793 /* Does what it says */
1794 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1795 walltime_accounting_destroy(walltime_accounting);
1797 // Ensure log file content is written
1798 if (logFileHandle)
1800 gmx_fio_flush(logFileHandle);
1803 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1804 * exceptions were enabled before function was called. */
1805 if (bEnableFPE)
1807 gmx_fedisableexcept();
1810 auto rc = static_cast<int>(gmx_get_stop_condition());
1812 #if GMX_THREAD_MPI
1813 /* we need to join all threads. The sub-threads join when they
1814 exit this function, but the master thread needs to be told to
1815 wait for that. */
1816 if (MASTER(cr))
1818 tMPI_Finalize();
1820 #endif
1821 return rc;
1822 } // namespace gmx
1824 Mdrunner::~Mdrunner()
1826 // Clean up of the Manager.
1827 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1828 // but okay as long as threads synchronize some time before adding or accessing
1829 // a new set of restraints.
1830 if (restraintManager_)
1832 restraintManager_->clear();
1833 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1834 "restraints added during runner life time should be cleared at runner "
1835 "destruction.");
1839 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1841 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1842 // Not sure if this should be logged through the md logger or something else,
1843 // but it is helpful to have some sort of INFO level message sent somewhere.
1844 // std::cout << "Registering restraint named " << name << std::endl;
1846 // When multiple restraints are used, it may be wasteful to register them separately.
1847 // Maybe instead register an entire Restraint Manager as a force provider.
1848 restraintManager_->addToSpec(std::move(puller), name);
1851 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1853 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1855 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1857 class Mdrunner::BuilderImplementation
1859 public:
1860 BuilderImplementation() = delete;
1861 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1862 ~BuilderImplementation();
1864 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1865 real forceWarningThreshold,
1866 StartingBehavior startingBehavior);
1868 void addDomdec(const DomdecOptions& options);
1870 void addVerletList(int nstlist);
1872 void addReplicaExchange(const ReplicaExchangeParameters& params);
1874 void addNonBonded(const char* nbpu_opt);
1876 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1878 void addBondedTaskAssignment(const char* bonded_opt);
1880 void addUpdateTaskAssignment(const char* update_opt);
1882 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1884 void addFilenames(ArrayRef<const t_filenm> filenames);
1886 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1888 void addLogFile(t_fileio* logFileHandle);
1890 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1892 Mdrunner build();
1894 private:
1895 // Default parameters copied from runner.h
1896 // \todo Clarify source(s) of default parameters.
1898 const char* nbpu_opt_ = nullptr;
1899 const char* pme_opt_ = nullptr;
1900 const char* pme_fft_opt_ = nullptr;
1901 const char* bonded_opt_ = nullptr;
1902 const char* update_opt_ = nullptr;
1904 MdrunOptions mdrunOptions_;
1906 DomdecOptions domdecOptions_;
1908 ReplicaExchangeParameters replicaExchangeParameters_;
1910 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1911 int nstlist_ = 0;
1913 //! Multisim communicator handle.
1914 gmx_multisim_t* multiSimulation_;
1916 //! mdrun communicator
1917 MPI_Comm communicator_ = MPI_COMM_NULL;
1919 //! Print a warning if any force is larger than this (in kJ/mol nm).
1920 real forceWarningThreshold_ = -1;
1922 //! Whether the simulation will start afresh, or restart with/without appending.
1923 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1925 //! The modules that comprise the functionality of mdrun.
1926 std::unique_ptr<MDModules> mdModules_;
1928 //! \brief Parallelism information.
1929 gmx_hw_opt_t hardwareOptions_;
1931 //! filename options for simulation.
1932 ArrayRef<const t_filenm> filenames_;
1934 /*! \brief Handle to output environment.
1936 * \todo gmx_output_env_t needs lifetime management.
1938 gmx_output_env_t* outputEnvironment_ = nullptr;
1940 /*! \brief Non-owning handle to MD log file.
1942 * \todo Context should own output facilities for client.
1943 * \todo Improve log file handle management.
1944 * \internal
1945 * Code managing the FILE* relies on the ability to set it to
1946 * nullptr to check whether the filehandle is valid.
1948 t_fileio* logFileHandle_ = nullptr;
1951 * \brief Builder for simulation stop signal handler.
1953 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1956 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1957 compat::not_null<SimulationContext*> context) :
1958 mdModules_(std::move(mdModules))
1960 communicator_ = context->communicator_;
1961 multiSimulation_ = context->multiSimulation_.get();
1964 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1966 Mdrunner::BuilderImplementation&
1967 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1968 const real forceWarningThreshold,
1969 const StartingBehavior startingBehavior)
1971 mdrunOptions_ = options;
1972 forceWarningThreshold_ = forceWarningThreshold;
1973 startingBehavior_ = startingBehavior;
1974 return *this;
1977 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1979 domdecOptions_ = options;
1982 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1984 nstlist_ = nstlist;
1987 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1989 replicaExchangeParameters_ = params;
1992 Mdrunner Mdrunner::BuilderImplementation::build()
1994 auto newRunner = Mdrunner(std::move(mdModules_));
1996 newRunner.mdrunOptions = mdrunOptions_;
1997 newRunner.pforce = forceWarningThreshold_;
1998 newRunner.startingBehavior = startingBehavior_;
1999 newRunner.domdecOptions = domdecOptions_;
2001 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2002 newRunner.hw_opt = hardwareOptions_;
2004 // No invariant to check. This parameter exists to optionally override other behavior.
2005 newRunner.nstlist_cmdline = nstlist_;
2007 newRunner.replExParams = replicaExchangeParameters_;
2009 newRunner.filenames = filenames_;
2011 newRunner.communicator = communicator_;
2013 // nullptr is a valid value for the multisim handle
2014 newRunner.ms = multiSimulation_;
2016 // \todo Clarify ownership and lifetime management for gmx_output_env_t
2017 // \todo Update sanity checking when output environment has clearly specified invariants.
2018 // Initialization and default values for oenv are not well specified in the current version.
2019 if (outputEnvironment_)
2021 newRunner.oenv = outputEnvironment_;
2023 else
2025 GMX_THROW(gmx::APIError(
2026 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2029 newRunner.logFileHandle = logFileHandle_;
2031 if (nbpu_opt_)
2033 newRunner.nbpu_opt = nbpu_opt_;
2035 else
2037 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2040 if (pme_opt_ && pme_fft_opt_)
2042 newRunner.pme_opt = pme_opt_;
2043 newRunner.pme_fft_opt = pme_fft_opt_;
2045 else
2047 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2050 if (bonded_opt_)
2052 newRunner.bonded_opt = bonded_opt_;
2054 else
2056 GMX_THROW(gmx::APIError(
2057 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2060 if (update_opt_)
2062 newRunner.update_opt = update_opt_;
2064 else
2066 GMX_THROW(gmx::APIError(
2067 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2071 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2073 if (stopHandlerBuilder_)
2075 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2077 else
2079 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2082 return newRunner;
2085 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2087 nbpu_opt_ = nbpu_opt;
2090 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2092 pme_opt_ = pme_opt;
2093 pme_fft_opt_ = pme_fft_opt;
2096 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2098 bonded_opt_ = bonded_opt;
2101 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2103 update_opt_ = update_opt;
2106 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2108 hardwareOptions_ = hardwareOptions;
2111 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2113 filenames_ = filenames;
2116 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2118 outputEnvironment_ = outputEnvironment;
2121 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2123 logFileHandle_ = logFileHandle;
2126 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2128 stopHandlerBuilder_ = std::move(builder);
2131 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2132 compat::not_null<SimulationContext*> context) :
2133 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2137 MdrunnerBuilder::~MdrunnerBuilder() = default;
2139 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2140 real forceWarningThreshold,
2141 const StartingBehavior startingBehavior)
2143 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2144 return *this;
2147 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2149 impl_->addDomdec(options);
2150 return *this;
2153 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2155 impl_->addVerletList(nstlist);
2156 return *this;
2159 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2161 impl_->addReplicaExchange(params);
2162 return *this;
2165 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2167 impl_->addNonBonded(nbpu_opt);
2168 return *this;
2171 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2173 // The builder method may become more general in the future, but in this version,
2174 // parameters for PME electrostatics are both required and the only parameters
2175 // available.
2176 if (pme_opt && pme_fft_opt)
2178 impl_->addPME(pme_opt, pme_fft_opt);
2180 else
2182 GMX_THROW(
2183 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2185 return *this;
2188 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2190 impl_->addBondedTaskAssignment(bonded_opt);
2191 return *this;
2194 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2196 impl_->addUpdateTaskAssignment(update_opt);
2197 return *this;
2200 Mdrunner MdrunnerBuilder::build()
2202 return impl_->build();
2205 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2207 impl_->addHardwareOptions(hardwareOptions);
2208 return *this;
2211 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2213 impl_->addFilenames(filenames);
2214 return *this;
2217 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2219 impl_->addOutputEnvironment(outputEnvironment);
2220 return *this;
2223 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2225 impl_->addLogFile(logFileHandle);
2226 return *this;
2229 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2231 impl_->addStopHandlerBuilder(std::move(builder));
2232 return *this;
2235 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2237 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;
2239 } // namespace gmx