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