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