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