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