2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/localatomsetmanager.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/ewald_utils.h"
65 #include "gromacs/ewald/pme.h"
66 #include "gromacs/ewald/pme_gpu_program.h"
67 #include "gromacs/fileio/checkpoint.h"
68 #include "gromacs/fileio/gmxfio.h"
69 #include "gromacs/fileio/oenv.h"
70 #include "gromacs/fileio/tpxio.h"
71 #include "gromacs/gmxlib/network.h"
72 #include "gromacs/gmxlib/nrnb.h"
73 #include "gromacs/gpu_utils/clfftinitializer.h"
74 #include "gromacs/gpu_utils/gpu_utils.h"
75 #include "gromacs/hardware/cpuinfo.h"
76 #include "gromacs/hardware/detecthardware.h"
77 #include "gromacs/hardware/printhardware.h"
78 #include "gromacs/listed_forces/disre.h"
79 #include "gromacs/listed_forces/gpubonded.h"
80 #include "gromacs/listed_forces/orires.h"
81 #include "gromacs/math/functions.h"
82 #include "gromacs/math/utilities.h"
83 #include "gromacs/math/vec.h"
84 #include "gromacs/mdlib/boxdeformation.h"
85 #include "gromacs/mdlib/broadcaststructs.h"
86 #include "gromacs/mdlib/calc_verletbuf.h"
87 #include "gromacs/mdlib/force.h"
88 #include "gromacs/mdlib/forcerec.h"
89 #include "gromacs/mdlib/gmx_omp_nthreads.h"
90 #include "gromacs/mdlib/makeconstraints.h"
91 #include "gromacs/mdlib/md_support.h"
92 #include "gromacs/mdlib/mdatoms.h"
93 #include "gromacs/mdlib/membed.h"
94 #include "gromacs/mdlib/ppforceworkload.h"
95 #include "gromacs/mdlib/qmmm.h"
96 #include "gromacs/mdlib/sighandler.h"
97 #include "gromacs/mdlib/stophandler.h"
98 #include "gromacs/mdrun/logging.h"
99 #include "gromacs/mdrun/multisim.h"
100 #include "gromacs/mdrun/simulationcontext.h"
101 #include "gromacs/mdrunutility/mdmodules.h"
102 #include "gromacs/mdrunutility/printtime.h"
103 #include "gromacs/mdrunutility/threadaffinity.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/enerdata.h"
106 #include "gromacs/mdtypes/fcdata.h"
107 #include "gromacs/mdtypes/inputrec.h"
108 #include "gromacs/mdtypes/md_enums.h"
109 #include "gromacs/mdtypes/mdrunoptions.h"
110 #include "gromacs/mdtypes/observableshistory.h"
111 #include "gromacs/mdtypes/state.h"
112 #include "gromacs/nbnxm/gpu_data_mgmt.h"
113 #include "gromacs/nbnxm/nbnxm.h"
114 #include "gromacs/nbnxm/pairlist_tuning.h"
115 #include "gromacs/pbcutil/pbc.h"
116 #include "gromacs/pulling/output.h"
117 #include "gromacs/pulling/pull.h"
118 #include "gromacs/pulling/pull_rotation.h"
119 #include "gromacs/restraint/manager.h"
120 #include "gromacs/restraint/restraintmdmodule.h"
121 #include "gromacs/restraint/restraintpotential.h"
122 #include "gromacs/swap/swapcoords.h"
123 #include "gromacs/taskassignment/decidegpuusage.h"
124 #include "gromacs/taskassignment/resourcedivision.h"
125 #include "gromacs/taskassignment/taskassignment.h"
126 #include "gromacs/taskassignment/usergpuids.h"
127 #include "gromacs/timing/gpu_timing.h"
128 #include "gromacs/timing/wallcycle.h"
129 #include "gromacs/timing/wallcyclereporting.h"
130 #include "gromacs/topology/mtop_util.h"
131 #include "gromacs/trajectory/trajectoryframe.h"
132 #include "gromacs/utility/basenetwork.h"
133 #include "gromacs/utility/cstringutil.h"
134 #include "gromacs/utility/exceptions.h"
135 #include "gromacs/utility/fatalerror.h"
136 #include "gromacs/utility/filestream.h"
137 #include "gromacs/utility/gmxassert.h"
138 #include "gromacs/utility/gmxmpi.h"
139 #include "gromacs/utility/logger.h"
140 #include "gromacs/utility/loggerbuilder.h"
141 #include "gromacs/utility/physicalnodecommunicator.h"
142 #include "gromacs/utility/pleasecite.h"
143 #include "gromacs/utility/programcontext.h"
144 #include "gromacs/utility/smalloc.h"
145 #include "gromacs/utility/stringutil.h"
147 #include "integrator.h"
148 #include "replicaexchange.h"
151 #include "corewrap.h"
157 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
159 * Used to ensure that the master thread does not modify mdrunner during copy
160 * on the spawned threads. */
161 static void threadMpiMdrunnerAccessBarrier()
164 MPI_Barrier(MPI_COMM_WORLD);
168 Mdrunner Mdrunner::cloneOnSpawnedThread() const
170 auto newRunner = Mdrunner(std::make_unique<MDModules>());
172 // All runners in the same process share a restraint manager resource because it is
173 // part of the interface to the client code, which is associated only with the
174 // original thread. Handles to the same resources can be obtained by copy.
176 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
179 // Copy original cr pointer before master thread can pass the thread barrier
180 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
182 // Copy members of master runner.
183 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
184 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
185 newRunner.hw_opt = hw_opt;
186 newRunner.filenames = filenames;
188 newRunner.oenv = oenv;
189 newRunner.mdrunOptions = mdrunOptions;
190 newRunner.domdecOptions = domdecOptions;
191 newRunner.nbpu_opt = nbpu_opt;
192 newRunner.pme_opt = pme_opt;
193 newRunner.pme_fft_opt = pme_fft_opt;
194 newRunner.bonded_opt = bonded_opt;
195 newRunner.nstlist_cmdline = nstlist_cmdline;
196 newRunner.replExParams = replExParams;
197 newRunner.pforce = pforce;
199 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
201 threadMpiMdrunnerAccessBarrier();
203 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
208 /*! \brief The callback used for running on spawned threads.
210 * Obtains the pointer to the master mdrunner object from the one
211 * argument permitted to the thread-launch API call, copies it to make
212 * a new runner for this thread, reinitializes necessary data, and
213 * proceeds to the simulation. */
214 static void mdrunner_start_fn(const void *arg)
218 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
219 /* copy the arg list to make sure that it's thread-local. This
220 doesn't copy pointed-to items, of course; fnm, cr and fplog
221 are reset in the call below, all others should be const. */
222 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
225 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
229 /*! \brief Start thread-MPI threads.
231 * Called by mdrunner() to start a specific number of threads
232 * (including the main thread) for thread-parallel runs. This in turn
233 * calls mdrunner() for each thread. All options are the same as for
235 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
238 /* first check whether we even need to start tMPI */
239 if (numThreadsToLaunch < 2)
245 /* now spawn new threads that start mdrunner_start_fn(), while
246 the main thread returns, we set thread affinity later */
247 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
248 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
250 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
253 threadMpiMdrunnerAccessBarrier();
255 GMX_UNUSED_VALUE(mdrunner_start_fn);
258 return reinitialize_commrec_for_this_thread(cr);
263 /*! \brief Initialize variables for Verlet scheme simulation */
264 static void prepare_verlet_scheme(FILE *fplog,
268 const gmx_mtop_t *mtop,
270 bool makeGpuPairList,
271 const gmx::CpuInfo &cpuinfo)
273 /* For NVE simulations, we will retain the initial list buffer */
274 if (EI_DYNAMICS(ir->eI) &&
275 ir->verletbuf_tol > 0 &&
276 !(EI_MD(ir->eI) && ir->etc == etcNO))
278 /* Update the Verlet buffer size for the current run setup */
280 /* Here we assume SIMD-enabled kernels are being used. But as currently
281 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
282 * and 4x2 gives a larger buffer than 4x4, this is ok.
284 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
285 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
287 const real rlist_new =
288 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
290 if (rlist_new != ir->rlist)
292 if (fplog != nullptr)
294 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
295 ir->rlist, rlist_new,
296 listSetup.cluster_size_i, listSetup.cluster_size_j);
298 ir->rlist = rlist_new;
302 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
304 gmx_fatal(FARGS, "Can not set nstlist without %s",
305 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
308 if (EI_DYNAMICS(ir->eI))
310 /* Set or try nstlist values */
311 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
315 /*! \brief Override the nslist value in inputrec
317 * with value passed on the command line (if any)
319 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
320 int64_t nsteps_cmdline,
325 /* override with anything else than the default -2 */
326 if (nsteps_cmdline > -2)
328 char sbuf_steps[STEPSTRSIZE];
329 char sbuf_msg[STRLEN];
331 ir->nsteps = nsteps_cmdline;
332 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
334 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
335 gmx_step_str(nsteps_cmdline, sbuf_steps),
336 fabs(nsteps_cmdline*ir->delta_t));
340 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
341 gmx_step_str(nsteps_cmdline, sbuf_steps));
344 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
346 else if (nsteps_cmdline < -2)
348 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
351 /* Do nothing if nsteps_cmdline == -2 */
357 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
359 * If not, and if a warning may be issued, logs a warning about
360 * falling back to CPU code. With thread-MPI, only the first
361 * call to this function should have \c issueWarning true. */
362 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
363 const t_inputrec *ir,
366 if (ir->opts.ngener - ir->nwall > 1)
368 /* The GPU code does not support more than one energy group.
369 * If the user requested GPUs explicitly, a fatal error is given later.
373 GMX_LOG(mdlog.warning).asParagraph()
374 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
375 "For better performance, run on the GPU without energy groups and then do "
376 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
383 //! Initializes the logger for mdrun.
384 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
386 gmx::LoggerBuilder builder;
387 if (fplog != nullptr)
389 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
391 if (cr == nullptr || SIMMASTER(cr))
393 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
394 &gmx::TextOutputFile::standardError());
396 return builder.build();
399 //! Make a TaskTarget from an mdrun argument string.
400 static TaskTarget findTaskTarget(const char *optionString)
402 TaskTarget returnValue = TaskTarget::Auto;
404 if (strncmp(optionString, "auto", 3) == 0)
406 returnValue = TaskTarget::Auto;
408 else if (strncmp(optionString, "cpu", 3) == 0)
410 returnValue = TaskTarget::Cpu;
412 else if (strncmp(optionString, "gpu", 3) == 0)
414 returnValue = TaskTarget::Gpu;
418 GMX_ASSERT(false, "Option string should have been checked for sanity already");
424 //! Finish run, aggregate data to print performance info.
425 static void finish_run(FILE *fplog,
426 const gmx::MDLogger &mdlog,
428 const t_inputrec *inputrec,
429 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
430 gmx_walltime_accounting_t walltime_accounting,
431 nonbonded_verlet_t *nbv,
432 const gmx_pme_t *pme,
435 t_nrnb *nrnb_tot = nullptr;
437 double nbfs = 0, mflop = 0;
439 elapsed_time_over_all_ranks,
440 elapsed_time_over_all_threads,
441 elapsed_time_over_all_threads_over_all_ranks;
442 /* Control whether it is valid to print a report. Only the
443 simulation master may print, but it should not do so if the run
444 terminated e.g. before a scheduled reset step. This is
445 complicated by the fact that PME ranks are unaware of the
446 reason why they were sent a pmerecvqxFINISH. To avoid
447 communication deadlocks, we always do the communication for the
448 report, even if we've decided not to write the report, because
449 how long it takes to finish the run is not important when we've
450 decided not to report on the simulation performance.
452 Further, we only report performance for dynamical integrators,
453 because those are the only ones for which we plan to
454 consider doing any optimizations. */
455 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
457 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
459 GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
467 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
476 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
477 elapsed_time_over_all_threads = walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
481 /* reduce elapsed_time over all MPI ranks in the current simulation */
482 MPI_Allreduce(&elapsed_time,
483 &elapsed_time_over_all_ranks,
484 1, MPI_DOUBLE, MPI_SUM,
486 elapsed_time_over_all_ranks /= cr->nnodes;
487 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
488 * current simulation. */
489 MPI_Allreduce(&elapsed_time_over_all_threads,
490 &elapsed_time_over_all_threads_over_all_ranks,
491 1, MPI_DOUBLE, MPI_SUM,
497 elapsed_time_over_all_ranks = elapsed_time;
498 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
503 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
510 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
512 print_dd_statistics(cr, inputrec, fplog);
515 /* TODO Move the responsibility for any scaling by thread counts
516 * to the code that handled the thread region, so that there's a
517 * mechanism to keep cycle counting working during the transition
518 * to task parallelism. */
519 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
520 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
521 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
522 auto cycle_sum(wallcycle_sum(cr, wcycle));
526 auto nbnxn_gpu_timings = (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
527 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
529 if (pme_gpu_task_enabled(pme))
531 pme_gpu_get_timings(pme, &pme_gpu_timings);
533 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
534 elapsed_time_over_all_ranks,
539 if (EI_DYNAMICS(inputrec->eI))
541 delta_t = inputrec->delta_t;
546 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
547 elapsed_time_over_all_ranks,
548 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
549 delta_t, nbfs, mflop);
553 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
554 elapsed_time_over_all_ranks,
555 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
556 delta_t, nbfs, mflop);
561 int Mdrunner::mdrunner()
565 t_forcerec *fr = nullptr;
566 t_fcdata *fcd = nullptr;
567 real ewaldcoeff_q = 0;
568 real ewaldcoeff_lj = 0;
569 int nChargePerturbed = -1, nTypePerturbed = 0;
570 gmx_wallcycle_t wcycle;
571 gmx_walltime_accounting_t walltime_accounting = nullptr;
572 gmx_membed_t * membed = nullptr;
573 gmx_hw_info_t *hwinfo = nullptr;
575 /* CAUTION: threads may be started later on in this function, so
576 cr doesn't reflect the final parallel state right now */
577 t_inputrec inputrecInstance;
578 t_inputrec *inputrec = &inputrecInstance;
581 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
582 bool doRerun = mdrunOptions.rerun;
584 // Handle task-assignment related user options.
585 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
586 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
587 std::vector<int> gpuIdsAvailable;
590 gpuIdsAvailable = parseUserGpuIdString(hw_opt.gpuIdsAvailable);
592 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
594 std::vector<int> userGpuTaskAssignment;
597 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
599 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
600 auto nonbondedTarget = findTaskTarget(nbpu_opt);
601 auto pmeTarget = findTaskTarget(pme_opt);
602 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
603 auto bondedTarget = findTaskTarget(bonded_opt);
604 PmeRunMode pmeRunMode = PmeRunMode::None;
606 // Here we assume that SIMMASTER(cr) does not change even after the
607 // threads are started.
609 FILE *fplog = nullptr;
610 // If we are appending, we don't write log output because we need
611 // to check that the old log file matches what the checkpoint file
612 // expects. Otherwise, we should start to write log output now if
613 // there is a file ready for it.
614 if (logFileHandle != nullptr && !mdrunOptions.continuationOptions.appendFiles)
616 fplog = gmx_fio_getfp(logFileHandle);
618 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
619 gmx::MDLogger mdlog(logOwner.logger());
621 // TODO The thread-MPI master rank makes a working
622 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
623 // after the threads have been launched. This works because no use
624 // is made of that communicator until after the execution paths
625 // have rejoined. But it is likely that we can improve the way
626 // this is expressed, e.g. by expressly running detection only the
627 // master rank for thread-MPI, rather than relying on the mutex
628 // and reference count.
629 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
630 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
632 gmx_print_detected_hardware(fplog, cr, ms, mdlog, hwinfo);
634 std::vector<int> gpuIdsToUse;
635 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
636 if (gpuIdsAvailable.empty())
638 gpuIdsToUse = compatibleGpus;
642 for (const auto &availableGpuId : gpuIdsAvailable)
644 bool availableGpuIsCompatible = false;
645 for (const auto &compatibleGpuId : compatibleGpus)
647 if (availableGpuId == compatibleGpuId)
649 availableGpuIsCompatible = true;
653 if (!availableGpuIsCompatible)
655 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);
657 gpuIdsToUse.push_back(availableGpuId);
661 if (fplog != nullptr)
663 /* Print references after all software/hardware printing */
664 please_cite(fplog, "Abraham2015");
665 please_cite(fplog, "Pall2015");
666 please_cite(fplog, "Pronk2013");
667 please_cite(fplog, "Hess2008b");
668 please_cite(fplog, "Spoel2005a");
669 please_cite(fplog, "Lindahl2001a");
670 please_cite(fplog, "Berendsen95a");
671 writeSourceDoi(fplog);
674 std::unique_ptr<t_state> globalState;
678 /* Only the master rank has the global state */
679 globalState = std::make_unique<t_state>();
681 /* Read (nearly) all data required for the simulation */
682 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
684 /* In rerun, set velocities to zero if present */
685 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
687 // rerun does not use velocities
688 GMX_LOG(mdlog.info).asParagraph().appendText(
689 "Rerun trajectory contains velocities. Rerun does only evaluate "
690 "potential energy and forces. The velocities will be ignored.");
691 for (int i = 0; i < globalState->natoms; i++)
693 clear_rvec(globalState->v[i]);
695 globalState->flags &= ~(1 << estV);
698 if (inputrec->cutoff_scheme != ecutsVERLET)
700 if (nstlist_cmdline > 0)
702 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
705 if (!compatibleGpus.empty())
707 GMX_LOG(mdlog.warning).asParagraph().appendText(
708 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
709 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
714 /* Check and update the hardware options for internal consistency */
715 check_and_update_hw_opt_1(mdlog, &hw_opt, cr, domdecOptions.numPmeRanks);
717 /* Early check for externally set process affinity. */
718 gmx_check_thread_affinity_set(mdlog, cr,
719 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
721 if (GMX_THREAD_MPI && SIMMASTER(cr))
723 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
725 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
728 /* Since the master knows the cut-off scheme, update hw_opt for this.
729 * This is done later for normal MPI and also once more with tMPI
730 * for all tMPI ranks.
732 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
734 bool useGpuForNonbonded = false;
735 bool useGpuForPme = false;
738 // If the user specified the number of ranks, then we must
739 // respect that, but in default mode, we need to allow for
740 // the number of GPUs to choose the number of ranks.
741 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
742 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
743 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
744 canUseGpuForNonbonded,
745 inputrec->cutoff_scheme == ecutsVERLET,
746 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
747 hw_opt.nthreads_tmpi);
748 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
749 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
750 *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
753 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
755 /* Determine how many thread-MPI ranks to start.
757 * TODO Over-writing the user-supplied value here does
758 * prevent any possible subsequent checks from working
760 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
769 // Now start the threads for thread MPI.
770 cr = spawnThreads(hw_opt.nthreads_tmpi);
771 /* The main thread continues here with a new cr. We don't deallocate
772 the old cr because other threads may still be reading it. */
773 // TODO Both master and spawned threads call dup_tfn and
774 // reinitialize_commrec_for_this_thread. Find a way to express
776 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
778 // END OF CAUTION: cr and physicalNodeComm are now reliable
782 /* now broadcast everything to the non-master nodes/threads: */
783 init_parallel(cr, inputrec, &mtop);
786 // Now each rank knows the inputrec that SIMMASTER read and used,
787 // and (if applicable) cr->nnodes has been assigned the number of
788 // thread-MPI ranks that have been chosen. The ranks can now all
789 // run the task-deciding functions and will agree on the result
790 // without needing to communicate.
792 // TODO Should we do the communication in debug mode to support
793 // having an assertion?
795 // Note that these variables describe only their own node.
797 // Note that when bonded interactions run on a GPU they always run
798 // alongside a nonbonded task, so do not influence task assignment
799 // even though they affect the force calculation workload.
800 bool useGpuForNonbonded = false;
801 bool useGpuForPme = false;
802 bool useGpuForBonded = false;
805 // It's possible that there are different numbers of GPUs on
806 // different nodes, which is the user's responsibilty to
807 // handle. If unsuitable, we will notice that during task
809 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
810 bool usingVerletScheme = inputrec->cutoff_scheme == ecutsVERLET;
811 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
812 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
814 canUseGpuForNonbonded,
816 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
818 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
819 *hwinfo, *inputrec, mtop,
820 cr->nnodes, domdecOptions.numPmeRanks,
822 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
824 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
825 bondedTarget, canUseGpuForBonded,
826 EVDW_PME(inputrec->vdwtype),
827 EEL_PME_EWALD(inputrec->coulombtype),
828 domdecOptions.numPmeRanks, gpusWereDetected);
830 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
831 if (pmeRunMode == PmeRunMode::GPU)
833 if (pmeFftTarget == TaskTarget::Cpu)
835 pmeRunMode = PmeRunMode::Mixed;
838 else if (pmeFftTarget == TaskTarget::Gpu)
840 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.");
843 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
846 // TODO: hide restraint implementation details from Mdrunner.
847 // There is nothing unique about restraints at this point as far as the
848 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
849 // factory functions from the SimulationContext on which to call mdModules_->add().
850 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
851 for (auto && restraint : restraintManager_->getRestraints())
853 auto module = RestraintMDModule::create(restraint,
855 mdModules_->add(std::move(module));
858 // TODO: Error handling
859 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
861 if (fplog != nullptr)
863 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
864 fprintf(fplog, "\n");
869 /* now make sure the state is initialized and propagated */
870 set_state_entries(globalState.get(), inputrec);
873 /* NM and TPI parallelize over force/energy calculations, not atoms,
874 * so we need to initialize and broadcast the global state.
876 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
880 globalState = std::make_unique<t_state>();
882 broadcastStateWithoutDynamics(cr, globalState.get());
885 /* A parallel command line option consistency check that we can
886 only do after any threads have started. */
887 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
888 domdecOptions.numCells[YY] > 1 ||
889 domdecOptions.numCells[ZZ] > 1 ||
890 domdecOptions.numPmeRanks > 0))
893 "The -dd or -npme option request a parallel simulation, "
895 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
898 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
900 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
906 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
908 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
911 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
913 if (domdecOptions.numPmeRanks > 0)
915 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
916 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
919 domdecOptions.numPmeRanks = 0;
922 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
924 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
925 * improve performance with many threads per GPU, since our OpenMP
926 * scaling is bad, but it's difficult to automate the setup.
928 domdecOptions.numPmeRanks = 0;
932 if (domdecOptions.numPmeRanks < 0)
934 domdecOptions.numPmeRanks = 0;
935 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
939 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
946 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
950 /* NMR restraints must be initialized before load_checkpoint,
951 * since with time averaging the history is added to t_state.
952 * For proper consistency check we therefore need to extend
954 * So the PME-only nodes (if present) will also initialize
955 * the distance restraints.
959 /* This needs to be called before read_checkpoint to extend the state */
960 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
962 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
964 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
966 ObservablesHistory observablesHistory = {};
968 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
970 if (continuationOptions.startedFromCheckpoint)
972 /* Check if checkpoint file exists before doing continuation.
973 * This way we can use identical input options for the first and subsequent runs...
977 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
979 cr, domdecOptions.numCells,
980 inputrec, globalState.get(),
981 &bReadEkin, &observablesHistory,
982 continuationOptions.appendFiles,
983 continuationOptions.appendFilesOptionSet,
984 mdrunOptions.reproducible);
988 continuationOptions.haveReadEkin = true;
991 if (continuationOptions.appendFiles && logFileHandle)
993 // Now we can start normal logging to the truncated log file.
994 fplog = gmx_fio_getfp(logFileHandle);
995 prepareLogAppending(fplog);
996 logOwner = buildLogger(fplog, cr);
997 mdlog = logOwner.logger();
1001 if (mdrunOptions.numStepsCommandline > -2)
1003 GMX_LOG(mdlog.info).asParagraph().
1004 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
1005 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
1007 /* override nsteps with value set on the commamdline */
1008 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1012 copy_mat(globalState->box, box);
1017 gmx_bcast(sizeof(box), box, cr);
1020 /* Update rlist and nstlist. */
1021 if (inputrec->cutoff_scheme == ecutsVERLET)
1023 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1024 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
1027 LocalAtomSetManager atomSets;
1029 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1030 inputrec->eI == eiNM))
1032 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
1034 box, positionsFromStatePointer(globalState.get()),
1036 // Note that local state still does not exist yet.
1040 /* PME, if used, is done on all nodes with 1D decomposition */
1042 cr->duty = (DUTY_PP | DUTY_PME);
1044 if (inputrec->ePBC == epbcSCREW)
1047 "pbc=screw is only implemented with domain decomposition");
1053 /* After possible communicator splitting in make_dd_communicators.
1054 * we can set up the intra/inter node communication.
1056 gmx_setup_nodecomm(fplog, cr);
1062 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1063 "This is simulation %d out of %d running as a composite GROMACS\n"
1064 "multi-simulation job. Setup for this simulation:\n",
1067 GMX_LOG(mdlog.warning).appendTextFormatted(
1068 "Using %d MPI %s\n",
1071 cr->nnodes == 1 ? "thread" : "threads"
1073 cr->nnodes == 1 ? "process" : "processes"
1079 /* Check and update hw_opt for the cut-off scheme */
1080 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
1082 /* Check and update the number of OpenMP threads requested */
1083 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1086 gmx_omp_nthreads_init(mdlog, cr,
1087 hwinfo->nthreads_hw_avail,
1088 physicalNodeComm.size_,
1089 hw_opt.nthreads_omp,
1090 hw_opt.nthreads_omp_pme,
1091 !thisRankHasDuty(cr, DUTY_PP),
1092 inputrec->cutoff_scheme == ecutsVERLET);
1094 // Enable FP exception detection for the Verlet scheme, but not in
1095 // Release mode and not for compilers with known buggy FP
1096 // exception support (clang with any optimization) or suspected
1097 // buggy FP exception support (gcc 7.* with optimization).
1098 #if !defined NDEBUG && \
1099 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1100 && defined __OPTIMIZE__)
1101 const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
1103 const bool bEnableFPE = false;
1105 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1108 gmx_feenableexcept();
1111 // Build a data structure that expresses which kinds of non-bonded
1112 // task are handled by this rank.
1114 // TODO Later, this might become a loop over all registered modules
1115 // relevant to the mdp inputs, to find those that have such tasks.
1117 // TODO This could move before init_domain_decomposition() as part
1118 // of refactoring that separates the responsibility for duty
1119 // assignment from setup for communication between tasks, and
1120 // setup for tasks handled with a domain (ie including short-ranged
1121 // tasks, bonded tasks, etc.).
1123 // Note that in general useGpuForNonbonded, etc. can have a value
1124 // that is inconsistent with the presence of actual GPUs on any
1125 // rank, and that is not known to be a problem until the
1126 // duty of the ranks on a node become known.
1128 // TODO Later we might need the concept of computeTasksOnThisRank,
1129 // from which we construct gpuTasksOnThisRank.
1131 // Currently the DD code assigns duty to ranks that can
1132 // include PP work that currently can be executed on a single
1133 // GPU, if present and compatible. This has to be coordinated
1134 // across PP ranks on a node, with possible multiple devices
1135 // or sharing devices on a node, either from the user
1136 // selection, or automatically.
1137 auto haveGpus = !gpuIdsToUse.empty();
1138 std::vector<GpuTask> gpuTasksOnThisRank;
1139 if (thisRankHasDuty(cr, DUTY_PP))
1141 if (useGpuForNonbonded)
1143 // Note that any bonded tasks on a GPU always accompany a
1147 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1149 else if (nonbondedTarget == TaskTarget::Gpu)
1151 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1153 else if (bondedTarget == TaskTarget::Gpu)
1155 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because there is none detected.");
1159 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1160 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1166 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1168 else if (pmeTarget == TaskTarget::Gpu)
1170 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1175 GpuTaskAssignment gpuTaskAssignment;
1178 // Produce the task assignment for this rank.
1179 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1180 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1181 useGpuForBonded, pmeRunMode);
1183 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1185 /* Prevent other ranks from continuing after an issue was found
1186 * and reported as a fatal error.
1188 * TODO This function implements a barrier so that MPI runtimes
1189 * can organize an orderly shutdown if one of the ranks has had to
1190 * issue a fatal error in various code already run. When we have
1191 * MPI-aware error handling and reporting, this should be
1196 MPI_Barrier(cr->mpi_comm_mysim);
1202 MPI_Barrier(ms->mpi_comm_masters);
1204 /* We need another barrier to prevent non-master ranks from contiuing
1205 * when an error occured in a different simulation.
1207 MPI_Barrier(cr->mpi_comm_mysim);
1211 /* Now that we know the setup is consistent, check for efficiency */
1212 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1215 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1217 if (thisRankHasDuty(cr, DUTY_PP))
1219 // This works because only one task of each type is currently permitted.
1220 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1221 hasTaskType<GpuTask::Nonbonded>);
1222 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1224 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1225 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1226 init_gpu(nonbondedDeviceInfo);
1228 if (DOMAINDECOMP(cr))
1230 /* When we share GPUs over ranks, we need to know this for the DLB */
1231 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1237 std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1239 gmx_device_info_t *pmeDeviceInfo = nullptr;
1240 // Later, this program could contain kernels that might be later
1241 // re-used as auto-tuning progresses, or subsequent simulations
1243 PmeGpuProgramStorage pmeGpuProgram;
1244 // This works because only one task of each type is currently permitted.
1245 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1246 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1247 if (thisRankHasPmeGpuTask)
1249 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1250 init_gpu(pmeDeviceInfo);
1251 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1252 // TODO It would be nice to move this logic into the factory
1253 // function. See Redmine #2535.
1254 bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1255 if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1257 initializedClfftLibrary = initializeClfftLibrary();
1261 /* getting number of PP/PME threads on this MPI / tMPI rank.
1262 PME: env variable should be read only on one node to make sure it is
1263 identical everywhere;
1265 const int numThreadsOnThisRank =
1266 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1267 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1268 *hwinfo->hardwareTopology,
1269 physicalNodeComm, mdlog);
1271 if (hw_opt.thread_affinity != threadaffOFF)
1273 /* Before setting affinity, check whether the affinity has changed
1274 * - which indicates that probably the OpenMP library has changed it
1275 * since we first checked).
1277 gmx_check_thread_affinity_set(mdlog, cr,
1278 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1280 int numThreadsOnThisNode, intraNodeThreadOffset;
1281 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1282 &intraNodeThreadOffset);
1284 /* Set the CPU affinity */
1285 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1286 numThreadsOnThisRank, numThreadsOnThisNode,
1287 intraNodeThreadOffset, nullptr);
1290 if (mdrunOptions.timingOptions.resetStep > -1)
1292 GMX_LOG(mdlog.info).asParagraph().
1293 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1295 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1299 /* Master synchronizes its value of reset_counters with all nodes
1300 * including PME only nodes */
1301 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1302 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1303 wcycle_set_reset_counters(wcycle, reset_counters);
1306 // Membrane embedding must be initialized before we call init_forcerec()
1311 fprintf(stderr, "Initializing membed");
1313 /* Note that membed cannot work in parallel because mtop is
1314 * changed here. Fix this if we ever want to make it run with
1315 * multiple ranks. */
1316 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1318 .checkpointOptions.period);
1321 std::unique_ptr<MDAtoms> mdAtoms;
1322 std::unique_ptr<gmx_vsite_t> vsite;
1325 if (thisRankHasDuty(cr, DUTY_PP))
1327 /* Initiate forcerecord */
1328 fr = new t_forcerec;
1329 fr->forceProviders = mdModules_->initForceProviders();
1330 init_forcerec(fplog, mdlog, fr, fcd,
1331 inputrec, &mtop, cr, box,
1332 opt2fn("-table", filenames.size(), filenames.data()),
1333 opt2fn("-tablep", filenames.size(), filenames.data()),
1334 opt2fns("-tableb", filenames.size(), filenames.data()),
1335 *hwinfo, nonbondedDeviceInfo,
1340 /* Initialize the mdAtoms structure.
1341 * mdAtoms is not filled with atom data,
1342 * as this can not be done now with domain decomposition.
1344 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1345 if (globalState && thisRankHasPmeGpuTask)
1347 // The pinning of coordinates in the global state object works, because we only use
1348 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1349 // points to the global state object without DD.
1350 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1351 // which should also perform the pinning.
1352 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1355 /* Initialize the virtual site communication */
1356 vsite = initVsite(mtop, cr);
1358 calc_shifts(box, fr->shift_vec);
1360 /* With periodic molecules the charge groups should be whole at start up
1361 * and the virtual sites should not be far from their proper positions.
1363 if (!inputrec->bContinuation && MASTER(cr) &&
1364 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1366 /* Make molecules whole at start of run */
1367 if (fr->ePBC != epbcNONE)
1369 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1373 /* Correct initial vsite positions are required
1374 * for the initial distribution in the domain decomposition
1375 * and for the initial shell prediction.
1377 constructVsitesGlobal(mtop, globalState->x);
1381 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1383 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1384 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1389 /* This is a PME only node */
1391 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1393 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1394 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1397 gmx_pme_t *sepPmeData = nullptr;
1398 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1399 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1400 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1402 /* Initiate PME if necessary,
1403 * either on all nodes or on dedicated PME nodes only. */
1404 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1406 if (mdAtoms && mdAtoms->mdatoms())
1408 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1409 if (EVDW_PME(inputrec->vdwtype))
1411 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1414 if (cr->npmenodes > 0)
1416 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1417 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1418 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1421 if (thisRankHasDuty(cr, DUTY_PME))
1425 pmedata = gmx_pme_init(cr,
1426 getNumPmeDomains(cr->dd),
1428 mtop.natoms, nChargePerturbed != 0, nTypePerturbed != 0,
1429 mdrunOptions.reproducible,
1430 ewaldcoeff_q, ewaldcoeff_lj,
1431 gmx_omp_nthreads_get(emntPME),
1432 pmeRunMode, nullptr,
1433 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1435 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1440 if (EI_DYNAMICS(inputrec->eI))
1442 /* Turn on signal handling on all nodes */
1444 * (A user signal from the PME nodes (if any)
1445 * is communicated to the PP nodes.
1447 signal_handler_install();
1450 if (thisRankHasDuty(cr, DUTY_PP))
1452 /* Assumes uniform use of the number of OpenMP threads */
1453 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1455 if (inputrec->bPull)
1457 /* Initialize pull code */
1458 inputrec->pull_work =
1459 init_pull(fplog, inputrec->pull, inputrec,
1460 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1461 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1463 initPullHistory(inputrec->pull_work, &observablesHistory);
1465 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1467 init_pull_output_files(inputrec->pull_work,
1468 filenames.size(), filenames.data(), oenv,
1469 continuationOptions);
1473 std::unique_ptr<EnforcedRotation> enforcedRotation;
1476 /* Initialize enforced rotation code */
1477 enforcedRotation = init_rot(fplog,
1489 if (inputrec->eSwapCoords != eswapNO)
1491 /* Initialize ion swapping code */
1492 init_swapcoords(fplog, inputrec, opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1493 &mtop, globalState.get(), &observablesHistory,
1494 cr, &atomSets, oenv, mdrunOptions);
1497 /* Let makeConstraints know whether we have essential dynamics constraints.
1498 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1500 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1501 || observablesHistory.edsamHistory);
1502 auto constr = makeConstraints(mtop, *inputrec, doEssentialDynamics,
1503 fplog, *mdAtoms->mdatoms(),
1504 cr, ms, nrnb, wcycle, fr->bMolPBC);
1506 /* Energy terms and groups */
1507 gmx_enerdata_t *enerd;
1509 init_enerdata(mtop.groups.grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
1511 if (DOMAINDECOMP(cr))
1513 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1514 /* This call is not included in init_domain_decomposition mainly
1515 * because fr->cginfo_mb is set later.
1517 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1518 domdecOptions.checkBondedInteractions,
1522 // TODO This is not the right place to manage the lifetime of
1523 // this data structure, but currently it's the easiest way to
1524 // make it work. Later, it should probably be made/updated
1525 // after the workload for the lifetime of a PP domain is
1527 PpForceWorkload ppForceWorkload;
1529 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to integrator.");
1530 /* Now do whatever the user wants us to do (how flexible...) */
1531 Integrator integrator {
1532 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1535 vsite.get(), constr.get(),
1536 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1538 mdModules_->outputProvider(),
1542 &observablesHistory,
1543 mdAtoms.get(), nrnb, wcycle, fr,
1548 walltime_accounting,
1549 std::move(stopHandlerBuilder_)
1551 integrator.run(inputrec->eI, doRerun);
1553 if (inputrec->bPull)
1555 finish_pull(inputrec->pull_work);
1557 destroy_enerdata(enerd);
1562 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1564 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1565 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1568 wallcycle_stop(wcycle, ewcRUN);
1570 /* Finish up, write some stuff
1571 * if rerunMD, don't write last frame again
1573 finish_run(fplog, mdlog, cr,
1574 inputrec, nrnb, wcycle, walltime_accounting,
1575 fr ? fr->nbv.get() : nullptr,
1577 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1579 // clean up cycle counter
1580 wallcycle_destroy(wcycle);
1585 gmx_pme_destroy(pmedata);
1589 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1590 // before we destroy the GPU context(s) in free_gpu_resources().
1591 // Pinned buffers are associated with contexts in CUDA.
1592 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1593 mdAtoms.reset(nullptr);
1594 globalState.reset(nullptr);
1595 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1597 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1598 free_gpu_resources(fr, physicalNodeComm);
1599 free_gpu(nonbondedDeviceInfo);
1600 free_gpu(pmeDeviceInfo);
1601 done_forcerec(fr, mtop.molblock.size(), mtop.groups.grps[egcENER].nr);
1606 free_membed(membed);
1609 gmx_hardware_info_free();
1611 /* Does what it says */
1612 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1613 walltime_accounting_destroy(walltime_accounting);
1616 // Ensure log file content is written
1619 gmx_fio_flush(logFileHandle);
1622 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1623 * exceptions were enabled before function was called. */
1626 gmx_fedisableexcept();
1629 auto rc = static_cast<int>(gmx_get_stop_condition());
1632 /* we need to join all threads. The sub-threads join when they
1633 exit this function, but the master thread needs to be told to
1635 if (PAR(cr) && MASTER(cr))
1639 //TODO free commrec in MPI simulations
1645 Mdrunner::~Mdrunner()
1647 // Clean up of the Manager.
1648 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1649 // but okay as long as threads synchronize some time before adding or accessing
1650 // a new set of restraints.
1651 if (restraintManager_)
1653 restraintManager_->clear();
1654 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1655 "restraints added during runner life time should be cleared at runner destruction.");
1659 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1662 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1663 // Not sure if this should be logged through the md logger or something else,
1664 // but it is helpful to have some sort of INFO level message sent somewhere.
1665 // std::cout << "Registering restraint named " << name << std::endl;
1667 // When multiple restraints are used, it may be wasteful to register them separately.
1668 // Maybe instead register an entire Restraint Manager as a force provider.
1669 restraintManager_->addToSpec(std::move(puller),
1673 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1674 : mdModules_(std::move(mdModules))
1678 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1680 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1681 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1683 class Mdrunner::BuilderImplementation
1686 BuilderImplementation() = delete;
1687 BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1688 SimulationContext * context);
1689 ~BuilderImplementation();
1691 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1692 real forceWarningThreshold);
1694 void addDomdec(const DomdecOptions &options);
1696 void addVerletList(int nstlist);
1698 void addReplicaExchange(const ReplicaExchangeParameters ¶ms);
1700 void addMultiSim(gmx_multisim_t* multisim);
1702 void addNonBonded(const char* nbpu_opt);
1704 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1706 void addBondedTaskAssignment(const char* bonded_opt);
1708 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1710 void addFilenames(ArrayRef <const t_filenm> filenames);
1712 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1714 void addLogFile(t_fileio *logFileHandle);
1716 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1722 // Default parameters copied from runner.h
1723 // \todo Clarify source(s) of default parameters.
1725 const char* nbpu_opt_ = nullptr;
1726 const char* pme_opt_ = nullptr;
1727 const char* pme_fft_opt_ = nullptr;
1728 const char *bonded_opt_ = nullptr;
1730 MdrunOptions mdrunOptions_;
1732 DomdecOptions domdecOptions_;
1734 ReplicaExchangeParameters replicaExchangeParameters_;
1736 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1739 //! Non-owning multisim communicator handle.
1740 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1742 //! Print a warning if any force is larger than this (in kJ/mol nm).
1743 real forceWarningThreshold_ = -1;
1745 //! The modules that comprise the functionality of mdrun.
1746 std::unique_ptr<MDModules> mdModules_;
1748 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1751 * \todo Establish robust protocol to make sure resources remain valid.
1752 * SimulationContext will likely be separated into multiple layers for
1753 * different levels of access and different phases of execution. Ref
1754 * https://redmine.gromacs.org/issues/2375
1755 * https://redmine.gromacs.org/issues/2587
1757 SimulationContext* context_ = nullptr;
1759 //! \brief Parallelism information.
1760 gmx_hw_opt_t hardwareOptions_;
1762 //! filename options for simulation.
1763 ArrayRef<const t_filenm> filenames_;
1765 /*! \brief Handle to output environment.
1767 * \todo gmx_output_env_t needs lifetime management.
1769 gmx_output_env_t* outputEnvironment_ = nullptr;
1771 /*! \brief Non-owning handle to MD log file.
1773 * \todo Context should own output facilities for client.
1774 * \todo Improve log file handle management.
1776 * Code managing the FILE* relies on the ability to set it to
1777 * nullptr to check whether the filehandle is valid.
1779 t_fileio* logFileHandle_ = nullptr;
1782 * \brief Builder for simulation stop signal handler.
1784 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1787 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1788 SimulationContext * context) :
1789 mdModules_(std::move(mdModules)),
1792 GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1795 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1797 Mdrunner::BuilderImplementation &
1798 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1799 real forceWarningThreshold)
1801 mdrunOptions_ = options;
1802 forceWarningThreshold_ = forceWarningThreshold;
1806 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1808 domdecOptions_ = options;
1811 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1816 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1818 replicaExchangeParameters_ = params;
1821 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1823 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1826 Mdrunner Mdrunner::BuilderImplementation::build()
1828 auto newRunner = Mdrunner(std::move(mdModules_));
1830 GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1832 newRunner.mdrunOptions = mdrunOptions_;
1833 newRunner.domdecOptions = domdecOptions_;
1835 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1836 newRunner.hw_opt = hardwareOptions_;
1838 // No invariant to check. This parameter exists to optionally override other behavior.
1839 newRunner.nstlist_cmdline = nstlist_;
1841 newRunner.replExParams = replicaExchangeParameters_;
1843 newRunner.filenames = filenames_;
1845 GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1846 newRunner.cr = context_->communicationRecord_;
1850 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1851 newRunner.ms = *multisim_;
1855 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1858 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1859 // \todo Update sanity checking when output environment has clearly specified invariants.
1860 // Initialization and default values for oenv are not well specified in the current version.
1861 if (outputEnvironment_)
1863 newRunner.oenv = outputEnvironment_;
1867 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1870 newRunner.logFileHandle = logFileHandle_;
1874 newRunner.nbpu_opt = nbpu_opt_;
1878 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1881 if (pme_opt_ && pme_fft_opt_)
1883 newRunner.pme_opt = pme_opt_;
1884 newRunner.pme_fft_opt = pme_fft_opt_;
1888 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1893 newRunner.bonded_opt = bonded_opt_;
1897 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1900 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1902 if (stopHandlerBuilder_)
1904 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1908 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1914 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1916 nbpu_opt_ = nbpu_opt;
1919 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1920 const char* pme_fft_opt)
1923 pme_fft_opt_ = pme_fft_opt;
1926 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1928 bonded_opt_ = bonded_opt;
1931 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1933 hardwareOptions_ = hardwareOptions;
1936 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1938 filenames_ = filenames;
1941 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1943 outputEnvironment_ = outputEnvironment;
1946 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1948 logFileHandle_ = logFileHandle;
1951 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1953 stopHandlerBuilder_ = std::move(builder);
1956 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules,
1957 compat::not_null<SimulationContext*> context) :
1958 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context)}
1962 MdrunnerBuilder::~MdrunnerBuilder() = default;
1964 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1965 real forceWarningThreshold)
1967 impl_->setExtraMdrunOptions(options, forceWarningThreshold);
1971 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1973 impl_->addDomdec(options);
1977 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1979 impl_->addVerletList(nstlist);
1983 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters ¶ms)
1985 impl_->addReplicaExchange(params);
1989 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1991 impl_->addMultiSim(multisim);
1995 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1997 impl_->addNonBonded(nbpu_opt);
2001 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
2002 const char* pme_fft_opt)
2004 // The builder method may become more general in the future, but in this version,
2005 // parameters for PME electrostatics are both required and the only parameters
2007 if (pme_opt && pme_fft_opt)
2009 impl_->addPME(pme_opt, pme_fft_opt);
2013 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2018 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2020 impl_->addBondedTaskAssignment(bonded_opt);
2024 Mdrunner MdrunnerBuilder::build()
2026 return impl_->build();
2029 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2031 impl_->addHardwareOptions(hardwareOptions);
2035 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2037 impl_->addFilenames(filenames);
2041 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2043 impl_->addOutputEnvironment(outputEnvironment);
2047 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2049 impl_->addLogFile(logFileHandle);
2053 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2055 impl_->addStopHandlerBuilder(std::move(builder));
2059 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2061 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;