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, 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_mdlib
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/ewald/ewald-utils.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/checkpoint.h"
63 #include "gromacs/fileio/oenv.h"
64 #include "gromacs/fileio/tpxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gpu_utils/gpu_utils.h"
67 #include "gromacs/hardware/cpuinfo.h"
68 #include "gromacs/hardware/detecthardware.h"
69 #include "gromacs/hardware/printhardware.h"
70 #include "gromacs/listed-forces/disre.h"
71 #include "gromacs/listed-forces/orires.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/calc_verletbuf.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/gmx_omp_nthreads.h"
80 #include "gromacs/mdlib/integrator.h"
81 #include "gromacs/mdlib/main.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/mdatoms.h"
84 #include "gromacs/mdlib/mdrun.h"
85 #include "gromacs/mdlib/minimize.h"
86 #include "gromacs/mdlib/nb_verlet.h"
87 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
88 #include "gromacs/mdlib/nbnxn_search.h"
89 #include "gromacs/mdlib/nbnxn_tuning.h"
90 #include "gromacs/mdlib/qmmm.h"
91 #include "gromacs/mdlib/sighandler.h"
92 #include "gromacs/mdlib/sim_util.h"
93 #include "gromacs/mdlib/tpi.h"
94 #include "gromacs/mdrunutility/mdmodules.h"
95 #include "gromacs/mdrunutility/threadaffinity.h"
96 #include "gromacs/mdtypes/commrec.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/md_enums.h"
99 #include "gromacs/mdtypes/observableshistory.h"
100 #include "gromacs/mdtypes/state.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/pulling/pull.h"
103 #include "gromacs/pulling/pull_rotation.h"
104 #include "gromacs/taskassignment/decidegpuusage.h"
105 #include "gromacs/taskassignment/resourcedivision.h"
106 #include "gromacs/taskassignment/taskassignment.h"
107 #include "gromacs/taskassignment/usergpuids.h"
108 #include "gromacs/timing/wallcycle.h"
109 #include "gromacs/topology/mtop_util.h"
110 #include "gromacs/trajectory/trajectoryframe.h"
111 #include "gromacs/utility/cstringutil.h"
112 #include "gromacs/utility/exceptions.h"
113 #include "gromacs/utility/fatalerror.h"
114 #include "gromacs/utility/filestream.h"
115 #include "gromacs/utility/gmxassert.h"
116 #include "gromacs/utility/gmxmpi.h"
117 #include "gromacs/utility/logger.h"
118 #include "gromacs/utility/loggerbuilder.h"
119 #include "gromacs/utility/pleasecite.h"
120 #include "gromacs/utility/programcontext.h"
121 #include "gromacs/utility/smalloc.h"
122 #include "gromacs/utility/stringutil.h"
130 #include "corewrap.h"
133 //! First step used in pressure scaling
134 gmx_int64_t deform_init_init_step_tpx;
135 //! Initial box for pressure scaling
136 matrix deform_init_box_tpx;
137 //! MPI variable for use in pressure scaling
138 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
143 void Mdrunner::reinitializeOnSpawnedThread()
145 // TODO This duplication is formally necessary if any thread might
146 // modify any memory in fnm or the pointers it contains. If the
147 // contents are ever provably const, then we can remove this
148 // allocation (and memory leak).
149 // TODO This should probably become part of a copy constructor for
151 fnm = dup_tfn(nfile, fnm);
153 cr = reinitialize_commrec_for_this_thread(cr);
157 // Only the master rank writes to the log files
162 /*! \brief The callback used for running on spawned threads.
164 * Obtains the pointer to the master mdrunner object from the one
165 * argument permitted to the thread-launch API call, copies it to make
166 * a new runner for this thread, reinitializes necessary data, and
167 * proceeds to the simulation. */
168 static void mdrunner_start_fn(void *arg)
172 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
173 /* copy the arg list to make sure that it's thread-local. This
174 doesn't copy pointed-to items, of course, but those are all
176 gmx::Mdrunner mdrunner = *masterMdrunner;
177 mdrunner.reinitializeOnSpawnedThread();
180 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
184 /*! \brief Start thread-MPI threads.
186 * Called by mdrunner() to start a specific number of threads
187 * (including the main thread) for thread-parallel runs. This in turn
188 * calls mdrunner() for each thread. All options are the same as for
190 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch)
193 /* first check whether we even need to start tMPI */
194 if (numThreadsToLaunch < 2)
199 gmx::Mdrunner spawnedMdrunner = *this;
200 // TODO This duplication is formally necessary if any thread might
201 // modify any memory in fnm or the pointers it contains. If the
202 // contents are ever provably const, then we can remove this
203 // allocation (and memory leak).
204 // TODO This should probably become part of a copy constructor for
206 spawnedMdrunner.fnm = dup_tfn(this->nfile, fnm);
209 /* now spawn new threads that start mdrunner_start_fn(), while
210 the main thread returns, we set thread affinity later */
211 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
212 mdrunner_start_fn, static_cast<void*>(&spawnedMdrunner)) != TMPI_SUCCESS)
214 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
217 GMX_UNUSED_VALUE(mdrunner_start_fn);
220 return reinitialize_commrec_for_this_thread(cr);
225 /*! \brief Initialize variables for Verlet scheme simulation */
226 static void prepare_verlet_scheme(FILE *fplog,
230 const gmx_mtop_t *mtop,
232 bool makeGpuPairList,
233 const gmx::CpuInfo &cpuinfo)
235 /* For NVE simulations, we will retain the initial list buffer */
236 if (EI_DYNAMICS(ir->eI) &&
237 ir->verletbuf_tol > 0 &&
238 !(EI_MD(ir->eI) && ir->etc == etcNO))
240 /* Update the Verlet buffer size for the current run setup */
242 /* Here we assume SIMD-enabled kernels are being used. But as currently
243 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
244 * and 4x2 gives a larger buffer than 4x4, this is ok.
246 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
247 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
250 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &listSetup, nullptr, &rlist_new);
252 if (rlist_new != ir->rlist)
254 if (fplog != nullptr)
256 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
257 ir->rlist, rlist_new,
258 listSetup.cluster_size_i, listSetup.cluster_size_j);
260 ir->rlist = rlist_new;
264 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
266 gmx_fatal(FARGS, "Can not set nstlist without %s",
267 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
270 if (EI_DYNAMICS(ir->eI))
272 /* Set or try nstlist values */
273 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
277 /*! \brief Override the nslist value in inputrec
279 * with value passed on the command line (if any)
281 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
282 gmx_int64_t nsteps_cmdline,
287 /* override with anything else than the default -2 */
288 if (nsteps_cmdline > -2)
290 char sbuf_steps[STEPSTRSIZE];
291 char sbuf_msg[STRLEN];
293 ir->nsteps = nsteps_cmdline;
294 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
296 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
297 gmx_step_str(nsteps_cmdline, sbuf_steps),
298 fabs(nsteps_cmdline*ir->delta_t));
302 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
303 gmx_step_str(nsteps_cmdline, sbuf_steps));
306 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
308 else if (nsteps_cmdline < -2)
310 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
313 /* Do nothing if nsteps_cmdline == -2 */
319 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
321 * If not, and if a warning may be issued, logs a warning about
322 * falling back to CPU code. With thread-MPI, only the first
323 * call to this function should have \c issueWarning true. */
324 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
325 const t_inputrec *ir,
328 if (ir->opts.ngener - ir->nwall > 1)
330 /* The GPU code does not support more than one energy group.
331 * If the user requested GPUs explicitly, a fatal error is given later.
335 GMX_LOG(mdlog.warning).asParagraph()
336 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
337 "For better performance, run on the GPU without energy groups and then do "
338 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
345 //! \brief Return the correct integrator function.
346 static integrator_t *my_integrator(unsigned int ei)
355 if (!EI_DYNAMICS(ei))
357 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
372 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
376 GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
378 GMX_THROW(APIError("Non existing integrator selected"));
382 //! Initializes the logger for mdrun.
383 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
385 gmx::LoggerBuilder builder;
386 if (fplog != nullptr)
388 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
390 if (cr == nullptr || SIMMASTER(cr))
392 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
393 &gmx::TextOutputFile::standardError());
395 return builder.build();
398 //! Make a TaskTarget from an mdrun argument string.
399 static TaskTarget findTaskTarget(const char *optionString)
401 TaskTarget returnValue = TaskTarget::Auto;
403 if (strncmp(optionString, "auto", 3) == 0)
405 returnValue = TaskTarget::Auto;
407 else if (strncmp(optionString, "cpu", 3) == 0)
409 returnValue = TaskTarget::Cpu;
411 else if (strncmp(optionString, "gpu", 3) == 0)
413 returnValue = TaskTarget::Gpu;
417 GMX_ASSERT(false, "Option string should have been checked for sanity already");
423 int Mdrunner::mdrunner()
426 gmx_ddbox_t ddbox = {0};
427 int npme_major, npme_minor;
429 gmx_mtop_t *mtop = nullptr;
430 t_forcerec *fr = nullptr;
431 t_fcdata *fcd = nullptr;
432 real ewaldcoeff_q = 0;
433 real ewaldcoeff_lj = 0;
434 gmx_vsite_t *vsite = nullptr;
436 int nChargePerturbed = -1, nTypePerturbed = 0;
437 gmx_wallcycle_t wcycle;
438 gmx_walltime_accounting_t walltime_accounting = nullptr;
440 gmx_int64_t reset_counters;
441 int nthreads_pme = 1;
442 gmx_membed_t * membed = nullptr;
443 gmx_hw_info_t *hwinfo = nullptr;
445 /* CAUTION: threads may be started later on in this function, so
446 cr doesn't reflect the final parallel state right now */
447 gmx::MDModules mdModules;
448 t_inputrec inputrecInstance;
449 t_inputrec *inputrec = &inputrecInstance;
452 if (mdrunOptions.continuationOptions.appendFiles)
457 bool doMembed = opt2bSet("-membed", nfile, fnm);
458 bool doRerun = mdrunOptions.rerun;
460 // Handle task-assignment related user options.
461 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
462 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
463 std::vector<int> gpuIdsAvailable;
466 gpuIdsAvailable = parseUserGpuIds(hw_opt.gpuIdsAvailable);
467 // TODO We could put the GPU IDs into a std::map to find
468 // duplicates, but for the small numbers of IDs involved, this
469 // code is simple and fast.
470 for (size_t i = 0; i != gpuIdsAvailable.size(); ++i)
472 for (size_t j = i+1; j != gpuIdsAvailable.size(); ++j)
474 if (gpuIdsAvailable[i] == gpuIdsAvailable[j])
476 GMX_THROW(InvalidInputError(formatString("The string of available GPU device IDs '%s' may not contain duplicate device IDs", hw_opt.gpuIdsAvailable.c_str())));
481 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
483 std::vector<int> userGpuTaskAssignment;
486 userGpuTaskAssignment = parseUserGpuIds(hw_opt.userGpuTaskAssignment);
488 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
489 auto nonbondedTarget = findTaskTarget(nbpu_opt);
490 auto pmeTarget = findTaskTarget(pme_opt);
491 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
492 PmeRunMode pmeRunMode = PmeRunMode::None;
494 // Here we assume that SIMMASTER(cr) does not change even after the
495 // threads are started.
496 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
497 gmx::MDLogger mdlog(logOwner.logger());
499 hwinfo = gmx_detect_hardware(mdlog, cr);
501 gmx_print_detected_hardware(fplog, cr, mdlog, hwinfo);
503 std::vector<int> gpuIdsToUse;
504 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
505 if (gpuIdsAvailable.empty())
507 gpuIdsToUse = compatibleGpus;
511 for (const auto &availableGpuId : gpuIdsAvailable)
513 bool availableGpuIsCompatible = false;
514 for (const auto &compatibleGpuId : compatibleGpus)
516 if (availableGpuId == compatibleGpuId)
518 availableGpuIsCompatible = true;
522 if (!availableGpuIsCompatible)
524 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);
526 gpuIdsToUse.push_back(availableGpuId);
530 if (fplog != nullptr)
532 /* Print references after all software/hardware printing */
533 please_cite(fplog, "Abraham2015");
534 please_cite(fplog, "Pall2015");
535 please_cite(fplog, "Pronk2013");
536 please_cite(fplog, "Hess2008b");
537 please_cite(fplog, "Spoel2005a");
538 please_cite(fplog, "Lindahl2001a");
539 please_cite(fplog, "Berendsen95a");
542 std::unique_ptr<t_state> globalState;
546 /* Only the master rank has the global state */
547 globalState = std::unique_ptr<t_state>(new t_state);
549 /* Read (nearly) all data required for the simulation */
550 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), mtop);
552 if (inputrec->cutoff_scheme != ecutsVERLET)
554 if (nstlist_cmdline > 0)
556 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
559 if (!compatibleGpus.empty())
561 GMX_LOG(mdlog.warning).asParagraph().appendText(
562 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
563 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
568 /* Check and update the hardware options for internal consistency */
569 check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
571 /* Early check for externally set process affinity. */
572 gmx_check_thread_affinity_set(mdlog, cr,
573 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
575 if (GMX_THREAD_MPI && SIMMASTER(cr))
577 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
579 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
582 /* Since the master knows the cut-off scheme, update hw_opt for this.
583 * This is done later for normal MPI and also once more with tMPI
584 * for all tMPI ranks.
586 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
588 bool useGpuForNonbonded = false;
589 bool useGpuForPme = false;
592 // If the user specified the number of ranks, then we must
593 // respect that, but in default mode, we need to allow for
594 // the number of GPUs to choose the number of ranks.
596 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
597 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
598 inputrec->cutoff_scheme == ecutsVERLET,
599 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
600 hw_opt.nthreads_tmpi);
601 auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
602 auto canUseGpuForPme = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
603 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
604 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
605 canUseGpuForPme, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
608 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
610 /* Determine how many thread-MPI ranks to start.
612 * TODO Over-writing the user-supplied value here does
613 * prevent any possible subsequent checks from working
615 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
624 // Now start the threads for thread MPI.
625 cr = spawnThreads(hw_opt.nthreads_tmpi);
626 /* The main thread continues here with a new cr. We don't deallocate
627 the old cr because other threads may still be reading it. */
628 // TODO Both master and spawned threads call dup_tfn and
629 // reinitialize_commrec_for_this_thread. Find a way to express
632 /* END OF CAUTION: cr is now reliable */
636 /* now broadcast everything to the non-master nodes/threads: */
637 init_parallel(cr, inputrec, mtop);
640 // Now each rank knows the inputrec that SIMMASTER read and used,
641 // and (if applicable) cr->nnodes has been assigned the number of
642 // thread-MPI ranks that have been chosen. The ranks can now all
643 // run the task-deciding functions and will agree on the result
644 // without needing to communicate.
646 // TODO Should we do the communication in debug mode to support
647 // having an assertion?
649 // Note that these variables describe only their own node.
650 bool useGpuForNonbonded = false;
651 bool useGpuForPme = false;
654 // It's possible that there are different numbers of GPUs on
655 // different nodes, which is the user's responsibilty to
656 // handle. If unsuitable, we will notice that during task
658 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
659 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
660 emulateGpuNonbonded, inputrec->cutoff_scheme == ecutsVERLET,
661 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
663 auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
664 auto canUseGpuForPme = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
665 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
666 canUseGpuForPme, cr->nnodes, domdecOptions.numPmeRanks,
669 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
670 if (pmeRunMode == PmeRunMode::GPU)
672 if (pmeFftTarget == TaskTarget::Cpu)
674 pmeRunMode = PmeRunMode::Mixed;
677 else if (pmeFftTarget == TaskTarget::Gpu)
679 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.");
682 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
684 // TODO: Error handling
685 mdModules.assignOptionsToModules(*inputrec->params, nullptr);
687 if (fplog != nullptr)
689 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
690 fprintf(fplog, "\n");
695 /* now make sure the state is initialized and propagated */
696 set_state_entries(globalState.get(), inputrec);
699 /* NM and TPI parallelize over force/energy calculations, not atoms,
700 * so we need to initialize and broadcast the global state.
702 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
706 globalState = std::unique_ptr<t_state>(new t_state);
708 broadcastStateWithoutDynamics(cr, globalState.get());
711 /* A parallel command line option consistency check that we can
712 only do after any threads have started. */
713 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
714 domdecOptions.numCells[YY] > 1 ||
715 domdecOptions.numCells[ZZ] > 1 ||
716 domdecOptions.numPmeRanks > 0))
719 "The -dd or -npme option request a parallel simulation, "
721 "but %s was compiled without threads or MPI enabled"
724 "but the number of MPI-threads (option -ntmpi) is not set or is 1"
726 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
729 , output_env_get_program_display_name(oenv)
734 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
736 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
739 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
741 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
744 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
746 if (domdecOptions.numPmeRanks > 0)
748 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
749 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
752 domdecOptions.numPmeRanks = 0;
755 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
757 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
758 * improve performance with many threads per GPU, since our OpenMP
759 * scaling is bad, but it's difficult to automate the setup.
761 domdecOptions.numPmeRanks = 0;
765 if (domdecOptions.numPmeRanks < 0)
767 domdecOptions.numPmeRanks = 0;
768 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
772 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
779 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
783 /* NMR restraints must be initialized before load_checkpoint,
784 * since with time averaging the history is added to t_state.
785 * For proper consistency check we therefore need to extend
787 * So the PME-only nodes (if present) will also initialize
788 * the distance restraints.
792 /* This needs to be called before read_checkpoint to extend the state */
793 init_disres(fplog, mtop, inputrec, cr, fcd, globalState.get(), replExParams.exchangeInterval > 0);
795 init_orires(fplog, mtop, inputrec, cr, globalState.get(), &(fcd->orires));
797 if (inputrecDeform(inputrec))
799 /* Store the deform reference box before reading the checkpoint */
802 copy_mat(globalState->box, box);
806 gmx_bcast(sizeof(box), box, cr);
808 /* Because we do not have the update struct available yet
809 * in which the reference values should be stored,
810 * we store them temporarily in static variables.
811 * This should be thread safe, since they are only written once
812 * and with identical values.
814 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
815 deform_init_init_step_tpx = inputrec->init_step;
816 copy_mat(box, deform_init_box_tpx);
817 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
820 ObservablesHistory observablesHistory = {};
822 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
824 if (continuationOptions.startedFromCheckpoint)
826 /* Check if checkpoint file exists before doing continuation.
827 * This way we can use identical input options for the first and subsequent runs...
831 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
832 cr, domdecOptions.numCells,
833 inputrec, globalState.get(),
834 &bReadEkin, &observablesHistory,
835 continuationOptions.appendFiles,
836 continuationOptions.appendFilesOptionSet,
837 mdrunOptions.reproducible);
841 continuationOptions.haveReadEkin = true;
845 if (SIMMASTER(cr) && continuationOptions.appendFiles)
847 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
848 continuationOptions.appendFiles, &fplog);
849 logOwner = buildLogger(fplog, nullptr);
850 mdlog = logOwner.logger();
853 if (mdrunOptions.numStepsCommandline > -2)
855 GMX_LOG(mdlog.info).asParagraph().
856 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
857 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
859 /* override nsteps with value set on the commamdline */
860 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
864 copy_mat(globalState->box, box);
869 gmx_bcast(sizeof(box), box, cr);
872 /* Update rlist and nstlist. */
873 if (inputrec->cutoff_scheme == ecutsVERLET)
875 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, mtop, box,
876 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
879 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
880 inputrec->eI == eiNM))
882 const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
884 cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
887 &ddbox, &npme_major, &npme_minor);
888 // Note that local state still does not exist yet.
892 /* PME, if used, is done on all nodes with 1D decomposition */
894 cr->duty = (DUTY_PP | DUTY_PME);
898 if (inputrec->ePBC == epbcSCREW)
901 "pbc=%s is only implemented with domain decomposition",
902 epbc_names[inputrec->ePBC]);
908 /* After possible communicator splitting in make_dd_communicators.
909 * we can set up the intra/inter node communication.
911 gmx_setup_nodecomm(fplog, cr);
914 /* Initialize per-physical-node MPI process/thread ID and counters. */
915 gmx_init_intranode_counters(cr);
916 if (cr->ms && cr->ms->nsim > 1 && !opt2bSet("-multidir", nfile, fnm))
918 GMX_LOG(mdlog.info).asParagraph().
919 appendText("The -multi flag is deprecated, and may be removed in a future version. Please "
920 "update your workflows to use -multidir instead.");
925 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
926 "This is simulation %d out of %d running as a composite GROMACS\n"
927 "multi-simulation job. Setup for this simulation:\n",
928 cr->ms->sim, cr->ms->nsim);
930 GMX_LOG(mdlog.warning).appendTextFormatted(
934 cr->nnodes == 1 ? "thread" : "threads"
936 cr->nnodes == 1 ? "process" : "processes"
942 /* Check and update hw_opt for the cut-off scheme */
943 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
945 /* Check and update the number of OpenMP threads requested */
946 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, pmeRunMode, *mtop);
948 gmx_omp_nthreads_init(mdlog, cr,
949 hwinfo->nthreads_hw_avail,
951 hw_opt.nthreads_omp_pme,
952 !thisRankHasDuty(cr, DUTY_PP),
953 inputrec->cutoff_scheme == ecutsVERLET);
956 if (EI_TPI(inputrec->eI) &&
957 inputrec->cutoff_scheme == ecutsVERLET)
959 gmx_feenableexcept();
963 // Build a data structure that expresses which kinds of non-bonded
964 // task are handled by this rank.
966 // TODO Later, this might become a loop over all registered modules
967 // relevant to the mdp inputs, to find those that have such tasks.
969 // TODO This could move before init_domain_decomposition() as part
970 // of refactoring that separates the responsibility for duty
971 // assignment from setup for communication between tasks, and
972 // setup for tasks handled with a domain (ie including short-ranged
973 // tasks, bonded tasks, etc.).
975 // Note that in general useGpuForNonbonded, etc. can have a value
976 // that is inconsistent with the presence of actual GPUs on any
977 // rank, and that is not known to be a problem until the
978 // duty of the ranks on a node become node.
980 // TODO Later we might need the concept of computeTasksOnThisRank,
981 // from which we construct gpuTasksOnThisRank.
983 // Currently the DD code assigns duty to ranks that can
984 // include PP work that currently can be executed on a single
985 // GPU, if present and compatible. This has to be coordinated
986 // across PP ranks on a node, with possible multiple devices
987 // or sharing devices on a node, either from the user
988 // selection, or automatically.
989 auto haveGpus = !gpuIdsToUse.empty();
990 std::vector<GpuTask> gpuTasksOnThisRank;
991 if (thisRankHasDuty(cr, DUTY_PP))
993 if (useGpuForNonbonded)
997 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
999 else if (nonbondedTarget == TaskTarget::Gpu)
1001 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1005 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1006 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1012 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1014 else if (pmeTarget == TaskTarget::Gpu)
1016 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1021 GpuTaskAssignment gpuTaskAssignment;
1024 // Produce the task assignment for this rank.
1025 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1026 mdlog, cr, gpuTasksOnThisRank);
1028 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1030 /* Prevent other ranks from continuing after an issue was found
1031 * and reported as a fatal error.
1033 * TODO This function implements a barrier so that MPI runtimes
1034 * can organize an orderly shutdown if one of the ranks has had to
1035 * issue a fatal error in various code already run. When we have
1036 * MPI-aware error handling and reporting, this should be
1041 MPI_Barrier(cr->mpi_comm_mysim);
1047 MPI_Barrier(cr->ms->mpi_comm_masters);
1049 /* We need another barrier to prevent non-master ranks from contiuing
1050 * when an error occured in a different simulation.
1052 MPI_Barrier(cr->mpi_comm_mysim);
1056 /* Now that we know the setup is consistent, check for efficiency */
1057 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1060 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1062 if (thisRankHasDuty(cr, DUTY_PP))
1064 // This works because only one task of each type is currently permitted.
1065 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1066 hasTaskType<GpuTask::Nonbonded>);
1067 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1069 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1070 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1071 init_gpu(mdlog, nonbondedDeviceInfo);
1073 if (DOMAINDECOMP(cr))
1075 /* When we share GPUs over ranks, we need to know this for the DLB */
1076 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1082 gmx_device_info_t *pmeDeviceInfo = nullptr;
1083 // This works because only one task of each type is currently permitted.
1084 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1085 if (pmeGpuTaskMapping != gpuTaskAssignment.end())
1087 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1088 init_gpu(mdlog, pmeDeviceInfo);
1091 /* getting number of PP/PME threads
1092 PME: env variable should be read only on one node to make sure it is
1093 identical everywhere;
1095 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1097 int numThreadsOnThisRank;
1098 /* threads on this MPI process or TMPI thread */
1099 if (thisRankHasDuty(cr, DUTY_PP))
1101 numThreadsOnThisRank = gmx_omp_nthreads_get(emntNonbonded);
1105 numThreadsOnThisRank = nthreads_pme;
1108 checkHardwareOversubscription(numThreadsOnThisRank,
1109 *hwinfo->hardwareTopology,
1112 if (hw_opt.thread_affinity != threadaffOFF)
1114 /* Before setting affinity, check whether the affinity has changed
1115 * - which indicates that probably the OpenMP library has changed it
1116 * since we first checked).
1118 gmx_check_thread_affinity_set(mdlog, cr,
1119 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1121 /* Set the CPU affinity */
1122 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1123 numThreadsOnThisRank, nullptr);
1126 if (mdrunOptions.timingOptions.resetStep > -1)
1128 GMX_LOG(mdlog.info).asParagraph().
1129 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1131 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1135 /* Master synchronizes its value of reset_counters with all nodes
1136 * including PME only nodes */
1137 reset_counters = wcycle_get_reset_counters(wcycle);
1138 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1139 wcycle_set_reset_counters(wcycle, reset_counters);
1142 // Membrane embedding must be initialized before we call init_forcerec()
1147 fprintf(stderr, "Initializing membed");
1149 /* Note that membed cannot work in parallel because mtop is
1150 * changed here. Fix this if we ever want to make it run with
1151 * multiple ranks. */
1152 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1155 std::unique_ptr<MDAtoms> mdAtoms;
1158 if (thisRankHasDuty(cr, DUTY_PP))
1160 /* Initiate forcerecord */
1162 fr->forceProviders = mdModules.initForceProviders();
1163 init_forcerec(fplog, mdlog, fr, fcd,
1164 inputrec, mtop, cr, box,
1165 opt2fn("-table", nfile, fnm),
1166 opt2fn("-tablep", nfile, fnm),
1167 getFilenm("-tableb", nfile, fnm),
1168 *hwinfo, nonbondedDeviceInfo,
1172 /* Initialize QM-MM */
1175 GMX_LOG(mdlog.info).asParagraph().
1176 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
1177 "version. Please get in touch with the developers if you find the support useful, "
1178 "as help is needed if the functionality is to continue to be available.");
1179 init_QMMMrec(cr, mtop, inputrec, fr);
1182 /* Initialize the mdAtoms structure.
1183 * mdAtoms is not filled with atom data,
1184 * as this can not be done now with domain decomposition.
1186 const bool useGpuForPme = (pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed);
1187 mdAtoms = makeMDAtoms(fplog, *mtop, *inputrec, useGpuForPme && thisRankHasDuty(cr, DUTY_PME));
1190 // The pinning of coordinates in the global state object works, because we only use
1191 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1192 // points to the global state object without DD.
1193 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1194 // which should also perform the pinning.
1195 changePinningPolicy(&globalState->x, useGpuForPme ? PinningPolicy::CanBePinned : PinningPolicy::CannotBePinned);
1198 /* Initialize the virtual site communication */
1199 vsite = initVsite(*mtop, cr);
1201 calc_shifts(box, fr->shift_vec);
1203 /* With periodic molecules the charge groups should be whole at start up
1204 * and the virtual sites should not be far from their proper positions.
1206 if (!inputrec->bContinuation && MASTER(cr) &&
1207 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1209 /* Make molecules whole at start of run */
1210 if (fr->ePBC != epbcNONE)
1212 rvec *xGlobal = as_rvec_array(globalState->x.data());
1213 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, xGlobal);
1217 /* Correct initial vsite positions are required
1218 * for the initial distribution in the domain decomposition
1219 * and for the initial shell prediction.
1221 constructVsitesGlobal(*mtop, globalState->x);
1225 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1227 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1228 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1233 /* This is a PME only node */
1235 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1237 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1238 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1241 gmx_pme_t *sepPmeData = nullptr;
1242 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1243 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1244 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1246 /* Initiate PME if necessary,
1247 * either on all nodes or on dedicated PME nodes only. */
1248 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1250 if (mdAtoms && mdAtoms->mdatoms())
1252 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1253 if (EVDW_PME(inputrec->vdwtype))
1255 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1258 if (cr->npmenodes > 0)
1260 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1261 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1262 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1265 if (thisRankHasDuty(cr, DUTY_PME))
1269 pmedata = gmx_pme_init(cr, npme_major, npme_minor, inputrec,
1270 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1271 mdrunOptions.reproducible,
1272 ewaldcoeff_q, ewaldcoeff_lj,
1274 pmeRunMode, nullptr, pmeDeviceInfo, mdlog);
1276 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1281 if (EI_DYNAMICS(inputrec->eI))
1283 /* Turn on signal handling on all nodes */
1285 * (A user signal from the PME nodes (if any)
1286 * is communicated to the PP nodes.
1288 signal_handler_install();
1291 if (thisRankHasDuty(cr, DUTY_PP))
1293 /* Assumes uniform use of the number of OpenMP threads */
1294 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1296 if (inputrec->bPull)
1298 /* Initialize pull code */
1299 inputrec->pull_work =
1300 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1301 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1302 EI_DYNAMICS(inputrec->eI) && MASTER(cr),
1303 continuationOptions);
1308 /* Initialize enforced rotation code */
1309 init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), mtop, oenv, mdrunOptions);
1312 /* Let init_constraints know whether we have essential dynamics constraints.
1313 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1315 bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1317 constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1319 if (DOMAINDECOMP(cr))
1321 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1322 /* This call is not included in init_domain_decomposition mainly
1323 * because fr->cginfo_mb is set later.
1325 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1326 domdecOptions.checkBondedInteractions,
1330 /* Now do whatever the user wants us to do (how flexible...) */
1331 my_integrator(inputrec->eI) (fplog, cr, mdlog, nfile, fnm,
1335 mdModules.outputProvider(),
1339 &observablesHistory,
1340 mdAtoms.get(), nrnb, wcycle, fr,
1343 walltime_accounting);
1347 finish_rot(inputrec->rot);
1350 if (inputrec->bPull)
1352 finish_pull(inputrec->pull_work);
1358 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1360 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1361 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1364 wallcycle_stop(wcycle, ewcRUN);
1366 /* Finish up, write some stuff
1367 * if rerunMD, don't write last frame again
1369 finish_run(fplog, mdlog, cr,
1370 inputrec, nrnb, wcycle, walltime_accounting,
1371 fr ? fr->nbv : nullptr,
1373 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1378 gmx_pme_destroy(pmedata);
1382 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1383 // before we destroy the GPU context(s) in free_gpu_resources().
1384 // Pinned buffers are associated with contexts in CUDA.
1385 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1386 mdAtoms.reset(nullptr);
1387 globalState.reset(nullptr);
1389 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1390 free_gpu_resources(fr, cr);
1391 free_gpu(nonbondedDeviceInfo);
1392 free_gpu(pmeDeviceInfo);
1396 free_membed(membed);
1399 gmx_hardware_info_free();
1401 /* Does what it says */
1402 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1403 walltime_accounting_destroy(walltime_accounting);
1405 /* Close logfile already here if we were appending to it */
1406 if (MASTER(cr) && continuationOptions.appendFiles)
1408 gmx_log_close(fplog);
1411 rc = (int)gmx_get_stop_condition();
1414 /* we need to join all threads. The sub-threads join when they
1415 exit this function, but the master thread needs to be told to
1417 if (PAR(cr) && MASTER(cr))