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