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