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