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/domdec/partition.h"
64 #include "gromacs/ewald/ewald_utils.h"
65 #include "gromacs/ewald/pme.h"
66 #include "gromacs/ewald/pme_gpu_program.h"
67 #include "gromacs/fileio/checkpoint.h"
68 #include "gromacs/fileio/gmxfio.h"
69 #include "gromacs/fileio/oenv.h"
70 #include "gromacs/fileio/tpxio.h"
71 #include "gromacs/gmxlib/network.h"
72 #include "gromacs/gmxlib/nrnb.h"
73 #include "gromacs/gpu_utils/clfftinitializer.h"
74 #include "gromacs/gpu_utils/gpu_utils.h"
75 #include "gromacs/hardware/cpuinfo.h"
76 #include "gromacs/hardware/detecthardware.h"
77 #include "gromacs/hardware/printhardware.h"
78 #include "gromacs/imd/imd.h"
79 #include "gromacs/listed_forces/disre.h"
80 #include "gromacs/listed_forces/gpubonded.h"
81 #include "gromacs/listed_forces/orires.h"
82 #include "gromacs/math/functions.h"
83 #include "gromacs/math/utilities.h"
84 #include "gromacs/math/vec.h"
85 #include "gromacs/mdlib/boxdeformation.h"
86 #include "gromacs/mdlib/broadcaststructs.h"
87 #include "gromacs/mdlib/calc_verletbuf.h"
88 #include "gromacs/mdlib/dispersioncorrection.h"
89 #include "gromacs/mdlib/enerdata_utils.h"
90 #include "gromacs/mdlib/force.h"
91 #include "gromacs/mdlib/forcerec.h"
92 #include "gromacs/mdlib/gmx_omp_nthreads.h"
93 #include "gromacs/mdlib/makeconstraints.h"
94 #include "gromacs/mdlib/md_support.h"
95 #include "gromacs/mdlib/mdatoms.h"
96 #include "gromacs/mdlib/membed.h"
97 #include "gromacs/mdlib/ppforceworkload.h"
98 #include "gromacs/mdlib/qmmm.h"
99 #include "gromacs/mdlib/sighandler.h"
100 #include "gromacs/mdlib/stophandler.h"
101 #include "gromacs/mdrun/mdmodules.h"
102 #include "gromacs/mdrun/simulationcontext.h"
103 #include "gromacs/mdrunutility/handlerestart.h"
104 #include "gromacs/mdrunutility/logging.h"
105 #include "gromacs/mdrunutility/mdmodulenotification.h"
106 #include "gromacs/mdrunutility/multisim.h"
107 #include "gromacs/mdrunutility/printtime.h"
108 #include "gromacs/mdrunutility/threadaffinity.h"
109 #include "gromacs/mdtypes/commrec.h"
110 #include "gromacs/mdtypes/enerdata.h"
111 #include "gromacs/mdtypes/fcdata.h"
112 #include "gromacs/mdtypes/group.h"
113 #include "gromacs/mdtypes/inputrec.h"
114 #include "gromacs/mdtypes/md_enums.h"
115 #include "gromacs/mdtypes/mdrunoptions.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/state.h"
118 #include "gromacs/nbnxm/gpu_data_mgmt.h"
119 #include "gromacs/nbnxm/nbnxm.h"
120 #include "gromacs/nbnxm/pairlist_tuning.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/output.h"
123 #include "gromacs/pulling/pull.h"
124 #include "gromacs/pulling/pull_rotation.h"
125 #include "gromacs/restraint/manager.h"
126 #include "gromacs/restraint/restraintmdmodule.h"
127 #include "gromacs/restraint/restraintpotential.h"
128 #include "gromacs/swap/swapcoords.h"
129 #include "gromacs/taskassignment/decidegpuusage.h"
130 #include "gromacs/taskassignment/resourcedivision.h"
131 #include "gromacs/taskassignment/taskassignment.h"
132 #include "gromacs/taskassignment/usergpuids.h"
133 #include "gromacs/timing/gpu_timing.h"
134 #include "gromacs/timing/wallcycle.h"
135 #include "gromacs/timing/wallcyclereporting.h"
136 #include "gromacs/topology/mtop_util.h"
137 #include "gromacs/trajectory/trajectoryframe.h"
138 #include "gromacs/utility/basenetwork.h"
139 #include "gromacs/utility/cstringutil.h"
140 #include "gromacs/utility/exceptions.h"
141 #include "gromacs/utility/fatalerror.h"
142 #include "gromacs/utility/filestream.h"
143 #include "gromacs/utility/gmxassert.h"
144 #include "gromacs/utility/gmxmpi.h"
145 #include "gromacs/utility/keyvaluetree.h"
146 #include "gromacs/utility/logger.h"
147 #include "gromacs/utility/loggerbuilder.h"
148 #include "gromacs/utility/physicalnodecommunicator.h"
149 #include "gromacs/utility/pleasecite.h"
150 #include "gromacs/utility/programcontext.h"
151 #include "gromacs/utility/smalloc.h"
152 #include "gromacs/utility/stringutil.h"
154 #include "isimulator.h"
155 #include "replicaexchange.h"
156 #include "simulatorbuilder.h"
159 #include "corewrap.h"
165 /*! \brief Log if development feature flags are encountered
167 * The use of dev features indicated by environment variables is logged
168 * in order to ensure that runs with such featrues enabled can be identified
169 * from their log and standard output.
171 * \param[in] mdlog Logger object.
173 static void reportDevelopmentFeatures(const gmx::MDLogger &mdlog)
175 const bool enableGpuBufOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
176 const bool useGpuUpdateConstrain = (getenv("GMX_UPDATE_CONSTRAIN_GPU") != nullptr);
180 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
181 "NOTE: This run uses the 'GPU buffer ops' feature, enabled by the GMX_USE_GPU_BUFFER_OPS environment variable.");
184 if (useGpuUpdateConstrain)
186 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
187 "NOTE: This run uses the 'GPU update/constraints' feature, enabled by the GMX_UPDATE_CONSTRAIN_GPU environment variable.");
191 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
193 * Used to ensure that the master thread does not modify mdrunner during copy
194 * on the spawned threads. */
195 static void threadMpiMdrunnerAccessBarrier()
198 MPI_Barrier(MPI_COMM_WORLD);
202 Mdrunner Mdrunner::cloneOnSpawnedThread() const
204 auto newRunner = Mdrunner(std::make_unique<MDModules>());
206 // All runners in the same process share a restraint manager resource because it is
207 // part of the interface to the client code, which is associated only with the
208 // original thread. Handles to the same resources can be obtained by copy.
210 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
213 // Copy original cr pointer before master thread can pass the thread barrier
214 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
216 // Copy members of master runner.
217 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
218 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
219 newRunner.hw_opt = hw_opt;
220 newRunner.filenames = filenames;
222 newRunner.oenv = oenv;
223 newRunner.mdrunOptions = mdrunOptions;
224 newRunner.domdecOptions = domdecOptions;
225 newRunner.nbpu_opt = nbpu_opt;
226 newRunner.pme_opt = pme_opt;
227 newRunner.pme_fft_opt = pme_fft_opt;
228 newRunner.bonded_opt = bonded_opt;
229 newRunner.nstlist_cmdline = nstlist_cmdline;
230 newRunner.replExParams = replExParams;
231 newRunner.pforce = pforce;
233 newRunner.startingBehavior = startingBehavior;
234 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
236 threadMpiMdrunnerAccessBarrier();
238 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
243 /*! \brief The callback used for running on spawned threads.
245 * Obtains the pointer to the master mdrunner object from the one
246 * argument permitted to the thread-launch API call, copies it to make
247 * a new runner for this thread, reinitializes necessary data, and
248 * proceeds to the simulation. */
249 static void mdrunner_start_fn(const void *arg)
253 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
254 /* copy the arg list to make sure that it's thread-local. This
255 doesn't copy pointed-to items, of course; fnm, cr and fplog
256 are reset in the call below, all others should be const. */
257 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
260 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
264 /*! \brief Start thread-MPI threads.
266 * Called by mdrunner() to start a specific number of threads
267 * (including the main thread) for thread-parallel runs. This in turn
268 * calls mdrunner() for each thread. All options are the same as for
270 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
273 /* first check whether we even need to start tMPI */
274 if (numThreadsToLaunch < 2)
280 /* now spawn new threads that start mdrunner_start_fn(), while
281 the main thread returns. Thread affinity is handled later. */
282 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
283 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
285 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
288 threadMpiMdrunnerAccessBarrier();
290 GMX_UNUSED_VALUE(mdrunner_start_fn);
293 return reinitialize_commrec_for_this_thread(cr);
298 /*! \brief Initialize variables for Verlet scheme simulation */
299 static void prepare_verlet_scheme(FILE *fplog,
303 const gmx_mtop_t *mtop,
305 bool makeGpuPairList,
306 const gmx::CpuInfo &cpuinfo)
308 /* For NVE simulations, we will retain the initial list buffer */
309 if (EI_DYNAMICS(ir->eI) &&
310 ir->verletbuf_tol > 0 &&
311 !(EI_MD(ir->eI) && ir->etc == etcNO))
313 /* Update the Verlet buffer size for the current run setup */
315 /* Here we assume SIMD-enabled kernels are being used. But as currently
316 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
317 * and 4x2 gives a larger buffer than 4x4, this is ok.
319 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
320 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
322 const real rlist_new =
323 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
325 if (rlist_new != ir->rlist)
327 if (fplog != nullptr)
329 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
330 ir->rlist, rlist_new,
331 listSetup.cluster_size_i, listSetup.cluster_size_j);
333 ir->rlist = rlist_new;
337 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
339 gmx_fatal(FARGS, "Can not set nstlist without %s",
340 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
343 if (EI_DYNAMICS(ir->eI))
345 /* Set or try nstlist values */
346 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
350 /*! \brief Override the nslist value in inputrec
352 * with value passed on the command line (if any)
354 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
355 int64_t nsteps_cmdline,
360 /* override with anything else than the default -2 */
361 if (nsteps_cmdline > -2)
363 char sbuf_steps[STEPSTRSIZE];
364 char sbuf_msg[STRLEN];
366 ir->nsteps = nsteps_cmdline;
367 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
369 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
370 gmx_step_str(nsteps_cmdline, sbuf_steps),
371 fabs(nsteps_cmdline*ir->delta_t));
375 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
376 gmx_step_str(nsteps_cmdline, sbuf_steps));
379 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
381 else if (nsteps_cmdline < -2)
383 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
386 /* Do nothing if nsteps_cmdline == -2 */
392 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
394 * If not, and if a warning may be issued, logs a warning about
395 * falling back to CPU code. With thread-MPI, only the first
396 * call to this function should have \c issueWarning true. */
397 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
398 const t_inputrec *ir,
401 if (ir->opts.ngener - ir->nwall > 1)
403 /* The GPU code does not support more than one energy group.
404 * If the user requested GPUs explicitly, a fatal error is given later.
408 GMX_LOG(mdlog.warning).asParagraph()
409 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
410 "For better performance, run on the GPU without energy groups and then do "
411 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
418 //! Initializes the logger for mdrun.
419 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
421 gmx::LoggerBuilder builder;
422 if (fplog != nullptr)
424 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
426 if (cr == nullptr || SIMMASTER(cr))
428 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
429 &gmx::TextOutputFile::standardError());
431 return builder.build();
434 //! Make a TaskTarget from an mdrun argument string.
435 static TaskTarget findTaskTarget(const char *optionString)
437 TaskTarget returnValue = TaskTarget::Auto;
439 if (strncmp(optionString, "auto", 3) == 0)
441 returnValue = TaskTarget::Auto;
443 else if (strncmp(optionString, "cpu", 3) == 0)
445 returnValue = TaskTarget::Cpu;
447 else if (strncmp(optionString, "gpu", 3) == 0)
449 returnValue = TaskTarget::Gpu;
453 GMX_ASSERT(false, "Option string should have been checked for sanity already");
459 //! Finish run, aggregate data to print performance info.
460 static void finish_run(FILE *fplog,
461 const gmx::MDLogger &mdlog,
463 const t_inputrec *inputrec,
464 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
465 gmx_walltime_accounting_t walltime_accounting,
466 nonbonded_verlet_t *nbv,
467 const gmx_pme_t *pme,
471 double nbfs = 0, mflop = 0;
473 elapsed_time_over_all_ranks,
474 elapsed_time_over_all_threads,
475 elapsed_time_over_all_threads_over_all_ranks;
476 /* Control whether it is valid to print a report. Only the
477 simulation master may print, but it should not do so if the run
478 terminated e.g. before a scheduled reset step. This is
479 complicated by the fact that PME ranks are unaware of the
480 reason why they were sent a pmerecvqxFINISH. To avoid
481 communication deadlocks, we always do the communication for the
482 report, even if we've decided not to write the report, because
483 how long it takes to finish the run is not important when we've
484 decided not to report on the simulation performance.
486 Further, we only report performance for dynamical integrators,
487 because those are the only ones for which we plan to
488 consider doing any optimizations. */
489 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
491 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
493 GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
498 std::unique_ptr<t_nrnb> nrnbTotalStorage;
501 nrnbTotalStorage = std::make_unique<t_nrnb>();
502 nrnb_tot = nrnbTotalStorage.get();
504 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
513 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
514 elapsed_time_over_all_threads = walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
518 /* reduce elapsed_time over all MPI ranks in the current simulation */
519 MPI_Allreduce(&elapsed_time,
520 &elapsed_time_over_all_ranks,
521 1, MPI_DOUBLE, MPI_SUM,
523 elapsed_time_over_all_ranks /= cr->nnodes;
524 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
525 * current simulation. */
526 MPI_Allreduce(&elapsed_time_over_all_threads,
527 &elapsed_time_over_all_threads_over_all_ranks,
528 1, MPI_DOUBLE, MPI_SUM,
534 elapsed_time_over_all_ranks = elapsed_time;
535 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
540 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
543 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
545 print_dd_statistics(cr, inputrec, fplog);
548 /* TODO Move the responsibility for any scaling by thread counts
549 * to the code that handled the thread region, so that there's a
550 * mechanism to keep cycle counting working during the transition
551 * to task parallelism. */
552 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
553 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
554 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
555 auto cycle_sum(wallcycle_sum(cr, wcycle));
559 auto nbnxn_gpu_timings = (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
560 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
562 if (pme_gpu_task_enabled(pme))
564 pme_gpu_get_timings(pme, &pme_gpu_timings);
566 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
567 elapsed_time_over_all_ranks,
572 if (EI_DYNAMICS(inputrec->eI))
574 delta_t = inputrec->delta_t;
579 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
580 elapsed_time_over_all_ranks,
581 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
582 delta_t, nbfs, mflop);
586 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
587 elapsed_time_over_all_ranks,
588 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
589 delta_t, nbfs, mflop);
594 int Mdrunner::mdrunner()
597 t_forcerec *fr = nullptr;
598 t_fcdata *fcd = nullptr;
599 real ewaldcoeff_q = 0;
600 real ewaldcoeff_lj = 0;
601 int nChargePerturbed = -1, nTypePerturbed = 0;
602 gmx_wallcycle_t wcycle;
603 gmx_walltime_accounting_t walltime_accounting = nullptr;
604 gmx_membed_t * membed = nullptr;
605 gmx_hw_info_t *hwinfo = nullptr;
607 /* CAUTION: threads may be started later on in this function, so
608 cr doesn't reflect the final parallel state right now */
609 t_inputrec inputrecInstance;
610 t_inputrec *inputrec = &inputrecInstance;
613 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
614 bool doRerun = mdrunOptions.rerun;
616 // Handle task-assignment related user options.
617 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
618 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
620 std::vector<int> userGpuTaskAssignment;
623 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
625 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
626 auto nonbondedTarget = findTaskTarget(nbpu_opt);
627 auto pmeTarget = findTaskTarget(pme_opt);
628 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
629 auto bondedTarget = findTaskTarget(bonded_opt);
630 PmeRunMode pmeRunMode = PmeRunMode::None;
632 // Here we assume that SIMMASTER(cr) does not change even after the
633 // threads are started.
635 FILE *fplog = nullptr;
636 // If we are appending, we don't write log output because we need
637 // to check that the old log file matches what the checkpoint file
638 // expects. Otherwise, we should start to write log output now if
639 // there is a file ready for it.
640 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
642 fplog = gmx_fio_getfp(logFileHandle);
644 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
645 gmx::MDLogger mdlog(logOwner.logger());
647 // report any development features that may be enabled by environment variables
648 reportDevelopmentFeatures(mdlog);
650 // With thread-MPI, the communicator changes after threads are
651 // launched, so this is rebuilt for the master rank at that
652 // time. The non-master ranks are fine to keep the one made here.
653 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
654 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
656 gmx_print_detected_hardware(fplog, isMasterSimMasterRank(ms, MASTER(cr)), mdlog, hwinfo);
658 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
660 // Print citation requests after all software/hardware printing
661 pleaseCiteGromacs(fplog);
663 std::unique_ptr<t_state> globalState;
667 /* Only the master rank has the global state */
668 globalState = std::make_unique<t_state>();
670 /* Read (nearly) all data required for the simulation */
671 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
674 /* Check and update the hardware options for internal consistency */
675 checkAndUpdateHardwareOptions(mdlog, &hw_opt, SIMMASTER(cr), domdecOptions.numPmeRanks);
677 if (GMX_THREAD_MPI && SIMMASTER(cr))
679 bool useGpuForNonbonded = false;
680 bool useGpuForPme = false;
683 // If the user specified the number of ranks, then we must
684 // respect that, but in default mode, we need to allow for
685 // the number of GPUs to choose the number of ranks.
686 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
687 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
688 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
689 canUseGpuForNonbonded,
690 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
691 hw_opt.nthreads_tmpi);
692 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
693 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
694 *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
697 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
699 /* Determine how many thread-MPI ranks to start.
701 * TODO Over-writing the user-supplied value here does
702 * prevent any possible subsequent checks from working
704 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
713 // Now start the threads for thread MPI.
714 cr = spawnThreads(hw_opt.nthreads_tmpi);
715 /* The main thread continues here with a new cr. We don't deallocate
716 the old cr because other threads may still be reading it. */
717 // TODO Both master and spawned threads call dup_tfn and
718 // reinitialize_commrec_for_this_thread. Find a way to express
720 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
722 // END OF CAUTION: cr and physicalNodeComm are now reliable
726 /* now broadcast everything to the non-master nodes/threads: */
727 init_parallel(cr, inputrec, &mtop);
730 // Now each rank knows the inputrec that SIMMASTER read and used,
731 // and (if applicable) cr->nnodes has been assigned the number of
732 // thread-MPI ranks that have been chosen. The ranks can now all
733 // run the task-deciding functions and will agree on the result
734 // without needing to communicate.
736 // TODO Should we do the communication in debug mode to support
737 // having an assertion?
739 // Note that these variables describe only their own node.
741 // Note that when bonded interactions run on a GPU they always run
742 // alongside a nonbonded task, so do not influence task assignment
743 // even though they affect the force calculation workload.
744 bool useGpuForNonbonded = false;
745 bool useGpuForPme = false;
746 bool useGpuForBonded = false;
749 // It's possible that there are different numbers of GPUs on
750 // different nodes, which is the user's responsibilty to
751 // handle. If unsuitable, we will notice that during task
753 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
754 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
755 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
757 canUseGpuForNonbonded,
758 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
760 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
761 *hwinfo, *inputrec, mtop,
762 cr->nnodes, domdecOptions.numPmeRanks,
764 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
766 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme,
767 bondedTarget, canUseGpuForBonded,
768 EVDW_PME(inputrec->vdwtype),
769 EEL_PME_EWALD(inputrec->coulombtype),
770 domdecOptions.numPmeRanks, gpusWereDetected);
772 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
773 if (pmeRunMode == PmeRunMode::GPU)
775 if (pmeFftTarget == TaskTarget::Cpu)
777 pmeRunMode = PmeRunMode::Mixed;
780 else if (pmeFftTarget == TaskTarget::Gpu)
782 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.");
785 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
788 // TODO: hide restraint implementation details from Mdrunner.
789 // There is nothing unique about restraints at this point as far as the
790 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
791 // factory functions from the SimulationContext on which to call mdModules_->add().
792 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
793 for (auto && restraint : restraintManager_->getRestraints())
795 auto module = RestraintMDModule::create(restraint,
797 mdModules_->add(std::move(module));
800 // TODO: Error handling
801 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
802 const auto &mdModulesNotifier = mdModules_->notifier().notifier_;
804 if (inputrec->internalParameters != nullptr)
806 mdModulesNotifier.notify(*inputrec->internalParameters);
809 if (fplog != nullptr)
811 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
812 fprintf(fplog, "\n");
817 /* In rerun, set velocities to zero if present */
818 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
820 // rerun does not use velocities
821 GMX_LOG(mdlog.info).asParagraph().appendText(
822 "Rerun trajectory contains velocities. Rerun does only evaluate "
823 "potential energy and forces. The velocities will be ignored.");
824 for (int i = 0; i < globalState->natoms; i++)
826 clear_rvec(globalState->v[i]);
828 globalState->flags &= ~(1 << estV);
831 /* now make sure the state is initialized and propagated */
832 set_state_entries(globalState.get(), inputrec);
835 /* NM and TPI parallelize over force/energy calculations, not atoms,
836 * so we need to initialize and broadcast the global state.
838 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
842 globalState = std::make_unique<t_state>();
844 broadcastStateWithoutDynamics(cr, globalState.get());
847 /* A parallel command line option consistency check that we can
848 only do after any threads have started. */
849 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
850 domdecOptions.numCells[YY] > 1 ||
851 domdecOptions.numCells[ZZ] > 1 ||
852 domdecOptions.numPmeRanks > 0))
855 "The -dd or -npme option request a parallel simulation, "
857 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
860 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
862 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
868 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
870 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
873 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
875 if (domdecOptions.numPmeRanks > 0)
877 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
878 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
881 domdecOptions.numPmeRanks = 0;
884 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
886 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
887 * improve performance with many threads per GPU, since our OpenMP
888 * scaling is bad, but it's difficult to automate the setup.
890 domdecOptions.numPmeRanks = 0;
894 if (domdecOptions.numPmeRanks < 0)
896 domdecOptions.numPmeRanks = 0;
897 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
901 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
908 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
912 /* NMR restraints must be initialized before load_checkpoint,
913 * since with time averaging the history is added to t_state.
914 * For proper consistency check we therefore need to extend
916 * So the PME-only nodes (if present) will also initialize
917 * the distance restraints.
921 /* This needs to be called before read_checkpoint to extend the state */
922 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
924 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
926 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
928 ObservablesHistory observablesHistory = {};
930 if (startingBehavior != StartingBehavior::NewSimulation)
932 /* Check if checkpoint file exists before doing continuation.
933 * This way we can use identical input options for the first and subsequent runs...
935 if (mdrunOptions.numStepsCommandline > -2)
937 /* Temporarily set the number of steps to unmlimited to avoid
938 * triggering the nsteps check in load_checkpoint().
939 * This hack will go away soon when the -nsteps option is removed.
941 inputrec->nsteps = -1;
944 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
946 cr, domdecOptions.numCells,
947 inputrec, globalState.get(),
949 mdrunOptions.reproducible);
951 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
953 // Now we can start normal logging to the truncated log file.
954 fplog = gmx_fio_getfp(logFileHandle);
955 prepareLogAppending(fplog);
956 logOwner = buildLogger(fplog, cr);
957 mdlog = logOwner.logger();
961 if (mdrunOptions.numStepsCommandline > -2)
963 GMX_LOG(mdlog.info).asParagraph().
964 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
965 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
967 /* override nsteps with value set on the commamdline */
968 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
972 copy_mat(globalState->box, box);
977 gmx_bcast(sizeof(box), box, cr);
980 if (inputrec->cutoff_scheme != ecutsVERLET)
982 gmx_fatal(FARGS, "This group-scheme .tpr file can no longer be run by mdrun. Please update to the Verlet scheme, or use an earlier version of GROMACS if necessary.");
984 /* Update rlist and nstlist. */
985 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
986 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
988 LocalAtomSetManager atomSets;
989 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
990 inputrec->eI == eiNM))
992 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
994 box, positionsFromStatePointer(globalState.get()),
996 mdModulesNotifier.notify(&atomSets);
997 // Note that local state still does not exist yet.
1001 /* PME, if used, is done on all nodes with 1D decomposition */
1003 cr->duty = (DUTY_PP | DUTY_PME);
1005 if (inputrec->ePBC == epbcSCREW)
1008 "pbc=screw is only implemented with domain decomposition");
1014 /* After possible communicator splitting in make_dd_communicators.
1015 * we can set up the intra/inter node communication.
1017 gmx_setup_nodecomm(fplog, cr);
1023 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1024 "This is simulation %d out of %d running as a composite GROMACS\n"
1025 "multi-simulation job. Setup for this simulation:\n",
1028 GMX_LOG(mdlog.warning).appendTextFormatted(
1029 "Using %d MPI %s\n",
1032 cr->nnodes == 1 ? "thread" : "threads"
1034 cr->nnodes == 1 ? "process" : "processes"
1040 // If mdrun -pin auto honors any affinity setting that already
1041 // exists. If so, it is nice to provide feedback about whether
1042 // that existing affinity setting was from OpenMP or something
1043 // else, so we run this code both before and after we initialize
1044 // the OpenMP support.
1045 gmx_check_thread_affinity_set(mdlog,
1046 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1047 /* Check and update the number of OpenMP threads requested */
1048 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1051 gmx_omp_nthreads_init(mdlog, cr,
1052 hwinfo->nthreads_hw_avail,
1053 physicalNodeComm.size_,
1054 hw_opt.nthreads_omp,
1055 hw_opt.nthreads_omp_pme,
1056 !thisRankHasDuty(cr, DUTY_PP));
1058 // Enable FP exception detection, but not in
1059 // Release mode and not for compilers with known buggy FP
1060 // exception support (clang with any optimization) or suspected
1061 // buggy FP exception support (gcc 7.* with optimization).
1062 #if !defined NDEBUG && \
1063 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1064 && defined __OPTIMIZE__)
1065 const bool bEnableFPE = true;
1067 const bool bEnableFPE = false;
1069 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1072 gmx_feenableexcept();
1075 // Build a data structure that expresses which kinds of non-bonded
1076 // task are handled by this rank.
1078 // TODO Later, this might become a loop over all registered modules
1079 // relevant to the mdp inputs, to find those that have such tasks.
1081 // TODO This could move before init_domain_decomposition() as part
1082 // of refactoring that separates the responsibility for duty
1083 // assignment from setup for communication between tasks, and
1084 // setup for tasks handled with a domain (ie including short-ranged
1085 // tasks, bonded tasks, etc.).
1087 // Note that in general useGpuForNonbonded, etc. can have a value
1088 // that is inconsistent with the presence of actual GPUs on any
1089 // rank, and that is not known to be a problem until the
1090 // duty of the ranks on a node become known.
1092 // TODO Later we might need the concept of computeTasksOnThisRank,
1093 // from which we construct gpuTasksOnThisRank.
1095 // Currently the DD code assigns duty to ranks that can
1096 // include PP work that currently can be executed on a single
1097 // GPU, if present and compatible. This has to be coordinated
1098 // across PP ranks on a node, with possible multiple devices
1099 // or sharing devices on a node, either from the user
1100 // selection, or automatically.
1101 auto haveGpus = !gpuIdsToUse.empty();
1102 std::vector<GpuTask> gpuTasksOnThisRank;
1103 if (thisRankHasDuty(cr, DUTY_PP))
1105 if (useGpuForNonbonded)
1107 // Note that any bonded tasks on a GPU always accompany a
1111 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1113 else if (nonbondedTarget == TaskTarget::Gpu)
1115 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because no GPU is detected.");
1117 else if (bondedTarget == TaskTarget::Gpu)
1119 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because no GPU is detected.");
1123 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1124 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1130 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1132 else if (pmeTarget == TaskTarget::Gpu)
1134 gmx_fatal(FARGS, "Cannot run PME on a GPU because no GPU is detected.");
1139 GpuTaskAssignment gpuTaskAssignment;
1142 // Produce the task assignment for this rank.
1143 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1144 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1145 useGpuForBonded, pmeRunMode);
1147 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1149 /* Prevent other ranks from continuing after an issue was found
1150 * and reported as a fatal error.
1152 * TODO This function implements a barrier so that MPI runtimes
1153 * can organize an orderly shutdown if one of the ranks has had to
1154 * issue a fatal error in various code already run. When we have
1155 * MPI-aware error handling and reporting, this should be
1160 MPI_Barrier(cr->mpi_comm_mysim);
1166 MPI_Barrier(ms->mpi_comm_masters);
1168 /* We need another barrier to prevent non-master ranks from contiuing
1169 * when an error occured in a different simulation.
1171 MPI_Barrier(cr->mpi_comm_mysim);
1175 /* Now that we know the setup is consistent, check for efficiency */
1176 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1179 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1181 if (thisRankHasDuty(cr, DUTY_PP))
1183 // This works because only one task of each type is currently permitted.
1184 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1185 hasTaskType<GpuTask::Nonbonded>);
1186 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1188 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1189 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1190 init_gpu(nonbondedDeviceInfo);
1192 if (DOMAINDECOMP(cr))
1194 /* When we share GPUs over ranks, we need to know this for the DLB */
1195 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1201 std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1203 gmx_device_info_t *pmeDeviceInfo = nullptr;
1204 // Later, this program could contain kernels that might be later
1205 // re-used as auto-tuning progresses, or subsequent simulations
1207 PmeGpuProgramStorage pmeGpuProgram;
1208 // This works because only one task of each type is currently permitted.
1209 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1210 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1211 if (thisRankHasPmeGpuTask)
1213 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1214 init_gpu(pmeDeviceInfo);
1215 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1216 // TODO It would be nice to move this logic into the factory
1217 // function. See Redmine #2535.
1218 bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1219 if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1221 initializedClfftLibrary = initializeClfftLibrary();
1225 /* getting number of PP/PME threads on this MPI / tMPI rank.
1226 PME: env variable should be read only on one node to make sure it is
1227 identical everywhere;
1229 const int numThreadsOnThisRank =
1230 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1231 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1232 *hwinfo->hardwareTopology,
1233 physicalNodeComm, mdlog);
1235 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1237 /* Before setting affinity, check whether the affinity has changed
1238 * - which indicates that probably the OpenMP library has changed it
1239 * since we first checked).
1241 gmx_check_thread_affinity_set(mdlog,
1242 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1244 int numThreadsOnThisNode, intraNodeThreadOffset;
1245 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1246 &intraNodeThreadOffset);
1248 /* Set the CPU affinity */
1249 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1250 numThreadsOnThisRank, numThreadsOnThisNode,
1251 intraNodeThreadOffset, nullptr);
1254 if (mdrunOptions.timingOptions.resetStep > -1)
1256 GMX_LOG(mdlog.info).asParagraph().
1257 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1259 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1263 /* Master synchronizes its value of reset_counters with all nodes
1264 * including PME only nodes */
1265 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1266 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1267 wcycle_set_reset_counters(wcycle, reset_counters);
1270 // Membrane embedding must be initialized before we call init_forcerec()
1275 fprintf(stderr, "Initializing membed");
1277 /* Note that membed cannot work in parallel because mtop is
1278 * changed here. Fix this if we ever want to make it run with
1279 * multiple ranks. */
1280 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1282 .checkpointOptions.period);
1285 std::unique_ptr<MDAtoms> mdAtoms;
1286 std::unique_ptr<gmx_vsite_t> vsite;
1289 if (thisRankHasDuty(cr, DUTY_PP))
1291 mdModulesNotifier.notify(*cr);
1292 /* Initiate forcerecord */
1293 fr = new t_forcerec;
1294 fr->forceProviders = mdModules_->initForceProviders();
1295 init_forcerec(fplog, mdlog, fr, fcd,
1296 inputrec, &mtop, cr, box,
1297 opt2fn("-table", filenames.size(), filenames.data()),
1298 opt2fn("-tablep", filenames.size(), filenames.data()),
1299 opt2fns("-tableb", filenames.size(), filenames.data()),
1300 *hwinfo, nonbondedDeviceInfo,
1306 /* Initialize the mdAtoms structure.
1307 * mdAtoms is not filled with atom data,
1308 * as this can not be done now with domain decomposition.
1310 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1311 if (globalState && thisRankHasPmeGpuTask)
1313 // The pinning of coordinates in the global state object works, because we only use
1314 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1315 // points to the global state object without DD.
1316 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1317 // which should also perform the pinning.
1318 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1321 /* Initialize the virtual site communication */
1322 vsite = initVsite(mtop, cr);
1324 calc_shifts(box, fr->shift_vec);
1326 /* With periodic molecules the charge groups should be whole at start up
1327 * and the virtual sites should not be far from their proper positions.
1329 if (!inputrec->bContinuation && MASTER(cr) &&
1330 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1332 /* Make molecules whole at start of run */
1333 if (fr->ePBC != epbcNONE)
1335 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1339 /* Correct initial vsite positions are required
1340 * for the initial distribution in the domain decomposition
1341 * and for the initial shell prediction.
1343 constructVsitesGlobal(mtop, globalState->x);
1347 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1349 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1350 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1355 /* This is a PME only node */
1357 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1359 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1360 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1363 gmx_pme_t *sepPmeData = nullptr;
1364 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1365 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1366 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1368 /* Initiate PME if necessary,
1369 * either on all nodes or on dedicated PME nodes only. */
1370 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1372 if (mdAtoms && mdAtoms->mdatoms())
1374 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1375 if (EVDW_PME(inputrec->vdwtype))
1377 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1380 if (cr->npmenodes > 0)
1382 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1383 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1384 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1387 if (thisRankHasDuty(cr, DUTY_PME))
1391 pmedata = gmx_pme_init(cr,
1392 getNumPmeDomains(cr->dd),
1394 nChargePerturbed != 0, nTypePerturbed != 0,
1395 mdrunOptions.reproducible,
1396 ewaldcoeff_q, ewaldcoeff_lj,
1397 gmx_omp_nthreads_get(emntPME),
1398 pmeRunMode, nullptr,
1399 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1401 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1406 if (EI_DYNAMICS(inputrec->eI))
1408 /* Turn on signal handling on all nodes */
1410 * (A user signal from the PME nodes (if any)
1411 * is communicated to the PP nodes.
1413 signal_handler_install();
1416 pull_t *pull_work = nullptr;
1417 if (thisRankHasDuty(cr, DUTY_PP))
1419 /* Assumes uniform use of the number of OpenMP threads */
1420 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1422 if (inputrec->bPull)
1424 /* Initialize pull code */
1426 init_pull(fplog, inputrec->pull, inputrec,
1427 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1428 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1430 initPullHistory(pull_work, &observablesHistory);
1432 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1434 init_pull_output_files(pull_work,
1435 filenames.size(), filenames.data(), oenv,
1440 std::unique_ptr<EnforcedRotation> enforcedRotation;
1443 /* Initialize enforced rotation code */
1444 enforcedRotation = init_rot(fplog,
1457 t_swap *swap = nullptr;
1458 if (inputrec->eSwapCoords != eswapNO)
1460 /* Initialize ion swapping code */
1461 swap = init_swapcoords(fplog, inputrec,
1462 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1463 &mtop, globalState.get(), &observablesHistory,
1464 cr, &atomSets, oenv, mdrunOptions,
1468 /* Let makeConstraints know whether we have essential dynamics constraints.
1469 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1471 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1472 || observablesHistory.edsamHistory);
1473 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics,
1474 fplog, *mdAtoms->mdatoms(),
1475 cr, ms, &nrnb, wcycle, fr->bMolPBC);
1477 /* Energy terms and groups */
1478 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(), inputrec->fepvals->n_lambda);
1480 /* Kinetic energy data */
1481 gmx_ekindata_t ekind;
1482 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1484 /* Set up interactive MD (IMD) */
1485 auto imdSession = makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1486 MASTER(cr) ? globalState->x.rvec_array() : nullptr,
1487 filenames.size(), filenames.data(), oenv, mdrunOptions.imdOptions,
1490 if (DOMAINDECOMP(cr))
1492 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1493 /* This call is not included in init_domain_decomposition mainly
1494 * because fr->cginfo_mb is set later.
1496 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1497 domdecOptions.checkBondedInteractions,
1501 // TODO This is not the right place to manage the lifetime of
1502 // this data structure, but currently it's the easiest way to
1503 // make it work. Later, it should probably be made/updated
1504 // after the workload for the lifetime of a PP domain is
1506 PpForceWorkload ppForceWorkload;
1508 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1509 SimulatorBuilder simulatorBuilder;
1511 // build and run simulator object based on user-input
1512 auto simulator = simulatorBuilder.build(
1513 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1517 vsite.get(), constr.get(),
1518 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1520 mdModules_->outputProvider(),
1521 inputrec, imdSession.get(), pull_work, swap, &mtop,
1524 &observablesHistory,
1525 mdAtoms.get(), &nrnb, wcycle, fr,
1531 walltime_accounting,
1532 std::move(stopHandlerBuilder_),
1536 if (inputrec->bPull)
1538 finish_pull(pull_work);
1540 finish_swapcoords(swap);
1544 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1546 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1547 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1550 wallcycle_stop(wcycle, ewcRUN);
1552 /* Finish up, write some stuff
1553 * if rerunMD, don't write last frame again
1555 finish_run(fplog, mdlog, cr,
1556 inputrec, &nrnb, wcycle, walltime_accounting,
1557 fr ? fr->nbv.get() : nullptr,
1559 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1561 // clean up cycle counter
1562 wallcycle_destroy(wcycle);
1567 gmx_pme_destroy(pmedata);
1571 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1572 // before we destroy the GPU context(s) in free_gpu_resources().
1573 // Pinned buffers are associated with contexts in CUDA.
1574 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1575 mdAtoms.reset(nullptr);
1576 globalState.reset(nullptr);
1577 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1579 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1580 free_gpu_resources(fr, physicalNodeComm, hwinfo->gpu_info);
1581 free_gpu(nonbondedDeviceInfo);
1582 free_gpu(pmeDeviceInfo);
1583 done_forcerec(fr, mtop.molblock.size());
1588 free_membed(membed);
1591 /* Does what it says */
1592 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1593 walltime_accounting_destroy(walltime_accounting);
1595 // Ensure log file content is written
1598 gmx_fio_flush(logFileHandle);
1601 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1602 * exceptions were enabled before function was called. */
1605 gmx_fedisableexcept();
1608 auto rc = static_cast<int>(gmx_get_stop_condition());
1611 /* we need to join all threads. The sub-threads join when they
1612 exit this function, but the master thread needs to be told to
1614 if (PAR(cr) && MASTER(cr))
1618 //TODO free commrec in MPI simulations
1624 Mdrunner::~Mdrunner()
1626 // Clean up of the Manager.
1627 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1628 // but okay as long as threads synchronize some time before adding or accessing
1629 // a new set of restraints.
1630 if (restraintManager_)
1632 restraintManager_->clear();
1633 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1634 "restraints added during runner life time should be cleared at runner destruction.");
1638 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1639 const std::string &name)
1641 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1642 // Not sure if this should be logged through the md logger or something else,
1643 // but it is helpful to have some sort of INFO level message sent somewhere.
1644 // std::cout << "Registering restraint named " << name << std::endl;
1646 // When multiple restraints are used, it may be wasteful to register them separately.
1647 // Maybe instead register an entire Restraint Manager as a force provider.
1648 restraintManager_->addToSpec(std::move(puller),
1652 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1653 : mdModules_(std::move(mdModules))
1657 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1659 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1660 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1662 class Mdrunner::BuilderImplementation
1665 BuilderImplementation() = delete;
1666 BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1667 SimulationContext * context);
1668 ~BuilderImplementation();
1670 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1671 real forceWarningThreshold,
1672 StartingBehavior startingBehavior);
1674 void addDomdec(const DomdecOptions &options);
1676 void addVerletList(int nstlist);
1678 void addReplicaExchange(const ReplicaExchangeParameters ¶ms);
1680 void addMultiSim(gmx_multisim_t* multisim);
1682 void addNonBonded(const char* nbpu_opt);
1684 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1686 void addBondedTaskAssignment(const char* bonded_opt);
1688 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1690 void addFilenames(ArrayRef <const t_filenm> filenames);
1692 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1694 void addLogFile(t_fileio *logFileHandle);
1696 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1702 // Default parameters copied from runner.h
1703 // \todo Clarify source(s) of default parameters.
1705 const char* nbpu_opt_ = nullptr;
1706 const char* pme_opt_ = nullptr;
1707 const char* pme_fft_opt_ = nullptr;
1708 const char *bonded_opt_ = nullptr;
1710 MdrunOptions mdrunOptions_;
1712 DomdecOptions domdecOptions_;
1714 ReplicaExchangeParameters replicaExchangeParameters_;
1716 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1719 //! Non-owning multisim communicator handle.
1720 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1722 //! Print a warning if any force is larger than this (in kJ/mol nm).
1723 real forceWarningThreshold_ = -1;
1725 //! Whether the simulation will start afresh, or restart with/without appending.
1726 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1728 //! The modules that comprise the functionality of mdrun.
1729 std::unique_ptr<MDModules> mdModules_;
1731 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1734 * \todo Establish robust protocol to make sure resources remain valid.
1735 * SimulationContext will likely be separated into multiple layers for
1736 * different levels of access and different phases of execution. Ref
1737 * https://redmine.gromacs.org/issues/2375
1738 * https://redmine.gromacs.org/issues/2587
1740 SimulationContext* context_ = nullptr;
1742 //! \brief Parallelism information.
1743 gmx_hw_opt_t hardwareOptions_;
1745 //! filename options for simulation.
1746 ArrayRef<const t_filenm> filenames_;
1748 /*! \brief Handle to output environment.
1750 * \todo gmx_output_env_t needs lifetime management.
1752 gmx_output_env_t* outputEnvironment_ = nullptr;
1754 /*! \brief Non-owning handle to MD log file.
1756 * \todo Context should own output facilities for client.
1757 * \todo Improve log file handle management.
1759 * Code managing the FILE* relies on the ability to set it to
1760 * nullptr to check whether the filehandle is valid.
1762 t_fileio* logFileHandle_ = nullptr;
1765 * \brief Builder for simulation stop signal handler.
1767 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1770 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1771 SimulationContext * context) :
1772 mdModules_(std::move(mdModules)),
1775 GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1778 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1780 Mdrunner::BuilderImplementation &
1781 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1782 const real forceWarningThreshold,
1783 const StartingBehavior startingBehavior)
1785 mdrunOptions_ = options;
1786 forceWarningThreshold_ = forceWarningThreshold;
1787 startingBehavior_ = startingBehavior;
1791 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1793 domdecOptions_ = options;
1796 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1801 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1803 replicaExchangeParameters_ = params;
1806 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1808 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1811 Mdrunner Mdrunner::BuilderImplementation::build()
1813 auto newRunner = Mdrunner(std::move(mdModules_));
1815 GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1817 newRunner.mdrunOptions = mdrunOptions_;
1818 newRunner.pforce = forceWarningThreshold_;
1819 newRunner.startingBehavior = startingBehavior_;
1820 newRunner.domdecOptions = domdecOptions_;
1822 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1823 newRunner.hw_opt = hardwareOptions_;
1825 // No invariant to check. This parameter exists to optionally override other behavior.
1826 newRunner.nstlist_cmdline = nstlist_;
1828 newRunner.replExParams = replicaExchangeParameters_;
1830 newRunner.filenames = filenames_;
1832 GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1833 newRunner.cr = context_->communicationRecord_;
1837 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1838 newRunner.ms = *multisim_;
1842 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1845 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1846 // \todo Update sanity checking when output environment has clearly specified invariants.
1847 // Initialization and default values for oenv are not well specified in the current version.
1848 if (outputEnvironment_)
1850 newRunner.oenv = outputEnvironment_;
1854 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1857 newRunner.logFileHandle = logFileHandle_;
1861 newRunner.nbpu_opt = nbpu_opt_;
1865 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1868 if (pme_opt_ && pme_fft_opt_)
1870 newRunner.pme_opt = pme_opt_;
1871 newRunner.pme_fft_opt = pme_fft_opt_;
1875 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1880 newRunner.bonded_opt = bonded_opt_;
1884 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1887 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1889 if (stopHandlerBuilder_)
1891 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1895 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1901 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1903 nbpu_opt_ = nbpu_opt;
1906 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1907 const char* pme_fft_opt)
1910 pme_fft_opt_ = pme_fft_opt;
1913 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1915 bonded_opt_ = bonded_opt;
1918 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1920 hardwareOptions_ = hardwareOptions;
1923 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1925 filenames_ = filenames;
1928 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1930 outputEnvironment_ = outputEnvironment;
1933 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1935 logFileHandle_ = logFileHandle;
1938 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1940 stopHandlerBuilder_ = std::move(builder);
1943 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
1944 compat::not_null<SimulationContext*> context) :
1945 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context)}
1949 MdrunnerBuilder::~MdrunnerBuilder() = default;
1951 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1952 real forceWarningThreshold,
1953 const StartingBehavior startingBehavior)
1955 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
1959 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1961 impl_->addDomdec(options);
1965 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1967 impl_->addVerletList(nstlist);
1971 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1973 impl_->addReplicaExchange(params);
1977 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1979 impl_->addMultiSim(multisim);
1983 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1985 impl_->addNonBonded(nbpu_opt);
1989 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
1990 const char* pme_fft_opt)
1992 // The builder method may become more general in the future, but in this version,
1993 // parameters for PME electrostatics are both required and the only parameters
1995 if (pme_opt && pme_fft_opt)
1997 impl_->addPME(pme_opt, pme_fft_opt);
2001 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2006 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2008 impl_->addBondedTaskAssignment(bonded_opt);
2012 Mdrunner MdrunnerBuilder::build()
2014 return impl_->build();
2017 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2019 impl_->addHardwareOptions(hardwareOptions);
2023 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2025 impl_->addFilenames(filenames);
2029 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2031 impl_->addOutputEnvironment(outputEnvironment);
2035 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2037 impl_->addLogFile(logFileHandle);
2041 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2043 impl_->addStopHandlerBuilder(std::move(builder));
2047 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2049 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;