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,2012,2013,2014,2015,2016,2017,2018,2019, 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/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/localatomsetmanager.h"
63 #include "gromacs/ewald/ewald_utils.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme_gpu_program.h"
66 #include "gromacs/fileio/checkpoint.h"
67 #include "gromacs/fileio/gmxfio.h"
68 #include "gromacs/fileio/oenv.h"
69 #include "gromacs/fileio/tpxio.h"
70 #include "gromacs/gmxlib/network.h"
71 #include "gromacs/gmxlib/nrnb.h"
72 #include "gromacs/gpu_utils/clfftinitializer.h"
73 #include "gromacs/gpu_utils/gpu_utils.h"
74 #include "gromacs/hardware/cpuinfo.h"
75 #include "gromacs/hardware/detecthardware.h"
76 #include "gromacs/hardware/printhardware.h"
77 #include "gromacs/listed_forces/disre.h"
78 #include "gromacs/listed_forces/gpubonded.h"
79 #include "gromacs/listed_forces/orires.h"
80 #include "gromacs/math/functions.h"
81 #include "gromacs/math/utilities.h"
82 #include "gromacs/math/vec.h"
83 #include "gromacs/mdlib/boxdeformation.h"
84 #include "gromacs/mdlib/calc_verletbuf.h"
85 #include "gromacs/mdlib/forcerec.h"
86 #include "gromacs/mdlib/gmx_omp_nthreads.h"
87 #include "gromacs/mdlib/makeconstraints.h"
88 #include "gromacs/mdlib/md_support.h"
89 #include "gromacs/mdlib/mdatoms.h"
90 #include "gromacs/mdlib/mdrun.h"
91 #include "gromacs/mdlib/membed.h"
92 #include "gromacs/mdlib/nb_verlet.h"
93 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
94 #include "gromacs/mdlib/nbnxn_search.h"
95 #include "gromacs/mdlib/nbnxn_tuning.h"
96 #include "gromacs/mdlib/ppforceworkload.h"
97 #include "gromacs/mdlib/qmmm.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/sim_util.h"
100 #include "gromacs/mdlib/stophandler.h"
101 #include "gromacs/mdrun/legacymdrunoptions.h"
102 #include "gromacs/mdrun/logging.h"
103 #include "gromacs/mdrun/multisim.h"
104 #include "gromacs/mdrun/simulationcontext.h"
105 #include "gromacs/mdrunutility/mdmodules.h"
106 #include "gromacs/mdrunutility/threadaffinity.h"
107 #include "gromacs/mdtypes/commrec.h"
108 #include "gromacs/mdtypes/fcdata.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/md_enums.h"
111 #include "gromacs/mdtypes/observableshistory.h"
112 #include "gromacs/mdtypes/state.h"
113 #include "gromacs/pbcutil/pbc.h"
114 #include "gromacs/pulling/output.h"
115 #include "gromacs/pulling/pull.h"
116 #include "gromacs/pulling/pull_rotation.h"
117 #include "gromacs/restraint/manager.h"
118 #include "gromacs/restraint/restraintmdmodule.h"
119 #include "gromacs/restraint/restraintpotential.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/taskassignment/decidegpuusage.h"
122 #include "gromacs/taskassignment/resourcedivision.h"
123 #include "gromacs/taskassignment/taskassignment.h"
124 #include "gromacs/taskassignment/usergpuids.h"
125 #include "gromacs/timing/wallcycle.h"
126 #include "gromacs/topology/mtop_util.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basenetwork.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/exceptions.h"
131 #include "gromacs/utility/fatalerror.h"
132 #include "gromacs/utility/filestream.h"
133 #include "gromacs/utility/gmxassert.h"
134 #include "gromacs/utility/gmxmpi.h"
135 #include "gromacs/utility/logger.h"
136 #include "gromacs/utility/loggerbuilder.h"
137 #include "gromacs/utility/physicalnodecommunicator.h"
138 #include "gromacs/utility/pleasecite.h"
139 #include "gromacs/utility/programcontext.h"
140 #include "gromacs/utility/smalloc.h"
141 #include "gromacs/utility/stringutil.h"
143 #include "integrator.h"
144 #include "replicaexchange.h"
147 #include "corewrap.h"
153 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
155 * Used to ensure that the master thread does not modify mdrunner during copy
156 * on the spawned threads. */
157 static void threadMpiMdrunnerAccessBarrier()
160 MPI_Barrier(MPI_COMM_WORLD);
164 Mdrunner Mdrunner::cloneOnSpawnedThread() const
166 auto newRunner = Mdrunner();
168 // All runners in the same process share a restraint manager resource because it is
169 // part of the interface to the client code, which is associated only with the
170 // original thread. Handles to the same resources can be obtained by copy.
172 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
175 // Copy original cr pointer before master thread can pass the thread barrier
176 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
178 // Copy members of master runner.
179 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
180 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
181 newRunner.hw_opt = hw_opt;
182 newRunner.filenames = filenames;
184 newRunner.oenv = oenv;
185 newRunner.mdrunOptions = mdrunOptions;
186 newRunner.domdecOptions = domdecOptions;
187 newRunner.nbpu_opt = nbpu_opt;
188 newRunner.pme_opt = pme_opt;
189 newRunner.pme_fft_opt = pme_fft_opt;
190 newRunner.bonded_opt = bonded_opt;
191 newRunner.nstlist_cmdline = nstlist_cmdline;
192 newRunner.replExParams = replExParams;
193 newRunner.pforce = pforce;
195 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
197 threadMpiMdrunnerAccessBarrier();
199 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
204 /*! \brief The callback used for running on spawned threads.
206 * Obtains the pointer to the master mdrunner object from the one
207 * argument permitted to the thread-launch API call, copies it to make
208 * a new runner for this thread, reinitializes necessary data, and
209 * proceeds to the simulation. */
210 static void mdrunner_start_fn(const void *arg)
214 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
215 /* copy the arg list to make sure that it's thread-local. This
216 doesn't copy pointed-to items, of course; fnm, cr and fplog
217 are reset in the call below, all others should be const. */
218 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
221 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
225 /*! \brief Start thread-MPI threads.
227 * Called by mdrunner() to start a specific number of threads
228 * (including the main thread) for thread-parallel runs. This in turn
229 * calls mdrunner() for each thread. All options are the same as for
231 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
234 /* first check whether we even need to start tMPI */
235 if (numThreadsToLaunch < 2)
241 /* now spawn new threads that start mdrunner_start_fn(), while
242 the main thread returns, we set thread affinity later */
243 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
244 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
246 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
249 threadMpiMdrunnerAccessBarrier();
251 GMX_UNUSED_VALUE(mdrunner_start_fn);
254 return reinitialize_commrec_for_this_thread(cr);
259 /*! \brief Initialize variables for Verlet scheme simulation */
260 static void prepare_verlet_scheme(FILE *fplog,
264 const gmx_mtop_t *mtop,
266 bool makeGpuPairList,
267 const gmx::CpuInfo &cpuinfo)
269 /* For NVE simulations, we will retain the initial list buffer */
270 if (EI_DYNAMICS(ir->eI) &&
271 ir->verletbuf_tol > 0 &&
272 !(EI_MD(ir->eI) && ir->etc == etcNO))
274 /* Update the Verlet buffer size for the current run setup */
276 /* Here we assume SIMD-enabled kernels are being used. But as currently
277 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
278 * and 4x2 gives a larger buffer than 4x4, this is ok.
280 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
281 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
284 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &listSetup, nullptr, &rlist_new);
286 if (rlist_new != ir->rlist)
288 if (fplog != nullptr)
290 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
291 ir->rlist, rlist_new,
292 listSetup.cluster_size_i, listSetup.cluster_size_j);
294 ir->rlist = rlist_new;
298 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
300 gmx_fatal(FARGS, "Can not set nstlist without %s",
301 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
304 if (EI_DYNAMICS(ir->eI))
306 /* Set or try nstlist values */
307 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
311 /*! \brief Override the nslist value in inputrec
313 * with value passed on the command line (if any)
315 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
316 int64_t nsteps_cmdline,
321 /* override with anything else than the default -2 */
322 if (nsteps_cmdline > -2)
324 char sbuf_steps[STEPSTRSIZE];
325 char sbuf_msg[STRLEN];
327 ir->nsteps = nsteps_cmdline;
328 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
330 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
331 gmx_step_str(nsteps_cmdline, sbuf_steps),
332 fabs(nsteps_cmdline*ir->delta_t));
336 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
337 gmx_step_str(nsteps_cmdline, sbuf_steps));
340 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
342 else if (nsteps_cmdline < -2)
344 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
347 /* Do nothing if nsteps_cmdline == -2 */
353 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
355 * If not, and if a warning may be issued, logs a warning about
356 * falling back to CPU code. With thread-MPI, only the first
357 * call to this function should have \c issueWarning true. */
358 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
359 const t_inputrec *ir,
362 if (ir->opts.ngener - ir->nwall > 1)
364 /* The GPU code does not support more than one energy group.
365 * If the user requested GPUs explicitly, a fatal error is given later.
369 GMX_LOG(mdlog.warning).asParagraph()
370 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
371 "For better performance, run on the GPU without energy groups and then do "
372 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
379 //! Initializes the logger for mdrun.
380 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
382 gmx::LoggerBuilder builder;
383 if (fplog != nullptr)
385 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
387 if (cr == nullptr || SIMMASTER(cr))
389 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
390 &gmx::TextOutputFile::standardError());
392 return builder.build();
395 //! Make a TaskTarget from an mdrun argument string.
396 static TaskTarget findTaskTarget(const char *optionString)
398 TaskTarget returnValue = TaskTarget::Auto;
400 if (strncmp(optionString, "auto", 3) == 0)
402 returnValue = TaskTarget::Auto;
404 else if (strncmp(optionString, "cpu", 3) == 0)
406 returnValue = TaskTarget::Cpu;
408 else if (strncmp(optionString, "gpu", 3) == 0)
410 returnValue = TaskTarget::Gpu;
414 GMX_ASSERT(false, "Option string should have been checked for sanity already");
420 int Mdrunner::mdrunner()
424 t_forcerec *fr = nullptr;
425 t_fcdata *fcd = nullptr;
426 real ewaldcoeff_q = 0;
427 real ewaldcoeff_lj = 0;
428 int nChargePerturbed = -1, nTypePerturbed = 0;
429 gmx_wallcycle_t wcycle;
430 gmx_walltime_accounting_t walltime_accounting = nullptr;
431 gmx_membed_t * membed = nullptr;
432 gmx_hw_info_t *hwinfo = nullptr;
434 /* CAUTION: threads may be started later on in this function, so
435 cr doesn't reflect the final parallel state right now */
436 std::unique_ptr<gmx::MDModules> mdModules(new gmx::MDModules);
437 t_inputrec inputrecInstance;
438 t_inputrec *inputrec = &inputrecInstance;
441 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
442 bool doRerun = mdrunOptions.rerun;
444 // Handle task-assignment related user options.
445 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
446 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
447 std::vector<int> gpuIdsAvailable;
450 gpuIdsAvailable = parseUserGpuIdString(hw_opt.gpuIdsAvailable);
452 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
454 std::vector<int> userGpuTaskAssignment;
457 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
459 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
460 auto nonbondedTarget = findTaskTarget(nbpu_opt);
461 auto pmeTarget = findTaskTarget(pme_opt);
462 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
463 auto bondedTarget = findTaskTarget(bonded_opt);
464 PmeRunMode pmeRunMode = PmeRunMode::None;
466 // Here we assume that SIMMASTER(cr) does not change even after the
467 // threads are started.
469 FILE *fplog = nullptr;
470 // If we are appending, we don't write log output because we need
471 // to check that the old log file matches what the checkpoint file
472 // expects. Otherwise, we should start to write log output now if
473 // there is a file ready for it.
474 if (logFileHandle != nullptr && !mdrunOptions.continuationOptions.appendFiles)
476 fplog = gmx_fio_getfp(logFileHandle);
478 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
479 gmx::MDLogger mdlog(logOwner.logger());
481 // TODO The thread-MPI master rank makes a working
482 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
483 // after the threads have been launched. This works because no use
484 // is made of that communicator until after the execution paths
485 // have rejoined. But it is likely that we can improve the way
486 // this is expressed, e.g. by expressly running detection only the
487 // master rank for thread-MPI, rather than relying on the mutex
488 // and reference count.
489 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
490 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
492 gmx_print_detected_hardware(fplog, cr, ms, mdlog, hwinfo);
494 std::vector<int> gpuIdsToUse;
495 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
496 if (gpuIdsAvailable.empty())
498 gpuIdsToUse = compatibleGpus;
502 for (const auto &availableGpuId : gpuIdsAvailable)
504 bool availableGpuIsCompatible = false;
505 for (const auto &compatibleGpuId : compatibleGpus)
507 if (availableGpuId == compatibleGpuId)
509 availableGpuIsCompatible = true;
513 if (!availableGpuIsCompatible)
515 gmx_fatal(FARGS, "You limited the set of compatible GPUs to a set that included ID #%d, but that ID is not for a compatible GPU. List only compatible GPUs.", availableGpuId);
517 gpuIdsToUse.push_back(availableGpuId);
521 if (fplog != nullptr)
523 /* Print references after all software/hardware printing */
524 please_cite(fplog, "Abraham2015");
525 please_cite(fplog, "Pall2015");
526 please_cite(fplog, "Pronk2013");
527 please_cite(fplog, "Hess2008b");
528 please_cite(fplog, "Spoel2005a");
529 please_cite(fplog, "Lindahl2001a");
530 please_cite(fplog, "Berendsen95a");
531 writeSourceDoi(fplog);
534 std::unique_ptr<t_state> globalState;
538 /* Only the master rank has the global state */
539 globalState = std::make_unique<t_state>();
541 /* Read (nearly) all data required for the simulation */
542 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
544 /* In rerun, set velocities to zero if present */
545 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
547 // rerun does not use velocities
548 GMX_LOG(mdlog.info).asParagraph().appendText(
549 "Rerun trajectory contains velocities. Rerun does only evaluate "
550 "potential energy and forces. The velocities will be ignored.");
551 for (int i = 0; i < globalState->natoms; i++)
553 clear_rvec(globalState->v[i]);
555 globalState->flags &= ~(1 << estV);
558 if (inputrec->cutoff_scheme != ecutsVERLET)
560 if (nstlist_cmdline > 0)
562 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
565 if (!compatibleGpus.empty())
567 GMX_LOG(mdlog.warning).asParagraph().appendText(
568 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
569 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
574 /* Check and update the hardware options for internal consistency */
575 check_and_update_hw_opt_1(mdlog, &hw_opt, cr, domdecOptions.numPmeRanks);
577 /* Early check for externally set process affinity. */
578 gmx_check_thread_affinity_set(mdlog, cr,
579 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
581 if (GMX_THREAD_MPI && SIMMASTER(cr))
583 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
585 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
588 /* Since the master knows the cut-off scheme, update hw_opt for this.
589 * This is done later for normal MPI and also once more with tMPI
590 * for all tMPI ranks.
592 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
594 bool useGpuForNonbonded = false;
595 bool useGpuForPme = false;
598 // If the user specified the number of ranks, then we must
599 // respect that, but in default mode, we need to allow for
600 // the number of GPUs to choose the number of ranks.
601 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
602 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
603 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
604 canUseGpuForNonbonded,
605 inputrec->cutoff_scheme == ecutsVERLET,
606 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
607 hw_opt.nthreads_tmpi);
608 auto canUseGpuForPme = pme_gpu_supports_build(*hwinfo, nullptr) && pme_gpu_supports_input(*inputrec, mtop, nullptr);
609 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
610 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
611 canUseGpuForPme, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
614 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
616 /* Determine how many thread-MPI ranks to start.
618 * TODO Over-writing the user-supplied value here does
619 * prevent any possible subsequent checks from working
621 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
630 // Now start the threads for thread MPI.
631 cr = spawnThreads(hw_opt.nthreads_tmpi);
632 /* The main thread continues here with a new cr. We don't deallocate
633 the old cr because other threads may still be reading it. */
634 // TODO Both master and spawned threads call dup_tfn and
635 // reinitialize_commrec_for_this_thread. Find a way to express
637 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
639 // END OF CAUTION: cr and physicalNodeComm are now reliable
643 /* now broadcast everything to the non-master nodes/threads: */
644 init_parallel(cr, inputrec, &mtop);
647 // Now each rank knows the inputrec that SIMMASTER read and used,
648 // and (if applicable) cr->nnodes has been assigned the number of
649 // thread-MPI ranks that have been chosen. The ranks can now all
650 // run the task-deciding functions and will agree on the result
651 // without needing to communicate.
653 // TODO Should we do the communication in debug mode to support
654 // having an assertion?
656 // Note that these variables describe only their own node.
658 // Note that when bonded interactions run on a GPU they always run
659 // alongside a nonbonded task, so do not influence task assignment
660 // even though they affect the force calculation workload.
661 bool useGpuForNonbonded = false;
662 bool useGpuForPme = false;
663 bool useGpuForBonded = false;
666 // It's possible that there are different numbers of GPUs on
667 // different nodes, which is the user's responsibilty to
668 // handle. If unsuitable, we will notice that during task
670 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
671 bool usingVerletScheme = inputrec->cutoff_scheme == ecutsVERLET;
672 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
673 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
675 canUseGpuForNonbonded,
677 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
679 auto canUseGpuForPme = pme_gpu_supports_build(*hwinfo, nullptr) && pme_gpu_supports_input(*inputrec, mtop, nullptr);
680 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
681 canUseGpuForPme, cr->nnodes, domdecOptions.numPmeRanks,
683 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
685 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
686 bondedTarget, canUseGpuForBonded,
687 EVDW_PME(inputrec->vdwtype),
688 EEL_PME_EWALD(inputrec->coulombtype),
689 domdecOptions.numPmeRanks, gpusWereDetected);
691 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
692 if (pmeRunMode == PmeRunMode::GPU)
694 if (pmeFftTarget == TaskTarget::Cpu)
696 pmeRunMode = PmeRunMode::Mixed;
699 else if (pmeFftTarget == TaskTarget::Gpu)
701 gmx_fatal(FARGS, "Assigning FFTs to GPU requires PME to be assigned to GPU as well. With PME on CPU you should not be using -pmefft.");
704 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
707 // TODO: hide restraint implementation details from Mdrunner.
708 // There is nothing unique about restraints at this point as far as the
709 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
710 // factory functions from the SimulationContext on which to call mdModules->add().
711 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
712 for (auto && restraint : restraintManager_->getRestraints())
714 auto module = RestraintMDModule::create(restraint,
716 mdModules->add(std::move(module));
719 // TODO: Error handling
720 mdModules->assignOptionsToModules(*inputrec->params, nullptr);
722 if (fplog != nullptr)
724 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
725 fprintf(fplog, "\n");
730 /* now make sure the state is initialized and propagated */
731 set_state_entries(globalState.get(), inputrec);
734 /* NM and TPI parallelize over force/energy calculations, not atoms,
735 * so we need to initialize and broadcast the global state.
737 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
741 globalState = std::make_unique<t_state>();
743 broadcastStateWithoutDynamics(cr, globalState.get());
746 /* A parallel command line option consistency check that we can
747 only do after any threads have started. */
748 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
749 domdecOptions.numCells[YY] > 1 ||
750 domdecOptions.numCells[ZZ] > 1 ||
751 domdecOptions.numPmeRanks > 0))
754 "The -dd or -npme option request a parallel simulation, "
756 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
759 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
761 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
767 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
769 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
772 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
774 if (domdecOptions.numPmeRanks > 0)
776 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
777 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
780 domdecOptions.numPmeRanks = 0;
783 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
785 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
786 * improve performance with many threads per GPU, since our OpenMP
787 * scaling is bad, but it's difficult to automate the setup.
789 domdecOptions.numPmeRanks = 0;
793 if (domdecOptions.numPmeRanks < 0)
795 domdecOptions.numPmeRanks = 0;
796 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
800 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
807 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
811 /* NMR restraints must be initialized before load_checkpoint,
812 * since with time averaging the history is added to t_state.
813 * For proper consistency check we therefore need to extend
815 * So the PME-only nodes (if present) will also initialize
816 * the distance restraints.
820 /* This needs to be called before read_checkpoint to extend the state */
821 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
823 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
825 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
827 ObservablesHistory observablesHistory = {};
829 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
831 if (continuationOptions.startedFromCheckpoint)
833 /* Check if checkpoint file exists before doing continuation.
834 * This way we can use identical input options for the first and subsequent runs...
838 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
840 cr, domdecOptions.numCells,
841 inputrec, globalState.get(),
842 &bReadEkin, &observablesHistory,
843 continuationOptions.appendFiles,
844 continuationOptions.appendFilesOptionSet,
845 mdrunOptions.reproducible);
849 continuationOptions.haveReadEkin = true;
852 if (continuationOptions.appendFiles && logFileHandle)
854 // Now we can start normal logging to the truncated log file.
855 fplog = gmx_fio_getfp(logFileHandle);
856 prepareLogAppending(fplog);
857 logOwner = buildLogger(fplog, cr);
858 mdlog = logOwner.logger();
862 if (mdrunOptions.numStepsCommandline > -2)
864 GMX_LOG(mdlog.info).asParagraph().
865 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
866 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
868 /* override nsteps with value set on the commamdline */
869 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
873 copy_mat(globalState->box, box);
878 gmx_bcast(sizeof(box), box, cr);
881 /* Update rlist and nstlist. */
882 if (inputrec->cutoff_scheme == ecutsVERLET)
884 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
885 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
888 LocalAtomSetManager atomSets;
890 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
891 inputrec->eI == eiNM))
893 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
895 box, positionsFromStatePointer(globalState.get()),
897 // Note that local state still does not exist yet.
901 /* PME, if used, is done on all nodes with 1D decomposition */
903 cr->duty = (DUTY_PP | DUTY_PME);
905 if (inputrec->ePBC == epbcSCREW)
908 "pbc=screw is only implemented with domain decomposition");
914 /* After possible communicator splitting in make_dd_communicators.
915 * we can set up the intra/inter node communication.
917 gmx_setup_nodecomm(fplog, cr);
923 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
924 "This is simulation %d out of %d running as a composite GROMACS\n"
925 "multi-simulation job. Setup for this simulation:\n",
928 GMX_LOG(mdlog.warning).appendTextFormatted(
932 cr->nnodes == 1 ? "thread" : "threads"
934 cr->nnodes == 1 ? "process" : "processes"
940 /* Check and update hw_opt for the cut-off scheme */
941 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
943 /* Check and update the number of OpenMP threads requested */
944 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
947 gmx_omp_nthreads_init(mdlog, cr,
948 hwinfo->nthreads_hw_avail,
949 physicalNodeComm.size_,
951 hw_opt.nthreads_omp_pme,
952 !thisRankHasDuty(cr, DUTY_PP),
953 inputrec->cutoff_scheme == ecutsVERLET);
955 // Enable FP exception detection for the Verlet scheme, but not in
956 // Release mode and not for compilers with known buggy FP
957 // exception support (clang with any optimization) or suspected
958 // buggy FP exception support (gcc 7.* with optimization).
959 #if !defined NDEBUG && \
960 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
961 && defined __OPTIMIZE__)
962 const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
964 const bool bEnableFPE = false;
966 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
969 gmx_feenableexcept();
972 // Build a data structure that expresses which kinds of non-bonded
973 // task are handled by this rank.
975 // TODO Later, this might become a loop over all registered modules
976 // relevant to the mdp inputs, to find those that have such tasks.
978 // TODO This could move before init_domain_decomposition() as part
979 // of refactoring that separates the responsibility for duty
980 // assignment from setup for communication between tasks, and
981 // setup for tasks handled with a domain (ie including short-ranged
982 // tasks, bonded tasks, etc.).
984 // Note that in general useGpuForNonbonded, etc. can have a value
985 // that is inconsistent with the presence of actual GPUs on any
986 // rank, and that is not known to be a problem until the
987 // duty of the ranks on a node become known.
989 // TODO Later we might need the concept of computeTasksOnThisRank,
990 // from which we construct gpuTasksOnThisRank.
992 // Currently the DD code assigns duty to ranks that can
993 // include PP work that currently can be executed on a single
994 // GPU, if present and compatible. This has to be coordinated
995 // across PP ranks on a node, with possible multiple devices
996 // or sharing devices on a node, either from the user
997 // selection, or automatically.
998 auto haveGpus = !gpuIdsToUse.empty();
999 std::vector<GpuTask> gpuTasksOnThisRank;
1000 if (thisRankHasDuty(cr, DUTY_PP))
1002 if (useGpuForNonbonded)
1004 // Note that any bonded tasks on a GPU always accompany a
1008 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1010 else if (nonbondedTarget == TaskTarget::Gpu)
1012 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1014 else if (bondedTarget == TaskTarget::Gpu)
1016 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because there is none detected.");
1020 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1021 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1027 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1029 else if (pmeTarget == TaskTarget::Gpu)
1031 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1036 GpuTaskAssignment gpuTaskAssignment;
1039 // Produce the task assignment for this rank.
1040 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1041 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1042 useGpuForBonded, pmeRunMode);
1044 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1046 /* Prevent other ranks from continuing after an issue was found
1047 * and reported as a fatal error.
1049 * TODO This function implements a barrier so that MPI runtimes
1050 * can organize an orderly shutdown if one of the ranks has had to
1051 * issue a fatal error in various code already run. When we have
1052 * MPI-aware error handling and reporting, this should be
1057 MPI_Barrier(cr->mpi_comm_mysim);
1063 MPI_Barrier(ms->mpi_comm_masters);
1065 /* We need another barrier to prevent non-master ranks from contiuing
1066 * when an error occured in a different simulation.
1068 MPI_Barrier(cr->mpi_comm_mysim);
1072 /* Now that we know the setup is consistent, check for efficiency */
1073 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1076 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1078 if (thisRankHasDuty(cr, DUTY_PP))
1080 // This works because only one task of each type is currently permitted.
1081 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1082 hasTaskType<GpuTask::Nonbonded>);
1083 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1085 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1086 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1087 init_gpu(nonbondedDeviceInfo);
1089 if (DOMAINDECOMP(cr))
1091 /* When we share GPUs over ranks, we need to know this for the DLB */
1092 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1098 std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1100 gmx_device_info_t *pmeDeviceInfo = nullptr;
1101 // Later, this program could contain kernels that might be later
1102 // re-used as auto-tuning progresses, or subsequent simulations
1104 PmeGpuProgramStorage pmeGpuProgram;
1105 // This works because only one task of each type is currently permitted.
1106 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1107 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1108 if (thisRankHasPmeGpuTask)
1110 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1111 init_gpu(pmeDeviceInfo);
1112 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1113 // TODO It would be nice to move this logic into the factory
1114 // function. See Redmine #2535.
1115 bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1116 if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1118 initializedClfftLibrary = initializeClfftLibrary();
1122 /* getting number of PP/PME threads on this MPI / tMPI rank.
1123 PME: env variable should be read only on one node to make sure it is
1124 identical everywhere;
1126 const int numThreadsOnThisRank =
1127 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1128 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1129 *hwinfo->hardwareTopology,
1130 physicalNodeComm, mdlog);
1132 if (hw_opt.thread_affinity != threadaffOFF)
1134 /* Before setting affinity, check whether the affinity has changed
1135 * - which indicates that probably the OpenMP library has changed it
1136 * since we first checked).
1138 gmx_check_thread_affinity_set(mdlog, cr,
1139 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1141 int numThreadsOnThisNode, intraNodeThreadOffset;
1142 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1143 &intraNodeThreadOffset);
1145 /* Set the CPU affinity */
1146 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1147 numThreadsOnThisRank, numThreadsOnThisNode,
1148 intraNodeThreadOffset, nullptr);
1151 if (mdrunOptions.timingOptions.resetStep > -1)
1153 GMX_LOG(mdlog.info).asParagraph().
1154 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1156 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1160 /* Master synchronizes its value of reset_counters with all nodes
1161 * including PME only nodes */
1162 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1163 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1164 wcycle_set_reset_counters(wcycle, reset_counters);
1167 // Membrane embedding must be initialized before we call init_forcerec()
1172 fprintf(stderr, "Initializing membed");
1174 /* Note that membed cannot work in parallel because mtop is
1175 * changed here. Fix this if we ever want to make it run with
1176 * multiple ranks. */
1177 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1179 .checkpointOptions.period);
1182 std::unique_ptr<MDAtoms> mdAtoms;
1183 std::unique_ptr<gmx_vsite_t> vsite;
1186 if (thisRankHasDuty(cr, DUTY_PP))
1188 /* Initiate forcerecord */
1190 fr->forceProviders = mdModules->initForceProviders();
1191 init_forcerec(fplog, mdlog, fr, fcd,
1192 inputrec, &mtop, cr, box,
1193 opt2fn("-table", filenames.size(), filenames.data()),
1194 opt2fn("-tablep", filenames.size(), filenames.data()),
1195 opt2fns("-tableb", filenames.size(), filenames.data()),
1196 *hwinfo, nonbondedDeviceInfo,
1201 /* Initialize the mdAtoms structure.
1202 * mdAtoms is not filled with atom data,
1203 * as this can not be done now with domain decomposition.
1205 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1206 if (globalState && thisRankHasPmeGpuTask)
1208 // The pinning of coordinates in the global state object works, because we only use
1209 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1210 // points to the global state object without DD.
1211 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1212 // which should also perform the pinning.
1213 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1216 /* Initialize the virtual site communication */
1217 vsite = initVsite(mtop, cr);
1219 calc_shifts(box, fr->shift_vec);
1221 /* With periodic molecules the charge groups should be whole at start up
1222 * and the virtual sites should not be far from their proper positions.
1224 if (!inputrec->bContinuation && MASTER(cr) &&
1225 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1227 /* Make molecules whole at start of run */
1228 if (fr->ePBC != epbcNONE)
1230 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1234 /* Correct initial vsite positions are required
1235 * for the initial distribution in the domain decomposition
1236 * and for the initial shell prediction.
1238 constructVsitesGlobal(mtop, globalState->x);
1242 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1244 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1245 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1250 /* This is a PME only node */
1252 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1254 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1255 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1258 gmx_pme_t *sepPmeData = nullptr;
1259 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1260 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1261 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1263 /* Initiate PME if necessary,
1264 * either on all nodes or on dedicated PME nodes only. */
1265 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1267 if (mdAtoms && mdAtoms->mdatoms())
1269 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1270 if (EVDW_PME(inputrec->vdwtype))
1272 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1275 if (cr->npmenodes > 0)
1277 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1278 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1279 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1282 if (thisRankHasDuty(cr, DUTY_PME))
1286 pmedata = gmx_pme_init(cr,
1287 getNumPmeDomains(cr->dd),
1289 mtop.natoms, nChargePerturbed != 0, nTypePerturbed != 0,
1290 mdrunOptions.reproducible,
1291 ewaldcoeff_q, ewaldcoeff_lj,
1292 gmx_omp_nthreads_get(emntPME),
1293 pmeRunMode, nullptr,
1294 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1296 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1301 if (EI_DYNAMICS(inputrec->eI))
1303 /* Turn on signal handling on all nodes */
1305 * (A user signal from the PME nodes (if any)
1306 * is communicated to the PP nodes.
1308 signal_handler_install();
1311 if (thisRankHasDuty(cr, DUTY_PP))
1313 /* Assumes uniform use of the number of OpenMP threads */
1314 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1316 if (inputrec->bPull)
1318 /* Initialize pull code */
1319 inputrec->pull_work =
1320 init_pull(fplog, inputrec->pull, inputrec,
1321 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1322 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1324 initPullHistory(inputrec->pull_work, &observablesHistory);
1326 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1328 init_pull_output_files(inputrec->pull_work,
1329 filenames.size(), filenames.data(), oenv,
1330 continuationOptions);
1334 std::unique_ptr<EnforcedRotation> enforcedRotation;
1337 /* Initialize enforced rotation code */
1338 enforcedRotation = init_rot(fplog,
1350 if (inputrec->eSwapCoords != eswapNO)
1352 /* Initialize ion swapping code */
1353 init_swapcoords(fplog, inputrec, opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1354 &mtop, globalState.get(), &observablesHistory,
1355 cr, &atomSets, oenv, mdrunOptions);
1358 /* Let makeConstraints know whether we have essential dynamics constraints.
1359 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1361 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1362 || observablesHistory.edsamHistory);
1363 auto constr = makeConstraints(mtop, *inputrec, doEssentialDynamics,
1364 fplog, *mdAtoms->mdatoms(),
1365 cr, ms, nrnb, wcycle, fr->bMolPBC);
1367 if (DOMAINDECOMP(cr))
1369 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1370 /* This call is not included in init_domain_decomposition mainly
1371 * because fr->cginfo_mb is set later.
1373 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1374 domdecOptions.checkBondedInteractions,
1378 // TODO This is not the right place to manage the lifetime of
1379 // this data structure, but currently it's the easiest way to
1380 // make it work. Later, it should probably be made/updated
1381 // after the workload for the lifetime of a PP domain is
1383 PpForceWorkload ppForceWorkload;
1385 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to integrator.");
1386 /* Now do whatever the user wants us to do (how flexible...) */
1387 Integrator integrator {
1388 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1391 vsite.get(), constr.get(),
1392 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1394 mdModules->outputProvider(),
1398 &observablesHistory,
1399 mdAtoms.get(), nrnb, wcycle, fr,
1403 walltime_accounting,
1404 std::move(stopHandlerBuilder_)
1406 integrator.run(inputrec->eI, doRerun);
1408 if (inputrec->bPull)
1410 finish_pull(inputrec->pull_work);
1416 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1418 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1419 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1422 wallcycle_stop(wcycle, ewcRUN);
1424 /* Finish up, write some stuff
1425 * if rerunMD, don't write last frame again
1427 finish_run(fplog, mdlog, cr,
1428 inputrec, nrnb, wcycle, walltime_accounting,
1429 fr ? fr->nbv : nullptr,
1431 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1433 // clean up cycle counter
1434 wallcycle_destroy(wcycle);
1439 gmx_pme_destroy(pmedata);
1443 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1444 // before we destroy the GPU context(s) in free_gpu_resources().
1445 // Pinned buffers are associated with contexts in CUDA.
1446 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1447 mdAtoms.reset(nullptr);
1448 globalState.reset(nullptr);
1449 mdModules.reset(nullptr); // destruct force providers here as they might also use the GPU
1451 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1452 free_gpu_resources(fr, physicalNodeComm);
1453 free_gpu(nonbondedDeviceInfo);
1454 free_gpu(pmeDeviceInfo);
1455 done_forcerec(fr, mtop.molblock.size(), mtop.groups.grps[egcENER].nr);
1460 free_membed(membed);
1463 gmx_hardware_info_free();
1465 /* Does what it says */
1466 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1467 walltime_accounting_destroy(walltime_accounting);
1470 // Ensure log file content is written
1473 gmx_fio_flush(logFileHandle);
1476 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1477 * exceptions were enabled before function was called. */
1480 gmx_fedisableexcept();
1483 auto rc = static_cast<int>(gmx_get_stop_condition());
1486 /* we need to join all threads. The sub-threads join when they
1487 exit this function, but the master thread needs to be told to
1489 if (PAR(cr) && MASTER(cr))
1493 //TODO free commrec in MPI simulations
1499 Mdrunner::~Mdrunner()
1501 // Clean up of the Manager.
1502 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1503 // but okay as long as threads synchronize some time before adding or accessing
1504 // a new set of restraints.
1505 if (restraintManager_)
1507 restraintManager_->clear();
1508 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1509 "restraints added during runner life time should be cleared at runner destruction.");
1513 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1516 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1517 // Not sure if this should be logged through the md logger or something else,
1518 // but it is helpful to have some sort of INFO level message sent somewhere.
1519 // std::cout << "Registering restraint named " << name << std::endl;
1521 // When multiple restraints are used, it may be wasteful to register them separately.
1522 // Maybe instead register an entire Restraint Manager as a force provider.
1523 restraintManager_->addToSpec(std::move(puller),
1527 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1529 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1530 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1532 class Mdrunner::BuilderImplementation
1535 BuilderImplementation() = delete;
1536 explicit BuilderImplementation(SimulationContext* context);
1537 ~BuilderImplementation();
1539 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1540 real forceWarningThreshold);
1542 void addDomdec(const DomdecOptions &options);
1544 void addVerletList(int nstlist);
1546 void addReplicaExchange(const ReplicaExchangeParameters ¶ms);
1548 void addMultiSim(gmx_multisim_t* multisim);
1550 void addNonBonded(const char* nbpu_opt);
1552 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1554 void addBondedTaskAssignment(const char* bonded_opt);
1556 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1558 void addFilenames(ArrayRef <const t_filenm> filenames);
1560 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1562 void addLogFile(t_fileio *logFileHandle);
1564 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1569 // Default parameters copied from runner.h
1570 // \todo Clarify source(s) of default parameters.
1572 const char* nbpu_opt_ = nullptr;
1573 const char* pme_opt_ = nullptr;
1574 const char* pme_fft_opt_ = nullptr;
1575 const char *bonded_opt_ = nullptr;
1577 MdrunOptions mdrunOptions_;
1579 DomdecOptions domdecOptions_;
1581 ReplicaExchangeParameters replicaExchangeParameters_;
1583 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1586 //! Non-owning multisim communicator handle.
1587 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1589 //! Print a warning if any force is larger than this (in kJ/mol nm).
1590 real forceWarningThreshold_ = -1;
1592 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1595 * \todo Establish robust protocol to make sure resources remain valid.
1596 * SimulationContext will likely be separated into multiple layers for
1597 * different levels of access and different phases of execution. Ref
1598 * https://redmine.gromacs.org/issues/2375
1599 * https://redmine.gromacs.org/issues/2587
1601 SimulationContext* context_ = nullptr;
1603 //! \brief Parallelism information.
1604 gmx_hw_opt_t hardwareOptions_;
1606 //! filename options for simulation.
1607 ArrayRef<const t_filenm> filenames_;
1609 /*! \brief Handle to output environment.
1611 * \todo gmx_output_env_t needs lifetime management.
1613 gmx_output_env_t* outputEnvironment_ = nullptr;
1615 /*! \brief Non-owning handle to MD log file.
1617 * \todo Context should own output facilities for client.
1618 * \todo Improve log file handle management.
1620 * Code managing the FILE* relies on the ability to set it to
1621 * nullptr to check whether the filehandle is valid.
1623 t_fileio* logFileHandle_ = nullptr;
1626 * \brief Builder for simulation stop signal handler.
1628 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1631 Mdrunner::BuilderImplementation::BuilderImplementation(SimulationContext* context) :
1634 GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1637 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1639 Mdrunner::BuilderImplementation &
1640 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1641 real forceWarningThreshold)
1643 mdrunOptions_ = options;
1644 forceWarningThreshold_ = forceWarningThreshold;
1648 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1650 domdecOptions_ = options;
1653 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1658 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1660 replicaExchangeParameters_ = params;
1663 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1665 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1668 Mdrunner Mdrunner::BuilderImplementation::build()
1670 auto newRunner = Mdrunner();
1672 GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1674 newRunner.mdrunOptions = mdrunOptions_;
1675 newRunner.domdecOptions = domdecOptions_;
1677 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1678 newRunner.hw_opt = hardwareOptions_;
1680 // No invariant to check. This parameter exists to optionally override other behavior.
1681 newRunner.nstlist_cmdline = nstlist_;
1683 newRunner.replExParams = replicaExchangeParameters_;
1685 newRunner.filenames = filenames_;
1687 GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1688 newRunner.cr = context_->communicationRecord_;
1692 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1693 newRunner.ms = *multisim_;
1697 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1700 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1701 // \todo Update sanity checking when output environment has clearly specified invariants.
1702 // Initialization and default values for oenv are not well specified in the current version.
1703 if (outputEnvironment_)
1705 newRunner.oenv = outputEnvironment_;
1709 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1712 newRunner.logFileHandle = logFileHandle_;
1716 newRunner.nbpu_opt = nbpu_opt_;
1720 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1723 if (pme_opt_ && pme_fft_opt_)
1725 newRunner.pme_opt = pme_opt_;
1726 newRunner.pme_fft_opt = pme_fft_opt_;
1730 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1735 newRunner.bonded_opt = bonded_opt_;
1739 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1742 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1744 if (stopHandlerBuilder_)
1746 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1750 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1756 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1758 nbpu_opt_ = nbpu_opt;
1761 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1762 const char* pme_fft_opt)
1765 pme_fft_opt_ = pme_fft_opt;
1768 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1770 bonded_opt_ = bonded_opt;
1773 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1775 hardwareOptions_ = hardwareOptions;
1778 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1780 filenames_ = filenames;
1783 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1785 outputEnvironment_ = outputEnvironment;
1788 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1790 logFileHandle_ = logFileHandle;
1793 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1795 stopHandlerBuilder_ = std::move(builder);
1798 MdrunnerBuilder::MdrunnerBuilder(compat::not_null<SimulationContext*> context) :
1799 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(context)}
1803 MdrunnerBuilder::~MdrunnerBuilder() = default;
1805 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1806 real forceWarningThreshold)
1808 impl_->setExtraMdrunOptions(options, forceWarningThreshold);
1812 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1814 impl_->addDomdec(options);
1818 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1820 impl_->addVerletList(nstlist);
1824 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1826 impl_->addReplicaExchange(params);
1830 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1832 impl_->addMultiSim(multisim);
1836 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1838 impl_->addNonBonded(nbpu_opt);
1842 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
1843 const char* pme_fft_opt)
1845 // The builder method may become more general in the future, but in this version,
1846 // parameters for PME electrostatics are both required and the only parameters
1848 if (pme_opt && pme_fft_opt)
1850 impl_->addPME(pme_opt, pme_fft_opt);
1854 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
1859 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
1861 impl_->addBondedTaskAssignment(bonded_opt);
1865 Mdrunner MdrunnerBuilder::build()
1867 return impl_->build();
1870 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1872 impl_->addHardwareOptions(hardwareOptions);
1876 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
1878 impl_->addFilenames(filenames);
1882 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1884 impl_->addOutputEnvironment(outputEnvironment);
1888 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
1890 impl_->addLogFile(logFileHandle);
1894 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1896 impl_->addStopHandlerBuilder(std::move(builder));
1900 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
1902 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;