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