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