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/gpu_utils.h"
77 #include "gromacs/hardware/cpuinfo.h"
78 #include "gromacs/hardware/detecthardware.h"
79 #include "gromacs/hardware/printhardware.h"
80 #include "gromacs/imd/imd.h"
81 #include "gromacs/listed_forces/disre.h"
82 #include "gromacs/listed_forces/gpubonded.h"
83 #include "gromacs/listed_forces/orires.h"
84 #include "gromacs/math/functions.h"
85 #include "gromacs/math/utilities.h"
86 #include "gromacs/math/vec.h"
87 #include "gromacs/mdlib/boxdeformation.h"
88 #include "gromacs/mdlib/broadcaststructs.h"
89 #include "gromacs/mdlib/calc_verletbuf.h"
90 #include "gromacs/mdlib/dispersioncorrection.h"
91 #include "gromacs/mdlib/enerdata_utils.h"
92 #include "gromacs/mdlib/force.h"
93 #include "gromacs/mdlib/forcerec.h"
94 #include "gromacs/mdlib/gmx_omp_nthreads.h"
95 #include "gromacs/mdlib/makeconstraints.h"
96 #include "gromacs/mdlib/md_support.h"
97 #include "gromacs/mdlib/mdatoms.h"
98 #include "gromacs/mdlib/membed.h"
99 #include "gromacs/mdlib/qmmm.h"
100 #include "gromacs/mdlib/sighandler.h"
101 #include "gromacs/mdlib/stophandler.h"
102 #include "gromacs/mdlib/updategroups.h"
103 #include "gromacs/mdrun/mdmodules.h"
104 #include "gromacs/mdrun/simulationcontext.h"
105 #include "gromacs/mdrunutility/handlerestart.h"
106 #include "gromacs/mdrunutility/logging.h"
107 #include "gromacs/mdrunutility/multisim.h"
108 #include "gromacs/mdrunutility/printtime.h"
109 #include "gromacs/mdrunutility/threadaffinity.h"
110 #include "gromacs/mdtypes/commrec.h"
111 #include "gromacs/mdtypes/enerdata.h"
112 #include "gromacs/mdtypes/fcdata.h"
113 #include "gromacs/mdtypes/group.h"
114 #include "gromacs/mdtypes/inputrec.h"
115 #include "gromacs/mdtypes/md_enums.h"
116 #include "gromacs/mdtypes/mdrunoptions.h"
117 #include "gromacs/mdtypes/observableshistory.h"
118 #include "gromacs/mdtypes/simulation_workload.h"
119 #include "gromacs/mdtypes/state.h"
120 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
121 #include "gromacs/nbnxm/gpu_data_mgmt.h"
122 #include "gromacs/nbnxm/nbnxm.h"
123 #include "gromacs/nbnxm/pairlist_tuning.h"
124 #include "gromacs/pbcutil/pbc.h"
125 #include "gromacs/pulling/output.h"
126 #include "gromacs/pulling/pull.h"
127 #include "gromacs/pulling/pull_rotation.h"
128 #include "gromacs/restraint/manager.h"
129 #include "gromacs/restraint/restraintmdmodule.h"
130 #include "gromacs/restraint/restraintpotential.h"
131 #include "gromacs/swap/swapcoords.h"
132 #include "gromacs/taskassignment/decidegpuusage.h"
133 #include "gromacs/taskassignment/decidesimulationworkload.h"
134 #include "gromacs/taskassignment/resourcedivision.h"
135 #include "gromacs/taskassignment/taskassignment.h"
136 #include "gromacs/taskassignment/usergpuids.h"
137 #include "gromacs/timing/gpu_timing.h"
138 #include "gromacs/timing/wallcycle.h"
139 #include "gromacs/timing/wallcyclereporting.h"
140 #include "gromacs/topology/mtop_util.h"
141 #include "gromacs/trajectory/trajectoryframe.h"
142 #include "gromacs/utility/basenetwork.h"
143 #include "gromacs/utility/cstringutil.h"
144 #include "gromacs/utility/exceptions.h"
145 #include "gromacs/utility/fatalerror.h"
146 #include "gromacs/utility/filestream.h"
147 #include "gromacs/utility/gmxassert.h"
148 #include "gromacs/utility/gmxmpi.h"
149 #include "gromacs/utility/keyvaluetree.h"
150 #include "gromacs/utility/logger.h"
151 #include "gromacs/utility/loggerbuilder.h"
152 #include "gromacs/utility/mdmodulenotification.h"
153 #include "gromacs/utility/physicalnodecommunicator.h"
154 #include "gromacs/utility/pleasecite.h"
155 #include "gromacs/utility/programcontext.h"
156 #include "gromacs/utility/smalloc.h"
157 #include "gromacs/utility/stringutil.h"
159 #include "isimulator.h"
160 #include "replicaexchange.h"
161 #include "simulatorbuilder.h"
164 # include "corewrap.h"
170 /*! \brief Structure that holds boolean flags corresponding to the development
171 * features present enabled through environment variables.
174 struct DevelopmentFeatureFlags
176 //! True if the Buffer ops development feature is enabled
177 // TODO: when the trigger of the buffer ops offload is fully automated this should go away
178 bool enableGpuBufferOps = false;
179 //! If true, forces 'mdrun -update auto' default to 'gpu'
180 bool forceGpuUpdateDefault = false;
181 //! True if the GPU halo exchange development feature is enabled
182 bool enableGpuHaloExchange = false;
183 //! True if the PME PP direct communication GPU development feature is enabled
184 bool enableGpuPmePPComm = false;
187 /*! \brief Manage any development feature flag variables encountered
189 * The use of dev features indicated by environment variables is
190 * logged in order to ensure that runs with such features enabled can
191 * be identified from their log and standard output. Any cross
192 * dependencies are also checked, and if unsatisfied, a fatal error
195 * Note that some development features overrides are applied already here:
196 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
198 * \param[in] mdlog Logger object.
199 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
200 * \param[in] pmeRunMode The PME run mode for this run
201 * \returns The object populated with development feature flags.
203 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
204 const bool useGpuForNonbonded,
205 const PmeRunMode pmeRunMode)
207 DevelopmentFeatureFlags devFlags;
209 // Some builds of GCC 5 give false positive warnings that these
210 // getenv results are ignored when clearly they are used.
211 #pragma GCC diagnostic push
212 #pragma GCC diagnostic ignored "-Wunused-result"
213 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
214 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
215 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
216 devFlags.enableGpuHaloExchange =
217 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
218 devFlags.enableGpuPmePPComm =
219 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
220 #pragma GCC diagnostic pop
222 if (devFlags.enableGpuBufferOps)
224 GMX_LOG(mdlog.warning)
226 .appendTextFormatted(
227 "This run uses the 'GPU buffer ops' feature, enabled by the "
228 "GMX_USE_GPU_BUFFER_OPS environment variable.");
231 if (devFlags.forceGpuUpdateDefault)
233 GMX_LOG(mdlog.warning)
235 .appendTextFormatted(
236 "This run will default to '-update gpu' as requested by the "
237 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
238 "decomposition lacks substantial testing and should be used with caution.");
241 if (devFlags.enableGpuHaloExchange)
243 if (useGpuForNonbonded)
245 if (!devFlags.enableGpuBufferOps)
247 GMX_LOG(mdlog.warning)
249 .appendTextFormatted(
250 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
251 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
252 devFlags.enableGpuBufferOps = true;
254 GMX_LOG(mdlog.warning)
256 .appendTextFormatted(
257 "This run uses the 'GPU halo exchange' feature, enabled by the "
258 "GMX_GPU_DD_COMMS environment variable.");
262 GMX_LOG(mdlog.warning)
264 .appendTextFormatted(
265 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
266 "halo exchange' feature will not be enabled as nonbonded interactions "
267 "are not offloaded.");
268 devFlags.enableGpuHaloExchange = false;
272 if (devFlags.enableGpuPmePPComm)
274 if (pmeRunMode == PmeRunMode::GPU)
276 GMX_LOG(mdlog.warning)
278 .appendTextFormatted(
279 "This run uses the 'GPU PME-PP communications' feature, enabled "
280 "by the GMX_GPU_PME_PP_COMMS environment variable.");
284 std::string clarification;
285 if (pmeRunMode == PmeRunMode::Mixed)
288 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
293 clarification = "PME is not offloaded to the GPU.";
295 GMX_LOG(mdlog.warning)
298 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
299 "'GPU PME-PP communications' feature was not enabled as "
301 devFlags.enableGpuPmePPComm = false;
308 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
310 * Used to ensure that the master thread does not modify mdrunner during copy
311 * on the spawned threads. */
312 static void threadMpiMdrunnerAccessBarrier()
315 MPI_Barrier(MPI_COMM_WORLD);
319 Mdrunner Mdrunner::cloneOnSpawnedThread() const
321 auto newRunner = Mdrunner(std::make_unique<MDModules>());
323 // All runners in the same process share a restraint manager resource because it is
324 // part of the interface to the client code, which is associated only with the
325 // original thread. Handles to the same resources can be obtained by copy.
327 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
330 // Copy members of master runner.
331 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
332 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
333 newRunner.hw_opt = hw_opt;
334 newRunner.filenames = filenames;
336 newRunner.oenv = oenv;
337 newRunner.mdrunOptions = mdrunOptions;
338 newRunner.domdecOptions = domdecOptions;
339 newRunner.nbpu_opt = nbpu_opt;
340 newRunner.pme_opt = pme_opt;
341 newRunner.pme_fft_opt = pme_fft_opt;
342 newRunner.bonded_opt = bonded_opt;
343 newRunner.update_opt = update_opt;
344 newRunner.nstlist_cmdline = nstlist_cmdline;
345 newRunner.replExParams = replExParams;
346 newRunner.pforce = pforce;
347 // Give the spawned thread the newly created valid communicator
348 // for the simulation.
349 newRunner.communicator = MPI_COMM_WORLD;
351 newRunner.startingBehavior = startingBehavior;
352 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
354 threadMpiMdrunnerAccessBarrier();
359 /*! \brief The callback used for running on spawned threads.
361 * Obtains the pointer to the master mdrunner object from the one
362 * argument permitted to the thread-launch API call, copies it to make
363 * a new runner for this thread, reinitializes necessary data, and
364 * proceeds to the simulation. */
365 static void mdrunner_start_fn(const void* arg)
369 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
370 /* copy the arg list to make sure that it's thread-local. This
371 doesn't copy pointed-to items, of course; fnm, cr and fplog
372 are reset in the call below, all others should be const. */
373 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
376 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
380 void Mdrunner::spawnThreads(int numThreadsToLaunch)
383 /* now spawn new threads that start mdrunner_start_fn(), while
384 the main thread returns. Thread affinity is handled later. */
385 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
386 static_cast<const void*>(this))
389 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
392 // Give the master thread the newly created valid communicator for
394 communicator = MPI_COMM_WORLD;
395 threadMpiMdrunnerAccessBarrier();
397 GMX_UNUSED_VALUE(numThreadsToLaunch);
398 GMX_UNUSED_VALUE(mdrunner_start_fn);
404 /*! \brief Initialize variables for Verlet scheme simulation */
405 static void prepare_verlet_scheme(FILE* fplog,
409 const gmx_mtop_t* mtop,
411 bool makeGpuPairList,
412 const gmx::CpuInfo& cpuinfo)
414 // We checked the cut-offs in grompp, but double-check here.
415 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
416 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
418 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
419 "With Verlet lists and PME we should have rcoulomb>=rvdw");
423 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
424 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
426 /* For NVE simulations, we will retain the initial list buffer */
427 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
429 /* Update the Verlet buffer size for the current run setup */
431 /* Here we assume SIMD-enabled kernels are being used. But as currently
432 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
433 * and 4x2 gives a larger buffer than 4x4, this is ok.
435 ListSetupType listType =
436 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
437 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
439 const real rlist_new =
440 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
442 if (rlist_new != ir->rlist)
444 if (fplog != nullptr)
447 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
448 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
450 ir->rlist = rlist_new;
454 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
456 gmx_fatal(FARGS, "Can not set nstlist without %s",
457 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
460 if (EI_DYNAMICS(ir->eI))
462 /* Set or try nstlist values */
463 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
467 /*! \brief Override the nslist value in inputrec
469 * with value passed on the command line (if any)
471 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
475 /* override with anything else than the default -2 */
476 if (nsteps_cmdline > -2)
478 char sbuf_steps[STEPSTRSIZE];
479 char sbuf_msg[STRLEN];
481 ir->nsteps = nsteps_cmdline;
482 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
485 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
486 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
490 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
491 gmx_step_str(nsteps_cmdline, sbuf_steps));
494 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
496 else if (nsteps_cmdline < -2)
498 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
500 /* Do nothing if nsteps_cmdline == -2 */
506 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
508 * If not, and if a warning may be issued, logs a warning about
509 * falling back to CPU code. With thread-MPI, only the first
510 * call to this function should have \c issueWarning true. */
511 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
513 bool gpuIsUseful = true;
516 if (ir.opts.ngener - ir.nwall > 1)
518 /* The GPU code does not support more than one energy group.
519 * If the user requested GPUs explicitly, a fatal error is given later.
523 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
524 "For better performance, run on the GPU without energy groups and then do "
525 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
531 warning = "TPI is not implemented for GPUs.";
534 if (!gpuIsUseful && issueWarning)
536 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
542 //! Initializes the logger for mdrun.
543 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
545 gmx::LoggerBuilder builder;
546 if (fplog != nullptr)
548 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
550 if (isSimulationMasterRank)
552 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
554 return builder.build();
557 //! Make a TaskTarget from an mdrun argument string.
558 static TaskTarget findTaskTarget(const char* optionString)
560 TaskTarget returnValue = TaskTarget::Auto;
562 if (strncmp(optionString, "auto", 3) == 0)
564 returnValue = TaskTarget::Auto;
566 else if (strncmp(optionString, "cpu", 3) == 0)
568 returnValue = TaskTarget::Cpu;
570 else if (strncmp(optionString, "gpu", 3) == 0)
572 returnValue = TaskTarget::Gpu;
576 GMX_ASSERT(false, "Option string should have been checked for sanity already");
582 //! Finish run, aggregate data to print performance info.
583 static void finish_run(FILE* fplog,
584 const gmx::MDLogger& mdlog,
586 const t_inputrec* inputrec,
588 gmx_wallcycle_t wcycle,
589 gmx_walltime_accounting_t walltime_accounting,
590 nonbonded_verlet_t* nbv,
591 const gmx_pme_t* pme,
595 double nbfs = 0, mflop = 0;
596 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
597 elapsed_time_over_all_threads_over_all_ranks;
598 /* Control whether it is valid to print a report. Only the
599 simulation master may print, but it should not do so if the run
600 terminated e.g. before a scheduled reset step. This is
601 complicated by the fact that PME ranks are unaware of the
602 reason why they were sent a pmerecvqxFINISH. To avoid
603 communication deadlocks, we always do the communication for the
604 report, even if we've decided not to write the report, because
605 how long it takes to finish the run is not important when we've
606 decided not to report on the simulation performance.
608 Further, we only report performance for dynamical integrators,
609 because those are the only ones for which we plan to
610 consider doing any optimizations. */
611 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
613 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
615 GMX_LOG(mdlog.warning)
617 .appendText("Simulation ended prematurely, no performance report will be written.");
622 std::unique_ptr<t_nrnb> nrnbTotalStorage;
625 nrnbTotalStorage = std::make_unique<t_nrnb>();
626 nrnb_tot = nrnbTotalStorage.get();
628 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
636 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
637 elapsed_time_over_all_threads =
638 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
642 /* reduce elapsed_time over all MPI ranks in the current simulation */
643 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
645 elapsed_time_over_all_ranks /= cr->nnodes;
646 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
647 * current simulation. */
648 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
649 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
654 elapsed_time_over_all_ranks = elapsed_time;
655 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
660 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
663 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
665 print_dd_statistics(cr, inputrec, fplog);
668 /* TODO Move the responsibility for any scaling by thread counts
669 * to the code that handled the thread region, so that there's a
670 * mechanism to keep cycle counting working during the transition
671 * to task parallelism. */
672 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
673 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
674 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
675 nthreads_pp, nthreads_pme);
676 auto cycle_sum(wallcycle_sum(cr, wcycle));
680 auto nbnxn_gpu_timings =
681 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
682 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
684 if (pme_gpu_task_enabled(pme))
686 pme_gpu_get_timings(pme, &pme_gpu_timings);
688 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
689 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
692 if (EI_DYNAMICS(inputrec->eI))
694 delta_t = inputrec->delta_t;
699 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
700 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
701 delta_t, nbfs, mflop);
705 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
706 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
707 delta_t, nbfs, mflop);
712 int Mdrunner::mdrunner()
715 t_forcerec* fr = nullptr;
716 t_fcdata* fcd = nullptr;
717 real ewaldcoeff_q = 0;
718 real ewaldcoeff_lj = 0;
719 int nChargePerturbed = -1, nTypePerturbed = 0;
720 gmx_wallcycle_t wcycle;
721 gmx_walltime_accounting_t walltime_accounting = nullptr;
722 gmx_membed_t* membed = nullptr;
723 gmx_hw_info_t* hwinfo = nullptr;
725 /* CAUTION: threads may be started later on in this function, so
726 cr doesn't reflect the final parallel state right now */
729 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
730 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
731 const bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
732 const bool doRerun = mdrunOptions.rerun;
734 // Handle task-assignment related user options.
735 EmulateGpuNonbonded emulateGpuNonbonded =
736 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
738 std::vector<int> userGpuTaskAssignment;
741 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
743 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
744 auto nonbondedTarget = findTaskTarget(nbpu_opt);
745 auto pmeTarget = findTaskTarget(pme_opt);
746 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
747 auto bondedTarget = findTaskTarget(bonded_opt);
748 auto updateTarget = findTaskTarget(update_opt);
750 FILE* fplog = nullptr;
751 // If we are appending, we don't write log output because we need
752 // to check that the old log file matches what the checkpoint file
753 // expects. Otherwise, we should start to write log output now if
754 // there is a file ready for it.
755 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
757 fplog = gmx_fio_getfp(logFileHandle);
759 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
760 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
761 gmx::MDLogger mdlog(logOwner.logger());
763 // TODO The thread-MPI master rank makes a working
764 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
765 // after the threads have been launched. This works because no use
766 // is made of that communicator until after the execution paths
767 // have rejoined. But it is likely that we can improve the way
768 // this is expressed, e.g. by expressly running detection only the
769 // master rank for thread-MPI, rather than relying on the mutex
770 // and reference count.
771 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
772 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
774 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
776 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
778 // Print citation requests after all software/hardware printing
779 pleaseCiteGromacs(fplog);
781 // TODO Replace this by unique_ptr once t_inputrec is C++
782 t_inputrec inputrecInstance;
783 t_inputrec* inputrec = nullptr;
784 std::unique_ptr<t_state> globalState;
786 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
788 if (isSimulationMasterRank)
790 /* Only the master rank has the global state */
791 globalState = std::make_unique<t_state>();
793 /* Read (nearly) all data required for the simulation
794 * and keep the partly serialized tpr contents to send to other ranks later
796 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
797 &inputrecInstance, globalState.get(), &mtop);
798 inputrec = &inputrecInstance;
801 /* Check and update the hardware options for internal consistency */
802 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
805 if (GMX_THREAD_MPI && isSimulationMasterRank)
807 bool useGpuForNonbonded = false;
808 bool useGpuForPme = false;
811 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
813 // If the user specified the number of ranks, then we must
814 // respect that, but in default mode, we need to allow for
815 // the number of GPUs to choose the number of ranks.
816 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
817 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
818 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
819 canUseGpuForNonbonded,
820 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
821 hw_opt.nthreads_tmpi);
822 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
823 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
824 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
826 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
828 /* Determine how many thread-MPI ranks to start.
830 * TODO Over-writing the user-supplied value here does
831 * prevent any possible subsequent checks from working
833 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
834 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
836 // Now start the threads for thread MPI.
837 spawnThreads(hw_opt.nthreads_tmpi);
838 // The spawned threads enter mdrunner() and execution of
839 // master and spawned threads joins at the end of this block.
840 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
843 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
844 CommrecHandle crHandle = init_commrec(communicator, ms);
845 t_commrec* cr = crHandle.get();
846 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
850 /* now broadcast everything to the non-master nodes/threads: */
851 if (!isSimulationMasterRank)
853 inputrec = &inputrecInstance;
855 init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
857 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
858 partialDeserializedTpr.reset(nullptr);
860 // Now the number of ranks is known to all ranks, and each knows
861 // the inputrec read by the master rank. The ranks can now all run
862 // the task-deciding functions and will agree on the result
863 // without needing to communicate.
864 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
866 // Note that these variables describe only their own node.
868 // Note that when bonded interactions run on a GPU they always run
869 // alongside a nonbonded task, so do not influence task assignment
870 // even though they affect the force calculation workload.
871 bool useGpuForNonbonded = false;
872 bool useGpuForPme = false;
873 bool useGpuForBonded = false;
874 bool useGpuForUpdate = false;
875 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
878 // It's possible that there are different numbers of GPUs on
879 // different nodes, which is the user's responsibility to
880 // handle. If unsuitable, we will notice that during task
882 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
883 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
884 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
885 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
886 useGpuForPme = decideWhetherToUseGpusForPme(
887 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
888 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
889 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
890 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
891 useGpuForBonded = decideWhetherToUseGpusForBonded(
892 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
893 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
894 domdecOptions.numPmeRanks, gpusWereDetected);
896 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
898 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
900 // Initialize development feature flags that enabled by environment variable
901 // and report those features that are enabled.
902 const DevelopmentFeatureFlags devFlags =
903 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
906 // TODO: hide restraint implementation details from Mdrunner.
907 // There is nothing unique about restraints at this point as far as the
908 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
909 // factory functions from the SimulationContext on which to call mdModules_->add().
910 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
911 for (auto&& restraint : restraintManager_->getRestraints())
913 auto module = RestraintMDModule::create(restraint, restraint->sites());
914 mdModules_->add(std::move(module));
917 // TODO: Error handling
918 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
919 const auto& mdModulesNotifier = mdModules_->notifier().notifier_;
921 if (inputrec->internalParameters != nullptr)
923 mdModulesNotifier.notify(*inputrec->internalParameters);
926 if (fplog != nullptr)
928 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
929 fprintf(fplog, "\n");
934 /* In rerun, set velocities to zero if present */
935 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
937 // rerun does not use velocities
941 "Rerun trajectory contains velocities. Rerun does only evaluate "
942 "potential energy and forces. The velocities will be ignored.");
943 for (int i = 0; i < globalState->natoms; i++)
945 clear_rvec(globalState->v[i]);
947 globalState->flags &= ~(1 << estV);
950 /* now make sure the state is initialized and propagated */
951 set_state_entries(globalState.get(), inputrec);
954 /* NM and TPI parallelize over force/energy calculations, not atoms,
955 * so we need to initialize and broadcast the global state.
957 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
961 globalState = std::make_unique<t_state>();
963 broadcastStateWithoutDynamics(cr, globalState.get());
966 /* A parallel command line option consistency check that we can
967 only do after any threads have started. */
969 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
970 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
973 "The -dd or -npme option request a parallel simulation, "
975 "but %s was compiled without threads or MPI enabled",
976 output_env_get_program_display_name(oenv));
978 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
980 "but %s was not started through mpirun/mpiexec or only one rank was requested "
981 "through mpirun/mpiexec",
982 output_env_get_program_display_name(oenv));
986 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
989 "The .mdp file specified an energy mininization or normal mode algorithm, and "
990 "these are not compatible with mdrun -rerun");
993 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
995 if (domdecOptions.numPmeRanks > 0)
997 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
998 "PME-only ranks are requested, but the system does not use PME "
999 "for electrostatics or LJ");
1002 domdecOptions.numPmeRanks = 0;
1005 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1007 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1008 * improve performance with many threads per GPU, since our OpenMP
1009 * scaling is bad, but it's difficult to automate the setup.
1011 domdecOptions.numPmeRanks = 0;
1015 if (domdecOptions.numPmeRanks < 0)
1017 domdecOptions.numPmeRanks = 0;
1018 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1022 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1023 "PME GPU decomposition is not supported");
1030 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1034 /* NMR restraints must be initialized before load_checkpoint,
1035 * since with time averaging the history is added to t_state.
1036 * For proper consistency check we therefore need to extend
1038 * So the PME-only nodes (if present) will also initialize
1039 * the distance restraints.
1043 /* This needs to be called before read_checkpoint to extend the state */
1044 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1046 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1048 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1050 ObservablesHistory observablesHistory = {};
1052 if (startingBehavior != StartingBehavior::NewSimulation)
1054 /* Check if checkpoint file exists before doing continuation.
1055 * This way we can use identical input options for the first and subsequent runs...
1057 if (mdrunOptions.numStepsCommandline > -2)
1059 /* Temporarily set the number of steps to unlimited to avoid
1060 * triggering the nsteps check in load_checkpoint().
1061 * This hack will go away soon when the -nsteps option is removed.
1063 inputrec->nsteps = -1;
1066 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1067 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1068 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1070 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1072 // Now we can start normal logging to the truncated log file.
1073 fplog = gmx_fio_getfp(logFileHandle);
1074 prepareLogAppending(fplog);
1075 logOwner = buildLogger(fplog, MASTER(cr));
1076 mdlog = logOwner.logger();
1080 if (mdrunOptions.numStepsCommandline > -2)
1085 "The -nsteps functionality is deprecated, and may be removed in a future "
1087 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1090 /* override nsteps with value set on the commandline */
1091 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1095 copy_mat(globalState->box, box);
1100 gmx_bcast(sizeof(box), box, cr);
1103 if (inputrec->cutoff_scheme != ecutsVERLET)
1106 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1107 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1109 /* Update rlist and nstlist. */
1110 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1111 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1114 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1115 // This builder is necessary while we have multi-part construction
1116 // of DD. Before DD is constructed, we use the existence of
1117 // the builder object to indicate that further construction of DD
1119 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1120 if (useDomainDecomposition)
1122 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1123 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1124 positionsFromStatePointer(globalState.get()));
1128 /* PME, if used, is done on all nodes with 1D decomposition */
1130 cr->duty = (DUTY_PP | DUTY_PME);
1132 if (inputrec->pbcType == PbcType::Screw)
1134 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1138 // Produce the task assignment for this rank.
1139 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1140 GpuTaskAssignments gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1141 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1142 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1143 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1144 // TODO cr->duty & DUTY_PME should imply that a PME
1145 // algorithm is active, but currently does not.
1146 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1148 // Get the device handles for the modules, nullptr when no task is assigned.
1149 gmx_device_info_t* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1150 gmx_device_info_t* pmeDeviceInfo = gpuTaskAssignments.initPmeDevice();
1152 // TODO Initialize GPU streams here.
1154 // TODO Currently this is always built, yet DD partition code
1155 // checks if it is built before using it. Probably it should
1156 // become an MDModule that is made only when another module
1157 // requires it (e.g. pull, CompEl, density fitting), so that we
1158 // don't update the local atom sets unilaterally every step.
1159 LocalAtomSetManager atomSets;
1162 // TODO Pass the GPU streams to ddBuilder to use in buffer
1163 // transfers (e.g. halo exchange)
1164 cr->dd = ddBuilder->build(&atomSets);
1165 // The builder's job is done, so destruct it
1166 ddBuilder.reset(nullptr);
1167 // Note that local state still does not exist yet.
1170 // The GPU update is decided here because we need to know whether the constraints or
1171 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1172 // defined). This is only known after DD is initialized, hence decision on using GPU
1173 // update is done so late.
1176 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1178 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1179 devFlags.forceGpuUpdateDefault, useDomainDecomposition, useUpdateGroups, pmeRunMode,
1180 domdecOptions.numPmeRanks > 0, useGpuForNonbonded, updateTarget, gpusWereDetected,
1181 *inputrec, mtop, doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1182 replExParams.exchangeInterval > 0, doRerun, mdlog);
1184 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1186 const bool printHostName = (cr->nnodes > 1);
1187 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1189 // If the user chose a task assignment, give them some hints
1190 // where appropriate.
1191 if (!userGpuTaskAssignment.empty())
1193 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1198 /* After possible communicator splitting in make_dd_communicators.
1199 * we can set up the intra/inter node communication.
1201 gmx_setup_nodecomm(fplog, cr);
1207 GMX_LOG(mdlog.warning)
1209 .appendTextFormatted(
1210 "This is simulation %d out of %d running as a composite GROMACS\n"
1211 "multi-simulation job. Setup for this simulation:\n",
1214 GMX_LOG(mdlog.warning)
1215 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1217 cr->nnodes == 1 ? "thread" : "threads"
1219 cr->nnodes == 1 ? "process" : "processes"
1225 // If mdrun -pin auto honors any affinity setting that already
1226 // exists. If so, it is nice to provide feedback about whether
1227 // that existing affinity setting was from OpenMP or something
1228 // else, so we run this code both before and after we initialize
1229 // the OpenMP support.
1230 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1231 /* Check and update the number of OpenMP threads requested */
1232 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1233 pmeRunMode, mtop, *inputrec);
1235 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1236 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1238 // Enable FP exception detection, but not in
1239 // Release mode and not for compilers with known buggy FP
1240 // exception support (clang with any optimization) or suspected
1241 // buggy FP exception support (gcc 7.* with optimization).
1242 #if !defined NDEBUG \
1243 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1244 && defined __OPTIMIZE__)
1245 const bool bEnableFPE = true;
1247 const bool bEnableFPE = false;
1249 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1252 gmx_feenableexcept();
1255 /* Now that we know the setup is consistent, check for efficiency */
1256 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1257 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1259 /* getting number of PP/PME threads on this MPI / tMPI rank.
1260 PME: env variable should be read only on one node to make sure it is
1261 identical everywhere;
1263 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1264 : gmx_omp_nthreads_get(emntPME);
1265 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1266 physicalNodeComm, mdlog);
1268 // Enable Peer access between GPUs where available
1269 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1270 // any of the GPU communication features are active.
1271 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1272 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1274 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1277 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1279 /* Before setting affinity, check whether the affinity has changed
1280 * - which indicates that probably the OpenMP library has changed it
1281 * since we first checked).
1283 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1285 int numThreadsOnThisNode, intraNodeThreadOffset;
1286 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1287 &intraNodeThreadOffset);
1289 /* Set the CPU affinity */
1290 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1291 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1294 if (mdrunOptions.timingOptions.resetStep > -1)
1299 "The -resetstep functionality is deprecated, and may be removed in a "
1302 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1306 /* Master synchronizes its value of reset_counters with all nodes
1307 * including PME only nodes */
1308 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1309 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1310 wcycle_set_reset_counters(wcycle, reset_counters);
1313 // Membrane embedding must be initialized before we call init_forcerec()
1318 fprintf(stderr, "Initializing membed");
1320 /* Note that membed cannot work in parallel because mtop is
1321 * changed here. Fix this if we ever want to make it run with
1322 * multiple ranks. */
1323 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1324 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1327 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1328 std::unique_ptr<MDAtoms> mdAtoms;
1329 std::unique_ptr<gmx_vsite_t> vsite;
1330 std::unique_ptr<GpuBonded> gpuBonded;
1333 if (thisRankHasDuty(cr, DUTY_PP))
1335 mdModulesNotifier.notify(*cr);
1336 mdModulesNotifier.notify(&atomSets);
1337 mdModulesNotifier.notify(PeriodicBoundaryConditionType{ inputrec->pbcType });
1338 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1339 /* Initiate forcerecord */
1340 fr = new t_forcerec;
1341 fr->forceProviders = mdModules_->initForceProviders();
1342 init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1343 opt2fn("-table", filenames.size(), filenames.data()),
1344 opt2fn("-tablep", filenames.size(), filenames.data()),
1345 opt2fns("-tableb", filenames.size(), filenames.data()),
1346 pmeRunMode == PmeRunMode::GPU && !thisRankHasDuty(cr, DUTY_PME), pforce);
1348 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, nonbondedDeviceInfo,
1349 &mtop, box, wcycle);
1350 if (useGpuForBonded)
1352 auto stream = havePPDomainDecomposition(cr)
1353 ? Nbnxm::gpu_get_command_stream(
1354 fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1355 : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1356 gmx::InteractionLocality::Local);
1357 gpuBonded = std::make_unique<GpuBonded>(mtop.ffparams, stream, wcycle);
1358 fr->gpuBonded = gpuBonded.get();
1361 // TODO Move this to happen during domain decomposition setup,
1362 // once stream and event handling works well with that.
1363 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
1364 if (havePPDomainDecomposition(cr) && prefer1DAnd1PulseDD && is1DAnd1PulseDD(*cr->dd))
1366 GMX_RELEASE_ASSERT(devFlags.enableGpuBufferOps,
1367 "Must use GMX_USE_GPU_BUFFER_OPS=1 to use GMX_GPU_DD_COMMS=1");
1369 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local);
1370 void* streamNonLocal =
1371 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal);
1372 GMX_LOG(mdlog.warning)
1374 .appendTextFormatted(
1375 "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the "
1376 "GMX_GPU_DD_COMMS environment variable.");
1377 cr->dd->gpuHaloExchange = std::make_unique<GpuHaloExchange>(
1378 cr->dd, cr->mpi_comm_mysim, streamLocal, streamNonLocal);
1381 /* Initialize the mdAtoms structure.
1382 * mdAtoms is not filled with atom data,
1383 * as this can not be done now with domain decomposition.
1385 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1386 if (globalState && thisRankHasPmeGpuTask)
1388 // The pinning of coordinates in the global state object works, because we only use
1389 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1390 // points to the global state object without DD.
1391 // FIXME: MD and EM separately set up the local state - this should happen in the same
1392 // function, which should also perform the pinning.
1393 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1396 /* Initialize the virtual site communication */
1397 vsite = initVsite(mtop, cr);
1399 calc_shifts(box, fr->shift_vec);
1401 /* With periodic molecules the charge groups should be whole at start up
1402 * and the virtual sites should not be far from their proper positions.
1404 if (!inputrec->bContinuation && MASTER(cr)
1405 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1407 /* Make molecules whole at start of run */
1408 if (fr->pbcType != PbcType::No)
1410 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1414 /* Correct initial vsite positions are required
1415 * for the initial distribution in the domain decomposition
1416 * and for the initial shell prediction.
1418 constructVsitesGlobal(mtop, globalState->x);
1422 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1424 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1425 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1430 /* This is a PME only node */
1432 GMX_ASSERT(globalState == nullptr,
1433 "We don't need the state on a PME only rank and expect it to be unitialized");
1435 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1436 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1439 gmx_pme_t* sepPmeData = nullptr;
1440 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1441 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1442 "Double-checking that only PME-only ranks have no forcerec");
1443 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1445 // TODO should live in ewald module once its testing is improved
1447 // Later, this program could contain kernels that might be later
1448 // re-used as auto-tuning progresses, or subsequent simulations
1450 PmeGpuProgramStorage pmeGpuProgram;
1451 if (thisRankHasPmeGpuTask)
1453 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1456 /* Initiate PME if necessary,
1457 * either on all nodes or on dedicated PME nodes only. */
1458 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1460 if (mdAtoms && mdAtoms->mdatoms())
1462 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1463 if (EVDW_PME(inputrec->vdwtype))
1465 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1468 if (cr->npmenodes > 0)
1470 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1471 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1472 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1475 if (thisRankHasDuty(cr, DUTY_PME))
1479 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1480 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1481 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1482 nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1484 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1489 if (EI_DYNAMICS(inputrec->eI))
1491 /* Turn on signal handling on all nodes */
1493 * (A user signal from the PME nodes (if any)
1494 * is communicated to the PP nodes.
1496 signal_handler_install();
1499 pull_t* pull_work = nullptr;
1500 if (thisRankHasDuty(cr, DUTY_PP))
1502 /* Assumes uniform use of the number of OpenMP threads */
1503 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1505 if (inputrec->bPull)
1507 /* Initialize pull code */
1508 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1509 inputrec->fepvals->init_lambda);
1510 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1512 initPullHistory(pull_work, &observablesHistory);
1514 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1516 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1520 std::unique_ptr<EnforcedRotation> enforcedRotation;
1523 /* Initialize enforced rotation code */
1525 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1526 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1529 t_swap* swap = nullptr;
1530 if (inputrec->eSwapCoords != eswapNO)
1532 /* Initialize ion swapping code */
1533 swap = init_swapcoords(fplog, inputrec,
1534 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1535 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1536 oenv, mdrunOptions, startingBehavior);
1539 /* Let makeConstraints know whether we have essential dynamics constraints. */
1540 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1541 *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1543 /* Energy terms and groups */
1544 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1545 inputrec->fepvals->n_lambda);
1547 /* Kinetic energy data */
1548 gmx_ekindata_t ekind;
1549 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1551 /* Set up interactive MD (IMD) */
1553 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1554 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1555 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1557 if (DOMAINDECOMP(cr))
1559 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1560 /* This call is not included in init_domain_decomposition mainly
1561 * because fr->cginfo_mb is set later.
1563 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1564 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1567 const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
1568 false, inputrec, doRerun, vsite.get(), ms, replExParams, fcd,
1569 static_cast<int>(filenames.size()), filenames.data(), &observablesHistory, membed);
1571 const bool useModularSimulator = inputIsCompatibleWithModularSimulator
1572 && !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
1574 // TODO This is not the right place to manage the lifetime of
1575 // this data structure, but currently it's the easiest way to
1577 MdrunScheduleWorkload runScheduleWork;
1578 // Also populates the simulation constant workload description.
1579 runScheduleWork.simulationWork = createSimulationWorkload(
1580 useGpuForNonbonded, pmeRunMode, useGpuForBonded, useGpuForUpdate,
1581 devFlags.enableGpuBufferOps, devFlags.enableGpuHaloExchange,
1582 devFlags.enableGpuPmePPComm, haveEwaldSurfaceContribution(*inputrec));
1584 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1585 if (gpusWereDetected
1586 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1587 || runScheduleWork.simulationWork.useGpuBufferOps))
1589 const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1590 const void* localStream =
1591 fr->nbv->gpu_nbv != nullptr
1592 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1594 const void* nonLocalStream =
1595 fr->nbv->gpu_nbv != nullptr
1596 ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1598 const void* deviceContext = pme_gpu_get_device_context(fr->pmedata);
1599 const int paddingSize = pme_gpu_get_padding_size(fr->pmedata);
1600 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1601 ? GpuApiCallBehavior::Async
1602 : GpuApiCallBehavior::Sync;
1604 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1605 pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
1606 fr->stateGpu = stateGpu.get();
1609 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1610 SimulatorBuilder simulatorBuilder;
1612 // build and run simulator object based on user-input
1613 auto simulator = simulatorBuilder.build(
1614 inputIsCompatibleWithModularSimulator, fplog, cr, ms, mdlog,
1615 static_cast<int>(filenames.size()), filenames.data(), oenv, mdrunOptions,
1616 startingBehavior, vsite.get(), constr.get(),
1617 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1618 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1619 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1620 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1621 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1624 if (fr->pmePpCommGpu)
1626 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1627 fr->pmePpCommGpu.reset();
1630 if (inputrec->bPull)
1632 finish_pull(pull_work);
1634 finish_swapcoords(swap);
1638 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1640 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1641 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1644 wallcycle_stop(wcycle, ewcRUN);
1646 /* Finish up, write some stuff
1647 * if rerunMD, don't write last frame again
1649 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1650 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1652 // clean up cycle counter
1653 wallcycle_destroy(wcycle);
1658 gmx_pme_destroy(pmedata);
1662 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1663 // before we destroy the GPU context(s) in free_gpu().
1664 // Pinned buffers are associated with contexts in CUDA.
1665 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1666 mdAtoms.reset(nullptr);
1667 globalState.reset(nullptr);
1668 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1669 gpuBonded.reset(nullptr);
1670 /* Free pinned buffers in *fr */
1674 if (hwinfo->gpu_info.n_dev > 0)
1676 /* stop the GPU profiler (only CUDA) */
1680 /* With tMPI we need to wait for all ranks to finish deallocation before
1681 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1684 * This is not a concern in OpenCL where we use one context per rank which
1685 * is freed in nbnxn_gpu_free().
1687 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1688 * but it is easier and more futureproof to call it on the whole node.
1690 * Note that this function needs to be called even if GPUs are not used
1691 * in this run because the PME ranks have no knowledge of whether GPUs
1692 * are used or not, but all ranks need to enter the barrier below.
1693 * \todo Remove this physical node barrier after making sure
1694 * that it's not needed anymore (with a shared GPU run).
1698 physicalNodeComm.barrier();
1701 free_gpu(nonbondedDeviceInfo);
1702 free_gpu(pmeDeviceInfo);
1707 free_membed(membed);
1710 /* Does what it says */
1711 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1712 walltime_accounting_destroy(walltime_accounting);
1714 // Ensure log file content is written
1717 gmx_fio_flush(logFileHandle);
1720 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1721 * exceptions were enabled before function was called. */
1724 gmx_fedisableexcept();
1727 auto rc = static_cast<int>(gmx_get_stop_condition());
1730 /* we need to join all threads. The sub-threads join when they
1731 exit this function, but the master thread needs to be told to
1733 if (PAR(cr) && MASTER(cr))
1741 Mdrunner::~Mdrunner()
1743 // Clean up of the Manager.
1744 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1745 // but okay as long as threads synchronize some time before adding or accessing
1746 // a new set of restraints.
1747 if (restraintManager_)
1749 restraintManager_->clear();
1750 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1751 "restraints added during runner life time should be cleared at runner "
1756 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1758 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1759 // Not sure if this should be logged through the md logger or something else,
1760 // but it is helpful to have some sort of INFO level message sent somewhere.
1761 // std::cout << "Registering restraint named " << name << std::endl;
1763 // When multiple restraints are used, it may be wasteful to register them separately.
1764 // Maybe instead register an entire Restraint Manager as a force provider.
1765 restraintManager_->addToSpec(std::move(puller), name);
1768 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1770 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1772 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1773 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1775 class Mdrunner::BuilderImplementation
1778 BuilderImplementation() = delete;
1779 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1780 ~BuilderImplementation();
1782 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1783 real forceWarningThreshold,
1784 StartingBehavior startingBehavior);
1786 void addDomdec(const DomdecOptions& options);
1788 void addVerletList(int nstlist);
1790 void addReplicaExchange(const ReplicaExchangeParameters& params);
1792 void addNonBonded(const char* nbpu_opt);
1794 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1796 void addBondedTaskAssignment(const char* bonded_opt);
1798 void addUpdateTaskAssignment(const char* update_opt);
1800 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1802 void addFilenames(ArrayRef<const t_filenm> filenames);
1804 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1806 void addLogFile(t_fileio* logFileHandle);
1808 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1813 // Default parameters copied from runner.h
1814 // \todo Clarify source(s) of default parameters.
1816 const char* nbpu_opt_ = nullptr;
1817 const char* pme_opt_ = nullptr;
1818 const char* pme_fft_opt_ = nullptr;
1819 const char* bonded_opt_ = nullptr;
1820 const char* update_opt_ = nullptr;
1822 MdrunOptions mdrunOptions_;
1824 DomdecOptions domdecOptions_;
1826 ReplicaExchangeParameters replicaExchangeParameters_;
1828 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1831 //! Multisim communicator handle.
1832 gmx_multisim_t* multiSimulation_;
1834 //! mdrun communicator
1835 MPI_Comm communicator_ = MPI_COMM_NULL;
1837 //! Print a warning if any force is larger than this (in kJ/mol nm).
1838 real forceWarningThreshold_ = -1;
1840 //! Whether the simulation will start afresh, or restart with/without appending.
1841 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1843 //! The modules that comprise the functionality of mdrun.
1844 std::unique_ptr<MDModules> mdModules_;
1846 //! \brief Parallelism information.
1847 gmx_hw_opt_t hardwareOptions_;
1849 //! filename options for simulation.
1850 ArrayRef<const t_filenm> filenames_;
1852 /*! \brief Handle to output environment.
1854 * \todo gmx_output_env_t needs lifetime management.
1856 gmx_output_env_t* outputEnvironment_ = nullptr;
1858 /*! \brief Non-owning handle to MD log file.
1860 * \todo Context should own output facilities for client.
1861 * \todo Improve log file handle management.
1863 * Code managing the FILE* relies on the ability to set it to
1864 * nullptr to check whether the filehandle is valid.
1866 t_fileio* logFileHandle_ = nullptr;
1869 * \brief Builder for simulation stop signal handler.
1871 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1874 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1875 compat::not_null<SimulationContext*> context) :
1876 mdModules_(std::move(mdModules))
1878 communicator_ = context->communicator_;
1879 multiSimulation_ = context->multiSimulation_.get();
1882 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1884 Mdrunner::BuilderImplementation&
1885 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1886 const real forceWarningThreshold,
1887 const StartingBehavior startingBehavior)
1889 mdrunOptions_ = options;
1890 forceWarningThreshold_ = forceWarningThreshold;
1891 startingBehavior_ = startingBehavior;
1895 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1897 domdecOptions_ = options;
1900 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1905 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1907 replicaExchangeParameters_ = params;
1910 Mdrunner Mdrunner::BuilderImplementation::build()
1912 auto newRunner = Mdrunner(std::move(mdModules_));
1914 newRunner.mdrunOptions = mdrunOptions_;
1915 newRunner.pforce = forceWarningThreshold_;
1916 newRunner.startingBehavior = startingBehavior_;
1917 newRunner.domdecOptions = domdecOptions_;
1919 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1920 newRunner.hw_opt = hardwareOptions_;
1922 // No invariant to check. This parameter exists to optionally override other behavior.
1923 newRunner.nstlist_cmdline = nstlist_;
1925 newRunner.replExParams = replicaExchangeParameters_;
1927 newRunner.filenames = filenames_;
1929 newRunner.communicator = communicator_;
1931 // nullptr is a valid value for the multisim handle
1932 newRunner.ms = multiSimulation_;
1934 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1935 // \todo Update sanity checking when output environment has clearly specified invariants.
1936 // Initialization and default values for oenv are not well specified in the current version.
1937 if (outputEnvironment_)
1939 newRunner.oenv = outputEnvironment_;
1943 GMX_THROW(gmx::APIError(
1944 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1947 newRunner.logFileHandle = logFileHandle_;
1951 newRunner.nbpu_opt = nbpu_opt_;
1955 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1958 if (pme_opt_ && pme_fft_opt_)
1960 newRunner.pme_opt = pme_opt_;
1961 newRunner.pme_fft_opt = pme_fft_opt_;
1965 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1970 newRunner.bonded_opt = bonded_opt_;
1974 GMX_THROW(gmx::APIError(
1975 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1980 newRunner.update_opt = update_opt_;
1984 GMX_THROW(gmx::APIError(
1985 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1989 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1991 if (stopHandlerBuilder_)
1993 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1997 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2003 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2005 nbpu_opt_ = nbpu_opt;
2008 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2011 pme_fft_opt_ = pme_fft_opt;
2014 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2016 bonded_opt_ = bonded_opt;
2019 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2021 update_opt_ = update_opt;
2024 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2026 hardwareOptions_ = hardwareOptions;
2029 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2031 filenames_ = filenames;
2034 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2036 outputEnvironment_ = outputEnvironment;
2039 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2041 logFileHandle_ = logFileHandle;
2044 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2046 stopHandlerBuilder_ = std::move(builder);
2049 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2050 compat::not_null<SimulationContext*> context) :
2051 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2055 MdrunnerBuilder::~MdrunnerBuilder() = default;
2057 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2058 real forceWarningThreshold,
2059 const StartingBehavior startingBehavior)
2061 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2065 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2067 impl_->addDomdec(options);
2071 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2073 impl_->addVerletList(nstlist);
2077 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2079 impl_->addReplicaExchange(params);
2083 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2085 impl_->addNonBonded(nbpu_opt);
2089 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2091 // The builder method may become more general in the future, but in this version,
2092 // parameters for PME electrostatics are both required and the only parameters
2094 if (pme_opt && pme_fft_opt)
2096 impl_->addPME(pme_opt, pme_fft_opt);
2101 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2106 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2108 impl_->addBondedTaskAssignment(bonded_opt);
2112 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2114 impl_->addUpdateTaskAssignment(update_opt);
2118 Mdrunner MdrunnerBuilder::build()
2120 return impl_->build();
2123 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2125 impl_->addHardwareOptions(hardwareOptions);
2129 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2131 impl_->addFilenames(filenames);
2135 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2137 impl_->addOutputEnvironment(outputEnvironment);
2141 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2143 impl_->addLogFile(logFileHandle);
2147 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2149 impl_->addStopHandlerBuilder(std::move(builder));
2153 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2155 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;