Replace DOMAINDECOMP macro by a renamed function
[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                                                                      numAvailableDevices,
911                                                                      userGpuTaskAssignment,
912                                                                      *hwinfo_,
913                                                                      *inputrec,
914                                                                      hw_opt.nthreads_tmpi,
915                                                                      domdecOptions.numPmeRanks);
916         }
917         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
918
919         /* Determine how many thread-MPI ranks to start.
920          *
921          * TODO Over-writing the user-supplied value here does
922          * prevent any possible subsequent checks from working
923          * correctly. */
924         hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo_,
925                                                 &hw_opt,
926                                                 numAvailableDevices,
927                                                 useGpuForNonbonded,
928                                                 useGpuForPme,
929                                                 inputrec.get(),
930                                                 mtop,
931                                                 mdlog,
932                                                 membedHolder.doMembed());
933
934         // Now start the threads for thread MPI.
935         spawnThreads(hw_opt.nthreads_tmpi);
936         // The spawned threads enter mdrunner() and execution of
937         // master and spawned threads joins at the end of this block.
938     }
939
940     GMX_RELEASE_ASSERT(ms || simulationCommunicator != MPI_COMM_NULL,
941                        "Must have valid communicator unless running a multi-simulation");
942     CommrecHandle crHandle = init_commrec(simulationCommunicator);
943     t_commrec*    cr       = crHandle.get();
944     GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
945
946     PhysicalNodeCommunicator physicalNodeComm(libraryWorldCommunicator, gmx_physicalnode_id_hash());
947
948     // If we detected the topology on this system, double-check that it makes sense
949     if (hwinfo_->hardwareTopology->isThisSystem())
950     {
951         hardwareTopologyDoubleCheckDetection(mdlog, *hwinfo_->hardwareTopology);
952     }
953
954     if (PAR(cr))
955     {
956         /* now broadcast everything to the non-master nodes/threads: */
957         if (!isSimulationMasterRank)
958         {
959             // Until now, only the master rank has a non-null pointer.
960             // On non-master ranks, allocate the object that will receive data in the following call.
961             inputrec = std::make_unique<t_inputrec>();
962         }
963         init_parallel(cr->mpiDefaultCommunicator,
964                       MASTER(cr),
965                       inputrec.get(),
966                       &mtop,
967                       partialDeserializedTpr.get());
968     }
969     GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
970     partialDeserializedTpr.reset(nullptr);
971
972     GMX_RELEASE_ASSERT(
973             !inputrec->useConstantAcceleration,
974             "Linear acceleration has been removed in GROMACS 2022, and was broken for many years "
975             "before that. Use GROMACS 4.5 or earlier if you need this feature.");
976
977     // Now we decide whether to use the domain decomposition machinery.
978     // Note that this does not necessarily imply actually using multiple domains.
979     // Now the number of ranks is known to all ranks, and each knows
980     // the inputrec read by the master rank. The ranks can now all run
981     // the task-deciding functions and will agree on the result
982     // without needing to communicate.
983     // The LBFGS minimizer, test-particle insertion, normal modes and shell dynamics don't support DD
984     const bool useDomainDecomposition =
985             !(inputrec->eI == IntegrationAlgorithm::LBFGS || EI_TPI(inputrec->eI)
986               || inputrec->eI == IntegrationAlgorithm::NM
987               || gmx_mtop_particletype_count(mtop)[ParticleType::Shell] > 0);
988
989     // Note that these variables describe only their own node.
990     //
991     // Note that when bonded interactions run on a GPU they always run
992     // alongside a nonbonded task, so do not influence task assignment
993     // even though they affect the force calculation workload.
994     bool useGpuForNonbonded = false;
995     bool useGpuForPme       = false;
996     bool useGpuForBonded    = false;
997     bool useGpuForUpdate    = false;
998     bool gpusWereDetected   = hwinfo_->ngpu_compatible_tot > 0;
999     try
1000     {
1001         // It's possible that there are different numbers of GPUs on
1002         // different nodes, which is the user's responsibility to
1003         // handle. If unsuitable, we will notice that during task
1004         // assignment.
1005         auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
1006         useGpuForNonbonded         = decideWhetherToUseGpusForNonbonded(
1007                 nonbondedTarget,
1008                 userGpuTaskAssignment,
1009                 emulateGpuNonbonded,
1010                 canUseGpuForNonbonded,
1011                 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI),
1012                 gpusWereDetected);
1013         useGpuForPme    = decideWhetherToUseGpusForPme(useGpuForNonbonded,
1014                                                     pmeTarget,
1015                                                     userGpuTaskAssignment,
1016                                                     *hwinfo_,
1017                                                     *inputrec,
1018                                                     cr->sizeOfDefaultCommunicator,
1019                                                     domdecOptions.numPmeRanks,
1020                                                     gpusWereDetected);
1021         useGpuForBonded = decideWhetherToUseGpusForBonded(
1022                 useGpuForNonbonded, useGpuForPme, bondedTarget, *inputrec, mtop, domdecOptions.numPmeRanks, gpusWereDetected);
1023     }
1024     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1025
1026     const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
1027
1028     // Initialize development feature flags that enabled by environment variable
1029     // and report those features that are enabled.
1030     const DevelopmentFeatureFlags devFlags =
1031             manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
1032
1033     const bool useModularSimulator = checkUseModularSimulator(false,
1034                                                               inputrec.get(),
1035                                                               doRerun,
1036                                                               mtop,
1037                                                               ms,
1038                                                               replExParams,
1039                                                               nullptr,
1040                                                               doEssentialDynamics,
1041                                                               membedHolder.doMembed());
1042
1043     ObservablesReducerBuilder observablesReducerBuilder;
1044
1045     // Build restraints.
1046     // TODO: hide restraint implementation details from Mdrunner.
1047     // There is nothing unique about restraints at this point as far as the
1048     // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
1049     // factory functions from the SimulationContext on which to call mdModules_->add().
1050     // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
1051     for (auto&& restraint : restraintManager_->getRestraints())
1052     {
1053         auto module = RestraintMDModule::create(restraint, restraint->sites());
1054         mdModules_->add(std::move(module));
1055     }
1056
1057     // TODO: Error handling
1058     mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
1059     // now that the MDModules know their options, they know which callbacks to sign up to
1060     mdModules_->subscribeToSimulationSetupNotifications();
1061     const auto& setupNotifier = mdModules_->notifiers().simulationSetupNotifier_;
1062
1063     // Notify MdModules of existing logger
1064     setupNotifier.notify(mdlog);
1065
1066     // Notify MdModules of internal parameters, saved into KVT
1067     if (inputrec->internalParameters != nullptr)
1068     {
1069         setupNotifier.notify(*inputrec->internalParameters);
1070     }
1071
1072     // Let MdModules know the .tpr filename
1073     {
1074         gmx::MdRunInputFilename mdRunInputFilename = { ftp2fn(efTPR, filenames.size(), filenames.data()) };
1075         setupNotifier.notify(mdRunInputFilename);
1076     }
1077
1078     if (fplog != nullptr)
1079     {
1080         pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
1081         fprintf(fplog, "\n");
1082     }
1083
1084     if (SIMMASTER(cr))
1085     {
1086         /* In rerun, set velocities to zero if present */
1087         if (doRerun && ((globalState->flags & enumValueToBitMask(StateEntry::V)) != 0))
1088         {
1089             // rerun does not use velocities
1090             GMX_LOG(mdlog.info)
1091                     .asParagraph()
1092                     .appendText(
1093                             "Rerun trajectory contains velocities. Rerun does only evaluate "
1094                             "potential energy and forces. The velocities will be ignored.");
1095             for (int i = 0; i < globalState->natoms; i++)
1096             {
1097                 clear_rvec(globalState->v[i]);
1098             }
1099             globalState->flags &= ~enumValueToBitMask(StateEntry::V);
1100         }
1101
1102         /* now make sure the state is initialized and propagated */
1103         set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
1104     }
1105
1106     /* NM and TPI parallelize over force/energy calculations, not atoms,
1107      * so we need to initialize and broadcast the global state.
1108      */
1109     if (inputrec->eI == IntegrationAlgorithm::NM || inputrec->eI == IntegrationAlgorithm::TPI)
1110     {
1111         if (!MASTER(cr))
1112         {
1113             globalState = std::make_unique<t_state>();
1114         }
1115         broadcastStateWithoutDynamics(
1116                 cr->mpiDefaultCommunicator, haveDDAtomOrdering(*cr), PAR(cr), globalState.get());
1117     }
1118
1119     /* A parallel command line option consistency check that we can
1120        only do after any threads have started. */
1121     if (!PAR(cr)
1122         && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
1123             || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
1124     {
1125         gmx_fatal(FARGS,
1126                   "The -dd or -npme option request a parallel simulation, "
1127 #if !GMX_MPI
1128                   "but %s was compiled without threads or MPI enabled",
1129                   output_env_get_program_display_name(oenv));
1130 #elif GMX_THREAD_MPI
1131                   "but the number of MPI-threads (option -ntmpi) is not set or is 1");
1132 #else
1133                   "but %s was not started through mpirun/mpiexec or only one rank was requested "
1134                   "through mpirun/mpiexec",
1135                   output_env_get_program_display_name(oenv));
1136 #endif
1137     }
1138
1139     if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || IntegrationAlgorithm::NM == inputrec->eI))
1140     {
1141         gmx_fatal(FARGS,
1142                   "The .mdp file specified an energy mininization or normal mode algorithm, and "
1143                   "these are not compatible with mdrun -rerun");
1144     }
1145
1146     /* NMR restraints must be initialized before load_checkpoint,
1147      * since with time averaging the history is added to t_state.
1148      * For proper consistency check we therefore need to extend
1149      * t_state here.
1150      * So the PME-only nodes (if present) will also initialize
1151      * the distance restraints.
1152      */
1153
1154     /* This needs to be called before read_checkpoint to extend the state */
1155     t_disresdata* disresdata;
1156     snew(disresdata, 1);
1157     init_disres(fplog,
1158                 mtop,
1159                 inputrec.get(),
1160                 DisResRunMode::MDRun,
1161                 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1162                 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1163                 cr->mpi_comm_mysim,
1164                 ms,
1165                 disresdata,
1166                 globalState.get(),
1167                 replExParams.exchangeInterval > 0);
1168
1169     if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0 && isSimulationMasterRank)
1170     {
1171         extendStateWithOriresHistory(mtop, *inputrec, globalState.get());
1172     }
1173
1174     auto deform = prepareBoxDeformation(globalState != nullptr ? globalState->box : box,
1175                                         MASTER(cr) ? DDRole::Master : DDRole::Agent,
1176                                         PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1177                                         cr->mpi_comm_mygroup,
1178                                         *inputrec);
1179
1180 #if GMX_FAHCORE
1181     /* We have to remember the generation's first step before reading checkpoint.
1182        This way, we can report to the F@H core both the generation's first step
1183        and the restored first step, thus making it able to distinguish between
1184        an interruption/resume and start of the n-th generation simulation.
1185        Having this information, the F@H core can correctly calculate and report
1186        the progress.
1187      */
1188     int gen_first_step = 0;
1189     if (MASTER(cr))
1190     {
1191         gen_first_step = inputrec->init_step;
1192     }
1193 #endif
1194
1195     ObservablesHistory observablesHistory = {};
1196
1197     auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1198     if (startingBehavior != StartingBehavior::NewSimulation)
1199     {
1200         /* Check if checkpoint file exists before doing continuation.
1201          * This way we can use identical input options for the first and subsequent runs...
1202          */
1203         if (mdrunOptions.numStepsCommandline > -2)
1204         {
1205             /* Temporarily set the number of steps to unlimited to avoid
1206              * triggering the nsteps check in load_checkpoint().
1207              * This hack will go away soon when the -nsteps option is removed.
1208              */
1209             inputrec->nsteps = -1;
1210         }
1211
1212         // Finish applying initial simulation state information from external sources on all ranks.
1213         // Reconcile checkpoint file data with Mdrunner state established up to this point.
1214         applyLocalState(*inputHolder_.get(),
1215                         logFileHandle,
1216                         cr,
1217                         domdecOptions.numCells,
1218                         inputrec.get(),
1219                         globalState.get(),
1220                         &observablesHistory,
1221                         mdrunOptions.reproducible,
1222                         mdModules_->notifiers(),
1223                         modularSimulatorCheckpointData.get(),
1224                         useModularSimulator);
1225         // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1226         //  invariants
1227         //  on all code paths.
1228         // Write checkpoint or provide hook to update SimulationInput.
1229         // If there was a checkpoint file, SimulationInput contains more information
1230         // than if there wasn't. At this point, we have synchronized the in-memory
1231         // state with the filesystem state only for restarted simulations. We should
1232         // be calling applyLocalState unconditionally and expect that the completeness
1233         // of SimulationInput is not dependent on its creation method.
1234
1235         if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1236         {
1237             // Now we can start normal logging to the truncated log file.
1238             fplog = gmx_fio_getfp(logFileHandle);
1239             prepareLogAppending(fplog);
1240             logOwner = buildLogger(fplog, MASTER(cr));
1241             mdlog    = logOwner.logger();
1242         }
1243     }
1244
1245 #if GMX_FAHCORE
1246     if (MASTER(cr))
1247     {
1248         fcRegisterSteps(inputrec->nsteps + inputrec->init_step, gen_first_step);
1249     }
1250 #endif
1251
1252     if (mdrunOptions.numStepsCommandline > -2)
1253     {
1254         GMX_LOG(mdlog.info)
1255                 .asParagraph()
1256                 .appendText(
1257                         "The -nsteps functionality is deprecated, and may be removed in a future "
1258                         "version. "
1259                         "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1260                         "file field.");
1261     }
1262     /* override nsteps with value set on the commandline */
1263     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1264
1265     if (isSimulationMasterRank)
1266     {
1267         copy_mat(globalState->box, box);
1268     }
1269
1270     if (PAR(cr))
1271     {
1272         gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1273     }
1274
1275     if (inputrec->cutoff_scheme != CutoffScheme::Verlet)
1276     {
1277         gmx_fatal(FARGS,
1278                   "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1279                   "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1280     }
1281     /* Update rlist and nstlist. */
1282     /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1283      * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1284      * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1285      */
1286     prepare_verlet_scheme(fplog,
1287                           cr,
1288                           inputrec.get(),
1289                           nstlist_cmdline,
1290                           mtop,
1291                           box,
1292                           useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1293                           *hwinfo_->cpuInfo);
1294
1295     // We need to decide on update groups early, as this affects
1296     // inter-domain communication distances.
1297     auto       updateGroupingsPerMoleculeType = makeUpdateGroupingsPerMoleculeType(mtop);
1298     const real maxUpdateGroupRadius           = computeMaxUpdateGroupRadius(
1299             mtop, updateGroupingsPerMoleculeType, maxReferenceTemperature(*inputrec));
1300     const real   cutoffMargin = std::sqrt(max_cutoff2(inputrec->pbcType, box)) - inputrec->rlist;
1301     UpdateGroups updateGroups = makeUpdateGroups(mdlog,
1302                                                  std::move(updateGroupingsPerMoleculeType),
1303                                                  maxUpdateGroupRadius,
1304                                                  useDomainDecomposition,
1305                                                  systemHasConstraintsOrVsites(mtop),
1306                                                  cutoffMargin);
1307
1308     try
1309     {
1310         const bool haveFrozenAtoms = inputrecFrozenAtoms(inputrec.get());
1311
1312         useGpuForUpdate = decideWhetherToUseGpuForUpdate(useDomainDecomposition,
1313                                                          updateGroups.useUpdateGroups(),
1314                                                          pmeRunMode,
1315                                                          domdecOptions.numPmeRanks > 0,
1316                                                          useGpuForNonbonded,
1317                                                          updateTarget,
1318                                                          gpusWereDetected,
1319                                                          *inputrec,
1320                                                          mtop,
1321                                                          doEssentialDynamics,
1322                                                          gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1323                                                          haveFrozenAtoms,
1324                                                          doRerun,
1325                                                          devFlags,
1326                                                          mdlog);
1327     }
1328     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1329
1330     bool useGpuDirectHalo = false;
1331
1332     if (useGpuForNonbonded)
1333     {
1334         // cr->npmenodes is not yet initialized.
1335         // domdecOptions.numPmeRanks == -1 results in 0 separate PME ranks when useGpuForNonbonded is true.
1336         // Todo: remove this assumption later once auto mode has support for separate PME rank
1337         const int numPmeRanks = domdecOptions.numPmeRanks > 0 ? domdecOptions.numPmeRanks : 0;
1338         bool      havePPDomainDecomposition = (cr->sizeOfDefaultCommunicator - numPmeRanks) > 1;
1339         useGpuDirectHalo                    = decideWhetherToUseGpuForHalo(devFlags,
1340                                                         havePPDomainDecomposition,
1341                                                         useGpuForNonbonded,
1342                                                         useModularSimulator,
1343                                                         doRerun,
1344                                                         EI_ENERGY_MINIMIZATION(inputrec->eI));
1345     }
1346
1347     // This builder is necessary while we have multi-part construction
1348     // of DD. Before DD is constructed, we use the existence of
1349     // the builder object to indicate that further construction of DD
1350     // is needed.
1351     std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1352     if (useDomainDecomposition)
1353     {
1354         // P2P GPU comm + GPU update leads to case in which we enqueue async work for multiple
1355         // timesteps. DLB needs to be disabled in that case
1356         const bool directGpuCommUsedWithGpuUpdate = GMX_THREAD_MPI && useGpuDirectHalo && useGpuForUpdate;
1357         ddBuilder                                 = std::make_unique<DomainDecompositionBuilder>(
1358                 mdlog,
1359                 cr,
1360                 domdecOptions,
1361                 mdrunOptions,
1362                 mtop,
1363                 *inputrec,
1364                 mdModules_->notifiers(),
1365                 box,
1366                 updateGroups.updateGroupingPerMoleculeType(),
1367                 updateGroups.useUpdateGroups(),
1368                 updateGroups.maxUpdateGroupRadius(),
1369                 positionsFromStatePointer(globalState.get()),
1370                 useGpuForNonbonded,
1371                 useGpuForPme,
1372                 directGpuCommUsedWithGpuUpdate);
1373     }
1374     else
1375     {
1376         /* PME, if used, is done on all nodes with 1D decomposition */
1377         cr->nnodes     = cr->sizeOfDefaultCommunicator;
1378         cr->sim_nodeid = cr->rankInDefaultCommunicator;
1379         cr->nodeid     = cr->rankInDefaultCommunicator;
1380         cr->npmenodes  = 0;
1381         cr->duty       = (DUTY_PP | DUTY_PME);
1382
1383         if (inputrec->pbcType == PbcType::Screw)
1384         {
1385             gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1386         }
1387     }
1388
1389     // Produce the task assignment for this rank - done after DD is constructed
1390     GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1391             availableDevices,
1392             userGpuTaskAssignment,
1393             *hwinfo_,
1394             simulationCommunicator,
1395             physicalNodeComm,
1396             nonbondedTarget,
1397             pmeTarget,
1398             bondedTarget,
1399             updateTarget,
1400             useGpuForNonbonded,
1401             useGpuForPme,
1402             thisRankHasDuty(cr, DUTY_PP),
1403             // TODO cr->duty & DUTY_PME should imply that a PME
1404             // algorithm is active, but currently does not.
1405             EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1406
1407     // Get the device handles for the modules, nullptr when no task is assigned.
1408     int                deviceId   = -1;
1409     DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1410
1411     // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?)
1412     bool useTiming = true;
1413
1414     if (GMX_GPU_CUDA)
1415     {
1416         /* WARNING: CUDA timings are incorrect with multiple streams.
1417          *          This is the main reason why they are disabled by default.
1418          */
1419         // TODO: Consider turning on by default when we can detect nr of streams.
1420         useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
1421     }
1422     else if (GMX_GPU_OPENCL)
1423     {
1424         useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
1425     }
1426
1427     // TODO Currently this is always built, yet DD partition code
1428     // checks if it is built before using it. Probably it should
1429     // become an MDModule that is made only when another module
1430     // requires it (e.g. pull, CompEl, density fitting), so that we
1431     // don't update the local atom sets unilaterally every step.
1432     LocalAtomSetManager atomSets;
1433
1434     // Local state and topology are declared (and perhaps constructed)
1435     // now, because DD needs them for the LocalTopologyChecker, but
1436     // they do not contain valid data until after the first DD
1437     // partition.
1438     std::unique_ptr<t_state> localStateInstance;
1439     t_state*                 localState;
1440     gmx_localtop_t           localTopology(mtop.ffparams);
1441
1442     if (ddBuilder)
1443     {
1444         localStateInstance = std::make_unique<t_state>();
1445         localState         = localStateInstance.get();
1446         // TODO Pass the GPU streams to ddBuilder to use in buffer
1447         // transfers (e.g. halo exchange)
1448         cr->dd = ddBuilder->build(&atomSets, localTopology, *localState, &observablesReducerBuilder);
1449         // The builder's job is done, so destruct it
1450         ddBuilder.reset(nullptr);
1451         // Note that local state still does not exist yet.
1452     }
1453     else
1454     {
1455         // Without DD, the local state is merely an alias to the global state,
1456         // so we don't need to allocate anything.
1457         localState = globalState.get();
1458     }
1459
1460     // Ensure that all atoms within the same update group are in the
1461     // same periodic image. Otherwise, a simulation that did not use
1462     // update groups (e.g. a single-rank simulation) cannot always be
1463     // correctly restarted in a way that does use update groups
1464     // (e.g. a multi-rank simulation).
1465     if (isSimulationMasterRank)
1466     {
1467         const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1468         if (useUpdateGroups)
1469         {
1470             putUpdateGroupAtomsInSamePeriodicImage(*cr->dd, mtop, globalState->box, globalState->x);
1471         }
1472     }
1473
1474     const bool printHostName = (cr->nnodes > 1);
1475     gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1476
1477     const bool disableNonbondedCalculation = (getenv("GMX_NO_NONBONDED") != nullptr);
1478     if (disableNonbondedCalculation)
1479     {
1480         /* turn off non-bonded calculations */
1481         GMX_LOG(mdlog.warning)
1482                 .asParagraph()
1483                 .appendText(
1484                         "Found environment variable GMX_NO_NONBONDED.\n"
1485                         "Disabling nonbonded calculations.");
1486     }
1487
1488     MdrunScheduleWorkload runScheduleWork;
1489
1490     // Also populates the simulation constant workload description.
1491     // Note: currently the default duty is DUTY_PP | DUTY_PME for all simulations, including those without PME,
1492     // so this boolean is sufficient on all ranks to determine whether separate PME ranks are used,
1493     // but this will no longer be the case if cr->duty is changed for !EEL_PME(fr->ic->eeltype).
1494     const bool haveSeparatePmeRank = (!thisRankHasDuty(cr, DUTY_PP) || !thisRankHasDuty(cr, DUTY_PME));
1495     runScheduleWork.simulationWork = createSimulationWorkload(*inputrec,
1496                                                               disableNonbondedCalculation,
1497                                                               devFlags,
1498                                                               havePPDomainDecomposition(cr),
1499                                                               haveSeparatePmeRank,
1500                                                               useGpuForNonbonded,
1501                                                               pmeRunMode,
1502                                                               useGpuForBonded,
1503                                                               useGpuForUpdate,
1504                                                               useGpuDirectHalo);
1505
1506     std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1507
1508     if (deviceInfo != nullptr)
1509     {
1510         if (runScheduleWork.simulationWork.havePpDomainDecomposition && thisRankHasDuty(cr, DUTY_PP))
1511         {
1512             dd_setup_dlb_resource_sharing(cr, deviceId);
1513         }
1514         deviceStreamManager = std::make_unique<DeviceStreamManager>(
1515                 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, useTiming);
1516     }
1517
1518     // If the user chose a task assignment, give them some hints
1519     // where appropriate.
1520     if (!userGpuTaskAssignment.empty())
1521     {
1522         gpuTaskAssignments.logPerformanceHints(mdlog, numAvailableDevices);
1523     }
1524
1525     if (PAR(cr))
1526     {
1527         /* After possible communicator splitting in make_dd_communicators.
1528          * we can set up the intra/inter node communication.
1529          */
1530         gmx_setup_nodecomm(fplog, cr);
1531     }
1532
1533 #if GMX_MPI
1534     if (isMultiSim(ms))
1535     {
1536         GMX_LOG(mdlog.warning)
1537                 .asParagraph()
1538                 .appendTextFormatted(
1539                         "This is simulation %d out of %d running as a composite GROMACS\n"
1540                         "multi-simulation job. Setup for this simulation:\n",
1541                         ms->simulationIndex_,
1542                         ms->numSimulations_);
1543     }
1544     GMX_LOG(mdlog.warning)
1545             .appendTextFormatted("Using %d MPI %s\n",
1546                                  cr->nnodes,
1547 #    if GMX_THREAD_MPI
1548                                  cr->nnodes == 1 ? "thread" : "threads"
1549 #    else
1550                                  cr->nnodes == 1 ? "process" : "processes"
1551 #    endif
1552             );
1553     fflush(stderr);
1554 #endif
1555
1556     // If mdrun -pin auto honors any affinity setting that already
1557     // exists. If so, it is nice to provide feedback about whether
1558     // that existing affinity setting was from OpenMP or something
1559     // else, so we run this code both before and after we initialize
1560     // the OpenMP support.
1561     gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, FALSE);
1562     /* Check and update the number of OpenMP threads requested */
1563     checkAndUpdateRequestedNumOpenmpThreads(
1564             &hw_opt, *hwinfo_, cr, ms, physicalNodeComm.size_, pmeRunMode, mtop, *inputrec);
1565
1566     gmx_omp_nthreads_init(mdlog,
1567                           cr,
1568                           hwinfo_->nthreads_hw_avail,
1569                           physicalNodeComm.size_,
1570                           hw_opt.nthreads_omp,
1571                           hw_opt.nthreads_omp_pme,
1572                           !thisRankHasDuty(cr, DUTY_PP));
1573
1574     const bool bEnableFPE = gmxShouldEnableFPExceptions();
1575     // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1576     if (bEnableFPE)
1577     {
1578         gmx_feenableexcept();
1579     }
1580
1581     /* Now that we know the setup is consistent, check for efficiency */
1582     check_resource_division_efficiency(
1583             hwinfo_, gpuTaskAssignments.thisRankHasAnyGpuTask(), mdrunOptions.ntompOptionIsSet, cr, mdlog);
1584
1585     /* getting number of PP/PME threads on this MPI / tMPI rank.
1586        PME: env variable should be read only on one node to make sure it is
1587        identical everywhere;
1588      */
1589     const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP)
1590                                              ? gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded)
1591                                              : gmx_omp_nthreads_get(ModuleMultiThread::Pme);
1592     checkHardwareOversubscription(
1593             numThreadsOnThisRank, cr->nodeid, *hwinfo_->hardwareTopology, physicalNodeComm, mdlog);
1594
1595     // Enable Peer access between GPUs where available
1596     // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1597     // any of the GPU communication features are active.
1598     if (haveDDAtomOrdering(*cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1599         && (runScheduleWork.simulationWork.useGpuHaloExchange
1600             || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1601     {
1602         setupGpuDevicePeerAccess(gpuTaskAssignments.deviceIdsAssigned(), mdlog);
1603     }
1604
1605     if (hw_opt.threadAffinity != ThreadAffinity::Off)
1606     {
1607         /* Before setting affinity, check whether the affinity has changed
1608          * - which indicates that probably the OpenMP library has changed it
1609          * since we first checked).
1610          */
1611         gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, TRUE);
1612
1613         int numThreadsOnThisNode, intraNodeThreadOffset;
1614         analyzeThreadsOnThisNode(
1615                 physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode, &intraNodeThreadOffset);
1616
1617         /* Set the CPU affinity */
1618         gmx_set_thread_affinity(mdlog,
1619                                 cr,
1620                                 &hw_opt,
1621                                 *hwinfo_->hardwareTopology,
1622                                 numThreadsOnThisRank,
1623                                 numThreadsOnThisNode,
1624                                 intraNodeThreadOffset,
1625                                 nullptr);
1626     }
1627
1628     if (mdrunOptions.timingOptions.resetStep > -1)
1629     {
1630         GMX_LOG(mdlog.info)
1631                 .asParagraph()
1632                 .appendText(
1633                         "The -resetstep functionality is deprecated, and may be removed in a "
1634                         "future version.");
1635     }
1636     std::unique_ptr<gmx_wallcycle> wcycle =
1637             wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1638
1639     if (PAR(cr))
1640     {
1641         /* Master synchronizes its value of reset_counters with all nodes
1642          * including PME only nodes */
1643         int64_t reset_counters = wcycle_get_reset_counters(wcycle.get());
1644         gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1645         wcycle_set_reset_counters(wcycle.get(), reset_counters);
1646     }
1647
1648     // Membrane embedding must be initialized before we call init_forcerec()
1649     membedHolder.initializeMembed(fplog,
1650                                   filenames.size(),
1651                                   filenames.data(),
1652                                   &mtop,
1653                                   inputrec.get(),
1654                                   globalState.get(),
1655                                   cr,
1656                                   &mdrunOptions.checkpointOptions.period);
1657
1658     const bool               thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1659     std::unique_ptr<MDAtoms> mdAtoms;
1660     std::unique_ptr<VirtualSitesHandler> vsite;
1661
1662     t_nrnb nrnb;
1663     if (thisRankHasDuty(cr, DUTY_PP))
1664     {
1665         setupNotifier.notify(*cr);
1666         setupNotifier.notify(&atomSets);
1667         setupNotifier.notify(mtop);
1668         setupNotifier.notify(inputrec->pbcType);
1669         setupNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1670         /* Initiate forcerecord */
1671         fr                 = std::make_unique<t_forcerec>();
1672         fr->forceProviders = mdModules_->initForceProviders();
1673         init_forcerec(fplog,
1674                       mdlog,
1675                       runScheduleWork.simulationWork,
1676                       fr.get(),
1677                       *inputrec,
1678                       mtop,
1679                       cr,
1680                       box,
1681                       opt2fn("-table", filenames.size(), filenames.data()),
1682                       opt2fn("-tablep", filenames.size(), filenames.data()),
1683                       opt2fns("-tableb", filenames.size(), filenames.data()),
1684                       pforce);
1685         // Dirty hack, for fixing disres and orires should be made mdmodules
1686         fr->fcdata->disres = disresdata;
1687         if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1688         {
1689             fr->fcdata->orires = std::make_unique<t_oriresdata>(
1690                     fplog, mtop, *inputrec, ms, globalState.get(), &atomSets);
1691         }
1692
1693         // Save a handle to device stream manager to use elsewhere in the code
1694         // TODO: Forcerec is not a correct place to store it.
1695         fr->deviceStreamManager = deviceStreamManager.get();
1696
1697         if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1698         {
1699             GMX_RELEASE_ASSERT(
1700                     deviceStreamManager != nullptr,
1701                     "GPU device stream manager should be valid in order to use PME-PP direct "
1702                     "communications.");
1703             GMX_RELEASE_ASSERT(
1704                     deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1705                     "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1706                     "communications.");
1707             fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1708                     cr->mpi_comm_mysim,
1709                     cr->dd->pme_nodeid,
1710                     &cr->dd->pmeForceReceiveBuffer,
1711                     deviceStreamManager->context(),
1712                     deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1713         }
1714
1715         fr->nbv = Nbnxm::init_nb_verlet(mdlog,
1716                                         *inputrec,
1717                                         *fr,
1718                                         cr,
1719                                         *hwinfo_,
1720                                         runScheduleWork.simulationWork.useGpuNonbonded,
1721                                         deviceStreamManager.get(),
1722                                         mtop,
1723                                         box,
1724                                         wcycle.get());
1725         // TODO: Move the logic below to a GPU bonded builder
1726         if (runScheduleWork.simulationWork.useGpuBonded)
1727         {
1728             GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1729                                "GPU device stream manager should be valid in order to use GPU "
1730                                "version of bonded forces.");
1731             fr->listedForcesGpu = std::make_unique<ListedForcesGpu>(
1732                     mtop.ffparams,
1733                     fr->ic->epsfac * fr->fudgeQQ,
1734                     deviceStreamManager->context(),
1735                     deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)),
1736                     wcycle.get());
1737         }
1738         fr->longRangeNonbondeds = std::make_unique<CpuPpLongRangeNonbondeds>(fr->n_tpi,
1739                                                                              fr->ic->ewaldcoeff_q,
1740                                                                              fr->ic->epsilon_r,
1741                                                                              fr->qsum,
1742                                                                              fr->ic->eeltype,
1743                                                                              fr->ic->vdwtype,
1744                                                                              *inputrec,
1745                                                                              &nrnb,
1746                                                                              wcycle.get(),
1747                                                                              fplog);
1748
1749         /* Initialize the mdAtoms structure.
1750          * mdAtoms is not filled with atom data,
1751          * as this can not be done now with domain decomposition.
1752          */
1753         mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1754         if (globalState && thisRankHasPmeGpuTask)
1755         {
1756             // The pinning of coordinates in the global state object works, because we only use
1757             // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1758             // points to the global state object without DD.
1759             // FIXME: MD and EM separately set up the local state - this should happen in the same
1760             // function, which should also perform the pinning.
1761             changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1762         }
1763
1764         /* Initialize the virtual site communication */
1765         vsite = makeVirtualSitesHandler(
1766                 mtop, cr, fr->pbcType, updateGroups.updateGroupingPerMoleculeType());
1767
1768         calc_shifts(box, fr->shift_vec);
1769
1770         /* With periodic molecules the charge groups should be whole at start up
1771          * and the virtual sites should not be far from their proper positions.
1772          */
1773         if (!inputrec->bContinuation && MASTER(cr)
1774             && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1775         {
1776             /* Make molecules whole at start of run */
1777             if (fr->pbcType != PbcType::No)
1778             {
1779                 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1780             }
1781             if (vsite)
1782             {
1783                 /* Correct initial vsite positions are required
1784                  * for the initial distribution in the domain decomposition
1785                  * and for the initial shell prediction.
1786                  */
1787                 constructVirtualSitesGlobal(mtop, globalState->x);
1788             }
1789         }
1790         // Make the DD reverse topology, now that any vsites that are present are available
1791         if (haveDDAtomOrdering(*cr))
1792         {
1793             dd_make_reverse_top(fplog, cr->dd, mtop, vsite.get(), *inputrec, domdecOptions.ddBondedChecking);
1794         }
1795
1796         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1797         {
1798             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1799             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1800         }
1801     }
1802     else
1803     {
1804         /* This is a PME only node */
1805
1806         GMX_ASSERT(globalState == nullptr,
1807                    "We don't need the state on a PME only rank and expect it to be unitialized");
1808
1809         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1810         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1811     }
1812
1813     gmx_pme_t* sepPmeData = nullptr;
1814     // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1815     GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1816                "Double-checking that only PME-only ranks have no forcerec");
1817     gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1818
1819     // TODO should live in ewald module once its testing is improved
1820     //
1821     // Later, this program could contain kernels that might be later
1822     // re-used as auto-tuning progresses, or subsequent simulations
1823     // are invoked.
1824     PmeGpuProgramStorage pmeGpuProgram;
1825     if (thisRankHasPmeGpuTask)
1826     {
1827         GMX_RELEASE_ASSERT(
1828                 (deviceStreamManager != nullptr),
1829                 "GPU device stream manager should be initialized in order to use GPU for PME.");
1830         GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1831                            "GPU device should be initialized in order to use GPU for PME.");
1832         pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1833     }
1834
1835     /* Initiate PME if necessary,
1836      * either on all nodes or on dedicated PME nodes only. */
1837     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1838     {
1839         if (mdAtoms && mdAtoms->mdatoms())
1840         {
1841             nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1842             if (EVDW_PME(inputrec->vdwtype))
1843             {
1844                 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1845             }
1846         }
1847         if (cr->npmenodes > 0)
1848         {
1849             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1850             gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1851             gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1852         }
1853
1854         if (thisRankHasDuty(cr, DUTY_PME))
1855         {
1856             try
1857             {
1858                 // TODO: This should be in the builder.
1859                 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1860                                            || (deviceStreamManager != nullptr),
1861                                    "Device stream manager should be valid in order to use GPU "
1862                                    "version of PME.");
1863                 GMX_RELEASE_ASSERT(
1864                         !runScheduleWork.simulationWork.useGpuPme
1865                                 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1866                         "GPU PME stream should be valid in order to use GPU version of PME.");
1867
1868                 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1869                                                              ? &deviceStreamManager->context()
1870                                                              : nullptr;
1871                 const DeviceStream*  pmeStream =
1872                         runScheduleWork.simulationWork.useGpuPme
1873                                  ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1874                                  : nullptr;
1875
1876                 pmedata = gmx_pme_init(cr,
1877                                        getNumPmeDomains(cr->dd),
1878                                        inputrec.get(),
1879                                        nChargePerturbed != 0,
1880                                        nTypePerturbed != 0,
1881                                        mdrunOptions.reproducible,
1882                                        ewaldcoeff_q,
1883                                        ewaldcoeff_lj,
1884                                        gmx_omp_nthreads_get(ModuleMultiThread::Pme),
1885                                        pmeRunMode,
1886                                        nullptr,
1887                                        deviceContext,
1888                                        pmeStream,
1889                                        pmeGpuProgram.get(),
1890                                        mdlog);
1891             }
1892             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1893         }
1894     }
1895
1896
1897     if (EI_DYNAMICS(inputrec->eI))
1898     {
1899         /* Turn on signal handling on all nodes */
1900         /*
1901          * (A user signal from the PME nodes (if any)
1902          * is communicated to the PP nodes.
1903          */
1904         signal_handler_install();
1905     }
1906
1907     pull_t* pull_work = nullptr;
1908     if (thisRankHasDuty(cr, DUTY_PP))
1909     {
1910         /* Assumes uniform use of the number of OpenMP threads */
1911         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(ModuleMultiThread::Default));
1912
1913         if (inputrec->bPull)
1914         {
1915             /* Initialize pull code */
1916             pull_work = init_pull(fplog,
1917                                   inputrec->pull.get(),
1918                                   inputrec.get(),
1919                                   mtop,
1920                                   cr,
1921                                   &atomSets,
1922                                   inputrec->fepvals->init_lambda);
1923             if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1924             {
1925                 initPullHistory(pull_work, &observablesHistory);
1926             }
1927             if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1928             {
1929                 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1930             }
1931         }
1932
1933         std::unique_ptr<EnforcedRotation> enforcedRotation;
1934         if (inputrec->bRot)
1935         {
1936             /* Initialize enforced rotation code */
1937             enforcedRotation = init_rot(fplog,
1938                                         inputrec.get(),
1939                                         filenames.size(),
1940                                         filenames.data(),
1941                                         cr,
1942                                         &atomSets,
1943                                         globalState.get(),
1944                                         mtop,
1945                                         oenv,
1946                                         mdrunOptions,
1947                                         startingBehavior);
1948         }
1949
1950         t_swap* swap = nullptr;
1951         if (inputrec->eSwapCoords != SwapType::No)
1952         {
1953             /* Initialize ion swapping code */
1954             swap = init_swapcoords(fplog,
1955                                    inputrec.get(),
1956                                    opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1957                                    mtop,
1958                                    globalState.get(),
1959                                    &observablesHistory,
1960                                    cr,
1961                                    &atomSets,
1962                                    oenv,
1963                                    mdrunOptions,
1964                                    startingBehavior);
1965         }
1966
1967         /* Let makeConstraints know whether we have essential dynamics constraints. */
1968         auto constr = makeConstraints(mtop,
1969                                       *inputrec,
1970                                       pull_work,
1971                                       doEssentialDynamics,
1972                                       fplog,
1973                                       cr,
1974                                       updateGroups.useUpdateGroups(),
1975                                       ms,
1976                                       &nrnb,
1977                                       wcycle.get(),
1978                                       fr->bMolPBC);
1979
1980         /* Energy terms and groups */
1981         gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1982                              inputrec->fepvals->n_lambda);
1983
1984         // cos acceleration is only supported by md, but older tpr
1985         // files might still combine it with other integrators
1986         GMX_RELEASE_ASSERT(inputrec->cos_accel == 0.0 || inputrec->eI == IntegrationAlgorithm::MD,
1987                            "cos_acceleration is only supported by integrator=md");
1988
1989         /* Kinetic energy data */
1990         gmx_ekindata_t ekind(inputrec->opts.ngtc,
1991                              inputrec->cos_accel,
1992                              gmx_omp_nthreads_get(ModuleMultiThread::Update));
1993
1994         /* Set up interactive MD (IMD) */
1995         auto imdSession = makeImdSession(inputrec.get(),
1996                                          cr,
1997                                          wcycle.get(),
1998                                          &enerd,
1999                                          ms,
2000                                          mtop,
2001                                          mdlog,
2002                                          MASTER(cr) ? globalState->x : gmx::ArrayRef<gmx::RVec>(),
2003                                          filenames.size(),
2004                                          filenames.data(),
2005                                          oenv,
2006                                          mdrunOptions.imdOptions,
2007                                          startingBehavior);
2008
2009         if (haveDDAtomOrdering(*cr))
2010         {
2011             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
2012             /* This call is not included in init_domain_decomposition
2013              * because fr->atomInfoForEachMoleculeBlock is set later.
2014              */
2015             makeBondedLinks(cr->dd, mtop, fr->atomInfoForEachMoleculeBlock);
2016         }
2017
2018         if (runScheduleWork.simulationWork.useGpuBufferOps)
2019         {
2020             fr->gpuForceReduction[gmx::AtomLocality::Local] = std::make_unique<gmx::GpuForceReduction>(
2021                     deviceStreamManager->context(),
2022                     deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedLocal),
2023                     wcycle.get());
2024             fr->gpuForceReduction[gmx::AtomLocality::NonLocal] = std::make_unique<gmx::GpuForceReduction>(
2025                     deviceStreamManager->context(),
2026                     deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedNonLocal),
2027                     wcycle.get());
2028         }
2029
2030         std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
2031         if (gpusWereDetected
2032             && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
2033                 || runScheduleWork.simulationWork.useGpuBufferOps))
2034         {
2035             GpuApiCallBehavior transferKind =
2036                     (inputrec->eI == IntegrationAlgorithm::MD && !doRerun && !useModularSimulator)
2037                             ? GpuApiCallBehavior::Async
2038                             : GpuApiCallBehavior::Sync;
2039             GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
2040                                "GPU device stream manager should be initialized to use GPU.");
2041             stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
2042                     *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle.get());
2043             fr->stateGpu = stateGpu.get();
2044         }
2045
2046         GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
2047         SimulatorBuilder simulatorBuilder;
2048
2049         simulatorBuilder.add(SimulatorStateData(
2050                 globalState.get(), localState, &observablesHistory, &enerd, &ekind));
2051         simulatorBuilder.add(std::move(membedHolder));
2052         simulatorBuilder.add(std::move(stopHandlerBuilder_));
2053         simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
2054
2055
2056         simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv, &observablesReducerBuilder));
2057         simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle.get()));
2058         simulatorBuilder.add(ConstraintsParam(
2059                 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, vsite.get()));
2060         // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
2061         simulatorBuilder.add(LegacyInput(
2062                 static_cast<int>(filenames.size()), filenames.data(), inputrec.get(), fr.get()));
2063         simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
2064         simulatorBuilder.add(InteractiveMD(imdSession.get()));
2065         simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifiers()));
2066         simulatorBuilder.add(CenterOfMassPulling(pull_work));
2067         // Todo move to an MDModule
2068         simulatorBuilder.add(IonSwapping(swap));
2069         simulatorBuilder.add(TopologyData(mtop, &localTopology, mdAtoms.get()));
2070         simulatorBuilder.add(BoxDeformationHandle(deform.get()));
2071         simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
2072
2073         // build and run simulator object based on user-input
2074         auto simulator = simulatorBuilder.build(useModularSimulator);
2075         simulator->run();
2076
2077         if (fr->pmePpCommGpu)
2078         {
2079             // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
2080             fr->pmePpCommGpu.reset();
2081         }
2082
2083         if (inputrec->bPull)
2084         {
2085             finish_pull(pull_work);
2086         }
2087         finish_swapcoords(swap);
2088     }
2089     else
2090     {
2091         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
2092         /* do PME only */
2093         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(ModuleMultiThread::Pme));
2094         gmx_pmeonly(pmedata,
2095                     cr,
2096                     &nrnb,
2097                     wcycle.get(),
2098                     walltime_accounting,
2099                     inputrec.get(),
2100                     pmeRunMode,
2101                     runScheduleWork.simulationWork.useGpuPmePpCommunication,
2102                     deviceStreamManager.get());
2103     }
2104
2105     wallcycle_stop(wcycle.get(), WallCycleCounter::Run);
2106
2107     /* Finish up, write some stuff
2108      * if rerunMD, don't write last frame again
2109      */
2110     finish_run(fplog,
2111                mdlog,
2112                cr,
2113                *inputrec,
2114                &nrnb,
2115                wcycle.get(),
2116                walltime_accounting,
2117                fr ? fr->nbv.get() : nullptr,
2118                pmedata,
2119                EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
2120
2121
2122     deviceStreamManager.reset(nullptr);
2123     // Free PME data
2124     if (pmedata)
2125     {
2126         gmx_pme_destroy(pmedata);
2127         pmedata = nullptr;
2128     }
2129
2130     // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
2131     // before we destroy the GPU context(s)
2132     // Pinned buffers are associated with contexts in CUDA.
2133     // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
2134     mdAtoms.reset(nullptr);
2135     globalState.reset(nullptr);
2136     localStateInstance.reset(nullptr);
2137     mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
2138     fr.reset(nullptr);         // destruct forcerec before gpu
2139     // TODO convert to C++ so we can get rid of these frees
2140     sfree(disresdata);
2141
2142     if (!hwinfo_->deviceInfoList.empty())
2143     {
2144         /* stop the GPU profiler (only CUDA) */
2145         stopGpuProfiler();
2146     }
2147
2148     /* With tMPI we need to wait for all ranks to finish deallocation before
2149      * destroying the CUDA context as some tMPI ranks may be sharing
2150      * GPU and context.
2151      *
2152      * This is not a concern in OpenCL where we use one context per rank.
2153      *
2154      * Note: it is safe to not call the barrier on the ranks which do not use GPU,
2155      * but it is easier and more futureproof to call it on the whole node.
2156      *
2157      * Note that this function needs to be called even if GPUs are not used
2158      * in this run because the PME ranks have no knowledge of whether GPUs
2159      * are used or not, but all ranks need to enter the barrier below.
2160      * \todo Remove this physical node barrier after making sure
2161      * that it's not needed anymore (with a shared GPU run).
2162      */
2163     if (GMX_THREAD_MPI)
2164     {
2165         physicalNodeComm.barrier();
2166     }
2167
2168     if (!devFlags.usingCudaAwareMpi)
2169     {
2170         // Don't reset GPU in case of CUDA-AWARE MPI
2171         // UCX creates CUDA buffers which are cleaned-up as part of MPI_Finalize()
2172         // resetting the device before MPI_Finalize() results in crashes inside UCX
2173         releaseDevice(deviceInfo);
2174     }
2175
2176     /* Does what it says */
2177     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
2178     walltime_accounting_destroy(walltime_accounting);
2179
2180     // Ensure log file content is written
2181     if (logFileHandle)
2182     {
2183         gmx_fio_flush(logFileHandle);
2184     }
2185
2186     /* Reset FPEs (important for unit tests) by disabling them. Assumes no
2187      * exceptions were enabled before function was called. */
2188     if (bEnableFPE)
2189     {
2190         gmx_fedisableexcept();
2191     }
2192
2193     auto rc = static_cast<int>(gmx_get_stop_condition());
2194
2195 #if GMX_THREAD_MPI
2196     /* we need to join all threads. The sub-threads join when they
2197        exit this function, but the master thread needs to be told to
2198        wait for that. */
2199     if (MASTER(cr))
2200     {
2201         tMPI_Finalize();
2202     }
2203 #endif
2204     return rc;
2205 } // namespace gmx
2206
2207 Mdrunner::~Mdrunner()
2208 {
2209     // Clean up of the Manager.
2210     // This will end up getting called on every thread-MPI rank, which is unnecessary,
2211     // but okay as long as threads synchronize some time before adding or accessing
2212     // a new set of restraints.
2213     if (restraintManager_)
2214     {
2215         restraintManager_->clear();
2216         GMX_ASSERT(restraintManager_->countRestraints() == 0,
2217                    "restraints added during runner life time should be cleared at runner "
2218                    "destruction.");
2219     }
2220 };
2221
2222 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
2223 {
2224     GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
2225     // Not sure if this should be logged through the md logger or something else,
2226     // but it is helpful to have some sort of INFO level message sent somewhere.
2227     //    std::cout << "Registering restraint named " << name << std::endl;
2228
2229     // When multiple restraints are used, it may be wasteful to register them separately.
2230     // Maybe instead register an entire Restraint Manager as a force provider.
2231     restraintManager_->addToSpec(std::move(puller), name);
2232 }
2233
2234 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
2235
2236 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
2237
2238 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265 in CentOS 7
2239 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
2240
2241 class Mdrunner::BuilderImplementation
2242 {
2243 public:
2244     BuilderImplementation() = delete;
2245     BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
2246     ~BuilderImplementation();
2247
2248     BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
2249                                                 real                forceWarningThreshold,
2250                                                 StartingBehavior    startingBehavior);
2251
2252     void addHardwareDetectionResult(const gmx_hw_info_t* hwinfo);
2253
2254     void addDomdec(const DomdecOptions& options);
2255
2256     void addInput(SimulationInputHandle inputHolder);
2257
2258     void addVerletList(int nstlist);
2259
2260     void addReplicaExchange(const ReplicaExchangeParameters& params);
2261
2262     void addNonBonded(const char* nbpu_opt);
2263
2264     void addPME(const char* pme_opt_, const char* pme_fft_opt_);
2265
2266     void addBondedTaskAssignment(const char* bonded_opt);
2267
2268     void addUpdateTaskAssignment(const char* update_opt);
2269
2270     void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
2271
2272     void addFilenames(ArrayRef<const t_filenm> filenames);
2273
2274     void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
2275
2276     void addLogFile(t_fileio* logFileHandle);
2277
2278     void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
2279
2280     Mdrunner build();
2281
2282 private:
2283     // Default parameters copied from runner.h
2284     // \todo Clarify source(s) of default parameters.
2285
2286     const char* nbpu_opt_    = nullptr;
2287     const char* pme_opt_     = nullptr;
2288     const char* pme_fft_opt_ = nullptr;
2289     const char* bonded_opt_  = nullptr;
2290     const char* update_opt_  = nullptr;
2291
2292     MdrunOptions mdrunOptions_;
2293
2294     DomdecOptions domdecOptions_;
2295
2296     ReplicaExchangeParameters replicaExchangeParameters_;
2297
2298     //! Command-line override for the duration of a neighbor list with the Verlet scheme.
2299     int nstlist_ = 0;
2300
2301     //! World communicator, used for hardware detection and task assignment
2302     MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
2303
2304     //! Multisim communicator handle.
2305     gmx_multisim_t* multiSimulation_;
2306
2307     //! mdrun communicator
2308     MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
2309
2310     //! Print a warning if any force is larger than this (in kJ/mol nm).
2311     real forceWarningThreshold_ = -1;
2312
2313     //! Whether the simulation will start afresh, or restart with/without appending.
2314     StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
2315
2316     //! The modules that comprise the functionality of mdrun.
2317     std::unique_ptr<MDModules> mdModules_;
2318
2319     //! Detected hardware.
2320     const gmx_hw_info_t* hwinfo_ = nullptr;
2321
2322     //! \brief Parallelism information.
2323     gmx_hw_opt_t hardwareOptions_;
2324
2325     //! filename options for simulation.
2326     ArrayRef<const t_filenm> filenames_;
2327
2328     /*! \brief Handle to output environment.
2329      *
2330      * \todo gmx_output_env_t needs lifetime management.
2331      */
2332     gmx_output_env_t* outputEnvironment_ = nullptr;
2333
2334     /*! \brief Non-owning handle to MD log file.
2335      *
2336      * \todo Context should own output facilities for client.
2337      * \todo Improve log file handle management.
2338      * \internal
2339      * Code managing the FILE* relies on the ability to set it to
2340      * nullptr to check whether the filehandle is valid.
2341      */
2342     t_fileio* logFileHandle_ = nullptr;
2343
2344     /*!
2345      * \brief Builder for simulation stop signal handler.
2346      */
2347     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
2348
2349     /*!
2350      * \brief Sources for initial simulation state.
2351      *
2352      * See issue #3652 for near-term refinements to the SimulationInput interface.
2353      *
2354      * See issue #3379 for broader discussion on API aspects of simulation inputs and outputs.
2355      */
2356     SimulationInputHandle inputHolder_;
2357 };
2358
2359 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
2360                                                        compat::not_null<SimulationContext*> context) :
2361     mdModules_(std::move(mdModules))
2362 {
2363     libraryWorldCommunicator_ = context->libraryWorldCommunicator_;
2364     simulationCommunicator_   = context->simulationCommunicator_;
2365     multiSimulation_          = context->multiSimulation_.get();
2366 }
2367
2368 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
2369
2370 Mdrunner::BuilderImplementation&
2371 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions&    options,
2372                                                       const real             forceWarningThreshold,
2373                                                       const StartingBehavior startingBehavior)
2374 {
2375     mdrunOptions_          = options;
2376     forceWarningThreshold_ = forceWarningThreshold;
2377     startingBehavior_      = startingBehavior;
2378     return *this;
2379 }
2380
2381 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
2382 {
2383     domdecOptions_ = options;
2384 }
2385
2386 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
2387 {
2388     nstlist_ = nstlist;
2389 }
2390
2391 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
2392 {
2393     replicaExchangeParameters_ = params;
2394 }
2395
2396 Mdrunner Mdrunner::BuilderImplementation::build()
2397 {
2398     auto newRunner = Mdrunner(std::move(mdModules_));
2399
2400     newRunner.mdrunOptions     = mdrunOptions_;
2401     newRunner.pforce           = forceWarningThreshold_;
2402     newRunner.startingBehavior = startingBehavior_;
2403     newRunner.domdecOptions    = domdecOptions_;
2404
2405     // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2406     newRunner.hw_opt = hardwareOptions_;
2407
2408     // No invariant to check. This parameter exists to optionally override other behavior.
2409     newRunner.nstlist_cmdline = nstlist_;
2410
2411     newRunner.replExParams = replicaExchangeParameters_;
2412
2413     newRunner.filenames = filenames_;
2414
2415     newRunner.libraryWorldCommunicator = libraryWorldCommunicator_;
2416
2417     newRunner.simulationCommunicator = simulationCommunicator_;
2418
2419     // nullptr is a valid value for the multisim handle
2420     newRunner.ms = multiSimulation_;
2421
2422     if (hwinfo_)
2423     {
2424         newRunner.hwinfo_ = hwinfo_;
2425     }
2426     else
2427     {
2428         GMX_THROW(gmx::APIError(
2429                 "MdrunnerBuilder::addHardwareDetectionResult() is required before build()"));
2430     }
2431
2432     if (inputHolder_)
2433     {
2434         newRunner.inputHolder_ = std::move(inputHolder_);
2435     }
2436     else
2437     {
2438         GMX_THROW(gmx::APIError("MdrunnerBuilder::addInput() is required before build()."));
2439     }
2440
2441     // \todo Clarify ownership and lifetime management for gmx_output_env_t
2442     // \todo Update sanity checking when output environment has clearly specified invariants.
2443     // Initialization and default values for oenv are not well specified in the current version.
2444     if (outputEnvironment_)
2445     {
2446         newRunner.oenv = outputEnvironment_;
2447     }
2448     else
2449     {
2450         GMX_THROW(gmx::APIError(
2451                 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2452     }
2453
2454     newRunner.logFileHandle = logFileHandle_;
2455
2456     if (nbpu_opt_)
2457     {
2458         newRunner.nbpu_opt = nbpu_opt_;
2459     }
2460     else
2461     {
2462         GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2463     }
2464
2465     if (pme_opt_ && pme_fft_opt_)
2466     {
2467         newRunner.pme_opt     = pme_opt_;
2468         newRunner.pme_fft_opt = pme_fft_opt_;
2469     }
2470     else
2471     {
2472         GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2473     }
2474
2475     if (bonded_opt_)
2476     {
2477         newRunner.bonded_opt = bonded_opt_;
2478     }
2479     else
2480     {
2481         GMX_THROW(gmx::APIError(
2482                 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2483     }
2484
2485     if (update_opt_)
2486     {
2487         newRunner.update_opt = update_opt_;
2488     }
2489     else
2490     {
2491         GMX_THROW(gmx::APIError(
2492                 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build()  "));
2493     }
2494
2495
2496     newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2497
2498     if (stopHandlerBuilder_)
2499     {
2500         newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2501     }
2502     else
2503     {
2504         newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2505     }
2506
2507     return newRunner;
2508 }
2509
2510 void Mdrunner::BuilderImplementation::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2511 {
2512     hwinfo_ = hwinfo;
2513 }
2514
2515 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2516 {
2517     nbpu_opt_ = nbpu_opt;
2518 }
2519
2520 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2521 {
2522     pme_opt_     = pme_opt;
2523     pme_fft_opt_ = pme_fft_opt;
2524 }
2525
2526 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2527 {
2528     bonded_opt_ = bonded_opt;
2529 }
2530
2531 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2532 {
2533     update_opt_ = update_opt;
2534 }
2535
2536 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2537 {
2538     hardwareOptions_ = hardwareOptions;
2539 }
2540
2541 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2542 {
2543     filenames_ = filenames;
2544 }
2545
2546 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2547 {
2548     outputEnvironment_ = outputEnvironment;
2549 }
2550
2551 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2552 {
2553     logFileHandle_ = logFileHandle;
2554 }
2555
2556 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2557 {
2558     stopHandlerBuilder_ = std::move(builder);
2559 }
2560
2561 void Mdrunner::BuilderImplementation::addInput(SimulationInputHandle inputHolder)
2562 {
2563     inputHolder_ = std::move(inputHolder);
2564 }
2565
2566 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules>           mdModules,
2567                                  compat::not_null<SimulationContext*> context) :
2568     impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2569 {
2570 }
2571
2572 MdrunnerBuilder::~MdrunnerBuilder() = default;
2573
2574 MdrunnerBuilder& MdrunnerBuilder::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2575 {
2576     impl_->addHardwareDetectionResult(hwinfo);
2577     return *this;
2578 }
2579
2580 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions&    options,
2581                                                       real                   forceWarningThreshold,
2582                                                       const StartingBehavior startingBehavior)
2583 {
2584     impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2585     return *this;
2586 }
2587
2588 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2589 {
2590     impl_->addDomdec(options);
2591     return *this;
2592 }
2593
2594 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2595 {
2596     impl_->addVerletList(nstlist);
2597     return *this;
2598 }
2599
2600 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2601 {
2602     impl_->addReplicaExchange(params);
2603     return *this;
2604 }
2605
2606 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2607 {
2608     impl_->addNonBonded(nbpu_opt);
2609     return *this;
2610 }
2611
2612 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2613 {
2614     // The builder method may become more general in the future, but in this version,
2615     // parameters for PME electrostatics are both required and the only parameters
2616     // available.
2617     if (pme_opt && pme_fft_opt)
2618     {
2619         impl_->addPME(pme_opt, pme_fft_opt);
2620     }
2621     else
2622     {
2623         GMX_THROW(
2624                 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2625     }
2626     return *this;
2627 }
2628
2629 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2630 {
2631     impl_->addBondedTaskAssignment(bonded_opt);
2632     return *this;
2633 }
2634
2635 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2636 {
2637     impl_->addUpdateTaskAssignment(update_opt);
2638     return *this;
2639 }
2640
2641 Mdrunner MdrunnerBuilder::build()
2642 {
2643     return impl_->build();
2644 }
2645
2646 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2647 {
2648     impl_->addHardwareOptions(hardwareOptions);
2649     return *this;
2650 }
2651
2652 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2653 {
2654     impl_->addFilenames(filenames);
2655     return *this;
2656 }
2657
2658 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2659 {
2660     impl_->addOutputEnvironment(outputEnvironment);
2661     return *this;
2662 }
2663
2664 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2665 {
2666     impl_->addLogFile(logFileHandle);
2667     return *this;
2668 }
2669
2670 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2671 {
2672     impl_->addStopHandlerBuilder(std::move(builder));
2673     return *this;
2674 }
2675
2676 MdrunnerBuilder& MdrunnerBuilder::addInput(SimulationInputHandle input)
2677 {
2678     impl_->addInput(std::move(input));
2679     return *this;
2680 }
2681
2682 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2683
2684 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;
2685
2686 } // namespace gmx