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