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/force.h"
89 #include "gromacs/mdlib/forcerec.h"
90 #include "gromacs/mdlib/gmx_omp_nthreads.h"
91 #include "gromacs/mdlib/makeconstraints.h"
92 #include "gromacs/mdlib/md_support.h"
93 #include "gromacs/mdlib/mdatoms.h"
94 #include "gromacs/mdlib/membed.h"
95 #include "gromacs/mdlib/ppforceworkload.h"
96 #include "gromacs/mdlib/qmmm.h"
97 #include "gromacs/mdlib/sighandler.h"
98 #include "gromacs/mdlib/stophandler.h"
99 #include "gromacs/mdrun/logging.h"
100 #include "gromacs/mdrun/multisim.h"
101 #include "gromacs/mdrun/simulationcontext.h"
102 #include "gromacs/mdrunutility/mdmodules.h"
103 #include "gromacs/mdrunutility/printtime.h"
104 #include "gromacs/mdrunutility/threadaffinity.h"
105 #include "gromacs/mdtypes/commrec.h"
106 #include "gromacs/mdtypes/enerdata.h"
107 #include "gromacs/mdtypes/fcdata.h"
108 #include "gromacs/mdtypes/inputrec.h"
109 #include "gromacs/mdtypes/md_enums.h"
110 #include "gromacs/mdtypes/mdrunoptions.h"
111 #include "gromacs/mdtypes/observableshistory.h"
112 #include "gromacs/mdtypes/state.h"
113 #include "gromacs/nbnxm/gpu_data_mgmt.h"
114 #include "gromacs/nbnxm/nbnxm.h"
115 #include "gromacs/nbnxm/pairlist_tuning.h"
116 #include "gromacs/pbcutil/pbc.h"
117 #include "gromacs/pulling/output.h"
118 #include "gromacs/pulling/pull.h"
119 #include "gromacs/pulling/pull_rotation.h"
120 #include "gromacs/restraint/manager.h"
121 #include "gromacs/restraint/restraintmdmodule.h"
122 #include "gromacs/restraint/restraintpotential.h"
123 #include "gromacs/swap/swapcoords.h"
124 #include "gromacs/taskassignment/decidegpuusage.h"
125 #include "gromacs/taskassignment/resourcedivision.h"
126 #include "gromacs/taskassignment/taskassignment.h"
127 #include "gromacs/taskassignment/usergpuids.h"
128 #include "gromacs/timing/gpu_timing.h"
129 #include "gromacs/timing/wallcycle.h"
130 #include "gromacs/timing/wallcyclereporting.h"
131 #include "gromacs/topology/mtop_util.h"
132 #include "gromacs/trajectory/trajectoryframe.h"
133 #include "gromacs/utility/basenetwork.h"
134 #include "gromacs/utility/cstringutil.h"
135 #include "gromacs/utility/exceptions.h"
136 #include "gromacs/utility/fatalerror.h"
137 #include "gromacs/utility/filestream.h"
138 #include "gromacs/utility/gmxassert.h"
139 #include "gromacs/utility/gmxmpi.h"
140 #include "gromacs/utility/logger.h"
141 #include "gromacs/utility/loggerbuilder.h"
142 #include "gromacs/utility/physicalnodecommunicator.h"
143 #include "gromacs/utility/pleasecite.h"
144 #include "gromacs/utility/programcontext.h"
145 #include "gromacs/utility/smalloc.h"
146 #include "gromacs/utility/stringutil.h"
148 #include "integrator.h"
149 #include "replicaexchange.h"
152 #include "corewrap.h"
158 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
160 * Used to ensure that the master thread does not modify mdrunner during copy
161 * on the spawned threads. */
162 static void threadMpiMdrunnerAccessBarrier()
165 MPI_Barrier(MPI_COMM_WORLD);
169 Mdrunner Mdrunner::cloneOnSpawnedThread() const
171 auto newRunner = Mdrunner(std::make_unique<MDModules>());
173 // All runners in the same process share a restraint manager resource because it is
174 // part of the interface to the client code, which is associated only with the
175 // original thread. Handles to the same resources can be obtained by copy.
177 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
180 // Copy original cr pointer before master thread can pass the thread barrier
181 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
183 // Copy members of master runner.
184 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
185 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
186 newRunner.hw_opt = hw_opt;
187 newRunner.filenames = filenames;
189 newRunner.oenv = oenv;
190 newRunner.mdrunOptions = mdrunOptions;
191 newRunner.domdecOptions = domdecOptions;
192 newRunner.nbpu_opt = nbpu_opt;
193 newRunner.pme_opt = pme_opt;
194 newRunner.pme_fft_opt = pme_fft_opt;
195 newRunner.bonded_opt = bonded_opt;
196 newRunner.nstlist_cmdline = nstlist_cmdline;
197 newRunner.replExParams = replExParams;
198 newRunner.pforce = pforce;
200 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
202 threadMpiMdrunnerAccessBarrier();
204 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
209 /*! \brief The callback used for running on spawned threads.
211 * Obtains the pointer to the master mdrunner object from the one
212 * argument permitted to the thread-launch API call, copies it to make
213 * a new runner for this thread, reinitializes necessary data, and
214 * proceeds to the simulation. */
215 static void mdrunner_start_fn(const void *arg)
219 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
220 /* copy the arg list to make sure that it's thread-local. This
221 doesn't copy pointed-to items, of course; fnm, cr and fplog
222 are reset in the call below, all others should be const. */
223 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
226 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
230 /*! \brief Start thread-MPI threads.
232 * Called by mdrunner() to start a specific number of threads
233 * (including the main thread) for thread-parallel runs. This in turn
234 * calls mdrunner() for each thread. All options are the same as for
236 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
239 /* first check whether we even need to start tMPI */
240 if (numThreadsToLaunch < 2)
246 /* now spawn new threads that start mdrunner_start_fn(), while
247 the main thread returns, we set thread affinity later */
248 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
249 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
251 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
254 threadMpiMdrunnerAccessBarrier();
256 GMX_UNUSED_VALUE(mdrunner_start_fn);
259 return reinitialize_commrec_for_this_thread(cr);
264 /*! \brief Initialize variables for Verlet scheme simulation */
265 static void prepare_verlet_scheme(FILE *fplog,
269 const gmx_mtop_t *mtop,
271 bool makeGpuPairList,
272 const gmx::CpuInfo &cpuinfo)
274 /* For NVE simulations, we will retain the initial list buffer */
275 if (EI_DYNAMICS(ir->eI) &&
276 ir->verletbuf_tol > 0 &&
277 !(EI_MD(ir->eI) && ir->etc == etcNO))
279 /* Update the Verlet buffer size for the current run setup */
281 /* Here we assume SIMD-enabled kernels are being used. But as currently
282 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
283 * and 4x2 gives a larger buffer than 4x4, this is ok.
285 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
286 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
288 const real rlist_new =
289 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
291 if (rlist_new != ir->rlist)
293 if (fplog != nullptr)
295 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
296 ir->rlist, rlist_new,
297 listSetup.cluster_size_i, listSetup.cluster_size_j);
299 ir->rlist = rlist_new;
303 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
305 gmx_fatal(FARGS, "Can not set nstlist without %s",
306 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
309 if (EI_DYNAMICS(ir->eI))
311 /* Set or try nstlist values */
312 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
316 /*! \brief Override the nslist value in inputrec
318 * with value passed on the command line (if any)
320 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
321 int64_t nsteps_cmdline,
326 /* override with anything else than the default -2 */
327 if (nsteps_cmdline > -2)
329 char sbuf_steps[STEPSTRSIZE];
330 char sbuf_msg[STRLEN];
332 ir->nsteps = nsteps_cmdline;
333 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
335 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
336 gmx_step_str(nsteps_cmdline, sbuf_steps),
337 fabs(nsteps_cmdline*ir->delta_t));
341 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
342 gmx_step_str(nsteps_cmdline, sbuf_steps));
345 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
347 else if (nsteps_cmdline < -2)
349 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
352 /* Do nothing if nsteps_cmdline == -2 */
358 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
360 * If not, and if a warning may be issued, logs a warning about
361 * falling back to CPU code. With thread-MPI, only the first
362 * call to this function should have \c issueWarning true. */
363 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
364 const t_inputrec *ir,
367 if (ir->opts.ngener - ir->nwall > 1)
369 /* The GPU code does not support more than one energy group.
370 * If the user requested GPUs explicitly, a fatal error is given later.
374 GMX_LOG(mdlog.warning).asParagraph()
375 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
376 "For better performance, run on the GPU without energy groups and then do "
377 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
384 //! Initializes the logger for mdrun.
385 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
387 gmx::LoggerBuilder builder;
388 if (fplog != nullptr)
390 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
392 if (cr == nullptr || SIMMASTER(cr))
394 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
395 &gmx::TextOutputFile::standardError());
397 return builder.build();
400 //! Make a TaskTarget from an mdrun argument string.
401 static TaskTarget findTaskTarget(const char *optionString)
403 TaskTarget returnValue = TaskTarget::Auto;
405 if (strncmp(optionString, "auto", 3) == 0)
407 returnValue = TaskTarget::Auto;
409 else if (strncmp(optionString, "cpu", 3) == 0)
411 returnValue = TaskTarget::Cpu;
413 else if (strncmp(optionString, "gpu", 3) == 0)
415 returnValue = TaskTarget::Gpu;
419 GMX_ASSERT(false, "Option string should have been checked for sanity already");
425 //! Finish run, aggregate data to print performance info.
426 static void finish_run(FILE *fplog,
427 const gmx::MDLogger &mdlog,
429 const t_inputrec *inputrec,
430 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
431 gmx_walltime_accounting_t walltime_accounting,
432 nonbonded_verlet_t *nbv,
433 const gmx_pme_t *pme,
436 t_nrnb *nrnb_tot = nullptr;
438 double nbfs = 0, mflop = 0;
440 elapsed_time_over_all_ranks,
441 elapsed_time_over_all_threads,
442 elapsed_time_over_all_threads_over_all_ranks;
443 /* Control whether it is valid to print a report. Only the
444 simulation master may print, but it should not do so if the run
445 terminated e.g. before a scheduled reset step. This is
446 complicated by the fact that PME ranks are unaware of the
447 reason why they were sent a pmerecvqxFINISH. To avoid
448 communication deadlocks, we always do the communication for the
449 report, even if we've decided not to write the report, because
450 how long it takes to finish the run is not important when we've
451 decided not to report on the simulation performance.
453 Further, we only report performance for dynamical integrators,
454 because those are the only ones for which we plan to
455 consider doing any optimizations. */
456 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
458 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
460 GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
468 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
477 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
478 elapsed_time_over_all_threads = walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
482 /* reduce elapsed_time over all MPI ranks in the current simulation */
483 MPI_Allreduce(&elapsed_time,
484 &elapsed_time_over_all_ranks,
485 1, MPI_DOUBLE, MPI_SUM,
487 elapsed_time_over_all_ranks /= cr->nnodes;
488 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
489 * current simulation. */
490 MPI_Allreduce(&elapsed_time_over_all_threads,
491 &elapsed_time_over_all_threads_over_all_ranks,
492 1, MPI_DOUBLE, MPI_SUM,
498 elapsed_time_over_all_ranks = elapsed_time;
499 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
504 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
511 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
513 print_dd_statistics(cr, inputrec, fplog);
516 /* TODO Move the responsibility for any scaling by thread counts
517 * to the code that handled the thread region, so that there's a
518 * mechanism to keep cycle counting working during the transition
519 * to task parallelism. */
520 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
521 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
522 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
523 auto cycle_sum(wallcycle_sum(cr, wcycle));
527 auto nbnxn_gpu_timings = (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
528 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
530 if (pme_gpu_task_enabled(pme))
532 pme_gpu_get_timings(pme, &pme_gpu_timings);
534 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
535 elapsed_time_over_all_ranks,
540 if (EI_DYNAMICS(inputrec->eI))
542 delta_t = inputrec->delta_t;
547 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
548 elapsed_time_over_all_ranks,
549 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
550 delta_t, nbfs, mflop);
554 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
555 elapsed_time_over_all_ranks,
556 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
557 delta_t, nbfs, mflop);
562 int Mdrunner::mdrunner()
566 t_forcerec *fr = nullptr;
567 t_fcdata *fcd = nullptr;
568 real ewaldcoeff_q = 0;
569 real ewaldcoeff_lj = 0;
570 int nChargePerturbed = -1, nTypePerturbed = 0;
571 gmx_wallcycle_t wcycle;
572 gmx_walltime_accounting_t walltime_accounting = nullptr;
573 gmx_membed_t * membed = nullptr;
574 gmx_hw_info_t *hwinfo = nullptr;
576 /* CAUTION: threads may be started later on in this function, so
577 cr doesn't reflect the final parallel state right now */
578 t_inputrec inputrecInstance;
579 t_inputrec *inputrec = &inputrecInstance;
582 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
583 bool doRerun = mdrunOptions.rerun;
585 // Handle task-assignment related user options.
586 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
587 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
588 std::vector<int> gpuIdsAvailable;
591 gpuIdsAvailable = parseUserGpuIdString(hw_opt.gpuIdsAvailable);
593 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
595 std::vector<int> userGpuTaskAssignment;
598 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
600 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
601 auto nonbondedTarget = findTaskTarget(nbpu_opt);
602 auto pmeTarget = findTaskTarget(pme_opt);
603 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
604 auto bondedTarget = findTaskTarget(bonded_opt);
605 PmeRunMode pmeRunMode = PmeRunMode::None;
607 // Here we assume that SIMMASTER(cr) does not change even after the
608 // threads are started.
610 FILE *fplog = nullptr;
611 // If we are appending, we don't write log output because we need
612 // to check that the old log file matches what the checkpoint file
613 // expects. Otherwise, we should start to write log output now if
614 // there is a file ready for it.
615 if (logFileHandle != nullptr && !mdrunOptions.continuationOptions.appendFiles)
617 fplog = gmx_fio_getfp(logFileHandle);
619 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
620 gmx::MDLogger mdlog(logOwner.logger());
622 // TODO The thread-MPI master rank makes a working
623 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
624 // after the threads have been launched. This works because no use
625 // is made of that communicator until after the execution paths
626 // have rejoined. But it is likely that we can improve the way
627 // this is expressed, e.g. by expressly running detection only the
628 // master rank for thread-MPI, rather than relying on the mutex
629 // and reference count.
630 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
631 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
633 gmx_print_detected_hardware(fplog, cr, ms, mdlog, hwinfo);
635 std::vector<int> gpuIdsToUse;
636 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
637 if (gpuIdsAvailable.empty())
639 gpuIdsToUse = compatibleGpus;
643 for (const auto &availableGpuId : gpuIdsAvailable)
645 bool availableGpuIsCompatible = false;
646 for (const auto &compatibleGpuId : compatibleGpus)
648 if (availableGpuId == compatibleGpuId)
650 availableGpuIsCompatible = true;
654 if (!availableGpuIsCompatible)
656 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);
658 gpuIdsToUse.push_back(availableGpuId);
662 if (fplog != nullptr)
664 /* Print references after all software/hardware printing */
665 please_cite(fplog, "Abraham2015");
666 please_cite(fplog, "Pall2015");
667 please_cite(fplog, "Pronk2013");
668 please_cite(fplog, "Hess2008b");
669 please_cite(fplog, "Spoel2005a");
670 please_cite(fplog, "Lindahl2001a");
671 please_cite(fplog, "Berendsen95a");
672 writeSourceDoi(fplog);
675 std::unique_ptr<t_state> globalState;
679 /* Only the master rank has the global state */
680 globalState = std::make_unique<t_state>();
682 /* Read (nearly) all data required for the simulation */
683 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
685 /* In rerun, set velocities to zero if present */
686 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
688 // rerun does not use velocities
689 GMX_LOG(mdlog.info).asParagraph().appendText(
690 "Rerun trajectory contains velocities. Rerun does only evaluate "
691 "potential energy and forces. The velocities will be ignored.");
692 for (int i = 0; i < globalState->natoms; i++)
694 clear_rvec(globalState->v[i]);
696 globalState->flags &= ~(1 << estV);
699 if (inputrec->cutoff_scheme != ecutsVERLET)
701 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.");
705 /* Check and update the hardware options for internal consistency */
706 check_and_update_hw_opt_1(mdlog, &hw_opt, cr, domdecOptions.numPmeRanks);
708 /* Early check for externally set process affinity. */
709 gmx_check_thread_affinity_set(mdlog, cr,
710 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
712 if (GMX_THREAD_MPI && SIMMASTER(cr))
714 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
716 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
719 /* Since the master knows the cut-off scheme, update hw_opt for this.
720 * This is done later for normal MPI and also once more with tMPI
721 * for all tMPI ranks.
723 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
725 bool useGpuForNonbonded = false;
726 bool useGpuForPme = false;
729 // If the user specified the number of ranks, then we must
730 // respect that, but in default mode, we need to allow for
731 // the number of GPUs to choose the number of ranks.
732 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
733 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
734 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
735 canUseGpuForNonbonded,
736 inputrec->cutoff_scheme == ecutsVERLET,
737 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
738 hw_opt.nthreads_tmpi);
739 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
740 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
741 *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
744 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
746 /* Determine how many thread-MPI ranks to start.
748 * TODO Over-writing the user-supplied value here does
749 * prevent any possible subsequent checks from working
751 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
760 // Now start the threads for thread MPI.
761 cr = spawnThreads(hw_opt.nthreads_tmpi);
762 /* The main thread continues here with a new cr. We don't deallocate
763 the old cr because other threads may still be reading it. */
764 // TODO Both master and spawned threads call dup_tfn and
765 // reinitialize_commrec_for_this_thread. Find a way to express
767 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
769 // END OF CAUTION: cr and physicalNodeComm are now reliable
773 /* now broadcast everything to the non-master nodes/threads: */
774 init_parallel(cr, inputrec, &mtop);
777 // Now each rank knows the inputrec that SIMMASTER read and used,
778 // and (if applicable) cr->nnodes has been assigned the number of
779 // thread-MPI ranks that have been chosen. The ranks can now all
780 // run the task-deciding functions and will agree on the result
781 // without needing to communicate.
783 // TODO Should we do the communication in debug mode to support
784 // having an assertion?
786 // Note that these variables describe only their own node.
788 // Note that when bonded interactions run on a GPU they always run
789 // alongside a nonbonded task, so do not influence task assignment
790 // even though they affect the force calculation workload.
791 bool useGpuForNonbonded = false;
792 bool useGpuForPme = false;
793 bool useGpuForBonded = false;
796 // It's possible that there are different numbers of GPUs on
797 // different nodes, which is the user's responsibilty to
798 // handle. If unsuitable, we will notice that during task
800 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
801 bool usingVerletScheme = inputrec->cutoff_scheme == ecutsVERLET;
802 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
803 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
805 canUseGpuForNonbonded,
807 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
809 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
810 *hwinfo, *inputrec, mtop,
811 cr->nnodes, domdecOptions.numPmeRanks,
813 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
815 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
816 bondedTarget, canUseGpuForBonded,
817 EVDW_PME(inputrec->vdwtype),
818 EEL_PME_EWALD(inputrec->coulombtype),
819 domdecOptions.numPmeRanks, gpusWereDetected);
821 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
822 if (pmeRunMode == PmeRunMode::GPU)
824 if (pmeFftTarget == TaskTarget::Cpu)
826 pmeRunMode = PmeRunMode::Mixed;
829 else if (pmeFftTarget == TaskTarget::Gpu)
831 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.");
834 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
837 // TODO: hide restraint implementation details from Mdrunner.
838 // There is nothing unique about restraints at this point as far as the
839 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
840 // factory functions from the SimulationContext on which to call mdModules_->add().
841 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
842 for (auto && restraint : restraintManager_->getRestraints())
844 auto module = RestraintMDModule::create(restraint,
846 mdModules_->add(std::move(module));
849 // TODO: Error handling
850 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
852 if (fplog != nullptr)
854 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
855 fprintf(fplog, "\n");
860 /* now make sure the state is initialized and propagated */
861 set_state_entries(globalState.get(), inputrec);
864 /* NM and TPI parallelize over force/energy calculations, not atoms,
865 * so we need to initialize and broadcast the global state.
867 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
871 globalState = std::make_unique<t_state>();
873 broadcastStateWithoutDynamics(cr, globalState.get());
876 /* A parallel command line option consistency check that we can
877 only do after any threads have started. */
878 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
879 domdecOptions.numCells[YY] > 1 ||
880 domdecOptions.numCells[ZZ] > 1 ||
881 domdecOptions.numPmeRanks > 0))
884 "The -dd or -npme option request a parallel simulation, "
886 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
889 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
891 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
897 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
899 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
902 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
904 if (domdecOptions.numPmeRanks > 0)
906 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
907 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
910 domdecOptions.numPmeRanks = 0;
913 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
915 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
916 * improve performance with many threads per GPU, since our OpenMP
917 * scaling is bad, but it's difficult to automate the setup.
919 domdecOptions.numPmeRanks = 0;
923 if (domdecOptions.numPmeRanks < 0)
925 domdecOptions.numPmeRanks = 0;
926 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
930 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
937 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
941 /* NMR restraints must be initialized before load_checkpoint,
942 * since with time averaging the history is added to t_state.
943 * For proper consistency check we therefore need to extend
945 * So the PME-only nodes (if present) will also initialize
946 * the distance restraints.
950 /* This needs to be called before read_checkpoint to extend the state */
951 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
953 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
955 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
957 ObservablesHistory observablesHistory = {};
959 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
961 if (continuationOptions.startedFromCheckpoint)
963 /* Check if checkpoint file exists before doing continuation.
964 * This way we can use identical input options for the first and subsequent runs...
968 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
970 cr, domdecOptions.numCells,
971 inputrec, globalState.get(),
972 &bReadEkin, &observablesHistory,
973 continuationOptions.appendFiles,
974 continuationOptions.appendFilesOptionSet,
975 mdrunOptions.reproducible);
979 continuationOptions.haveReadEkin = true;
982 if (continuationOptions.appendFiles && logFileHandle)
984 // Now we can start normal logging to the truncated log file.
985 fplog = gmx_fio_getfp(logFileHandle);
986 prepareLogAppending(fplog);
987 logOwner = buildLogger(fplog, cr);
988 mdlog = logOwner.logger();
992 if (mdrunOptions.numStepsCommandline > -2)
994 GMX_LOG(mdlog.info).asParagraph().
995 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
996 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
998 /* override nsteps with value set on the commamdline */
999 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1003 copy_mat(globalState->box, box);
1008 gmx_bcast(sizeof(box), box, cr);
1011 /* Update rlist and nstlist. */
1012 if (inputrec->cutoff_scheme == ecutsVERLET)
1014 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1015 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
1018 LocalAtomSetManager atomSets;
1020 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1021 inputrec->eI == eiNM))
1023 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
1025 box, positionsFromStatePointer(globalState.get()),
1027 // Note that local state still does not exist yet.
1031 /* PME, if used, is done on all nodes with 1D decomposition */
1033 cr->duty = (DUTY_PP | DUTY_PME);
1035 if (inputrec->ePBC == epbcSCREW)
1038 "pbc=screw is only implemented with domain decomposition");
1044 /* After possible communicator splitting in make_dd_communicators.
1045 * we can set up the intra/inter node communication.
1047 gmx_setup_nodecomm(fplog, cr);
1053 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1054 "This is simulation %d out of %d running as a composite GROMACS\n"
1055 "multi-simulation job. Setup for this simulation:\n",
1058 GMX_LOG(mdlog.warning).appendTextFormatted(
1059 "Using %d MPI %s\n",
1062 cr->nnodes == 1 ? "thread" : "threads"
1064 cr->nnodes == 1 ? "process" : "processes"
1070 /* Check and update hw_opt for the cut-off scheme */
1071 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
1073 /* Check and update the number of OpenMP threads requested */
1074 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1077 gmx_omp_nthreads_init(mdlog, cr,
1078 hwinfo->nthreads_hw_avail,
1079 physicalNodeComm.size_,
1080 hw_opt.nthreads_omp,
1081 hw_opt.nthreads_omp_pme,
1082 !thisRankHasDuty(cr, DUTY_PP),
1083 inputrec->cutoff_scheme == ecutsVERLET);
1085 // Enable FP exception detection for the Verlet scheme, but not in
1086 // Release mode and not for compilers with known buggy FP
1087 // exception support (clang with any optimization) or suspected
1088 // buggy FP exception support (gcc 7.* with optimization).
1089 #if !defined NDEBUG && \
1090 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1091 && defined __OPTIMIZE__)
1092 const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
1094 const bool bEnableFPE = false;
1096 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1099 gmx_feenableexcept();
1102 // Build a data structure that expresses which kinds of non-bonded
1103 // task are handled by this rank.
1105 // TODO Later, this might become a loop over all registered modules
1106 // relevant to the mdp inputs, to find those that have such tasks.
1108 // TODO This could move before init_domain_decomposition() as part
1109 // of refactoring that separates the responsibility for duty
1110 // assignment from setup for communication between tasks, and
1111 // setup for tasks handled with a domain (ie including short-ranged
1112 // tasks, bonded tasks, etc.).
1114 // Note that in general useGpuForNonbonded, etc. can have a value
1115 // that is inconsistent with the presence of actual GPUs on any
1116 // rank, and that is not known to be a problem until the
1117 // duty of the ranks on a node become known.
1119 // TODO Later we might need the concept of computeTasksOnThisRank,
1120 // from which we construct gpuTasksOnThisRank.
1122 // Currently the DD code assigns duty to ranks that can
1123 // include PP work that currently can be executed on a single
1124 // GPU, if present and compatible. This has to be coordinated
1125 // across PP ranks on a node, with possible multiple devices
1126 // or sharing devices on a node, either from the user
1127 // selection, or automatically.
1128 auto haveGpus = !gpuIdsToUse.empty();
1129 std::vector<GpuTask> gpuTasksOnThisRank;
1130 if (thisRankHasDuty(cr, DUTY_PP))
1132 if (useGpuForNonbonded)
1134 // Note that any bonded tasks on a GPU always accompany a
1138 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1140 else if (nonbondedTarget == TaskTarget::Gpu)
1142 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1144 else if (bondedTarget == TaskTarget::Gpu)
1146 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because there is none detected.");
1150 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1151 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1157 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1159 else if (pmeTarget == TaskTarget::Gpu)
1161 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1166 GpuTaskAssignment gpuTaskAssignment;
1169 // Produce the task assignment for this rank.
1170 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1171 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1172 useGpuForBonded, pmeRunMode);
1174 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1176 /* Prevent other ranks from continuing after an issue was found
1177 * and reported as a fatal error.
1179 * TODO This function implements a barrier so that MPI runtimes
1180 * can organize an orderly shutdown if one of the ranks has had to
1181 * issue a fatal error in various code already run. When we have
1182 * MPI-aware error handling and reporting, this should be
1187 MPI_Barrier(cr->mpi_comm_mysim);
1193 MPI_Barrier(ms->mpi_comm_masters);
1195 /* We need another barrier to prevent non-master ranks from contiuing
1196 * when an error occured in a different simulation.
1198 MPI_Barrier(cr->mpi_comm_mysim);
1202 /* Now that we know the setup is consistent, check for efficiency */
1203 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1206 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1208 if (thisRankHasDuty(cr, DUTY_PP))
1210 // This works because only one task of each type is currently permitted.
1211 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1212 hasTaskType<GpuTask::Nonbonded>);
1213 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1215 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1216 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1217 init_gpu(nonbondedDeviceInfo);
1219 if (DOMAINDECOMP(cr))
1221 /* When we share GPUs over ranks, we need to know this for the DLB */
1222 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1228 std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1230 gmx_device_info_t *pmeDeviceInfo = nullptr;
1231 // Later, this program could contain kernels that might be later
1232 // re-used as auto-tuning progresses, or subsequent simulations
1234 PmeGpuProgramStorage pmeGpuProgram;
1235 // This works because only one task of each type is currently permitted.
1236 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1237 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1238 if (thisRankHasPmeGpuTask)
1240 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1241 init_gpu(pmeDeviceInfo);
1242 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1243 // TODO It would be nice to move this logic into the factory
1244 // function. See Redmine #2535.
1245 bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1246 if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1248 initializedClfftLibrary = initializeClfftLibrary();
1252 /* getting number of PP/PME threads on this MPI / tMPI rank.
1253 PME: env variable should be read only on one node to make sure it is
1254 identical everywhere;
1256 const int numThreadsOnThisRank =
1257 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1258 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1259 *hwinfo->hardwareTopology,
1260 physicalNodeComm, mdlog);
1262 if (hw_opt.thread_affinity != threadaffOFF)
1264 /* Before setting affinity, check whether the affinity has changed
1265 * - which indicates that probably the OpenMP library has changed it
1266 * since we first checked).
1268 gmx_check_thread_affinity_set(mdlog, cr,
1269 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1271 int numThreadsOnThisNode, intraNodeThreadOffset;
1272 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1273 &intraNodeThreadOffset);
1275 /* Set the CPU affinity */
1276 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1277 numThreadsOnThisRank, numThreadsOnThisNode,
1278 intraNodeThreadOffset, nullptr);
1281 if (mdrunOptions.timingOptions.resetStep > -1)
1283 GMX_LOG(mdlog.info).asParagraph().
1284 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1286 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1290 /* Master synchronizes its value of reset_counters with all nodes
1291 * including PME only nodes */
1292 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1293 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1294 wcycle_set_reset_counters(wcycle, reset_counters);
1297 // Membrane embedding must be initialized before we call init_forcerec()
1302 fprintf(stderr, "Initializing membed");
1304 /* Note that membed cannot work in parallel because mtop is
1305 * changed here. Fix this if we ever want to make it run with
1306 * multiple ranks. */
1307 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1309 .checkpointOptions.period);
1312 std::unique_ptr<MDAtoms> mdAtoms;
1313 std::unique_ptr<gmx_vsite_t> vsite;
1316 if (thisRankHasDuty(cr, DUTY_PP))
1318 /* Initiate forcerecord */
1319 fr = new t_forcerec;
1320 fr->forceProviders = mdModules_->initForceProviders();
1321 init_forcerec(fplog, mdlog, fr, fcd,
1322 inputrec, &mtop, cr, box,
1323 opt2fn("-table", filenames.size(), filenames.data()),
1324 opt2fn("-tablep", filenames.size(), filenames.data()),
1325 opt2fns("-tableb", filenames.size(), filenames.data()),
1326 *hwinfo, nonbondedDeviceInfo,
1331 /* Initialize the mdAtoms structure.
1332 * mdAtoms is not filled with atom data,
1333 * as this can not be done now with domain decomposition.
1335 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1336 if (globalState && thisRankHasPmeGpuTask)
1338 // The pinning of coordinates in the global state object works, because we only use
1339 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1340 // points to the global state object without DD.
1341 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1342 // which should also perform the pinning.
1343 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1346 /* Initialize the virtual site communication */
1347 vsite = initVsite(mtop, cr);
1349 calc_shifts(box, fr->shift_vec);
1351 /* With periodic molecules the charge groups should be whole at start up
1352 * and the virtual sites should not be far from their proper positions.
1354 if (!inputrec->bContinuation && MASTER(cr) &&
1355 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1357 /* Make molecules whole at start of run */
1358 if (fr->ePBC != epbcNONE)
1360 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1364 /* Correct initial vsite positions are required
1365 * for the initial distribution in the domain decomposition
1366 * and for the initial shell prediction.
1368 constructVsitesGlobal(mtop, globalState->x);
1372 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1374 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1375 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1380 /* This is a PME only node */
1382 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1384 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1385 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1388 gmx_pme_t *sepPmeData = nullptr;
1389 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1390 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1391 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1393 /* Initiate PME if necessary,
1394 * either on all nodes or on dedicated PME nodes only. */
1395 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1397 if (mdAtoms && mdAtoms->mdatoms())
1399 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1400 if (EVDW_PME(inputrec->vdwtype))
1402 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1405 if (cr->npmenodes > 0)
1407 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1408 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1409 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1412 if (thisRankHasDuty(cr, DUTY_PME))
1416 pmedata = gmx_pme_init(cr,
1417 getNumPmeDomains(cr->dd),
1419 mtop.natoms, nChargePerturbed != 0, nTypePerturbed != 0,
1420 mdrunOptions.reproducible,
1421 ewaldcoeff_q, ewaldcoeff_lj,
1422 gmx_omp_nthreads_get(emntPME),
1423 pmeRunMode, nullptr,
1424 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1426 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1431 if (EI_DYNAMICS(inputrec->eI))
1433 /* Turn on signal handling on all nodes */
1435 * (A user signal from the PME nodes (if any)
1436 * is communicated to the PP nodes.
1438 signal_handler_install();
1441 if (thisRankHasDuty(cr, DUTY_PP))
1443 /* Assumes uniform use of the number of OpenMP threads */
1444 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1446 if (inputrec->bPull)
1448 /* Initialize pull code */
1449 inputrec->pull_work =
1450 init_pull(fplog, inputrec->pull, inputrec,
1451 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1452 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1454 initPullHistory(inputrec->pull_work, &observablesHistory);
1456 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1458 init_pull_output_files(inputrec->pull_work,
1459 filenames.size(), filenames.data(), oenv,
1460 continuationOptions);
1464 std::unique_ptr<EnforcedRotation> enforcedRotation;
1467 /* Initialize enforced rotation code */
1468 enforcedRotation = init_rot(fplog,
1480 if (inputrec->eSwapCoords != eswapNO)
1482 /* Initialize ion swapping code */
1483 init_swapcoords(fplog, inputrec, opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1484 &mtop, globalState.get(), &observablesHistory,
1485 cr, &atomSets, oenv, mdrunOptions);
1488 /* Let makeConstraints know whether we have essential dynamics constraints.
1489 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1491 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1492 || observablesHistory.edsamHistory);
1493 auto constr = makeConstraints(mtop, *inputrec, doEssentialDynamics,
1494 fplog, *mdAtoms->mdatoms(),
1495 cr, ms, nrnb, wcycle, fr->bMolPBC);
1497 /* Energy terms and groups */
1498 gmx_enerdata_t *enerd;
1500 init_enerdata(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].nr, inputrec->fepvals->n_lambda, enerd);
1502 /* Set up interactive MD (IMD) */
1503 auto imdSession = makeImdSession(inputrec, cr, wcycle, enerd, ms, &mtop, mdlog,
1504 MASTER(cr) ? globalState->x.rvec_array() : nullptr,
1505 filenames.size(), filenames.data(), oenv, mdrunOptions);
1507 if (DOMAINDECOMP(cr))
1509 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1510 /* This call is not included in init_domain_decomposition mainly
1511 * because fr->cginfo_mb is set later.
1513 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1514 domdecOptions.checkBondedInteractions,
1518 // TODO This is not the right place to manage the lifetime of
1519 // this data structure, but currently it's the easiest way to
1520 // make it work. Later, it should probably be made/updated
1521 // after the workload for the lifetime of a PP domain is
1523 PpForceWorkload ppForceWorkload;
1525 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to integrator.");
1526 /* Now do whatever the user wants us to do (how flexible...) */
1527 Integrator integrator {
1528 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1531 vsite.get(), constr.get(),
1532 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1534 mdModules_->outputProvider(),
1535 inputrec, imdSession.get(), &mtop,
1538 &observablesHistory,
1539 mdAtoms.get(), nrnb, wcycle, fr,
1544 walltime_accounting,
1545 std::move(stopHandlerBuilder_)
1547 integrator.run(inputrec->eI, doRerun);
1549 if (inputrec->bPull)
1551 finish_pull(inputrec->pull_work);
1553 destroy_enerdata(enerd);
1558 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1560 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1561 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1564 wallcycle_stop(wcycle, ewcRUN);
1566 /* Finish up, write some stuff
1567 * if rerunMD, don't write last frame again
1569 finish_run(fplog, mdlog, cr,
1570 inputrec, nrnb, wcycle, walltime_accounting,
1571 fr ? fr->nbv.get() : nullptr,
1573 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1575 // clean up cycle counter
1576 wallcycle_destroy(wcycle);
1581 gmx_pme_destroy(pmedata);
1585 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1586 // before we destroy the GPU context(s) in free_gpu_resources().
1587 // Pinned buffers are associated with contexts in CUDA.
1588 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1589 mdAtoms.reset(nullptr);
1590 globalState.reset(nullptr);
1591 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1593 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1594 free_gpu_resources(fr, physicalNodeComm);
1595 free_gpu(nonbondedDeviceInfo);
1596 free_gpu(pmeDeviceInfo);
1597 done_forcerec(fr, mtop.molblock.size(), mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].nr);
1602 free_membed(membed);
1605 gmx_hardware_info_free();
1607 /* Does what it says */
1608 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1609 walltime_accounting_destroy(walltime_accounting);
1612 // Ensure log file content is written
1615 gmx_fio_flush(logFileHandle);
1618 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1619 * exceptions were enabled before function was called. */
1622 gmx_fedisableexcept();
1625 auto rc = static_cast<int>(gmx_get_stop_condition());
1628 /* we need to join all threads. The sub-threads join when they
1629 exit this function, but the master thread needs to be told to
1631 if (PAR(cr) && MASTER(cr))
1635 //TODO free commrec in MPI simulations
1641 Mdrunner::~Mdrunner()
1643 // Clean up of the Manager.
1644 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1645 // but okay as long as threads synchronize some time before adding or accessing
1646 // a new set of restraints.
1647 if (restraintManager_)
1649 restraintManager_->clear();
1650 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1651 "restraints added during runner life time should be cleared at runner destruction.");
1655 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1658 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1659 // Not sure if this should be logged through the md logger or something else,
1660 // but it is helpful to have some sort of INFO level message sent somewhere.
1661 // std::cout << "Registering restraint named " << name << std::endl;
1663 // When multiple restraints are used, it may be wasteful to register them separately.
1664 // Maybe instead register an entire Restraint Manager as a force provider.
1665 restraintManager_->addToSpec(std::move(puller),
1669 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1670 : mdModules_(std::move(mdModules))
1674 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1676 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1677 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1679 class Mdrunner::BuilderImplementation
1682 BuilderImplementation() = delete;
1683 BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1684 SimulationContext * context);
1685 ~BuilderImplementation();
1687 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1688 real forceWarningThreshold);
1690 void addDomdec(const DomdecOptions &options);
1692 void addVerletList(int nstlist);
1694 void addReplicaExchange(const ReplicaExchangeParameters ¶ms);
1696 void addMultiSim(gmx_multisim_t* multisim);
1698 void addNonBonded(const char* nbpu_opt);
1700 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1702 void addBondedTaskAssignment(const char* bonded_opt);
1704 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1706 void addFilenames(ArrayRef <const t_filenm> filenames);
1708 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1710 void addLogFile(t_fileio *logFileHandle);
1712 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1718 // Default parameters copied from runner.h
1719 // \todo Clarify source(s) of default parameters.
1721 const char* nbpu_opt_ = nullptr;
1722 const char* pme_opt_ = nullptr;
1723 const char* pme_fft_opt_ = nullptr;
1724 const char *bonded_opt_ = nullptr;
1726 MdrunOptions mdrunOptions_;
1728 DomdecOptions domdecOptions_;
1730 ReplicaExchangeParameters replicaExchangeParameters_;
1732 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1735 //! Non-owning multisim communicator handle.
1736 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1738 //! Print a warning if any force is larger than this (in kJ/mol nm).
1739 real forceWarningThreshold_ = -1;
1741 //! The modules that comprise the functionality of mdrun.
1742 std::unique_ptr<MDModules> mdModules_;
1744 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1747 * \todo Establish robust protocol to make sure resources remain valid.
1748 * SimulationContext will likely be separated into multiple layers for
1749 * different levels of access and different phases of execution. Ref
1750 * https://redmine.gromacs.org/issues/2375
1751 * https://redmine.gromacs.org/issues/2587
1753 SimulationContext* context_ = nullptr;
1755 //! \brief Parallelism information.
1756 gmx_hw_opt_t hardwareOptions_;
1758 //! filename options for simulation.
1759 ArrayRef<const t_filenm> filenames_;
1761 /*! \brief Handle to output environment.
1763 * \todo gmx_output_env_t needs lifetime management.
1765 gmx_output_env_t* outputEnvironment_ = nullptr;
1767 /*! \brief Non-owning handle to MD log file.
1769 * \todo Context should own output facilities for client.
1770 * \todo Improve log file handle management.
1772 * Code managing the FILE* relies on the ability to set it to
1773 * nullptr to check whether the filehandle is valid.
1775 t_fileio* logFileHandle_ = nullptr;
1778 * \brief Builder for simulation stop signal handler.
1780 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1783 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1784 SimulationContext * context) :
1785 mdModules_(std::move(mdModules)),
1788 GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1791 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1793 Mdrunner::BuilderImplementation &
1794 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1795 real forceWarningThreshold)
1797 mdrunOptions_ = options;
1798 forceWarningThreshold_ = forceWarningThreshold;
1802 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1804 domdecOptions_ = options;
1807 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1812 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1814 replicaExchangeParameters_ = params;
1817 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1819 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1822 Mdrunner Mdrunner::BuilderImplementation::build()
1824 auto newRunner = Mdrunner(std::move(mdModules_));
1826 GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1828 newRunner.mdrunOptions = mdrunOptions_;
1829 newRunner.domdecOptions = domdecOptions_;
1831 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1832 newRunner.hw_opt = hardwareOptions_;
1834 // No invariant to check. This parameter exists to optionally override other behavior.
1835 newRunner.nstlist_cmdline = nstlist_;
1837 newRunner.replExParams = replicaExchangeParameters_;
1839 newRunner.filenames = filenames_;
1841 GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1842 newRunner.cr = context_->communicationRecord_;
1846 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1847 newRunner.ms = *multisim_;
1851 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1854 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1855 // \todo Update sanity checking when output environment has clearly specified invariants.
1856 // Initialization and default values for oenv are not well specified in the current version.
1857 if (outputEnvironment_)
1859 newRunner.oenv = outputEnvironment_;
1863 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1866 newRunner.logFileHandle = logFileHandle_;
1870 newRunner.nbpu_opt = nbpu_opt_;
1874 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1877 if (pme_opt_ && pme_fft_opt_)
1879 newRunner.pme_opt = pme_opt_;
1880 newRunner.pme_fft_opt = pme_fft_opt_;
1884 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1889 newRunner.bonded_opt = bonded_opt_;
1893 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1896 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1898 if (stopHandlerBuilder_)
1900 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1904 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1910 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1912 nbpu_opt_ = nbpu_opt;
1915 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1916 const char* pme_fft_opt)
1919 pme_fft_opt_ = pme_fft_opt;
1922 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1924 bonded_opt_ = bonded_opt;
1927 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1929 hardwareOptions_ = hardwareOptions;
1932 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1934 filenames_ = filenames;
1937 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1939 outputEnvironment_ = outputEnvironment;
1942 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1944 logFileHandle_ = logFileHandle;
1947 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1949 stopHandlerBuilder_ = std::move(builder);
1952 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
1953 compat::not_null<SimulationContext*> context) :
1954 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context)}
1958 MdrunnerBuilder::~MdrunnerBuilder() = default;
1960 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1961 real forceWarningThreshold)
1963 impl_->setExtraMdrunOptions(options, forceWarningThreshold);
1967 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1969 impl_->addDomdec(options);
1973 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1975 impl_->addVerletList(nstlist);
1979 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1981 impl_->addReplicaExchange(params);
1985 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1987 impl_->addMultiSim(multisim);
1991 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1993 impl_->addNonBonded(nbpu_opt);
1997 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
1998 const char* pme_fft_opt)
2000 // The builder method may become more general in the future, but in this version,
2001 // parameters for PME electrostatics are both required and the only parameters
2003 if (pme_opt && pme_fft_opt)
2005 impl_->addPME(pme_opt, pme_fft_opt);
2009 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2014 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2016 impl_->addBondedTaskAssignment(bonded_opt);
2020 Mdrunner MdrunnerBuilder::build()
2022 return impl_->build();
2025 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2027 impl_->addHardwareOptions(hardwareOptions);
2031 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2033 impl_->addFilenames(filenames);
2037 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2039 impl_->addOutputEnvironment(outputEnvironment);
2043 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2045 impl_->addLogFile(logFileHandle);
2049 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2051 impl_->addStopHandlerBuilder(std::move(builder));
2055 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2057 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;