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