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