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