Use CenterOfMassPulling
[gromacs.git] / src / gromacs / mdrun / runner.cpp
blobe0c24925fa59e3e5a3f8e435181b9cc950fc2059
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \internal \file
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
44 #include "gmxpre.h"
46 #include "runner.h"
48 #include "config.h"
50 #include <cassert>
51 #include <cinttypes>
52 #include <csignal>
53 #include <cstdlib>
54 #include <cstring>
56 #include <algorithm>
57 #include <memory>
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/device_context.h"
77 #include "gromacs/gpu_utils/device_stream_manager.h"
78 #include "gromacs/gpu_utils/gpu_utils.h"
79 #include "gromacs/hardware/cpuinfo.h"
80 #include "gromacs/hardware/detecthardware.h"
81 #include "gromacs/hardware/printhardware.h"
82 #include "gromacs/imd/imd.h"
83 #include "gromacs/listed_forces/disre.h"
84 #include "gromacs/listed_forces/gpubonded.h"
85 #include "gromacs/listed_forces/listed_forces.h"
86 #include "gromacs/listed_forces/orires.h"
87 #include "gromacs/math/functions.h"
88 #include "gromacs/math/utilities.h"
89 #include "gromacs/math/vec.h"
90 #include "gromacs/mdlib/boxdeformation.h"
91 #include "gromacs/mdlib/broadcaststructs.h"
92 #include "gromacs/mdlib/calc_verletbuf.h"
93 #include "gromacs/mdlib/dispersioncorrection.h"
94 #include "gromacs/mdlib/enerdata_utils.h"
95 #include "gromacs/mdlib/force.h"
96 #include "gromacs/mdlib/forcerec.h"
97 #include "gromacs/mdlib/gmx_omp_nthreads.h"
98 #include "gromacs/mdlib/makeconstraints.h"
99 #include "gromacs/mdlib/md_support.h"
100 #include "gromacs/mdlib/mdatoms.h"
101 #include "gromacs/mdlib/sighandler.h"
102 #include "gromacs/mdlib/stophandler.h"
103 #include "gromacs/mdlib/tgroup.h"
104 #include "gromacs/mdlib/updategroups.h"
105 #include "gromacs/mdlib/vsite.h"
106 #include "gromacs/mdrun/mdmodules.h"
107 #include "gromacs/mdrun/simulationcontext.h"
108 #include "gromacs/mdrunutility/handlerestart.h"
109 #include "gromacs/mdrunutility/logging.h"
110 #include "gromacs/mdrunutility/multisim.h"
111 #include "gromacs/mdrunutility/printtime.h"
112 #include "gromacs/mdrunutility/threadaffinity.h"
113 #include "gromacs/mdtypes/commrec.h"
114 #include "gromacs/mdtypes/enerdata.h"
115 #include "gromacs/mdtypes/fcdata.h"
116 #include "gromacs/mdtypes/forcerec.h"
117 #include "gromacs/mdtypes/group.h"
118 #include "gromacs/mdtypes/inputrec.h"
119 #include "gromacs/mdtypes/interaction_const.h"
120 #include "gromacs/mdtypes/md_enums.h"
121 #include "gromacs/mdtypes/mdatom.h"
122 #include "gromacs/mdtypes/mdrunoptions.h"
123 #include "gromacs/mdtypes/observableshistory.h"
124 #include "gromacs/mdtypes/simulation_workload.h"
125 #include "gromacs/mdtypes/state.h"
126 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
127 #include "gromacs/modularsimulator/modularsimulator.h"
128 #include "gromacs/nbnxm/gpu_data_mgmt.h"
129 #include "gromacs/nbnxm/nbnxm.h"
130 #include "gromacs/nbnxm/pairlist_tuning.h"
131 #include "gromacs/pbcutil/pbc.h"
132 #include "gromacs/pulling/output.h"
133 #include "gromacs/pulling/pull.h"
134 #include "gromacs/pulling/pull_rotation.h"
135 #include "gromacs/restraint/manager.h"
136 #include "gromacs/restraint/restraintmdmodule.h"
137 #include "gromacs/restraint/restraintpotential.h"
138 #include "gromacs/swap/swapcoords.h"
139 #include "gromacs/taskassignment/decidegpuusage.h"
140 #include "gromacs/taskassignment/decidesimulationworkload.h"
141 #include "gromacs/taskassignment/resourcedivision.h"
142 #include "gromacs/taskassignment/taskassignment.h"
143 #include "gromacs/taskassignment/usergpuids.h"
144 #include "gromacs/timing/gpu_timing.h"
145 #include "gromacs/timing/wallcycle.h"
146 #include "gromacs/timing/wallcyclereporting.h"
147 #include "gromacs/topology/mtop_util.h"
148 #include "gromacs/trajectory/trajectoryframe.h"
149 #include "gromacs/utility/basenetwork.h"
150 #include "gromacs/utility/cstringutil.h"
151 #include "gromacs/utility/exceptions.h"
152 #include "gromacs/utility/fatalerror.h"
153 #include "gromacs/utility/filestream.h"
154 #include "gromacs/utility/gmxassert.h"
155 #include "gromacs/utility/gmxmpi.h"
156 #include "gromacs/utility/keyvaluetree.h"
157 #include "gromacs/utility/logger.h"
158 #include "gromacs/utility/loggerbuilder.h"
159 #include "gromacs/utility/mdmodulenotification.h"
160 #include "gromacs/utility/physicalnodecommunicator.h"
161 #include "gromacs/utility/pleasecite.h"
162 #include "gromacs/utility/programcontext.h"
163 #include "gromacs/utility/smalloc.h"
164 #include "gromacs/utility/stringutil.h"
166 #include "isimulator.h"
167 #include "membedholder.h"
168 #include "replicaexchange.h"
169 #include "simulatorbuilder.h"
171 #if GMX_FAHCORE
172 # include "corewrap.h"
173 #endif
175 namespace gmx
179 /*! \brief Manage any development feature flag variables encountered
181 * The use of dev features indicated by environment variables is
182 * logged in order to ensure that runs with such features enabled can
183 * be identified from their log and standard output. Any cross
184 * dependencies are also checked, and if unsatisfied, a fatal error
185 * issued.
187 * Note that some development features overrides are applied already here:
188 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
190 * \param[in] mdlog Logger object.
191 * \param[in] useGpuForNonbonded True if the nonbonded task is offloaded in this run.
192 * \param[in] pmeRunMode The PME run mode for this run
193 * \returns The object populated with development feature flags.
195 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
196 const bool useGpuForNonbonded,
197 const PmeRunMode pmeRunMode)
199 DevelopmentFeatureFlags devFlags;
201 // Some builds of GCC 5 give false positive warnings that these
202 // getenv results are ignored when clearly they are used.
203 #pragma GCC diagnostic push
204 #pragma GCC diagnostic ignored "-Wunused-result"
205 devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
206 && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
207 devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
208 devFlags.enableGpuHaloExchange =
209 (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
210 devFlags.enableGpuPmePPComm =
211 (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
212 #pragma GCC diagnostic pop
214 if (devFlags.enableGpuBufferOps)
216 GMX_LOG(mdlog.warning)
217 .asParagraph()
218 .appendTextFormatted(
219 "This run uses the 'GPU buffer ops' feature, enabled by the "
220 "GMX_USE_GPU_BUFFER_OPS environment variable.");
223 if (devFlags.forceGpuUpdateDefault)
225 GMX_LOG(mdlog.warning)
226 .asParagraph()
227 .appendTextFormatted(
228 "This run will default to '-update gpu' as requested by the "
229 "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
230 "decomposition lacks substantial testing and should be used with caution.");
233 if (devFlags.enableGpuHaloExchange)
235 if (useGpuForNonbonded)
237 if (!devFlags.enableGpuBufferOps)
239 GMX_LOG(mdlog.warning)
240 .asParagraph()
241 .appendTextFormatted(
242 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
243 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
244 devFlags.enableGpuBufferOps = true;
246 GMX_LOG(mdlog.warning)
247 .asParagraph()
248 .appendTextFormatted(
249 "This run has requested the 'GPU halo exchange' feature, enabled by "
250 "the "
251 "GMX_GPU_DD_COMMS environment variable.");
253 else
255 GMX_LOG(mdlog.warning)
256 .asParagraph()
257 .appendTextFormatted(
258 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
259 "halo exchange' feature will not be enabled as nonbonded interactions "
260 "are not offloaded.");
261 devFlags.enableGpuHaloExchange = false;
265 if (devFlags.enableGpuPmePPComm)
267 if (pmeRunMode == PmeRunMode::GPU)
269 GMX_LOG(mdlog.warning)
270 .asParagraph()
271 .appendTextFormatted(
272 "This run uses the 'GPU PME-PP communications' feature, enabled "
273 "by the GMX_GPU_PME_PP_COMMS environment variable.");
275 else
277 std::string clarification;
278 if (pmeRunMode == PmeRunMode::Mixed)
280 clarification =
281 "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
282 "mode).";
284 else
286 clarification = "PME is not offloaded to the GPU.";
288 GMX_LOG(mdlog.warning)
289 .asParagraph()
290 .appendText(
291 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
292 "'GPU PME-PP communications' feature was not enabled as "
293 + clarification);
294 devFlags.enableGpuPmePPComm = false;
298 return devFlags;
301 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
303 * Used to ensure that the master thread does not modify mdrunner during copy
304 * on the spawned threads. */
305 static void threadMpiMdrunnerAccessBarrier()
307 #if GMX_THREAD_MPI
308 MPI_Barrier(MPI_COMM_WORLD);
309 #endif
312 Mdrunner Mdrunner::cloneOnSpawnedThread() const
314 auto newRunner = Mdrunner(std::make_unique<MDModules>());
316 // All runners in the same process share a restraint manager resource because it is
317 // part of the interface to the client code, which is associated only with the
318 // original thread. Handles to the same resources can be obtained by copy.
320 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
323 // Copy members of master runner.
324 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
325 // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
326 newRunner.hw_opt = hw_opt;
327 newRunner.filenames = filenames;
329 newRunner.oenv = oenv;
330 newRunner.mdrunOptions = mdrunOptions;
331 newRunner.domdecOptions = domdecOptions;
332 newRunner.nbpu_opt = nbpu_opt;
333 newRunner.pme_opt = pme_opt;
334 newRunner.pme_fft_opt = pme_fft_opt;
335 newRunner.bonded_opt = bonded_opt;
336 newRunner.update_opt = update_opt;
337 newRunner.nstlist_cmdline = nstlist_cmdline;
338 newRunner.replExParams = replExParams;
339 newRunner.pforce = pforce;
340 // Give the spawned thread the newly created valid communicator
341 // for the simulation.
342 newRunner.communicator = MPI_COMM_WORLD;
343 newRunner.ms = ms;
344 newRunner.startingBehavior = startingBehavior;
345 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
347 threadMpiMdrunnerAccessBarrier();
349 return newRunner;
352 /*! \brief The callback used for running on spawned threads.
354 * Obtains the pointer to the master mdrunner object from the one
355 * argument permitted to the thread-launch API call, copies it to make
356 * a new runner for this thread, reinitializes necessary data, and
357 * proceeds to the simulation. */
358 static void mdrunner_start_fn(const void* arg)
362 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
363 /* copy the arg list to make sure that it's thread-local. This
364 doesn't copy pointed-to items, of course; fnm, cr and fplog
365 are reset in the call below, all others should be const. */
366 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
367 mdrunner.mdrunner();
369 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
373 void Mdrunner::spawnThreads(int numThreadsToLaunch)
375 #if GMX_THREAD_MPI
376 /* now spawn new threads that start mdrunner_start_fn(), while
377 the main thread returns. Thread affinity is handled later. */
378 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
379 static_cast<const void*>(this))
380 != TMPI_SUCCESS)
382 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
385 // Give the master thread the newly created valid communicator for
386 // the simulation.
387 communicator = MPI_COMM_WORLD;
388 threadMpiMdrunnerAccessBarrier();
389 #else
390 GMX_UNUSED_VALUE(numThreadsToLaunch);
391 GMX_UNUSED_VALUE(mdrunner_start_fn);
392 #endif
395 } // namespace gmx
397 /*! \brief Initialize variables for Verlet scheme simulation */
398 static void prepare_verlet_scheme(FILE* fplog,
399 t_commrec* cr,
400 t_inputrec* ir,
401 int nstlist_cmdline,
402 const gmx_mtop_t* mtop,
403 const matrix box,
404 bool makeGpuPairList,
405 const gmx::CpuInfo& cpuinfo)
407 // We checked the cut-offs in grompp, but double-check here.
408 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
409 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
411 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
412 "With Verlet lists and PME we should have rcoulomb>=rvdw");
414 else
416 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
417 "With Verlet lists and no PME rcoulomb and rvdw should be identical");
419 /* For NVE simulations, we will retain the initial list buffer */
420 if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
422 /* Update the Verlet buffer size for the current run setup */
424 /* Here we assume SIMD-enabled kernels are being used. But as currently
425 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
426 * and 4x2 gives a larger buffer than 4x4, this is ok.
428 ListSetupType listType =
429 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
430 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
432 const real rlist_new =
433 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
435 if (rlist_new != ir->rlist)
437 if (fplog != nullptr)
439 fprintf(fplog,
440 "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
441 ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
443 ir->rlist = rlist_new;
447 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
449 gmx_fatal(FARGS, "Can not set nstlist without %s",
450 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
453 if (EI_DYNAMICS(ir->eI))
455 /* Set or try nstlist values */
456 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
460 /*! \brief Override the nslist value in inputrec
462 * with value passed on the command line (if any)
464 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
466 assert(ir);
468 /* override with anything else than the default -2 */
469 if (nsteps_cmdline > -2)
471 char sbuf_steps[STEPSTRSIZE];
472 char sbuf_msg[STRLEN];
474 ir->nsteps = nsteps_cmdline;
475 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
477 sprintf(sbuf_msg,
478 "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
479 gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
481 else
483 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
484 gmx_step_str(nsteps_cmdline, sbuf_steps));
487 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
489 else if (nsteps_cmdline < -2)
491 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
493 /* Do nothing if nsteps_cmdline == -2 */
496 namespace gmx
499 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
501 * If not, and if a warning may be issued, logs a warning about
502 * falling back to CPU code. With thread-MPI, only the first
503 * call to this function should have \c issueWarning true. */
504 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
506 bool gpuIsUseful = true;
507 std::string warning;
509 if (ir.opts.ngener - ir.nwall > 1)
511 /* The GPU code does not support more than one energy group.
512 * If the user requested GPUs explicitly, a fatal error is given later.
514 gpuIsUseful = false;
515 warning =
516 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
517 "For better performance, run on the GPU without energy groups and then do "
518 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
521 if (EI_TPI(ir.eI))
523 gpuIsUseful = false;
524 warning = "TPI is not implemented for GPUs.";
527 if (!gpuIsUseful && issueWarning)
529 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
532 return gpuIsUseful;
535 //! Initializes the logger for mdrun.
536 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
538 gmx::LoggerBuilder builder;
539 if (fplog != nullptr)
541 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
543 if (isSimulationMasterRank)
545 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
547 return builder.build();
550 //! Make a TaskTarget from an mdrun argument string.
551 static TaskTarget findTaskTarget(const char* optionString)
553 TaskTarget returnValue = TaskTarget::Auto;
555 if (strncmp(optionString, "auto", 3) == 0)
557 returnValue = TaskTarget::Auto;
559 else if (strncmp(optionString, "cpu", 3) == 0)
561 returnValue = TaskTarget::Cpu;
563 else if (strncmp(optionString, "gpu", 3) == 0)
565 returnValue = TaskTarget::Gpu;
567 else
569 GMX_ASSERT(false, "Option string should have been checked for sanity already");
572 return returnValue;
575 //! Finish run, aggregate data to print performance info.
576 static void finish_run(FILE* fplog,
577 const gmx::MDLogger& mdlog,
578 const t_commrec* cr,
579 const t_inputrec* inputrec,
580 t_nrnb nrnb[],
581 gmx_wallcycle_t wcycle,
582 gmx_walltime_accounting_t walltime_accounting,
583 nonbonded_verlet_t* nbv,
584 const gmx_pme_t* pme,
585 gmx_bool bWriteStat)
587 double delta_t = 0;
588 double nbfs = 0, mflop = 0;
589 double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
590 elapsed_time_over_all_threads_over_all_ranks;
591 /* Control whether it is valid to print a report. Only the
592 simulation master may print, but it should not do so if the run
593 terminated e.g. before a scheduled reset step. This is
594 complicated by the fact that PME ranks are unaware of the
595 reason why they were sent a pmerecvqxFINISH. To avoid
596 communication deadlocks, we always do the communication for the
597 report, even if we've decided not to write the report, because
598 how long it takes to finish the run is not important when we've
599 decided not to report on the simulation performance.
601 Further, we only report performance for dynamical integrators,
602 because those are the only ones for which we plan to
603 consider doing any optimizations. */
604 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
606 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
608 GMX_LOG(mdlog.warning)
609 .asParagraph()
610 .appendText("Simulation ended prematurely, no performance report will be written.");
611 printReport = false;
614 t_nrnb* nrnb_tot;
615 std::unique_ptr<t_nrnb> nrnbTotalStorage;
616 if (cr->nnodes > 1)
618 nrnbTotalStorage = std::make_unique<t_nrnb>();
619 nrnb_tot = nrnbTotalStorage.get();
620 #if GMX_MPI
621 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
622 #endif
624 else
626 nrnb_tot = nrnb;
629 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
630 elapsed_time_over_all_threads =
631 walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
632 if (cr->nnodes > 1)
634 #if GMX_MPI
635 /* reduce elapsed_time over all MPI ranks in the current simulation */
636 MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
637 cr->mpi_comm_mysim);
638 elapsed_time_over_all_ranks /= cr->nnodes;
639 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
640 * current simulation. */
641 MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
642 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
643 #endif
645 else
647 elapsed_time_over_all_ranks = elapsed_time;
648 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
651 if (printReport)
653 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
656 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
658 print_dd_statistics(cr, inputrec, fplog);
661 /* TODO Move the responsibility for any scaling by thread counts
662 * to the code that handled the thread region, so that there's a
663 * mechanism to keep cycle counting working during the transition
664 * to task parallelism. */
665 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
666 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
667 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
668 nthreads_pp, nthreads_pme);
669 auto cycle_sum(wallcycle_sum(cr, wcycle));
671 if (printReport)
673 auto nbnxn_gpu_timings =
674 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
675 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
677 if (pme_gpu_task_enabled(pme))
679 pme_gpu_get_timings(pme, &pme_gpu_timings);
681 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
682 elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
683 &pme_gpu_timings);
685 if (EI_DYNAMICS(inputrec->eI))
687 delta_t = inputrec->delta_t;
690 if (fplog)
692 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
693 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
694 delta_t, nbfs, mflop);
696 if (bWriteStat)
698 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
699 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
700 delta_t, nbfs, mflop);
705 int Mdrunner::mdrunner()
707 matrix box;
708 t_forcerec* fr = nullptr;
709 real ewaldcoeff_q = 0;
710 real ewaldcoeff_lj = 0;
711 int nChargePerturbed = -1, nTypePerturbed = 0;
712 gmx_wallcycle_t wcycle;
713 gmx_walltime_accounting_t walltime_accounting = nullptr;
714 MembedHolder membedHolder(filenames.size(), filenames.data());
715 gmx_hw_info_t* hwinfo = nullptr;
717 /* CAUTION: threads may be started later on in this function, so
718 cr doesn't reflect the final parallel state right now */
719 gmx_mtop_t mtop;
721 /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
722 const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
723 const bool doRerun = mdrunOptions.rerun;
725 // Handle task-assignment related user options.
726 EmulateGpuNonbonded emulateGpuNonbonded =
727 (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
729 std::vector<int> userGpuTaskAssignment;
732 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
734 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
735 auto nonbondedTarget = findTaskTarget(nbpu_opt);
736 auto pmeTarget = findTaskTarget(pme_opt);
737 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
738 auto bondedTarget = findTaskTarget(bonded_opt);
739 auto updateTarget = findTaskTarget(update_opt);
741 FILE* fplog = nullptr;
742 // If we are appending, we don't write log output because we need
743 // to check that the old log file matches what the checkpoint file
744 // expects. Otherwise, we should start to write log output now if
745 // there is a file ready for it.
746 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
748 fplog = gmx_fio_getfp(logFileHandle);
750 const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
751 gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
752 gmx::MDLogger mdlog(logOwner.logger());
754 // TODO The thread-MPI master rank makes a working
755 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
756 // after the threads have been launched. This works because no use
757 // is made of that communicator until after the execution paths
758 // have rejoined. But it is likely that we can improve the way
759 // this is expressed, e.g. by expressly running detection only the
760 // master rank for thread-MPI, rather than relying on the mutex
761 // and reference count.
762 PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
763 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
765 gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
767 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
769 // Print citation requests after all software/hardware printing
770 pleaseCiteGromacs(fplog);
772 // TODO Replace this by unique_ptr once t_inputrec is C++
773 t_inputrec inputrecInstance;
774 t_inputrec* inputrec = nullptr;
775 std::unique_ptr<t_state> globalState;
777 auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
779 if (isSimulationMasterRank)
781 /* Only the master rank has the global state */
782 globalState = std::make_unique<t_state>();
784 /* Read (nearly) all data required for the simulation
785 * and keep the partly serialized tpr contents to send to other ranks later
787 *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
788 &inputrecInstance, globalState.get(), &mtop);
789 inputrec = &inputrecInstance;
792 /* Check and update the hardware options for internal consistency */
793 checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
794 inputrec);
796 if (GMX_THREAD_MPI && isSimulationMasterRank)
798 bool useGpuForNonbonded = false;
799 bool useGpuForPme = false;
802 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
804 // If the user specified the number of ranks, then we must
805 // respect that, but in default mode, we need to allow for
806 // the number of GPUs to choose the number of ranks.
807 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
808 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi(
809 nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
810 canUseGpuForNonbonded,
811 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
812 hw_opt.nthreads_tmpi);
813 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
814 useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
815 *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
817 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
819 /* Determine how many thread-MPI ranks to start.
821 * TODO Over-writing the user-supplied value here does
822 * prevent any possible subsequent checks from working
823 * correctly. */
824 hw_opt.nthreads_tmpi =
825 get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded, useGpuForPme,
826 inputrec, &mtop, mdlog, membedHolder.doMembed());
828 // Now start the threads for thread MPI.
829 spawnThreads(hw_opt.nthreads_tmpi);
830 // The spawned threads enter mdrunner() and execution of
831 // master and spawned threads joins at the end of this block.
832 physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
835 GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
836 CommrecHandle crHandle = init_commrec(communicator, ms);
837 t_commrec* cr = crHandle.get();
838 GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
840 if (PAR(cr))
842 /* now broadcast everything to the non-master nodes/threads: */
843 if (!isSimulationMasterRank)
845 inputrec = &inputrecInstance;
847 init_parallel(cr->mpi_comm_mygroup, MASTER(cr), inputrec, &mtop, partialDeserializedTpr.get());
849 GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
850 partialDeserializedTpr.reset(nullptr);
852 // Now the number of ranks is known to all ranks, and each knows
853 // the inputrec read by the master rank. The ranks can now all run
854 // the task-deciding functions and will agree on the result
855 // without needing to communicate.
856 const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
858 // Note that these variables describe only their own node.
860 // Note that when bonded interactions run on a GPU they always run
861 // alongside a nonbonded task, so do not influence task assignment
862 // even though they affect the force calculation workload.
863 bool useGpuForNonbonded = false;
864 bool useGpuForPme = false;
865 bool useGpuForBonded = false;
866 bool useGpuForUpdate = false;
867 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
870 // It's possible that there are different numbers of GPUs on
871 // different nodes, which is the user's responsibility to
872 // handle. If unsuitable, we will notice that during task
873 // assignment.
874 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
875 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(
876 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
877 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
878 useGpuForPme = decideWhetherToUseGpusForPme(
879 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
880 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
881 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
882 && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
883 useGpuForBonded = decideWhetherToUseGpusForBonded(
884 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
885 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
886 domdecOptions.numPmeRanks, gpusWereDetected);
888 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
890 const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
892 // Initialize development feature flags that enabled by environment variable
893 // and report those features that are enabled.
894 const DevelopmentFeatureFlags devFlags =
895 manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
897 const bool useModularSimulator =
898 checkUseModularSimulator(false, inputrec, doRerun, mtop, ms, replExParams, nullptr,
899 doEssentialDynamics, membedHolder.doMembed());
901 // Build restraints.
902 // TODO: hide restraint implementation details from Mdrunner.
903 // There is nothing unique about restraints at this point as far as the
904 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
905 // factory functions from the SimulationContext on which to call mdModules_->add().
906 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
907 for (auto&& restraint : restraintManager_->getRestraints())
909 auto module = RestraintMDModule::create(restraint, restraint->sites());
910 mdModules_->add(std::move(module));
913 // TODO: Error handling
914 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
915 // now that the MdModules know their options, they know which callbacks to sign up to
916 mdModules_->subscribeToSimulationSetupNotifications();
917 const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
919 if (inputrec->internalParameters != nullptr)
921 mdModulesNotifier.notify(*inputrec->internalParameters);
924 if (fplog != nullptr)
926 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
927 fprintf(fplog, "\n");
930 if (SIMMASTER(cr))
932 /* In rerun, set velocities to zero if present */
933 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
935 // rerun does not use velocities
936 GMX_LOG(mdlog.info)
937 .asParagraph()
938 .appendText(
939 "Rerun trajectory contains velocities. Rerun does only evaluate "
940 "potential energy and forces. The velocities will be ignored.");
941 for (int i = 0; i < globalState->natoms; i++)
943 clear_rvec(globalState->v[i]);
945 globalState->flags &= ~(1 << estV);
948 /* now make sure the state is initialized and propagated */
949 set_state_entries(globalState.get(), inputrec, useModularSimulator);
952 /* NM and TPI parallelize over force/energy calculations, not atoms,
953 * so we need to initialize and broadcast the global state.
955 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
957 if (!MASTER(cr))
959 globalState = std::make_unique<t_state>();
961 broadcastStateWithoutDynamics(cr->mpi_comm_mygroup, DOMAINDECOMP(cr), PAR(cr), globalState.get());
964 /* A parallel command line option consistency check that we can
965 only do after any threads have started. */
966 if (!PAR(cr)
967 && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
968 || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
970 gmx_fatal(FARGS,
971 "The -dd or -npme option request a parallel simulation, "
972 #if !GMX_MPI
973 "but %s was compiled without threads or MPI enabled",
974 output_env_get_program_display_name(oenv));
975 #elif GMX_THREAD_MPI
976 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
977 #else
978 "but %s was not started through mpirun/mpiexec or only one rank was requested "
979 "through mpirun/mpiexec",
980 output_env_get_program_display_name(oenv));
981 #endif
984 if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
986 gmx_fatal(FARGS,
987 "The .mdp file specified an energy mininization or normal mode algorithm, and "
988 "these are not compatible with mdrun -rerun");
991 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
993 if (domdecOptions.numPmeRanks > 0)
995 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
996 "PME-only ranks are requested, but the system does not use PME "
997 "for electrostatics or LJ");
1000 domdecOptions.numPmeRanks = 0;
1003 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1005 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1006 * improve performance with many threads per GPU, since our OpenMP
1007 * scaling is bad, but it's difficult to automate the setup.
1009 domdecOptions.numPmeRanks = 0;
1011 if (useGpuForPme)
1013 if (domdecOptions.numPmeRanks < 0)
1015 domdecOptions.numPmeRanks = 0;
1016 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1018 else
1020 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1021 "PME GPU decomposition is not supported");
1025 #if GMX_FAHCORE
1026 if (MASTER(cr))
1028 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1030 #endif
1032 /* NMR restraints must be initialized before load_checkpoint,
1033 * since with time averaging the history is added to t_state.
1034 * For proper consistency check we therefore need to extend
1035 * t_state here.
1036 * So the PME-only nodes (if present) will also initialize
1037 * the distance restraints.
1040 /* This needs to be called before read_checkpoint to extend the state */
1041 t_disresdata* disresdata;
1042 snew(disresdata, 1);
1043 init_disres(fplog, &mtop, inputrec, DisResRunMode::MDRun, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1044 PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, disresdata,
1045 globalState.get(), replExParams.exchangeInterval > 0);
1047 t_oriresdata* oriresdata;
1048 snew(oriresdata, 1);
1049 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), oriresdata);
1051 auto deform = prepareBoxDeformation(globalState->box, MASTER(cr) ? DDRole::Master : DDRole::Agent,
1052 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1053 cr->mpi_comm_mygroup, *inputrec);
1055 ObservablesHistory observablesHistory = {};
1057 if (startingBehavior != StartingBehavior::NewSimulation)
1059 /* Check if checkpoint file exists before doing continuation.
1060 * This way we can use identical input options for the first and subsequent runs...
1062 if (mdrunOptions.numStepsCommandline > -2)
1064 /* Temporarily set the number of steps to unlimited to avoid
1065 * triggering the nsteps check in load_checkpoint().
1066 * This hack will go away soon when the -nsteps option is removed.
1068 inputrec->nsteps = -1;
1071 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1072 logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1073 &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1075 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1077 // Now we can start normal logging to the truncated log file.
1078 fplog = gmx_fio_getfp(logFileHandle);
1079 prepareLogAppending(fplog);
1080 logOwner = buildLogger(fplog, MASTER(cr));
1081 mdlog = logOwner.logger();
1085 if (mdrunOptions.numStepsCommandline > -2)
1087 GMX_LOG(mdlog.info)
1088 .asParagraph()
1089 .appendText(
1090 "The -nsteps functionality is deprecated, and may be removed in a future "
1091 "version. "
1092 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1093 "file field.");
1095 /* override nsteps with value set on the commandline */
1096 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1098 if (SIMMASTER(cr))
1100 copy_mat(globalState->box, box);
1103 if (PAR(cr))
1105 gmx_bcast(sizeof(box), box, cr->mpi_comm_mygroup);
1108 if (inputrec->cutoff_scheme != ecutsVERLET)
1110 gmx_fatal(FARGS,
1111 "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1112 "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1114 /* Update rlist and nstlist. */
1115 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1116 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1117 *hwinfo->cpuInfo);
1119 const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1120 // This builder is necessary while we have multi-part construction
1121 // of DD. Before DD is constructed, we use the existence of
1122 // the builder object to indicate that further construction of DD
1123 // is needed.
1124 std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1125 if (useDomainDecomposition)
1127 ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1128 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1129 positionsFromStatePointer(globalState.get()));
1131 else
1133 /* PME, if used, is done on all nodes with 1D decomposition */
1134 cr->npmenodes = 0;
1135 cr->duty = (DUTY_PP | DUTY_PME);
1137 if (inputrec->pbcType == PbcType::Screw)
1139 gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1143 // Produce the task assignment for this rank - done after DD is constructed
1144 GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1145 gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1146 nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1147 useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1148 // TODO cr->duty & DUTY_PME should imply that a PME
1149 // algorithm is active, but currently does not.
1150 EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1152 // Get the device handles for the modules, nullptr when no task is assigned.
1153 int deviceId = -1;
1154 DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1156 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1157 bool useTiming = true;
1158 if (GMX_GPU == GMX_GPU_CUDA)
1160 /* WARNING: CUDA timings are incorrect with multiple streams.
1161 * This is the main reason why they are disabled by default.
1163 // TODO: Consider turning on by default when we can detect nr of streams.
1164 useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1166 else if (GMX_GPU == GMX_GPU_OPENCL)
1168 useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1171 // TODO Currently this is always built, yet DD partition code
1172 // checks if it is built before using it. Probably it should
1173 // become an MDModule that is made only when another module
1174 // requires it (e.g. pull, CompEl, density fitting), so that we
1175 // don't update the local atom sets unilaterally every step.
1176 LocalAtomSetManager atomSets;
1177 if (ddBuilder)
1179 // TODO Pass the GPU streams to ddBuilder to use in buffer
1180 // transfers (e.g. halo exchange)
1181 cr->dd = ddBuilder->build(&atomSets);
1182 // The builder's job is done, so destruct it
1183 ddBuilder.reset(nullptr);
1184 // Note that local state still does not exist yet.
1187 // The GPU update is decided here because we need to know whether the constraints or
1188 // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1189 // defined). This is only known after DD is initialized, hence decision on using GPU
1190 // update is done so late.
1193 const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1195 useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1196 useDomainDecomposition, useUpdateGroups, pmeRunMode, domdecOptions.numPmeRanks > 0,
1197 useGpuForNonbonded, updateTarget, gpusWereDetected, *inputrec, mtop,
1198 doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1199 replExParams.exchangeInterval > 0, doRerun, devFlags, mdlog);
1201 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1203 const bool printHostName = (cr->nnodes > 1);
1204 gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1206 std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1208 if (deviceInfo != nullptr)
1210 if (DOMAINDECOMP(cr) && thisRankHasDuty(cr, DUTY_PP))
1212 dd_setup_dlb_resource_sharing(cr, deviceId);
1214 deviceStreamManager = std::make_unique<DeviceStreamManager>(
1215 *deviceInfo, useGpuForPme, useGpuForNonbonded, havePPDomainDecomposition(cr),
1216 useGpuForUpdate, useTiming);
1219 // If the user chose a task assignment, give them some hints
1220 // where appropriate.
1221 if (!userGpuTaskAssignment.empty())
1223 gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1226 if (PAR(cr))
1228 /* After possible communicator splitting in make_dd_communicators.
1229 * we can set up the intra/inter node communication.
1231 gmx_setup_nodecomm(fplog, cr);
1234 #if GMX_MPI
1235 if (isMultiSim(ms))
1237 GMX_LOG(mdlog.warning)
1238 .asParagraph()
1239 .appendTextFormatted(
1240 "This is simulation %d out of %d running as a composite GROMACS\n"
1241 "multi-simulation job. Setup for this simulation:\n",
1242 ms->sim, ms->nsim);
1244 GMX_LOG(mdlog.warning)
1245 .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1246 # if GMX_THREAD_MPI
1247 cr->nnodes == 1 ? "thread" : "threads"
1248 # else
1249 cr->nnodes == 1 ? "process" : "processes"
1250 # endif
1252 fflush(stderr);
1253 #endif
1255 // If mdrun -pin auto honors any affinity setting that already
1256 // exists. If so, it is nice to provide feedback about whether
1257 // that existing affinity setting was from OpenMP or something
1258 // else, so we run this code both before and after we initialize
1259 // the OpenMP support.
1260 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1261 /* Check and update the number of OpenMP threads requested */
1262 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1263 pmeRunMode, mtop, *inputrec);
1265 gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1266 hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1268 // Enable FP exception detection, but not in
1269 // Release mode and not for compilers with known buggy FP
1270 // exception support (clang with any optimization) or suspected
1271 // buggy FP exception support (gcc 7.* with optimization).
1272 #if !defined NDEBUG \
1273 && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1274 && defined __OPTIMIZE__)
1275 const bool bEnableFPE = true;
1276 #else
1277 const bool bEnableFPE = false;
1278 #endif
1279 // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1280 if (bEnableFPE)
1282 gmx_feenableexcept();
1285 /* Now that we know the setup is consistent, check for efficiency */
1286 check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1287 mdrunOptions.ntompOptionIsSet, cr, mdlog);
1289 /* getting number of PP/PME threads on this MPI / tMPI rank.
1290 PME: env variable should be read only on one node to make sure it is
1291 identical everywhere;
1293 const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1294 : gmx_omp_nthreads_get(emntPME);
1295 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1296 physicalNodeComm, mdlog);
1298 // Enable Peer access between GPUs where available
1299 // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1300 // any of the GPU communication features are active.
1301 if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1302 && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1304 setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1307 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1309 /* Before setting affinity, check whether the affinity has changed
1310 * - which indicates that probably the OpenMP library has changed it
1311 * since we first checked).
1313 gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1315 int numThreadsOnThisNode, intraNodeThreadOffset;
1316 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1317 &intraNodeThreadOffset);
1319 /* Set the CPU affinity */
1320 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1321 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1324 if (mdrunOptions.timingOptions.resetStep > -1)
1326 GMX_LOG(mdlog.info)
1327 .asParagraph()
1328 .appendText(
1329 "The -resetstep functionality is deprecated, and may be removed in a "
1330 "future version.");
1332 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1334 if (PAR(cr))
1336 /* Master synchronizes its value of reset_counters with all nodes
1337 * including PME only nodes */
1338 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1339 gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1340 wcycle_set_reset_counters(wcycle, reset_counters);
1343 // Membrane embedding must be initialized before we call init_forcerec()
1344 membedHolder.initializeMembed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1345 globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1347 const bool thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1348 std::unique_ptr<MDAtoms> mdAtoms;
1349 std::unique_ptr<VirtualSitesHandler> vsite;
1350 std::unique_ptr<GpuBonded> gpuBonded;
1352 t_nrnb nrnb;
1353 if (thisRankHasDuty(cr, DUTY_PP))
1355 mdModulesNotifier.notify(*cr);
1356 mdModulesNotifier.notify(&atomSets);
1357 mdModulesNotifier.notify(inputrec->pbcType);
1358 mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1359 /* Initiate forcerecord */
1360 fr = new t_forcerec;
1361 fr->forceProviders = mdModules_->initForceProviders();
1362 init_forcerec(fplog, mdlog, fr, inputrec, &mtop, cr, box,
1363 opt2fn("-table", filenames.size(), filenames.data()),
1364 opt2fn("-tablep", filenames.size(), filenames.data()),
1365 opt2fns("-tableb", filenames.size(), filenames.data()), pforce);
1366 // Dirty hack, for fixing disres and orires should be made mdmodules
1367 fr->listedForces->fcdata().disres = disresdata;
1368 fr->listedForces->fcdata().orires = oriresdata;
1370 // Save a handle to device stream manager to use elsewhere in the code
1371 // TODO: Forcerec is not a correct place to store it.
1372 fr->deviceStreamManager = deviceStreamManager.get();
1374 if (devFlags.enableGpuPmePPComm && !thisRankHasDuty(cr, DUTY_PME))
1376 GMX_RELEASE_ASSERT(
1377 deviceStreamManager != nullptr,
1378 "GPU device stream manager should be valid in order to use PME-PP direct "
1379 "communications.");
1380 GMX_RELEASE_ASSERT(
1381 deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1382 "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1383 "communications.");
1384 fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1385 cr->mpi_comm_mysim, cr->dd->pme_nodeid, deviceStreamManager->context(),
1386 deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1389 fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, useGpuForNonbonded,
1390 deviceStreamManager.get(), &mtop, box, wcycle);
1391 // TODO: Move the logic below to a GPU bonded builder
1392 if (useGpuForBonded)
1394 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1395 "GPU device stream manager should be valid in order to use GPU "
1396 "version of bonded forces.");
1397 gpuBonded = std::make_unique<GpuBonded>(
1398 mtop.ffparams, fr->ic->epsfac * fr->fudgeQQ, deviceStreamManager->context(),
1399 deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)), wcycle);
1400 fr->gpuBonded = gpuBonded.get();
1403 /* Initialize the mdAtoms structure.
1404 * mdAtoms is not filled with atom data,
1405 * as this can not be done now with domain decomposition.
1407 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1408 if (globalState && thisRankHasPmeGpuTask)
1410 // The pinning of coordinates in the global state object works, because we only use
1411 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1412 // points to the global state object without DD.
1413 // FIXME: MD and EM separately set up the local state - this should happen in the same
1414 // function, which should also perform the pinning.
1415 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1418 /* Initialize the virtual site communication */
1419 vsite = makeVirtualSitesHandler(mtop, cr, fr->pbcType);
1421 calc_shifts(box, fr->shift_vec);
1423 /* With periodic molecules the charge groups should be whole at start up
1424 * and the virtual sites should not be far from their proper positions.
1426 if (!inputrec->bContinuation && MASTER(cr)
1427 && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1429 /* Make molecules whole at start of run */
1430 if (fr->pbcType != PbcType::No)
1432 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1434 if (vsite)
1436 /* Correct initial vsite positions are required
1437 * for the initial distribution in the domain decomposition
1438 * and for the initial shell prediction.
1440 constructVirtualSitesGlobal(mtop, globalState->x);
1444 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1446 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1447 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1450 else
1452 /* This is a PME only node */
1454 GMX_ASSERT(globalState == nullptr,
1455 "We don't need the state on a PME only rank and expect it to be unitialized");
1457 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1458 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1461 gmx_pme_t* sepPmeData = nullptr;
1462 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1463 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1464 "Double-checking that only PME-only ranks have no forcerec");
1465 gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1467 // TODO should live in ewald module once its testing is improved
1469 // Later, this program could contain kernels that might be later
1470 // re-used as auto-tuning progresses, or subsequent simulations
1471 // are invoked.
1472 PmeGpuProgramStorage pmeGpuProgram;
1473 if (thisRankHasPmeGpuTask)
1475 GMX_RELEASE_ASSERT(
1476 (deviceStreamManager != nullptr),
1477 "GPU device stream manager should be initialized in order to use GPU for PME.");
1478 GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1479 "GPU device should be initialized in order to use GPU for PME.");
1480 pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1483 /* Initiate PME if necessary,
1484 * either on all nodes or on dedicated PME nodes only. */
1485 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1487 if (mdAtoms && mdAtoms->mdatoms())
1489 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1490 if (EVDW_PME(inputrec->vdwtype))
1492 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1495 if (cr->npmenodes > 0)
1497 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1498 gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1499 gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1502 if (thisRankHasDuty(cr, DUTY_PME))
1506 // TODO: This should be in the builder.
1507 GMX_RELEASE_ASSERT(!useGpuForPme || (deviceStreamManager != nullptr),
1508 "Device stream manager should be valid in order to use GPU "
1509 "version of PME.");
1510 GMX_RELEASE_ASSERT(
1511 !useGpuForPme || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1512 "GPU PME stream should be valid in order to use GPU version of PME.");
1514 const DeviceContext* deviceContext =
1515 useGpuForPme ? &deviceStreamManager->context() : nullptr;
1516 const DeviceStream* pmeStream =
1517 useGpuForPme ? &deviceStreamManager->stream(DeviceStreamType::Pme) : nullptr;
1519 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1520 nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1521 ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1522 nullptr, deviceContext, pmeStream, pmeGpuProgram.get(), mdlog);
1524 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1529 if (EI_DYNAMICS(inputrec->eI))
1531 /* Turn on signal handling on all nodes */
1533 * (A user signal from the PME nodes (if any)
1534 * is communicated to the PP nodes.
1536 signal_handler_install();
1539 pull_t* pull_work = nullptr;
1540 if (thisRankHasDuty(cr, DUTY_PP))
1542 /* Assumes uniform use of the number of OpenMP threads */
1543 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1545 if (inputrec->bPull)
1547 /* Initialize pull code */
1548 pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1549 inputrec->fepvals->init_lambda);
1550 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1552 initPullHistory(pull_work, &observablesHistory);
1554 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1556 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1560 std::unique_ptr<EnforcedRotation> enforcedRotation;
1561 if (inputrec->bRot)
1563 /* Initialize enforced rotation code */
1564 enforcedRotation =
1565 init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1566 globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1569 t_swap* swap = nullptr;
1570 if (inputrec->eSwapCoords != eswapNO)
1572 /* Initialize ion swapping code */
1573 swap = init_swapcoords(fplog, inputrec,
1574 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1575 &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1576 oenv, mdrunOptions, startingBehavior);
1579 /* Let makeConstraints know whether we have essential dynamics constraints. */
1580 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog, cr,
1581 ms, &nrnb, wcycle, fr->bMolPBC);
1583 /* Energy terms and groups */
1584 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1585 inputrec->fepvals->n_lambda);
1587 /* Kinetic energy data */
1588 gmx_ekindata_t ekind;
1589 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1591 /* Set up interactive MD (IMD) */
1592 auto imdSession =
1593 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1594 MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1595 filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1597 if (DOMAINDECOMP(cr))
1599 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1600 /* This call is not included in init_domain_decomposition mainly
1601 * because fr->cginfo_mb is set later.
1603 dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1604 domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1607 // TODO This is not the right place to manage the lifetime of
1608 // this data structure, but currently it's the easiest way to
1609 // make it work.
1610 MdrunScheduleWorkload runScheduleWork;
1611 // Also populates the simulation constant workload description.
1612 runScheduleWork.simulationWork =
1613 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1614 useGpuForUpdate, devFlags.enableGpuBufferOps,
1615 devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1617 std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1618 if (gpusWereDetected
1619 && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1620 || runScheduleWork.simulationWork.useGpuBufferOps))
1622 GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1623 ? GpuApiCallBehavior::Async
1624 : GpuApiCallBehavior::Sync;
1625 GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1626 "GPU device stream manager should be initialized to use GPU.");
1627 stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1628 *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle);
1629 fr->stateGpu = stateGpu.get();
1632 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1633 SimulatorBuilder simulatorBuilder;
1635 simulatorBuilder.add(SimulatorStateData(globalState.get(), &observablesHistory, &enerd, &ekind));
1636 simulatorBuilder.add(std::move(membedHolder));
1637 simulatorBuilder.add(std::move(stopHandlerBuilder_));
1638 simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
1641 simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv));
1642 simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle));
1643 simulatorBuilder.add(ConstraintsParam(
1644 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1645 vsite.get()));
1646 // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
1647 simulatorBuilder.add(
1648 LegacyInput(static_cast<int>(filenames.size()), filenames.data(), inputrec, fr));
1649 simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
1650 simulatorBuilder.add(InteractiveMD(imdSession.get()));
1651 simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifier()));
1652 simulatorBuilder.add(CenterOfMassPulling(pull_work));
1653 // Todo move to an MDModule
1654 simulatorBuilder.add(IonSwapping(swap));
1655 simulatorBuilder.add(TopologyData(&mtop, mdAtoms.get()));
1657 // build and run simulator object based on user-input
1658 auto simulator =
1659 simulatorBuilder.build(useModularSimulator, deform.get(), swap, &mtop, mdAtoms.get());
1660 simulator->run();
1662 if (fr->pmePpCommGpu)
1664 // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1665 fr->pmePpCommGpu.reset();
1668 if (inputrec->bPull)
1670 finish_pull(pull_work);
1672 finish_swapcoords(swap);
1674 else
1676 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1677 /* do PME only */
1678 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1679 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode,
1680 deviceStreamManager.get());
1683 wallcycle_stop(wcycle, ewcRUN);
1685 /* Finish up, write some stuff
1686 * if rerunMD, don't write last frame again
1688 finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1689 fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1691 // clean up cycle counter
1692 wallcycle_destroy(wcycle);
1694 deviceStreamManager.reset(nullptr);
1695 // Free PME data
1696 if (pmedata)
1698 gmx_pme_destroy(pmedata);
1699 pmedata = nullptr;
1702 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1703 // before we destroy the GPU context(s) in free_gpu().
1704 // Pinned buffers are associated with contexts in CUDA.
1705 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1706 mdAtoms.reset(nullptr);
1707 globalState.reset(nullptr);
1708 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1709 gpuBonded.reset(nullptr);
1710 /* Free pinned buffers in *fr */
1711 delete fr;
1712 fr = nullptr;
1713 // TODO convert to C++ so we can get rid of these frees
1714 sfree(disresdata);
1715 sfree(oriresdata);
1717 if (hwinfo->gpu_info.n_dev > 0)
1719 /* stop the GPU profiler (only CUDA) */
1720 stopGpuProfiler();
1723 /* With tMPI we need to wait for all ranks to finish deallocation before
1724 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1725 * GPU and context.
1727 * This is not a concern in OpenCL where we use one context per rank which
1728 * is freed in nbnxn_gpu_free().
1730 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1731 * but it is easier and more futureproof to call it on the whole node.
1733 * Note that this function needs to be called even if GPUs are not used
1734 * in this run because the PME ranks have no knowledge of whether GPUs
1735 * are used or not, but all ranks need to enter the barrier below.
1736 * \todo Remove this physical node barrier after making sure
1737 * that it's not needed anymore (with a shared GPU run).
1739 if (GMX_THREAD_MPI)
1741 physicalNodeComm.barrier();
1744 free_gpu(deviceInfo);
1746 /* Does what it says */
1747 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1748 walltime_accounting_destroy(walltime_accounting);
1750 // Ensure log file content is written
1751 if (logFileHandle)
1753 gmx_fio_flush(logFileHandle);
1756 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1757 * exceptions were enabled before function was called. */
1758 if (bEnableFPE)
1760 gmx_fedisableexcept();
1763 auto rc = static_cast<int>(gmx_get_stop_condition());
1765 #if GMX_THREAD_MPI
1766 /* we need to join all threads. The sub-threads join when they
1767 exit this function, but the master thread needs to be told to
1768 wait for that. */
1769 if (PAR(cr) && MASTER(cr))
1771 tMPI_Finalize();
1773 #endif
1774 return rc;
1775 } // namespace gmx
1777 Mdrunner::~Mdrunner()
1779 // Clean up of the Manager.
1780 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1781 // but okay as long as threads synchronize some time before adding or accessing
1782 // a new set of restraints.
1783 if (restraintManager_)
1785 restraintManager_->clear();
1786 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1787 "restraints added during runner life time should be cleared at runner "
1788 "destruction.");
1792 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1794 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1795 // Not sure if this should be logged through the md logger or something else,
1796 // but it is helpful to have some sort of INFO level message sent somewhere.
1797 // std::cout << "Registering restraint named " << name << std::endl;
1799 // When multiple restraints are used, it may be wasteful to register them separately.
1800 // Maybe instead register an entire Restraint Manager as a force provider.
1801 restraintManager_->addToSpec(std::move(puller), name);
1804 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1806 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1808 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept = default;
1810 class Mdrunner::BuilderImplementation
1812 public:
1813 BuilderImplementation() = delete;
1814 BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1815 ~BuilderImplementation();
1817 BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1818 real forceWarningThreshold,
1819 StartingBehavior startingBehavior);
1821 void addDomdec(const DomdecOptions& options);
1823 void addVerletList(int nstlist);
1825 void addReplicaExchange(const ReplicaExchangeParameters& params);
1827 void addNonBonded(const char* nbpu_opt);
1829 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1831 void addBondedTaskAssignment(const char* bonded_opt);
1833 void addUpdateTaskAssignment(const char* update_opt);
1835 void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1837 void addFilenames(ArrayRef<const t_filenm> filenames);
1839 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1841 void addLogFile(t_fileio* logFileHandle);
1843 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1845 Mdrunner build();
1847 private:
1848 // Default parameters copied from runner.h
1849 // \todo Clarify source(s) of default parameters.
1851 const char* nbpu_opt_ = nullptr;
1852 const char* pme_opt_ = nullptr;
1853 const char* pme_fft_opt_ = nullptr;
1854 const char* bonded_opt_ = nullptr;
1855 const char* update_opt_ = nullptr;
1857 MdrunOptions mdrunOptions_;
1859 DomdecOptions domdecOptions_;
1861 ReplicaExchangeParameters replicaExchangeParameters_;
1863 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1864 int nstlist_ = 0;
1866 //! Multisim communicator handle.
1867 gmx_multisim_t* multiSimulation_;
1869 //! mdrun communicator
1870 MPI_Comm communicator_ = MPI_COMM_NULL;
1872 //! Print a warning if any force is larger than this (in kJ/mol nm).
1873 real forceWarningThreshold_ = -1;
1875 //! Whether the simulation will start afresh, or restart with/without appending.
1876 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1878 //! The modules that comprise the functionality of mdrun.
1879 std::unique_ptr<MDModules> mdModules_;
1881 //! \brief Parallelism information.
1882 gmx_hw_opt_t hardwareOptions_;
1884 //! filename options for simulation.
1885 ArrayRef<const t_filenm> filenames_;
1887 /*! \brief Handle to output environment.
1889 * \todo gmx_output_env_t needs lifetime management.
1891 gmx_output_env_t* outputEnvironment_ = nullptr;
1893 /*! \brief Non-owning handle to MD log file.
1895 * \todo Context should own output facilities for client.
1896 * \todo Improve log file handle management.
1897 * \internal
1898 * Code managing the FILE* relies on the ability to set it to
1899 * nullptr to check whether the filehandle is valid.
1901 t_fileio* logFileHandle_ = nullptr;
1904 * \brief Builder for simulation stop signal handler.
1906 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1909 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1910 compat::not_null<SimulationContext*> context) :
1911 mdModules_(std::move(mdModules))
1913 communicator_ = context->communicator_;
1914 multiSimulation_ = context->multiSimulation_.get();
1917 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1919 Mdrunner::BuilderImplementation&
1920 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions& options,
1921 const real forceWarningThreshold,
1922 const StartingBehavior startingBehavior)
1924 mdrunOptions_ = options;
1925 forceWarningThreshold_ = forceWarningThreshold;
1926 startingBehavior_ = startingBehavior;
1927 return *this;
1930 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1932 domdecOptions_ = options;
1935 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1937 nstlist_ = nstlist;
1940 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1942 replicaExchangeParameters_ = params;
1945 Mdrunner Mdrunner::BuilderImplementation::build()
1947 auto newRunner = Mdrunner(std::move(mdModules_));
1949 newRunner.mdrunOptions = mdrunOptions_;
1950 newRunner.pforce = forceWarningThreshold_;
1951 newRunner.startingBehavior = startingBehavior_;
1952 newRunner.domdecOptions = domdecOptions_;
1954 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1955 newRunner.hw_opt = hardwareOptions_;
1957 // No invariant to check. This parameter exists to optionally override other behavior.
1958 newRunner.nstlist_cmdline = nstlist_;
1960 newRunner.replExParams = replicaExchangeParameters_;
1962 newRunner.filenames = filenames_;
1964 newRunner.communicator = communicator_;
1966 // nullptr is a valid value for the multisim handle
1967 newRunner.ms = multiSimulation_;
1969 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1970 // \todo Update sanity checking when output environment has clearly specified invariants.
1971 // Initialization and default values for oenv are not well specified in the current version.
1972 if (outputEnvironment_)
1974 newRunner.oenv = outputEnvironment_;
1976 else
1978 GMX_THROW(gmx::APIError(
1979 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1982 newRunner.logFileHandle = logFileHandle_;
1984 if (nbpu_opt_)
1986 newRunner.nbpu_opt = nbpu_opt_;
1988 else
1990 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1993 if (pme_opt_ && pme_fft_opt_)
1995 newRunner.pme_opt = pme_opt_;
1996 newRunner.pme_fft_opt = pme_fft_opt_;
1998 else
2000 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2003 if (bonded_opt_)
2005 newRunner.bonded_opt = bonded_opt_;
2007 else
2009 GMX_THROW(gmx::APIError(
2010 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2013 if (update_opt_)
2015 newRunner.update_opt = update_opt_;
2017 else
2019 GMX_THROW(gmx::APIError(
2020 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
2024 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2026 if (stopHandlerBuilder_)
2028 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2030 else
2032 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2035 return newRunner;
2038 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2040 nbpu_opt_ = nbpu_opt;
2043 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2045 pme_opt_ = pme_opt;
2046 pme_fft_opt_ = pme_fft_opt;
2049 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2051 bonded_opt_ = bonded_opt;
2054 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2056 update_opt_ = update_opt;
2059 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2061 hardwareOptions_ = hardwareOptions;
2064 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2066 filenames_ = filenames;
2069 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2071 outputEnvironment_ = outputEnvironment;
2074 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2076 logFileHandle_ = logFileHandle;
2079 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2081 stopHandlerBuilder_ = std::move(builder);
2084 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
2085 compat::not_null<SimulationContext*> context) :
2086 impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2090 MdrunnerBuilder::~MdrunnerBuilder() = default;
2092 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions& options,
2093 real forceWarningThreshold,
2094 const StartingBehavior startingBehavior)
2096 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2097 return *this;
2100 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2102 impl_->addDomdec(options);
2103 return *this;
2106 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2108 impl_->addVerletList(nstlist);
2109 return *this;
2112 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2114 impl_->addReplicaExchange(params);
2115 return *this;
2118 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2120 impl_->addNonBonded(nbpu_opt);
2121 return *this;
2124 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2126 // The builder method may become more general in the future, but in this version,
2127 // parameters for PME electrostatics are both required and the only parameters
2128 // available.
2129 if (pme_opt && pme_fft_opt)
2131 impl_->addPME(pme_opt, pme_fft_opt);
2133 else
2135 GMX_THROW(
2136 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2138 return *this;
2141 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2143 impl_->addBondedTaskAssignment(bonded_opt);
2144 return *this;
2147 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2149 impl_->addUpdateTaskAssignment(update_opt);
2150 return *this;
2153 Mdrunner MdrunnerBuilder::build()
2155 return impl_->build();
2158 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2160 impl_->addHardwareOptions(hardwareOptions);
2161 return *this;
2164 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2166 impl_->addFilenames(filenames);
2167 return *this;
2170 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2172 impl_->addOutputEnvironment(outputEnvironment);
2173 return *this;
2176 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2178 impl_->addLogFile(logFileHandle);
2179 return *this;
2182 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2184 impl_->addStopHandlerBuilder(std::move(builder));
2185 return *this;
2188 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2190 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;
2192 } // namespace gmx