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/timing/wallcycle.h"
106 #include "gromacs/topology/mtop_util.h"
107 #include "gromacs/trajectory/trajectoryframe.h"
108 #include "gromacs/utility/cstringutil.h"
109 #include "gromacs/utility/exceptions.h"
110 #include "gromacs/utility/fatalerror.h"
111 #include "gromacs/utility/filestream.h"
112 #include "gromacs/utility/gmxassert.h"
113 #include "gromacs/utility/gmxmpi.h"
114 #include "gromacs/utility/logger.h"
115 #include "gromacs/utility/loggerbuilder.h"
116 #include "gromacs/utility/pleasecite.h"
117 #include "gromacs/utility/programcontext.h"
118 #include "gromacs/utility/smalloc.h"
126 #include "corewrap.h"
129 //! First step used in pressure scaling
130 gmx_int64_t deform_init_init_step_tpx;
131 //! Initial box for pressure scaling
132 matrix deform_init_box_tpx;
133 //! MPI variable for use in pressure scaling
134 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
137 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
138 * the number of threads will get lowered.
140 #define MIN_ATOMS_PER_MPI_THREAD 90
141 #define MIN_ATOMS_PER_GPU 900
146 void Mdrunner::reinitializeOnSpawnedThread()
148 // TODO This duplication is formally necessary if any thread might
149 // modify any memory in fnm or the pointers it contains. If the
150 // contents are ever provably const, then we can remove this
151 // allocation (and memory leak).
152 // TODO This should probably become part of a copy constructor for
154 fnm = dup_tfn(nfile, fnm);
156 cr = reinitialize_commrec_for_this_thread(cr);
160 // Only the master rank writes to the log files
165 /*! \brief The callback used for running on spawned threads.
167 * Obtains the pointer to the master mdrunner object from the one
168 * argument permitted to the thread-launch API call, copies it to make
169 * a new runner for this thread, reinitializes necessary data, and
170 * proceeds to the simulation. */
171 static void mdrunner_start_fn(void *arg)
175 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
176 /* copy the arg list to make sure that it's thread-local. This
177 doesn't copy pointed-to items, of course, but those are all
179 gmx::Mdrunner mdrunner = *masterMdrunner;
180 mdrunner.reinitializeOnSpawnedThread();
183 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
187 /*! \brief Start thread-MPI threads.
189 * Called by mdrunner() to start a specific number of threads
190 * (including the main thread) for thread-parallel runs. This in turn
191 * calls mdrunner() for each thread. All options are the same as for
193 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch)
196 /* first check whether we even need to start tMPI */
197 if (numThreadsToLaunch < 2)
202 gmx::Mdrunner spawnedMdrunner = *this;
203 // TODO This duplication is formally necessary if any thread might
204 // modify any memory in fnm or the pointers it contains. If the
205 // contents are ever provably const, then we can remove this
206 // allocation (and memory leak).
207 // TODO This should probably become part of a copy constructor for
209 spawnedMdrunner.fnm = dup_tfn(this->nfile, fnm);
211 /* now spawn new threads that start mdrunner_start_fn(), while
212 the main thread returns, we set thread affinity later */
213 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
214 mdrunner_start_fn, static_cast<void*>(&spawnedMdrunner)) != TMPI_SUCCESS)
216 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
219 return reinitialize_commrec_for_this_thread(cr);
224 #endif /* GMX_THREAD_MPI */
226 /*! \brief Initialize variables for Verlet scheme simulation */
227 static void prepare_verlet_scheme(FILE *fplog,
231 const gmx_mtop_t *mtop,
233 bool makeGpuPairList,
234 const gmx::CpuInfo &cpuinfo)
236 /* For NVE simulations, we will retain the initial list buffer */
237 if (EI_DYNAMICS(ir->eI) &&
238 ir->verletbuf_tol > 0 &&
239 !(EI_MD(ir->eI) && ir->etc == etcNO))
241 /* Update the Verlet buffer size for the current run setup */
242 verletbuf_list_setup_t ls;
245 /* Here we assume SIMD-enabled kernels are being used. But as currently
246 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
247 * and 4x2 gives a larger buffer than 4x4, this is ok.
249 verletbuf_get_list_setup(true, makeGpuPairList, &ls);
251 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &ls, nullptr, &rlist_new);
253 if (rlist_new != ir->rlist)
255 if (fplog != nullptr)
257 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
258 ir->rlist, rlist_new,
259 ls.cluster_size_i, ls.cluster_size_j);
261 ir->rlist = rlist_new;
265 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
267 gmx_fatal(FARGS, "Can not set nstlist without %s",
268 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
271 if (EI_DYNAMICS(ir->eI))
273 /* Set or try nstlist values */
274 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
278 /*! \brief Override the nslist value in inputrec
280 * with value passed on the command line (if any)
282 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
283 gmx_int64_t nsteps_cmdline,
288 /* override with anything else than the default -2 */
289 if (nsteps_cmdline > -2)
291 char sbuf_steps[STEPSTRSIZE];
292 char sbuf_msg[STRLEN];
294 ir->nsteps = nsteps_cmdline;
295 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
297 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
298 gmx_step_str(nsteps_cmdline, sbuf_steps),
299 fabs(nsteps_cmdline*ir->delta_t));
303 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
304 gmx_step_str(nsteps_cmdline, sbuf_steps));
307 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
309 else if (nsteps_cmdline < -2)
311 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
314 /* Do nothing if nsteps_cmdline == -2 */
320 //! Halt the run if there are inconsistences between user choices to run with GPUs and/or hardware detection.
321 static void exitIfCannotForceGpuRun(bool requirePhysicalGpu,
322 EmulateGpuNonbonded emulateGpuNonbonded,
323 bool useVerletScheme,
324 bool compatibleGpusFound)
326 /* Was GPU acceleration either explicitly (-nb gpu) or implicitly
327 * (gpu ID passed) requested? */
328 if (!requirePhysicalGpu)
333 if (GMX_GPU == GMX_GPU_NONE)
335 gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!",
336 gmx::getProgramContext().displayName());
339 if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
341 gmx_fatal(FARGS, "GPU emulation cannot be requested together with GPU acceleration!");
344 if (!useVerletScheme)
346 gmx_fatal(FARGS, "GPU acceleration requested, but can't be used without cutoff-scheme=Verlet");
349 if (!compatibleGpusFound)
351 gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
355 /*! \brief Return whether GPU acceleration is useful with the given settings.
357 * If not, logs a message about falling back to CPU code. */
358 static bool gpuAccelerationIsUseful(const MDLogger &mdlog,
359 const t_inputrec *ir,
362 if (doRerun && ir->opts.ngener > 1)
364 /* Rerun execution time is dominated by I/O and pair search,
365 * so GPUs are not very useful, plus they do not support more
366 * than one energy group. If the user requested GPUs
367 * explicitly, a fatal error is given later. With non-reruns,
368 * we fall back to a single whole-of system energy group
369 * (which runs much faster than a multiple-energy-groups
370 * implementation would), and issue a note in the .log
371 * file. Users can re-run if they want the information. */
372 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");
379 //! \brief Return the correct integrator function.
380 static integrator_t *my_integrator(unsigned int ei)
389 if (!EI_DYNAMICS(ei))
391 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
406 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
410 GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
412 GMX_THROW(APIError("Non existing integrator selected"));
416 //! Initializes the logger for mdrun.
417 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
419 gmx::LoggerBuilder builder;
420 if (fplog != nullptr)
422 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
424 if (cr == nullptr || SIMMASTER(cr))
426 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
427 &gmx::TextOutputFile::standardError());
429 return builder.build();
432 int Mdrunner::mdrunner()
435 gmx_ddbox_t ddbox = {0};
436 int npme_major, npme_minor;
438 gmx_mtop_t *mtop = nullptr;
439 t_mdatoms *mdatoms = 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 bool forceUseCpu = (strncmp(nbpu_opt, "cpu", 3) == 0);
476 if (!hw_opt.gpuIdTaskAssignment.empty() && forceUseCpu)
478 gmx_fatal(FARGS, "GPU IDs were specified, and short-ranged interactions were assigned to the CPU. Make no more than one of these choices.");
480 bool forceUsePhysicalGpu = (strncmp(nbpu_opt, "gpu", 3) == 0) || !hw_opt.gpuIdTaskAssignment.empty();
481 bool tryUsePhysicalGpu = (strncmp(nbpu_opt, "auto", 4) == 0) && hw_opt.gpuIdTaskAssignment.empty() && (emulateGpuNonbonded == EmulateGpuNonbonded::No);
482 GMX_RELEASE_ASSERT(!(forceUsePhysicalGpu && tryUsePhysicalGpu), "Must either force use of "
483 "GPUs for short-ranged interactions, or try to use them, not both.");
484 const PmeRunMode pmeRunMode = PmeRunMode::CPU;
485 //TODO this is a placeholder as PME on GPU is not permitted yet
486 //TODO should there exist a PmeRunMode::None value for consistency?
488 // Here we assume that SIMMASTER(cr) does not change even after the
489 // threads are started.
490 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
491 gmx::MDLogger mdlog(logOwner.logger());
493 hwinfo = gmx_detect_hardware(mdlog, cr);
495 gmx_print_detected_hardware(fplog, cr, mdlog, hwinfo);
497 if (fplog != nullptr)
499 /* Print references after all software/hardware printing */
500 please_cite(fplog, "Abraham2015");
501 please_cite(fplog, "Pall2015");
502 please_cite(fplog, "Pronk2013");
503 please_cite(fplog, "Hess2008b");
504 please_cite(fplog, "Spoel2005a");
505 please_cite(fplog, "Lindahl2001a");
506 please_cite(fplog, "Berendsen95a");
509 std::unique_ptr<t_state> globalState;
513 /* Only the master rank has the global state */
514 globalState = std::unique_ptr<t_state>(new t_state);
516 /* Read (nearly) all data required for the simulation */
517 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), mtop);
519 exitIfCannotForceGpuRun(forceUsePhysicalGpu,
521 inputrec->cutoff_scheme == ecutsVERLET,
522 compatibleGpusFound(hwinfo->gpu_info));
524 if (inputrec->cutoff_scheme == ecutsVERLET)
526 /* TODO This logic could run later, e.g. before -npme -1
527 is handled. If inputrec has already been communicated,
528 then the resulting tryUsePhysicalGpu does not need to
530 if ((tryUsePhysicalGpu || forceUsePhysicalGpu) &&
531 !gpuAccelerationIsUseful(mdlog, inputrec, doRerun))
533 /* Fallback message printed by nbnxn_acceleration_supported */
534 if (forceUsePhysicalGpu)
536 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
538 tryUsePhysicalGpu = false;
543 if (nstlist_cmdline > 0)
545 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
548 if (compatibleGpusFound(hwinfo->gpu_info))
550 GMX_LOG(mdlog.warning).asParagraph().appendText(
551 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
552 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
554 tryUsePhysicalGpu = false;
557 bool nonbondedOnGpu = (tryUsePhysicalGpu || forceUsePhysicalGpu) && compatibleGpusFound(hwinfo->gpu_info);
559 /* Check and update the hardware options for internal consistency */
560 check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
562 /* Early check for externally set process affinity. */
563 gmx_check_thread_affinity_set(mdlog, cr,
564 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
569 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
571 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
574 /* Since the master knows the cut-off scheme, update hw_opt for this.
575 * This is done later for normal MPI and also once more with tMPI
576 * for all tMPI ranks.
578 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
580 /* Determine how many thread-MPI ranks to start.
582 * TODO Over-writing the user-supplied value here does
583 * prevent any possible subsequent checks from working
585 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
587 domdecOptions.numPmeRanks,
593 // Now start the threads for thread MPI.
594 cr = spawnThreads(hw_opt.nthreads_tmpi);
595 /* The main thread continues here with a new cr. We don't deallocate
596 the old cr because other threads may still be reading it. */
597 // TODO Both master and spawned threads call dup_tfn and
598 // reinitialize_commrec_for_this_thread. Find a way to express
602 /* END OF CAUTION: cr is now reliable */
606 /* now broadcast everything to the non-master nodes/threads: */
607 init_parallel(cr, inputrec, mtop);
609 gmx_bcast_sim(sizeof(nonbondedOnGpu), &nonbondedOnGpu, cr);
611 // TODO: Error handling
612 mdModules.assignOptionsToModules(*inputrec->params, nullptr);
614 if (fplog != nullptr)
616 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
617 fprintf(fplog, "\n");
622 /* now make sure the state is initialized and propagated */
623 set_state_entries(globalState.get(), inputrec);
626 /* NM and TPI parallelize over force/energy calculations, not atoms,
627 * so we need to initialize and broadcast the global state.
629 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
633 globalState = std::unique_ptr<t_state>(new t_state);
635 broadcastStateWithoutDynamics(cr, globalState.get());
638 /* A parallel command line option consistency check that we can
639 only do after any threads have started. */
640 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
641 domdecOptions.numCells[YY] > 1 ||
642 domdecOptions.numCells[ZZ] > 1 ||
643 domdecOptions.numPmeRanks > 0))
646 "The -dd or -npme option request a parallel simulation, "
648 "but %s was compiled without threads or MPI enabled"
651 "but the number of MPI-threads (option -ntmpi) is not set or is 1"
653 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
656 , output_env_get_program_display_name(oenv)
661 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
663 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
666 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
668 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
671 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
673 if (domdecOptions.numPmeRanks > 0)
675 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
676 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
679 domdecOptions.numPmeRanks = 0;
682 if (nonbondedOnGpu && domdecOptions.numPmeRanks < 0)
684 /* With GPUs we don't automatically use PME-only ranks. PME ranks can
685 * improve performance with many threads per GPU, since our OpenMP
686 * scaling is bad, but it's difficult to automate the setup.
688 domdecOptions.numPmeRanks = 0;
694 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
698 /* NMR restraints must be initialized before load_checkpoint,
699 * since with time averaging the history is added to t_state.
700 * For proper consistency check we therefore need to extend
702 * So the PME-only nodes (if present) will also initialize
703 * the distance restraints.
707 /* This needs to be called before read_checkpoint to extend the state */
708 init_disres(fplog, mtop, inputrec, cr, fcd, globalState.get(), replExParams.exchangeInterval > 0);
710 init_orires(fplog, mtop, inputrec, cr, globalState.get(), &(fcd->orires));
712 if (inputrecDeform(inputrec))
714 /* Store the deform reference box before reading the checkpoint */
717 copy_mat(globalState->box, box);
721 gmx_bcast(sizeof(box), box, cr);
723 /* Because we do not have the update struct available yet
724 * in which the reference values should be stored,
725 * we store them temporarily in static variables.
726 * This should be thread safe, since they are only written once
727 * and with identical values.
729 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
730 deform_init_init_step_tpx = inputrec->init_step;
731 copy_mat(box, deform_init_box_tpx);
732 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
735 ObservablesHistory observablesHistory = {};
737 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
739 if (continuationOptions.startedFromCheckpoint)
741 /* Check if checkpoint file exists before doing continuation.
742 * This way we can use identical input options for the first and subsequent runs...
746 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
747 cr, domdecOptions.numCells,
748 inputrec, globalState.get(),
749 &bReadEkin, &observablesHistory,
750 continuationOptions.appendFiles,
751 continuationOptions.appendFilesOptionSet,
752 mdrunOptions.reproducible);
756 continuationOptions.haveReadEkin = true;
760 if (SIMMASTER(cr) && continuationOptions.appendFiles)
762 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
763 continuationOptions.appendFiles, &fplog);
764 logOwner = buildLogger(fplog, nullptr);
765 mdlog = logOwner.logger();
768 /* override nsteps with value set on the commamdline */
769 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
773 copy_mat(globalState->box, box);
778 gmx_bcast(sizeof(box), box, cr);
781 /* Update rlist and nstlist. */
782 if (inputrec->cutoff_scheme == ecutsVERLET)
784 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, mtop, box,
785 nonbondedOnGpu || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
788 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
789 inputrec->eI == eiNM))
791 const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
793 cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
796 &ddbox, &npme_major, &npme_minor);
797 // Note that local state still does not exist yet.
801 /* PME, if used, is done on all nodes with 1D decomposition */
803 cr->duty = (DUTY_PP | DUTY_PME);
807 if (inputrec->ePBC == epbcSCREW)
810 "pbc=%s is only implemented with domain decomposition",
811 epbc_names[inputrec->ePBC]);
817 /* After possible communicator splitting in make_dd_communicators.
818 * we can set up the intra/inter node communication.
820 gmx_setup_nodecomm(fplog, cr);
823 /* Initialize per-physical-node MPI process/thread ID and counters. */
824 gmx_init_intranode_counters(cr);
828 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
829 "This is simulation %d out of %d running as a composite GROMACS\n"
830 "multi-simulation job. Setup for this simulation:\n",
831 cr->ms->sim, cr->ms->nsim);
833 GMX_LOG(mdlog.warning).appendTextFormatted(
837 cr->nnodes == 1 ? "thread" : "threads"
839 cr->nnodes == 1 ? "process" : "processes"
845 /* Check and update hw_opt for the cut-off scheme */
846 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
848 /* Check and update hw_opt for the number of MPI ranks */
849 check_and_update_hw_opt_3(&hw_opt);
851 gmx_omp_nthreads_init(mdlog, cr,
852 hwinfo->nthreads_hw_avail,
854 hw_opt.nthreads_omp_pme,
855 !thisRankHasDuty(cr, DUTY_PP),
856 inputrec->cutoff_scheme == ecutsVERLET);
859 if (EI_TPI(inputrec->eI) &&
860 inputrec->cutoff_scheme == ecutsVERLET)
862 gmx_feenableexcept();
866 // Contains the ID of the GPU used by each PP rank on this node,
867 // indexed by that rank. Empty if no GPUs are selected for use on
869 std::vector<int> gpuTaskAssignment;
872 /* Currently the DD code assigns duty to ranks that can
873 * include PP work that currently can be executed on a single
874 * GPU, if present and compatible. This has to be coordinated
875 * across PP ranks on a node, with possible multiple devices
876 * or sharing devices on a node, either from the user
877 * selection, or automatically. */
878 bool rankCanUseGpu = thisRankHasDuty(cr, DUTY_PP);
879 gpuTaskAssignment = mapPpRanksToGpus(rankCanUseGpu, cr, hwinfo->gpu_info, hw_opt);
882 reportGpuUsage(mdlog, hwinfo->gpu_info, !hw_opt.gpuIdTaskAssignment.empty(),
883 gpuTaskAssignment, cr->nrank_pp_intranode, cr->nnodes > 1);
885 if (!gpuTaskAssignment.empty())
887 GMX_RELEASE_ASSERT(cr->nrank_pp_intranode == static_cast<int>(gpuTaskAssignment.size()),
888 "The number of PP ranks on each node must equal the number of GPU tasks used on each node");
891 /* Prevent other ranks from continuing after an issue was found
892 * and reported as a fatal error.
894 * TODO This function implements a barrier so that MPI runtimes
895 * can organize an orderly shutdown if one of the ranks has had to
896 * issue a fatal error in various code already run. When we have
897 * MPI-aware error handling and reporting, this should be
902 MPI_Barrier(cr->mpi_comm_mysim);
906 /* Now that we know the setup is consistent, check for efficiency */
907 check_resource_division_efficiency(hwinfo, hw_opt.nthreads_tot, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
910 gmx_device_info_t *shortRangedDeviceInfo = nullptr;
911 int shortRangedDeviceId = -1;
912 if (thisRankHasDuty(cr, DUTY_PP))
914 if (!gpuTaskAssignment.empty())
916 shortRangedDeviceId = gpuTaskAssignment[cr->rank_pp_intranode];
917 shortRangedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, shortRangedDeviceId);
921 if (DOMAINDECOMP(cr))
923 /* When we share GPUs over ranks, we need to know this for the DLB */
924 dd_setup_dlb_resource_sharing(cr, shortRangedDeviceId);
927 /* getting number of PP/PME threads
928 PME: env variable should be read only on one node to make sure it is
929 identical everywhere;
931 nthreads_pme = gmx_omp_nthreads_get(emntPME);
933 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
937 /* Master synchronizes its value of reset_counters with all nodes
938 * including PME only nodes */
939 reset_counters = wcycle_get_reset_counters(wcycle);
940 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
941 wcycle_set_reset_counters(wcycle, reset_counters);
944 // Membrane embedding must be initialized before we call init_forcerec()
949 fprintf(stderr, "Initializing membed");
951 /* Note that membed cannot work in parallel because mtop is
952 * changed here. Fix this if we ever want to make it run with
954 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
958 if (thisRankHasDuty(cr, DUTY_PP))
960 /* Initiate forcerecord */
962 fr->forceProviders = mdModules.initForceProviders();
963 init_forcerec(fplog, mdlog, fr, fcd,
964 inputrec, mtop, cr, box,
965 opt2fn("-table", nfile, fnm),
966 opt2fn("-tablep", nfile, fnm),
967 getFilenm("-tableb", nfile, fnm),
968 shortRangedDeviceInfo,
972 /* Initialize QM-MM */
975 init_QMMMrec(cr, mtop, inputrec, fr);
978 /* Initialize the mdatoms structure.
979 * mdatoms is not filled with atom data,
980 * as this can not be done now with domain decomposition.
982 mdatoms = init_mdatoms(fplog, *mtop, *inputrec);
984 /* Initialize the virtual site communication */
985 vsite = initVsite(*mtop, cr);
987 calc_shifts(box, fr->shift_vec);
989 /* With periodic molecules the charge groups should be whole at start up
990 * and the virtual sites should not be far from their proper positions.
992 if (!inputrec->bContinuation && MASTER(cr) &&
993 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
995 /* Make molecules whole at start of run */
996 if (fr->ePBC != epbcNONE)
998 rvec *xGlobal = as_rvec_array(globalState->x.data());
999 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, xGlobal);
1003 /* Correct initial vsite positions are required
1004 * for the initial distribution in the domain decomposition
1005 * and for the initial shell prediction.
1007 constructVsitesGlobal(*mtop, globalState->x);
1011 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1013 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1014 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1019 /* This is a PME only node */
1021 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1023 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1024 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1027 gmx_pme_t *sepPmeData = nullptr;
1028 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1029 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1030 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1032 if (hw_opt.thread_affinity != threadaffOFF)
1034 /* Before setting affinity, check whether the affinity has changed
1035 * - which indicates that probably the OpenMP library has changed it
1036 * since we first checked).
1038 gmx_check_thread_affinity_set(mdlog, cr,
1039 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1042 /* threads on this MPI process or TMPI thread */
1043 if (thisRankHasDuty(cr, DUTY_PP))
1045 nthread_local = gmx_omp_nthreads_get(emntNonbonded);
1049 nthread_local = gmx_omp_nthreads_get(emntPME);
1052 /* Set the CPU affinity */
1053 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1054 nthread_local, nullptr);
1057 /* Initiate PME if necessary,
1058 * either on all nodes or on dedicated PME nodes only. */
1059 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1063 nChargePerturbed = mdatoms->nChargePerturbed;
1064 if (EVDW_PME(inputrec->vdwtype))
1066 nTypePerturbed = mdatoms->nTypePerturbed;
1069 if (cr->npmenodes > 0)
1071 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1072 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1073 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1076 if (thisRankHasDuty(cr, DUTY_PME))
1080 gmx_device_info_t *pmeGpuInfo = nullptr;
1081 pmedata = gmx_pme_init(cr, npme_major, npme_minor, inputrec,
1082 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1083 mdrunOptions.reproducible,
1084 ewaldcoeff_q, ewaldcoeff_lj,
1086 pmeRunMode, nullptr, pmeGpuInfo, mdlog);
1088 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1093 if (EI_DYNAMICS(inputrec->eI))
1095 /* Turn on signal handling on all nodes */
1097 * (A user signal from the PME nodes (if any)
1098 * is communicated to the PP nodes.
1100 signal_handler_install();
1103 if (thisRankHasDuty(cr, DUTY_PP))
1105 /* Assumes uniform use of the number of OpenMP threads */
1106 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1108 if (inputrec->bPull)
1110 /* Initialize pull code */
1111 inputrec->pull_work =
1112 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1113 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1114 EI_DYNAMICS(inputrec->eI) && MASTER(cr),
1115 continuationOptions);
1120 /* Initialize enforced rotation code */
1121 init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), mtop, oenv, mdrunOptions);
1124 /* Let init_constraints know whether we have essential dynamics constraints.
1125 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1127 bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1129 constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1131 if (DOMAINDECOMP(cr))
1133 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1134 /* This call is not included in init_domain_decomposition mainly
1135 * because fr->cginfo_mb is set later.
1137 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1138 domdecOptions.checkBondedInteractions,
1142 /* Now do whatever the user wants us to do (how flexible...) */
1143 my_integrator(inputrec->eI) (fplog, cr, mdlog, nfile, fnm,
1147 mdModules.outputProvider(),
1151 &observablesHistory,
1152 mdatoms, nrnb, wcycle, fr,
1155 walltime_accounting);
1159 finish_rot(inputrec->rot);
1162 if (inputrec->bPull)
1164 finish_pull(inputrec->pull_work);
1170 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1172 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1173 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1176 wallcycle_stop(wcycle, ewcRUN);
1178 /* Finish up, write some stuff
1179 * if rerunMD, don't write last frame again
1181 finish_run(fplog, mdlog, cr,
1182 inputrec, nrnb, wcycle, walltime_accounting,
1183 fr ? fr->nbv : nullptr,
1185 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1190 gmx_pme_destroy(pmedata);
1194 /* Free GPU memory and context */
1195 free_gpu_resources(fr, cr, shortRangedDeviceInfo);
1199 free_membed(membed);
1202 gmx_hardware_info_free();
1204 /* Does what it says */
1205 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1206 walltime_accounting_destroy(walltime_accounting);
1208 /* Close logfile already here if we were appending to it */
1209 if (MASTER(cr) && continuationOptions.appendFiles)
1211 gmx_log_close(fplog);
1214 rc = (int)gmx_get_stop_condition();
1217 /* we need to join all threads. The sub-threads join when they
1218 exit this function, but the master thread needs to be told to
1220 if (PAR(cr) && MASTER(cr))