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