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