2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/device_context.h"
77 #include "gromacs/gpu_utils/device_stream_manager.h"
78 #include "gromacs/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/mdrun/simulationinput.h"
109 #include "gromacs/mdrun/simulationinputhandle.h"
110 #include "gromacs/mdrunutility/handlerestart.h"
111 #include "gromacs/mdrunutility/logging.h"
112 #include "gromacs/mdrunutility/multisim.h"
113 #include "gromacs/mdrunutility/printtime.h"
114 #include "gromacs/mdrunutility/threadaffinity.h"
115 #include "gromacs/mdtypes/checkpointdata.h"
116 #include "gromacs/mdtypes/commrec.h"
117 #include "gromacs/mdtypes/enerdata.h"
118 #include "gromacs/mdtypes/fcdata.h"
119 #include "gromacs/mdtypes/forcerec.h"
120 #include "gromacs/mdtypes/group.h"
121 #include "gromacs/mdtypes/inputrec.h"
122 #include "gromacs/mdtypes/interaction_const.h"
123 #include "gromacs/mdtypes/md_enums.h"
124 #include "gromacs/mdtypes/mdatom.h"
125 #include "gromacs/mdtypes/mdrunoptions.h"
126 #include "gromacs/mdtypes/observableshistory.h"
127 #include "gromacs/mdtypes/simulation_workload.h"
128 #include "gromacs/mdtypes/state.h"
129 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
130 #include "gromacs/modularsimulator/modularsimulator.h"
131 #include "gromacs/nbnxm/gpu_data_mgmt.h"
132 #include "gromacs/nbnxm/nbnxm.h"
133 #include "gromacs/nbnxm/pairlist_tuning.h"
134 #include "gromacs/pbcutil/pbc.h"
135 #include "gromacs/pulling/output.h"
136 #include "gromacs/pulling/pull.h"
137 #include "gromacs/pulling/pull_rotation.h"
138 #include "gromacs/restraint/manager.h"
139 #include "gromacs/restraint/restraintmdmodule.h"
140 #include "gromacs/restraint/restraintpotential.h"
141 #include "gromacs/swap/swapcoords.h"
142 #include "gromacs/taskassignment/decidegpuusage.h"
143 #include "gromacs/taskassignment/decidesimulationworkload.h"
144 #include "gromacs/taskassignment/resourcedivision.h"
145 #include "gromacs/taskassignment/taskassignment.h"
146 #include "gromacs/taskassignment/usergpuids.h"
147 #include "gromacs/timing/gpu_timing.h"
148 #include "gromacs/timing/wallcycle.h"
149 #include "gromacs/timing/wallcyclereporting.h"
150 #include "gromacs/topology/mtop_util.h"
151 #include "gromacs/trajectory/trajectoryframe.h"
152 #include "gromacs/utility/basenetwork.h"
153 #include "gromacs/utility/cstringutil.h"
154 #include "gromacs/utility/exceptions.h"
155 #include "gromacs/utility/fatalerror.h"
156 #include "gromacs/utility/filestream.h"
157 #include "gromacs/utility/gmxassert.h"
158 #include "gromacs/utility/gmxmpi.h"
159 #include "gromacs/utility/keyvaluetree.h"
160 #include "gromacs/utility/logger.h"
161 #include "gromacs/utility/loggerbuilder.h"
162 #include "gromacs/utility/mdmodulenotification.h"
163 #include "gromacs/utility/physicalnodecommunicator.h"
164 #include "gromacs/utility/pleasecite.h"
165 #include "gromacs/utility/programcontext.h"
166 #include "gromacs/utility/smalloc.h"
167 #include "gromacs/utility/stringutil.h"
169 #include "isimulator.h"
170 #include "membedholder.h"
171 #include "replicaexchange.h"
172 #include "simulatorbuilder.h"
175 # include "corewrap.h"
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
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.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
213 devFlags.enableGpuPmePPComm =
214 GMX_GPU_CUDA && GMX_THREAD_MPI && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
216 #pragma GCC diagnostic pop
218 if (devFlags.enableGpuBufferOps)
220 GMX_LOG(mdlog.warning)
222 .appendTextFormatted(
223 "This run uses the 'GPU buffer ops' feature, enabled by the "
224 "GMX_USE_GPU_BUFFER_OPS environment variable.");
227 if (devFlags.forceGpuUpdateDefault)
229 GMX_LOG(mdlog.warning)
231 .appendTextFormatted(
232 "This run will default to '-update gpu' as requested by the "
233 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
234 "decomposition lacks substantial testing and should be used with caution.");
237 if (devFlags.enableGpuHaloExchange)
239 if (useGpuForNonbonded)
241 if (!devFlags.enableGpuBufferOps)
243 GMX_LOG(mdlog.warning)
245 .appendTextFormatted(
246 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
247 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
248 devFlags.enableGpuBufferOps = true;
250 GMX_LOG(mdlog.warning)
252 .appendTextFormatted(
253 "This run has requested the 'GPU halo exchange' feature, enabled by "
255 "GMX_GPU_DD_COMMS environment variable.");
259 GMX_LOG(mdlog.warning)
261 .appendTextFormatted(
262 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
263 "halo exchange' feature will not be enabled as nonbonded interactions "
264 "are not offloaded.");
265 devFlags.enableGpuHaloExchange = false;
269 if (devFlags.enableGpuPmePPComm)
271 if (pmeRunMode == PmeRunMode::GPU)
273 if (!devFlags.enableGpuBufferOps)
275 GMX_LOG(mdlog.warning)
277 .appendTextFormatted(
278 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
279 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
280 devFlags.enableGpuBufferOps = true;
282 GMX_LOG(mdlog.warning)
284 .appendTextFormatted(
285 "This run uses the 'GPU PME-PP communications' feature, enabled "
286 "by the GMX_GPU_PME_PP_COMMS environment variable.");
290 std::string clarification;
291 if (pmeRunMode == PmeRunMode::Mixed)
294 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
299 clarification = "PME is not offloaded to the GPU.";
301 GMX_LOG(mdlog.warning)
304 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
305 "'GPU PME-PP communications' feature was not enabled as "
307 devFlags.enableGpuPmePPComm = false;
314 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
316 * Used to ensure that the master thread does not modify mdrunner during copy
317 * on the spawned threads. */
318 static void threadMpiMdrunnerAccessBarrier()
321 MPI_Barrier(MPI_COMM_WORLD);
325 Mdrunner Mdrunner::cloneOnSpawnedThread() const
327 auto newRunner = Mdrunner(std::make_unique<MDModules>());
329 // All runners in the same process share a restraint manager resource because it is
330 // part of the interface to the client code, which is associated only with the
331 // original thread. Handles to the same resources can be obtained by copy.
333 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
336 // Copy members of master runner.
337 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
338 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
339 newRunner.hw_opt = hw_opt;
340 newRunner.filenames = filenames;
342 newRunner.oenv = oenv;
343 newRunner.mdrunOptions = mdrunOptions;
344 newRunner.domdecOptions = domdecOptions;
345 newRunner.nbpu_opt = nbpu_opt;
346 newRunner.pme_opt = pme_opt;
347 newRunner.pme_fft_opt = pme_fft_opt;
348 newRunner.bonded_opt = bonded_opt;
349 newRunner.update_opt = update_opt;
350 newRunner.nstlist_cmdline = nstlist_cmdline;
351 newRunner.replExParams = replExParams;
352 newRunner.pforce = pforce;
353 // Give the spawned thread the newly created valid communicator
354 // for the simulation.
355 newRunner.communicator = MPI_COMM_WORLD;
357 newRunner.startingBehavior = startingBehavior;
358 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
359 newRunner.inputHolder_ = inputHolder_;
361 threadMpiMdrunnerAccessBarrier();
366 /*! \brief The callback used for running on spawned threads.
368 * Obtains the pointer to the master mdrunner object from the one
369 * argument permitted to the thread-launch API call, copies it to make
370 * a new runner for this thread, reinitializes necessary data, and
371 * proceeds to the simulation. */
372 static void mdrunner_start_fn(const void* arg)
376 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
377 /* copy the arg list to make sure that it's thread-local. This
378 doesn't copy pointed-to items, of course; fnm, cr and fplog
379 are reset in the call below, all others should be const. */
380 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
383 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
387 void Mdrunner::spawnThreads(int numThreadsToLaunch)
390 /* now spawn new threads that start mdrunner_start_fn(), while
391 the main thread returns. Thread affinity is handled later. */
392 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
393 static_cast<const void*>(this))
396 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
399 // Give the master thread the newly created valid communicator for
401 communicator = MPI_COMM_WORLD;
402 threadMpiMdrunnerAccessBarrier();
404 GMX_UNUSED_VALUE(numThreadsToLaunch);
405 GMX_UNUSED_VALUE(mdrunner_start_fn);
411 /*! \brief Initialize variables for Verlet scheme simulation */
412 static void prepare_verlet_scheme(FILE* fplog,
416 const gmx_mtop_t* mtop,
418 bool makeGpuPairList,
419 const gmx::CpuInfo& cpuinfo)
421 // We checked the cut-offs in grompp, but double-check here.
422 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
423 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
425 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
426 "With Verlet lists and PME we should have rcoulomb>=rvdw");
430 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
431 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
433 /* For NVE simulations, we will retain the initial list buffer */
434 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
436 /* Update the Verlet buffer size for the current run setup */
438 /* Here we assume SIMD-enabled kernels are being used. But as currently
439 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
440 * and 4x2 gives a larger buffer than 4x4, this is ok.
442 ListSetupType listType =
443 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
444 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
446 const real rlist_new =
447 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
449 if (rlist_new != ir->rlist)
451 if (fplog != nullptr)
454 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
455 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
457 ir->rlist = rlist_new;
461 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
463 gmx_fatal(FARGS, "Can not set nstlist without %s",
464 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
467 if (EI_DYNAMICS(ir->eI))
469 /* Set or try nstlist values */
470 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
474 /*! \brief Override the nslist value in inputrec
476 * with value passed on the command line (if any)
478 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
482 /* override with anything else than the default -2 */
483 if (nsteps_cmdline > -2)
485 char sbuf_steps[STEPSTRSIZE];
486 char sbuf_msg[STRLEN];
488 ir->nsteps = nsteps_cmdline;
489 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
492 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
493 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
497 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
498 gmx_step_str(nsteps_cmdline, sbuf_steps));
501 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
503 else if (nsteps_cmdline < -2)
505 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
507 /* Do nothing if nsteps_cmdline == -2 */
513 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
515 * If not, and if a warning may be issued, logs a warning about
516 * falling back to CPU code. With thread-MPI, only the first
517 * call to this function should have \c issueWarning true. */
518 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
520 bool gpuIsUseful = true;
523 if (ir.opts.ngener - ir.nwall > 1)
525 /* The GPU code does not support more than one energy group.
526 * If the user requested GPUs explicitly, a fatal error is given later.
530 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
531 "For better performance, run on the GPU without energy groups and then do "
532 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
538 warning = "TPI is not implemented for GPUs.";
541 if (!gpuIsUseful && issueWarning)
543 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
549 //! Initializes the logger for mdrun.
550 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
552 gmx::LoggerBuilder builder;
553 if (fplog != nullptr)
555 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
557 if (isSimulationMasterRank)
559 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
561 return builder.build();
564 //! Make a TaskTarget from an mdrun argument string.
565 static TaskTarget findTaskTarget(const char* optionString)
567 TaskTarget returnValue = TaskTarget::Auto;
569 if (strncmp(optionString, "auto", 3) == 0)
571 returnValue = TaskTarget::Auto;
573 else if (strncmp(optionString, "cpu", 3) == 0)
575 returnValue = TaskTarget::Cpu;
577 else if (strncmp(optionString, "gpu", 3) == 0)
579 returnValue = TaskTarget::Gpu;
583 GMX_ASSERT(false, "Option string should have been checked for sanity already");
589 //! Finish run, aggregate data to print performance info.
590 static void finish_run(FILE* fplog,
591 const gmx::MDLogger& mdlog,
593 const t_inputrec* inputrec,
595 gmx_wallcycle_t wcycle,
596 gmx_walltime_accounting_t walltime_accounting,
597 nonbonded_verlet_t* nbv,
598 const gmx_pme_t* pme,
602 double nbfs = 0, mflop = 0;
603 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
604 elapsed_time_over_all_threads_over_all_ranks;
605 /* Control whether it is valid to print a report. Only the
606 simulation master may print, but it should not do so if the run
607 terminated e.g. before a scheduled reset step. This is
608 complicated by the fact that PME ranks are unaware of the
609 reason why they were sent a pmerecvqxFINISH. To avoid
610 communication deadlocks, we always do the communication for the
611 report, even if we've decided not to write the report, because
612 how long it takes to finish the run is not important when we've
613 decided not to report on the simulation performance.
615 Further, we only report performance for dynamical integrators,
616 because those are the only ones for which we plan to
617 consider doing any optimizations. */
618 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
620 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
622 GMX_LOG(mdlog.warning)
624 .appendText("Simulation ended prematurely, no performance report will be written.");
629 std::unique_ptr<t_nrnb> nrnbTotalStorage;
632 nrnbTotalStorage = std::make_unique<t_nrnb>();
633 nrnb_tot = nrnbTotalStorage.get();
635 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
643 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
644 elapsed_time_over_all_threads =
645 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
649 /* reduce elapsed_time over all MPI ranks in the current simulation */
650 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
652 elapsed_time_over_all_ranks /= cr->nnodes;
653 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
654 * current simulation. */
655 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
656 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
661 elapsed_time_over_all_ranks = elapsed_time;
662 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
667 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
670 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
672 print_dd_statistics(cr, inputrec, fplog);
675 /* TODO Move the responsibility for any scaling by thread counts
676 * to the code that handled the thread region, so that there's a
677 * mechanism to keep cycle counting working during the transition
678 * to task parallelism. */
679 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
680 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
681 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
682 nthreads_pp, nthreads_pme);
683 auto cycle_sum(wallcycle_sum(cr, wcycle));
687 auto nbnxn_gpu_timings =
688 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
689 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
691 if (pme_gpu_task_enabled(pme))
693 pme_gpu_get_timings(pme, &pme_gpu_timings);
695 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
696 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
699 if (EI_DYNAMICS(inputrec->eI))
701 delta_t = inputrec->delta_t;
706 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
707 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
708 delta_t, nbfs, mflop);
712 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
713 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
714 delta_t, nbfs, mflop);
719 int Mdrunner::mdrunner()
722 t_forcerec* fr = nullptr;
723 real ewaldcoeff_q = 0;
724 real ewaldcoeff_lj = 0;
725 int nChargePerturbed = -1, nTypePerturbed = 0;
726 gmx_wallcycle_t wcycle;
727 gmx_walltime_accounting_t walltime_accounting = nullptr;
728 MembedHolder membedHolder(filenames.size(), filenames.data());
729 gmx_hw_info_t* hwinfo = nullptr;
731 /* CAUTION: threads may be started later on in this function, so
732 cr doesn't reflect the final parallel state right now */
735 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
736 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
737 const bool doRerun = mdrunOptions.rerun;
739 // Handle task-assignment related user options.
740 EmulateGpuNonbonded emulateGpuNonbonded =
741 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
743 std::vector<int> userGpuTaskAssignment;
746 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
748 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
749 auto nonbondedTarget = findTaskTarget(nbpu_opt);
750 auto pmeTarget = findTaskTarget(pme_opt);
751 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
752 auto bondedTarget = findTaskTarget(bonded_opt);
753 auto updateTarget = findTaskTarget(update_opt);
755 FILE* fplog = nullptr;
756 // If we are appending, we don't write log output because we need
757 // to check that the old log file matches what the checkpoint file
758 // expects. Otherwise, we should start to write log output now if
759 // there is a file ready for it.
760 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
762 fplog = gmx_fio_getfp(logFileHandle);
764 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
765 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
766 gmx::MDLogger mdlog(logOwner.logger());
768 // TODO The thread-MPI master rank makes a working
769 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
770 // after the threads have been launched. This works because no use
771 // is made of that communicator until after the execution paths
772 // have rejoined. But it is likely that we can improve the way
773 // this is expressed, e.g. by expressly running detection only the
774 // master rank for thread-MPI, rather than relying on the mutex
775 // and reference count.
776 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
777 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
779 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
781 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->deviceInfoList, hw_opt.gpuIdsAvailable);
783 // Print citation requests after all software/hardware printing
784 pleaseCiteGromacs(fplog);
786 // Note: legacy program logic relies on checking whether these pointers are assigned.
787 // Objects may or may not be allocated later.
788 std::unique_ptr<t_inputrec> inputrec;
789 std::unique_ptr<t_state> globalState;
791 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
793 if (isSimulationMasterRank)
795 // Allocate objects to be initialized by later function calls.
796 /* Only the master rank has the global state */
797 globalState = std::make_unique<t_state>();
798 inputrec = std::make_unique<t_inputrec>();
800 /* Read (nearly) all data required for the simulation
801 * and keep the partly serialized tpr contents to send to other ranks later
803 applyGlobalSimulationState(*inputHolder_.get(), partialDeserializedTpr.get(),
804 globalState.get(), inputrec.get(), &mtop);
807 /* Check and update the hardware options for internal consistency */
808 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
811 if (GMX_THREAD_MPI && isSimulationMasterRank)
813 bool useGpuForNonbonded = false;
814 bool useGpuForPme = false;
817 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
819 // If the user specified the number of ranks, then we must
820 // respect that, but in default mode, we need to allow for
821 // the number of GPUs to choose the number of ranks.
822 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
823 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
824 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
825 canUseGpuForNonbonded,
826 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
827 hw_opt.nthreads_tmpi);
828 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
829 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
830 *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
832 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
834 /* Determine how many thread-MPI ranks to start.
836 * TODO Over-writing the user-supplied value here does
837 * prevent any possible subsequent checks from working
839 hw_opt.nthreads_tmpi =
840 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
841 inputrec.get(), &mtop, mdlog, membedHolder.doMembed());
843 // Now start the threads for thread MPI.
844 spawnThreads(hw_opt.nthreads_tmpi);
845 // The spawned threads enter mdrunner() and execution of
846 // master and spawned threads joins at the end of this block.
847 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
850 GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
851 "Must have valid world communicator unless running a multi-simulation");
852 CommrecHandle crHandle = init_commrec(communicator);
853 t_commrec* cr = crHandle.get();
854 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
858 /* now broadcast everything to the non-master nodes/threads: */
859 if (!isSimulationMasterRank)
861 // Until now, only the master rank has a non-null pointer.
862 // On non-master ranks, allocate the object that will receive data in the following call.
863 inputrec = std::make_unique<t_inputrec>();
865 init_parallel(cr->mpiDefaultCommunicator, MASTER(cr), inputrec.get(), &mtop,
866 partialDeserializedTpr.get());
868 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
869 partialDeserializedTpr.reset(nullptr);
871 // Now the number of ranks is known to all ranks, and each knows
872 // the inputrec read by the master rank. The ranks can now all run
873 // the task-deciding functions and will agree on the result
874 // without needing to communicate.
875 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
877 // Note that these variables describe only their own node.
879 // Note that when bonded interactions run on a GPU they always run
880 // alongside a nonbonded task, so do not influence task assignment
881 // even though they affect the force calculation workload.
882 bool useGpuForNonbonded = false;
883 bool useGpuForPme = false;
884 bool useGpuForBonded = false;
885 bool useGpuForUpdate = false;
886 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
889 // It's possible that there are different numbers of GPUs on
890 // different nodes, which is the user's responsibility to
891 // handle. If unsuitable, we will notice that during task
893 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
894 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
895 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
896 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
897 useGpuForPme = decideWhetherToUseGpusForPme(
898 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec,
899 cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
900 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
901 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
902 useGpuForBonded = decideWhetherToUseGpusForBonded(
903 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
904 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
905 domdecOptions.numPmeRanks, gpusWereDetected);
907 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
909 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
911 // Initialize development feature flags that enabled by environment variable
912 // and report those features that are enabled.
913 const DevelopmentFeatureFlags devFlags =
914 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
916 const bool useModularSimulator =
917 checkUseModularSimulator(false, inputrec.get(), doRerun, mtop, ms, replExParams,
918 nullptr, doEssentialDynamics, membedHolder.doMembed());
921 // TODO: hide restraint implementation details from Mdrunner.
922 // There is nothing unique about restraints at this point as far as the
923 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
924 // factory functions from the SimulationContext on which to call mdModules_->add().
925 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
926 for (auto&& restraint : restraintManager_->getRestraints())
928 auto module = RestraintMDModule::create(restraint, restraint->sites());
929 mdModules_->add(std::move(module));
932 // TODO: Error handling
933 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
934 // now that the MdModules know their options, they know which callbacks to sign up to
935 mdModules_->subscribeToSimulationSetupNotifications();
936 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
938 if (inputrec->internalParameters != nullptr)
940 mdModulesNotifier.notify(*inputrec->internalParameters);
943 if (fplog != nullptr)
945 pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
946 fprintf(fplog, "\n");
951 /* In rerun, set velocities to zero if present */
952 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
954 // rerun does not use velocities
958 "Rerun trajectory contains velocities. Rerun does only evaluate "
959 "potential energy and forces. The velocities will be ignored.");
960 for (int i = 0; i < globalState->natoms; i++)
962 clear_rvec(globalState->v[i]);
964 globalState->flags &= ~(1 << estV);
967 /* now make sure the state is initialized and propagated */
968 set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
971 /* NM and TPI parallelize over force/energy calculations, not atoms,
972 * so we need to initialize and broadcast the global state.
974 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
978 globalState = std::make_unique<t_state>();
980 broadcastStateWithoutDynamics(cr->mpiDefaultCommunicator, DOMAINDECOMP(cr), PAR(cr),
984 /* A parallel command line option consistency check that we can
985 only do after any threads have started. */
987 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
988 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
991 "The -dd or -npme option request a parallel simulation, "
993 "but %s was compiled without threads or MPI enabled",
994 output_env_get_program_display_name(oenv));
996 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
998 "but %s was not started through mpirun/mpiexec or only one rank was requested "
999 "through mpirun/mpiexec",
1000 output_env_get_program_display_name(oenv));
1004 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1007 "The .mdp file specified an energy mininization or normal mode algorithm, and "
1008 "these are not compatible with mdrun -rerun");
1011 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1013 if (domdecOptions.numPmeRanks > 0)
1015 gmx_fatal_collective(FARGS, cr->mpiDefaultCommunicator, MASTER(cr),
1016 "PME-only ranks are requested, but the system does not use PME "
1017 "for electrostatics or LJ");
1020 domdecOptions.numPmeRanks = 0;
1023 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1025 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1026 * improve performance with many threads per GPU, since our OpenMP
1027 * scaling is bad, but it's difficult to automate the setup.
1029 domdecOptions.numPmeRanks = 0;
1033 if (domdecOptions.numPmeRanks < 0)
1035 domdecOptions.numPmeRanks = 0;
1036 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1040 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1041 "PME GPU decomposition is not supported");
1048 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1052 /* NMR restraints must be initialized before load_checkpoint,
1053 * since with time averaging the history is added to t_state.
1054 * For proper consistency check we therefore need to extend
1056 * So the PME-only nodes (if present) will also initialize
1057 * the distance restraints.
1060 /* This needs to be called before read_checkpoint to extend the state */
1061 t_disresdata* disresdata;
1062 snew(disresdata, 1);
1063 init_disres(fplog, &mtop, inputrec.get(), DisResRunMode::MDRun,
1064 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1065 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1066 globalState.get(), replExParams.exchangeInterval > 0);
1068 t_oriresdata* oriresdata;
1069 snew(oriresdata, 1);
1070 init_orires(fplog, &mtop, inputrec.get(), cr, ms, globalState.get(), oriresdata);
1072 auto deform = prepareBoxDeformation(
1073 globalState != nullptr ? globalState->box : box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1074 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mygroup, *inputrec);
1076 ObservablesHistory observablesHistory = {};
1078 auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1079 if (startingBehavior != StartingBehavior::NewSimulation)
1081 /* Check if checkpoint file exists before doing continuation.
1082 * This way we can use identical input options for the first and subsequent runs...
1084 if (mdrunOptions.numStepsCommandline > -2)
1086 /* Temporarily set the number of steps to unlimited to avoid
1087 * triggering the nsteps check in load_checkpoint().
1088 * This hack will go away soon when the -nsteps option is removed.
1090 inputrec->nsteps = -1;
1093 // Finish applying initial simulation state information from external sources on all ranks.
1094 // Reconcile checkpoint file data with Mdrunner state established up to this point.
1095 applyLocalState(*inputHolder_.get(), logFileHandle, cr, domdecOptions.numCells,
1096 inputrec.get(), globalState.get(), &observablesHistory,
1097 mdrunOptions.reproducible, mdModules_->notifier(),
1098 modularSimulatorCheckpointData.get(), useModularSimulator);
1099 // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1101 // on all code paths.
1102 // Write checkpoint or provide hook to update SimulationInput.
1103 // If there was a checkpoint file, SimulationInput contains more information
1104 // than if there wasn't. At this point, we have synchronized the in-memory
1105 // state with the filesystem state only for restarted simulations. We should
1106 // be calling applyLocalState unconditionally and expect that the completeness
1107 // of SimulationInput is not dependent on its creation method.
1109 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1111 // Now we can start normal logging to the truncated log file.
1112 fplog = gmx_fio_getfp(logFileHandle);
1113 prepareLogAppending(fplog);
1114 logOwner = buildLogger(fplog, MASTER(cr));
1115 mdlog = logOwner.logger();
1119 if (mdrunOptions.numStepsCommandline > -2)
1124 "The -nsteps functionality is deprecated, and may be removed in a future "
1126 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1129 /* override nsteps with value set on the commandline */
1130 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1132 if (isSimulationMasterRank)
1134 copy_mat(globalState->box, box);
1139 gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1142 if (inputrec->cutoff_scheme != ecutsVERLET)
1145 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1146 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1148 /* Update rlist and nstlist. */
1149 /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1150 * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1151 * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1153 prepare_verlet_scheme(fplog, cr, inputrec.get(), nstlist_cmdline, &mtop, box,
1154 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1157 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1158 // This builder is necessary while we have multi-part construction
1159 // of DD. Before DD is constructed, we use the existence of
1160 // the builder object to indicate that further construction of DD
1162 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1163 if (useDomainDecomposition)
1165 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1166 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1167 positionsFromStatePointer(globalState.get()));
1171 /* PME, if used, is done on all nodes with 1D decomposition */
1172 cr->nnodes = cr->sizeOfDefaultCommunicator;
1173 cr->sim_nodeid = cr->rankInDefaultCommunicator;
1174 cr->nodeid = cr->rankInDefaultCommunicator;
1176 cr->duty = (DUTY_PP | DUTY_PME);
1178 if (inputrec->pbcType == PbcType::Screw)
1180 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1184 // Produce the task assignment for this rank - done after DD is constructed
1185 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1186 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1187 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1188 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1189 // TODO cr->duty & DUTY_PME should imply that a PME
1190 // algorithm is active, but currently does not.
1191 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1193 // Get the device handles for the modules, nullptr when no task is assigned.
1195 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1197 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1198 bool useTiming = true;
1202 /* WARNING: CUDA timings are incorrect with multiple streams.
1203 * This is the main reason why they are disabled by default.
1205 // TODO: Consider turning on by default when we can detect nr of streams.
1206 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1208 else if (GMX_GPU_OPENCL)
1210 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1213 // TODO Currently this is always built, yet DD partition code
1214 // checks if it is built before using it. Probably it should
1215 // become an MDModule that is made only when another module
1216 // requires it (e.g. pull, CompEl, density fitting), so that we
1217 // don't update the local atom sets unilaterally every step.
1218 LocalAtomSetManager atomSets;
1221 // TODO Pass the GPU streams to ddBuilder to use in buffer
1222 // transfers (e.g. halo exchange)
1223 cr->dd = ddBuilder->build(&atomSets);
1224 // The builder's job is done, so destruct it
1225 ddBuilder.reset(nullptr);
1226 // Note that local state still does not exist yet.
1229 // The GPU update is decided here because we need to know whether the constraints or
1230 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1231 // defined). This is only known after DD is initialized, hence decision on using GPU
1232 // update is done so late.
1235 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1237 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1238 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1239 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1240 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1241 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1243 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1245 const bool printHostName = (cr->nnodes > 1);
1246 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1248 MdrunScheduleWorkload runScheduleWork;
1249 // Also populates the simulation constant workload description.
1250 runScheduleWork.simulationWork = createSimulationWorkload(
1251 *inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded, useGpuForUpdate,
1252 devFlags.enableGpuBufferOps, devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1254 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1256 if (deviceInfo != nullptr)
1258 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1260 dd_setup_dlb_resource_sharing(cr, deviceId);
1262 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1263 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, 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));
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);
1284 GMX_LOG(mdlog.warning)
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,
1294 cr->nnodes == 1 ? "thread" : "threads"
1296 cr->nnodes == 1 ? "process" : "processes"
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;
1324 const bool bEnableFPE = false;
1326 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
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 && (runScheduleWork.simulationWork.useGpuHaloExchange
1350 || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1352 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1355 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1357 /* Before setting affinity, check whether the affinity has changed
1358 * - which indicates that probably the OpenMP library has changed it
1359 * since we first checked).
1361 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1363 int numThreadsOnThisNode, intraNodeThreadOffset;
1364 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1365 &intraNodeThreadOffset);
1367 /* Set the CPU affinity */
1368 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1369 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1372 if (mdrunOptions.timingOptions.resetStep > -1)
1377 "The -resetstep functionality is deprecated, and may be removed in a "
1380 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1384 /* Master synchronizes its value of reset_counters with all nodes
1385 * including PME only nodes */
1386 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1387 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1388 wcycle_set_reset_counters(wcycle, reset_counters);
1391 // Membrane embedding must be initialized before we call init_forcerec()
1392 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec.get(),
1393 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1395 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1396 std::unique_ptr<MDAtoms> mdAtoms;
1397 std::unique_ptr<VirtualSitesHandler> vsite;
1398 std::unique_ptr<GpuBonded> gpuBonded;
1401 if (thisRankHasDuty(cr, DUTY_PP))
1403 mdModulesNotifier.notify(*cr);
1404 mdModulesNotifier.notify(&atomSets);
1405 mdModulesNotifier.notify(inputrec->pbcType);
1406 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1407 /* Initiate forcerecord */
1408 fr = new t_forcerec;
1409 fr->forceProviders = mdModules_->initForceProviders();
1410 init_forcerec(fplog, mdlog, fr, inputrec.get(), &mtop, cr, box,
1411 opt2fn("-table", filenames.size(), filenames.data()),
1412 opt2fn("-tablep", filenames.size(), filenames.data()),
1413 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1414 // Dirty hack, for fixing disres and orires should be made mdmodules
1415 fr->fcdata->disres = disresdata;
1416 fr->fcdata->orires = oriresdata;
1418 // Save a handle to device stream manager to use elsewhere in the code
1419 // TODO: Forcerec is not a correct place to store it.
1420 fr->deviceStreamManager = deviceStreamManager.get();
1422 if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1425 deviceStreamManager != nullptr,
1426 "GPU device stream manager should be valid in order to use PME-PP direct "
1429 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1430 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1432 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1433 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1434 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1437 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec.get(), fr, cr, *hwinfo,
1438 runScheduleWork.simulationWork.useGpuNonbonded,
1439 deviceStreamManager.get(), &mtop, box, wcycle);
1440 // TODO: Move the logic below to a GPU bonded builder
1441 if (runScheduleWork.simulationWork.useGpuBonded)
1443 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1444 "GPU device stream manager should be valid in order to use GPU "
1445 "version of bonded forces.");
1446 gpuBonded = std::make_unique<GpuBonded>(
1447 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1448 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1449 fr->gpuBonded = gpuBonded.get();
1452 /* Initialize the mdAtoms structure.
1453 * mdAtoms is not filled with atom data,
1454 * as this can not be done now with domain decomposition.
1456 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1457 if (globalState && thisRankHasPmeGpuTask)
1459 // The pinning of coordinates in the global state object works, because we only use
1460 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1461 // points to the global state object without DD.
1462 // FIXME: MD and EM separately set up the local state - this should happen in the same
1463 // function, which should also perform the pinning.
1464 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1467 /* Initialize the virtual site communication */
1468 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1470 calc_shifts(box, fr->shift_vec);
1472 /* With periodic molecules the charge groups should be whole at start up
1473 * and the virtual sites should not be far from their proper positions.
1475 if (!inputrec->bContinuation && MASTER(cr)
1476 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1478 /* Make molecules whole at start of run */
1479 if (fr->pbcType != PbcType::No)
1481 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1485 /* Correct initial vsite positions are required
1486 * for the initial distribution in the domain decomposition
1487 * and for the initial shell prediction.
1489 constructVirtualSitesGlobal(mtop, globalState->x);
1493 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1495 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1496 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1501 /* This is a PME only node */
1503 GMX_ASSERT(globalState == nullptr,
1504 "We don't need the state on a PME only rank and expect it to be unitialized");
1506 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1507 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1510 gmx_pme_t* sepPmeData = nullptr;
1511 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1512 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1513 "Double-checking that only PME-only ranks have no forcerec");
1514 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1516 // TODO should live in ewald module once its testing is improved
1518 // Later, this program could contain kernels that might be later
1519 // re-used as auto-tuning progresses, or subsequent simulations
1521 PmeGpuProgramStorage pmeGpuProgram;
1522 if (thisRankHasPmeGpuTask)
1525 (deviceStreamManager != nullptr),
1526 "GPU device stream manager should be initialized in order to use GPU for PME.");
1527 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1528 "GPU device should be initialized in order to use GPU for PME.");
1529 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1532 /* Initiate PME if necessary,
1533 * either on all nodes or on dedicated PME nodes only. */
1534 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1536 if (mdAtoms && mdAtoms->mdatoms())
1538 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1539 if (EVDW_PME(inputrec->vdwtype))
1541 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1544 if (cr->npmenodes > 0)
1546 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1547 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1548 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1551 if (thisRankHasDuty(cr, DUTY_PME))
1555 // TODO: This should be in the builder.
1556 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1557 || (deviceStreamManager != nullptr),
1558 "Device stream manager should be valid in order to use GPU "
1561 !runScheduleWork.simulationWork.useGpuPme
1562 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1563 "GPU PME stream should be valid in order to use GPU version of PME.");
1565 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1566 ? &deviceStreamManager->context()
1568 const DeviceStream* pmeStream =
1569 runScheduleWork.simulationWork.useGpuPme
1570 ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1573 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec.get(),
1574 nChargePerturbed != 0, nTypePerturbed != 0,
1575 mdrunOptions.reproducible, ewaldcoeff_q, ewaldcoeff_lj,
1576 gmx_omp_nthreads_get(emntPME), pmeRunMode, nullptr,
1577 deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1579 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1584 if (EI_DYNAMICS(inputrec->eI))
1586 /* Turn on signal handling on all nodes */
1588 * (A user signal from the PME nodes (if any)
1589 * is communicated to the PP nodes.
1591 signal_handler_install();
1594 pull_t* pull_work = nullptr;
1595 if (thisRankHasDuty(cr, DUTY_PP))
1597 /* Assumes uniform use of the number of OpenMP threads */
1598 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1600 if (inputrec->bPull)
1602 /* Initialize pull code */
1603 pull_work = init_pull(fplog, inputrec->pull, inputrec.get(), &mtop, cr, &atomSets,
1604 inputrec->fepvals->init_lambda);
1605 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1607 initPullHistory(pull_work, &observablesHistory);
1609 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1611 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1615 std::unique_ptr<EnforcedRotation> enforcedRotation;
1618 /* Initialize enforced rotation code */
1619 enforcedRotation = init_rot(fplog, inputrec.get(), filenames.size(), filenames.data(),
1620 cr, &atomSets, globalState.get(), &mtop, oenv, mdrunOptions,
1624 t_swap* swap = nullptr;
1625 if (inputrec->eSwapCoords != eswapNO)
1627 /* Initialize ion swapping code */
1628 swap = init_swapcoords(fplog, inputrec.get(),
1629 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1630 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1631 oenv, mdrunOptions, startingBehavior);
1634 /* Let makeConstraints know whether we have essential dynamics constraints. */
1635 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1636 ms, &nrnb, wcycle, fr->bMolPBC);
1638 /* Energy terms and groups */
1639 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1640 inputrec->fepvals->n_lambda);
1642 /* Kinetic energy data */
1643 gmx_ekindata_t ekind;
1644 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1646 /* Set up interactive MD (IMD) */
1648 makeImdSession(inputrec.get(), cr, wcycle, &enerd, ms, &mtop, mdlog,
1649 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1650 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1652 if (DOMAINDECOMP(cr))
1654 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1655 /* This call is not included in init_domain_decomposition mainly
1656 * because fr->cginfo_mb is set later.
1658 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec.get(),
1659 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1662 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1663 if (gpusWereDetected
1664 && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
1665 || runScheduleWork.simulationWork.useGpuBufferOps))
1667 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1668 ? GpuApiCallBehavior::Async
1669 : GpuApiCallBehavior::Sync;
1670 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1671 "GPU device stream manager should be initialized to use GPU.");
1672 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1673 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1674 fr->stateGpu = stateGpu.get();
1677 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1678 SimulatorBuilder simulatorBuilder;
1680 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1681 simulatorBuilder.add(std::move(membedHolder));
1682 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1683 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1686 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1687 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1688 simulatorBuilder.add(ConstraintsParam(
1689 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1691 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1692 simulatorBuilder.add(LegacyInput(static_cast<int>(filenames.size()), filenames.data(),
1693 inputrec.get(), fr));
1694 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1695 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1696 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1697 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1698 // Todo move to an MDModule
1699 simulatorBuilder.add(IonSwapping(swap));
1700 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1701 simulatorBuilder.add(BoxDeformationHandle(deform.get()));
1702 simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
1704 // build and run simulator object based on user-input
1705 auto simulator = simulatorBuilder.build(useModularSimulator);
1708 if (fr->pmePpCommGpu)
1710 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1711 fr->pmePpCommGpu.reset();
1714 if (inputrec->bPull)
1716 finish_pull(pull_work);
1718 finish_swapcoords(swap);
1722 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1724 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1725 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec.get(), pmeRunMode,
1726 deviceStreamManager.get());
1729 wallcycle_stop(wcycle, ewcRUN);
1731 /* Finish up, write some stuff
1732 * if rerunMD, don't write last frame again
1734 finish_run(fplog, mdlog, cr, inputrec.get(), &nrnb, wcycle, walltime_accounting,
1735 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1737 // clean up cycle counter
1738 wallcycle_destroy(wcycle);
1740 deviceStreamManager.reset(nullptr);
1744 gmx_pme_destroy(pmedata);
1748 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1749 // before we destroy the GPU context(s)
1750 // Pinned buffers are associated with contexts in CUDA.
1751 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1752 mdAtoms.reset(nullptr);
1753 globalState.reset(nullptr);
1754 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1755 gpuBonded.reset(nullptr);
1756 /* Free pinned buffers in *fr */
1759 // TODO convert to C++ so we can get rid of these frees
1763 if (!hwinfo->deviceInfoList.empty())
1765 /* stop the GPU profiler (only CUDA) */
1769 /* With tMPI we need to wait for all ranks to finish deallocation before
1770 * destroying the CUDA context as some tMPI ranks may be sharing
1773 * This is not a concern in OpenCL where we use one context per rank.
1775 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1776 * but it is easier and more futureproof to call it on the whole node.
1778 * Note that this function needs to be called even if GPUs are not used
1779 * in this run because the PME ranks have no knowledge of whether GPUs
1780 * are used or not, but all ranks need to enter the barrier below.
1781 * \todo Remove this physical node barrier after making sure
1782 * that it's not needed anymore (with a shared GPU run).
1786 physicalNodeComm.barrier();
1788 releaseDevice(deviceInfo);
1790 /* Does what it says */
1791 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1792 walltime_accounting_destroy(walltime_accounting);
1794 // Ensure log file content is written
1797 gmx_fio_flush(logFileHandle);
1800 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1801 * exceptions were enabled before function was called. */
1804 gmx_fedisableexcept();
1807 auto rc = static_cast<int>(gmx_get_stop_condition());
1810 /* we need to join all threads. The sub-threads join when they
1811 exit this function, but the master thread needs to be told to
1821 Mdrunner::~Mdrunner()
1823 // Clean up of the Manager.
1824 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1825 // but okay as long as threads synchronize some time before adding or accessing
1826 // a new set of restraints.
1827 if (restraintManager_)
1829 restraintManager_->clear();
1830 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1831 "restraints added during runner life time should be cleared at runner "
1836 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1838 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1839 // Not sure if this should be logged through the md logger or something else,
1840 // but it is helpful to have some sort of INFO level message sent somewhere.
1841 // std::cout << "Registering restraint named " << name << std::endl;
1843 // When multiple restraints are used, it may be wasteful to register them separately.
1844 // Maybe instead register an entire Restraint Manager as a force provider.
1845 restraintManager_->addToSpec(std::move(puller), name);
1848 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1850 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1852 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1854 class Mdrunner::BuilderImplementation
1857 BuilderImplementation() = delete;
1858 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1859 ~BuilderImplementation();
1861 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1862 real forceWarningThreshold,
1863 StartingBehavior startingBehavior);
1865 void addDomdec(const DomdecOptions& options);
1867 void addInput(SimulationInputHandle inputHolder);
1869 void addVerletList(int nstlist);
1871 void addReplicaExchange(const ReplicaExchangeParameters& params);
1873 void addNonBonded(const char* nbpu_opt);
1875 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1877 void addBondedTaskAssignment(const char* bonded_opt);
1879 void addUpdateTaskAssignment(const char* update_opt);
1881 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1883 void addFilenames(ArrayRef<const t_filenm> filenames);
1885 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1887 void addLogFile(t_fileio* logFileHandle);
1889 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1894 // Default parameters copied from runner.h
1895 // \todo Clarify source(s) of default parameters.
1897 const char* nbpu_opt_ = nullptr;
1898 const char* pme_opt_ = nullptr;
1899 const char* pme_fft_opt_ = nullptr;
1900 const char* bonded_opt_ = nullptr;
1901 const char* update_opt_ = nullptr;
1903 MdrunOptions mdrunOptions_;
1905 DomdecOptions domdecOptions_;
1907 ReplicaExchangeParameters replicaExchangeParameters_;
1909 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1912 //! Multisim communicator handle.
1913 gmx_multisim_t* multiSimulation_;
1915 //! mdrun communicator
1916 MPI_Comm communicator_ = MPI_COMM_NULL;
1918 //! Print a warning if any force is larger than this (in kJ/mol nm).
1919 real forceWarningThreshold_ = -1;
1921 //! Whether the simulation will start afresh, or restart with/without appending.
1922 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1924 //! The modules that comprise the functionality of mdrun.
1925 std::unique_ptr<MDModules> mdModules_;
1927 //! \brief Parallelism information.
1928 gmx_hw_opt_t hardwareOptions_;
1930 //! filename options for simulation.
1931 ArrayRef<const t_filenm> filenames_;
1933 /*! \brief Handle to output environment.
1935 * \todo gmx_output_env_t needs lifetime management.
1937 gmx_output_env_t* outputEnvironment_ = nullptr;
1939 /*! \brief Non-owning handle to MD log file.
1941 * \todo Context should own output facilities for client.
1942 * \todo Improve log file handle management.
1944 * Code managing the FILE* relies on the ability to set it to
1945 * nullptr to check whether the filehandle is valid.
1947 t_fileio* logFileHandle_ = nullptr;
1950 * \brief Builder for simulation stop signal handler.
1952 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1954 SimulationInputHandle inputHolder_;
1957 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1958 compat::not_null<SimulationContext*> context) :
1959 mdModules_(std::move(mdModules))
1961 communicator_ = context->communicator_;
1962 multiSimulation_ = context->multiSimulation_.get();
1965 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1967 Mdrunner::BuilderImplementation&
1968 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1969 const real forceWarningThreshold,
1970 const StartingBehavior startingBehavior)
1972 mdrunOptions_ = options;
1973 forceWarningThreshold_ = forceWarningThreshold;
1974 startingBehavior_ = startingBehavior;
1978 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1980 domdecOptions_ = options;
1983 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1988 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1990 replicaExchangeParameters_ = params;
1993 Mdrunner Mdrunner::BuilderImplementation::build()
1995 auto newRunner = Mdrunner(std::move(mdModules_));
1997 newRunner.mdrunOptions = mdrunOptions_;
1998 newRunner.pforce = forceWarningThreshold_;
1999 newRunner.startingBehavior = startingBehavior_;
2000 newRunner.domdecOptions = domdecOptions_;
2002 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2003 newRunner.hw_opt = hardwareOptions_;
2005 // No invariant to check. This parameter exists to optionally override other behavior.
2006 newRunner.nstlist_cmdline = nstlist_;
2008 newRunner.replExParams = replicaExchangeParameters_;
2010 newRunner.filenames = filenames_;
2012 newRunner.communicator = communicator_;
2014 // nullptr is a valid value for the multisim handle
2015 newRunner.ms = multiSimulation_;
2019 newRunner.inputHolder_ = std::move(inputHolder_);
2023 GMX_THROW(gmx::APIError("MdrunnerBuilder::addInput() is required before build()."));
2026 // \todo Clarify ownership and lifetime management for gmx_output_env_t
2027 // \todo Update sanity checking when output environment has clearly specified invariants.
2028 // Initialization and default values for oenv are not well specified in the current version.
2029 if (outputEnvironment_)
2031 newRunner.oenv = outputEnvironment_;
2035 GMX_THROW(gmx::APIError(
2036 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2039 newRunner.logFileHandle = logFileHandle_;
2043 newRunner.nbpu_opt = nbpu_opt_;
2047 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2050 if (pme_opt_ && pme_fft_opt_)
2052 newRunner.pme_opt = pme_opt_;
2053 newRunner.pme_fft_opt = pme_fft_opt_;
2057 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2062 newRunner.bonded_opt = bonded_opt_;
2066 GMX_THROW(gmx::APIError(
2067 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2072 newRunner.update_opt = update_opt_;
2076 GMX_THROW(gmx::APIError(
2077 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2081 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2083 if (stopHandlerBuilder_)
2085 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2089 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2095 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2097 nbpu_opt_ = nbpu_opt;
2100 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2103 pme_fft_opt_ = pme_fft_opt;
2106 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2108 bonded_opt_ = bonded_opt;
2111 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2113 update_opt_ = update_opt;
2116 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2118 hardwareOptions_ = hardwareOptions;
2121 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2123 filenames_ = filenames;
2126 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2128 outputEnvironment_ = outputEnvironment;
2131 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2133 logFileHandle_ = logFileHandle;
2136 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2138 stopHandlerBuilder_ = std::move(builder);
2141 void Mdrunner::BuilderImplementation::addInput(SimulationInputHandle inputHolder)
2143 inputHolder_ = std::move(inputHolder);
2146 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2147 compat::not_null<SimulationContext*> context) :
2148 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2152 MdrunnerBuilder::~MdrunnerBuilder() = default;
2154 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2155 real forceWarningThreshold,
2156 const StartingBehavior startingBehavior)
2158 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2162 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2164 impl_->addDomdec(options);
2168 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2170 impl_->addVerletList(nstlist);
2174 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2176 impl_->addReplicaExchange(params);
2180 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2182 impl_->addNonBonded(nbpu_opt);
2186 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2188 // The builder method may become more general in the future, but in this version,
2189 // parameters for PME electrostatics are both required and the only parameters
2191 if (pme_opt && pme_fft_opt)
2193 impl_->addPME(pme_opt, pme_fft_opt);
2198 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2203 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2205 impl_->addBondedTaskAssignment(bonded_opt);
2209 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2211 impl_->addUpdateTaskAssignment(update_opt);
2215 Mdrunner MdrunnerBuilder::build()
2217 return impl_->build();
2220 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2222 impl_->addHardwareOptions(hardwareOptions);
2226 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2228 impl_->addFilenames(filenames);
2232 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2234 impl_->addOutputEnvironment(outputEnvironment);
2238 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2240 impl_->addLogFile(logFileHandle);
2244 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2246 impl_->addStopHandlerBuilder(std::move(builder));
2250 MdrunnerBuilder& MdrunnerBuilder::addInput(SimulationInputHandle input)
2252 impl_->addInput(std::move(input));
2256 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2258 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;