Use more forward declarations to removed header dependencies
[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-2019,2020, 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 <cinttypes>
52 #include <csignal>
53 #include <cstdlib>
54 #include <cstring>
55
56 #include <algorithm>
57 #include <memory>
58
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/ewald/pme_only.h"
69 #include "gromacs/ewald/pme_pp_comm_gpu.h"
70 #include "gromacs/fileio/checkpoint.h"
71 #include "gromacs/fileio/gmxfio.h"
72 #include "gromacs/fileio/oenv.h"
73 #include "gromacs/fileio/tpxio.h"
74 #include "gromacs/gmxlib/network.h"
75 #include "gromacs/gmxlib/nrnb.h"
76 #include "gromacs/gpu_utils/gpu_utils.h"
77 #include "gromacs/hardware/cpuinfo.h"
78 #include "gromacs/hardware/detecthardware.h"
79 #include "gromacs/hardware/printhardware.h"
80 #include "gromacs/imd/imd.h"
81 #include "gromacs/listed_forces/disre.h"
82 #include "gromacs/listed_forces/gpubonded.h"
83 #include "gromacs/listed_forces/orires.h"
84 #include "gromacs/math/functions.h"
85 #include "gromacs/math/utilities.h"
86 #include "gromacs/math/vec.h"
87 #include "gromacs/mdlib/boxdeformation.h"
88 #include "gromacs/mdlib/broadcaststructs.h"
89 #include "gromacs/mdlib/calc_verletbuf.h"
90 #include "gromacs/mdlib/dispersioncorrection.h"
91 #include "gromacs/mdlib/enerdata_utils.h"
92 #include "gromacs/mdlib/force.h"
93 #include "gromacs/mdlib/forcerec.h"
94 #include "gromacs/mdlib/gmx_omp_nthreads.h"
95 #include "gromacs/mdlib/makeconstraints.h"
96 #include "gromacs/mdlib/md_support.h"
97 #include "gromacs/mdlib/mdatoms.h"
98 #include "gromacs/mdlib/membed.h"
99 #include "gromacs/mdlib/qmmm.h"
100 #include "gromacs/mdlib/sighandler.h"
101 #include "gromacs/mdlib/stophandler.h"
102 #include "gromacs/mdlib/updategroups.h"
103 #include "gromacs/mdlib/vsite.h"
104 #include "gromacs/mdrun/mdmodules.h"
105 #include "gromacs/mdrun/simulationcontext.h"
106 #include "gromacs/mdrunutility/handlerestart.h"
107 #include "gromacs/mdrunutility/logging.h"
108 #include "gromacs/mdrunutility/multisim.h"
109 #include "gromacs/mdrunutility/printtime.h"
110 #include "gromacs/mdrunutility/threadaffinity.h"
111 #include "gromacs/mdtypes/commrec.h"
112 #include "gromacs/mdtypes/enerdata.h"
113 #include "gromacs/mdtypes/fcdata.h"
114 #include "gromacs/mdtypes/forcerec.h"
115 #include "gromacs/mdtypes/group.h"
116 #include "gromacs/mdtypes/inputrec.h"
117 #include "gromacs/mdtypes/interaction_const.h"
118 #include "gromacs/mdtypes/md_enums.h"
119 #include "gromacs/mdtypes/mdatom.h"
120 #include "gromacs/mdtypes/mdrunoptions.h"
121 #include "gromacs/mdtypes/observableshistory.h"
122 #include "gromacs/mdtypes/simulation_workload.h"
123 #include "gromacs/mdtypes/state.h"
124 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
125 #include "gromacs/nbnxm/gpu_data_mgmt.h"
126 #include "gromacs/nbnxm/nbnxm.h"
127 #include "gromacs/nbnxm/pairlist_tuning.h"
128 #include "gromacs/pbcutil/pbc.h"
129 #include "gromacs/pulling/output.h"
130 #include "gromacs/pulling/pull.h"
131 #include "gromacs/pulling/pull_rotation.h"
132 #include "gromacs/restraint/manager.h"
133 #include "gromacs/restraint/restraintmdmodule.h"
134 #include "gromacs/restraint/restraintpotential.h"
135 #include "gromacs/swap/swapcoords.h"
136 #include "gromacs/taskassignment/decidegpuusage.h"
137 #include "gromacs/taskassignment/decidesimulationworkload.h"
138 #include "gromacs/taskassignment/resourcedivision.h"
139 #include "gromacs/taskassignment/taskassignment.h"
140 #include "gromacs/taskassignment/usergpuids.h"
141 #include "gromacs/timing/gpu_timing.h"
142 #include "gromacs/timing/wallcycle.h"
143 #include "gromacs/timing/wallcyclereporting.h"
144 #include "gromacs/topology/mtop_util.h"
145 #include "gromacs/trajectory/trajectoryframe.h"
146 #include "gromacs/utility/basenetwork.h"
147 #include "gromacs/utility/cstringutil.h"
148 #include "gromacs/utility/exceptions.h"
149 #include "gromacs/utility/fatalerror.h"
150 #include "gromacs/utility/filestream.h"
151 #include "gromacs/utility/gmxassert.h"
152 #include "gromacs/utility/gmxmpi.h"
153 #include "gromacs/utility/keyvaluetree.h"
154 #include "gromacs/utility/logger.h"
155 #include "gromacs/utility/loggerbuilder.h"
156 #include "gromacs/utility/mdmodulenotification.h"
157 #include "gromacs/utility/physicalnodecommunicator.h"
158 #include "gromacs/utility/pleasecite.h"
159 #include "gromacs/utility/programcontext.h"
160 #include "gromacs/utility/smalloc.h"
161 #include "gromacs/utility/stringutil.h"
162
163 #include "isimulator.h"
164 #include "replicaexchange.h"
165 #include "simulatorbuilder.h"
166
167 #if GMX_FAHCORE
168 #    include "corewrap.h"
169 #endif
170
171 namespace gmx
172 {
173
174 /*! \brief Structure that holds boolean flags corresponding to the development
175  *        features present enabled through environment variables.
176  *
177  */
178 struct DevelopmentFeatureFlags
179 {
180     //! True if the Buffer ops development feature is enabled
181     // TODO: when the trigger of the buffer ops offload is fully automated this should go away
182     bool enableGpuBufferOps = false;
183     //! If true, forces 'mdrun -update auto' default to 'gpu'
184     bool forceGpuUpdateDefault = false;
185     //! True if the GPU halo exchange development feature is enabled
186     bool enableGpuHaloExchange = false;
187     //! True if the PME PP direct communication GPU development feature is enabled
188     bool enableGpuPmePPComm = false;
189 };
190
191 /*! \brief Manage any development feature flag variables encountered
192  *
193  * The use of dev features indicated by environment variables is
194  * logged in order to ensure that runs with such features enabled can
195  * be identified from their log and standard output. Any cross
196  * dependencies are also checked, and if unsatisfied, a fatal error
197  * issued.
198  *
199  * Note that some development features overrides are applied already here:
200  * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
201  *
202  * \param[in]  mdlog                Logger object.
203  * \param[in]  useGpuForNonbonded   True if the nonbonded task is offloaded in this run.
204  * \param[in]  pmeRunMode           The PME run mode for this run
205  * \returns                         The object populated with development feature flags.
206  */
207 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
208                                                          const bool           useGpuForNonbonded,
209                                                          const PmeRunMode     pmeRunMode)
210 {
211     DevelopmentFeatureFlags devFlags;
212
213     // Some builds of GCC 5 give false positive warnings that these
214     // getenv results are ignored when clearly they are used.
215 #pragma GCC diagnostic push
216 #pragma GCC diagnostic ignored "-Wunused-result"
217     devFlags.enableGpuBufferOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr)
218                                   && (GMX_GPU == GMX_GPU_CUDA) && useGpuForNonbonded;
219     devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr);
220     devFlags.enableGpuHaloExchange =
221             (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
222     devFlags.enableGpuPmePPComm =
223             (getenv("GMX_GPU_PME_PP_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
224 #pragma GCC diagnostic pop
225
226     if (devFlags.enableGpuBufferOps)
227     {
228         GMX_LOG(mdlog.warning)
229                 .asParagraph()
230                 .appendTextFormatted(
231                         "This run uses the 'GPU buffer ops' feature, enabled by the "
232                         "GMX_USE_GPU_BUFFER_OPS environment variable.");
233     }
234
235     if (devFlags.forceGpuUpdateDefault)
236     {
237         GMX_LOG(mdlog.warning)
238                 .asParagraph()
239                 .appendTextFormatted(
240                         "This run will default to '-update gpu' as requested by the "
241                         "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
242                         "decomposition lacks substantial testing and should be used with caution.");
243     }
244
245     if (devFlags.enableGpuHaloExchange)
246     {
247         if (useGpuForNonbonded)
248         {
249             if (!devFlags.enableGpuBufferOps)
250             {
251                 GMX_LOG(mdlog.warning)
252                         .asParagraph()
253                         .appendTextFormatted(
254                                 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
255                                 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
256                 devFlags.enableGpuBufferOps = true;
257             }
258             GMX_LOG(mdlog.warning)
259                     .asParagraph()
260                     .appendTextFormatted(
261                             "This run has requested the 'GPU halo exchange' feature, enabled by "
262                             "the "
263                             "GMX_GPU_DD_COMMS environment variable.");
264         }
265         else
266         {
267             GMX_LOG(mdlog.warning)
268                     .asParagraph()
269                     .appendTextFormatted(
270                             "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
271                             "halo exchange' feature will not be enabled as nonbonded interactions "
272                             "are not offloaded.");
273             devFlags.enableGpuHaloExchange = false;
274         }
275     }
276
277     if (devFlags.enableGpuPmePPComm)
278     {
279         if (pmeRunMode == PmeRunMode::GPU)
280         {
281             GMX_LOG(mdlog.warning)
282                     .asParagraph()
283                     .appendTextFormatted(
284                             "This run uses the 'GPU PME-PP communications' feature, enabled "
285                             "by the GMX_GPU_PME_PP_COMMS environment variable.");
286         }
287         else
288         {
289             std::string clarification;
290             if (pmeRunMode == PmeRunMode::Mixed)
291             {
292                 clarification =
293                         "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
294                         "mode).";
295             }
296             else
297             {
298                 clarification = "PME is not offloaded to the GPU.";
299             }
300             GMX_LOG(mdlog.warning)
301                     .asParagraph()
302                     .appendText(
303                             "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
304                             "'GPU PME-PP communications' feature was not enabled as "
305                             + clarification);
306             devFlags.enableGpuPmePPComm = false;
307         }
308     }
309
310     return devFlags;
311 }
312
313 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
314  *
315  * Used to ensure that the master thread does not modify mdrunner during copy
316  * on the spawned threads. */
317 static void threadMpiMdrunnerAccessBarrier()
318 {
319 #if GMX_THREAD_MPI
320     MPI_Barrier(MPI_COMM_WORLD);
321 #endif
322 }
323
324 Mdrunner Mdrunner::cloneOnSpawnedThread() const
325 {
326     auto newRunner = Mdrunner(std::make_unique<MDModules>());
327
328     // All runners in the same process share a restraint manager resource because it is
329     // part of the interface to the client code, which is associated only with the
330     // original thread. Handles to the same resources can be obtained by copy.
331     {
332         newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
333     }
334
335     // Copy members of master runner.
336     // \todo Replace with builder when Simulation context and/or runner phases are better defined.
337     // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
338     newRunner.hw_opt    = hw_opt;
339     newRunner.filenames = filenames;
340
341     newRunner.oenv            = oenv;
342     newRunner.mdrunOptions    = mdrunOptions;
343     newRunner.domdecOptions   = domdecOptions;
344     newRunner.nbpu_opt        = nbpu_opt;
345     newRunner.pme_opt         = pme_opt;
346     newRunner.pme_fft_opt     = pme_fft_opt;
347     newRunner.bonded_opt      = bonded_opt;
348     newRunner.update_opt      = update_opt;
349     newRunner.nstlist_cmdline = nstlist_cmdline;
350     newRunner.replExParams    = replExParams;
351     newRunner.pforce          = pforce;
352     // Give the spawned thread the newly created valid communicator
353     // for the simulation.
354     newRunner.communicator        = MPI_COMM_WORLD;
355     newRunner.ms                  = ms;
356     newRunner.startingBehavior    = startingBehavior;
357     newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
358
359     threadMpiMdrunnerAccessBarrier();
360
361     return newRunner;
362 }
363
364 /*! \brief The callback used for running on spawned threads.
365  *
366  * Obtains the pointer to the master mdrunner object from the one
367  * argument permitted to the thread-launch API call, copies it to make
368  * a new runner for this thread, reinitializes necessary data, and
369  * proceeds to the simulation. */
370 static void mdrunner_start_fn(const void* arg)
371 {
372     try
373     {
374         auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
375         /* copy the arg list to make sure that it's thread-local. This
376            doesn't copy pointed-to items, of course; fnm, cr and fplog
377            are reset in the call below, all others should be const. */
378         gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
379         mdrunner.mdrunner();
380     }
381     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
382 }
383
384
385 void Mdrunner::spawnThreads(int numThreadsToLaunch)
386 {
387 #if GMX_THREAD_MPI
388     /* now spawn new threads that start mdrunner_start_fn(), while
389        the main thread returns. Thread affinity is handled later. */
390     if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn,
391                      static_cast<const void*>(this))
392         != TMPI_SUCCESS)
393     {
394         GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
395     }
396
397     // Give the master thread the newly created valid communicator for
398     // the simulation.
399     communicator = MPI_COMM_WORLD;
400     threadMpiMdrunnerAccessBarrier();
401 #else
402     GMX_UNUSED_VALUE(numThreadsToLaunch);
403     GMX_UNUSED_VALUE(mdrunner_start_fn);
404 #endif
405 }
406
407 } // namespace gmx
408
409 /*! \brief Initialize variables for Verlet scheme simulation */
410 static void prepare_verlet_scheme(FILE*               fplog,
411                                   t_commrec*          cr,
412                                   t_inputrec*         ir,
413                                   int                 nstlist_cmdline,
414                                   const gmx_mtop_t*   mtop,
415                                   const matrix        box,
416                                   bool                makeGpuPairList,
417                                   const gmx::CpuInfo& cpuinfo)
418 {
419     // We checked the cut-offs in grompp, but double-check here.
420     // We have PME+LJcutoff kernels for rcoulomb>rvdw.
421     if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
422     {
423         GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
424                            "With Verlet lists and PME we should have rcoulomb>=rvdw");
425     }
426     else
427     {
428         GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
429                            "With Verlet lists and no PME rcoulomb and rvdw should be identical");
430     }
431     /* For NVE simulations, we will retain the initial list buffer */
432     if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
433     {
434         /* Update the Verlet buffer size for the current run setup */
435
436         /* Here we assume SIMD-enabled kernels are being used. But as currently
437          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
438          * and 4x2 gives a larger buffer than 4x4, this is ok.
439          */
440         ListSetupType listType =
441                 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
442         VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
443
444         const real rlist_new =
445                 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
446
447         if (rlist_new != ir->rlist)
448         {
449             if (fplog != nullptr)
450             {
451                 fprintf(fplog,
452                         "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
453                         ir->rlist, rlist_new, listSetup.cluster_size_i, listSetup.cluster_size_j);
454             }
455             ir->rlist = rlist_new;
456         }
457     }
458
459     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
460     {
461         gmx_fatal(FARGS, "Can not set nstlist without %s",
462                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
463     }
464
465     if (EI_DYNAMICS(ir->eI))
466     {
467         /* Set or try nstlist values */
468         increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
469     }
470 }
471
472 /*! \brief Override the nslist value in inputrec
473  *
474  * with value passed on the command line (if any)
475  */
476 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
477 {
478     assert(ir);
479
480     /* override with anything else than the default -2 */
481     if (nsteps_cmdline > -2)
482     {
483         char sbuf_steps[STEPSTRSIZE];
484         char sbuf_msg[STRLEN];
485
486         ir->nsteps = nsteps_cmdline;
487         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
488         {
489             sprintf(sbuf_msg,
490                     "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
491                     gmx_step_str(nsteps_cmdline, sbuf_steps), fabs(nsteps_cmdline * ir->delta_t));
492         }
493         else
494         {
495             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
496                     gmx_step_str(nsteps_cmdline, sbuf_steps));
497         }
498
499         GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
500     }
501     else if (nsteps_cmdline < -2)
502     {
503         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
504     }
505     /* Do nothing if nsteps_cmdline == -2 */
506 }
507
508 namespace gmx
509 {
510
511 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
512  *
513  * If not, and if a warning may be issued, logs a warning about
514  * falling back to CPU code. With thread-MPI, only the first
515  * call to this function should have \c issueWarning true. */
516 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
517 {
518     bool        gpuIsUseful = true;
519     std::string warning;
520
521     if (ir.opts.ngener - ir.nwall > 1)
522     {
523         /* The GPU code does not support more than one energy group.
524          * If the user requested GPUs explicitly, a fatal error is given later.
525          */
526         gpuIsUseful = false;
527         warning =
528                 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
529                 "For better performance, run on the GPU without energy groups and then do "
530                 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
531     }
532
533     if (EI_TPI(ir.eI))
534     {
535         gpuIsUseful = false;
536         warning     = "TPI is not implemented for GPUs.";
537     }
538
539     if (!gpuIsUseful && issueWarning)
540     {
541         GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
542     }
543
544     return gpuIsUseful;
545 }
546
547 //! Initializes the logger for mdrun.
548 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
549 {
550     gmx::LoggerBuilder builder;
551     if (fplog != nullptr)
552     {
553         builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
554     }
555     if (isSimulationMasterRank)
556     {
557         builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
558     }
559     return builder.build();
560 }
561
562 //! Make a TaskTarget from an mdrun argument string.
563 static TaskTarget findTaskTarget(const char* optionString)
564 {
565     TaskTarget returnValue = TaskTarget::Auto;
566
567     if (strncmp(optionString, "auto", 3) == 0)
568     {
569         returnValue = TaskTarget::Auto;
570     }
571     else if (strncmp(optionString, "cpu", 3) == 0)
572     {
573         returnValue = TaskTarget::Cpu;
574     }
575     else if (strncmp(optionString, "gpu", 3) == 0)
576     {
577         returnValue = TaskTarget::Gpu;
578     }
579     else
580     {
581         GMX_ASSERT(false, "Option string should have been checked for sanity already");
582     }
583
584     return returnValue;
585 }
586
587 //! Finish run, aggregate data to print performance info.
588 static void finish_run(FILE*                     fplog,
589                        const gmx::MDLogger&      mdlog,
590                        const t_commrec*          cr,
591                        const t_inputrec*         inputrec,
592                        t_nrnb                    nrnb[],
593                        gmx_wallcycle_t           wcycle,
594                        gmx_walltime_accounting_t walltime_accounting,
595                        nonbonded_verlet_t*       nbv,
596                        const gmx_pme_t*          pme,
597                        gmx_bool                  bWriteStat)
598 {
599     double delta_t = 0;
600     double nbfs = 0, mflop = 0;
601     double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
602             elapsed_time_over_all_threads_over_all_ranks;
603     /* Control whether it is valid to print a report. Only the
604        simulation master may print, but it should not do so if the run
605        terminated e.g. before a scheduled reset step. This is
606        complicated by the fact that PME ranks are unaware of the
607        reason why they were sent a pmerecvqxFINISH. To avoid
608        communication deadlocks, we always do the communication for the
609        report, even if we've decided not to write the report, because
610        how long it takes to finish the run is not important when we've
611        decided not to report on the simulation performance.
612
613        Further, we only report performance for dynamical integrators,
614        because those are the only ones for which we plan to
615        consider doing any optimizations. */
616     bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
617
618     if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
619     {
620         GMX_LOG(mdlog.warning)
621                 .asParagraph()
622                 .appendText("Simulation ended prematurely, no performance report will be written.");
623         printReport = false;
624     }
625
626     t_nrnb*                 nrnb_tot;
627     std::unique_ptr<t_nrnb> nrnbTotalStorage;
628     if (cr->nnodes > 1)
629     {
630         nrnbTotalStorage = std::make_unique<t_nrnb>();
631         nrnb_tot         = nrnbTotalStorage.get();
632 #if GMX_MPI
633         MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
634 #endif
635     }
636     else
637     {
638         nrnb_tot = nrnb;
639     }
640
641     elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
642     elapsed_time_over_all_threads =
643             walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
644     if (cr->nnodes > 1)
645     {
646 #if GMX_MPI
647         /* reduce elapsed_time over all MPI ranks in the current simulation */
648         MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM,
649                       cr->mpi_comm_mysim);
650         elapsed_time_over_all_ranks /= cr->nnodes;
651         /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
652          * current simulation. */
653         MPI_Allreduce(&elapsed_time_over_all_threads, &elapsed_time_over_all_threads_over_all_ranks,
654                       1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
655 #endif
656     }
657     else
658     {
659         elapsed_time_over_all_ranks                  = elapsed_time;
660         elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
661     }
662
663     if (printReport)
664     {
665         print_flop(fplog, nrnb_tot, &nbfs, &mflop);
666     }
667
668     if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
669     {
670         print_dd_statistics(cr, inputrec, fplog);
671     }
672
673     /* TODO Move the responsibility for any scaling by thread counts
674      * to the code that handled the thread region, so that there's a
675      * mechanism to keep cycle counting working during the transition
676      * to task parallelism. */
677     int nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
678     int nthreads_pme = gmx_omp_nthreads_get(emntPME);
679     wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP),
680                                    nthreads_pp, nthreads_pme);
681     auto cycle_sum(wallcycle_sum(cr, wcycle));
682
683     if (printReport)
684     {
685         auto nbnxn_gpu_timings =
686                 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
687         gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
688
689         if (pme_gpu_task_enabled(pme))
690         {
691             pme_gpu_get_timings(pme, &pme_gpu_timings);
692         }
693         wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
694                         elapsed_time_over_all_ranks, wcycle, cycle_sum, nbnxn_gpu_timings,
695                         &pme_gpu_timings);
696
697         if (EI_DYNAMICS(inputrec->eI))
698         {
699             delta_t = inputrec->delta_t;
700         }
701
702         if (fplog)
703         {
704             print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
705                        walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
706                        delta_t, nbfs, mflop);
707         }
708         if (bWriteStat)
709         {
710             print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks, elapsed_time_over_all_ranks,
711                        walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
712                        delta_t, nbfs, mflop);
713         }
714     }
715 }
716
717 int Mdrunner::mdrunner()
718 {
719     matrix                    box;
720     t_forcerec*               fr               = nullptr;
721     t_fcdata*                 fcd              = nullptr;
722     real                      ewaldcoeff_q     = 0;
723     real                      ewaldcoeff_lj    = 0;
724     int                       nChargePerturbed = -1, nTypePerturbed = 0;
725     gmx_wallcycle_t           wcycle;
726     gmx_walltime_accounting_t walltime_accounting = nullptr;
727     gmx_membed_t*             membed              = nullptr;
728     gmx_hw_info_t*            hwinfo              = nullptr;
729
730     /* CAUTION: threads may be started later on in this function, so
731        cr doesn't reflect the final parallel state right now */
732     gmx_mtop_t mtop;
733
734     /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
735     const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
736     const bool doMembed            = opt2bSet("-membed", filenames.size(), filenames.data());
737     const bool doRerun             = mdrunOptions.rerun;
738
739     // Handle task-assignment related user options.
740     EmulateGpuNonbonded emulateGpuNonbonded =
741             (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
742
743     std::vector<int> userGpuTaskAssignment;
744     try
745     {
746         userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
747     }
748     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
749     auto nonbondedTarget = findTaskTarget(nbpu_opt);
750     auto pmeTarget       = findTaskTarget(pme_opt);
751     auto pmeFftTarget    = findTaskTarget(pme_fft_opt);
752     auto bondedTarget    = findTaskTarget(bonded_opt);
753     auto updateTarget    = findTaskTarget(update_opt);
754
755     FILE* fplog = nullptr;
756     // If we are appending, we don't write log output because we need
757     // to check that the old log file matches what the checkpoint file
758     // expects. Otherwise, we should start to write log output now if
759     // there is a file ready for it.
760     if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
761     {
762         fplog = gmx_fio_getfp(logFileHandle);
763     }
764     const bool       isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
765     gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
766     gmx::MDLogger    mdlog(logOwner.logger());
767
768     // TODO The thread-MPI master rank makes a working
769     // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
770     // after the threads have been launched. This works because no use
771     // is made of that communicator until after the execution paths
772     // have rejoined. But it is likely that we can improve the way
773     // this is expressed, e.g. by expressly running detection only the
774     // master rank for thread-MPI, rather than relying on the mutex
775     // and reference count.
776     PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
777     hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
778
779     gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
780
781     std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
782
783     // Print citation requests after all software/hardware printing
784     pleaseCiteGromacs(fplog);
785
786     // TODO Replace this by unique_ptr once t_inputrec is C++
787     t_inputrec               inputrecInstance;
788     t_inputrec*              inputrec = nullptr;
789     std::unique_ptr<t_state> globalState;
790
791     auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
792
793     if (isSimulationMasterRank)
794     {
795         /* Only the master rank has the global state */
796         globalState = std::make_unique<t_state>();
797
798         /* Read (nearly) all data required for the simulation
799          * and keep the partly serialized tpr contents to send to other ranks later
800          */
801         *partialDeserializedTpr = read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
802                                                  &inputrecInstance, globalState.get(), &mtop);
803         inputrec                = &inputrecInstance;
804     }
805
806     /* Check and update the hardware options for internal consistency */
807     checkAndUpdateHardwareOptions(mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks,
808                                   inputrec);
809
810     if (GMX_THREAD_MPI && isSimulationMasterRank)
811     {
812         bool useGpuForNonbonded = false;
813         bool useGpuForPme       = false;
814         try
815         {
816             GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
817
818             // If the user specified the number of ranks, then we must
819             // respect that, but in default mode, we need to allow for
820             // the number of GPUs to choose the number of ranks.
821             auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
822             useGpuForNonbonded         = decideWhetherToUseGpusForNonbondedWithThreadMpi(
823                     nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
824                     canUseGpuForNonbonded,
825                     gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
826                     hw_opt.nthreads_tmpi);
827             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
828                     useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
829                     *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
830         }
831         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
832
833         /* Determine how many thread-MPI ranks to start.
834          *
835          * TODO Over-writing the user-supplied value here does
836          * prevent any possible subsequent checks from working
837          * correctly. */
838         hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo, &hw_opt, gpuIdsToUse, useGpuForNonbonded,
839                                                 useGpuForPme, inputrec, &mtop, mdlog, doMembed);
840
841         // Now start the threads for thread MPI.
842         spawnThreads(hw_opt.nthreads_tmpi);
843         // The spawned threads enter mdrunner() and execution of
844         // master and spawned threads joins at the end of this block.
845         physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
846     }
847
848     GMX_RELEASE_ASSERT(communicator == MPI_COMM_WORLD, "Must have valid world communicator");
849     CommrecHandle crHandle = init_commrec(communicator, ms);
850     t_commrec*    cr       = crHandle.get();
851     GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
852
853     if (PAR(cr))
854     {
855         /* now broadcast everything to the non-master nodes/threads: */
856         if (!isSimulationMasterRank)
857         {
858             inputrec = &inputrecInstance;
859         }
860         init_parallel(cr, inputrec, &mtop, partialDeserializedTpr.get());
861     }
862     GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
863     partialDeserializedTpr.reset(nullptr);
864
865     // Now the number of ranks is known to all ranks, and each knows
866     // the inputrec read by the master rank. The ranks can now all run
867     // the task-deciding functions and will agree on the result
868     // without needing to communicate.
869     const bool useDomainDecomposition = (PAR(cr) && !(EI_TPI(inputrec->eI) || inputrec->eI == eiNM));
870
871     // Note that these variables describe only their own node.
872     //
873     // Note that when bonded interactions run on a GPU they always run
874     // alongside a nonbonded task, so do not influence task assignment
875     // even though they affect the force calculation workload.
876     bool useGpuForNonbonded = false;
877     bool useGpuForPme       = false;
878     bool useGpuForBonded    = false;
879     bool useGpuForUpdate    = false;
880     bool gpusWereDetected   = hwinfo->ngpu_compatible_tot > 0;
881     try
882     {
883         // It's possible that there are different numbers of GPUs on
884         // different nodes, which is the user's responsibility to
885         // handle. If unsuitable, we will notice that during task
886         // assignment.
887         auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
888         useGpuForNonbonded         = decideWhetherToUseGpusForNonbonded(
889                 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
890                 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
891         useGpuForPme = decideWhetherToUseGpusForPme(
892                 useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo, *inputrec, mtop,
893                 cr->nnodes, domdecOptions.numPmeRanks, gpusWereDetected);
894         auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr)
895                                   && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
896         useGpuForBonded = decideWhetherToUseGpusForBonded(
897                 useGpuForNonbonded, useGpuForPme, bondedTarget, canUseGpuForBonded,
898                 EVDW_PME(inputrec->vdwtype), EEL_PME_EWALD(inputrec->coulombtype),
899                 domdecOptions.numPmeRanks, gpusWereDetected);
900     }
901     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
902
903     const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
904
905     // Initialize development feature flags that enabled by environment variable
906     // and report those features that are enabled.
907     const DevelopmentFeatureFlags devFlags =
908             manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
909
910     // Build restraints.
911     // TODO: hide restraint implementation details from Mdrunner.
912     // There is nothing unique about restraints at this point as far as the
913     // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
914     // factory functions from the SimulationContext on which to call mdModules_->add().
915     // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
916     for (auto&& restraint : restraintManager_->getRestraints())
917     {
918         auto module = RestraintMDModule::create(restraint, restraint->sites());
919         mdModules_->add(std::move(module));
920     }
921
922     // TODO: Error handling
923     mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
924     const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
925
926     if (inputrec->internalParameters != nullptr)
927     {
928         mdModulesNotifier.notify(*inputrec->internalParameters);
929     }
930
931     if (fplog != nullptr)
932     {
933         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
934         fprintf(fplog, "\n");
935     }
936
937     if (SIMMASTER(cr))
938     {
939         /* In rerun, set velocities to zero if present */
940         if (doRerun && ((globalState->flags & (1 << estV)) != 0))
941         {
942             // rerun does not use velocities
943             GMX_LOG(mdlog.info)
944                     .asParagraph()
945                     .appendText(
946                             "Rerun trajectory contains velocities. Rerun does only evaluate "
947                             "potential energy and forces. The velocities will be ignored.");
948             for (int i = 0; i < globalState->natoms; i++)
949             {
950                 clear_rvec(globalState->v[i]);
951             }
952             globalState->flags &= ~(1 << estV);
953         }
954
955         /* now make sure the state is initialized and propagated */
956         set_state_entries(globalState.get(), inputrec);
957     }
958
959     /* NM and TPI parallelize over force/energy calculations, not atoms,
960      * so we need to initialize and broadcast the global state.
961      */
962     if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
963     {
964         if (!MASTER(cr))
965         {
966             globalState = std::make_unique<t_state>();
967         }
968         broadcastStateWithoutDynamics(cr, globalState.get());
969     }
970
971     /* A parallel command line option consistency check that we can
972        only do after any threads have started. */
973     if (!PAR(cr)
974         && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
975             || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
976     {
977         gmx_fatal(FARGS,
978                   "The -dd or -npme option request a parallel simulation, "
979 #if !GMX_MPI
980                   "but %s was compiled without threads or MPI enabled",
981                   output_env_get_program_display_name(oenv));
982 #elif GMX_THREAD_MPI
983                   "but the number of MPI-threads (option -ntmpi) is not set or is 1");
984 #else
985                   "but %s was not started through mpirun/mpiexec or only one rank was requested "
986                   "through mpirun/mpiexec",
987                   output_env_get_program_display_name(oenv));
988 #endif
989     }
990
991     if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
992     {
993         gmx_fatal(FARGS,
994                   "The .mdp file specified an energy mininization or normal mode algorithm, and "
995                   "these are not compatible with mdrun -rerun");
996     }
997
998     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
999     {
1000         if (domdecOptions.numPmeRanks > 0)
1001         {
1002             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
1003                                  "PME-only ranks are requested, but the system does not use PME "
1004                                  "for electrostatics or LJ");
1005         }
1006
1007         domdecOptions.numPmeRanks = 0;
1008     }
1009
1010     if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1011     {
1012         /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1013          * improve performance with many threads per GPU, since our OpenMP
1014          * scaling is bad, but it's difficult to automate the setup.
1015          */
1016         domdecOptions.numPmeRanks = 0;
1017     }
1018     if (useGpuForPme)
1019     {
1020         if (domdecOptions.numPmeRanks < 0)
1021         {
1022             domdecOptions.numPmeRanks = 0;
1023             // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1024         }
1025         else
1026         {
1027             GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
1028                                "PME GPU decomposition is not supported");
1029         }
1030     }
1031
1032 #if GMX_FAHCORE
1033     if (MASTER(cr))
1034     {
1035         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1036     }
1037 #endif
1038
1039     /* NMR restraints must be initialized before load_checkpoint,
1040      * since with time averaging the history is added to t_state.
1041      * For proper consistency check we therefore need to extend
1042      * t_state here.
1043      * So the PME-only nodes (if present) will also initialize
1044      * the distance restraints.
1045      */
1046     snew(fcd, 1);
1047
1048     /* This needs to be called before read_checkpoint to extend the state */
1049     init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1050
1051     init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1052
1053     auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1054
1055     ObservablesHistory observablesHistory = {};
1056
1057     if (startingBehavior != StartingBehavior::NewSimulation)
1058     {
1059         /* Check if checkpoint file exists before doing continuation.
1060          * This way we can use identical input options for the first and subsequent runs...
1061          */
1062         if (mdrunOptions.numStepsCommandline > -2)
1063         {
1064             /* Temporarily set the number of steps to unlimited to avoid
1065              * triggering the nsteps check in load_checkpoint().
1066              * This hack will go away soon when the -nsteps option is removed.
1067              */
1068             inputrec->nsteps = -1;
1069         }
1070
1071         load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1072                         logFileHandle, cr, domdecOptions.numCells, inputrec, globalState.get(),
1073                         &observablesHistory, mdrunOptions.reproducible, mdModules_->notifier());
1074
1075         if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1076         {
1077             // Now we can start normal logging to the truncated log file.
1078             fplog = gmx_fio_getfp(logFileHandle);
1079             prepareLogAppending(fplog);
1080             logOwner = buildLogger(fplog, MASTER(cr));
1081             mdlog    = logOwner.logger();
1082         }
1083     }
1084
1085     if (mdrunOptions.numStepsCommandline > -2)
1086     {
1087         GMX_LOG(mdlog.info)
1088                 .asParagraph()
1089                 .appendText(
1090                         "The -nsteps functionality is deprecated, and may be removed in a future "
1091                         "version. "
1092                         "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1093                         "file field.");
1094     }
1095     /* override nsteps with value set on the commandline */
1096     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1097
1098     if (SIMMASTER(cr))
1099     {
1100         copy_mat(globalState->box, box);
1101     }
1102
1103     if (PAR(cr))
1104     {
1105         gmx_bcast(sizeof(box), box, cr);
1106     }
1107
1108     if (inputrec->cutoff_scheme != ecutsVERLET)
1109     {
1110         gmx_fatal(FARGS,
1111                   "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1112                   "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1113     }
1114     /* Update rlist and nstlist. */
1115     prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1116                           useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1117                           *hwinfo->cpuInfo);
1118
1119     const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1120     // This builder is necessary while we have multi-part construction
1121     // of DD. Before DD is constructed, we use the existence of
1122     // the builder object to indicate that further construction of DD
1123     // is needed.
1124     std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1125     if (useDomainDecomposition)
1126     {
1127         ddBuilder = std::make_unique<DomainDecompositionBuilder>(
1128                 mdlog, cr, domdecOptions, mdrunOptions, prefer1DAnd1PulseDD, mtop, *inputrec, box,
1129                 positionsFromStatePointer(globalState.get()));
1130     }
1131     else
1132     {
1133         /* PME, if used, is done on all nodes with 1D decomposition */
1134         cr->npmenodes = 0;
1135         cr->duty      = (DUTY_PP | DUTY_PME);
1136
1137         if (inputrec->pbcType == PbcType::Screw)
1138         {
1139             gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1140         }
1141     }
1142
1143     // Produce the task assignment for this rank.
1144     GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1145     GpuTaskAssignments        gpuTaskAssignments = gpuTaskAssignmentsBuilder.build(
1146             gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
1147             nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
1148             useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
1149             // TODO cr->duty & DUTY_PME should imply that a PME
1150             // algorithm is active, but currently does not.
1151             EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1152
1153     // Get the device handles for the modules, nullptr when no task is assigned.
1154     DeviceInformation* nonbondedDeviceInfo = gpuTaskAssignments.initNonbondedDevice(cr);
1155     DeviceInformation* pmeDeviceInfo       = gpuTaskAssignments.initPmeDevice();
1156
1157     // TODO Initialize GPU streams here.
1158
1159     // TODO Currently this is always built, yet DD partition code
1160     // checks if it is built before using it. Probably it should
1161     // become an MDModule that is made only when another module
1162     // requires it (e.g. pull, CompEl, density fitting), so that we
1163     // don't update the local atom sets unilaterally every step.
1164     LocalAtomSetManager atomSets;
1165     if (ddBuilder)
1166     {
1167         // TODO Pass the GPU streams to ddBuilder to use in buffer
1168         // transfers (e.g. halo exchange)
1169         cr->dd = ddBuilder->build(&atomSets);
1170         // The builder's job is done, so destruct it
1171         ddBuilder.reset(nullptr);
1172         // Note that local state still does not exist yet.
1173     }
1174
1175     // The GPU update is decided here because we need to know whether the constraints or
1176     // SETTLEs can span accross the domain borders (i.e. whether or not update groups are
1177     // defined). This is only known after DD is initialized, hence decision on using GPU
1178     // update is done so late.
1179     try
1180     {
1181         const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1182
1183         useGpuForUpdate = decideWhetherToUseGpuForUpdate(
1184                 devFlags.forceGpuUpdateDefault, useDomainDecomposition, useUpdateGroups, pmeRunMode,
1185                 domdecOptions.numPmeRanks > 0, useGpuForNonbonded, updateTarget, gpusWereDetected,
1186                 *inputrec, mtop, doEssentialDynamics, gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1187                 replExParams.exchangeInterval > 0, doRerun, mdlog);
1188     }
1189     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1190
1191     const bool printHostName = (cr->nnodes > 1);
1192     gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1193
1194     // If the user chose a task assignment, give them some hints
1195     // where appropriate.
1196     if (!userGpuTaskAssignment.empty())
1197     {
1198         gpuTaskAssignments.logPerformanceHints(mdlog, ssize(gpuIdsToUse));
1199     }
1200
1201     if (PAR(cr))
1202     {
1203         /* After possible communicator splitting in make_dd_communicators.
1204          * we can set up the intra/inter node communication.
1205          */
1206         gmx_setup_nodecomm(fplog, cr);
1207     }
1208
1209 #if GMX_MPI
1210     if (isMultiSim(ms))
1211     {
1212         GMX_LOG(mdlog.warning)
1213                 .asParagraph()
1214                 .appendTextFormatted(
1215                         "This is simulation %d out of %d running as a composite GROMACS\n"
1216                         "multi-simulation job. Setup for this simulation:\n",
1217                         ms->sim, ms->nsim);
1218     }
1219     GMX_LOG(mdlog.warning)
1220             .appendTextFormatted("Using %d MPI %s\n", cr->nnodes,
1221 #    if GMX_THREAD_MPI
1222                                  cr->nnodes == 1 ? "thread" : "threads"
1223 #    else
1224                                  cr->nnodes == 1 ? "process" : "processes"
1225 #    endif
1226             );
1227     fflush(stderr);
1228 #endif
1229
1230     // If mdrun -pin auto honors any affinity setting that already
1231     // exists. If so, it is nice to provide feedback about whether
1232     // that existing affinity setting was from OpenMP or something
1233     // else, so we run this code both before and after we initialize
1234     // the OpenMP support.
1235     gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1236     /* Check and update the number of OpenMP threads requested */
1237     checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1238                                             pmeRunMode, mtop, *inputrec);
1239
1240     gmx_omp_nthreads_init(mdlog, cr, hwinfo->nthreads_hw_avail, physicalNodeComm.size_,
1241                           hw_opt.nthreads_omp, hw_opt.nthreads_omp_pme, !thisRankHasDuty(cr, DUTY_PP));
1242
1243     // Enable FP exception detection, but not in
1244     // Release mode and not for compilers with known buggy FP
1245     // exception support (clang with any optimization) or suspected
1246     // buggy FP exception support (gcc 7.* with optimization).
1247 #if !defined NDEBUG                                                                         \
1248         && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1249              && defined __OPTIMIZE__)
1250     const bool bEnableFPE = true;
1251 #else
1252     const bool bEnableFPE = false;
1253 #endif
1254     // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1255     if (bEnableFPE)
1256     {
1257         gmx_feenableexcept();
1258     }
1259
1260     /* Now that we know the setup is consistent, check for efficiency */
1261     check_resource_division_efficiency(hwinfo, gpuTaskAssignments.thisRankHasAnyGpuTask(),
1262                                        mdrunOptions.ntompOptionIsSet, cr, mdlog);
1263
1264     /* getting number of PP/PME threads on this MPI / tMPI rank.
1265        PME: env variable should be read only on one node to make sure it is
1266        identical everywhere;
1267      */
1268     const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded)
1269                                                                   : gmx_omp_nthreads_get(emntPME);
1270     checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid, *hwinfo->hardwareTopology,
1271                                   physicalNodeComm, mdlog);
1272
1273     // Enable Peer access between GPUs where available
1274     // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1275     // any of the GPU communication features are active.
1276     if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1277         && (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1278     {
1279         setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1280     }
1281
1282     if (hw_opt.threadAffinity != ThreadAffinity::Off)
1283     {
1284         /* Before setting affinity, check whether the affinity has changed
1285          * - which indicates that probably the OpenMP library has changed it
1286          * since we first checked).
1287          */
1288         gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1289
1290         int numThreadsOnThisNode, intraNodeThreadOffset;
1291         analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1292                                  &intraNodeThreadOffset);
1293
1294         /* Set the CPU affinity */
1295         gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology, numThreadsOnThisRank,
1296                                 numThreadsOnThisNode, intraNodeThreadOffset, nullptr);
1297     }
1298
1299     if (mdrunOptions.timingOptions.resetStep > -1)
1300     {
1301         GMX_LOG(mdlog.info)
1302                 .asParagraph()
1303                 .appendText(
1304                         "The -resetstep functionality is deprecated, and may be removed in a "
1305                         "future version.");
1306     }
1307     wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1308
1309     if (PAR(cr))
1310     {
1311         /* Master synchronizes its value of reset_counters with all nodes
1312          * including PME only nodes */
1313         int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1314         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1315         wcycle_set_reset_counters(wcycle, reset_counters);
1316     }
1317
1318     // Membrane embedding must be initialized before we call init_forcerec()
1319     if (doMembed)
1320     {
1321         if (MASTER(cr))
1322         {
1323             fprintf(stderr, "Initializing membed");
1324         }
1325         /* Note that membed cannot work in parallel because mtop is
1326          * changed here. Fix this if we ever want to make it run with
1327          * multiple ranks. */
1328         membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec,
1329                              globalState.get(), cr, &mdrunOptions.checkpointOptions.period);
1330     }
1331
1332     const bool                   thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1333     std::unique_ptr<MDAtoms>     mdAtoms;
1334     std::unique_ptr<gmx_vsite_t> vsite;
1335     std::unique_ptr<GpuBonded>   gpuBonded;
1336
1337     t_nrnb nrnb;
1338     if (thisRankHasDuty(cr, DUTY_PP))
1339     {
1340         mdModulesNotifier.notify(*cr);
1341         mdModulesNotifier.notify(&atomSets);
1342         mdModulesNotifier.notify(inputrec->pbcType);
1343         mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1344         /* Initiate forcerecord */
1345         fr                 = new t_forcerec;
1346         fr->forceProviders = mdModules_->initForceProviders();
1347         init_forcerec(fplog, mdlog, fr, fcd, inputrec, &mtop, cr, box,
1348                       opt2fn("-table", filenames.size(), filenames.data()),
1349                       opt2fn("-tablep", filenames.size(), filenames.data()),
1350                       opt2fns("-tableb", filenames.size(), filenames.data()),
1351                       pmeRunMode == PmeRunMode::GPU && !thisRankHasDuty(cr, DUTY_PME), pforce);
1352
1353         fr->nbv = Nbnxm::init_nb_verlet(mdlog, inputrec, fr, cr, *hwinfo, nonbondedDeviceInfo,
1354                                         &mtop, box, wcycle);
1355         if (useGpuForBonded)
1356         {
1357             auto stream = havePPDomainDecomposition(cr)
1358                                   ? Nbnxm::gpu_get_command_stream(
1359                                             fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1360                                   : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1361                                                                   gmx::InteractionLocality::Local);
1362             gpuBonded     = std::make_unique<GpuBonded>(mtop.ffparams, stream, wcycle);
1363             fr->gpuBonded = gpuBonded.get();
1364         }
1365
1366         /* Initialize the mdAtoms structure.
1367          * mdAtoms is not filled with atom data,
1368          * as this can not be done now with domain decomposition.
1369          */
1370         mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1371         if (globalState && thisRankHasPmeGpuTask)
1372         {
1373             // The pinning of coordinates in the global state object works, because we only use
1374             // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1375             // points to the global state object without DD.
1376             // FIXME: MD and EM separately set up the local state - this should happen in the same
1377             // function, which should also perform the pinning.
1378             changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1379         }
1380
1381         /* Initialize the virtual site communication */
1382         vsite = initVsite(mtop, cr);
1383
1384         calc_shifts(box, fr->shift_vec);
1385
1386         /* With periodic molecules the charge groups should be whole at start up
1387          * and the virtual sites should not be far from their proper positions.
1388          */
1389         if (!inputrec->bContinuation && MASTER(cr)
1390             && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1391         {
1392             /* Make molecules whole at start of run */
1393             if (fr->pbcType != PbcType::No)
1394             {
1395                 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1396             }
1397             if (vsite)
1398             {
1399                 /* Correct initial vsite positions are required
1400                  * for the initial distribution in the domain decomposition
1401                  * and for the initial shell prediction.
1402                  */
1403                 constructVsitesGlobal(mtop, globalState->x);
1404             }
1405         }
1406
1407         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1408         {
1409             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1410             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1411         }
1412     }
1413     else
1414     {
1415         /* This is a PME only node */
1416
1417         GMX_ASSERT(globalState == nullptr,
1418                    "We don't need the state on a PME only rank and expect it to be unitialized");
1419
1420         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1421         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1422     }
1423
1424     gmx_pme_t* sepPmeData = nullptr;
1425     // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1426     GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1427                "Double-checking that only PME-only ranks have no forcerec");
1428     gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1429
1430     // TODO should live in ewald module once its testing is improved
1431     //
1432     // Later, this program could contain kernels that might be later
1433     // re-used as auto-tuning progresses, or subsequent simulations
1434     // are invoked.
1435     PmeGpuProgramStorage pmeGpuProgram;
1436     if (thisRankHasPmeGpuTask)
1437     {
1438         pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1439     }
1440
1441     /* Initiate PME if necessary,
1442      * either on all nodes or on dedicated PME nodes only. */
1443     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1444     {
1445         if (mdAtoms && mdAtoms->mdatoms())
1446         {
1447             nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1448             if (EVDW_PME(inputrec->vdwtype))
1449             {
1450                 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1451             }
1452         }
1453         if (cr->npmenodes > 0)
1454         {
1455             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1456             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1457             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1458         }
1459
1460         if (thisRankHasDuty(cr, DUTY_PME))
1461         {
1462             try
1463             {
1464                 pmedata = gmx_pme_init(cr, getNumPmeDomains(cr->dd), inputrec, nChargePerturbed != 0,
1465                                        nTypePerturbed != 0, mdrunOptions.reproducible, ewaldcoeff_q,
1466                                        ewaldcoeff_lj, gmx_omp_nthreads_get(emntPME), pmeRunMode,
1467                                        nullptr, pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1468             }
1469             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1470         }
1471     }
1472
1473
1474     if (EI_DYNAMICS(inputrec->eI))
1475     {
1476         /* Turn on signal handling on all nodes */
1477         /*
1478          * (A user signal from the PME nodes (if any)
1479          * is communicated to the PP nodes.
1480          */
1481         signal_handler_install();
1482     }
1483
1484     pull_t* pull_work = nullptr;
1485     if (thisRankHasDuty(cr, DUTY_PP))
1486     {
1487         /* Assumes uniform use of the number of OpenMP threads */
1488         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1489
1490         if (inputrec->bPull)
1491         {
1492             /* Initialize pull code */
1493             pull_work = init_pull(fplog, inputrec->pull, inputrec, &mtop, cr, &atomSets,
1494                                   inputrec->fepvals->init_lambda);
1495             if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1496             {
1497                 initPullHistory(pull_work, &observablesHistory);
1498             }
1499             if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1500             {
1501                 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1502             }
1503         }
1504
1505         std::unique_ptr<EnforcedRotation> enforcedRotation;
1506         if (inputrec->bRot)
1507         {
1508             /* Initialize enforced rotation code */
1509             enforcedRotation =
1510                     init_rot(fplog, inputrec, filenames.size(), filenames.data(), cr, &atomSets,
1511                              globalState.get(), &mtop, oenv, mdrunOptions, startingBehavior);
1512         }
1513
1514         t_swap* swap = nullptr;
1515         if (inputrec->eSwapCoords != eswapNO)
1516         {
1517             /* Initialize ion swapping code */
1518             swap = init_swapcoords(fplog, inputrec,
1519                                    opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1520                                    &mtop, globalState.get(), &observablesHistory, cr, &atomSets,
1521                                    oenv, mdrunOptions, startingBehavior);
1522         }
1523
1524         /* Let makeConstraints know whether we have essential dynamics constraints. */
1525         auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics, fplog,
1526                                       *mdAtoms->mdatoms(), cr, ms, &nrnb, wcycle, fr->bMolPBC);
1527
1528         /* Energy terms and groups */
1529         gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1530                              inputrec->fepvals->n_lambda);
1531
1532         /* Kinetic energy data */
1533         gmx_ekindata_t ekind;
1534         init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1535
1536         /* Set up interactive MD (IMD) */
1537         auto imdSession =
1538                 makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1539                                MASTER(cr) ? globalState->x.rvec_array() : nullptr, filenames.size(),
1540                                filenames.data(), oenv, mdrunOptions.imdOptions, startingBehavior);
1541
1542         if (DOMAINDECOMP(cr))
1543         {
1544             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1545             /* This call is not included in init_domain_decomposition mainly
1546              * because fr->cginfo_mb is set later.
1547              */
1548             dd_init_bondeds(fplog, cr->dd, mtop, vsite.get(), inputrec,
1549                             domdecOptions.checkBondedInteractions, fr->cginfo_mb);
1550         }
1551
1552         const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
1553                 false, inputrec, doRerun, vsite.get(), ms, replExParams, fcd,
1554                 static_cast<int>(filenames.size()), filenames.data(), &observablesHistory, membed);
1555
1556         const bool useModularSimulator = inputIsCompatibleWithModularSimulator
1557                                          && !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
1558
1559         // TODO This is not the right place to manage the lifetime of
1560         // this data structure, but currently it's the easiest way to
1561         // make it work.
1562         MdrunScheduleWorkload runScheduleWork;
1563         // Also populates the simulation constant workload description.
1564         runScheduleWork.simulationWork =
1565                 createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
1566                                          useGpuForUpdate, devFlags.enableGpuBufferOps,
1567                                          devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
1568
1569         std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1570         if (gpusWereDetected
1571             && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME))
1572                 || runScheduleWork.simulationWork.useGpuBufferOps))
1573         {
1574             const void* pmeStream = pme_gpu_get_device_stream(fr->pmedata);
1575             const void* localStream =
1576                     fr->nbv->gpu_nbv != nullptr
1577                             ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local)
1578                             : nullptr;
1579             const void* nonLocalStream =
1580                     fr->nbv->gpu_nbv != nullptr
1581                             ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal)
1582                             : nullptr;
1583             const void*        deviceContext = pme_gpu_get_device_context(fr->pmedata);
1584             const int          paddingSize   = pme_gpu_get_padding_size(fr->pmedata);
1585             GpuApiCallBehavior transferKind = (inputrec->eI == eiMD && !doRerun && !useModularSimulator)
1586                                                       ? GpuApiCallBehavior::Async
1587                                                       : GpuApiCallBehavior::Sync;
1588
1589             stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
1590                     pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize, wcycle);
1591             fr->stateGpu = stateGpu.get();
1592         }
1593
1594         GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1595         SimulatorBuilder simulatorBuilder;
1596
1597         // build and run simulator object based on user-input
1598         auto simulator = simulatorBuilder.build(
1599                 inputIsCompatibleWithModularSimulator, fplog, cr, ms, mdlog,
1600                 static_cast<int>(filenames.size()), filenames.data(), oenv, mdrunOptions,
1601                 startingBehavior, vsite.get(), constr.get(),
1602                 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, deform.get(),
1603                 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
1604                 pull_work, swap, &mtop, fcd, globalState.get(), &observablesHistory, mdAtoms.get(),
1605                 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
1606                 walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
1607         simulator->run();
1608
1609         if (fr->pmePpCommGpu)
1610         {
1611             // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
1612             fr->pmePpCommGpu.reset();
1613         }
1614
1615         if (inputrec->bPull)
1616         {
1617             finish_pull(pull_work);
1618         }
1619         finish_swapcoords(swap);
1620     }
1621     else
1622     {
1623         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1624         /* do PME only */
1625         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1626         gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1627     }
1628
1629     wallcycle_stop(wcycle, ewcRUN);
1630
1631     /* Finish up, write some stuff
1632      * if rerunMD, don't write last frame again
1633      */
1634     finish_run(fplog, mdlog, cr, inputrec, &nrnb, wcycle, walltime_accounting,
1635                fr ? fr->nbv.get() : nullptr, pmedata, EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1636
1637     // clean up cycle counter
1638     wallcycle_destroy(wcycle);
1639
1640     // Free PME data
1641     if (pmedata)
1642     {
1643         gmx_pme_destroy(pmedata);
1644         pmedata = nullptr;
1645     }
1646
1647     // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1648     // before we destroy the GPU context(s) in free_gpu().
1649     // Pinned buffers are associated with contexts in CUDA.
1650     // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1651     mdAtoms.reset(nullptr);
1652     globalState.reset(nullptr);
1653     mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1654     gpuBonded.reset(nullptr);
1655     /* Free pinned buffers in *fr */
1656     delete fr;
1657     fr = nullptr;
1658
1659     if (hwinfo->gpu_info.n_dev > 0)
1660     {
1661         /* stop the GPU profiler (only CUDA) */
1662         stopGpuProfiler();
1663     }
1664
1665     /* With tMPI we need to wait for all ranks to finish deallocation before
1666      * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1667      * GPU and context.
1668      *
1669      * This is not a concern in OpenCL where we use one context per rank which
1670      * is freed in nbnxn_gpu_free().
1671      *
1672      * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1673      * but it is easier and more futureproof to call it on the whole node.
1674      *
1675      * Note that this function needs to be called even if GPUs are not used
1676      * in this run because the PME ranks have no knowledge of whether GPUs
1677      * are used or not, but all ranks need to enter the barrier below.
1678      * \todo Remove this physical node barrier after making sure
1679      * that it's not needed anymore (with a shared GPU run).
1680      */
1681     if (GMX_THREAD_MPI)
1682     {
1683         physicalNodeComm.barrier();
1684     }
1685
1686     free_gpu(nonbondedDeviceInfo);
1687     free_gpu(pmeDeviceInfo);
1688     sfree(fcd);
1689
1690     if (doMembed)
1691     {
1692         free_membed(membed);
1693     }
1694
1695     /* Does what it says */
1696     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1697     walltime_accounting_destroy(walltime_accounting);
1698
1699     // Ensure log file content is written
1700     if (logFileHandle)
1701     {
1702         gmx_fio_flush(logFileHandle);
1703     }
1704
1705     /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1706      * exceptions were enabled before function was called. */
1707     if (bEnableFPE)
1708     {
1709         gmx_fedisableexcept();
1710     }
1711
1712     auto rc = static_cast<int>(gmx_get_stop_condition());
1713
1714 #if GMX_THREAD_MPI
1715     /* we need to join all threads. The sub-threads join when they
1716        exit this function, but the master thread needs to be told to
1717        wait for that. */
1718     if (PAR(cr) && MASTER(cr))
1719     {
1720         tMPI_Finalize();
1721     }
1722 #endif
1723     return rc;
1724 }
1725
1726 Mdrunner::~Mdrunner()
1727 {
1728     // Clean up of the Manager.
1729     // This will end up getting called on every thread-MPI rank, which is unnecessary,
1730     // but okay as long as threads synchronize some time before adding or accessing
1731     // a new set of restraints.
1732     if (restraintManager_)
1733     {
1734         restraintManager_->clear();
1735         GMX_ASSERT(restraintManager_->countRestraints() == 0,
1736                    "restraints added during runner life time should be cleared at runner "
1737                    "destruction.");
1738     }
1739 };
1740
1741 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
1742 {
1743     GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1744     // Not sure if this should be logged through the md logger or something else,
1745     // but it is helpful to have some sort of INFO level message sent somewhere.
1746     //    std::cout << "Registering restraint named " << name << std::endl;
1747
1748     // When multiple restraints are used, it may be wasteful to register them separately.
1749     // Maybe instead register an entire Restraint Manager as a force provider.
1750     restraintManager_->addToSpec(std::move(puller), name);
1751 }
1752
1753 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
1754
1755 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
1756
1757 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1758 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1759
1760 class Mdrunner::BuilderImplementation
1761 {
1762 public:
1763     BuilderImplementation() = delete;
1764     BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
1765     ~BuilderImplementation();
1766
1767     BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
1768                                                 real                forceWarningThreshold,
1769                                                 StartingBehavior    startingBehavior);
1770
1771     void addDomdec(const DomdecOptions& options);
1772
1773     void addVerletList(int nstlist);
1774
1775     void addReplicaExchange(const ReplicaExchangeParameters& params);
1776
1777     void addNonBonded(const char* nbpu_opt);
1778
1779     void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1780
1781     void addBondedTaskAssignment(const char* bonded_opt);
1782
1783     void addUpdateTaskAssignment(const char* update_opt);
1784
1785     void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
1786
1787     void addFilenames(ArrayRef<const t_filenm> filenames);
1788
1789     void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1790
1791     void addLogFile(t_fileio* logFileHandle);
1792
1793     void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1794
1795     Mdrunner build();
1796
1797 private:
1798     // Default parameters copied from runner.h
1799     // \todo Clarify source(s) of default parameters.
1800
1801     const char* nbpu_opt_    = nullptr;
1802     const char* pme_opt_     = nullptr;
1803     const char* pme_fft_opt_ = nullptr;
1804     const char* bonded_opt_  = nullptr;
1805     const char* update_opt_  = nullptr;
1806
1807     MdrunOptions mdrunOptions_;
1808
1809     DomdecOptions domdecOptions_;
1810
1811     ReplicaExchangeParameters replicaExchangeParameters_;
1812
1813     //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1814     int nstlist_ = 0;
1815
1816     //! Multisim communicator handle.
1817     gmx_multisim_t* multiSimulation_;
1818
1819     //! mdrun communicator
1820     MPI_Comm communicator_ = MPI_COMM_NULL;
1821
1822     //! Print a warning if any force is larger than this (in kJ/mol nm).
1823     real forceWarningThreshold_ = -1;
1824
1825     //! Whether the simulation will start afresh, or restart with/without appending.
1826     StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1827
1828     //! The modules that comprise the functionality of mdrun.
1829     std::unique_ptr<MDModules> mdModules_;
1830
1831     //! \brief Parallelism information.
1832     gmx_hw_opt_t hardwareOptions_;
1833
1834     //! filename options for simulation.
1835     ArrayRef<const t_filenm> filenames_;
1836
1837     /*! \brief Handle to output environment.
1838      *
1839      * \todo gmx_output_env_t needs lifetime management.
1840      */
1841     gmx_output_env_t* outputEnvironment_ = nullptr;
1842
1843     /*! \brief Non-owning handle to MD log file.
1844      *
1845      * \todo Context should own output facilities for client.
1846      * \todo Improve log file handle management.
1847      * \internal
1848      * Code managing the FILE* relies on the ability to set it to
1849      * nullptr to check whether the filehandle is valid.
1850      */
1851     t_fileio* logFileHandle_ = nullptr;
1852
1853     /*!
1854      * \brief Builder for simulation stop signal handler.
1855      */
1856     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1857 };
1858
1859 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
1860                                                        compat::not_null<SimulationContext*> context) :
1861     mdModules_(std::move(mdModules))
1862 {
1863     communicator_    = context->communicator_;
1864     multiSimulation_ = context->multiSimulation_.get();
1865 }
1866
1867 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1868
1869 Mdrunner::BuilderImplementation&
1870 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions&    options,
1871                                                       const real             forceWarningThreshold,
1872                                                       const StartingBehavior startingBehavior)
1873 {
1874     mdrunOptions_          = options;
1875     forceWarningThreshold_ = forceWarningThreshold;
1876     startingBehavior_      = startingBehavior;
1877     return *this;
1878 }
1879
1880 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
1881 {
1882     domdecOptions_ = options;
1883 }
1884
1885 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1886 {
1887     nstlist_ = nstlist;
1888 }
1889
1890 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
1891 {
1892     replicaExchangeParameters_ = params;
1893 }
1894
1895 Mdrunner Mdrunner::BuilderImplementation::build()
1896 {
1897     auto newRunner = Mdrunner(std::move(mdModules_));
1898
1899     newRunner.mdrunOptions     = mdrunOptions_;
1900     newRunner.pforce           = forceWarningThreshold_;
1901     newRunner.startingBehavior = startingBehavior_;
1902     newRunner.domdecOptions    = domdecOptions_;
1903
1904     // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1905     newRunner.hw_opt = hardwareOptions_;
1906
1907     // No invariant to check. This parameter exists to optionally override other behavior.
1908     newRunner.nstlist_cmdline = nstlist_;
1909
1910     newRunner.replExParams = replicaExchangeParameters_;
1911
1912     newRunner.filenames = filenames_;
1913
1914     newRunner.communicator = communicator_;
1915
1916     // nullptr is a valid value for the multisim handle
1917     newRunner.ms = multiSimulation_;
1918
1919     // \todo Clarify ownership and lifetime management for gmx_output_env_t
1920     // \todo Update sanity checking when output environment has clearly specified invariants.
1921     // Initialization and default values for oenv are not well specified in the current version.
1922     if (outputEnvironment_)
1923     {
1924         newRunner.oenv = outputEnvironment_;
1925     }
1926     else
1927     {
1928         GMX_THROW(gmx::APIError(
1929                 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1930     }
1931
1932     newRunner.logFileHandle = logFileHandle_;
1933
1934     if (nbpu_opt_)
1935     {
1936         newRunner.nbpu_opt = nbpu_opt_;
1937     }
1938     else
1939     {
1940         GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1941     }
1942
1943     if (pme_opt_ && pme_fft_opt_)
1944     {
1945         newRunner.pme_opt     = pme_opt_;
1946         newRunner.pme_fft_opt = pme_fft_opt_;
1947     }
1948     else
1949     {
1950         GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1951     }
1952
1953     if (bonded_opt_)
1954     {
1955         newRunner.bonded_opt = bonded_opt_;
1956     }
1957     else
1958     {
1959         GMX_THROW(gmx::APIError(
1960                 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1961     }
1962
1963     if (update_opt_)
1964     {
1965         newRunner.update_opt = update_opt_;
1966     }
1967     else
1968     {
1969         GMX_THROW(gmx::APIError(
1970                 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build()  "));
1971     }
1972
1973
1974     newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1975
1976     if (stopHandlerBuilder_)
1977     {
1978         newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1979     }
1980     else
1981     {
1982         newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1983     }
1984
1985     return newRunner;
1986 }
1987
1988 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1989 {
1990     nbpu_opt_ = nbpu_opt;
1991 }
1992
1993 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
1994 {
1995     pme_opt_     = pme_opt;
1996     pme_fft_opt_ = pme_fft_opt;
1997 }
1998
1999 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2000 {
2001     bonded_opt_ = bonded_opt;
2002 }
2003
2004 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2005 {
2006     update_opt_ = update_opt;
2007 }
2008
2009 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2010 {
2011     hardwareOptions_ = hardwareOptions;
2012 }
2013
2014 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2015 {
2016     filenames_ = filenames;
2017 }
2018
2019 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2020 {
2021     outputEnvironment_ = outputEnvironment;
2022 }
2023
2024 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2025 {
2026     logFileHandle_ = logFileHandle;
2027 }
2028
2029 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2030 {
2031     stopHandlerBuilder_ = std::move(builder);
2032 }
2033
2034 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules>           mdModules,
2035                                  compat::not_null<SimulationContext*> context) :
2036     impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2037 {
2038 }
2039
2040 MdrunnerBuilder::~MdrunnerBuilder() = default;
2041
2042 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions&    options,
2043                                                       real                   forceWarningThreshold,
2044                                                       const StartingBehavior startingBehavior)
2045 {
2046     impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2047     return *this;
2048 }
2049
2050 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2051 {
2052     impl_->addDomdec(options);
2053     return *this;
2054 }
2055
2056 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2057 {
2058     impl_->addVerletList(nstlist);
2059     return *this;
2060 }
2061
2062 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2063 {
2064     impl_->addReplicaExchange(params);
2065     return *this;
2066 }
2067
2068 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2069 {
2070     impl_->addNonBonded(nbpu_opt);
2071     return *this;
2072 }
2073
2074 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2075 {
2076     // The builder method may become more general in the future, but in this version,
2077     // parameters for PME electrostatics are both required and the only parameters
2078     // available.
2079     if (pme_opt && pme_fft_opt)
2080     {
2081         impl_->addPME(pme_opt, pme_fft_opt);
2082     }
2083     else
2084     {
2085         GMX_THROW(
2086                 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2087     }
2088     return *this;
2089 }
2090
2091 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2092 {
2093     impl_->addBondedTaskAssignment(bonded_opt);
2094     return *this;
2095 }
2096
2097 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2098 {
2099     impl_->addUpdateTaskAssignment(update_opt);
2100     return *this;
2101 }
2102
2103 Mdrunner MdrunnerBuilder::build()
2104 {
2105     return impl_->build();
2106 }
2107
2108 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2109 {
2110     impl_->addHardwareOptions(hardwareOptions);
2111     return *this;
2112 }
2113
2114 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2115 {
2116     impl_->addFilenames(filenames);
2117     return *this;
2118 }
2119
2120 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2121 {
2122     impl_->addOutputEnvironment(outputEnvironment);
2123     return *this;
2124 }
2125
2126 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2127 {
2128     impl_->addLogFile(logFileHandle);
2129     return *this;
2130 }
2131
2132 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2133 {
2134     impl_->addStopHandlerBuilder(std::move(builder));
2135     return *this;
2136 }
2137
2138 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2139
2140 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;
2141
2142 } // namespace gmx