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/logging.h"
104 #include "gromacs/mdrunutility/multisim.h"
105 #include "gromacs/mdrunutility/printtime.h"
106 #include "gromacs/mdrunutility/threadaffinity.h"
107 #include "gromacs/mdtypes/commrec.h"
108 #include "gromacs/mdtypes/enerdata.h"
109 #include "gromacs/mdtypes/fcdata.h"
110 #include "gromacs/mdtypes/inputrec.h"
111 #include "gromacs/mdtypes/md_enums.h"
112 #include "gromacs/mdtypes/mdrunoptions.h"
113 #include "gromacs/mdtypes/observableshistory.h"
114 #include "gromacs/mdtypes/state.h"
115 #include "gromacs/nbnxm/gpu_data_mgmt.h"
116 #include "gromacs/nbnxm/nbnxm.h"
117 #include "gromacs/nbnxm/pairlist_tuning.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/output.h"
120 #include "gromacs/pulling/pull.h"
121 #include "gromacs/pulling/pull_rotation.h"
122 #include "gromacs/restraint/manager.h"
123 #include "gromacs/restraint/restraintmdmodule.h"
124 #include "gromacs/restraint/restraintpotential.h"
125 #include "gromacs/swap/swapcoords.h"
126 #include "gromacs/taskassignment/decidegpuusage.h"
127 #include "gromacs/taskassignment/resourcedivision.h"
128 #include "gromacs/taskassignment/taskassignment.h"
129 #include "gromacs/taskassignment/usergpuids.h"
130 #include "gromacs/timing/gpu_timing.h"
131 #include "gromacs/timing/wallcycle.h"
132 #include "gromacs/timing/wallcyclereporting.h"
133 #include "gromacs/topology/mtop_util.h"
134 #include "gromacs/trajectory/trajectoryframe.h"
135 #include "gromacs/utility/basenetwork.h"
136 #include "gromacs/utility/cstringutil.h"
137 #include "gromacs/utility/exceptions.h"
138 #include "gromacs/utility/fatalerror.h"
139 #include "gromacs/utility/filestream.h"
140 #include "gromacs/utility/gmxassert.h"
141 #include "gromacs/utility/gmxmpi.h"
142 #include "gromacs/utility/logger.h"
143 #include "gromacs/utility/loggerbuilder.h"
144 #include "gromacs/utility/physicalnodecommunicator.h"
145 #include "gromacs/utility/pleasecite.h"
146 #include "gromacs/utility/programcontext.h"
147 #include "gromacs/utility/smalloc.h"
148 #include "gromacs/utility/stringutil.h"
150 #include "integrator.h"
151 #include "replicaexchange.h"
154 #include "corewrap.h"
160 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
162 * Used to ensure that the master thread does not modify mdrunner during copy
163 * on the spawned threads. */
164 static void threadMpiMdrunnerAccessBarrier()
167 MPI_Barrier(MPI_COMM_WORLD);
171 Mdrunner Mdrunner::cloneOnSpawnedThread() const
173 auto newRunner = Mdrunner(std::make_unique<MDModules>());
175 // All runners in the same process share a restraint manager resource because it is
176 // part of the interface to the client code, which is associated only with the
177 // original thread. Handles to the same resources can be obtained by copy.
179 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
182 // Copy original cr pointer before master thread can pass the thread barrier
183 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
185 // Copy members of master runner.
186 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
187 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
188 newRunner.hw_opt = hw_opt;
189 newRunner.filenames = filenames;
191 newRunner.oenv = oenv;
192 newRunner.mdrunOptions = mdrunOptions;
193 newRunner.domdecOptions = domdecOptions;
194 newRunner.nbpu_opt = nbpu_opt;
195 newRunner.pme_opt = pme_opt;
196 newRunner.pme_fft_opt = pme_fft_opt;
197 newRunner.bonded_opt = bonded_opt;
198 newRunner.nstlist_cmdline = nstlist_cmdline;
199 newRunner.replExParams = replExParams;
200 newRunner.pforce = pforce;
202 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
204 threadMpiMdrunnerAccessBarrier();
206 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
211 /*! \brief The callback used for running on spawned threads.
213 * Obtains the pointer to the master mdrunner object from the one
214 * argument permitted to the thread-launch API call, copies it to make
215 * a new runner for this thread, reinitializes necessary data, and
216 * proceeds to the simulation. */
217 static void mdrunner_start_fn(const void *arg)
221 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
222 /* copy the arg list to make sure that it's thread-local. This
223 doesn't copy pointed-to items, of course; fnm, cr and fplog
224 are reset in the call below, all others should be const. */
225 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
228 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
232 /*! \brief Start thread-MPI threads.
234 * Called by mdrunner() to start a specific number of threads
235 * (including the main thread) for thread-parallel runs. This in turn
236 * calls mdrunner() for each thread. All options are the same as for
238 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
241 /* first check whether we even need to start tMPI */
242 if (numThreadsToLaunch < 2)
248 /* now spawn new threads that start mdrunner_start_fn(), while
249 the main thread returns, we set thread affinity later */
250 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
251 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
253 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
256 threadMpiMdrunnerAccessBarrier();
258 GMX_UNUSED_VALUE(mdrunner_start_fn);
261 return reinitialize_commrec_for_this_thread(cr);
266 /*! \brief Initialize variables for Verlet scheme simulation */
267 static void prepare_verlet_scheme(FILE *fplog,
271 const gmx_mtop_t *mtop,
273 bool makeGpuPairList,
274 const gmx::CpuInfo &cpuinfo)
276 /* For NVE simulations, we will retain the initial list buffer */
277 if (EI_DYNAMICS(ir->eI) &&
278 ir->verletbuf_tol > 0 &&
279 !(EI_MD(ir->eI) && ir->etc == etcNO))
281 /* Update the Verlet buffer size for the current run setup */
283 /* Here we assume SIMD-enabled kernels are being used. But as currently
284 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
285 * and 4x2 gives a larger buffer than 4x4, this is ok.
287 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
288 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
290 const real rlist_new =
291 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
293 if (rlist_new != ir->rlist)
295 if (fplog != nullptr)
297 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
298 ir->rlist, rlist_new,
299 listSetup.cluster_size_i, listSetup.cluster_size_j);
301 ir->rlist = rlist_new;
305 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
307 gmx_fatal(FARGS, "Can not set nstlist without %s",
308 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
311 if (EI_DYNAMICS(ir->eI))
313 /* Set or try nstlist values */
314 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
318 /*! \brief Override the nslist value in inputrec
320 * with value passed on the command line (if any)
322 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
323 int64_t nsteps_cmdline,
328 /* override with anything else than the default -2 */
329 if (nsteps_cmdline > -2)
331 char sbuf_steps[STEPSTRSIZE];
332 char sbuf_msg[STRLEN];
334 ir->nsteps = nsteps_cmdline;
335 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
337 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
338 gmx_step_str(nsteps_cmdline, sbuf_steps),
339 fabs(nsteps_cmdline*ir->delta_t));
343 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
344 gmx_step_str(nsteps_cmdline, sbuf_steps));
347 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
349 else if (nsteps_cmdline < -2)
351 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
354 /* Do nothing if nsteps_cmdline == -2 */
360 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
362 * If not, and if a warning may be issued, logs a warning about
363 * falling back to CPU code. With thread-MPI, only the first
364 * call to this function should have \c issueWarning true. */
365 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
366 const t_inputrec *ir,
369 if (ir->opts.ngener - ir->nwall > 1)
371 /* The GPU code does not support more than one energy group.
372 * If the user requested GPUs explicitly, a fatal error is given later.
376 GMX_LOG(mdlog.warning).asParagraph()
377 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
378 "For better performance, run on the GPU without energy groups and then do "
379 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
386 //! Initializes the logger for mdrun.
387 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
389 gmx::LoggerBuilder builder;
390 if (fplog != nullptr)
392 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
394 if (cr == nullptr || SIMMASTER(cr))
396 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
397 &gmx::TextOutputFile::standardError());
399 return builder.build();
402 //! Make a TaskTarget from an mdrun argument string.
403 static TaskTarget findTaskTarget(const char *optionString)
405 TaskTarget returnValue = TaskTarget::Auto;
407 if (strncmp(optionString, "auto", 3) == 0)
409 returnValue = TaskTarget::Auto;
411 else if (strncmp(optionString, "cpu", 3) == 0)
413 returnValue = TaskTarget::Cpu;
415 else if (strncmp(optionString, "gpu", 3) == 0)
417 returnValue = TaskTarget::Gpu;
421 GMX_ASSERT(false, "Option string should have been checked for sanity already");
427 //! Finish run, aggregate data to print performance info.
428 static void finish_run(FILE *fplog,
429 const gmx::MDLogger &mdlog,
431 const t_inputrec *inputrec,
432 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
433 gmx_walltime_accounting_t walltime_accounting,
434 nonbonded_verlet_t *nbv,
435 const gmx_pme_t *pme,
438 t_nrnb *nrnb_tot = nullptr;
440 double nbfs = 0, mflop = 0;
442 elapsed_time_over_all_ranks,
443 elapsed_time_over_all_threads,
444 elapsed_time_over_all_threads_over_all_ranks;
445 /* Control whether it is valid to print a report. Only the
446 simulation master may print, but it should not do so if the run
447 terminated e.g. before a scheduled reset step. This is
448 complicated by the fact that PME ranks are unaware of the
449 reason why they were sent a pmerecvqxFINISH. To avoid
450 communication deadlocks, we always do the communication for the
451 report, even if we've decided not to write the report, because
452 how long it takes to finish the run is not important when we've
453 decided not to report on the simulation performance.
455 Further, we only report performance for dynamical integrators,
456 because those are the only ones for which we plan to
457 consider doing any optimizations. */
458 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
460 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
462 GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
470 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
479 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
480 elapsed_time_over_all_threads = walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
484 /* reduce elapsed_time over all MPI ranks in the current simulation */
485 MPI_Allreduce(&elapsed_time,
486 &elapsed_time_over_all_ranks,
487 1, MPI_DOUBLE, MPI_SUM,
489 elapsed_time_over_all_ranks /= cr->nnodes;
490 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
491 * current simulation. */
492 MPI_Allreduce(&elapsed_time_over_all_threads,
493 &elapsed_time_over_all_threads_over_all_ranks,
494 1, MPI_DOUBLE, MPI_SUM,
500 elapsed_time_over_all_ranks = elapsed_time;
501 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
506 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
513 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
515 print_dd_statistics(cr, inputrec, fplog);
518 /* TODO Move the responsibility for any scaling by thread counts
519 * to the code that handled the thread region, so that there's a
520 * mechanism to keep cycle counting working during the transition
521 * to task parallelism. */
522 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
523 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
524 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
525 auto cycle_sum(wallcycle_sum(cr, wcycle));
529 auto nbnxn_gpu_timings = (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
530 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
532 if (pme_gpu_task_enabled(pme))
534 pme_gpu_get_timings(pme, &pme_gpu_timings);
536 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
537 elapsed_time_over_all_ranks,
542 if (EI_DYNAMICS(inputrec->eI))
544 delta_t = inputrec->delta_t;
549 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
550 elapsed_time_over_all_ranks,
551 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
552 delta_t, nbfs, mflop);
556 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
557 elapsed_time_over_all_ranks,
558 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
559 delta_t, nbfs, mflop);
564 int Mdrunner::mdrunner()
568 t_forcerec *fr = nullptr;
569 t_fcdata *fcd = nullptr;
570 real ewaldcoeff_q = 0;
571 real ewaldcoeff_lj = 0;
572 int nChargePerturbed = -1, nTypePerturbed = 0;
573 gmx_wallcycle_t wcycle;
574 gmx_walltime_accounting_t walltime_accounting = nullptr;
575 gmx_membed_t * membed = nullptr;
576 gmx_hw_info_t *hwinfo = nullptr;
578 /* CAUTION: threads may be started later on in this function, so
579 cr doesn't reflect the final parallel state right now */
580 t_inputrec inputrecInstance;
581 t_inputrec *inputrec = &inputrecInstance;
584 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
585 bool doRerun = mdrunOptions.rerun;
587 // Handle task-assignment related user options.
588 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
589 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
590 std::vector<int> gpuIdsAvailable;
593 gpuIdsAvailable = parseUserGpuIdString(hw_opt.gpuIdsAvailable);
595 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
597 std::vector<int> userGpuTaskAssignment;
600 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
602 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
603 auto nonbondedTarget = findTaskTarget(nbpu_opt);
604 auto pmeTarget = findTaskTarget(pme_opt);
605 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
606 auto bondedTarget = findTaskTarget(bonded_opt);
607 PmeRunMode pmeRunMode = PmeRunMode::None;
609 // Here we assume that SIMMASTER(cr) does not change even after the
610 // threads are started.
612 FILE *fplog = nullptr;
613 // If we are appending, we don't write log output because we need
614 // to check that the old log file matches what the checkpoint file
615 // expects. Otherwise, we should start to write log output now if
616 // there is a file ready for it.
617 if (logFileHandle != nullptr && !mdrunOptions.continuationOptions.appendFiles)
619 fplog = gmx_fio_getfp(logFileHandle);
621 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
622 gmx::MDLogger mdlog(logOwner.logger());
624 // TODO The thread-MPI master rank makes a working
625 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
626 // after the threads have been launched. This works because no use
627 // is made of that communicator until after the execution paths
628 // have rejoined. But it is likely that we can improve the way
629 // this is expressed, e.g. by expressly running detection only the
630 // master rank for thread-MPI, rather than relying on the mutex
631 // and reference count.
632 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
633 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
635 gmx_print_detected_hardware(fplog, isMasterSimMasterRank(ms, MASTER(cr)), mdlog, hwinfo);
637 std::vector<int> gpuIdsToUse;
638 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
639 if (gpuIdsAvailable.empty())
641 gpuIdsToUse = compatibleGpus;
645 for (const auto &availableGpuId : gpuIdsAvailable)
647 bool availableGpuIsCompatible = false;
648 for (const auto &compatibleGpuId : compatibleGpus)
650 if (availableGpuId == compatibleGpuId)
652 availableGpuIsCompatible = true;
656 if (!availableGpuIsCompatible)
658 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);
660 gpuIdsToUse.push_back(availableGpuId);
664 if (fplog != nullptr)
666 /* Print references after all software/hardware printing */
667 please_cite(fplog, "Abraham2015");
668 please_cite(fplog, "Pall2015");
669 please_cite(fplog, "Pronk2013");
670 please_cite(fplog, "Hess2008b");
671 please_cite(fplog, "Spoel2005a");
672 please_cite(fplog, "Lindahl2001a");
673 please_cite(fplog, "Berendsen95a");
674 writeSourceDoi(fplog);
677 std::unique_ptr<t_state> globalState;
681 /* Only the master rank has the global state */
682 globalState = std::make_unique<t_state>();
684 /* Read (nearly) all data required for the simulation */
685 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
687 /* In rerun, set velocities to zero if present */
688 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
690 // rerun does not use velocities
691 GMX_LOG(mdlog.info).asParagraph().appendText(
692 "Rerun trajectory contains velocities. Rerun does only evaluate "
693 "potential energy and forces. The velocities will be ignored.");
694 for (int i = 0; i < globalState->natoms; i++)
696 clear_rvec(globalState->v[i]);
698 globalState->flags &= ~(1 << estV);
701 if (inputrec->cutoff_scheme != ecutsVERLET)
703 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.");
707 /* Check and update the hardware options for internal consistency */
708 check_and_update_hw_opt_1(mdlog, &hw_opt, cr, domdecOptions.numPmeRanks);
710 /* Early check for externally set process affinity. */
711 gmx_check_thread_affinity_set(mdlog, cr,
712 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
714 if (GMX_THREAD_MPI && SIMMASTER(cr))
716 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
718 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
721 /* Since the master knows the cut-off scheme, update hw_opt for this.
722 * This is done later for normal MPI and also once more with tMPI
723 * for all tMPI ranks.
725 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
727 bool useGpuForNonbonded = false;
728 bool useGpuForPme = false;
731 // If the user specified the number of ranks, then we must
732 // respect that, but in default mode, we need to allow for
733 // the number of GPUs to choose the number of ranks.
734 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
735 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
736 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
737 canUseGpuForNonbonded,
738 inputrec->cutoff_scheme == ecutsVERLET,
739 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
740 hw_opt.nthreads_tmpi);
741 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
742 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
743 *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
746 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
748 /* Determine how many thread-MPI ranks to start.
750 * TODO Over-writing the user-supplied value here does
751 * prevent any possible subsequent checks from working
753 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
762 // Now start the threads for thread MPI.
763 cr = spawnThreads(hw_opt.nthreads_tmpi);
764 /* The main thread continues here with a new cr. We don't deallocate
765 the old cr because other threads may still be reading it. */
766 // TODO Both master and spawned threads call dup_tfn and
767 // reinitialize_commrec_for_this_thread. Find a way to express
769 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
771 // END OF CAUTION: cr and physicalNodeComm are now reliable
775 /* now broadcast everything to the non-master nodes/threads: */
776 init_parallel(cr, inputrec, &mtop);
779 // Now each rank knows the inputrec that SIMMASTER read and used,
780 // and (if applicable) cr->nnodes has been assigned the number of
781 // thread-MPI ranks that have been chosen. The ranks can now all
782 // run the task-deciding functions and will agree on the result
783 // without needing to communicate.
785 // TODO Should we do the communication in debug mode to support
786 // having an assertion?
788 // Note that these variables describe only their own node.
790 // Note that when bonded interactions run on a GPU they always run
791 // alongside a nonbonded task, so do not influence task assignment
792 // even though they affect the force calculation workload.
793 bool useGpuForNonbonded = false;
794 bool useGpuForPme = false;
795 bool useGpuForBonded = false;
798 // It's possible that there are different numbers of GPUs on
799 // different nodes, which is the user's responsibilty to
800 // handle. If unsuitable, we will notice that during task
802 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
803 bool usingVerletScheme = inputrec->cutoff_scheme == ecutsVERLET;
804 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
805 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
807 canUseGpuForNonbonded,
809 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
811 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
812 *hwinfo, *inputrec, mtop,
813 cr->nnodes, domdecOptions.numPmeRanks,
815 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
817 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
818 bondedTarget, canUseGpuForBonded,
819 EVDW_PME(inputrec->vdwtype),
820 EEL_PME_EWALD(inputrec->coulombtype),
821 domdecOptions.numPmeRanks, gpusWereDetected);
823 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
824 if (pmeRunMode == PmeRunMode::GPU)
826 if (pmeFftTarget == TaskTarget::Cpu)
828 pmeRunMode = PmeRunMode::Mixed;
831 else if (pmeFftTarget == TaskTarget::Gpu)
833 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.");
836 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
839 // TODO: hide restraint implementation details from Mdrunner.
840 // There is nothing unique about restraints at this point as far as the
841 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
842 // factory functions from the SimulationContext on which to call mdModules_->add().
843 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
844 for (auto && restraint : restraintManager_->getRestraints())
846 auto module = RestraintMDModule::create(restraint,
848 mdModules_->add(std::move(module));
851 // TODO: Error handling
852 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
854 if (fplog != nullptr)
856 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
857 fprintf(fplog, "\n");
862 /* now make sure the state is initialized and propagated */
863 set_state_entries(globalState.get(), inputrec);
866 /* NM and TPI parallelize over force/energy calculations, not atoms,
867 * so we need to initialize and broadcast the global state.
869 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
873 globalState = std::make_unique<t_state>();
875 broadcastStateWithoutDynamics(cr, globalState.get());
878 /* A parallel command line option consistency check that we can
879 only do after any threads have started. */
880 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
881 domdecOptions.numCells[YY] > 1 ||
882 domdecOptions.numCells[ZZ] > 1 ||
883 domdecOptions.numPmeRanks > 0))
886 "The -dd or -npme option request a parallel simulation, "
888 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
891 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
893 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
899 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
901 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
904 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
906 if (domdecOptions.numPmeRanks > 0)
908 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
909 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
912 domdecOptions.numPmeRanks = 0;
915 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
917 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
918 * improve performance with many threads per GPU, since our OpenMP
919 * scaling is bad, but it's difficult to automate the setup.
921 domdecOptions.numPmeRanks = 0;
925 if (domdecOptions.numPmeRanks < 0)
927 domdecOptions.numPmeRanks = 0;
928 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
932 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
939 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
943 /* NMR restraints must be initialized before load_checkpoint,
944 * since with time averaging the history is added to t_state.
945 * For proper consistency check we therefore need to extend
947 * So the PME-only nodes (if present) will also initialize
948 * the distance restraints.
952 /* This needs to be called before read_checkpoint to extend the state */
953 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
955 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
957 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
959 ObservablesHistory observablesHistory = {};
961 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
963 if (continuationOptions.startedFromCheckpoint)
965 /* Check if checkpoint file exists before doing continuation.
966 * This way we can use identical input options for the first and subsequent runs...
970 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
972 cr, domdecOptions.numCells,
973 inputrec, globalState.get(),
974 &bReadEkin, &observablesHistory,
975 continuationOptions.appendFiles,
976 continuationOptions.appendFilesOptionSet,
977 mdrunOptions.reproducible);
981 continuationOptions.haveReadEkin = true;
984 if (continuationOptions.appendFiles && logFileHandle)
986 // Now we can start normal logging to the truncated log file.
987 fplog = gmx_fio_getfp(logFileHandle);
988 prepareLogAppending(fplog);
989 logOwner = buildLogger(fplog, cr);
990 mdlog = logOwner.logger();
994 if (mdrunOptions.numStepsCommandline > -2)
996 GMX_LOG(mdlog.info).asParagraph().
997 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
998 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
1000 /* override nsteps with value set on the commamdline */
1001 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1005 copy_mat(globalState->box, box);
1010 gmx_bcast(sizeof(box), box, cr);
1013 /* Update rlist and nstlist. */
1014 if (inputrec->cutoff_scheme == ecutsVERLET)
1016 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1017 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
1020 LocalAtomSetManager atomSets;
1022 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1023 inputrec->eI == eiNM))
1025 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
1027 box, positionsFromStatePointer(globalState.get()),
1029 // Note that local state still does not exist yet.
1033 /* PME, if used, is done on all nodes with 1D decomposition */
1035 cr->duty = (DUTY_PP | DUTY_PME);
1037 if (inputrec->ePBC == epbcSCREW)
1040 "pbc=screw is only implemented with domain decomposition");
1046 /* After possible communicator splitting in make_dd_communicators.
1047 * we can set up the intra/inter node communication.
1049 gmx_setup_nodecomm(fplog, cr);
1055 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1056 "This is simulation %d out of %d running as a composite GROMACS\n"
1057 "multi-simulation job. Setup for this simulation:\n",
1060 GMX_LOG(mdlog.warning).appendTextFormatted(
1061 "Using %d MPI %s\n",
1064 cr->nnodes == 1 ? "thread" : "threads"
1066 cr->nnodes == 1 ? "process" : "processes"
1072 /* Check and update hw_opt for the cut-off scheme */
1073 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
1075 /* Check and update the number of OpenMP threads requested */
1076 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1079 gmx_omp_nthreads_init(mdlog, cr,
1080 hwinfo->nthreads_hw_avail,
1081 physicalNodeComm.size_,
1082 hw_opt.nthreads_omp,
1083 hw_opt.nthreads_omp_pme,
1084 !thisRankHasDuty(cr, DUTY_PP),
1085 inputrec->cutoff_scheme == ecutsVERLET);
1087 // Enable FP exception detection for the Verlet scheme, but not in
1088 // Release mode and not for compilers with known buggy FP
1089 // exception support (clang with any optimization) or suspected
1090 // buggy FP exception support (gcc 7.* with optimization).
1091 #if !defined NDEBUG && \
1092 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1093 && defined __OPTIMIZE__)
1094 const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
1096 const bool bEnableFPE = false;
1098 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1101 gmx_feenableexcept();
1104 // Build a data structure that expresses which kinds of non-bonded
1105 // task are handled by this rank.
1107 // TODO Later, this might become a loop over all registered modules
1108 // relevant to the mdp inputs, to find those that have such tasks.
1110 // TODO This could move before init_domain_decomposition() as part
1111 // of refactoring that separates the responsibility for duty
1112 // assignment from setup for communication between tasks, and
1113 // setup for tasks handled with a domain (ie including short-ranged
1114 // tasks, bonded tasks, etc.).
1116 // Note that in general useGpuForNonbonded, etc. can have a value
1117 // that is inconsistent with the presence of actual GPUs on any
1118 // rank, and that is not known to be a problem until the
1119 // duty of the ranks on a node become known.
1121 // TODO Later we might need the concept of computeTasksOnThisRank,
1122 // from which we construct gpuTasksOnThisRank.
1124 // Currently the DD code assigns duty to ranks that can
1125 // include PP work that currently can be executed on a single
1126 // GPU, if present and compatible. This has to be coordinated
1127 // across PP ranks on a node, with possible multiple devices
1128 // or sharing devices on a node, either from the user
1129 // selection, or automatically.
1130 auto haveGpus = !gpuIdsToUse.empty();
1131 std::vector<GpuTask> gpuTasksOnThisRank;
1132 if (thisRankHasDuty(cr, DUTY_PP))
1134 if (useGpuForNonbonded)
1136 // Note that any bonded tasks on a GPU always accompany a
1140 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1142 else if (nonbondedTarget == TaskTarget::Gpu)
1144 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1146 else if (bondedTarget == TaskTarget::Gpu)
1148 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because there is none detected.");
1152 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1153 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1159 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1161 else if (pmeTarget == TaskTarget::Gpu)
1163 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1168 GpuTaskAssignment gpuTaskAssignment;
1171 // Produce the task assignment for this rank.
1172 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1173 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1174 useGpuForBonded, pmeRunMode);
1176 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1178 /* Prevent other ranks from continuing after an issue was found
1179 * and reported as a fatal error.
1181 * TODO This function implements a barrier so that MPI runtimes
1182 * can organize an orderly shutdown if one of the ranks has had to
1183 * issue a fatal error in various code already run. When we have
1184 * MPI-aware error handling and reporting, this should be
1189 MPI_Barrier(cr->mpi_comm_mysim);
1195 MPI_Barrier(ms->mpi_comm_masters);
1197 /* We need another barrier to prevent non-master ranks from contiuing
1198 * when an error occured in a different simulation.
1200 MPI_Barrier(cr->mpi_comm_mysim);
1204 /* Now that we know the setup is consistent, check for efficiency */
1205 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1208 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1210 if (thisRankHasDuty(cr, DUTY_PP))
1212 // This works because only one task of each type is currently permitted.
1213 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1214 hasTaskType<GpuTask::Nonbonded>);
1215 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1217 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1218 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1219 init_gpu(nonbondedDeviceInfo);
1221 if (DOMAINDECOMP(cr))
1223 /* When we share GPUs over ranks, we need to know this for the DLB */
1224 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1230 std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1232 gmx_device_info_t *pmeDeviceInfo = nullptr;
1233 // Later, this program could contain kernels that might be later
1234 // re-used as auto-tuning progresses, or subsequent simulations
1236 PmeGpuProgramStorage pmeGpuProgram;
1237 // This works because only one task of each type is currently permitted.
1238 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1239 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1240 if (thisRankHasPmeGpuTask)
1242 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1243 init_gpu(pmeDeviceInfo);
1244 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1245 // TODO It would be nice to move this logic into the factory
1246 // function. See Redmine #2535.
1247 bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1248 if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1250 initializedClfftLibrary = initializeClfftLibrary();
1254 /* getting number of PP/PME threads on this MPI / tMPI rank.
1255 PME: env variable should be read only on one node to make sure it is
1256 identical everywhere;
1258 const int numThreadsOnThisRank =
1259 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1260 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1261 *hwinfo->hardwareTopology,
1262 physicalNodeComm, mdlog);
1264 if (hw_opt.thread_affinity != threadaffOFF)
1266 /* Before setting affinity, check whether the affinity has changed
1267 * - which indicates that probably the OpenMP library has changed it
1268 * since we first checked).
1270 gmx_check_thread_affinity_set(mdlog, cr,
1271 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1273 int numThreadsOnThisNode, intraNodeThreadOffset;
1274 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1275 &intraNodeThreadOffset);
1277 /* Set the CPU affinity */
1278 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1279 numThreadsOnThisRank, numThreadsOnThisNode,
1280 intraNodeThreadOffset, nullptr);
1283 if (mdrunOptions.timingOptions.resetStep > -1)
1285 GMX_LOG(mdlog.info).asParagraph().
1286 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1288 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1292 /* Master synchronizes its value of reset_counters with all nodes
1293 * including PME only nodes */
1294 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1295 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1296 wcycle_set_reset_counters(wcycle, reset_counters);
1299 // Membrane embedding must be initialized before we call init_forcerec()
1304 fprintf(stderr, "Initializing membed");
1306 /* Note that membed cannot work in parallel because mtop is
1307 * changed here. Fix this if we ever want to make it run with
1308 * multiple ranks. */
1309 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1311 .checkpointOptions.period);
1314 std::unique_ptr<MDAtoms> mdAtoms;
1315 std::unique_ptr<gmx_vsite_t> vsite;
1318 if (thisRankHasDuty(cr, DUTY_PP))
1320 /* Initiate forcerecord */
1321 fr = new t_forcerec;
1322 fr->forceProviders = mdModules_->initForceProviders();
1323 init_forcerec(fplog, mdlog, fr, fcd,
1324 inputrec, &mtop, cr, box,
1325 opt2fn("-table", filenames.size(), filenames.data()),
1326 opt2fn("-tablep", filenames.size(), filenames.data()),
1327 opt2fns("-tableb", filenames.size(), filenames.data()),
1328 *hwinfo, nonbondedDeviceInfo,
1333 /* Initialize the mdAtoms structure.
1334 * mdAtoms is not filled with atom data,
1335 * as this can not be done now with domain decomposition.
1337 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1338 if (globalState && thisRankHasPmeGpuTask)
1340 // The pinning of coordinates in the global state object works, because we only use
1341 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1342 // points to the global state object without DD.
1343 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1344 // which should also perform the pinning.
1345 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1348 /* Initialize the virtual site communication */
1349 vsite = initVsite(mtop, cr);
1351 calc_shifts(box, fr->shift_vec);
1353 /* With periodic molecules the charge groups should be whole at start up
1354 * and the virtual sites should not be far from their proper positions.
1356 if (!inputrec->bContinuation && MASTER(cr) &&
1357 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1359 /* Make molecules whole at start of run */
1360 if (fr->ePBC != epbcNONE)
1362 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1366 /* Correct initial vsite positions are required
1367 * for the initial distribution in the domain decomposition
1368 * and for the initial shell prediction.
1370 constructVsitesGlobal(mtop, globalState->x);
1374 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1376 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1377 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1382 /* This is a PME only node */
1384 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1386 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1387 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1390 gmx_pme_t *sepPmeData = nullptr;
1391 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1392 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1393 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1395 /* Initiate PME if necessary,
1396 * either on all nodes or on dedicated PME nodes only. */
1397 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1399 if (mdAtoms && mdAtoms->mdatoms())
1401 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1402 if (EVDW_PME(inputrec->vdwtype))
1404 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1407 if (cr->npmenodes > 0)
1409 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1410 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1411 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1414 if (thisRankHasDuty(cr, DUTY_PME))
1418 pmedata = gmx_pme_init(cr,
1419 getNumPmeDomains(cr->dd),
1421 mtop.natoms, nChargePerturbed != 0, nTypePerturbed != 0,
1422 mdrunOptions.reproducible,
1423 ewaldcoeff_q, ewaldcoeff_lj,
1424 gmx_omp_nthreads_get(emntPME),
1425 pmeRunMode, nullptr,
1426 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1428 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1433 if (EI_DYNAMICS(inputrec->eI))
1435 /* Turn on signal handling on all nodes */
1437 * (A user signal from the PME nodes (if any)
1438 * is communicated to the PP nodes.
1440 signal_handler_install();
1443 pull_t *pull_work = nullptr;
1444 if (thisRankHasDuty(cr, DUTY_PP))
1446 /* Assumes uniform use of the number of OpenMP threads */
1447 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1449 if (inputrec->bPull)
1451 /* Initialize pull code */
1453 init_pull(fplog, inputrec->pull, inputrec,
1454 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1455 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1457 initPullHistory(pull_work, &observablesHistory);
1459 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1461 init_pull_output_files(pull_work,
1462 filenames.size(), filenames.data(), oenv,
1463 continuationOptions);
1467 std::unique_ptr<EnforcedRotation> enforcedRotation;
1470 /* Initialize enforced rotation code */
1471 enforcedRotation = init_rot(fplog,
1483 t_swap *swap = nullptr;
1484 if (inputrec->eSwapCoords != eswapNO)
1486 /* Initialize ion swapping code */
1487 swap = init_swapcoords(fplog, inputrec,
1488 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1489 &mtop, globalState.get(), &observablesHistory,
1490 cr, &atomSets, oenv, mdrunOptions);
1493 /* Let makeConstraints know whether we have essential dynamics constraints.
1494 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1496 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1497 || observablesHistory.edsamHistory);
1498 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics,
1499 fplog, *mdAtoms->mdatoms(),
1500 cr, ms, nrnb, wcycle, fr->bMolPBC);
1502 /* Energy terms and groups */
1503 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(), inputrec->fepvals->n_lambda);
1505 /* Set up interactive MD (IMD) */
1506 auto imdSession = makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1507 MASTER(cr) ? globalState->x.rvec_array() : nullptr,
1508 filenames.size(), filenames.data(), oenv, mdrunOptions);
1510 if (DOMAINDECOMP(cr))
1512 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1513 /* This call is not included in init_domain_decomposition mainly
1514 * because fr->cginfo_mb is set later.
1516 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1517 domdecOptions.checkBondedInteractions,
1521 // TODO This is not the right place to manage the lifetime of
1522 // this data structure, but currently it's the easiest way to
1523 // make it work. Later, it should probably be made/updated
1524 // after the workload for the lifetime of a PP domain is
1526 PpForceWorkload ppForceWorkload;
1528 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to integrator.");
1529 /* Now do whatever the user wants us to do (how flexible...) */
1530 Integrator integrator {
1531 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1534 vsite.get(), constr.get(),
1535 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1537 mdModules_->outputProvider(),
1538 inputrec, imdSession.get(), pull_work, swap, &mtop,
1541 &observablesHistory,
1542 mdAtoms.get(), nrnb, wcycle, fr,
1547 walltime_accounting,
1548 std::move(stopHandlerBuilder_)
1550 integrator.run(inputrec->eI, doRerun);
1552 if (inputrec->bPull)
1554 finish_pull(pull_work);
1556 finish_swapcoords(swap);
1560 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1562 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1563 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1566 wallcycle_stop(wcycle, ewcRUN);
1568 /* Finish up, write some stuff
1569 * if rerunMD, don't write last frame again
1571 finish_run(fplog, mdlog, cr,
1572 inputrec, nrnb, wcycle, walltime_accounting,
1573 fr ? fr->nbv.get() : nullptr,
1575 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1577 // clean up cycle counter
1578 wallcycle_destroy(wcycle);
1583 gmx_pme_destroy(pmedata);
1587 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1588 // before we destroy the GPU context(s) in free_gpu_resources().
1589 // Pinned buffers are associated with contexts in CUDA.
1590 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1591 mdAtoms.reset(nullptr);
1592 globalState.reset(nullptr);
1593 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1595 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1596 free_gpu_resources(fr, physicalNodeComm);
1597 free_gpu(nonbondedDeviceInfo);
1598 free_gpu(pmeDeviceInfo);
1599 done_forcerec(fr, mtop.molblock.size());
1604 free_membed(membed);
1607 gmx_hardware_info_free();
1609 /* Does what it says */
1610 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1611 walltime_accounting_destroy(walltime_accounting);
1614 // Ensure log file content is written
1617 gmx_fio_flush(logFileHandle);
1620 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1621 * exceptions were enabled before function was called. */
1624 gmx_fedisableexcept();
1627 auto rc = static_cast<int>(gmx_get_stop_condition());
1630 /* we need to join all threads. The sub-threads join when they
1631 exit this function, but the master thread needs to be told to
1633 if (PAR(cr) && MASTER(cr))
1637 //TODO free commrec in MPI simulations
1643 Mdrunner::~Mdrunner()
1645 // Clean up of the Manager.
1646 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1647 // but okay as long as threads synchronize some time before adding or accessing
1648 // a new set of restraints.
1649 if (restraintManager_)
1651 restraintManager_->clear();
1652 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1653 "restraints added during runner life time should be cleared at runner destruction.");
1657 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1660 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1661 // Not sure if this should be logged through the md logger or something else,
1662 // but it is helpful to have some sort of INFO level message sent somewhere.
1663 // std::cout << "Registering restraint named " << name << std::endl;
1665 // When multiple restraints are used, it may be wasteful to register them separately.
1666 // Maybe instead register an entire Restraint Manager as a force provider.
1667 restraintManager_->addToSpec(std::move(puller),
1671 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1672 : mdModules_(std::move(mdModules))
1676 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1678 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1679 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1681 class Mdrunner::BuilderImplementation
1684 BuilderImplementation() = delete;
1685 BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1686 SimulationContext * context);
1687 ~BuilderImplementation();
1689 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1690 real forceWarningThreshold);
1692 void addDomdec(const DomdecOptions &options);
1694 void addVerletList(int nstlist);
1696 void addReplicaExchange(const ReplicaExchangeParameters ¶ms);
1698 void addMultiSim(gmx_multisim_t* multisim);
1700 void addNonBonded(const char* nbpu_opt);
1702 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1704 void addBondedTaskAssignment(const char* bonded_opt);
1706 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1708 void addFilenames(ArrayRef <const t_filenm> filenames);
1710 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1712 void addLogFile(t_fileio *logFileHandle);
1714 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1720 // Default parameters copied from runner.h
1721 // \todo Clarify source(s) of default parameters.
1723 const char* nbpu_opt_ = nullptr;
1724 const char* pme_opt_ = nullptr;
1725 const char* pme_fft_opt_ = nullptr;
1726 const char *bonded_opt_ = nullptr;
1728 MdrunOptions mdrunOptions_;
1730 DomdecOptions domdecOptions_;
1732 ReplicaExchangeParameters replicaExchangeParameters_;
1734 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1737 //! Non-owning multisim communicator handle.
1738 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1740 //! Print a warning if any force is larger than this (in kJ/mol nm).
1741 real forceWarningThreshold_ = -1;
1743 //! The modules that comprise the functionality of mdrun.
1744 std::unique_ptr<MDModules> mdModules_;
1746 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1749 * \todo Establish robust protocol to make sure resources remain valid.
1750 * SimulationContext will likely be separated into multiple layers for
1751 * different levels of access and different phases of execution. Ref
1752 * https://redmine.gromacs.org/issues/2375
1753 * https://redmine.gromacs.org/issues/2587
1755 SimulationContext* context_ = nullptr;
1757 //! \brief Parallelism information.
1758 gmx_hw_opt_t hardwareOptions_;
1760 //! filename options for simulation.
1761 ArrayRef<const t_filenm> filenames_;
1763 /*! \brief Handle to output environment.
1765 * \todo gmx_output_env_t needs lifetime management.
1767 gmx_output_env_t* outputEnvironment_ = nullptr;
1769 /*! \brief Non-owning handle to MD log file.
1771 * \todo Context should own output facilities for client.
1772 * \todo Improve log file handle management.
1774 * Code managing the FILE* relies on the ability to set it to
1775 * nullptr to check whether the filehandle is valid.
1777 t_fileio* logFileHandle_ = nullptr;
1780 * \brief Builder for simulation stop signal handler.
1782 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1785 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1786 SimulationContext * context) :
1787 mdModules_(std::move(mdModules)),
1790 GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1793 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1795 Mdrunner::BuilderImplementation &
1796 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1797 real forceWarningThreshold)
1799 mdrunOptions_ = options;
1800 forceWarningThreshold_ = forceWarningThreshold;
1804 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1806 domdecOptions_ = options;
1809 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1814 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1816 replicaExchangeParameters_ = params;
1819 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1821 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1824 Mdrunner Mdrunner::BuilderImplementation::build()
1826 auto newRunner = Mdrunner(std::move(mdModules_));
1828 GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1830 newRunner.mdrunOptions = mdrunOptions_;
1831 newRunner.domdecOptions = domdecOptions_;
1833 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1834 newRunner.hw_opt = hardwareOptions_;
1836 // No invariant to check. This parameter exists to optionally override other behavior.
1837 newRunner.nstlist_cmdline = nstlist_;
1839 newRunner.replExParams = replicaExchangeParameters_;
1841 newRunner.filenames = filenames_;
1843 GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1844 newRunner.cr = context_->communicationRecord_;
1848 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1849 newRunner.ms = *multisim_;
1853 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1856 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1857 // \todo Update sanity checking when output environment has clearly specified invariants.
1858 // Initialization and default values for oenv are not well specified in the current version.
1859 if (outputEnvironment_)
1861 newRunner.oenv = outputEnvironment_;
1865 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1868 newRunner.logFileHandle = logFileHandle_;
1872 newRunner.nbpu_opt = nbpu_opt_;
1876 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1879 if (pme_opt_ && pme_fft_opt_)
1881 newRunner.pme_opt = pme_opt_;
1882 newRunner.pme_fft_opt = pme_fft_opt_;
1886 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1891 newRunner.bonded_opt = bonded_opt_;
1895 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1898 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1900 if (stopHandlerBuilder_)
1902 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1906 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1912 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1914 nbpu_opt_ = nbpu_opt;
1917 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1918 const char* pme_fft_opt)
1921 pme_fft_opt_ = pme_fft_opt;
1924 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1926 bonded_opt_ = bonded_opt;
1929 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1931 hardwareOptions_ = hardwareOptions;
1934 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1936 filenames_ = filenames;
1939 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1941 outputEnvironment_ = outputEnvironment;
1944 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1946 logFileHandle_ = logFileHandle;
1949 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1951 stopHandlerBuilder_ = std::move(builder);
1954 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
1955 compat::not_null<SimulationContext*> context) :
1956 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context)}
1960 MdrunnerBuilder::~MdrunnerBuilder() = default;
1962 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1963 real forceWarningThreshold)
1965 impl_->setExtraMdrunOptions(options, forceWarningThreshold);
1969 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1971 impl_->addDomdec(options);
1975 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1977 impl_->addVerletList(nstlist);
1981 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1983 impl_->addReplicaExchange(params);
1987 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1989 impl_->addMultiSim(multisim);
1993 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1995 impl_->addNonBonded(nbpu_opt);
1999 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
2000 const char* pme_fft_opt)
2002 // The builder method may become more general in the future, but in this version,
2003 // parameters for PME electrostatics are both required and the only parameters
2005 if (pme_opt && pme_fft_opt)
2007 impl_->addPME(pme_opt, pme_fft_opt);
2011 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2016 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2018 impl_->addBondedTaskAssignment(bonded_opt);
2022 Mdrunner MdrunnerBuilder::build()
2024 return impl_->build();
2027 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2029 impl_->addHardwareOptions(hardwareOptions);
2033 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2035 impl_->addFilenames(filenames);
2039 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2041 impl_->addOutputEnvironment(outputEnvironment);
2045 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2047 impl_->addLogFile(logFileHandle);
2051 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2053 impl_->addStopHandlerBuilder(std::move(builder));
2057 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2059 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;