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