53c55ff09fe9d25ef69b513ab9e31902911f42d8
[alexxy/gromacs.git] / src / programs / mdrun / runner.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \internal \file
38  *
39  * \brief Implements the MD runner routine calling all integrators.
40  *
41  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42  * \ingroup module_mdlib
43  */
44 #include "gmxpre.h"
45
46 #include "runner.h"
47
48 #include "config.h"
49
50 #include <cassert>
51 #include <csignal>
52 #include <cstdlib>
53 #include <cstring>
54
55 #include <algorithm>
56
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/ewald/ewald-utils.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/checkpoint.h"
63 #include "gromacs/fileio/oenv.h"
64 #include "gromacs/fileio/tpxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gpu_utils/gpu_utils.h"
67 #include "gromacs/hardware/cpuinfo.h"
68 #include "gromacs/hardware/detecthardware.h"
69 #include "gromacs/hardware/printhardware.h"
70 #include "gromacs/listed-forces/disre.h"
71 #include "gromacs/listed-forces/orires.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/utilities.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/calc_verletbuf.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/gmx_omp_nthreads.h"
80 #include "gromacs/mdlib/integrator.h"
81 #include "gromacs/mdlib/main.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/mdatoms.h"
84 #include "gromacs/mdlib/mdrun.h"
85 #include "gromacs/mdlib/minimize.h"
86 #include "gromacs/mdlib/nb_verlet.h"
87 #include "gromacs/mdlib/nbnxn_search.h"
88 #include "gromacs/mdlib/nbnxn_tuning.h"
89 #include "gromacs/mdlib/qmmm.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/sim_util.h"
92 #include "gromacs/mdlib/tpi.h"
93 #include "gromacs/mdrunutility/mdmodules.h"
94 #include "gromacs/mdrunutility/threadaffinity.h"
95 #include "gromacs/mdtypes/commrec.h"
96 #include "gromacs/mdtypes/inputrec.h"
97 #include "gromacs/mdtypes/md_enums.h"
98 #include "gromacs/mdtypes/observableshistory.h"
99 #include "gromacs/mdtypes/state.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/pulling/pull_rotation.h"
103 #include "gromacs/taskassignment/decidegpuusage.h"
104 #include "gromacs/taskassignment/resourcedivision.h"
105 #include "gromacs/taskassignment/taskassignment.h"
106 #include "gromacs/taskassignment/usergpuids.h"
107 #include "gromacs/timing/wallcycle.h"
108 #include "gromacs/topology/mtop_util.h"
109 #include "gromacs/trajectory/trajectoryframe.h"
110 #include "gromacs/utility/cstringutil.h"
111 #include "gromacs/utility/exceptions.h"
112 #include "gromacs/utility/fatalerror.h"
113 #include "gromacs/utility/filestream.h"
114 #include "gromacs/utility/gmxassert.h"
115 #include "gromacs/utility/gmxmpi.h"
116 #include "gromacs/utility/logger.h"
117 #include "gromacs/utility/loggerbuilder.h"
118 #include "gromacs/utility/pleasecite.h"
119 #include "gromacs/utility/programcontext.h"
120 #include "gromacs/utility/smalloc.h"
121 #include "gromacs/utility/stringutil.h"
122
123 #include "deform.h"
124 #include "md.h"
125 #include "membed.h"
126 #include "repl_ex.h"
127
128 #ifdef GMX_FAHCORE
129 #include "corewrap.h"
130 #endif
131
132 //! First step used in pressure scaling
133 gmx_int64_t         deform_init_init_step_tpx;
134 //! Initial box for pressure scaling
135 matrix              deform_init_box_tpx;
136 //! MPI variable for use in pressure scaling
137 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
138
139 namespace gmx
140 {
141
142 void Mdrunner::reinitializeOnSpawnedThread()
143 {
144     // TODO This duplication is formally necessary if any thread might
145     // modify any memory in fnm or the pointers it contains. If the
146     // contents are ever provably const, then we can remove this
147     // allocation (and memory leak).
148     // TODO This should probably become part of a copy constructor for
149     // Mdrunner.
150     fnm = dup_tfn(nfile, fnm);
151
152     cr  = reinitialize_commrec_for_this_thread(cr);
153
154     if (!MASTER(cr))
155     {
156         // Only the master rank writes to the log files
157         fplog = nullptr;
158     }
159 }
160
161 /*! \brief The callback used for running on spawned threads.
162  *
163  * Obtains the pointer to the master mdrunner object from the one
164  * argument permitted to the thread-launch API call, copies it to make
165  * a new runner for this thread, reinitializes necessary data, and
166  * proceeds to the simulation. */
167 static void mdrunner_start_fn(void *arg)
168 {
169     try
170     {
171         auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
172         /* copy the arg list to make sure that it's thread-local. This
173            doesn't copy pointed-to items, of course, but those are all
174            const. */
175         gmx::Mdrunner mdrunner = *masterMdrunner;
176         mdrunner.reinitializeOnSpawnedThread();
177         mdrunner.mdrunner();
178     }
179     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
180 }
181
182
183 /*! \brief Start thread-MPI threads.
184  *
185  * Called by mdrunner() to start a specific number of threads
186  * (including the main thread) for thread-parallel runs. This in turn
187  * calls mdrunner() for each thread. All options are the same as for
188  * mdrunner(). */
189 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch)
190 {
191
192     /* first check whether we even need to start tMPI */
193     if (numThreadsToLaunch < 2)
194     {
195         return cr;
196     }
197
198     gmx::Mdrunner spawnedMdrunner = *this;
199     // TODO This duplication is formally necessary if any thread might
200     // modify any memory in fnm or the pointers it contains. If the
201     // contents are ever provably const, then we can remove this
202     // allocation (and memory leak).
203     // TODO This should probably become part of a copy constructor for
204     // Mdrunner.
205     spawnedMdrunner.fnm = dup_tfn(this->nfile, fnm);
206
207 #if GMX_THREAD_MPI
208     /* now spawn new threads that start mdrunner_start_fn(), while
209        the main thread returns, we set thread affinity later */
210     if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
211                      mdrunner_start_fn, static_cast<void*>(&spawnedMdrunner)) != TMPI_SUCCESS)
212     {
213         GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
214     }
215 #else
216     GMX_UNUSED_VALUE(mdrunner_start_fn);
217 #endif
218
219     return reinitialize_commrec_for_this_thread(cr);
220 }
221
222 }      // namespace
223
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)
233 {
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))
238     {
239         /* Update the Verlet buffer size for the current run setup */
240
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.
244          */
245         ListSetupType      listType  = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
246         VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
247
248         real               rlist_new;
249         calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1, -1, &listSetup, nullptr, &rlist_new);
250
251         if (rlist_new != ir->rlist)
252         {
253             if (fplog != nullptr)
254             {
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);
258             }
259             ir->rlist     = rlist_new;
260         }
261     }
262
263     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
264     {
265         gmx_fatal(FARGS, "Can not set nstlist without %s",
266                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
267     }
268
269     if (EI_DYNAMICS(ir->eI))
270     {
271         /* Set or try nstlist values */
272         increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
273     }
274 }
275
276 /*! \brief Override the nslist value in inputrec
277  *
278  * with value passed on the command line (if any)
279  */
280 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
281                                     gmx_int64_t          nsteps_cmdline,
282                                     t_inputrec          *ir)
283 {
284     assert(ir);
285
286     /* override with anything else than the default -2 */
287     if (nsteps_cmdline > -2)
288     {
289         char sbuf_steps[STEPSTRSIZE];
290         char sbuf_msg[STRLEN];
291
292         ir->nsteps = nsteps_cmdline;
293         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
294         {
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));
298         }
299         else
300         {
301             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
302                     gmx_step_str(nsteps_cmdline, sbuf_steps));
303         }
304
305         GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
306     }
307     else if (nsteps_cmdline < -2)
308     {
309         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
310                   nsteps_cmdline);
311     }
312     /* Do nothing if nsteps_cmdline == -2 */
313 }
314
315 namespace gmx
316 {
317
318 /*! \brief Return whether GPU acceleration of nonbondeds is useful with the given settings.
319  *
320  * If not, logs a message about falling back to CPU code. */
321 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger   &mdlog,
322                                                const t_inputrec *ir,
323                                                bool              doRerun)
324 {
325     if (doRerun && ir->opts.ngener > 1)
326     {
327         /* Rerun execution time is dominated by I/O and pair search,
328          * so GPUs are not very useful, plus they do not support more
329          * than one energy group. If the user requested GPUs
330          * explicitly, a fatal error is given later.  With non-reruns,
331          * we fall back to a single whole-of system energy group
332          * (which runs much faster than a multiple-energy-groups
333          * implementation would), and issue a note in the .log
334          * file. Users can re-run if they want the information. */
335         GMX_LOG(mdlog.warning).asParagraph().appendText("Multiple energy groups is not implemented for GPUs, so is not useful for this rerun, so falling back to the CPU");
336         return false;
337     }
338
339     return true;
340 }
341
342 //! \brief Return the correct integrator function.
343 static integrator_t *my_integrator(unsigned int ei)
344 {
345     switch (ei)
346     {
347         case eiMD:
348         case eiBD:
349         case eiSD1:
350         case eiVV:
351         case eiVVAK:
352             if (!EI_DYNAMICS(ei))
353             {
354                 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
355             }
356             return do_md;
357         case eiSteep:
358             return do_steep;
359         case eiCG:
360             return do_cg;
361         case eiNM:
362             return do_nm;
363         case eiLBFGS:
364             return do_lbfgs;
365         case eiTPI:
366         case eiTPIC:
367             if (!EI_TPI(ei))
368             {
369                 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
370             }
371             return do_tpi;
372         case eiSD2_REMOVED:
373             GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
374         default:
375             GMX_THROW(APIError("Non existing integrator selected"));
376     }
377 }
378
379 //! Initializes the logger for mdrun.
380 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
381 {
382     gmx::LoggerBuilder builder;
383     if (fplog != nullptr)
384     {
385         builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
386     }
387     if (cr == nullptr || SIMMASTER(cr))
388     {
389         builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
390                                 &gmx::TextOutputFile::standardError());
391     }
392     return builder.build();
393 }
394
395 //! Make a TaskTarget from an mdrun argument string.
396 static TaskTarget findTaskTarget(const char *optionString)
397 {
398     TaskTarget returnValue = TaskTarget::Auto;
399
400     if (strncmp(optionString, "auto", 3) == 0)
401     {
402         returnValue = TaskTarget::Auto;
403     }
404     else if (strncmp(optionString, "cpu", 3) == 0)
405     {
406         returnValue = TaskTarget::Cpu;
407     }
408     else if (strncmp(optionString, "gpu", 3) == 0)
409     {
410         returnValue = TaskTarget::Gpu;
411     }
412     else
413     {
414         GMX_ASSERT(false, "Option string should have been checked for sanity already");
415     }
416
417     return returnValue;
418 }
419
420 int Mdrunner::mdrunner()
421 {
422     matrix                    box;
423     gmx_ddbox_t               ddbox = {0};
424     int                       npme_major, npme_minor;
425     t_nrnb                   *nrnb;
426     gmx_mtop_t               *mtop          = nullptr;
427     t_forcerec               *fr            = nullptr;
428     t_fcdata                 *fcd           = nullptr;
429     real                      ewaldcoeff_q  = 0;
430     real                      ewaldcoeff_lj = 0;
431     gmx_vsite_t              *vsite         = nullptr;
432     gmx_constr_t              constr;
433     int                       nChargePerturbed = -1, nTypePerturbed = 0;
434     gmx_wallcycle_t           wcycle;
435     gmx_walltime_accounting_t walltime_accounting = nullptr;
436     int                       rc;
437     gmx_int64_t               reset_counters;
438     int                       nthreads_pme = 1;
439     gmx_membed_t *            membed       = nullptr;
440     gmx_hw_info_t            *hwinfo       = nullptr;
441
442     /* CAUTION: threads may be started later on in this function, so
443        cr doesn't reflect the final parallel state right now */
444     gmx::MDModules mdModules;
445     t_inputrec     inputrecInstance;
446     t_inputrec    *inputrec = &inputrecInstance;
447     snew(mtop, 1);
448
449     if (mdrunOptions.continuationOptions.appendFiles)
450     {
451         fplog = nullptr;
452     }
453
454     bool doMembed = opt2bSet("-membed", nfile, fnm);
455     bool doRerun  = mdrunOptions.rerun;
456
457     // Handle task-assignment related user options.
458     EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
459                                                EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
460     std::vector<int>    gpuIdsAvailable;
461     try
462     {
463         gpuIdsAvailable = parseUserGpuIds(hw_opt.gpuIdsAvailable);
464         // TODO We could put the GPU IDs into a std::map to find
465         // duplicates, but for the small numbers of IDs involved, this
466         // code is simple and fast.
467         for (size_t i = 0; i != gpuIdsAvailable.size(); ++i)
468         {
469             for (size_t j = i+1; j != gpuIdsAvailable.size(); ++j)
470             {
471                 if (gpuIdsAvailable[i] == gpuIdsAvailable[j])
472                 {
473                     GMX_THROW(InvalidInputError(formatString("The string of available GPU device IDs '%s' may not contain duplicate device IDs", hw_opt.gpuIdsAvailable.c_str())));
474                 }
475             }
476         }
477     }
478     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
479
480     std::vector<int> userGpuTaskAssignment;
481     try
482     {
483         userGpuTaskAssignment = parseUserGpuIds(hw_opt.userGpuTaskAssignment);
484     }
485     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
486     auto             nonbondedTarget = findTaskTarget(nbpu_opt);
487     // TODO Connect these to actual mdrun arguments and some functionality
488     const char      *pme_opt     = "cpu";
489     auto             pmeTarget   = findTaskTarget(pme_opt);
490
491     // TODO find a sensible home and behaviour for this
492     //const char      *pme_fft_opt = "auto";
493     //auto pmeFftTarget    = findTaskTarget(pme_fft_opt);
494
495     const PmeRunMode pmeRunMode = PmeRunMode::CPU;
496     //TODO this is a placeholder as PME on GPU is not permitted yet
497     //TODO should there exist a PmeRunMode::None value for consistency?
498
499     // Here we assume that SIMMASTER(cr) does not change even after the
500     // threads are started.
501     gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
502     gmx::MDLogger    mdlog(logOwner.logger());
503
504     hwinfo = gmx_detect_hardware(mdlog, cr);
505
506     gmx_print_detected_hardware(fplog, cr, mdlog, hwinfo);
507
508     std::vector<int> gpuIdsToUse;
509     auto             compatibleGpus = getCompatibleGpus(hwinfo->gpu_info);
510     if (gpuIdsAvailable.empty())
511     {
512         gpuIdsToUse = compatibleGpus;
513     }
514     else
515     {
516         for (const auto &availableGpuId : gpuIdsAvailable)
517         {
518             bool availableGpuIsCompatible = false;
519             for (const auto &compatibleGpuId : compatibleGpus)
520             {
521                 if (availableGpuId == compatibleGpuId)
522                 {
523                     availableGpuIsCompatible = true;
524                     break;
525                 }
526             }
527             if (!availableGpuIsCompatible)
528             {
529                 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);
530             }
531             gpuIdsToUse.push_back(availableGpuId);
532         }
533     }
534
535     if (fplog != nullptr)
536     {
537         /* Print references after all software/hardware printing */
538         please_cite(fplog, "Abraham2015");
539         please_cite(fplog, "Pall2015");
540         please_cite(fplog, "Pronk2013");
541         please_cite(fplog, "Hess2008b");
542         please_cite(fplog, "Spoel2005a");
543         please_cite(fplog, "Lindahl2001a");
544         please_cite(fplog, "Berendsen95a");
545     }
546
547     std::unique_ptr<t_state> globalState;
548
549     if (SIMMASTER(cr))
550     {
551         /* Only the master rank has the global state */
552         globalState = std::unique_ptr<t_state>(new t_state);
553
554         /* Read (nearly) all data required for the simulation */
555         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, globalState.get(), mtop);
556
557         if (inputrec->cutoff_scheme != ecutsVERLET)
558         {
559             if (nstlist_cmdline > 0)
560             {
561                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
562             }
563
564             if (!compatibleGpus.empty())
565             {
566                 GMX_LOG(mdlog.warning).asParagraph().appendText(
567                         "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
568                         "      To use a GPU, set the mdp option: cutoff-scheme = Verlet");
569             }
570         }
571     }
572
573     /* Check and update the hardware options for internal consistency */
574     check_and_update_hw_opt_1(&hw_opt, cr, domdecOptions.numPmeRanks);
575
576     /* Early check for externally set process affinity. */
577     gmx_check_thread_affinity_set(mdlog, cr,
578                                   &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
579
580     if (GMX_THREAD_MPI && SIMMASTER(cr))
581     {
582         if (domdecOptions.numPmeRanks > 0 && hw_opt.nthreads_tmpi <= 0)
583         {
584             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
585         }
586
587         /* Since the master knows the cut-off scheme, update hw_opt for this.
588          * This is done later for normal MPI and also once more with tMPI
589          * for all tMPI ranks.
590          */
591         check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
592
593         bool useGpuForNonbonded = false;
594         bool useGpuForPme       = false;
595         try
596         {
597             // If the user specified the number of ranks, then we must
598             // respect that, but in default mode, we need to allow for
599             // the number of GPUs to choose the number of ranks.
600
601             useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
602                     (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
603                     inputrec->cutoff_scheme == ecutsVERLET,
604                     gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, doRerun),
605                     hw_opt.nthreads_tmpi);
606             auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
607             auto canUseGpuForPme   = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
608             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
609                     (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
610                     canUseGpuForPme, hw_opt.nthreads_tmpi);
611         }
612         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
613         /* Determine how many thread-MPI ranks to start.
614          *
615          * TODO Over-writing the user-supplied value here does
616          * prevent any possible subsequent checks from working
617          * correctly. */
618         hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
619                                                 &hw_opt,
620                                                 gpuIdsToUse,
621                                                 useGpuForNonbonded,
622                                                 useGpuForPme,
623                                                 inputrec, mtop,
624                                                 mdlog,
625                                                 doMembed);
626
627         // Now start the threads for thread MPI.
628         cr = spawnThreads(hw_opt.nthreads_tmpi);
629         /* The main thread continues here with a new cr. We don't deallocate
630            the old cr because other threads may still be reading it. */
631         // TODO Both master and spawned threads call dup_tfn and
632         // reinitialize_commrec_for_this_thread. Find a way to express
633         // this better.
634     }
635     /* END OF CAUTION: cr is now reliable */
636
637     if (PAR(cr))
638     {
639         /* now broadcast everything to the non-master nodes/threads: */
640         init_parallel(cr, inputrec, mtop);
641     }
642
643     // Now each rank knows the inputrec that SIMMASTER read and used,
644     // and (if applicable) cr->nnodes has been assigned the number of
645     // thread-MPI ranks that have been chosen. The ranks can now all
646     // run the task-deciding functions and will agree on the result
647     // without needing to communicate.
648     //
649     // TODO Should we do the communication in debug mode to support
650     // having an assertion?
651     //
652     // Note that these variables describe only their own node.
653     bool useGpuForNonbonded = false;
654     bool useGpuForPme       = false;
655     try
656     {
657         useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment,
658                                                                 emulateGpuNonbonded, inputrec->cutoff_scheme == ecutsVERLET,
659                                                                 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, doRerun));
660         auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
661         auto canUseGpuForPme   = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
662         useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, canUseGpuForPme, cr->nnodes);
663     }
664     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
665
666     // TODO: Error handling
667     mdModules.assignOptionsToModules(*inputrec->params, nullptr);
668
669     if (fplog != nullptr)
670     {
671         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
672         fprintf(fplog, "\n");
673     }
674
675     if (SIMMASTER(cr))
676     {
677         /* now make sure the state is initialized and propagated */
678         set_state_entries(globalState.get(), inputrec);
679     }
680
681     /* NM and TPI parallelize over force/energy calculations, not atoms,
682      * so we need to initialize and broadcast the global state.
683      */
684     if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
685     {
686         if (!MASTER(cr))
687         {
688             globalState = std::unique_ptr<t_state>(new t_state);
689         }
690         broadcastStateWithoutDynamics(cr, globalState.get());
691     }
692
693     /* A parallel command line option consistency check that we can
694        only do after any threads have started. */
695     if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
696                      domdecOptions.numCells[YY] > 1 ||
697                      domdecOptions.numCells[ZZ] > 1 ||
698                      domdecOptions.numPmeRanks > 0))
699     {
700         gmx_fatal(FARGS,
701                   "The -dd or -npme option request a parallel simulation, "
702 #if !GMX_MPI
703                   "but %s was compiled without threads or MPI enabled"
704 #else
705 #if GMX_THREAD_MPI
706                   "but the number of MPI-threads (option -ntmpi) is not set or is 1"
707 #else
708                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
709 #endif
710 #endif
711                   , output_env_get_program_display_name(oenv)
712                   );
713     }
714
715     if (doRerun &&
716         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
717     {
718         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
719     }
720
721     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
722     {
723         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
724     }
725
726     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
727     {
728         if (domdecOptions.numPmeRanks > 0)
729         {
730             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
731                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
732         }
733
734         domdecOptions.numPmeRanks = 0;
735     }
736
737     if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
738     {
739         /* With GPUs we don't automatically use PME-only ranks. PME ranks can
740          * improve performance with many threads per GPU, since our OpenMP
741          * scaling is bad, but it's difficult to automate the setup.
742          */
743         domdecOptions.numPmeRanks = 0;
744     }
745
746 #ifdef GMX_FAHCORE
747     if (MASTER(cr))
748     {
749         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
750     }
751 #endif
752
753     /* NMR restraints must be initialized before load_checkpoint,
754      * since with time averaging the history is added to t_state.
755      * For proper consistency check we therefore need to extend
756      * t_state here.
757      * So the PME-only nodes (if present) will also initialize
758      * the distance restraints.
759      */
760     snew(fcd, 1);
761
762     /* This needs to be called before read_checkpoint to extend the state */
763     init_disres(fplog, mtop, inputrec, cr, fcd, globalState.get(), replExParams.exchangeInterval > 0);
764
765     init_orires(fplog, mtop, inputrec, cr, globalState.get(), &(fcd->orires));
766
767     if (inputrecDeform(inputrec))
768     {
769         /* Store the deform reference box before reading the checkpoint */
770         if (SIMMASTER(cr))
771         {
772             copy_mat(globalState->box, box);
773         }
774         if (PAR(cr))
775         {
776             gmx_bcast(sizeof(box), box, cr);
777         }
778         /* Because we do not have the update struct available yet
779          * in which the reference values should be stored,
780          * we store them temporarily in static variables.
781          * This should be thread safe, since they are only written once
782          * and with identical values.
783          */
784         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
785         deform_init_init_step_tpx = inputrec->init_step;
786         copy_mat(box, deform_init_box_tpx);
787         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
788     }
789
790     ObservablesHistory   observablesHistory = {};
791
792     ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
793
794     if (continuationOptions.startedFromCheckpoint)
795     {
796         /* Check if checkpoint file exists before doing continuation.
797          * This way we can use identical input options for the first and subsequent runs...
798          */
799         gmx_bool bReadEkin;
800
801         load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
802                         cr, domdecOptions.numCells,
803                         inputrec, globalState.get(),
804                         &bReadEkin, &observablesHistory,
805                         continuationOptions.appendFiles,
806                         continuationOptions.appendFilesOptionSet,
807                         mdrunOptions.reproducible);
808
809         if (bReadEkin)
810         {
811             continuationOptions.haveReadEkin = true;
812         }
813     }
814
815     if (SIMMASTER(cr) && continuationOptions.appendFiles)
816     {
817         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
818                      continuationOptions.appendFiles, &fplog);
819         logOwner = buildLogger(fplog, nullptr);
820         mdlog    = logOwner.logger();
821     }
822
823     /* override nsteps with value set on the commamdline */
824     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
825
826     if (SIMMASTER(cr))
827     {
828         copy_mat(globalState->box, box);
829     }
830
831     if (PAR(cr))
832     {
833         gmx_bcast(sizeof(box), box, cr);
834     }
835
836     /* Update rlist and nstlist. */
837     if (inputrec->cutoff_scheme == ecutsVERLET)
838     {
839         prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, mtop, box,
840                               useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
841     }
842
843     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
844                      inputrec->eI == eiNM))
845     {
846         const rvec *xOnMaster = (SIMMASTER(cr) ? as_rvec_array(globalState->x.data()) : nullptr);
847
848         cr->dd = init_domain_decomposition(fplog, cr, domdecOptions, mdrunOptions,
849                                            mtop, inputrec,
850                                            box, xOnMaster,
851                                            &ddbox, &npme_major, &npme_minor);
852         // Note that local state still does not exist yet.
853     }
854     else
855     {
856         /* PME, if used, is done on all nodes with 1D decomposition */
857         cr->npmenodes = 0;
858         cr->duty      = (DUTY_PP | DUTY_PME);
859         npme_major    = 1;
860         npme_minor    = 1;
861
862         if (inputrec->ePBC == epbcSCREW)
863         {
864             gmx_fatal(FARGS,
865                       "pbc=%s is only implemented with domain decomposition",
866                       epbc_names[inputrec->ePBC]);
867         }
868     }
869
870     if (PAR(cr))
871     {
872         /* After possible communicator splitting in make_dd_communicators.
873          * we can set up the intra/inter node communication.
874          */
875         gmx_setup_nodecomm(fplog, cr);
876     }
877
878     /* Initialize per-physical-node MPI process/thread ID and counters. */
879     gmx_init_intranode_counters(cr);
880 #if GMX_MPI
881     if (MULTISIM(cr))
882     {
883         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
884                 "This is simulation %d out of %d running as a composite GROMACS\n"
885                 "multi-simulation job. Setup for this simulation:\n",
886                 cr->ms->sim, cr->ms->nsim);
887     }
888     GMX_LOG(mdlog.warning).appendTextFormatted(
889             "Using %d MPI %s\n",
890             cr->nnodes,
891 #if GMX_THREAD_MPI
892             cr->nnodes == 1 ? "thread" : "threads"
893 #else
894             cr->nnodes == 1 ? "process" : "processes"
895 #endif
896             );
897     fflush(stderr);
898 #endif
899
900     /* Check and update hw_opt for the cut-off scheme */
901     check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
902
903     /* Check and update the number of OpenMP threads requested */
904     checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, pmeRunMode, *mtop);
905
906     gmx_omp_nthreads_init(mdlog, cr,
907                           hwinfo->nthreads_hw_avail,
908                           hw_opt.nthreads_omp,
909                           hw_opt.nthreads_omp_pme,
910                           !thisRankHasDuty(cr, DUTY_PP),
911                           inputrec->cutoff_scheme == ecutsVERLET);
912
913 #ifndef NDEBUG
914     if (EI_TPI(inputrec->eI) &&
915         inputrec->cutoff_scheme == ecutsVERLET)
916     {
917         gmx_feenableexcept();
918     }
919 #endif
920
921     // Build a data structure that expresses which kinds of non-bonded
922     // task are handled by this rank.
923     //
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.
926     //
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.).
932     //
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.
937     //
938     // TODO Later we might need the concept of computeTasksOnThisRank,
939     // from which we construct gpuTasksOnThisRank.
940     //
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))
950     {
951         if (useGpuForNonbonded)
952         {
953             if (haveGpus)
954             {
955                 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
956             }
957             else if (nonbondedTarget == TaskTarget::Gpu)
958             {
959                 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
960             }
961         }
962     }
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)))
965     {
966         if (useGpuForPme)
967         {
968             if (haveGpus)
969             {
970                 gpuTasksOnThisRank.push_back(GpuTask::Pme);
971             }
972             else if (pmeTarget == TaskTarget::Gpu)
973             {
974                 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
975             }
976         }
977     }
978
979     GpuTaskAssignment gpuTaskAssignment;
980     try
981     {
982         // Produce the task assignment for this rank.
983         gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
984                                               mdlog, cr, gpuTasksOnThisRank);
985     }
986     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
987
988     /* Prevent other ranks from continuing after an issue was found
989      * and reported as a fatal error.
990      *
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))
998     {
999         MPI_Barrier(cr->mpi_comm_mysim);
1000     }
1001     if (MULTISIM(cr))
1002     {
1003         MPI_Barrier(cr->ms->mpi_comm_masters);
1004     }
1005 #endif
1006
1007     /* Now that we know the setup is consistent, check for efficiency */
1008     check_resource_division_efficiency(hwinfo, hw_opt.nthreads_tot, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1009                                        cr, mdlog);
1010
1011     gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1012     int                nonbondedDeviceId   = -1;
1013     if (thisRankHasDuty(cr, DUTY_PP))
1014     {
1015         if (!gpuTaskAssignment.empty())
1016         {
1017             GMX_RELEASE_ASSERT(gpuTaskAssignment.size() == 1, "A valid GPU assignment can only have one task per rank");
1018             GMX_RELEASE_ASSERT(gpuTaskAssignment[0].task_ == gmx::GpuTask::Nonbonded, "A valid GPU assignment can only include short-ranged tasks");
1019             nonbondedDeviceId   = gpuTaskAssignment[0].deviceId_;
1020             nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1021         }
1022     }
1023
1024     if (DOMAINDECOMP(cr))
1025     {
1026         /* When we share GPUs over ranks, we need to know this for the DLB */
1027         dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1028     }
1029
1030     /* getting number of PP/PME threads
1031        PME: env variable should be read only on one node to make sure it is
1032        identical everywhere;
1033      */
1034     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1035
1036     wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1037
1038     if (PAR(cr))
1039     {
1040         /* Master synchronizes its value of reset_counters with all nodes
1041          * including PME only nodes */
1042         reset_counters = wcycle_get_reset_counters(wcycle);
1043         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1044         wcycle_set_reset_counters(wcycle, reset_counters);
1045     }
1046
1047     // Membrane embedding must be initialized before we call init_forcerec()
1048     if (doMembed)
1049     {
1050         if (MASTER(cr))
1051         {
1052             fprintf(stderr, "Initializing membed");
1053         }
1054         /* Note that membed cannot work in parallel because mtop is
1055          * changed here. Fix this if we ever want to make it run with
1056          * multiple ranks. */
1057         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1058     }
1059
1060     std::unique_ptr<MDAtoms> mdAtoms;
1061
1062     snew(nrnb, 1);
1063     if (thisRankHasDuty(cr, DUTY_PP))
1064     {
1065         /* Initiate forcerecord */
1066         fr                 = mk_forcerec();
1067         fr->forceProviders = mdModules.initForceProviders();
1068         init_forcerec(fplog, mdlog, fr, fcd,
1069                       inputrec, mtop, cr, box,
1070                       opt2fn("-table", nfile, fnm),
1071                       opt2fn("-tablep", nfile, fnm),
1072                       getFilenm("-tableb", nfile, fnm),
1073                       nonbondedDeviceInfo,
1074                       FALSE,
1075                       pforce);
1076
1077         /* Initialize QM-MM */
1078         if (fr->bQMMM)
1079         {
1080             init_QMMMrec(cr, mtop, inputrec, fr);
1081         }
1082
1083         /* Initialize the mdAtoms structure.
1084          * mdAtoms is not filled with atom data,
1085          * as this can not be done now with domain decomposition.
1086          */
1087         const bool useGpuForPme = (pmeRunMode == PmeRunMode::GPU) || (pmeRunMode == PmeRunMode::Hybrid);
1088         mdAtoms = makeMDAtoms(fplog, *mtop, *inputrec, useGpuForPme && thisRankHasDuty(cr, DUTY_PME));
1089         if (globalState)
1090         {
1091             // The pinning of coordinates in the global state object works, because we only use
1092             // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1093             // points to the global state object without DD.
1094             // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1095             // which should also perform the pinning.
1096             changePinningPolicy(&globalState->x, useGpuForPme ? PinningPolicy::CanBePinned : PinningPolicy::CannotBePinned);
1097         }
1098
1099         /* Initialize the virtual site communication */
1100         vsite = initVsite(*mtop, cr);
1101
1102         calc_shifts(box, fr->shift_vec);
1103
1104         /* With periodic molecules the charge groups should be whole at start up
1105          * and the virtual sites should not be far from their proper positions.
1106          */
1107         if (!inputrec->bContinuation && MASTER(cr) &&
1108             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1109         {
1110             /* Make molecules whole at start of run */
1111             if (fr->ePBC != epbcNONE)
1112             {
1113                 rvec *xGlobal = as_rvec_array(globalState->x.data());
1114                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, xGlobal);
1115             }
1116             if (vsite)
1117             {
1118                 /* Correct initial vsite positions are required
1119                  * for the initial distribution in the domain decomposition
1120                  * and for the initial shell prediction.
1121                  */
1122                 constructVsitesGlobal(*mtop, globalState->x);
1123             }
1124         }
1125
1126         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1127         {
1128             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1129             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1130         }
1131     }
1132     else
1133     {
1134         /* This is a PME only node */
1135
1136         GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1137
1138         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1139         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1140     }
1141
1142     gmx_pme_t *sepPmeData = nullptr;
1143     // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1144     GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1145     gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1146
1147     if (hw_opt.thread_affinity != threadaffOFF)
1148     {
1149         /* Before setting affinity, check whether the affinity has changed
1150          * - which indicates that probably the OpenMP library has changed it
1151          * since we first checked).
1152          */
1153         gmx_check_thread_affinity_set(mdlog, cr,
1154                                       &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1155
1156         int nthread_local;
1157         /* threads on this MPI process or TMPI thread */
1158         if (thisRankHasDuty(cr, DUTY_PP))
1159         {
1160             nthread_local = gmx_omp_nthreads_get(emntNonbonded);
1161         }
1162         else
1163         {
1164             nthread_local = gmx_omp_nthreads_get(emntPME);
1165         }
1166
1167         /* Set the CPU affinity */
1168         gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1169                                 nthread_local, nullptr);
1170     }
1171
1172     /* Initiate PME if necessary,
1173      * either on all nodes or on dedicated PME nodes only. */
1174     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1175     {
1176         if (mdAtoms && mdAtoms->mdatoms())
1177         {
1178             nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1179             if (EVDW_PME(inputrec->vdwtype))
1180             {
1181                 nTypePerturbed   = mdAtoms->mdatoms()->nTypePerturbed;
1182             }
1183         }
1184         if (cr->npmenodes > 0)
1185         {
1186             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1187             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1188             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1189         }
1190
1191         if (thisRankHasDuty(cr, DUTY_PME))
1192         {
1193             try
1194             {
1195                 gmx_device_info_t *pmeGpuInfo = nullptr;
1196                 pmedata = gmx_pme_init(cr, npme_major, npme_minor, inputrec,
1197                                        mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1198                                        mdrunOptions.reproducible,
1199                                        ewaldcoeff_q, ewaldcoeff_lj,
1200                                        nthreads_pme,
1201                                        pmeRunMode, nullptr, pmeGpuInfo, mdlog);
1202             }
1203             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1204         }
1205     }
1206
1207
1208     if (EI_DYNAMICS(inputrec->eI))
1209     {
1210         /* Turn on signal handling on all nodes */
1211         /*
1212          * (A user signal from the PME nodes (if any)
1213          * is communicated to the PP nodes.
1214          */
1215         signal_handler_install();
1216     }
1217
1218     if (thisRankHasDuty(cr, DUTY_PP))
1219     {
1220         /* Assumes uniform use of the number of OpenMP threads */
1221         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1222
1223         if (inputrec->bPull)
1224         {
1225             /* Initialize pull code */
1226             inputrec->pull_work =
1227                 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1228                           mtop, cr, oenv, inputrec->fepvals->init_lambda,
1229                           EI_DYNAMICS(inputrec->eI) && MASTER(cr),
1230                           continuationOptions);
1231         }
1232
1233         if (inputrec->bRot)
1234         {
1235             /* Initialize enforced rotation code */
1236             init_rot(fplog, inputrec, nfile, fnm, cr, globalState.get(), mtop, oenv, mdrunOptions);
1237         }
1238
1239         /* Let init_constraints know whether we have essential dynamics constraints.
1240          * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1241          */
1242         bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1243
1244         constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1245
1246         if (DOMAINDECOMP(cr))
1247         {
1248             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1249             /* This call is not included in init_domain_decomposition mainly
1250              * because fr->cginfo_mb is set later.
1251              */
1252             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1253                             domdecOptions.checkBondedInteractions,
1254                             fr->cginfo_mb);
1255         }
1256
1257         /* Now do whatever the user wants us to do (how flexible...) */
1258         my_integrator(inputrec->eI) (fplog, cr, mdlog, nfile, fnm,
1259                                      oenv,
1260                                      mdrunOptions,
1261                                      vsite, constr,
1262                                      mdModules.outputProvider(),
1263                                      inputrec, mtop,
1264                                      fcd,
1265                                      globalState.get(),
1266                                      &observablesHistory,
1267                                      mdAtoms.get(), nrnb, wcycle, fr,
1268                                      replExParams,
1269                                      membed,
1270                                      walltime_accounting);
1271
1272         if (inputrec->bRot)
1273         {
1274             finish_rot(inputrec->rot);
1275         }
1276
1277         if (inputrec->bPull)
1278         {
1279             finish_pull(inputrec->pull_work);
1280         }
1281
1282     }
1283     else
1284     {
1285         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1286         /* do PME only */
1287         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1288         gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1289     }
1290
1291     wallcycle_stop(wcycle, ewcRUN);
1292
1293     /* Finish up, write some stuff
1294      * if rerunMD, don't write last frame again
1295      */
1296     finish_run(fplog, mdlog, cr,
1297                inputrec, nrnb, wcycle, walltime_accounting,
1298                fr ? fr->nbv : nullptr,
1299                pmedata,
1300                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1301
1302     // Free PME data
1303     if (pmedata)
1304     {
1305         gmx_pme_destroy(pmedata);
1306         pmedata = nullptr;
1307     }
1308
1309     /* Free GPU memory and context */
1310     free_gpu_resources(fr, cr, nonbondedDeviceInfo);
1311
1312     if (doMembed)
1313     {
1314         free_membed(membed);
1315     }
1316
1317     gmx_hardware_info_free();
1318
1319     /* Does what it says */
1320     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1321     walltime_accounting_destroy(walltime_accounting);
1322
1323     /* Close logfile already here if we were appending to it */
1324     if (MASTER(cr) && continuationOptions.appendFiles)
1325     {
1326         gmx_log_close(fplog);
1327     }
1328
1329     rc = (int)gmx_get_stop_condition();
1330
1331 #if GMX_THREAD_MPI
1332     /* we need to join all threads. The sub-threads join when they
1333        exit this function, but the master thread needs to be told to
1334        wait for that. */
1335     if (PAR(cr) && MASTER(cr))
1336     {
1337         tMPI_Finalize();
1338     }
1339 #endif
1340
1341     return rc;
1342 }
1343
1344 } // namespace gmx