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