Clean up PME domain count passing
[gromacs.git] / src / gromacs / mdrun / runner.cpp
blob729a293626a517db580cce17db6b9ace28a40522
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \internal \file
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
44 #include "gmxpre.h"
46 #include "runner.h"
48 #include "config.h"
50 #include <cassert>
51 #include <csignal>
52 #include <cstdlib>
53 #include <cstring>
55 #include <algorithm>
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/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/hardware/cpuinfo.h"
69 #include "gromacs/hardware/detecthardware.h"
70 #include "gromacs/hardware/printhardware.h"
71 #include "gromacs/listed-forces/disre.h"
72 #include "gromacs/listed-forces/orires.h"
73 #include "gromacs/math/functions.h"
74 #include "gromacs/math/utilities.h"
75 #include "gromacs/math/vec.h"
76 #include "gromacs/mdlib/calc_verletbuf.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/deform.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/forcerec.h"
81 #include "gromacs/mdlib/gmx_omp_nthreads.h"
82 #include "gromacs/mdlib/main.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdrun.h"
86 #include "gromacs/mdlib/membed.h"
87 #include "gromacs/mdlib/nb_verlet.h"
88 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
89 #include "gromacs/mdlib/nbnxn_search.h"
90 #include "gromacs/mdlib/nbnxn_tuning.h"
91 #include "gromacs/mdlib/qmmm.h"
92 #include "gromacs/mdlib/repl_ex.h"
93 #include "gromacs/mdlib/sighandler.h"
94 #include "gromacs/mdlib/sim_util.h"
95 #include "gromacs/mdrunutility/mdmodules.h"
96 #include "gromacs/mdrunutility/threadaffinity.h"
97 #include "gromacs/mdtypes/commrec.h"
98 #include "gromacs/mdtypes/inputrec.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/mdtypes/observableshistory.h"
101 #include "gromacs/mdtypes/state.h"
102 #include "gromacs/pbcutil/pbc.h"
103 #include "gromacs/pulling/pull.h"
104 #include "gromacs/pulling/pull_rotation.h"
105 #include "gromacs/taskassignment/decidegpuusage.h"
106 #include "gromacs/taskassignment/resourcedivision.h"
107 #include "gromacs/taskassignment/taskassignment.h"
108 #include "gromacs/taskassignment/usergpuids.h"
109 #include "gromacs/timing/wallcycle.h"
110 #include "gromacs/topology/mtop_util.h"
111 #include "gromacs/trajectory/trajectoryframe.h"
112 #include "gromacs/utility/basenetwork.h"
113 #include "gromacs/utility/cstringutil.h"
114 #include "gromacs/utility/exceptions.h"
115 #include "gromacs/utility/fatalerror.h"
116 #include "gromacs/utility/filestream.h"
117 #include "gromacs/utility/gmxassert.h"
118 #include "gromacs/utility/gmxmpi.h"
119 #include "gromacs/utility/logger.h"
120 #include "gromacs/utility/loggerbuilder.h"
121 #include "gromacs/utility/physicalnodecommunicator.h"
122 #include "gromacs/utility/pleasecite.h"
123 #include "gromacs/utility/programcontext.h"
124 #include "gromacs/utility/smalloc.h"
125 #include "gromacs/utility/stringutil.h"
127 #include "integrator.h"
129 #ifdef GMX_FAHCORE
130 #include "corewrap.h"
131 #endif
133 //! First step used in pressure scaling
134 gmx_int64_t deform_init_init_step_tpx;
135 //! Initial box for pressure scaling
136 matrix deform_init_box_tpx;
137 //! MPI variable for use in pressure scaling
138 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
140 namespace gmx
143 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
145 * Used to ensure that the master thread does not modify mdrunner during copy
146 * on the spawned threads. */
147 static void threadMpiMdrunnerAccessBarrier()
149 #if GMX_THREAD_MPI
150 MPI_Barrier(MPI_COMM_WORLD);
151 #endif
154 void Mdrunner::reinitializeOnSpawnedThread()
156 threadMpiMdrunnerAccessBarrier();
158 cr = reinitialize_commrec_for_this_thread(cr);
160 GMX_RELEASE_ASSERT(!MASTER(cr), "reinitializeOnSpawnedThread should only be called on spawned threads");
162 // Only the master rank writes to the log file
163 fplog = nullptr;
166 /*! \brief The callback used for running on spawned threads.
168 * Obtains the pointer to the master mdrunner object from the one
169 * argument permitted to the thread-launch API call, copies it to make
170 * a new runner for this thread, reinitializes necessary data, and
171 * proceeds to the simulation. */
172 static void mdrunner_start_fn(const void *arg)
176 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
177 /* copy the arg list to make sure that it's thread-local. This
178 doesn't copy pointed-to items, of course; fnm, cr and fplog
179 are reset in the call below, all others should be const. */
180 gmx::Mdrunner mdrunner = *masterMdrunner;
181 mdrunner.reinitializeOnSpawnedThread();
182 mdrunner.mdrunner();
184 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
188 /*! \brief Start thread-MPI threads.
190 * Called by mdrunner() to start a specific number of threads
191 * (including the main thread) for thread-parallel runs. This in turn
192 * calls mdrunner() for each thread. All options are the same as for
193 * mdrunner(). */
194 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
197 /* first check whether we even need to start tMPI */
198 if (numThreadsToLaunch < 2)
200 return cr;
203 #if GMX_THREAD_MPI
204 /* now spawn new threads that start mdrunner_start_fn(), while
205 the main thread returns, we set thread affinity later */
206 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
207 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
209 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
212 threadMpiMdrunnerAccessBarrier();
213 #else
214 GMX_UNUSED_VALUE(mdrunner_start_fn);
215 #endif
217 return reinitialize_commrec_for_this_thread(cr);
220 } // namespace
222 /*! \brief Initialize variables for Verlet scheme simulation */
223 static void prepare_verlet_scheme(FILE *fplog,
224 t_commrec *cr,
225 t_inputrec *ir,
226 int nstlist_cmdline,
227 const gmx_mtop_t *mtop,
228 const matrix box,
229 bool makeGpuPairList,
230 const gmx::CpuInfo &cpuinfo)
232 /* For NVE simulations, we will retain the initial list buffer */
233 if (EI_DYNAMICS(ir->eI) &&
234 ir->verletbuf_tol > 0 &&
235 !(EI_MD(ir->eI) && ir->etc == etcNO))
237 /* Update the Verlet buffer size for the current run setup */
239 /* Here we assume SIMD-enabled kernels are being used. But as currently
240 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
241 * and 4x2 gives a larger buffer than 4x4, this is ok.
243 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
244 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
246 real rlist_new;
247 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &listSetup, nullptr, &rlist_new);
249 if (rlist_new != ir->rlist)
251 if (fplog != nullptr)
253 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
254 ir->rlist, rlist_new,
255 listSetup.cluster_size_i, listSetup.cluster_size_j);
257 ir->rlist = rlist_new;
261 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
263 gmx_fatal(FARGS, "Can not set nstlist without %s",
264 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
267 if (EI_DYNAMICS(ir->eI))
269 /* Set or try nstlist values */
270 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
274 /*! \brief Override the nslist value in inputrec
276 * with value passed on the command line (if any)
278 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
279 gmx_int64_t nsteps_cmdline,
280 t_inputrec *ir)
282 assert(ir);
284 /* override with anything else than the default -2 */
285 if (nsteps_cmdline > -2)
287 char sbuf_steps[STEPSTRSIZE];
288 char sbuf_msg[STRLEN];
290 ir->nsteps = nsteps_cmdline;
291 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
293 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
294 gmx_step_str(nsteps_cmdline, sbuf_steps),
295 fabs(nsteps_cmdline*ir->delta_t));
297 else
299 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
300 gmx_step_str(nsteps_cmdline, sbuf_steps));
303 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
305 else if (nsteps_cmdline < -2)
307 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
308 nsteps_cmdline);
310 /* Do nothing if nsteps_cmdline == -2 */
313 namespace gmx
316 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
318 * If not, and if a warning may be issued, logs a warning about
319 * falling back to CPU code. With thread-MPI, only the first
320 * call to this function should have \c issueWarning true. */
321 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
322 const t_inputrec *ir,
323 bool issueWarning)
325 if (ir->opts.ngener - ir->nwall > 1)
327 /* The GPU code does not support more than one energy group.
328 * If the user requested GPUs explicitly, a fatal error is given later.
330 if (issueWarning)
332 GMX_LOG(mdlog.warning).asParagraph()
333 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
334 "For better performance, run on the GPU without energy groups and then do "
335 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
337 return false;
339 return true;
342 //! Initializes the logger for mdrun.
343 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
345 gmx::LoggerBuilder builder;
346 if (fplog != nullptr)
348 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
350 if (cr == nullptr || SIMMASTER(cr))
352 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
353 &gmx::TextOutputFile::standardError());
355 return builder.build();
358 //! Make a TaskTarget from an mdrun argument string.
359 static TaskTarget findTaskTarget(const char *optionString)
361 TaskTarget returnValue = TaskTarget::Auto;
363 if (strncmp(optionString, "auto", 3) == 0)
365 returnValue = TaskTarget::Auto;
367 else if (strncmp(optionString, "cpu", 3) == 0)
369 returnValue = TaskTarget::Cpu;
371 else if (strncmp(optionString, "gpu", 3) == 0)
373 returnValue = TaskTarget::Gpu;
375 else
377 GMX_ASSERT(false, "Option string should have been checked for sanity already");
380 return returnValue;
383 int Mdrunner::mdrunner()
385 matrix box;
386 t_nrnb *nrnb;
387 t_forcerec *fr = nullptr;
388 t_fcdata *fcd = nullptr;
389 real ewaldcoeff_q = 0;
390 real ewaldcoeff_lj = 0;
391 gmx_vsite_t *vsite = nullptr;
392 int nChargePerturbed = -1, nTypePerturbed = 0;
393 gmx_wallcycle_t wcycle;
394 gmx_walltime_accounting_t walltime_accounting = nullptr;
395 int rc;
396 gmx_int64_t reset_counters;
397 int nthreads_pme = 1;
398 gmx_membed_t * membed = nullptr;
399 gmx_hw_info_t *hwinfo = nullptr;
401 /* CAUTION: threads may be started later on in this function, so
402 cr doesn't reflect the final parallel state right now */
403 std::unique_ptr<gmx::MDModules> mdModules(new gmx::MDModules);
404 t_inputrec inputrecInstance;
405 t_inputrec *inputrec = &inputrecInstance;
406 gmx_mtop_t mtop;
408 if (mdrunOptions.continuationOptions.appendFiles)
410 fplog = nullptr;
413 bool doMembed = opt2bSet("-membed", nfile, fnm);
414 bool doRerun = mdrunOptions.rerun;
416 // Handle task-assignment related user options.
417 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
418 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
419 std::vector<int> gpuIdsAvailable;
422 gpuIdsAvailable = parseUserGpuIds(hw_opt.gpuIdsAvailable);
423 // TODO We could put the GPU IDs into a std::map to find
424 // duplicates, but for the small numbers of IDs involved, this
425 // code is simple and fast.
426 for (size_t i = 0; i != gpuIdsAvailable.size(); ++i)
428 for (size_t j = i+1; j != gpuIdsAvailable.size(); ++j)
430 if (gpuIdsAvailable[i] == gpuIdsAvailable[j])
432 GMX_THROW(InvalidInputError(formatString("The string of available GPU device IDs '%s' may not contain duplicate device IDs", hw_opt.gpuIdsAvailable.c_str())));
437 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
439 std::vector<int> userGpuTaskAssignment;
442 userGpuTaskAssignment = parseUserGpuIds(hw_opt.userGpuTaskAssignment);
444 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
445 auto nonbondedTarget = findTaskTarget(nbpu_opt);
446 auto pmeTarget = findTaskTarget(pme_opt);
447 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
448 PmeRunMode pmeRunMode = PmeRunMode::None;
450 // Here we assume that SIMMASTER(cr) does not change even after the
451 // threads are started.
452 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
453 gmx::MDLogger mdlog(logOwner.logger());
455 // TODO The thread-MPI master rank makes a working
456 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
457 // after the threads have been launched. This works because no use
458 // is made of that communicator until after the execution paths
459 // have rejoined. But it is likely that we can improve the way
460 // this is expressed, e.g. by expressly running detection only the
461 // master rank for thread-MPI, rather than relying on the mutex
462 // and reference count.
463 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
464 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
466 gmx_print_detected_hardware(fplog, cr, ms, mdlog, hwinfo);
468 std::vector<int> gpuIdsToUse;
469 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
470 if (gpuIdsAvailable.empty())
472 gpuIdsToUse = compatibleGpus;
474 else
476 for (const auto &availableGpuId : gpuIdsAvailable)
478 bool availableGpuIsCompatible = false;
479 for (const auto &compatibleGpuId : compatibleGpus)
481 if (availableGpuId == compatibleGpuId)
483 availableGpuIsCompatible = true;
484 break;
487 if (!availableGpuIsCompatible)
489 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);
491 gpuIdsToUse.push_back(availableGpuId);
495 if (fplog != nullptr)
497 /* Print references after all software/hardware printing */
498 please_cite(fplog, "Abraham2015");
499 please_cite(fplog, "Pall2015");
500 please_cite(fplog, "Pronk2013");
501 please_cite(fplog, "Hess2008b");
502 please_cite(fplog, "Spoel2005a");
503 please_cite(fplog, "Lindahl2001a");
504 please_cite(fplog, "Berendsen95a");
507 std::unique_ptr<t_state> globalState;
509 if (SIMMASTER(cr))
511 /* Only the master rank has the global state */
512 globalState = std::unique_ptr<t_state>(new t_state);
514 /* Read (nearly) all data required for the simulation */
515 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), &mtop);
517 if (inputrec->cutoff_scheme != ecutsVERLET)
519 if (nstlist_cmdline > 0)
521 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
524 if (!compatibleGpus.empty())
526 GMX_LOG(mdlog.warning).asParagraph().appendText(
527 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
528 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
533 /* Check and update the hardware options for internal consistency */
534 check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
536 /* Early check for externally set process affinity. */
537 gmx_check_thread_affinity_set(mdlog, cr,
538 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
540 if (GMX_THREAD_MPI && SIMMASTER(cr))
542 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
544 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
547 /* Since the master knows the cut-off scheme, update hw_opt for this.
548 * This is done later for normal MPI and also once more with tMPI
549 * for all tMPI ranks.
551 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
553 bool useGpuForNonbonded = false;
554 bool useGpuForPme = false;
557 // If the user specified the number of ranks, then we must
558 // respect that, but in default mode, we need to allow for
559 // the number of GPUs to choose the number of ranks.
561 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
562 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
563 inputrec->cutoff_scheme == ecutsVERLET,
564 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
565 hw_opt.nthreads_tmpi);
566 auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
567 auto canUseGpuForPme = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
568 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
569 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
570 canUseGpuForPme, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
573 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
575 /* Determine how many thread-MPI ranks to start.
577 * TODO Over-writing the user-supplied value here does
578 * prevent any possible subsequent checks from working
579 * correctly. */
580 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
581 &hw_opt,
582 gpuIdsToUse,
583 useGpuForNonbonded,
584 useGpuForPme,
585 inputrec, &mtop,
586 mdlog,
587 doMembed);
589 // Now start the threads for thread MPI.
590 cr = spawnThreads(hw_opt.nthreads_tmpi);
591 /* The main thread continues here with a new cr. We don't deallocate
592 the old cr because other threads may still be reading it. */
593 // TODO Both master and spawned threads call dup_tfn and
594 // reinitialize_commrec_for_this_thread. Find a way to express
595 // this better.
596 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
598 // END OF CAUTION: cr and physicalNodeComm are now reliable
600 if (PAR(cr))
602 /* now broadcast everything to the non-master nodes/threads: */
603 init_parallel(cr, inputrec, &mtop);
606 // Now each rank knows the inputrec that SIMMASTER read and used,
607 // and (if applicable) cr->nnodes has been assigned the number of
608 // thread-MPI ranks that have been chosen. The ranks can now all
609 // run the task-deciding functions and will agree on the result
610 // without needing to communicate.
612 // TODO Should we do the communication in debug mode to support
613 // having an assertion?
615 // Note that these variables describe only their own node.
616 bool useGpuForNonbonded = false;
617 bool useGpuForPme = false;
620 // It's possible that there are different numbers of GPUs on
621 // different nodes, which is the user's responsibilty to
622 // handle. If unsuitable, we will notice that during task
623 // assignment.
624 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
625 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
626 emulateGpuNonbonded, inputrec->cutoff_scheme == ecutsVERLET,
627 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
628 gpusWereDetected);
629 auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
630 auto canUseGpuForPme = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
631 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
632 canUseGpuForPme, cr->nnodes, domdecOptions.numPmeRanks,
633 gpusWereDetected);
635 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
636 if (pmeRunMode == PmeRunMode::GPU)
638 if (pmeFftTarget == TaskTarget::Cpu)
640 pmeRunMode = PmeRunMode::Mixed;
643 else if (pmeFftTarget == TaskTarget::Gpu)
645 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.");
648 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
650 // TODO: Error handling
651 mdModules->assignOptionsToModules(*inputrec->params, nullptr);
653 if (fplog != nullptr)
655 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
656 fprintf(fplog, "\n");
659 if (SIMMASTER(cr))
661 /* now make sure the state is initialized and propagated */
662 set_state_entries(globalState.get(), inputrec);
665 /* NM and TPI parallelize over force/energy calculations, not atoms,
666 * so we need to initialize and broadcast the global state.
668 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
670 if (!MASTER(cr))
672 globalState = std::unique_ptr<t_state>(new t_state);
674 broadcastStateWithoutDynamics(cr, globalState.get());
677 /* A parallel command line option consistency check that we can
678 only do after any threads have started. */
679 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
680 domdecOptions.numCells[YY] > 1 ||
681 domdecOptions.numCells[ZZ] > 1 ||
682 domdecOptions.numPmeRanks > 0))
684 gmx_fatal(FARGS,
685 "The -dd or -npme option request a parallel simulation, "
686 #if !GMX_MPI
687 "but %s was compiled without threads or MPI enabled"
688 #else
689 #if GMX_THREAD_MPI
690 "but the number of MPI-threads (option -ntmpi) is not set or is 1"
691 #else
692 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
693 #endif
694 #endif
695 , output_env_get_program_display_name(oenv)
699 if (doRerun &&
700 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
702 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
705 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
707 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
710 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
712 if (domdecOptions.numPmeRanks > 0)
714 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
715 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
718 domdecOptions.numPmeRanks = 0;
721 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
723 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
724 * improve performance with many threads per GPU, since our OpenMP
725 * scaling is bad, but it's difficult to automate the setup.
727 domdecOptions.numPmeRanks = 0;
729 if (useGpuForPme)
731 if (domdecOptions.numPmeRanks < 0)
733 domdecOptions.numPmeRanks = 0;
734 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
736 else
738 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
742 #ifdef GMX_FAHCORE
743 if (MASTER(cr))
745 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
747 #endif
749 /* NMR restraints must be initialized before load_checkpoint,
750 * since with time averaging the history is added to t_state.
751 * For proper consistency check we therefore need to extend
752 * t_state here.
753 * So the PME-only nodes (if present) will also initialize
754 * the distance restraints.
756 snew(fcd, 1);
758 /* This needs to be called before read_checkpoint to extend the state */
759 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
761 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
763 if (inputrecDeform(inputrec))
765 /* Store the deform reference box before reading the checkpoint */
766 if (SIMMASTER(cr))
768 copy_mat(globalState->box, box);
770 if (PAR(cr))
772 gmx_bcast(sizeof(box), box, cr);
774 /* Because we do not have the update struct available yet
775 * in which the reference values should be stored,
776 * we store them temporarily in static variables.
777 * This should be thread safe, since they are only written once
778 * and with identical values.
780 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
781 deform_init_init_step_tpx = inputrec->init_step;
782 copy_mat(box, deform_init_box_tpx);
783 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
786 ObservablesHistory observablesHistory = {};
788 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
790 if (continuationOptions.startedFromCheckpoint)
792 /* Check if checkpoint file exists before doing continuation.
793 * This way we can use identical input options for the first and subsequent runs...
795 gmx_bool bReadEkin;
797 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
798 cr, domdecOptions.numCells,
799 inputrec, globalState.get(),
800 &bReadEkin, &observablesHistory,
801 continuationOptions.appendFiles,
802 continuationOptions.appendFilesOptionSet,
803 mdrunOptions.reproducible);
805 if (bReadEkin)
807 continuationOptions.haveReadEkin = true;
811 if (SIMMASTER(cr) && continuationOptions.appendFiles)
813 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
814 continuationOptions.appendFiles, &fplog);
815 logOwner = buildLogger(fplog, nullptr);
816 mdlog = logOwner.logger();
819 if (mdrunOptions.numStepsCommandline > -2)
821 GMX_LOG(mdlog.info).asParagraph().
822 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
823 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
825 /* override nsteps with value set on the commamdline */
826 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
828 if (SIMMASTER(cr))
830 copy_mat(globalState->box, box);
833 if (PAR(cr))
835 gmx_bcast(sizeof(box), box, cr);
838 /* Update rlist and nstlist. */
839 if (inputrec->cutoff_scheme == ecutsVERLET)
841 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
842 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
845 /* Initalize the domain decomposition */
846 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
847 inputrec->eI == eiNM))
849 const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
851 cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
852 &mtop, inputrec,
853 box, xOnMaster);
854 // Note that local state still does not exist yet.
856 else
858 /* PME, if used, is done on all nodes with 1D decomposition */
859 cr->npmenodes = 0;
860 cr->duty = (DUTY_PP | DUTY_PME);
862 if (inputrec->ePBC == epbcSCREW)
864 gmx_fatal(FARGS,
865 "pbc=%s is only implemented with domain decomposition",
866 epbc_names[inputrec->ePBC]);
870 if (PAR(cr))
872 /* After possible communicator splitting in make_dd_communicators.
873 * we can set up the intra/inter node communication.
875 gmx_setup_nodecomm(fplog, cr);
878 #if GMX_MPI
879 if (isMultiSim(ms))
881 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
882 "This is simulation %d out of %d running as a composite GROMACS\n"
883 "multi-simulation job. Setup for this simulation:\n",
884 ms->sim, ms->nsim);
886 GMX_LOG(mdlog.warning).appendTextFormatted(
887 "Using %d MPI %s\n",
888 cr->nnodes,
889 #if GMX_THREAD_MPI
890 cr->nnodes == 1 ? "thread" : "threads"
891 #else
892 cr->nnodes == 1 ? "process" : "processes"
893 #endif
895 fflush(stderr);
896 #endif
898 /* Check and update hw_opt for the cut-off scheme */
899 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
901 /* Check and update the number of OpenMP threads requested */
902 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
903 pmeRunMode, mtop);
905 gmx_omp_nthreads_init(mdlog, cr,
906 hwinfo->nthreads_hw_avail,
907 physicalNodeComm.size_,
908 hw_opt.nthreads_omp,
909 hw_opt.nthreads_omp_pme,
910 !thisRankHasDuty(cr, DUTY_PP),
911 inputrec->cutoff_scheme == ecutsVERLET);
913 #ifndef NDEBUG
914 if (EI_TPI(inputrec->eI) &&
915 inputrec->cutoff_scheme == ecutsVERLET)
917 gmx_feenableexcept();
919 #endif
921 // Build a data structure that expresses which kinds of non-bonded
922 // task are handled by this rank.
924 // TODO Later, this might become a loop over all registered modules
925 // relevant to the mdp inputs, to find those that have such tasks.
927 // TODO This could move before init_domain_decomposition() as part
928 // of refactoring that separates the responsibility for duty
929 // assignment from setup for communication between tasks, and
930 // setup for tasks handled with a domain (ie including short-ranged
931 // tasks, bonded tasks, etc.).
933 // Note that in general useGpuForNonbonded, etc. can have a value
934 // that is inconsistent with the presence of actual GPUs on any
935 // rank, and that is not known to be a problem until the
936 // duty of the ranks on a node become node.
938 // TODO Later we might need the concept of computeTasksOnThisRank,
939 // from which we construct gpuTasksOnThisRank.
941 // Currently the DD code assigns duty to ranks that can
942 // include PP work that currently can be executed on a single
943 // GPU, if present and compatible. This has to be coordinated
944 // across PP ranks on a node, with possible multiple devices
945 // or sharing devices on a node, either from the user
946 // selection, or automatically.
947 auto haveGpus = !gpuIdsToUse.empty();
948 std::vector<GpuTask> gpuTasksOnThisRank;
949 if (thisRankHasDuty(cr, DUTY_PP))
951 if (useGpuForNonbonded)
953 if (haveGpus)
955 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
957 else if (nonbondedTarget == TaskTarget::Gpu)
959 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
963 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
964 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
966 if (useGpuForPme)
968 if (haveGpus)
970 gpuTasksOnThisRank.push_back(GpuTask::Pme);
972 else if (pmeTarget == TaskTarget::Gpu)
974 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
979 GpuTaskAssignment gpuTaskAssignment;
982 // Produce the task assignment for this rank.
983 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
984 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank);
986 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
988 /* Prevent other ranks from continuing after an issue was found
989 * and reported as a fatal error.
991 * TODO This function implements a barrier so that MPI runtimes
992 * can organize an orderly shutdown if one of the ranks has had to
993 * issue a fatal error in various code already run. When we have
994 * MPI-aware error handling and reporting, this should be
995 * improved. */
996 #if GMX_MPI
997 if (PAR(cr))
999 MPI_Barrier(cr->mpi_comm_mysim);
1001 if (isMultiSim(ms))
1003 if (SIMMASTER(cr))
1005 MPI_Barrier(ms->mpi_comm_masters);
1007 /* We need another barrier to prevent non-master ranks from contiuing
1008 * when an error occured in a different simulation.
1010 MPI_Barrier(cr->mpi_comm_mysim);
1012 #endif
1014 /* Now that we know the setup is consistent, check for efficiency */
1015 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1016 cr, mdlog);
1018 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1020 if (thisRankHasDuty(cr, DUTY_PP))
1022 // This works because only one task of each type is currently permitted.
1023 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1024 hasTaskType<GpuTask::Nonbonded>);
1025 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1027 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1028 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1029 init_gpu(mdlog, nonbondedDeviceInfo);
1031 if (DOMAINDECOMP(cr))
1033 /* When we share GPUs over ranks, we need to know this for the DLB */
1034 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1040 gmx_device_info_t *pmeDeviceInfo = nullptr;
1041 // This works because only one task of each type is currently permitted.
1042 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1043 if (pmeGpuTaskMapping != gpuTaskAssignment.end())
1045 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1046 init_gpu(mdlog, pmeDeviceInfo);
1049 /* getting number of PP/PME threads
1050 PME: env variable should be read only on one node to make sure it is
1051 identical everywhere;
1053 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1055 int numThreadsOnThisRank;
1056 /* threads on this MPI process or TMPI thread */
1057 if (thisRankHasDuty(cr, DUTY_PP))
1059 numThreadsOnThisRank = gmx_omp_nthreads_get(emntNonbonded);
1061 else
1063 numThreadsOnThisRank = nthreads_pme;
1066 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1067 *hwinfo->hardwareTopology,
1068 physicalNodeComm, mdlog);
1070 if (hw_opt.thread_affinity != threadaffOFF)
1072 /* Before setting affinity, check whether the affinity has changed
1073 * - which indicates that probably the OpenMP library has changed it
1074 * since we first checked).
1076 gmx_check_thread_affinity_set(mdlog, cr,
1077 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1079 int numThreadsOnThisNode, intraNodeThreadOffset;
1080 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1081 &intraNodeThreadOffset);
1083 /* Set the CPU affinity */
1084 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1085 numThreadsOnThisRank, numThreadsOnThisNode,
1086 intraNodeThreadOffset, nullptr);
1089 if (mdrunOptions.timingOptions.resetStep > -1)
1091 GMX_LOG(mdlog.info).asParagraph().
1092 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1094 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1096 if (PAR(cr))
1098 /* Master synchronizes its value of reset_counters with all nodes
1099 * including PME only nodes */
1100 reset_counters = wcycle_get_reset_counters(wcycle);
1101 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1102 wcycle_set_reset_counters(wcycle, reset_counters);
1105 // Membrane embedding must be initialized before we call init_forcerec()
1106 if (doMembed)
1108 if (MASTER(cr))
1110 fprintf(stderr, "Initializing membed");
1112 /* Note that membed cannot work in parallel because mtop is
1113 * changed here. Fix this if we ever want to make it run with
1114 * multiple ranks. */
1115 membed = init_membed(fplog, nfile, fnm, &mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1118 std::unique_ptr<MDAtoms> mdAtoms;
1120 snew(nrnb, 1);
1121 if (thisRankHasDuty(cr, DUTY_PP))
1123 /* Initiate forcerecord */
1124 fr = mk_forcerec();
1125 fr->forceProviders = mdModules->initForceProviders();
1126 init_forcerec(fplog, mdlog, fr, fcd,
1127 inputrec, &mtop, cr, box,
1128 opt2fn("-table", nfile, fnm),
1129 opt2fn("-tablep", nfile, fnm),
1130 opt2fns("-tableb", nfile, fnm),
1131 *hwinfo, nonbondedDeviceInfo,
1132 FALSE,
1133 pforce);
1135 /* Initialize QM-MM */
1136 if (fr->bQMMM)
1138 GMX_LOG(mdlog.info).asParagraph().
1139 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
1140 "version. Please get in touch with the developers if you find the support useful, "
1141 "as help is needed if the functionality is to continue to be available.");
1142 init_QMMMrec(cr, &mtop, inputrec, fr);
1145 /* Initialize the mdAtoms structure.
1146 * mdAtoms is not filled with atom data,
1147 * as this can not be done now with domain decomposition.
1149 const bool useGpuForPme = (pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed);
1150 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, useGpuForPme && thisRankHasDuty(cr, DUTY_PME));
1151 if (globalState)
1153 // The pinning of coordinates in the global state object works, because we only use
1154 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1155 // points to the global state object without DD.
1156 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1157 // which should also perform the pinning.
1158 changePinningPolicy(&globalState->x, useGpuForPme ? PinningPolicy::CanBePinned : PinningPolicy::CannotBePinned);
1161 /* Initialize the virtual site communication */
1162 vsite = initVsite(mtop, cr);
1164 calc_shifts(box, fr->shift_vec);
1166 /* With periodic molecules the charge groups should be whole at start up
1167 * and the virtual sites should not be far from their proper positions.
1169 if (!inputrec->bContinuation && MASTER(cr) &&
1170 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1172 /* Make molecules whole at start of run */
1173 if (fr->ePBC != epbcNONE)
1175 rvec *xGlobal = as_rvec_array(globalState->x.data());
1176 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, xGlobal);
1178 if (vsite)
1180 /* Correct initial vsite positions are required
1181 * for the initial distribution in the domain decomposition
1182 * and for the initial shell prediction.
1184 constructVsitesGlobal(mtop, globalState->x);
1188 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1190 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1191 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1194 else
1196 /* This is a PME only node */
1198 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1200 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1201 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1204 gmx_pme_t *sepPmeData = nullptr;
1205 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1206 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1207 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1209 /* Initiate PME if necessary,
1210 * either on all nodes or on dedicated PME nodes only. */
1211 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1213 if (mdAtoms && mdAtoms->mdatoms())
1215 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1216 if (EVDW_PME(inputrec->vdwtype))
1218 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1221 if (cr->npmenodes > 0)
1223 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1224 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1225 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1228 if (thisRankHasDuty(cr, DUTY_PME))
1232 pmedata = gmx_pme_init(cr,
1233 getNumPmeDomains(cr->dd),
1234 inputrec,
1235 mtop.natoms, nChargePerturbed, nTypePerturbed,
1236 mdrunOptions.reproducible,
1237 ewaldcoeff_q, ewaldcoeff_lj,
1238 nthreads_pme,
1239 pmeRunMode, nullptr, pmeDeviceInfo, mdlog);
1241 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1246 if (EI_DYNAMICS(inputrec->eI))
1248 /* Turn on signal handling on all nodes */
1250 * (A user signal from the PME nodes (if any)
1251 * is communicated to the PP nodes.
1253 signal_handler_install();
1256 if (thisRankHasDuty(cr, DUTY_PP))
1258 /* Assumes uniform use of the number of OpenMP threads */
1259 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1261 if (inputrec->bPull)
1263 /* Initialize pull code */
1264 inputrec->pull_work =
1265 init_pull(fplog, inputrec->pull, inputrec,
1266 &mtop, cr, inputrec->fepvals->init_lambda);
1267 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1269 init_pull_output_files(inputrec->pull_work,
1270 nfile, fnm, oenv,
1271 continuationOptions);
1275 if (inputrec->bRot)
1277 /* Initialize enforced rotation code */
1278 init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), &mtop, oenv, mdrunOptions);
1281 /* Let init_constraints know whether we have essential dynamics constraints.
1282 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1284 bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1286 Constraints *constr = init_constraints(fplog, &mtop, inputrec, doEdsam, cr);
1288 if (DOMAINDECOMP(cr))
1290 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1291 /* This call is not included in init_domain_decomposition mainly
1292 * because fr->cginfo_mb is set later.
1294 dd_init_bondeds(fplog, cr->dd, &mtop, vsite, inputrec,
1295 domdecOptions.checkBondedInteractions,
1296 fr->cginfo_mb);
1299 /* Now do whatever the user wants us to do (how flexible...) */
1300 Integrator integrator {
1301 fplog, cr, ms, mdlog, nfile, fnm,
1302 oenv,
1303 mdrunOptions,
1304 vsite, constr,
1305 mdModules->outputProvider(),
1306 inputrec, &mtop,
1307 fcd,
1308 globalState.get(),
1309 &observablesHistory,
1310 mdAtoms.get(), nrnb, wcycle, fr,
1311 replExParams,
1312 membed,
1313 walltime_accounting
1315 integrator.run(inputrec->eI);
1316 if (inputrec->bRot)
1318 finish_rot(inputrec->rot);
1321 if (inputrec->bPull)
1323 finish_pull(inputrec->pull_work);
1327 else
1329 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1330 /* do PME only */
1331 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1332 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1335 wallcycle_stop(wcycle, ewcRUN);
1337 /* Finish up, write some stuff
1338 * if rerunMD, don't write last frame again
1340 finish_run(fplog, mdlog, cr,
1341 inputrec, nrnb, wcycle, walltime_accounting,
1342 fr ? fr->nbv : nullptr,
1343 pmedata,
1344 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1346 // Free PME data
1347 if (pmedata)
1349 gmx_pme_destroy(pmedata);
1350 pmedata = nullptr;
1353 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1354 // before we destroy the GPU context(s) in free_gpu_resources().
1355 // Pinned buffers are associated with contexts in CUDA.
1356 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1357 mdAtoms.reset(nullptr);
1358 globalState.reset(nullptr);
1359 mdModules.reset(nullptr); // destruct force providers here as they might also use the GPU
1361 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1362 free_gpu_resources(fr, physicalNodeComm);
1363 free_gpu(nonbondedDeviceInfo);
1364 free_gpu(pmeDeviceInfo);
1365 done_forcerec(fr, mtop.molblock.size(), mtop.groups.grps[egcENER].nr);
1366 sfree(fcd);
1368 if (doMembed)
1370 free_membed(membed);
1373 gmx_hardware_info_free();
1375 /* Does what it says */
1376 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1377 walltime_accounting_destroy(walltime_accounting);
1378 sfree(nrnb);
1380 /* Close logfile already here if we were appending to it */
1381 if (MASTER(cr) && continuationOptions.appendFiles)
1383 gmx_log_close(fplog);
1384 fplog = nullptr;
1387 rc = (int)gmx_get_stop_condition();
1389 #if GMX_THREAD_MPI
1390 /* we need to join all threads. The sub-threads join when they
1391 exit this function, but the master thread needs to be told to
1392 wait for that. */
1393 if (PAR(cr) && MASTER(cr))
1395 done_commrec(cr);
1396 tMPI_Finalize();
1398 #endif
1400 return rc;
1403 } // namespace gmx