Free more memory in grompp and mdrun
[gromacs.git] / src / programs / mdrun / runner.cpp
blob306b21777e267e1d6908f3c2ee7abed74fd583ce
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_mdlib
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/gpu_utils/gpu_utils.h"
67 #include "gromacs/hardware/cpuinfo.h"
68 #include "gromacs/hardware/detecthardware.h"
69 #include "gromacs/hardware/printhardware.h"
70 #include "gromacs/listed-forces/disre.h"
71 #include "gromacs/listed-forces/orires.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/calc_verletbuf.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/gmx_omp_nthreads.h"
80 #include "gromacs/mdlib/integrator.h"
81 #include "gromacs/mdlib/main.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/mdatoms.h"
84 #include "gromacs/mdlib/mdrun.h"
85 #include "gromacs/mdlib/minimize.h"
86 #include "gromacs/mdlib/nb_verlet.h"
87 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
88 #include "gromacs/mdlib/nbnxn_search.h"
89 #include "gromacs/mdlib/nbnxn_tuning.h"
90 #include "gromacs/mdlib/qmmm.h"
91 #include "gromacs/mdlib/sighandler.h"
92 #include "gromacs/mdlib/sim_util.h"
93 #include "gromacs/mdlib/tpi.h"
94 #include "gromacs/mdrunutility/mdmodules.h"
95 #include "gromacs/mdrunutility/threadaffinity.h"
96 #include "gromacs/mdtypes/commrec.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/md_enums.h"
99 #include "gromacs/mdtypes/observableshistory.h"
100 #include "gromacs/mdtypes/state.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/pulling/pull.h"
103 #include "gromacs/pulling/pull_rotation.h"
104 #include "gromacs/taskassignment/decidegpuusage.h"
105 #include "gromacs/taskassignment/resourcedivision.h"
106 #include "gromacs/taskassignment/taskassignment.h"
107 #include "gromacs/taskassignment/usergpuids.h"
108 #include "gromacs/timing/wallcycle.h"
109 #include "gromacs/topology/mtop_util.h"
110 #include "gromacs/trajectory/trajectoryframe.h"
111 #include "gromacs/utility/basenetwork.h"
112 #include "gromacs/utility/cstringutil.h"
113 #include "gromacs/utility/exceptions.h"
114 #include "gromacs/utility/fatalerror.h"
115 #include "gromacs/utility/filestream.h"
116 #include "gromacs/utility/gmxassert.h"
117 #include "gromacs/utility/gmxmpi.h"
118 #include "gromacs/utility/logger.h"
119 #include "gromacs/utility/loggerbuilder.h"
120 #include "gromacs/utility/physicalnodecommunicator.h"
121 #include "gromacs/utility/pleasecite.h"
122 #include "gromacs/utility/programcontext.h"
123 #include "gromacs/utility/smalloc.h"
124 #include "gromacs/utility/stringutil.h"
126 #include "deform.h"
127 #include "md.h"
128 #include "membed.h"
129 #include "repl_ex.h"
131 #ifdef GMX_FAHCORE
132 #include "corewrap.h"
133 #endif
135 //! First step used in pressure scaling
136 gmx_int64_t deform_init_init_step_tpx;
137 //! Initial box for pressure scaling
138 matrix deform_init_box_tpx;
139 //! MPI variable for use in pressure scaling
140 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
142 namespace gmx
145 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
147 * Used to ensure that the master thread does not modify mdrunner during copy
148 * on the spawned threads. */
149 static void threadMpiMdrunnerAccessBarrier()
151 #if GMX_THREAD_MPI
152 MPI_Barrier(MPI_COMM_WORLD);
153 #endif
156 void Mdrunner::reinitializeOnSpawnedThread()
158 threadMpiMdrunnerAccessBarrier();
160 cr = reinitialize_commrec_for_this_thread(cr);
162 GMX_RELEASE_ASSERT(!MASTER(cr), "reinitializeOnSpawnedThread should only be called on spawned threads");
164 // Only the master rank writes to the log file
165 fplog = nullptr;
168 /*! \brief The callback used for running on spawned threads.
170 * Obtains the pointer to the master mdrunner object from the one
171 * argument permitted to the thread-launch API call, copies it to make
172 * a new runner for this thread, reinitializes necessary data, and
173 * proceeds to the simulation. */
174 static void mdrunner_start_fn(const void *arg)
178 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
179 /* copy the arg list to make sure that it's thread-local. This
180 doesn't copy pointed-to items, of course; fnm, cr and fplog
181 are reset in the call below, all others should be const. */
182 gmx::Mdrunner mdrunner = *masterMdrunner;
183 mdrunner.reinitializeOnSpawnedThread();
184 mdrunner.mdrunner();
186 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
190 /*! \brief Start thread-MPI threads.
192 * Called by mdrunner() to start a specific number of threads
193 * (including the main thread) for thread-parallel runs. This in turn
194 * calls mdrunner() for each thread. All options are the same as for
195 * mdrunner(). */
196 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
199 /* first check whether we even need to start tMPI */
200 if (numThreadsToLaunch < 2)
202 return cr;
205 #if GMX_THREAD_MPI
206 /* now spawn new threads that start mdrunner_start_fn(), while
207 the main thread returns, we set thread affinity later */
208 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
209 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
211 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
214 threadMpiMdrunnerAccessBarrier();
215 #else
216 GMX_UNUSED_VALUE(mdrunner_start_fn);
217 #endif
219 return reinitialize_commrec_for_this_thread(cr);
222 } // namespace
224 /*! \brief Initialize variables for Verlet scheme simulation */
225 static void prepare_verlet_scheme(FILE *fplog,
226 t_commrec *cr,
227 t_inputrec *ir,
228 int nstlist_cmdline,
229 const gmx_mtop_t *mtop,
230 const matrix box,
231 bool makeGpuPairList,
232 const gmx::CpuInfo &cpuinfo)
234 /* For NVE simulations, we will retain the initial list buffer */
235 if (EI_DYNAMICS(ir->eI) &&
236 ir->verletbuf_tol > 0 &&
237 !(EI_MD(ir->eI) && ir->etc == etcNO))
239 /* Update the Verlet buffer size for the current run setup */
241 /* Here we assume SIMD-enabled kernels are being used. But as currently
242 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
243 * and 4x2 gives a larger buffer than 4x4, this is ok.
245 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
246 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
248 real rlist_new;
249 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &listSetup, nullptr, &rlist_new);
251 if (rlist_new != ir->rlist)
253 if (fplog != nullptr)
255 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
256 ir->rlist, rlist_new,
257 listSetup.cluster_size_i, listSetup.cluster_size_j);
259 ir->rlist = rlist_new;
263 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
265 gmx_fatal(FARGS, "Can not set nstlist without %s",
266 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
269 if (EI_DYNAMICS(ir->eI))
271 /* Set or try nstlist values */
272 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
276 /*! \brief Override the nslist value in inputrec
278 * with value passed on the command line (if any)
280 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
281 gmx_int64_t nsteps_cmdline,
282 t_inputrec *ir)
284 assert(ir);
286 /* override with anything else than the default -2 */
287 if (nsteps_cmdline > -2)
289 char sbuf_steps[STEPSTRSIZE];
290 char sbuf_msg[STRLEN];
292 ir->nsteps = nsteps_cmdline;
293 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
295 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
296 gmx_step_str(nsteps_cmdline, sbuf_steps),
297 fabs(nsteps_cmdline*ir->delta_t));
299 else
301 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
302 gmx_step_str(nsteps_cmdline, sbuf_steps));
305 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
307 else if (nsteps_cmdline < -2)
309 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
310 nsteps_cmdline);
312 /* Do nothing if nsteps_cmdline == -2 */
315 namespace gmx
318 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
320 * If not, and if a warning may be issued, logs a warning about
321 * falling back to CPU code. With thread-MPI, only the first
322 * call to this function should have \c issueWarning true. */
323 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
324 const t_inputrec *ir,
325 bool issueWarning)
327 if (ir->opts.ngener - ir->nwall > 1)
329 /* The GPU code does not support more than one energy group.
330 * If the user requested GPUs explicitly, a fatal error is given later.
332 if (issueWarning)
334 GMX_LOG(mdlog.warning).asParagraph()
335 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
336 "For better performance, run on the GPU without energy groups and then do "
337 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
339 return false;
341 return true;
344 //! \brief Return the correct integrator function.
345 static integrator_t *my_integrator(unsigned int ei)
347 switch (ei)
349 case eiMD:
350 case eiBD:
351 case eiSD1:
352 case eiVV:
353 case eiVVAK:
354 if (!EI_DYNAMICS(ei))
356 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
358 return do_md;
359 case eiSteep:
360 return do_steep;
361 case eiCG:
362 return do_cg;
363 case eiNM:
364 return do_nm;
365 case eiLBFGS:
366 return do_lbfgs;
367 case eiTPI:
368 case eiTPIC:
369 if (!EI_TPI(ei))
371 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
373 return do_tpi;
374 case eiSD2_REMOVED:
375 GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
376 default:
377 GMX_THROW(APIError("Non existing integrator selected"));
381 //! Initializes the logger for mdrun.
382 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
384 gmx::LoggerBuilder builder;
385 if (fplog != nullptr)
387 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
389 if (cr == nullptr || SIMMASTER(cr))
391 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
392 &gmx::TextOutputFile::standardError());
394 return builder.build();
397 //! Make a TaskTarget from an mdrun argument string.
398 static TaskTarget findTaskTarget(const char *optionString)
400 TaskTarget returnValue = TaskTarget::Auto;
402 if (strncmp(optionString, "auto", 3) == 0)
404 returnValue = TaskTarget::Auto;
406 else if (strncmp(optionString, "cpu", 3) == 0)
408 returnValue = TaskTarget::Cpu;
410 else if (strncmp(optionString, "gpu", 3) == 0)
412 returnValue = TaskTarget::Gpu;
414 else
416 GMX_ASSERT(false, "Option string should have been checked for sanity already");
419 return returnValue;
422 int Mdrunner::mdrunner()
424 matrix box;
425 gmx_ddbox_t ddbox = {0};
426 int npme_major, npme_minor;
427 t_nrnb *nrnb;
428 gmx_mtop_t *mtop = nullptr;
429 t_forcerec *fr = nullptr;
430 t_fcdata *fcd = nullptr;
431 real ewaldcoeff_q = 0;
432 real ewaldcoeff_lj = 0;
433 gmx_vsite_t *vsite = nullptr;
434 gmx_constr_t constr;
435 int nChargePerturbed = -1, nTypePerturbed = 0;
436 gmx_wallcycle_t wcycle;
437 gmx_walltime_accounting_t walltime_accounting = nullptr;
438 int rc;
439 gmx_int64_t reset_counters;
440 int nthreads_pme = 1;
441 gmx_membed_t * membed = nullptr;
442 gmx_hw_info_t *hwinfo = nullptr;
444 /* CAUTION: threads may be started later on in this function, so
445 cr doesn't reflect the final parallel state right now */
446 std::unique_ptr<gmx::MDModules> mdModules(new gmx::MDModules);
447 t_inputrec inputrecInstance;
448 t_inputrec *inputrec = &inputrecInstance;
449 snew(mtop, 1);
451 if (mdrunOptions.continuationOptions.appendFiles)
453 fplog = nullptr;
456 bool doMembed = opt2bSet("-membed", nfile, fnm);
457 bool doRerun = mdrunOptions.rerun;
459 // Handle task-assignment related user options.
460 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
461 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
462 std::vector<int> gpuIdsAvailable;
465 gpuIdsAvailable = parseUserGpuIds(hw_opt.gpuIdsAvailable);
466 // TODO We could put the GPU IDs into a std::map to find
467 // duplicates, but for the small numbers of IDs involved, this
468 // code is simple and fast.
469 for (size_t i = 0; i != gpuIdsAvailable.size(); ++i)
471 for (size_t j = i+1; j != gpuIdsAvailable.size(); ++j)
473 if (gpuIdsAvailable[i] == gpuIdsAvailable[j])
475 GMX_THROW(InvalidInputError(formatString("The string of available GPU device IDs '%s' may not contain duplicate device IDs", hw_opt.gpuIdsAvailable.c_str())));
480 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
482 std::vector<int> userGpuTaskAssignment;
485 userGpuTaskAssignment = parseUserGpuIds(hw_opt.userGpuTaskAssignment);
487 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
488 auto nonbondedTarget = findTaskTarget(nbpu_opt);
489 auto pmeTarget = findTaskTarget(pme_opt);
490 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
491 PmeRunMode pmeRunMode = PmeRunMode::None;
493 // Here we assume that SIMMASTER(cr) does not change even after the
494 // threads are started.
495 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
496 gmx::MDLogger mdlog(logOwner.logger());
498 // TODO The thread-MPI master rank makes a working
499 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
500 // after the threads have been launched. This works because no use
501 // is made of that communicator until after the execution paths
502 // have rejoined. But it is likely that we can improve the way
503 // this is expressed, e.g. by expressly running detection only the
504 // master rank for thread-MPI, rather than relying on the mutex
505 // and reference count.
506 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
507 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
509 gmx_print_detected_hardware(fplog, cr, ms, mdlog, hwinfo);
511 std::vector<int> gpuIdsToUse;
512 auto compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
513 if (gpuIdsAvailable.empty())
515 gpuIdsToUse = compatibleGpus;
517 else
519 for (const auto &availableGpuId : gpuIdsAvailable)
521 bool availableGpuIsCompatible = false;
522 for (const auto &compatibleGpuId : compatibleGpus)
524 if (availableGpuId == compatibleGpuId)
526 availableGpuIsCompatible = true;
527 break;
530 if (!availableGpuIsCompatible)
532 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);
534 gpuIdsToUse.push_back(availableGpuId);
538 if (fplog != nullptr)
540 /* Print references after all software/hardware printing */
541 please_cite(fplog, "Abraham2015");
542 please_cite(fplog, "Pall2015");
543 please_cite(fplog, "Pronk2013");
544 please_cite(fplog, "Hess2008b");
545 please_cite(fplog, "Spoel2005a");
546 please_cite(fplog, "Lindahl2001a");
547 please_cite(fplog, "Berendsen95a");
550 std::unique_ptr<t_state> globalState;
552 if (SIMMASTER(cr))
554 /* Only the master rank has the global state */
555 globalState = std::unique_ptr<t_state>(new t_state);
557 /* Read (nearly) all data required for the simulation */
558 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), mtop);
560 if (inputrec->cutoff_scheme != ecutsVERLET)
562 if (nstlist_cmdline > 0)
564 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
567 if (!compatibleGpus.empty())
569 GMX_LOG(mdlog.warning).asParagraph().appendText(
570 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
571 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
576 /* Check and update the hardware options for internal consistency */
577 check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
579 /* Early check for externally set process affinity. */
580 gmx_check_thread_affinity_set(mdlog, cr,
581 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
583 if (GMX_THREAD_MPI && SIMMASTER(cr))
585 if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
587 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
590 /* Since the master knows the cut-off scheme, update hw_opt for this.
591 * This is done later for normal MPI and also once more with tMPI
592 * for all tMPI ranks.
594 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
596 bool useGpuForNonbonded = false;
597 bool useGpuForPme = false;
600 // If the user specified the number of ranks, then we must
601 // respect that, but in default mode, we need to allow for
602 // the number of GPUs to choose the number of ranks.
604 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
605 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
606 inputrec->cutoff_scheme == ecutsVERLET,
607 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
608 hw_opt.nthreads_tmpi);
609 auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
610 auto canUseGpuForPme = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
611 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
612 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
613 canUseGpuForPme, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
616 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
618 /* Determine how many thread-MPI ranks to start.
620 * TODO Over-writing the user-supplied value here does
621 * prevent any possible subsequent checks from working
622 * correctly. */
623 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
624 &hw_opt,
625 gpuIdsToUse,
626 useGpuForNonbonded,
627 useGpuForPme,
628 inputrec, mtop,
629 mdlog,
630 doMembed);
632 // Now start the threads for thread MPI.
633 cr = spawnThreads(hw_opt.nthreads_tmpi);
634 /* The main thread continues here with a new cr. We don't deallocate
635 the old cr because other threads may still be reading it. */
636 // TODO Both master and spawned threads call dup_tfn and
637 // reinitialize_commrec_for_this_thread. Find a way to express
638 // this better.
639 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
641 // END OF CAUTION: cr and physicalNodeComm are now reliable
643 if (PAR(cr))
645 /* now broadcast everything to the non-master nodes/threads: */
646 init_parallel(cr, inputrec, mtop);
649 // Now each rank knows the inputrec that SIMMASTER read and used,
650 // and (if applicable) cr->nnodes has been assigned the number of
651 // thread-MPI ranks that have been chosen. The ranks can now all
652 // run the task-deciding functions and will agree on the result
653 // without needing to communicate.
655 // TODO Should we do the communication in debug mode to support
656 // having an assertion?
658 // Note that these variables describe only their own node.
659 bool useGpuForNonbonded = false;
660 bool useGpuForPme = false;
663 // It's possible that there are different numbers of GPUs on
664 // different nodes, which is the user's responsibilty to
665 // handle. If unsuitable, we will notice that during task
666 // assignment.
667 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
668 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
669 emulateGpuNonbonded, inputrec->cutoff_scheme == ecutsVERLET,
670 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
671 gpusWereDetected);
672 auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
673 auto canUseGpuForPme = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
674 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
675 canUseGpuForPme, cr->nnodes, domdecOptions.numPmeRanks,
676 gpusWereDetected);
678 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
679 if (pmeRunMode == PmeRunMode::GPU)
681 if (pmeFftTarget == TaskTarget::Cpu)
683 pmeRunMode = PmeRunMode::Mixed;
686 else if (pmeFftTarget == TaskTarget::Gpu)
688 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.");
691 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
693 // TODO: Error handling
694 mdModules->assignOptionsToModules(*inputrec->params, nullptr);
696 if (fplog != nullptr)
698 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
699 fprintf(fplog, "\n");
702 if (SIMMASTER(cr))
704 /* now make sure the state is initialized and propagated */
705 set_state_entries(globalState.get(), inputrec);
708 /* NM and TPI parallelize over force/energy calculations, not atoms,
709 * so we need to initialize and broadcast the global state.
711 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
713 if (!MASTER(cr))
715 globalState = std::unique_ptr<t_state>(new t_state);
717 broadcastStateWithoutDynamics(cr, globalState.get());
720 /* A parallel command line option consistency check that we can
721 only do after any threads have started. */
722 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
723 domdecOptions.numCells[YY] > 1 ||
724 domdecOptions.numCells[ZZ] > 1 ||
725 domdecOptions.numPmeRanks > 0))
727 gmx_fatal(FARGS,
728 "The -dd or -npme option request a parallel simulation, "
729 #if !GMX_MPI
730 "but %s was compiled without threads or MPI enabled"
731 #else
732 #if GMX_THREAD_MPI
733 "but the number of MPI-threads (option -ntmpi) is not set or is 1"
734 #else
735 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
736 #endif
737 #endif
738 , output_env_get_program_display_name(oenv)
742 if (doRerun &&
743 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
745 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
748 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
750 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
753 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
755 if (domdecOptions.numPmeRanks > 0)
757 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
758 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
761 domdecOptions.numPmeRanks = 0;
764 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
766 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
767 * improve performance with many threads per GPU, since our OpenMP
768 * scaling is bad, but it's difficult to automate the setup.
770 domdecOptions.numPmeRanks = 0;
772 if (useGpuForPme)
774 if (domdecOptions.numPmeRanks < 0)
776 domdecOptions.numPmeRanks = 0;
777 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
779 else
781 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
785 #ifdef GMX_FAHCORE
786 if (MASTER(cr))
788 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
790 #endif
792 /* NMR restraints must be initialized before load_checkpoint,
793 * since with time averaging the history is added to t_state.
794 * For proper consistency check we therefore need to extend
795 * t_state here.
796 * So the PME-only nodes (if present) will also initialize
797 * the distance restraints.
799 snew(fcd, 1);
801 /* This needs to be called before read_checkpoint to extend the state */
802 init_disres(fplog, mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
804 init_orires(fplog, mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
806 if (inputrecDeform(inputrec))
808 /* Store the deform reference box before reading the checkpoint */
809 if (SIMMASTER(cr))
811 copy_mat(globalState->box, box);
813 if (PAR(cr))
815 gmx_bcast(sizeof(box), box, cr);
817 /* Because we do not have the update struct available yet
818 * in which the reference values should be stored,
819 * we store them temporarily in static variables.
820 * This should be thread safe, since they are only written once
821 * and with identical values.
823 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
824 deform_init_init_step_tpx = inputrec->init_step;
825 copy_mat(box, deform_init_box_tpx);
826 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
829 ObservablesHistory observablesHistory = {};
831 ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
833 if (continuationOptions.startedFromCheckpoint)
835 /* Check if checkpoint file exists before doing continuation.
836 * This way we can use identical input options for the first and subsequent runs...
838 gmx_bool bReadEkin;
840 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
841 cr, domdecOptions.numCells,
842 inputrec, globalState.get(),
843 &bReadEkin, &observablesHistory,
844 continuationOptions.appendFiles,
845 continuationOptions.appendFilesOptionSet,
846 mdrunOptions.reproducible);
848 if (bReadEkin)
850 continuationOptions.haveReadEkin = true;
854 if (SIMMASTER(cr) && continuationOptions.appendFiles)
856 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
857 continuationOptions.appendFiles, &fplog);
858 logOwner = buildLogger(fplog, nullptr);
859 mdlog = logOwner.logger();
862 if (mdrunOptions.numStepsCommandline > -2)
864 GMX_LOG(mdlog.info).asParagraph().
865 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
866 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
868 /* override nsteps with value set on the commamdline */
869 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
871 if (SIMMASTER(cr))
873 copy_mat(globalState->box, box);
876 if (PAR(cr))
878 gmx_bcast(sizeof(box), box, cr);
881 /* Update rlist and nstlist. */
882 if (inputrec->cutoff_scheme == ecutsVERLET)
884 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, mtop, box,
885 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
888 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
889 inputrec->eI == eiNM))
891 const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
893 cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
894 mtop, inputrec,
895 box, xOnMaster,
896 &ddbox, &npme_major, &npme_minor);
897 // Note that local state still does not exist yet.
899 else
901 /* PME, if used, is done on all nodes with 1D decomposition */
902 cr->npmenodes = 0;
903 cr->duty = (DUTY_PP | DUTY_PME);
904 npme_major = 1;
905 npme_minor = 1;
907 if (inputrec->ePBC == epbcSCREW)
909 gmx_fatal(FARGS,
910 "pbc=%s is only implemented with domain decomposition",
911 epbc_names[inputrec->ePBC]);
915 if (PAR(cr))
917 /* After possible communicator splitting in make_dd_communicators.
918 * we can set up the intra/inter node communication.
920 gmx_setup_nodecomm(fplog, cr);
923 #if GMX_MPI
924 if (isMultiSim(ms))
926 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
927 "This is simulation %d out of %d running as a composite GROMACS\n"
928 "multi-simulation job. Setup for this simulation:\n",
929 ms->sim, ms->nsim);
931 GMX_LOG(mdlog.warning).appendTextFormatted(
932 "Using %d MPI %s\n",
933 cr->nnodes,
934 #if GMX_THREAD_MPI
935 cr->nnodes == 1 ? "thread" : "threads"
936 #else
937 cr->nnodes == 1 ? "process" : "processes"
938 #endif
940 fflush(stderr);
941 #endif
943 /* Check and update hw_opt for the cut-off scheme */
944 check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
946 /* Check and update the number of OpenMP threads requested */
947 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
948 pmeRunMode, *mtop);
950 gmx_omp_nthreads_init(mdlog, cr,
951 hwinfo->nthreads_hw_avail,
952 physicalNodeComm.size_,
953 hw_opt.nthreads_omp,
954 hw_opt.nthreads_omp_pme,
955 !thisRankHasDuty(cr, DUTY_PP),
956 inputrec->cutoff_scheme == ecutsVERLET);
958 #ifndef NDEBUG
959 if (EI_TPI(inputrec->eI) &&
960 inputrec->cutoff_scheme == ecutsVERLET)
962 gmx_feenableexcept();
964 #endif
966 // Build a data structure that expresses which kinds of non-bonded
967 // task are handled by this rank.
969 // TODO Later, this might become a loop over all registered modules
970 // relevant to the mdp inputs, to find those that have such tasks.
972 // TODO This could move before init_domain_decomposition() as part
973 // of refactoring that separates the responsibility for duty
974 // assignment from setup for communication between tasks, and
975 // setup for tasks handled with a domain (ie including short-ranged
976 // tasks, bonded tasks, etc.).
978 // Note that in general useGpuForNonbonded, etc. can have a value
979 // that is inconsistent with the presence of actual GPUs on any
980 // rank, and that is not known to be a problem until the
981 // duty of the ranks on a node become node.
983 // TODO Later we might need the concept of computeTasksOnThisRank,
984 // from which we construct gpuTasksOnThisRank.
986 // Currently the DD code assigns duty to ranks that can
987 // include PP work that currently can be executed on a single
988 // GPU, if present and compatible. This has to be coordinated
989 // across PP ranks on a node, with possible multiple devices
990 // or sharing devices on a node, either from the user
991 // selection, or automatically.
992 auto haveGpus = !gpuIdsToUse.empty();
993 std::vector<GpuTask> gpuTasksOnThisRank;
994 if (thisRankHasDuty(cr, DUTY_PP))
996 if (useGpuForNonbonded)
998 if (haveGpus)
1000 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1002 else if (nonbondedTarget == TaskTarget::Gpu)
1004 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1008 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1009 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1011 if (useGpuForPme)
1013 if (haveGpus)
1015 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1017 else if (pmeTarget == TaskTarget::Gpu)
1019 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1024 GpuTaskAssignment gpuTaskAssignment;
1027 // Produce the task assignment for this rank.
1028 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1029 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank);
1031 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1033 /* Prevent other ranks from continuing after an issue was found
1034 * and reported as a fatal error.
1036 * TODO This function implements a barrier so that MPI runtimes
1037 * can organize an orderly shutdown if one of the ranks has had to
1038 * issue a fatal error in various code already run. When we have
1039 * MPI-aware error handling and reporting, this should be
1040 * improved. */
1041 #if GMX_MPI
1042 if (PAR(cr))
1044 MPI_Barrier(cr->mpi_comm_mysim);
1046 if (isMultiSim(ms))
1048 if (SIMMASTER(cr))
1050 MPI_Barrier(ms->mpi_comm_masters);
1052 /* We need another barrier to prevent non-master ranks from contiuing
1053 * when an error occured in a different simulation.
1055 MPI_Barrier(cr->mpi_comm_mysim);
1057 #endif
1059 /* Now that we know the setup is consistent, check for efficiency */
1060 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1061 cr, mdlog);
1063 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1065 if (thisRankHasDuty(cr, DUTY_PP))
1067 // This works because only one task of each type is currently permitted.
1068 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1069 hasTaskType<GpuTask::Nonbonded>);
1070 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1072 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1073 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1074 init_gpu(mdlog, nonbondedDeviceInfo);
1076 if (DOMAINDECOMP(cr))
1078 /* When we share GPUs over ranks, we need to know this for the DLB */
1079 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1085 gmx_device_info_t *pmeDeviceInfo = nullptr;
1086 // This works because only one task of each type is currently permitted.
1087 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1088 if (pmeGpuTaskMapping != gpuTaskAssignment.end())
1090 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1091 init_gpu(mdlog, pmeDeviceInfo);
1094 /* getting number of PP/PME threads
1095 PME: env variable should be read only on one node to make sure it is
1096 identical everywhere;
1098 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1100 int numThreadsOnThisRank;
1101 /* threads on this MPI process or TMPI thread */
1102 if (thisRankHasDuty(cr, DUTY_PP))
1104 numThreadsOnThisRank = gmx_omp_nthreads_get(emntNonbonded);
1106 else
1108 numThreadsOnThisRank = nthreads_pme;
1111 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1112 *hwinfo->hardwareTopology,
1113 physicalNodeComm, mdlog);
1115 if (hw_opt.thread_affinity != threadaffOFF)
1117 /* Before setting affinity, check whether the affinity has changed
1118 * - which indicates that probably the OpenMP library has changed it
1119 * since we first checked).
1121 gmx_check_thread_affinity_set(mdlog, cr,
1122 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1124 int numThreadsOnThisNode, intraNodeThreadOffset;
1125 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1126 &intraNodeThreadOffset);
1128 /* Set the CPU affinity */
1129 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1130 numThreadsOnThisRank, numThreadsOnThisNode,
1131 intraNodeThreadOffset, nullptr);
1134 if (mdrunOptions.timingOptions.resetStep > -1)
1136 GMX_LOG(mdlog.info).asParagraph().
1137 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1139 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1141 if (PAR(cr))
1143 /* Master synchronizes its value of reset_counters with all nodes
1144 * including PME only nodes */
1145 reset_counters = wcycle_get_reset_counters(wcycle);
1146 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1147 wcycle_set_reset_counters(wcycle, reset_counters);
1150 // Membrane embedding must be initialized before we call init_forcerec()
1151 if (doMembed)
1153 if (MASTER(cr))
1155 fprintf(stderr, "Initializing membed");
1157 /* Note that membed cannot work in parallel because mtop is
1158 * changed here. Fix this if we ever want to make it run with
1159 * multiple ranks. */
1160 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1163 std::unique_ptr<MDAtoms> mdAtoms;
1165 snew(nrnb, 1);
1166 if (thisRankHasDuty(cr, DUTY_PP))
1168 /* Initiate forcerecord */
1169 fr = mk_forcerec();
1170 fr->forceProviders = mdModules->initForceProviders();
1171 init_forcerec(fplog, mdlog, fr, fcd,
1172 inputrec, mtop, cr, box,
1173 opt2fn("-table", nfile, fnm),
1174 opt2fn("-tablep", nfile, fnm),
1175 opt2fns("-tableb", nfile, fnm),
1176 *hwinfo, nonbondedDeviceInfo,
1177 FALSE,
1178 pforce);
1180 /* Initialize QM-MM */
1181 if (fr->bQMMM)
1183 GMX_LOG(mdlog.info).asParagraph().
1184 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
1185 "version. Please get in touch with the developers if you find the support useful, "
1186 "as help is needed if the functionality is to continue to be available.");
1187 init_QMMMrec(cr, mtop, inputrec, fr);
1190 /* Initialize the mdAtoms structure.
1191 * mdAtoms is not filled with atom data,
1192 * as this can not be done now with domain decomposition.
1194 const bool useGpuForPme = (pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Mixed);
1195 mdAtoms = makeMDAtoms(fplog, *mtop, *inputrec, useGpuForPme && thisRankHasDuty(cr, DUTY_PME));
1196 if (globalState)
1198 // The pinning of coordinates in the global state object works, because we only use
1199 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1200 // points to the global state object without DD.
1201 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1202 // which should also perform the pinning.
1203 changePinningPolicy(&globalState->x, useGpuForPme ? PinningPolicy::CanBePinned : PinningPolicy::CannotBePinned);
1206 /* Initialize the virtual site communication */
1207 vsite = initVsite(*mtop, cr);
1209 calc_shifts(box, fr->shift_vec);
1211 /* With periodic molecules the charge groups should be whole at start up
1212 * and the virtual sites should not be far from their proper positions.
1214 if (!inputrec->bContinuation && MASTER(cr) &&
1215 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1217 /* Make molecules whole at start of run */
1218 if (fr->ePBC != epbcNONE)
1220 rvec *xGlobal = as_rvec_array(globalState->x.data());
1221 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, xGlobal);
1223 if (vsite)
1225 /* Correct initial vsite positions are required
1226 * for the initial distribution in the domain decomposition
1227 * and for the initial shell prediction.
1229 constructVsitesGlobal(*mtop, globalState->x);
1233 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1235 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1236 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1239 else
1241 /* This is a PME only node */
1243 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1245 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1246 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1249 gmx_pme_t *sepPmeData = nullptr;
1250 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1251 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1252 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1254 /* Initiate PME if necessary,
1255 * either on all nodes or on dedicated PME nodes only. */
1256 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1258 if (mdAtoms && mdAtoms->mdatoms())
1260 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1261 if (EVDW_PME(inputrec->vdwtype))
1263 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1266 if (cr->npmenodes > 0)
1268 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1269 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1270 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1273 if (thisRankHasDuty(cr, DUTY_PME))
1277 pmedata = gmx_pme_init(cr, npme_major, npme_minor, inputrec,
1278 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1279 mdrunOptions.reproducible,
1280 ewaldcoeff_q, ewaldcoeff_lj,
1281 nthreads_pme,
1282 pmeRunMode, nullptr, pmeDeviceInfo, mdlog);
1284 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1289 if (EI_DYNAMICS(inputrec->eI))
1291 /* Turn on signal handling on all nodes */
1293 * (A user signal from the PME nodes (if any)
1294 * is communicated to the PP nodes.
1296 signal_handler_install();
1299 if (thisRankHasDuty(cr, DUTY_PP))
1301 /* Assumes uniform use of the number of OpenMP threads */
1302 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1304 if (inputrec->bPull)
1306 /* Initialize pull code */
1307 inputrec->pull_work =
1308 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1309 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1310 EI_DYNAMICS(inputrec->eI) && MASTER(cr),
1311 continuationOptions);
1314 if (inputrec->bRot)
1316 /* Initialize enforced rotation code */
1317 init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), mtop, oenv, mdrunOptions);
1320 /* Let init_constraints know whether we have essential dynamics constraints.
1321 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1323 bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1325 constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1327 if (DOMAINDECOMP(cr))
1329 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1330 /* This call is not included in init_domain_decomposition mainly
1331 * because fr->cginfo_mb is set later.
1333 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1334 domdecOptions.checkBondedInteractions,
1335 fr->cginfo_mb);
1338 /* Now do whatever the user wants us to do (how flexible...) */
1339 my_integrator(inputrec->eI) (fplog, cr, ms, mdlog, nfile, fnm,
1340 oenv,
1341 mdrunOptions,
1342 vsite, constr,
1343 mdModules->outputProvider(),
1344 inputrec, mtop,
1345 fcd,
1346 globalState.get(),
1347 &observablesHistory,
1348 mdAtoms.get(), nrnb, wcycle, fr,
1349 replExParams,
1350 membed,
1351 walltime_accounting);
1353 if (inputrec->bRot)
1355 finish_rot(inputrec->rot);
1358 if (inputrec->bPull)
1360 finish_pull(inputrec->pull_work);
1364 else
1366 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1367 /* do PME only */
1368 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1369 gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1372 wallcycle_stop(wcycle, ewcRUN);
1374 /* Finish up, write some stuff
1375 * if rerunMD, don't write last frame again
1377 finish_run(fplog, mdlog, cr,
1378 inputrec, nrnb, wcycle, walltime_accounting,
1379 fr ? fr->nbv : nullptr,
1380 pmedata,
1381 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1383 // Free PME data
1384 if (pmedata)
1386 gmx_pme_destroy(pmedata);
1387 pmedata = nullptr;
1390 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1391 // before we destroy the GPU context(s) in free_gpu_resources().
1392 // Pinned buffers are associated with contexts in CUDA.
1393 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1394 mdAtoms.reset(nullptr);
1395 globalState.reset(nullptr);
1396 mdModules.reset(nullptr); // destruct force providers here as they might also use the GPU
1398 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1399 free_gpu_resources(fr, physicalNodeComm);
1400 free_gpu(nonbondedDeviceInfo);
1401 free_gpu(pmeDeviceInfo);
1402 done_forcerec(fr, mtop->nmolblock, mtop->groups.grps[egcENER].nr);
1403 done_mtop(mtop);
1404 sfree(fcd);
1406 if (doMembed)
1408 free_membed(membed);
1411 gmx_hardware_info_free();
1413 /* Does what it says */
1414 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1415 walltime_accounting_destroy(walltime_accounting);
1416 sfree(nrnb);
1418 /* Close logfile already here if we were appending to it */
1419 if (MASTER(cr) && continuationOptions.appendFiles)
1421 gmx_log_close(fplog);
1422 fplog = nullptr;
1425 rc = (int)gmx_get_stop_condition();
1427 #if GMX_THREAD_MPI
1428 /* we need to join all threads. The sub-threads join when they
1429 exit this function, but the master thread needs to be told to
1430 wait for that. */
1431 if (PAR(cr) && MASTER(cr))
1433 done_commrec(cr);
1434 tMPI_Finalize();
1436 #endif
1438 return rc;
1441 } // namespace gmx