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