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