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