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