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