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/enerdata_utils.h"
89 #include "gromacs/mdlib/force.h"
90 #include "gromacs/mdlib/forcerec.h"
91 #include "gromacs/mdlib/gmx_omp_nthreads.h"
92 #include "gromacs/mdlib/makeconstraints.h"
93 #include "gromacs/mdlib/md_support.h"
94 #include "gromacs/mdlib/mdatoms.h"
95 #include "gromacs/mdlib/membed.h"
96 #include "gromacs/mdlib/ppforceworkload.h"
97 #include "gromacs/mdlib/qmmm.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/stophandler.h"
100 #include "gromacs/mdrun/logging.h"
101 #include "gromacs/mdrun/mdmodules.h"
102 #include "gromacs/mdrun/multisim.h"
103 #include "gromacs/mdrun/simulationcontext.h"
104 #include "gromacs/mdrunutility/printtime.h"
105 #include "gromacs/mdrunutility/threadaffinity.h"
106 #include "gromacs/mdtypes/commrec.h"
107 #include "gromacs/mdtypes/enerdata.h"
108 #include "gromacs/mdtypes/fcdata.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/md_enums.h"
111 #include "gromacs/mdtypes/mdrunoptions.h"
112 #include "gromacs/mdtypes/observableshistory.h"
113 #include "gromacs/mdtypes/state.h"
114 #include "gromacs/nbnxm/gpu_data_mgmt.h"
115 #include "gromacs/nbnxm/nbnxm.h"
116 #include "gromacs/nbnxm/pairlist_tuning.h"
117 #include "gromacs/pbcutil/pbc.h"
118 #include "gromacs/pulling/output.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/pulling/pull_rotation.h"
121 #include "gromacs/restraint/manager.h"
122 #include "gromacs/restraint/restraintmdmodule.h"
123 #include "gromacs/restraint/restraintpotential.h"
124 #include "gromacs/swap/swapcoords.h"
125 #include "gromacs/taskassignment/decidegpuusage.h"
126 #include "gromacs/taskassignment/resourcedivision.h"
127 #include "gromacs/taskassignment/taskassignment.h"
128 #include "gromacs/taskassignment/usergpuids.h"
129 #include "gromacs/timing/gpu_timing.h"
130 #include "gromacs/timing/wallcycle.h"
131 #include "gromacs/timing/wallcyclereporting.h"
132 #include "gromacs/topology/mtop_util.h"
133 #include "gromacs/trajectory/trajectoryframe.h"
134 #include "gromacs/utility/basenetwork.h"
135 #include "gromacs/utility/cstringutil.h"
136 #include "gromacs/utility/exceptions.h"
137 #include "gromacs/utility/fatalerror.h"
138 #include "gromacs/utility/filestream.h"
139 #include "gromacs/utility/gmxassert.h"
140 #include "gromacs/utility/gmxmpi.h"
141 #include "gromacs/utility/logger.h"
142 #include "gromacs/utility/loggerbuilder.h"
143 #include "gromacs/utility/physicalnodecommunicator.h"
144 #include "gromacs/utility/pleasecite.h"
145 #include "gromacs/utility/programcontext.h"
146 #include "gromacs/utility/smalloc.h"
147 #include "gromacs/utility/stringutil.h"
149 #include "integrator.h"
150 #include "replicaexchange.h"
153 #include "corewrap.h"
159 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
161 * Used to ensure that the master thread does not modify mdrunner during copy
162 * on the spawned threads. */
163 static void threadMpiMdrunnerAccessBarrier()
166 MPI_Barrier(MPI_COMM_WORLD);
170 Mdrunner Mdrunner::cloneOnSpawnedThread() const
172 auto newRunner = Mdrunner(std::make_unique<MDModules>());
174 // All runners in the same process share a restraint manager resource because it is
175 // part of the interface to the client code, which is associated only with the
176 // original thread. Handles to the same resources can be obtained by copy.
178 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
181 // Copy original cr pointer before master thread can pass the thread barrier
182 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
184 // Copy members of master runner.
185 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
186 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
187 newRunner.hw_opt = hw_opt;
188 newRunner.filenames = filenames;
190 newRunner.oenv = oenv;
191 newRunner.mdrunOptions = mdrunOptions;
192 newRunner.domdecOptions = domdecOptions;
193 newRunner.nbpu_opt = nbpu_opt;
194 newRunner.pme_opt = pme_opt;
195 newRunner.pme_fft_opt = pme_fft_opt;
196 newRunner.bonded_opt = bonded_opt;
197 newRunner.nstlist_cmdline = nstlist_cmdline;
198 newRunner.replExParams = replExParams;
199 newRunner.pforce = pforce;
201 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
203 threadMpiMdrunnerAccessBarrier();
205 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
210 /*! \brief The callback used for running on spawned threads.
212 * Obtains the pointer to the master mdrunner object from the one
213 * argument permitted to the thread-launch API call, copies it to make
214 * a new runner for this thread, reinitializes necessary data, and
215 * proceeds to the simulation. */
216 static void mdrunner_start_fn(const void *arg)
220 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
221 /* copy the arg list to make sure that it's thread-local. This
222 doesn't copy pointed-to items, of course; fnm, cr and fplog
223 are reset in the call below, all others should be const. */
224 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
227 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
231 /*! \brief Start thread-MPI threads.
233 * Called by mdrunner() to start a specific number of threads
234 * (including the main thread) for thread-parallel runs. This in turn
235 * calls mdrunner() for each thread. All options are the same as for
237 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
240 /* first check whether we even need to start tMPI */
241 if (numThreadsToLaunch < 2)
247 /* now spawn new threads that start mdrunner_start_fn(), while
248 the main thread returns, we set thread affinity later */
249 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
250 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
252 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
255 threadMpiMdrunnerAccessBarrier();
257 GMX_UNUSED_VALUE(mdrunner_start_fn);
260 return reinitialize_commrec_for_this_thread(cr);
265 /*! \brief Initialize variables for Verlet scheme simulation */
266 static void prepare_verlet_scheme(FILE *fplog,
270 const gmx_mtop_t *mtop,
272 bool makeGpuPairList,
273 const gmx::CpuInfo &cpuinfo)
275 /* For NVE simulations, we will retain the initial list buffer */
276 if (EI_DYNAMICS(ir->eI) &&
277 ir->verletbuf_tol > 0 &&
278 !(EI_MD(ir->eI) && ir->etc == etcNO))
280 /* Update the Verlet buffer size for the current run setup */
282 /* Here we assume SIMD-enabled kernels are being used. But as currently
283 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
284 * and 4x2 gives a larger buffer than 4x4, this is ok.
286 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
287 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
289 const real rlist_new =
290 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
292 if (rlist_new != ir->rlist)
294 if (fplog != nullptr)
296 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
297 ir->rlist, rlist_new,
298 listSetup.cluster_size_i, listSetup.cluster_size_j);
300 ir->rlist = rlist_new;
304 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
306 gmx_fatal(FARGS, "Can not set nstlist without %s",
307 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
310 if (EI_DYNAMICS(ir->eI))
312 /* Set or try nstlist values */
313 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
317 /*! \brief Override the nslist value in inputrec
319 * with value passed on the command line (if any)
321 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
322 int64_t nsteps_cmdline,
327 /* override with anything else than the default -2 */
328 if (nsteps_cmdline > -2)
330 char sbuf_steps[STEPSTRSIZE];
331 char sbuf_msg[STRLEN];
333 ir->nsteps = nsteps_cmdline;
334 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
336 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
337 gmx_step_str(nsteps_cmdline, sbuf_steps),
338 fabs(nsteps_cmdline*ir->delta_t));
342 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
343 gmx_step_str(nsteps_cmdline, sbuf_steps));
346 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
348 else if (nsteps_cmdline < -2)
350 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
353 /* Do nothing if nsteps_cmdline == -2 */
359 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
361 * If not, and if a warning may be issued, logs a warning about
362 * falling back to CPU code. With thread-MPI, only the first
363 * call to this function should have \c issueWarning true. */
364 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
365 const t_inputrec *ir,
368 if (ir->opts.ngener - ir->nwall > 1)
370 /* The GPU code does not support more than one energy group.
371 * If the user requested GPUs explicitly, a fatal error is given later.
375 GMX_LOG(mdlog.warning).asParagraph()
376 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
377 "For better performance, run on the GPU without energy groups and then do "
378 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
385 //! Initializes the logger for mdrun.
386 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
388 gmx::LoggerBuilder builder;
389 if (fplog != nullptr)
391 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
393 if (cr == nullptr || SIMMASTER(cr))
395 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
396 &gmx::TextOutputFile::standardError());
398 return builder.build();
401 //! Make a TaskTarget from an mdrun argument string.
402 static TaskTarget findTaskTarget(const char *optionString)
404 TaskTarget returnValue = TaskTarget::Auto;
406 if (strncmp(optionString, "auto", 3) == 0)
408 returnValue = TaskTarget::Auto;
410 else if (strncmp(optionString, "cpu", 3) == 0)
412 returnValue = TaskTarget::Cpu;
414 else if (strncmp(optionString, "gpu", 3) == 0)
416 returnValue = TaskTarget::Gpu;
420 GMX_ASSERT(false, "Option string should have been checked for sanity already");
426 //! Finish run, aggregate data to print performance info.
427 static void finish_run(FILE *fplog,
428 const gmx::MDLogger &mdlog,
430 const t_inputrec *inputrec,
431 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
432 gmx_walltime_accounting_t walltime_accounting,
433 nonbonded_verlet_t *nbv,
434 const gmx_pme_t *pme,
437 t_nrnb *nrnb_tot = nullptr;
439 double nbfs = 0, mflop = 0;
441 elapsed_time_over_all_ranks,
442 elapsed_time_over_all_threads,
443 elapsed_time_over_all_threads_over_all_ranks;
444 /* Control whether it is valid to print a report. Only the
445 simulation master may print, but it should not do so if the run
446 terminated e.g. before a scheduled reset step. This is
447 complicated by the fact that PME ranks are unaware of the
448 reason why they were sent a pmerecvqxFINISH. To avoid
449 communication deadlocks, we always do the communication for the
450 report, even if we've decided not to write the report, because
451 how long it takes to finish the run is not important when we've
452 decided not to report on the simulation performance.
454 Further, we only report performance for dynamical integrators,
455 because those are the only ones for which we plan to
456 consider doing any optimizations. */
457 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
459 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
461 GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
469 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
478 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
479 elapsed_time_over_all_threads = walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
483 /* reduce elapsed_time over all MPI ranks in the current simulation */
484 MPI_Allreduce(&elapsed_time,
485 &elapsed_time_over_all_ranks,
486 1, MPI_DOUBLE, MPI_SUM,
488 elapsed_time_over_all_ranks /= cr->nnodes;
489 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
490 * current simulation. */
491 MPI_Allreduce(&elapsed_time_over_all_threads,
492 &elapsed_time_over_all_threads_over_all_ranks,
493 1, MPI_DOUBLE, MPI_SUM,
499 elapsed_time_over_all_ranks = elapsed_time;
500 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
505 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
512 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
514 print_dd_statistics(cr, inputrec, fplog);
517 /* TODO Move the responsibility for any scaling by thread counts
518 * to the code that handled the thread region, so that there's a
519 * mechanism to keep cycle counting working during the transition
520 * to task parallelism. */
521 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
522 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
523 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
524 auto cycle_sum(wallcycle_sum(cr, wcycle));
528 auto nbnxn_gpu_timings = (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
529 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
531 if (pme_gpu_task_enabled(pme))
533 pme_gpu_get_timings(pme, &pme_gpu_timings);
535 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
536 elapsed_time_over_all_ranks,
541 if (EI_DYNAMICS(inputrec->eI))
543 delta_t = inputrec->delta_t;
548 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
549 elapsed_time_over_all_ranks,
550 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
551 delta_t, nbfs, mflop);
555 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
556 elapsed_time_over_all_ranks,
557 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
558 delta_t, nbfs, mflop);
563 int Mdrunner::mdrunner()
567 t_forcerec *fr = nullptr;
568 t_fcdata *fcd = nullptr;
569 real ewaldcoeff_q = 0;
570 real ewaldcoeff_lj = 0;
571 int nChargePerturbed = -1, nTypePerturbed = 0;
572 gmx_wallcycle_t wcycle;
573 gmx_walltime_accounting_t walltime_accounting = nullptr;
574 gmx_membed_t * membed = nullptr;
575 gmx_hw_info_t *hwinfo = nullptr;
577 /* CAUTION: threads may be started later on in this function, so
578 cr doesn't reflect the final parallel state right now */
579 t_inputrec inputrecInstance;
580 t_inputrec *inputrec = &inputrecInstance;
583 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
584 bool doRerun = mdrunOptions.rerun;
586 // Handle task-assignment related user options.
587 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
588 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
589 std::vector<int> gpuIdsAvailable;
592 gpuIdsAvailable = parseUserGpuIdString(hw_opt.gpuIdsAvailable);
594 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
596 std::vector<int> userGpuTaskAssignment;
599 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
601 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
602 auto nonbondedTarget = findTaskTarget(nbpu_opt);
603 auto pmeTarget = findTaskTarget(pme_opt);
604 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
605 auto bondedTarget = findTaskTarget(bonded_opt);
606 PmeRunMode pmeRunMode = PmeRunMode::None;
608 // Here we assume that SIMMASTER(cr) does not change even after the
609 // threads are started.
611 FILE *fplog = nullptr;
612 // If we are appending, we don't write log output because we need
613 // to check that the old log file matches what the checkpoint file
614 // expects. Otherwise, we should start to write log output now if
615 // there is a file ready for it.
616 if (logFileHandle != nullptr && !mdrunOptions.continuationOptions.appendFiles)
618 fplog = gmx_fio_getfp(logFileHandle);
620 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
621 gmx::MDLogger mdlog(logOwner.logger());
623 // TODO The thread-MPI master rank makes a working
624 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
625 // after the threads have been launched. This works because no use
626 // is made of that communicator until after the execution paths
627 // have rejoined. But it is likely that we can improve the way
628 // this is expressed, e.g. by expressly running detection only the
629 // master rank for thread-MPI, rather than relying on the mutex
630 // and reference count.
631 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
632 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
634 gmx_print_detected_hardware(fplog, cr, ms, mdlog, hwinfo);
636 std::vector<int> gpuIdsToUse;
637 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
638 if (gpuIdsAvailable.empty())
640 gpuIdsToUse = compatibleGpus;
644 for (const auto &availableGpuId : gpuIdsAvailable)
646 bool availableGpuIsCompatible = false;
647 for (const auto &compatibleGpuId : compatibleGpus)
649 if (availableGpuId == compatibleGpuId)
651 availableGpuIsCompatible = true;
655 if (!availableGpuIsCompatible)
657 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);
659 gpuIdsToUse.push_back(availableGpuId);
663 if (fplog != nullptr)
665 /* Print references after all software/hardware printing */
666 please_cite(fplog, "Abraham2015");
667 please_cite(fplog, "Pall2015");
668 please_cite(fplog, "Pronk2013");
669 please_cite(fplog, "Hess2008b");
670 please_cite(fplog, "Spoel2005a");
671 please_cite(fplog, "Lindahl2001a");
672 please_cite(fplog, "Berendsen95a");
673 writeSourceDoi(fplog);
676 std::unique_ptr<t_state> globalState;
680 /* Only the master rank has the global state */
681 globalState = std::make_unique<t_state>();
683 /* Read (nearly) all data required for the simulation */
684 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
686 /* In rerun, set velocities to zero if present */
687 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
689 // rerun does not use velocities
690 GMX_LOG(mdlog.info).asParagraph().appendText(
691 "Rerun trajectory contains velocities. Rerun does only evaluate "
692 "potential energy and forces. The velocities will be ignored.");
693 for (int i = 0; i < globalState->natoms; i++)
695 clear_rvec(globalState->v[i]);
697 globalState->flags &= ~(1 << estV);
700 if (inputrec->cutoff_scheme != ecutsVERLET)
702 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.");
706 /* Check and update the hardware options for internal consistency */
707 check_and_update_hw_opt_1(mdlog, &hw_opt, cr, domdecOptions.numPmeRanks);
709 /* Early check for externally set process affinity. */
710 gmx_check_thread_affinity_set(mdlog, cr,
711 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
713 if (GMX_THREAD_MPI && SIMMASTER(cr))
715 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
717 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
720 /* Since the master knows the cut-off scheme, update hw_opt for this.
721 * This is done later for normal MPI and also once more with tMPI
722 * for all tMPI ranks.
724 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
726 bool useGpuForNonbonded = false;
727 bool useGpuForPme = false;
730 // If the user specified the number of ranks, then we must
731 // respect that, but in default mode, we need to allow for
732 // the number of GPUs to choose the number of ranks.
733 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
734 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
735 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
736 canUseGpuForNonbonded,
737 inputrec->cutoff_scheme == ecutsVERLET,
738 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
739 hw_opt.nthreads_tmpi);
740 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
741 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
742 *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
745 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
747 /* Determine how many thread-MPI ranks to start.
749 * TODO Over-writing the user-supplied value here does
750 * prevent any possible subsequent checks from working
752 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
761 // Now start the threads for thread MPI.
762 cr = spawnThreads(hw_opt.nthreads_tmpi);
763 /* The main thread continues here with a new cr. We don't deallocate
764 the old cr because other threads may still be reading it. */
765 // TODO Both master and spawned threads call dup_tfn and
766 // reinitialize_commrec_for_this_thread. Find a way to express
768 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
770 // END OF CAUTION: cr and physicalNodeComm are now reliable
774 /* now broadcast everything to the non-master nodes/threads: */
775 init_parallel(cr, inputrec, &mtop);
778 // Now each rank knows the inputrec that SIMMASTER read and used,
779 // and (if applicable) cr->nnodes has been assigned the number of
780 // thread-MPI ranks that have been chosen. The ranks can now all
781 // run the task-deciding functions and will agree on the result
782 // without needing to communicate.
784 // TODO Should we do the communication in debug mode to support
785 // having an assertion?
787 // Note that these variables describe only their own node.
789 // Note that when bonded interactions run on a GPU they always run
790 // alongside a nonbonded task, so do not influence task assignment
791 // even though they affect the force calculation workload.
792 bool useGpuForNonbonded = false;
793 bool useGpuForPme = false;
794 bool useGpuForBonded = false;
797 // It's possible that there are different numbers of GPUs on
798 // different nodes, which is the user's responsibilty to
799 // handle. If unsuitable, we will notice that during task
801 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
802 bool usingVerletScheme = inputrec->cutoff_scheme == ecutsVERLET;
803 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
804 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
806 canUseGpuForNonbonded,
808 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
810 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
811 *hwinfo, *inputrec, mtop,
812 cr->nnodes, domdecOptions.numPmeRanks,
814 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
816 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
817 bondedTarget, canUseGpuForBonded,
818 EVDW_PME(inputrec->vdwtype),
819 EEL_PME_EWALD(inputrec->coulombtype),
820 domdecOptions.numPmeRanks, gpusWereDetected);
822 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
823 if (pmeRunMode == PmeRunMode::GPU)
825 if (pmeFftTarget == TaskTarget::Cpu)
827 pmeRunMode = PmeRunMode::Mixed;
830 else if (pmeFftTarget == TaskTarget::Gpu)
832 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.");
835 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
838 // TODO: hide restraint implementation details from Mdrunner.
839 // There is nothing unique about restraints at this point as far as the
840 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
841 // factory functions from the SimulationContext on which to call mdModules_->add().
842 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
843 for (auto && restraint : restraintManager_->getRestraints())
845 auto module = RestraintMDModule::create(restraint,
847 mdModules_->add(std::move(module));
850 // TODO: Error handling
851 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
853 if (fplog != nullptr)
855 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
856 fprintf(fplog, "\n");
861 /* now make sure the state is initialized and propagated */
862 set_state_entries(globalState.get(), inputrec);
865 /* NM and TPI parallelize over force/energy calculations, not atoms,
866 * so we need to initialize and broadcast the global state.
868 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
872 globalState = std::make_unique<t_state>();
874 broadcastStateWithoutDynamics(cr, globalState.get());
877 /* A parallel command line option consistency check that we can
878 only do after any threads have started. */
879 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
880 domdecOptions.numCells[YY] > 1 ||
881 domdecOptions.numCells[ZZ] > 1 ||
882 domdecOptions.numPmeRanks > 0))
885 "The -dd or -npme option request a parallel simulation, "
887 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
890 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
892 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
898 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
900 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
903 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
905 if (domdecOptions.numPmeRanks > 0)
907 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
908 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
911 domdecOptions.numPmeRanks = 0;
914 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
916 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
917 * improve performance with many threads per GPU, since our OpenMP
918 * scaling is bad, but it's difficult to automate the setup.
920 domdecOptions.numPmeRanks = 0;
924 if (domdecOptions.numPmeRanks < 0)
926 domdecOptions.numPmeRanks = 0;
927 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
931 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
938 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
942 /* NMR restraints must be initialized before load_checkpoint,
943 * since with time averaging the history is added to t_state.
944 * For proper consistency check we therefore need to extend
946 * So the PME-only nodes (if present) will also initialize
947 * the distance restraints.
951 /* This needs to be called before read_checkpoint to extend the state */
952 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
954 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
956 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
958 ObservablesHistory observablesHistory = {};
960 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
962 if (continuationOptions.startedFromCheckpoint)
964 /* Check if checkpoint file exists before doing continuation.
965 * This way we can use identical input options for the first and subsequent runs...
969 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
971 cr, domdecOptions.numCells,
972 inputrec, globalState.get(),
973 &bReadEkin, &observablesHistory,
974 continuationOptions.appendFiles,
975 continuationOptions.appendFilesOptionSet,
976 mdrunOptions.reproducible);
980 continuationOptions.haveReadEkin = true;
983 if (continuationOptions.appendFiles && logFileHandle)
985 // Now we can start normal logging to the truncated log file.
986 fplog = gmx_fio_getfp(logFileHandle);
987 prepareLogAppending(fplog);
988 logOwner = buildLogger(fplog, cr);
989 mdlog = logOwner.logger();
993 if (mdrunOptions.numStepsCommandline > -2)
995 GMX_LOG(mdlog.info).asParagraph().
996 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
997 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
999 /* override nsteps with value set on the commamdline */
1000 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1004 copy_mat(globalState->box, box);
1009 gmx_bcast(sizeof(box), box, cr);
1012 /* Update rlist and nstlist. */
1013 if (inputrec->cutoff_scheme == ecutsVERLET)
1015 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1016 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
1019 LocalAtomSetManager atomSets;
1021 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1022 inputrec->eI == eiNM))
1024 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
1026 box, positionsFromStatePointer(globalState.get()),
1028 // Note that local state still does not exist yet.
1032 /* PME, if used, is done on all nodes with 1D decomposition */
1034 cr->duty = (DUTY_PP | DUTY_PME);
1036 if (inputrec->ePBC == epbcSCREW)
1039 "pbc=screw is only implemented with domain decomposition");
1045 /* After possible communicator splitting in make_dd_communicators.
1046 * we can set up the intra/inter node communication.
1048 gmx_setup_nodecomm(fplog, cr);
1054 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1055 "This is simulation %d out of %d running as a composite GROMACS\n"
1056 "multi-simulation job. Setup for this simulation:\n",
1059 GMX_LOG(mdlog.warning).appendTextFormatted(
1060 "Using %d MPI %s\n",
1063 cr->nnodes == 1 ? "thread" : "threads"
1065 cr->nnodes == 1 ? "process" : "processes"
1071 /* Check and update hw_opt for the cut-off scheme */
1072 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
1074 /* Check and update the number of OpenMP threads requested */
1075 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1078 gmx_omp_nthreads_init(mdlog, cr,
1079 hwinfo->nthreads_hw_avail,
1080 physicalNodeComm.size_,
1081 hw_opt.nthreads_omp,
1082 hw_opt.nthreads_omp_pme,
1083 !thisRankHasDuty(cr, DUTY_PP),
1084 inputrec->cutoff_scheme == ecutsVERLET);
1086 // Enable FP exception detection for the Verlet scheme, but not in
1087 // Release mode and not for compilers with known buggy FP
1088 // exception support (clang with any optimization) or suspected
1089 // buggy FP exception support (gcc 7.* with optimization).
1090 #if !defined NDEBUG && \
1091 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1092 && defined __OPTIMIZE__)
1093 const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
1095 const bool bEnableFPE = false;
1097 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1100 gmx_feenableexcept();
1103 // Build a data structure that expresses which kinds of non-bonded
1104 // task are handled by this rank.
1106 // TODO Later, this might become a loop over all registered modules
1107 // relevant to the mdp inputs, to find those that have such tasks.
1109 // TODO This could move before init_domain_decomposition() as part
1110 // of refactoring that separates the responsibility for duty
1111 // assignment from setup for communication between tasks, and
1112 // setup for tasks handled with a domain (ie including short-ranged
1113 // tasks, bonded tasks, etc.).
1115 // Note that in general useGpuForNonbonded, etc. can have a value
1116 // that is inconsistent with the presence of actual GPUs on any
1117 // rank, and that is not known to be a problem until the
1118 // duty of the ranks on a node become known.
1120 // TODO Later we might need the concept of computeTasksOnThisRank,
1121 // from which we construct gpuTasksOnThisRank.
1123 // Currently the DD code assigns duty to ranks that can
1124 // include PP work that currently can be executed on a single
1125 // GPU, if present and compatible. This has to be coordinated
1126 // across PP ranks on a node, with possible multiple devices
1127 // or sharing devices on a node, either from the user
1128 // selection, or automatically.
1129 auto haveGpus = !gpuIdsToUse.empty();
1130 std::vector<GpuTask> gpuTasksOnThisRank;
1131 if (thisRankHasDuty(cr, DUTY_PP))
1133 if (useGpuForNonbonded)
1135 // Note that any bonded tasks on a GPU always accompany a
1139 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1141 else if (nonbondedTarget == TaskTarget::Gpu)
1143 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1145 else if (bondedTarget == TaskTarget::Gpu)
1147 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because there is none detected.");
1151 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1152 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1158 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1160 else if (pmeTarget == TaskTarget::Gpu)
1162 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1167 GpuTaskAssignment gpuTaskAssignment;
1170 // Produce the task assignment for this rank.
1171 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1172 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1173 useGpuForBonded, pmeRunMode);
1175 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1177 /* Prevent other ranks from continuing after an issue was found
1178 * and reported as a fatal error.
1180 * TODO This function implements a barrier so that MPI runtimes
1181 * can organize an orderly shutdown if one of the ranks has had to
1182 * issue a fatal error in various code already run. When we have
1183 * MPI-aware error handling and reporting, this should be
1188 MPI_Barrier(cr->mpi_comm_mysim);
1194 MPI_Barrier(ms->mpi_comm_masters);
1196 /* We need another barrier to prevent non-master ranks from contiuing
1197 * when an error occured in a different simulation.
1199 MPI_Barrier(cr->mpi_comm_mysim);
1203 /* Now that we know the setup is consistent, check for efficiency */
1204 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1207 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1209 if (thisRankHasDuty(cr, DUTY_PP))
1211 // This works because only one task of each type is currently permitted.
1212 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1213 hasTaskType<GpuTask::Nonbonded>);
1214 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1216 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1217 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1218 init_gpu(nonbondedDeviceInfo);
1220 if (DOMAINDECOMP(cr))
1222 /* When we share GPUs over ranks, we need to know this for the DLB */
1223 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1229 std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1231 gmx_device_info_t *pmeDeviceInfo = nullptr;
1232 // Later, this program could contain kernels that might be later
1233 // re-used as auto-tuning progresses, or subsequent simulations
1235 PmeGpuProgramStorage pmeGpuProgram;
1236 // This works because only one task of each type is currently permitted.
1237 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1238 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1239 if (thisRankHasPmeGpuTask)
1241 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1242 init_gpu(pmeDeviceInfo);
1243 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1244 // TODO It would be nice to move this logic into the factory
1245 // function. See Redmine #2535.
1246 bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1247 if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1249 initializedClfftLibrary = initializeClfftLibrary();
1253 /* getting number of PP/PME threads on this MPI / tMPI rank.
1254 PME: env variable should be read only on one node to make sure it is
1255 identical everywhere;
1257 const int numThreadsOnThisRank =
1258 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1259 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1260 *hwinfo->hardwareTopology,
1261 physicalNodeComm, mdlog);
1263 if (hw_opt.thread_affinity != threadaffOFF)
1265 /* Before setting affinity, check whether the affinity has changed
1266 * - which indicates that probably the OpenMP library has changed it
1267 * since we first checked).
1269 gmx_check_thread_affinity_set(mdlog, cr,
1270 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1272 int numThreadsOnThisNode, intraNodeThreadOffset;
1273 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1274 &intraNodeThreadOffset);
1276 /* Set the CPU affinity */
1277 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1278 numThreadsOnThisRank, numThreadsOnThisNode,
1279 intraNodeThreadOffset, nullptr);
1282 if (mdrunOptions.timingOptions.resetStep > -1)
1284 GMX_LOG(mdlog.info).asParagraph().
1285 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1287 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1291 /* Master synchronizes its value of reset_counters with all nodes
1292 * including PME only nodes */
1293 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1294 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1295 wcycle_set_reset_counters(wcycle, reset_counters);
1298 // Membrane embedding must be initialized before we call init_forcerec()
1303 fprintf(stderr, "Initializing membed");
1305 /* Note that membed cannot work in parallel because mtop is
1306 * changed here. Fix this if we ever want to make it run with
1307 * multiple ranks. */
1308 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1310 .checkpointOptions.period);
1313 std::unique_ptr<MDAtoms> mdAtoms;
1314 std::unique_ptr<gmx_vsite_t> vsite;
1317 if (thisRankHasDuty(cr, DUTY_PP))
1319 /* Initiate forcerecord */
1320 fr = new t_forcerec;
1321 fr->forceProviders = mdModules_->initForceProviders();
1322 init_forcerec(fplog, mdlog, fr, fcd,
1323 inputrec, &mtop, cr, box,
1324 opt2fn("-table", filenames.size(), filenames.data()),
1325 opt2fn("-tablep", filenames.size(), filenames.data()),
1326 opt2fns("-tableb", filenames.size(), filenames.data()),
1327 *hwinfo, nonbondedDeviceInfo,
1332 /* Initialize the mdAtoms structure.
1333 * mdAtoms is not filled with atom data,
1334 * as this can not be done now with domain decomposition.
1336 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1337 if (globalState && thisRankHasPmeGpuTask)
1339 // The pinning of coordinates in the global state object works, because we only use
1340 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1341 // points to the global state object without DD.
1342 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1343 // which should also perform the pinning.
1344 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1347 /* Initialize the virtual site communication */
1348 vsite = initVsite(mtop, cr);
1350 calc_shifts(box, fr->shift_vec);
1352 /* With periodic molecules the charge groups should be whole at start up
1353 * and the virtual sites should not be far from their proper positions.
1355 if (!inputrec->bContinuation && MASTER(cr) &&
1356 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1358 /* Make molecules whole at start of run */
1359 if (fr->ePBC != epbcNONE)
1361 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1365 /* Correct initial vsite positions are required
1366 * for the initial distribution in the domain decomposition
1367 * and for the initial shell prediction.
1369 constructVsitesGlobal(mtop, globalState->x);
1373 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1375 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1376 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1381 /* This is a PME only node */
1383 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1385 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1386 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1389 gmx_pme_t *sepPmeData = nullptr;
1390 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1391 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1392 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1394 /* Initiate PME if necessary,
1395 * either on all nodes or on dedicated PME nodes only. */
1396 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1398 if (mdAtoms && mdAtoms->mdatoms())
1400 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1401 if (EVDW_PME(inputrec->vdwtype))
1403 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1406 if (cr->npmenodes > 0)
1408 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1409 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1410 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1413 if (thisRankHasDuty(cr, DUTY_PME))
1417 pmedata = gmx_pme_init(cr,
1418 getNumPmeDomains(cr->dd),
1420 mtop.natoms, nChargePerturbed != 0, nTypePerturbed != 0,
1421 mdrunOptions.reproducible,
1422 ewaldcoeff_q, ewaldcoeff_lj,
1423 gmx_omp_nthreads_get(emntPME),
1424 pmeRunMode, nullptr,
1425 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1427 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1432 if (EI_DYNAMICS(inputrec->eI))
1434 /* Turn on signal handling on all nodes */
1436 * (A user signal from the PME nodes (if any)
1437 * is communicated to the PP nodes.
1439 signal_handler_install();
1442 if (thisRankHasDuty(cr, DUTY_PP))
1444 /* Assumes uniform use of the number of OpenMP threads */
1445 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1447 if (inputrec->bPull)
1449 /* Initialize pull code */
1450 inputrec->pull_work =
1451 init_pull(fplog, inputrec->pull, inputrec,
1452 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1453 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1455 initPullHistory(inputrec->pull_work, &observablesHistory);
1457 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1459 init_pull_output_files(inputrec->pull_work,
1460 filenames.size(), filenames.data(), oenv,
1461 continuationOptions);
1465 std::unique_ptr<EnforcedRotation> enforcedRotation;
1468 /* Initialize enforced rotation code */
1469 enforcedRotation = init_rot(fplog,
1481 if (inputrec->eSwapCoords != eswapNO)
1483 /* Initialize ion swapping code */
1484 init_swapcoords(fplog, inputrec, opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1485 &mtop, globalState.get(), &observablesHistory,
1486 cr, &atomSets, oenv, mdrunOptions);
1489 /* Let makeConstraints know whether we have essential dynamics constraints.
1490 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1492 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1493 || observablesHistory.edsamHistory);
1494 auto constr = makeConstraints(mtop, *inputrec, doEssentialDynamics,
1495 fplog, *mdAtoms->mdatoms(),
1496 cr, ms, nrnb, wcycle, fr->bMolPBC);
1498 /* Energy terms and groups */
1499 gmx_enerdata_t *enerd;
1501 init_enerdata(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].nr, inputrec->fepvals->n_lambda, enerd);
1503 /* Set up interactive MD (IMD) */
1504 auto imdSession = makeImdSession(inputrec, cr, wcycle, enerd, ms, &mtop, mdlog,
1505 MASTER(cr) ? globalState->x.rvec_array() : nullptr,
1506 filenames.size(), filenames.data(), oenv, mdrunOptions);
1508 if (DOMAINDECOMP(cr))
1510 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1511 /* This call is not included in init_domain_decomposition mainly
1512 * because fr->cginfo_mb is set later.
1514 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1515 domdecOptions.checkBondedInteractions,
1519 // TODO This is not the right place to manage the lifetime of
1520 // this data structure, but currently it's the easiest way to
1521 // make it work. Later, it should probably be made/updated
1522 // after the workload for the lifetime of a PP domain is
1524 PpForceWorkload ppForceWorkload;
1526 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to integrator.");
1527 /* Now do whatever the user wants us to do (how flexible...) */
1528 Integrator integrator {
1529 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1532 vsite.get(), constr.get(),
1533 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1535 mdModules_->outputProvider(),
1536 inputrec, imdSession.get(), &mtop,
1539 &observablesHistory,
1540 mdAtoms.get(), nrnb, wcycle, fr,
1545 walltime_accounting,
1546 std::move(stopHandlerBuilder_)
1548 integrator.run(inputrec->eI, doRerun);
1550 if (inputrec->bPull)
1552 finish_pull(inputrec->pull_work);
1554 destroy_enerdata(enerd);
1559 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1561 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1562 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1565 wallcycle_stop(wcycle, ewcRUN);
1567 /* Finish up, write some stuff
1568 * if rerunMD, don't write last frame again
1570 finish_run(fplog, mdlog, cr,
1571 inputrec, nrnb, wcycle, walltime_accounting,
1572 fr ? fr->nbv.get() : nullptr,
1574 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1576 // clean up cycle counter
1577 wallcycle_destroy(wcycle);
1582 gmx_pme_destroy(pmedata);
1586 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1587 // before we destroy the GPU context(s) in free_gpu_resources().
1588 // Pinned buffers are associated with contexts in CUDA.
1589 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1590 mdAtoms.reset(nullptr);
1591 globalState.reset(nullptr);
1592 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1594 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1595 free_gpu_resources(fr, physicalNodeComm);
1596 free_gpu(nonbondedDeviceInfo);
1597 free_gpu(pmeDeviceInfo);
1598 done_forcerec(fr, mtop.molblock.size());
1603 free_membed(membed);
1606 gmx_hardware_info_free();
1608 /* Does what it says */
1609 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1610 walltime_accounting_destroy(walltime_accounting);
1613 // Ensure log file content is written
1616 gmx_fio_flush(logFileHandle);
1619 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1620 * exceptions were enabled before function was called. */
1623 gmx_fedisableexcept();
1626 auto rc = static_cast<int>(gmx_get_stop_condition());
1629 /* we need to join all threads. The sub-threads join when they
1630 exit this function, but the master thread needs to be told to
1632 if (PAR(cr) && MASTER(cr))
1636 //TODO free commrec in MPI simulations
1642 Mdrunner::~Mdrunner()
1644 // Clean up of the Manager.
1645 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1646 // but okay as long as threads synchronize some time before adding or accessing
1647 // a new set of restraints.
1648 if (restraintManager_)
1650 restraintManager_->clear();
1651 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1652 "restraints added during runner life time should be cleared at runner destruction.");
1656 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1659 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1660 // Not sure if this should be logged through the md logger or something else,
1661 // but it is helpful to have some sort of INFO level message sent somewhere.
1662 // std::cout << "Registering restraint named " << name << std::endl;
1664 // When multiple restraints are used, it may be wasteful to register them separately.
1665 // Maybe instead register an entire Restraint Manager as a force provider.
1666 restraintManager_->addToSpec(std::move(puller),
1670 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1671 : mdModules_(std::move(mdModules))
1675 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1677 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1678 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1680 class Mdrunner::BuilderImplementation
1683 BuilderImplementation() = delete;
1684 BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1685 SimulationContext * context);
1686 ~BuilderImplementation();
1688 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1689 real forceWarningThreshold);
1691 void addDomdec(const DomdecOptions &options);
1693 void addVerletList(int nstlist);
1695 void addReplicaExchange(const ReplicaExchangeParameters ¶ms);
1697 void addMultiSim(gmx_multisim_t* multisim);
1699 void addNonBonded(const char* nbpu_opt);
1701 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1703 void addBondedTaskAssignment(const char* bonded_opt);
1705 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1707 void addFilenames(ArrayRef <const t_filenm> filenames);
1709 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1711 void addLogFile(t_fileio *logFileHandle);
1713 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1719 // Default parameters copied from runner.h
1720 // \todo Clarify source(s) of default parameters.
1722 const char* nbpu_opt_ = nullptr;
1723 const char* pme_opt_ = nullptr;
1724 const char* pme_fft_opt_ = nullptr;
1725 const char *bonded_opt_ = nullptr;
1727 MdrunOptions mdrunOptions_;
1729 DomdecOptions domdecOptions_;
1731 ReplicaExchangeParameters replicaExchangeParameters_;
1733 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1736 //! Non-owning multisim communicator handle.
1737 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1739 //! Print a warning if any force is larger than this (in kJ/mol nm).
1740 real forceWarningThreshold_ = -1;
1742 //! The modules that comprise the functionality of mdrun.
1743 std::unique_ptr<MDModules> mdModules_;
1745 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1748 * \todo Establish robust protocol to make sure resources remain valid.
1749 * SimulationContext will likely be separated into multiple layers for
1750 * different levels of access and different phases of execution. Ref
1751 * https://redmine.gromacs.org/issues/2375
1752 * https://redmine.gromacs.org/issues/2587
1754 SimulationContext* context_ = nullptr;
1756 //! \brief Parallelism information.
1757 gmx_hw_opt_t hardwareOptions_;
1759 //! filename options for simulation.
1760 ArrayRef<const t_filenm> filenames_;
1762 /*! \brief Handle to output environment.
1764 * \todo gmx_output_env_t needs lifetime management.
1766 gmx_output_env_t* outputEnvironment_ = nullptr;
1768 /*! \brief Non-owning handle to MD log file.
1770 * \todo Context should own output facilities for client.
1771 * \todo Improve log file handle management.
1773 * Code managing the FILE* relies on the ability to set it to
1774 * nullptr to check whether the filehandle is valid.
1776 t_fileio* logFileHandle_ = nullptr;
1779 * \brief Builder for simulation stop signal handler.
1781 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1784 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1785 SimulationContext * context) :
1786 mdModules_(std::move(mdModules)),
1789 GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1792 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1794 Mdrunner::BuilderImplementation &
1795 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1796 real forceWarningThreshold)
1798 mdrunOptions_ = options;
1799 forceWarningThreshold_ = forceWarningThreshold;
1803 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1805 domdecOptions_ = options;
1808 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1813 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1815 replicaExchangeParameters_ = params;
1818 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1820 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1823 Mdrunner Mdrunner::BuilderImplementation::build()
1825 auto newRunner = Mdrunner(std::move(mdModules_));
1827 GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1829 newRunner.mdrunOptions = mdrunOptions_;
1830 newRunner.domdecOptions = domdecOptions_;
1832 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1833 newRunner.hw_opt = hardwareOptions_;
1835 // No invariant to check. This parameter exists to optionally override other behavior.
1836 newRunner.nstlist_cmdline = nstlist_;
1838 newRunner.replExParams = replicaExchangeParameters_;
1840 newRunner.filenames = filenames_;
1842 GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1843 newRunner.cr = context_->communicationRecord_;
1847 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1848 newRunner.ms = *multisim_;
1852 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1855 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1856 // \todo Update sanity checking when output environment has clearly specified invariants.
1857 // Initialization and default values for oenv are not well specified in the current version.
1858 if (outputEnvironment_)
1860 newRunner.oenv = outputEnvironment_;
1864 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1867 newRunner.logFileHandle = logFileHandle_;
1871 newRunner.nbpu_opt = nbpu_opt_;
1875 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1878 if (pme_opt_ && pme_fft_opt_)
1880 newRunner.pme_opt = pme_opt_;
1881 newRunner.pme_fft_opt = pme_fft_opt_;
1885 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1890 newRunner.bonded_opt = bonded_opt_;
1894 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1897 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1899 if (stopHandlerBuilder_)
1901 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1905 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1911 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1913 nbpu_opt_ = nbpu_opt;
1916 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1917 const char* pme_fft_opt)
1920 pme_fft_opt_ = pme_fft_opt;
1923 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1925 bonded_opt_ = bonded_opt;
1928 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1930 hardwareOptions_ = hardwareOptions;
1933 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1935 filenames_ = filenames;
1938 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1940 outputEnvironment_ = outputEnvironment;
1943 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1945 logFileHandle_ = logFileHandle;
1948 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1950 stopHandlerBuilder_ = std::move(builder);
1953 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
1954 compat::not_null<SimulationContext*> context) :
1955 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context)}
1959 MdrunnerBuilder::~MdrunnerBuilder() = default;
1961 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1962 real forceWarningThreshold)
1964 impl_->setExtraMdrunOptions(options, forceWarningThreshold);
1968 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1970 impl_->addDomdec(options);
1974 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1976 impl_->addVerletList(nstlist);
1980 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1982 impl_->addReplicaExchange(params);
1986 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1988 impl_->addMultiSim(multisim);
1992 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1994 impl_->addNonBonded(nbpu_opt);
1998 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
1999 const char* pme_fft_opt)
2001 // The builder method may become more general in the future, but in this version,
2002 // parameters for PME electrostatics are both required and the only parameters
2004 if (pme_opt && pme_fft_opt)
2006 impl_->addPME(pme_opt, pme_fft_opt);
2010 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2015 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2017 impl_->addBondedTaskAssignment(bonded_opt);
2021 Mdrunner MdrunnerBuilder::build()
2023 return impl_->build();
2026 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2028 impl_->addHardwareOptions(hardwareOptions);
2032 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2034 impl_->addFilenames(filenames);
2038 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2040 impl_->addOutputEnvironment(outputEnvironment);
2044 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2046 impl_->addLogFile(logFileHandle);
2050 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2052 impl_->addStopHandlerBuilder(std::move(builder));
2056 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2058 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;