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