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