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, 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_search.h"
88 #include "gromacs/mdlib/nbnxn_tuning.h"
89 #include "gromacs/mdlib/qmmm.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/sim_util.h"
92 #include "gromacs/mdlib/tpi.h"
93 #include "gromacs/mdrunutility/mdmodules.h"
94 #include "gromacs/mdrunutility/threadaffinity.h"
95 #include "gromacs/mdtypes/commrec.h"
96 #include "gromacs/mdtypes/inputrec.h"
97 #include "gromacs/mdtypes/md_enums.h"
98 #include "gromacs/mdtypes/observableshistory.h"
99 #include "gromacs/mdtypes/state.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/pulling/pull_rotation.h"
103 #include "gromacs/taskassignment/hardwareassign.h"
104 #include "gromacs/taskassignment/resourcedivision.h"
105 #include "gromacs/taskassignment/usergpuids.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/topology/mtop_util.h"
108 #include "gromacs/trajectory/trajectoryframe.h"
109 #include "gromacs/utility/cstringutil.h"
110 #include "gromacs/utility/exceptions.h"
111 #include "gromacs/utility/fatalerror.h"
112 #include "gromacs/utility/filestream.h"
113 #include "gromacs/utility/gmxassert.h"
114 #include "gromacs/utility/gmxmpi.h"
115 #include "gromacs/utility/logger.h"
116 #include "gromacs/utility/loggerbuilder.h"
117 #include "gromacs/utility/pleasecite.h"
118 #include "gromacs/utility/programcontext.h"
119 #include "gromacs/utility/smalloc.h"
127 #include "corewrap.h"
130 //! First step used in pressure scaling
131 gmx_int64_t deform_init_init_step_tpx;
132 //! Initial box for pressure scaling
133 matrix deform_init_box_tpx;
134 //! MPI variable for use in pressure scaling
135 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
138 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
139 * the number of threads will get lowered.
141 #define MIN_ATOMS_PER_MPI_THREAD 90
142 #define MIN_ATOMS_PER_GPU 900
147 void Mdrunner::reinitializeOnSpawnedThread()
149 // TODO This duplication is formally necessary if any thread might
150 // modify any memory in fnm or the pointers it contains. If the
151 // contents are ever provably const, then we can remove this
152 // allocation (and memory leak).
153 // TODO This should probably become part of a copy constructor for
155 fnm = dup_tfn(nfile, fnm);
157 cr = reinitialize_commrec_for_this_thread(cr);
161 // Only the master rank writes to the log files
166 /*! \brief The callback used for running on spawned threads.
168 * Obtains the pointer to the master mdrunner object from the one
169 * argument permitted to the thread-launch API call, copies it to make
170 * a new runner for this thread, reinitializes necessary data, and
171 * proceeds to the simulation. */
172 static void mdrunner_start_fn(void *arg)
176 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
177 /* copy the arg list to make sure that it's thread-local. This
178 doesn't copy pointed-to items, of course, but those are all
180 gmx::Mdrunner mdrunner = *masterMdrunner;
181 mdrunner.reinitializeOnSpawnedThread();
184 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
188 /*! \brief Start thread-MPI threads.
190 * Called by mdrunner() to start a specific number of threads
191 * (including the main thread) for thread-parallel runs. This in turn
192 * calls mdrunner() for each thread. All options are the same as for
194 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch)
197 /* first check whether we even need to start tMPI */
198 if (numThreadsToLaunch < 2)
203 gmx::Mdrunner spawnedMdrunner = *this;
204 // TODO This duplication is formally necessary if any thread might
205 // modify any memory in fnm or the pointers it contains. If the
206 // contents are ever provably const, then we can remove this
207 // allocation (and memory leak).
208 // TODO This should probably become part of a copy constructor for
210 spawnedMdrunner.fnm = dup_tfn(this->nfile, fnm);
212 /* now spawn new threads that start mdrunner_start_fn(), while
213 the main thread returns, we set thread affinity later */
214 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
215 mdrunner_start_fn, static_cast<void*>(&spawnedMdrunner)) != TMPI_SUCCESS)
217 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
220 return reinitialize_commrec_for_this_thread(cr);
225 #endif /* GMX_THREAD_MPI */
227 /*! \brief Initialize variables for Verlet scheme simulation */
228 static void prepare_verlet_scheme(FILE *fplog,
232 const gmx_mtop_t *mtop,
234 bool makeGpuPairList,
235 const gmx::CpuInfo &cpuinfo)
237 /* For NVE simulations, we will retain the initial list buffer */
238 if (EI_DYNAMICS(ir->eI) &&
239 ir->verletbuf_tol > 0 &&
240 !(EI_MD(ir->eI) && ir->etc == etcNO))
242 /* Update the Verlet buffer size for the current run setup */
244 /* Here we assume SIMD-enabled kernels are being used. But as currently
245 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
246 * and 4x2 gives a larger buffer than 4x4, this is ok.
248 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
249 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
252 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &listSetup, nullptr, &rlist_new);
254 if (rlist_new != ir->rlist)
256 if (fplog != nullptr)
258 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
259 ir->rlist, rlist_new,
260 listSetup.cluster_size_i, listSetup.cluster_size_j);
262 ir->rlist = rlist_new;
266 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
268 gmx_fatal(FARGS, "Can not set nstlist without %s",
269 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
272 if (EI_DYNAMICS(ir->eI))
274 /* Set or try nstlist values */
275 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
279 /*! \brief Override the nslist value in inputrec
281 * with value passed on the command line (if any)
283 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
284 gmx_int64_t nsteps_cmdline,
289 /* override with anything else than the default -2 */
290 if (nsteps_cmdline > -2)
292 char sbuf_steps[STEPSTRSIZE];
293 char sbuf_msg[STRLEN];
295 ir->nsteps = nsteps_cmdline;
296 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
298 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
299 gmx_step_str(nsteps_cmdline, sbuf_steps),
300 fabs(nsteps_cmdline*ir->delta_t));
304 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
305 gmx_step_str(nsteps_cmdline, sbuf_steps));
308 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
310 else if (nsteps_cmdline < -2)
312 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
315 /* Do nothing if nsteps_cmdline == -2 */
321 //! Halt the run if there are inconsistences between user choices to run with GPUs and/or hardware detection.
322 static void exitIfCannotForceGpuRun(bool requirePhysicalGpu,
323 EmulateGpuNonbonded emulateGpuNonbonded,
324 bool useVerletScheme,
325 bool compatibleGpusFound)
327 /* Was GPU acceleration either explicitly (-nb gpu) or implicitly
328 * (gpu ID passed) requested? */
329 if (!requirePhysicalGpu)
334 if (GMX_GPU == GMX_GPU_NONE)
336 gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!",
337 gmx::getProgramContext().displayName());
340 if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
342 gmx_fatal(FARGS, "GPU emulation cannot be requested together with GPU acceleration!");
345 if (!useVerletScheme)
347 gmx_fatal(FARGS, "GPU acceleration requested, but can't be used without cutoff-scheme=Verlet");
350 if (!compatibleGpusFound)
352 gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
356 /*! \brief Return whether GPU acceleration is useful with the given settings.
358 * If not, logs a message about falling back to CPU code. */
359 static bool gpuAccelerationIsUseful(const MDLogger &mdlog,
360 const t_inputrec *ir,
363 if (doRerun && ir->opts.ngener > 1)
365 /* Rerun execution time is dominated by I/O and pair search,
366 * so GPUs are not very useful, plus they do not support more
367 * than one energy group. If the user requested GPUs
368 * explicitly, a fatal error is given later. With non-reruns,
369 * we fall back to a single whole-of system energy group
370 * (which runs much faster than a multiple-energy-groups
371 * implementation would), and issue a note in the .log
372 * file. Users can re-run if they want the information. */
373 GMX_LOG(mdlog.warning).asParagraph().appendText("Multiple energy groups is not implemented for GPUs, so is not useful for this rerun, so falling back to the CPU");
380 //! \brief Return the correct integrator function.
381 static integrator_t *my_integrator(unsigned int ei)
390 if (!EI_DYNAMICS(ei))
392 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
407 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
411 GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
413 GMX_THROW(APIError("Non existing integrator selected"));
417 //! Initializes the logger for mdrun.
418 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
420 gmx::LoggerBuilder builder;
421 if (fplog != nullptr)
423 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
425 if (cr == nullptr || SIMMASTER(cr))
427 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
428 &gmx::TextOutputFile::standardError());
430 return builder.build();
433 int Mdrunner::mdrunner()
436 gmx_ddbox_t ddbox = {0};
437 int npme_major, npme_minor;
439 gmx_mtop_t *mtop = nullptr;
440 t_forcerec *fr = nullptr;
441 t_fcdata *fcd = nullptr;
442 real ewaldcoeff_q = 0;
443 real ewaldcoeff_lj = 0;
444 gmx_vsite_t *vsite = nullptr;
446 int nChargePerturbed = -1, nTypePerturbed = 0;
447 gmx_wallcycle_t wcycle;
448 gmx_walltime_accounting_t walltime_accounting = nullptr;
450 gmx_int64_t reset_counters;
451 int nthreads_pme = 1;
452 gmx_membed_t * membed = nullptr;
453 gmx_hw_info_t *hwinfo = nullptr;
455 /* CAUTION: threads may be started later on in this function, so
456 cr doesn't reflect the final parallel state right now */
457 gmx::MDModules mdModules;
458 t_inputrec inputrecInstance;
459 t_inputrec *inputrec = &inputrecInstance;
462 if (mdrunOptions.continuationOptions.appendFiles)
467 bool doMembed = opt2bSet("-membed", nfile, fnm);
468 bool doRerun = mdrunOptions.rerun;
470 /* Handle GPU-related user options. Later, we check consistency
471 * with things like whether support is compiled, or tMPI thread
473 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
474 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
475 std::vector<int> userGpuIds;
478 userGpuIds = parseUserGpuIds(hw_opt.gpuIdTaskAssignment);
480 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
482 bool forceUseCpu = (strncmp(nbpu_opt, "cpu", 3) == 0);
483 if (!userGpuIds.empty() && forceUseCpu)
485 gmx_fatal(FARGS, "GPU IDs were specified, and short-ranged interactions were assigned to the CPU. Make no more than one of these choices.");
487 bool forceUsePhysicalGpu = (strncmp(nbpu_opt, "gpu", 3) == 0) || !userGpuIds.empty();
488 bool tryUsePhysicalGpu = (strncmp(nbpu_opt, "auto", 4) == 0) && userGpuIds.empty() && (emulateGpuNonbonded == EmulateGpuNonbonded::No);
489 GMX_RELEASE_ASSERT(!(forceUsePhysicalGpu && tryUsePhysicalGpu), "Must either force use of "
490 "GPUs for short-ranged interactions, or try to use them, not both.");
491 const PmeRunMode pmeRunMode = PmeRunMode::CPU;
492 //TODO this is a placeholder as PME on GPU is not permitted yet
493 //TODO should there exist a PmeRunMode::None value for consistency?
495 // Here we assume that SIMMASTER(cr) does not change even after the
496 // threads are started.
497 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
498 gmx::MDLogger mdlog(logOwner.logger());
500 hwinfo = gmx_detect_hardware(mdlog, cr);
502 gmx_print_detected_hardware(fplog, cr, mdlog, hwinfo);
504 if (fplog != nullptr)
506 /* Print references after all software/hardware printing */
507 please_cite(fplog, "Abraham2015");
508 please_cite(fplog, "Pall2015");
509 please_cite(fplog, "Pronk2013");
510 please_cite(fplog, "Hess2008b");
511 please_cite(fplog, "Spoel2005a");
512 please_cite(fplog, "Lindahl2001a");
513 please_cite(fplog, "Berendsen95a");
516 std::unique_ptr<t_state> globalState;
520 /* Only the master rank has the global state */
521 globalState = std::unique_ptr<t_state>(new t_state);
523 /* Read (nearly) all data required for the simulation */
524 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), mtop);
526 exitIfCannotForceGpuRun(forceUsePhysicalGpu,
528 inputrec->cutoff_scheme == ecutsVERLET,
529 compatibleGpusFound(hwinfo->gpu_info));
531 if (inputrec->cutoff_scheme == ecutsVERLET)
533 /* TODO This logic could run later, e.g. before -npme -1
534 is handled. If inputrec has already been communicated,
535 then the resulting tryUsePhysicalGpu does not need to
537 if ((tryUsePhysicalGpu || forceUsePhysicalGpu) &&
538 !gpuAccelerationIsUseful(mdlog, inputrec, doRerun))
540 /* Fallback message printed by nbnxn_acceleration_supported */
541 if (forceUsePhysicalGpu)
543 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
545 tryUsePhysicalGpu = false;
550 if (nstlist_cmdline > 0)
552 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
555 if (compatibleGpusFound(hwinfo->gpu_info))
557 GMX_LOG(mdlog.warning).asParagraph().appendText(
558 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
559 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
561 tryUsePhysicalGpu = false;
564 bool nonbondedOnGpu = (tryUsePhysicalGpu || forceUsePhysicalGpu) && compatibleGpusFound(hwinfo->gpu_info);
566 /* Check and update the hardware options for internal consistency */
567 check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
569 /* Early check for externally set process affinity. */
570 gmx_check_thread_affinity_set(mdlog, cr,
571 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
576 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
578 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
581 /* Since the master knows the cut-off scheme, update hw_opt for this.
582 * This is done later for normal MPI and also once more with tMPI
583 * for all tMPI ranks.
585 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
587 /* Determine how many thread-MPI ranks to start.
589 * TODO Over-writing the user-supplied value here does
590 * prevent any possible subsequent checks from working
592 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
595 domdecOptions.numPmeRanks,
601 // Now start the threads for thread MPI.
602 cr = spawnThreads(hw_opt.nthreads_tmpi);
603 /* The main thread continues here with a new cr. We don't deallocate
604 the old cr because other threads may still be reading it. */
605 // TODO Both master and spawned threads call dup_tfn and
606 // reinitialize_commrec_for_this_thread. Find a way to express
610 /* END OF CAUTION: cr is now reliable */
614 /* now broadcast everything to the non-master nodes/threads: */
615 init_parallel(cr, inputrec, mtop);
617 gmx_bcast_sim(sizeof(nonbondedOnGpu), &nonbondedOnGpu, cr);
619 // TODO: Error handling
620 mdModules.assignOptionsToModules(*inputrec->params, nullptr);
622 if (fplog != nullptr)
624 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
625 fprintf(fplog, "\n");
630 /* now make sure the state is initialized and propagated */
631 set_state_entries(globalState.get(), inputrec);
634 /* NM and TPI parallelize over force/energy calculations, not atoms,
635 * so we need to initialize and broadcast the global state.
637 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
641 globalState = std::unique_ptr<t_state>(new t_state);
643 broadcastStateWithoutDynamics(cr, globalState.get());
646 /* A parallel command line option consistency check that we can
647 only do after any threads have started. */
648 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
649 domdecOptions.numCells[YY] > 1 ||
650 domdecOptions.numCells[ZZ] > 1 ||
651 domdecOptions.numPmeRanks > 0))
654 "The -dd or -npme option request a parallel simulation, "
656 "but %s was compiled without threads or MPI enabled"
659 "but the number of MPI-threads (option -ntmpi) is not set or is 1"
661 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
664 , output_env_get_program_display_name(oenv)
669 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
671 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
674 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
676 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
679 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
681 if (domdecOptions.numPmeRanks > 0)
683 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
684 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
687 domdecOptions.numPmeRanks = 0;
690 if (nonbondedOnGpu && domdecOptions.numPmeRanks < 0)
692 /* With GPUs we don't automatically use PME-only ranks. PME ranks can
693 * improve performance with many threads per GPU, since our OpenMP
694 * scaling is bad, but it's difficult to automate the setup.
696 domdecOptions.numPmeRanks = 0;
702 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
706 /* NMR restraints must be initialized before load_checkpoint,
707 * since with time averaging the history is added to t_state.
708 * For proper consistency check we therefore need to extend
710 * So the PME-only nodes (if present) will also initialize
711 * the distance restraints.
715 /* This needs to be called before read_checkpoint to extend the state */
716 init_disres(fplog, mtop, inputrec, cr, fcd, globalState.get(), replExParams.exchangeInterval > 0);
718 init_orires(fplog, mtop, inputrec, cr, globalState.get(), &(fcd->orires));
720 if (inputrecDeform(inputrec))
722 /* Store the deform reference box before reading the checkpoint */
725 copy_mat(globalState->box, box);
729 gmx_bcast(sizeof(box), box, cr);
731 /* Because we do not have the update struct available yet
732 * in which the reference values should be stored,
733 * we store them temporarily in static variables.
734 * This should be thread safe, since they are only written once
735 * and with identical values.
737 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
738 deform_init_init_step_tpx = inputrec->init_step;
739 copy_mat(box, deform_init_box_tpx);
740 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
743 ObservablesHistory observablesHistory = {};
745 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
747 if (continuationOptions.startedFromCheckpoint)
749 /* Check if checkpoint file exists before doing continuation.
750 * This way we can use identical input options for the first and subsequent runs...
754 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
755 cr, domdecOptions.numCells,
756 inputrec, globalState.get(),
757 &bReadEkin, &observablesHistory,
758 continuationOptions.appendFiles,
759 continuationOptions.appendFilesOptionSet,
760 mdrunOptions.reproducible);
764 continuationOptions.haveReadEkin = true;
768 if (SIMMASTER(cr) && continuationOptions.appendFiles)
770 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
771 continuationOptions.appendFiles, &fplog);
772 logOwner = buildLogger(fplog, nullptr);
773 mdlog = logOwner.logger();
776 /* override nsteps with value set on the commamdline */
777 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
781 copy_mat(globalState->box, box);
786 gmx_bcast(sizeof(box), box, cr);
789 /* Update rlist and nstlist. */
790 if (inputrec->cutoff_scheme == ecutsVERLET)
792 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, mtop, box,
793 nonbondedOnGpu || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
796 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
797 inputrec->eI == eiNM))
799 const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
801 cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
804 &ddbox, &npme_major, &npme_minor);
805 // Note that local state still does not exist yet.
809 /* PME, if used, is done on all nodes with 1D decomposition */
811 cr->duty = (DUTY_PP | DUTY_PME);
815 if (inputrec->ePBC == epbcSCREW)
818 "pbc=%s is only implemented with domain decomposition",
819 epbc_names[inputrec->ePBC]);
825 /* After possible communicator splitting in make_dd_communicators.
826 * we can set up the intra/inter node communication.
828 gmx_setup_nodecomm(fplog, cr);
831 /* Initialize per-physical-node MPI process/thread ID and counters. */
832 gmx_init_intranode_counters(cr);
836 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
837 "This is simulation %d out of %d running as a composite GROMACS\n"
838 "multi-simulation job. Setup for this simulation:\n",
839 cr->ms->sim, cr->ms->nsim);
841 GMX_LOG(mdlog.warning).appendTextFormatted(
845 cr->nnodes == 1 ? "thread" : "threads"
847 cr->nnodes == 1 ? "process" : "processes"
853 /* Check and update hw_opt for the cut-off scheme */
854 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
856 /* Check and update hw_opt for the number of MPI ranks */
857 check_and_update_hw_opt_3(&hw_opt);
859 gmx_omp_nthreads_init(mdlog, cr,
860 hwinfo->nthreads_hw_avail,
862 hw_opt.nthreads_omp_pme,
863 !thisRankHasDuty(cr, DUTY_PP),
864 inputrec->cutoff_scheme == ecutsVERLET);
867 if (EI_TPI(inputrec->eI) &&
868 inputrec->cutoff_scheme == ecutsVERLET)
870 gmx_feenableexcept();
874 // Contains the ID of the GPU used by each PP rank on this node,
875 // indexed by that rank. Empty if no GPUs are selected for use on
877 std::vector<int> gpuTaskAssignment;
880 // Currently the DD code assigns duty to ranks that can
881 // include PP work that currently can be executed on a single
882 // GPU, if present and compatible. This has to be coordinated
883 // across PP ranks on a node, with possible multiple devices
884 // or sharing devices on a node, either from the user
885 // selection, or automatically.
887 // GPU ID assignment strings, if provided, cover all the ranks on
888 // a node. If nodes or the process placement on them are
889 // heterogeneous, then the GMX_GPU_ID environment variable must be
890 // set by a user who also wishes to direct GPU ID assignment.
891 // Thus the implementation of task assignment can assume it has a
892 // GPU ID assignment appropriate for the node upon which its
893 // process is running.
895 // Valid GPU ID assignments are an ordered set of digits that
896 // identify GPU device IDs (e.g. as understood by the GPU runtime,
897 // and subject to environment modification such as with
898 // CUDA_VISIBLE_DEVICES) that will be used for the GPU-suitable
899 // tasks on all of the ranks of that node.
900 bool rankCanUseGpu = thisRankHasDuty(cr, DUTY_PP);
901 gpuTaskAssignment = mapPpRanksToGpus(rankCanUseGpu, cr, hwinfo->gpu_info, hwinfo->compatibleGpus, userGpuIds);
904 reportGpuUsage(mdlog, hwinfo->gpu_info, !userGpuIds.empty(),
905 gpuTaskAssignment, cr->nrank_pp_intranode, cr->nnodes > 1);
907 if (!gpuTaskAssignment.empty())
909 GMX_RELEASE_ASSERT(cr->nrank_pp_intranode == static_cast<int>(gpuTaskAssignment.size()),
910 "The number of PP ranks on each node must equal the number of GPU tasks used on each node");
913 /* Prevent other ranks from continuing after an issue was found
914 * and reported as a fatal error.
916 * TODO This function implements a barrier so that MPI runtimes
917 * can organize an orderly shutdown if one of the ranks has had to
918 * issue a fatal error in various code already run. When we have
919 * MPI-aware error handling and reporting, this should be
924 MPI_Barrier(cr->mpi_comm_mysim);
928 /* Now that we know the setup is consistent, check for efficiency */
929 check_resource_division_efficiency(hwinfo, hw_opt.nthreads_tot, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
932 gmx_device_info_t *shortRangedDeviceInfo = nullptr;
933 int shortRangedDeviceId = -1;
934 if (thisRankHasDuty(cr, DUTY_PP))
936 if (!gpuTaskAssignment.empty())
938 shortRangedDeviceId = gpuTaskAssignment[cr->rank_pp_intranode];
939 shortRangedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, shortRangedDeviceId);
943 if (DOMAINDECOMP(cr))
945 /* When we share GPUs over ranks, we need to know this for the DLB */
946 dd_setup_dlb_resource_sharing(cr, shortRangedDeviceId);
949 /* getting number of PP/PME threads
950 PME: env variable should be read only on one node to make sure it is
951 identical everywhere;
953 nthreads_pme = gmx_omp_nthreads_get(emntPME);
955 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
959 /* Master synchronizes its value of reset_counters with all nodes
960 * including PME only nodes */
961 reset_counters = wcycle_get_reset_counters(wcycle);
962 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
963 wcycle_set_reset_counters(wcycle, reset_counters);
966 // Membrane embedding must be initialized before we call init_forcerec()
971 fprintf(stderr, "Initializing membed");
973 /* Note that membed cannot work in parallel because mtop is
974 * changed here. Fix this if we ever want to make it run with
976 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
979 std::unique_ptr<MDAtoms> mdAtoms;
982 if (thisRankHasDuty(cr, DUTY_PP))
984 /* Initiate forcerecord */
986 fr->forceProviders = mdModules.initForceProviders();
987 init_forcerec(fplog, mdlog, fr, fcd,
988 inputrec, mtop, cr, box,
989 opt2fn("-table", nfile, fnm),
990 opt2fn("-tablep", nfile, fnm),
991 getFilenm("-tableb", nfile, fnm),
992 shortRangedDeviceInfo,
996 /* Initialize QM-MM */
999 init_QMMMrec(cr, mtop, inputrec, fr);
1002 /* Initialize the mdAtoms structure.
1003 * mdAtoms is not filled with atom data,
1004 * as this can not be done now with domain decomposition.
1006 mdAtoms = makeMDAtoms(fplog, *mtop, *inputrec);
1008 /* Initialize the virtual site communication */
1009 vsite = initVsite(*mtop, cr);
1011 calc_shifts(box, fr->shift_vec);
1013 /* With periodic molecules the charge groups should be whole at start up
1014 * and the virtual sites should not be far from their proper positions.
1016 if (!inputrec->bContinuation && MASTER(cr) &&
1017 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1019 /* Make molecules whole at start of run */
1020 if (fr->ePBC != epbcNONE)
1022 rvec *xGlobal = as_rvec_array(globalState->x.data());
1023 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, xGlobal);
1027 /* Correct initial vsite positions are required
1028 * for the initial distribution in the domain decomposition
1029 * and for the initial shell prediction.
1031 constructVsitesGlobal(*mtop, globalState->x);
1035 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1037 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1038 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1043 /* This is a PME only node */
1045 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1047 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1048 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1051 gmx_pme_t *sepPmeData = nullptr;
1052 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1053 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1054 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1056 if (hw_opt.thread_affinity != threadaffOFF)
1058 /* Before setting affinity, check whether the affinity has changed
1059 * - which indicates that probably the OpenMP library has changed it
1060 * since we first checked).
1062 gmx_check_thread_affinity_set(mdlog, cr,
1063 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1066 /* threads on this MPI process or TMPI thread */
1067 if (thisRankHasDuty(cr, DUTY_PP))
1069 nthread_local = gmx_omp_nthreads_get(emntNonbonded);
1073 nthread_local = gmx_omp_nthreads_get(emntPME);
1076 /* Set the CPU affinity */
1077 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1078 nthread_local, nullptr);
1081 /* Initiate PME if necessary,
1082 * either on all nodes or on dedicated PME nodes only. */
1083 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1085 if (mdAtoms && mdAtoms->mdatoms())
1087 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1088 if (EVDW_PME(inputrec->vdwtype))
1090 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1093 if (cr->npmenodes > 0)
1095 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1096 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1097 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1100 if (thisRankHasDuty(cr, DUTY_PME))
1104 gmx_device_info_t *pmeGpuInfo = nullptr;
1105 pmedata = gmx_pme_init(cr, npme_major, npme_minor, inputrec,
1106 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1107 mdrunOptions.reproducible,
1108 ewaldcoeff_q, ewaldcoeff_lj,
1110 pmeRunMode, nullptr, pmeGpuInfo, mdlog);
1112 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1117 if (EI_DYNAMICS(inputrec->eI))
1119 /* Turn on signal handling on all nodes */
1121 * (A user signal from the PME nodes (if any)
1122 * is communicated to the PP nodes.
1124 signal_handler_install();
1127 if (thisRankHasDuty(cr, DUTY_PP))
1129 /* Assumes uniform use of the number of OpenMP threads */
1130 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1132 if (inputrec->bPull)
1134 /* Initialize pull code */
1135 inputrec->pull_work =
1136 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1137 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1138 EI_DYNAMICS(inputrec->eI) && MASTER(cr),
1139 continuationOptions);
1144 /* Initialize enforced rotation code */
1145 init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), mtop, oenv, mdrunOptions);
1148 /* Let init_constraints know whether we have essential dynamics constraints.
1149 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1151 bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1153 constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1155 if (DOMAINDECOMP(cr))
1157 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1158 /* This call is not included in init_domain_decomposition mainly
1159 * because fr->cginfo_mb is set later.
1161 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1162 domdecOptions.checkBondedInteractions,
1166 /* Now do whatever the user wants us to do (how flexible...) */
1167 my_integrator(inputrec->eI) (fplog, cr, mdlog, nfile, fnm,
1171 mdModules.outputProvider(),
1175 &observablesHistory,
1176 mdAtoms.get(), nrnb, wcycle, fr,
1179 walltime_accounting);
1183 finish_rot(inputrec->rot);
1186 if (inputrec->bPull)
1188 finish_pull(inputrec->pull_work);
1194 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1196 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1197 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1200 wallcycle_stop(wcycle, ewcRUN);
1202 /* Finish up, write some stuff
1203 * if rerunMD, don't write last frame again
1205 finish_run(fplog, mdlog, cr,
1206 inputrec, nrnb, wcycle, walltime_accounting,
1207 fr ? fr->nbv : nullptr,
1209 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1214 gmx_pme_destroy(pmedata);
1218 /* Free GPU memory and context */
1219 free_gpu_resources(fr, cr, shortRangedDeviceInfo);
1223 free_membed(membed);
1226 gmx_hardware_info_free();
1228 /* Does what it says */
1229 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1230 walltime_accounting_destroy(walltime_accounting);
1232 /* Close logfile already here if we were appending to it */
1233 if (MASTER(cr) && continuationOptions.appendFiles)
1235 gmx_log_close(fplog);
1238 rc = (int)gmx_get_stop_condition();
1241 /* we need to join all threads. The sub-threads join when they
1242 exit this function, but the master thread needs to be told to
1244 if (PAR(cr) && MASTER(cr))