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