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