88fa3ea8dacccdb905ba1571324d6dceeff2afa0
[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     // NOTE: The devFlags need decideWhetherToUseGpusForNonbonded(...) and decideWhetherToUseGpusForPme(...) for overrides,
889     //       decideWhetherToUseGpuForUpdate() needs devFlags for the '-update auto' override, hence the interleaving.
890     // NOTE: When the simulationWork is constructed, the useGpuForUpdate overrides the devFlags.enableGpuBufferOps.
891     try
892     {
893         useGpuForUpdate = decideWhetherToUseGpuForUpdate(devFlags.forceGpuUpdateDefaultOn,
894                                                          useDomainDecomposition,
895                                                          useGpuForPme,
896                                                          useGpuForNonbonded,
897                                                          updateTarget,
898                                                          gpusWereDetected,
899                                                          *inputrec,
900                                                          gmx_mtop_interaction_count(mtop, IF_VSITE) > 0,
901                                                          doEssentialDynamics,
902                                                          gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
903                                                          replExParams.exchangeInterval > 0);
904     }
905     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
906
907
908     // Build restraints.
909     // TODO: hide restraint implementation details from Mdrunner.
910     // There is nothing unique about restraints at this point as far as the
911     // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
912     // factory functions from the SimulationContext on which to call mdModules_->add().
913     // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
914     for (auto && restraint : restraintManager_->getRestraints())
915     {
916         auto module = RestraintMDModule::create(restraint,
917                                                 restraint->sites());
918         mdModules_->add(std::move(module));
919     }
920
921     // TODO: Error handling
922     mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
923     const auto &mdModulesNotifier = mdModules_->notifier().notifier_;
924
925     if (inputrec->internalParameters != nullptr)
926     {
927         mdModulesNotifier.notify(*inputrec->internalParameters);
928     }
929
930     if (fplog != nullptr)
931     {
932         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
933         fprintf(fplog, "\n");
934     }
935
936     if (SIMMASTER(cr))
937     {
938         /* In rerun, set velocities to zero if present */
939         if (doRerun && ((globalState->flags & (1 << estV)) != 0))
940         {
941             // rerun does not use velocities
942             GMX_LOG(mdlog.info).asParagraph().appendText(
943                     "Rerun trajectory contains velocities. Rerun does only evaluate "
944                     "potential energy and forces. The velocities will be ignored.");
945             for (int i = 0; i < globalState->natoms; i++)
946             {
947                 clear_rvec(globalState->v[i]);
948             }
949             globalState->flags &= ~(1 << estV);
950         }
951
952         /* now make sure the state is initialized and propagated */
953         set_state_entries(globalState.get(), inputrec);
954     }
955
956     /* NM and TPI parallelize over force/energy calculations, not atoms,
957      * so we need to initialize and broadcast the global state.
958      */
959     if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
960     {
961         if (!MASTER(cr))
962         {
963             globalState = std::make_unique<t_state>();
964         }
965         broadcastStateWithoutDynamics(cr, globalState.get());
966     }
967
968     /* A parallel command line option consistency check that we can
969        only do after any threads have started. */
970     if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
971                      domdecOptions.numCells[YY] > 1 ||
972                      domdecOptions.numCells[ZZ] > 1 ||
973                      domdecOptions.numPmeRanks > 0))
974     {
975         gmx_fatal(FARGS,
976                   "The -dd or -npme option request a parallel simulation, "
977 #if !GMX_MPI
978                   "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
979 #elif GMX_THREAD_MPI
980                   "but the number of MPI-threads (option -ntmpi) is not set or is 1");
981 #else
982                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
983 #endif
984     }
985
986     if (doRerun &&
987         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
988     {
989         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
990     }
991
992     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
993     {
994         if (domdecOptions.numPmeRanks > 0)
995         {
996             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
997                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
998         }
999
1000         domdecOptions.numPmeRanks = 0;
1001     }
1002
1003     if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
1004     {
1005         /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
1006          * improve performance with many threads per GPU, since our OpenMP
1007          * scaling is bad, but it's difficult to automate the setup.
1008          */
1009         domdecOptions.numPmeRanks = 0;
1010     }
1011     if (useGpuForPme)
1012     {
1013         if (domdecOptions.numPmeRanks < 0)
1014         {
1015             domdecOptions.numPmeRanks = 0;
1016             // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
1017         }
1018         else
1019         {
1020             GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
1021         }
1022     }
1023
1024 #if GMX_FAHCORE
1025     if (MASTER(cr))
1026     {
1027         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1028     }
1029 #endif
1030
1031     /* NMR restraints must be initialized before load_checkpoint,
1032      * since with time averaging the history is added to t_state.
1033      * For proper consistency check we therefore need to extend
1034      * t_state here.
1035      * So the PME-only nodes (if present) will also initialize
1036      * the distance restraints.
1037      */
1038     snew(fcd, 1);
1039
1040     /* This needs to be called before read_checkpoint to extend the state */
1041     init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
1042
1043     init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
1044
1045     auto                 deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
1046
1047     ObservablesHistory   observablesHistory = {};
1048
1049     if (startingBehavior != StartingBehavior::NewSimulation)
1050     {
1051         /* Check if checkpoint file exists before doing continuation.
1052          * This way we can use identical input options for the first and subsequent runs...
1053          */
1054         if (mdrunOptions.numStepsCommandline > -2)
1055         {
1056             /* Temporarily set the number of steps to unlimited to avoid
1057              * triggering the nsteps check in load_checkpoint().
1058              * This hack will go away soon when the -nsteps option is removed.
1059              */
1060             inputrec->nsteps = -1;
1061         }
1062
1063         load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
1064                         logFileHandle,
1065                         cr, domdecOptions.numCells,
1066                         inputrec, globalState.get(),
1067                         &observablesHistory,
1068                         mdrunOptions.reproducible, mdModules_->notifier());
1069
1070         if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1071         {
1072             // Now we can start normal logging to the truncated log file.
1073             fplog    = gmx_fio_getfp(logFileHandle);
1074             prepareLogAppending(fplog);
1075             logOwner = buildLogger(fplog, MASTER(cr));
1076             mdlog    = logOwner.logger();
1077         }
1078     }
1079
1080     if (mdrunOptions.numStepsCommandline > -2)
1081     {
1082         GMX_LOG(mdlog.info).asParagraph().
1083             appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
1084                        "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
1085     }
1086     /* override nsteps with value set on the commandline */
1087     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1088
1089     if (SIMMASTER(cr))
1090     {
1091         copy_mat(globalState->box, box);
1092     }
1093
1094     if (PAR(cr))
1095     {
1096         gmx_bcast(sizeof(box), box, cr);
1097     }
1098
1099     if (inputrec->cutoff_scheme != ecutsVERLET)
1100     {
1101         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.");
1102     }
1103     /* Update rlist and nstlist. */
1104     prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1105                           useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
1106
1107     const bool prefer1DAnd1PulseDD = (devFlags.enableGpuHaloExchange && useGpuForNonbonded);
1108     // This builder is necessary while we have multi-part construction
1109     // of DD. Before DD is constructed, we use the existence of
1110     // the builder object to indicate that further construction of DD
1111     // is needed.
1112     std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1113     if (useDomainDecomposition)
1114     {
1115         ddBuilder = std::make_unique<DomainDecompositionBuilder>
1116                 (mdlog, cr, domdecOptions, mdrunOptions,
1117                 prefer1DAnd1PulseDD, mtop, *inputrec,
1118                 box, positionsFromStatePointer(globalState.get()));
1119     }
1120     else
1121     {
1122         /* PME, if used, is done on all nodes with 1D decomposition */
1123         cr->npmenodes = 0;
1124         cr->duty      = (DUTY_PP | DUTY_PME);
1125
1126         if (inputrec->ePBC == epbcSCREW)
1127         {
1128             gmx_fatal(FARGS,
1129                       "pbc=screw is only implemented with domain decomposition");
1130         }
1131     }
1132
1133     // Produce the task assignment for this rank.
1134     GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder;
1135     GpuTaskAssignments        gpuTaskAssignments =
1136         gpuTaskAssignmentsBuilder.build(gpuIdsToUse,
1137                                         userGpuTaskAssignment,
1138                                         *hwinfo,
1139                                         cr,
1140                                         ms,
1141                                         physicalNodeComm,
1142                                         nonbondedTarget,
1143                                         pmeTarget,
1144                                         bondedTarget,
1145                                         updateTarget,
1146                                         useGpuForNonbonded,
1147                                         useGpuForPme,
1148                                         thisRankHasDuty(cr, DUTY_PP),
1149                                         // TODO cr->duty & DUTY_PME should imply that a PME
1150                                         // algorithm is active, but currently does not.
1151                                         EEL_PME(inputrec->coulombtype) &&
1152                                         thisRankHasDuty(cr, DUTY_PME));
1153
1154     const bool printHostName = (cr->nnodes > 1);
1155     gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode);
1156
1157     // If the user chose a task assignment, give them some hints
1158     // where appropriate.
1159     if (!userGpuTaskAssignment.empty())
1160     {
1161         gpuTaskAssignments.logPerformanceHints(mdlog,
1162                                                ssize(gpuIdsToUse));
1163     }
1164
1165     // Get the device handles for the modules, nullptr when no task is assigned.
1166     gmx_device_info_t *nonbondedDeviceInfo   = gpuTaskAssignments.initNonbondedDevice(cr);
1167     gmx_device_info_t *pmeDeviceInfo         = gpuTaskAssignments.initPmeDevice();
1168
1169     // TODO Initialize GPU streams here.
1170
1171     // TODO Currently this is always built, yet DD partition code
1172     // checks if it is built before using it. Probably it should
1173     // become an MDModule that is made only when another module
1174     // requires it (e.g. pull, CompEl, density fitting), so that we
1175     // don't update the local atom sets unilaterally every step.
1176     LocalAtomSetManager atomSets;
1177     if (ddBuilder)
1178     {
1179         // TODO Pass the GPU streams to ddBuilder to use in buffer
1180         // transfers (e.g. halo exchange)
1181         cr->dd = ddBuilder->build(&atomSets);
1182         // The builder's job is done, so destruct it
1183         ddBuilder.reset(nullptr);
1184         // Note that local state still does not exist yet.
1185     }
1186
1187     if (PAR(cr))
1188     {
1189         /* After possible communicator splitting in make_dd_communicators.
1190          * we can set up the intra/inter node communication.
1191          */
1192         gmx_setup_nodecomm(fplog, cr);
1193     }
1194
1195 #if GMX_MPI
1196     if (isMultiSim(ms))
1197     {
1198         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1199                 "This is simulation %d out of %d running as a composite GROMACS\n"
1200                 "multi-simulation job. Setup for this simulation:\n",
1201                 ms->sim, ms->nsim);
1202     }
1203     GMX_LOG(mdlog.warning).appendTextFormatted(
1204             "Using %d MPI %s\n",
1205             cr->nnodes,
1206 #if GMX_THREAD_MPI
1207             cr->nnodes == 1 ? "thread" : "threads"
1208 #else
1209             cr->nnodes == 1 ? "process" : "processes"
1210 #endif
1211             );
1212     fflush(stderr);
1213 #endif
1214
1215     // If mdrun -pin auto honors any affinity setting that already
1216     // exists. If so, it is nice to provide feedback about whether
1217     // that existing affinity setting was from OpenMP or something
1218     // else, so we run this code both before and after we initialize
1219     // the OpenMP support.
1220     gmx_check_thread_affinity_set(mdlog,
1221                                   &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1222     /* Check and update the number of OpenMP threads requested */
1223     checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1224                                             pmeRunMode, mtop, *inputrec);
1225
1226     gmx_omp_nthreads_init(mdlog, cr,
1227                           hwinfo->nthreads_hw_avail,
1228                           physicalNodeComm.size_,
1229                           hw_opt.nthreads_omp,
1230                           hw_opt.nthreads_omp_pme,
1231                           !thisRankHasDuty(cr, DUTY_PP));
1232
1233     // Enable FP exception detection, but not in
1234     // Release mode and not for compilers with known buggy FP
1235     // exception support (clang with any optimization) or suspected
1236     // buggy FP exception support (gcc 7.* with optimization).
1237 #if !defined NDEBUG && \
1238     !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1239     && defined __OPTIMIZE__)
1240     const bool bEnableFPE = true;
1241 #else
1242     const bool bEnableFPE = false;
1243 #endif
1244     //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1245     if (bEnableFPE)
1246     {
1247         gmx_feenableexcept();
1248     }
1249
1250     /* Now that we know the setup is consistent, check for efficiency */
1251     check_resource_division_efficiency(hwinfo,
1252                                        gpuTaskAssignments.thisRankHasAnyGpuTask(),
1253                                        mdrunOptions.ntompOptionIsSet,
1254                                        cr,
1255                                        mdlog);
1256
1257     /* getting number of PP/PME threads on this MPI / tMPI rank.
1258        PME: env variable should be read only on one node to make sure it is
1259        identical everywhere;
1260      */
1261     const int numThreadsOnThisRank =
1262         thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1263     checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1264                                   *hwinfo->hardwareTopology,
1265                                   physicalNodeComm, mdlog);
1266
1267     // Enable Peer access between GPUs where available
1268     // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1269     // any of the GPU communication features are active.
1270     if (DOMAINDECOMP(cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI &&
1271         (devFlags.enableGpuHaloExchange || devFlags.enableGpuPmePPComm))
1272     {
1273         setupGpuDevicePeerAccess(gpuIdsToUse, mdlog);
1274     }
1275
1276     if (hw_opt.threadAffinity != ThreadAffinity::Off)
1277     {
1278         /* Before setting affinity, check whether the affinity has changed
1279          * - which indicates that probably the OpenMP library has changed it
1280          * since we first checked).
1281          */
1282         gmx_check_thread_affinity_set(mdlog,
1283                                       &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1284
1285         int numThreadsOnThisNode, intraNodeThreadOffset;
1286         analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1287                                  &intraNodeThreadOffset);
1288
1289         /* Set the CPU affinity */
1290         gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1291                                 numThreadsOnThisRank, numThreadsOnThisNode,
1292                                 intraNodeThreadOffset, nullptr);
1293     }
1294
1295     if (mdrunOptions.timingOptions.resetStep > -1)
1296     {
1297         GMX_LOG(mdlog.info).asParagraph().
1298             appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1299     }
1300     wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1301
1302     if (PAR(cr))
1303     {
1304         /* Master synchronizes its value of reset_counters with all nodes
1305          * including PME only nodes */
1306         int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1307         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1308         wcycle_set_reset_counters(wcycle, reset_counters);
1309     }
1310
1311     // Membrane embedding must be initialized before we call init_forcerec()
1312     if (doMembed)
1313     {
1314         if (MASTER(cr))
1315         {
1316             fprintf(stderr, "Initializing membed");
1317         }
1318         /* Note that membed cannot work in parallel because mtop is
1319          * changed here. Fix this if we ever want to make it run with
1320          * multiple ranks. */
1321         membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1322                              &mdrunOptions
1323                                  .checkpointOptions.period);
1324     }
1325
1326     const bool                   thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1327     std::unique_ptr<MDAtoms>     mdAtoms;
1328     std::unique_ptr<gmx_vsite_t> vsite;
1329
1330     t_nrnb nrnb;
1331     if (thisRankHasDuty(cr, DUTY_PP))
1332     {
1333         mdModulesNotifier.notify(*cr);
1334         mdModulesNotifier.notify(&atomSets);
1335         mdModulesNotifier.notify(PeriodicBoundaryConditionType {inputrec->ePBC});
1336         mdModulesNotifier.notify(SimulationTimeStep { inputrec->delta_t });
1337         /* Initiate forcerecord */
1338         fr                 = new t_forcerec;
1339         fr->forceProviders = mdModules_->initForceProviders();
1340         init_forcerec(fplog, mdlog, fr, fcd,
1341                       inputrec, &mtop, cr, box,
1342                       opt2fn("-table", filenames.size(), filenames.data()),
1343                       opt2fn("-tablep", filenames.size(), filenames.data()),
1344                       opt2fns("-tableb", filenames.size(), filenames.data()),
1345                       *hwinfo, nonbondedDeviceInfo,
1346                       useGpuForBonded,
1347                       pmeRunMode == PmeRunMode::GPU && !thisRankHasDuty(cr, DUTY_PME),
1348                       pforce,
1349                       wcycle);
1350
1351         // TODO Move this to happen during domain decomposition setup,
1352         // once stream and event handling works well with that.
1353         // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
1354         if (havePPDomainDecomposition(cr) && prefer1DAnd1PulseDD && is1DAnd1PulseDD(*cr->dd))
1355         {
1356             GMX_RELEASE_ASSERT(devFlags.enableGpuBufferOps, "Must use GMX_USE_GPU_BUFFER_OPS=1 to use GMX_GPU_DD_COMMS=1");
1357             void *streamLocal              = Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local);
1358             void *streamNonLocal           = Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal);
1359             void *coordinatesOnDeviceEvent = fr->nbv->get_x_on_device_event();
1360             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1361                     "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the GMX_GPU_DD_COMMS environment variable.");
1362             cr->dd->gpuHaloExchange = std::make_unique<GpuHaloExchange>(cr->dd, cr->mpi_comm_mysim, streamLocal,
1363                                                                         streamNonLocal, coordinatesOnDeviceEvent);
1364         }
1365
1366         /* Initialize the mdAtoms structure.
1367          * mdAtoms is not filled with atom data,
1368          * as this can not be done now with domain decomposition.
1369          */
1370         mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1371         if (globalState && thisRankHasPmeGpuTask)
1372         {
1373             // The pinning of coordinates in the global state object works, because we only use
1374             // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1375             // points to the global state object without DD.
1376             // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1377             // which should also perform the pinning.
1378             changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1379         }
1380
1381         /* Initialize the virtual site communication */
1382         vsite = initVsite(mtop, cr);
1383
1384         calc_shifts(box, fr->shift_vec);
1385
1386         /* With periodic molecules the charge groups should be whole at start up
1387          * and the virtual sites should not be far from their proper positions.
1388          */
1389         if (!inputrec->bContinuation && MASTER(cr) &&
1390             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1391         {
1392             /* Make molecules whole at start of run */
1393             if (fr->ePBC != epbcNONE)
1394             {
1395                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1396             }
1397             if (vsite)
1398             {
1399                 /* Correct initial vsite positions are required
1400                  * for the initial distribution in the domain decomposition
1401                  * and for the initial shell prediction.
1402                  */
1403                 constructVsitesGlobal(mtop, globalState->x);
1404             }
1405         }
1406
1407         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1408         {
1409             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1410             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1411         }
1412     }
1413     else
1414     {
1415         /* This is a PME only node */
1416
1417         GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1418
1419         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1420         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1421     }
1422
1423     gmx_pme_t *sepPmeData = nullptr;
1424     // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1425     GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1426     gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1427
1428     // TODO should live in ewald module once its testing is improved
1429     //
1430     // Later, this program could contain kernels that might be later
1431     // re-used as auto-tuning progresses, or subsequent simulations
1432     // are invoked.
1433     PmeGpuProgramStorage pmeGpuProgram;
1434     if (thisRankHasPmeGpuTask)
1435     {
1436         pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1437     }
1438
1439     /* Initiate PME if necessary,
1440      * either on all nodes or on dedicated PME nodes only. */
1441     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1442     {
1443         if (mdAtoms && mdAtoms->mdatoms())
1444         {
1445             nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1446             if (EVDW_PME(inputrec->vdwtype))
1447             {
1448                 nTypePerturbed   = mdAtoms->mdatoms()->nTypePerturbed;
1449             }
1450         }
1451         if (cr->npmenodes > 0)
1452         {
1453             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1454             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1455             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1456         }
1457
1458         if (thisRankHasDuty(cr, DUTY_PME))
1459         {
1460             try
1461             {
1462                 pmedata = gmx_pme_init(cr,
1463                                        getNumPmeDomains(cr->dd),
1464                                        inputrec,
1465                                        nChargePerturbed != 0, nTypePerturbed != 0,
1466                                        mdrunOptions.reproducible,
1467                                        ewaldcoeff_q, ewaldcoeff_lj,
1468                                        gmx_omp_nthreads_get(emntPME),
1469                                        pmeRunMode, nullptr,
1470                                        pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1471             }
1472             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1473         }
1474     }
1475
1476
1477     if (EI_DYNAMICS(inputrec->eI))
1478     {
1479         /* Turn on signal handling on all nodes */
1480         /*
1481          * (A user signal from the PME nodes (if any)
1482          * is communicated to the PP nodes.
1483          */
1484         signal_handler_install();
1485     }
1486
1487     pull_t *pull_work = nullptr;
1488     if (thisRankHasDuty(cr, DUTY_PP))
1489     {
1490         /* Assumes uniform use of the number of OpenMP threads */
1491         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1492
1493         if (inputrec->bPull)
1494         {
1495             /* Initialize pull code */
1496             pull_work =
1497                 init_pull(fplog, inputrec->pull, inputrec,
1498                           &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1499             if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1500             {
1501                 initPullHistory(pull_work, &observablesHistory);
1502             }
1503             if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1504             {
1505                 init_pull_output_files(pull_work,
1506                                        filenames.size(), filenames.data(), oenv,
1507                                        startingBehavior);
1508             }
1509         }
1510
1511         std::unique_ptr<EnforcedRotation> enforcedRotation;
1512         if (inputrec->bRot)
1513         {
1514             /* Initialize enforced rotation code */
1515             enforcedRotation = init_rot(fplog,
1516                                         inputrec,
1517                                         filenames.size(),
1518                                         filenames.data(),
1519                                         cr,
1520                                         &atomSets,
1521                                         globalState.get(),
1522                                         &mtop,
1523                                         oenv,
1524                                         mdrunOptions,
1525                                         startingBehavior);
1526         }
1527
1528         t_swap *swap = nullptr;
1529         if (inputrec->eSwapCoords != eswapNO)
1530         {
1531             /* Initialize ion swapping code */
1532             swap = init_swapcoords(fplog, inputrec,
1533                                    opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1534                                    &mtop, globalState.get(), &observablesHistory,
1535                                    cr, &atomSets, oenv, mdrunOptions,
1536                                    startingBehavior);
1537         }
1538
1539         /* Let makeConstraints know whether we have essential dynamics constraints. */
1540         auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics,
1541                                       fplog, *mdAtoms->mdatoms(),
1542                                       cr, ms, &nrnb, wcycle, fr->bMolPBC);
1543
1544         /* Energy terms and groups */
1545         gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(), inputrec->fepvals->n_lambda);
1546
1547         /* Kinetic energy data */
1548         gmx_ekindata_t ekind;
1549         init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1550
1551         /* Set up interactive MD (IMD) */
1552         auto imdSession = makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1553                                          MASTER(cr) ? globalState->x.rvec_array() : nullptr,
1554                                          filenames.size(), filenames.data(), oenv, mdrunOptions.imdOptions,
1555                                          startingBehavior);
1556
1557         if (DOMAINDECOMP(cr))
1558         {
1559             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1560             /* This call is not included in init_domain_decomposition mainly
1561              * because fr->cginfo_mb is set later.
1562              */
1563             dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1564                             domdecOptions.checkBondedInteractions,
1565                             fr->cginfo_mb);
1566         }
1567
1568         const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
1569                     false,
1570                     inputrec, doRerun, vsite.get(), ms, replExParams,
1571                     fcd, static_cast<int>(filenames.size()), filenames.data(),
1572                     &observablesHistory, membed);
1573
1574         const bool useModularSimulator = inputIsCompatibleWithModularSimulator && !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
1575
1576         // TODO This is not the right place to manage the lifetime of
1577         // this data structure, but currently it's the easiest way to
1578         // make it work.
1579         MdrunScheduleWorkload runScheduleWork;
1580         // Also populates the simulation constant workload description.
1581         runScheduleWork.simulationWork = createSimulationWorkload(useGpuForNonbonded,
1582                                                                   pmeRunMode,
1583                                                                   useGpuForBonded,
1584                                                                   useGpuForUpdate,
1585                                                                   devFlags.enableGpuBufferOps,
1586                                                                   devFlags.enableGpuHaloExchange,
1587                                                                   devFlags.enableGpuPmePPComm);
1588
1589         std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
1590         if (gpusWereDetected && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME)) || runScheduleWork.simulationWork.useGpuBufferOps))
1591         {
1592             const void         *pmeStream      = pme_gpu_get_device_stream(fr->pmedata);
1593             const void         *localStream    = fr->nbv->gpu_nbv != nullptr ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::Local) : nullptr;
1594             const void         *nonLocalStream = fr->nbv->gpu_nbv != nullptr ? Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, InteractionLocality::NonLocal) : nullptr;
1595             const void         *deviceContext  = pme_gpu_get_device_context(fr->pmedata);
1596             const int           paddingSize    = pme_gpu_get_padding_size(fr->pmedata);
1597             GpuApiCallBehavior  transferKind   = (inputrec->eI == eiMD && !doRerun && !useModularSimulator) ? GpuApiCallBehavior::Async : GpuApiCallBehavior::Sync;
1598
1599             stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(pmeStream,
1600                                                                      localStream,
1601                                                                      nonLocalStream,
1602                                                                      deviceContext,
1603                                                                      transferKind,
1604                                                                      paddingSize);
1605             fr->stateGpu = stateGpu.get();
1606         }
1607
1608         GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1609         SimulatorBuilder simulatorBuilder;
1610
1611         // build and run simulator object based on user-input
1612         auto simulator = simulatorBuilder.build(
1613                     inputIsCompatibleWithModularSimulator,
1614                     fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1615                     oenv,
1616                     mdrunOptions,
1617                     startingBehavior,
1618                     vsite.get(), constr.get(),
1619                     enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1620                     deform.get(),
1621                     mdModules_->outputProvider(),
1622                     mdModules_->notifier(),
1623                     inputrec, imdSession.get(), pull_work, swap, &mtop,
1624                     fcd,
1625                     globalState.get(),
1626                     &observablesHistory,
1627                     mdAtoms.get(), &nrnb, wcycle, fr,
1628                     &enerd,
1629                     &ekind,
1630                     &runScheduleWork,
1631                     replExParams,
1632                     membed,
1633                     walltime_accounting,
1634                     std::move(stopHandlerBuilder_),
1635                     doRerun);
1636         simulator->run();
1637
1638         if (inputrec->bPull)
1639         {
1640             finish_pull(pull_work);
1641         }
1642         finish_swapcoords(swap);
1643     }
1644     else
1645     {
1646         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1647         /* do PME only */
1648         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1649         gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1650     }
1651
1652     wallcycle_stop(wcycle, ewcRUN);
1653
1654     /* Finish up, write some stuff
1655      * if rerunMD, don't write last frame again
1656      */
1657     finish_run(fplog, mdlog, cr,
1658                inputrec, &nrnb, wcycle, walltime_accounting,
1659                fr ? fr->nbv.get() : nullptr,
1660                pmedata,
1661                EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1662
1663     // clean up cycle counter
1664     wallcycle_destroy(wcycle);
1665
1666 // Free PME data
1667     if (pmedata)
1668     {
1669         gmx_pme_destroy(pmedata);
1670         pmedata = nullptr;
1671     }
1672
1673     // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1674     // before we destroy the GPU context(s) in free_gpu_resources().
1675     // Pinned buffers are associated with contexts in CUDA.
1676     // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1677     mdAtoms.reset(nullptr);
1678     globalState.reset(nullptr);
1679     mdModules_.reset(nullptr);   // destruct force providers here as they might also use the GPU
1680
1681     /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1682     free_gpu_resources(fr, physicalNodeComm, hwinfo->gpu_info);
1683     free_gpu(nonbondedDeviceInfo);
1684     free_gpu(pmeDeviceInfo);
1685     done_forcerec(fr, mtop.molblock.size());
1686     sfree(fcd);
1687
1688     if (doMembed)
1689     {
1690         free_membed(membed);
1691     }
1692
1693     /* Does what it says */
1694     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1695     walltime_accounting_destroy(walltime_accounting);
1696
1697     // Ensure log file content is written
1698     if (logFileHandle)
1699     {
1700         gmx_fio_flush(logFileHandle);
1701     }
1702
1703     /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1704      * exceptions were enabled before function was called. */
1705     if (bEnableFPE)
1706     {
1707         gmx_fedisableexcept();
1708     }
1709
1710     auto rc = static_cast<int>(gmx_get_stop_condition());
1711
1712 #if GMX_THREAD_MPI
1713     /* we need to join all threads. The sub-threads join when they
1714        exit this function, but the master thread needs to be told to
1715        wait for that. */
1716     if (PAR(cr) && MASTER(cr))
1717     {
1718         tMPI_Finalize();
1719     }
1720 #endif
1721     return rc;
1722 }
1723
1724 Mdrunner::~Mdrunner()
1725 {
1726     // Clean up of the Manager.
1727     // This will end up getting called on every thread-MPI rank, which is unnecessary,
1728     // but okay as long as threads synchronize some time before adding or accessing
1729     // a new set of restraints.
1730     if (restraintManager_)
1731     {
1732         restraintManager_->clear();
1733         GMX_ASSERT(restraintManager_->countRestraints() == 0,
1734                    "restraints added during runner life time should be cleared at runner destruction.");
1735     }
1736 };
1737
1738 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1739                             const std::string                        &name)
1740 {
1741     GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1742     // Not sure if this should be logged through the md logger or something else,
1743     // but it is helpful to have some sort of INFO level message sent somewhere.
1744     //    std::cout << "Registering restraint named " << name << std::endl;
1745
1746     // When multiple restraints are used, it may be wasteful to register them separately.
1747     // Maybe instead register an entire Restraint Manager as a force provider.
1748     restraintManager_->addToSpec(std::move(puller),
1749                                  name);
1750 }
1751
1752 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1753     : mdModules_(std::move(mdModules))
1754 {
1755 }
1756
1757 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1758
1759 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1760 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1761
1762 class Mdrunner::BuilderImplementation
1763 {
1764     public:
1765         BuilderImplementation() = delete;
1766         BuilderImplementation(std::unique_ptr<MDModules>           mdModules,
1767                               compat::not_null<SimulationContext*> context);
1768         ~BuilderImplementation();
1769
1770         BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1771                                                     real                forceWarningThreshold,
1772                                                     StartingBehavior    startingBehavior);
1773
1774         void addDomdec(const DomdecOptions &options);
1775
1776         void addVerletList(int nstlist);
1777
1778         void addReplicaExchange(const ReplicaExchangeParameters &params);
1779
1780         void addNonBonded(const char* nbpu_opt);
1781
1782         void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1783
1784         void addBondedTaskAssignment(const char* bonded_opt);
1785
1786         void addUpdateTaskAssignment(const char* update_opt);
1787
1788         void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1789
1790         void addFilenames(ArrayRef <const t_filenm> filenames);
1791
1792         void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1793
1794         void addLogFile(t_fileio *logFileHandle);
1795
1796         void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1797
1798         Mdrunner build();
1799
1800     private:
1801
1802         // Default parameters copied from runner.h
1803         // \todo Clarify source(s) of default parameters.
1804
1805         const char* nbpu_opt_          = nullptr;
1806         const char* pme_opt_           = nullptr;
1807         const char* pme_fft_opt_       = nullptr;
1808         const char *bonded_opt_        = nullptr;
1809         const char *update_opt_        = nullptr;
1810
1811         MdrunOptions                          mdrunOptions_;
1812
1813         DomdecOptions                         domdecOptions_;
1814
1815         ReplicaExchangeParameters             replicaExchangeParameters_;
1816
1817         //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1818         int         nstlist_ = 0;
1819
1820         //! Multisim communicator handle.
1821         gmx_multisim_t *multiSimulation_;
1822
1823         //! mdrun communicator
1824         MPI_Comm communicator_ = MPI_COMM_NULL;
1825
1826         //! Print a warning if any force is larger than this (in kJ/mol nm).
1827         real forceWarningThreshold_ = -1;
1828
1829         //! Whether the simulation will start afresh, or restart with/without appending.
1830         StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1831
1832         //! The modules that comprise the functionality of mdrun.
1833         std::unique_ptr<MDModules> mdModules_;
1834
1835         //! \brief Parallelism information.
1836         gmx_hw_opt_t hardwareOptions_;
1837
1838         //! filename options for simulation.
1839         ArrayRef<const t_filenm> filenames_;
1840
1841         /*! \brief Handle to output environment.
1842          *
1843          * \todo gmx_output_env_t needs lifetime management.
1844          */
1845         gmx_output_env_t*    outputEnvironment_ = nullptr;
1846
1847         /*! \brief Non-owning handle to MD log file.
1848          *
1849          * \todo Context should own output facilities for client.
1850          * \todo Improve log file handle management.
1851          * \internal
1852          * Code managing the FILE* relies on the ability to set it to
1853          * nullptr to check whether the filehandle is valid.
1854          */
1855         t_fileio* logFileHandle_ = nullptr;
1856
1857         /*!
1858          * \brief Builder for simulation stop signal handler.
1859          */
1860         std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1861 };
1862
1863 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules>           mdModules,
1864                                                        compat::not_null<SimulationContext*> context) :
1865     mdModules_(std::move(mdModules))
1866 {
1867     communicator_    = context->communicator_;
1868     multiSimulation_ = context->multiSimulation_.get();
1869 }
1870
1871 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1872
1873 Mdrunner::BuilderImplementation &
1874 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions    &options,
1875                                                       const real             forceWarningThreshold,
1876                                                       const StartingBehavior startingBehavior)
1877 {
1878     mdrunOptions_          = options;
1879     forceWarningThreshold_ = forceWarningThreshold;
1880     startingBehavior_      = startingBehavior;
1881     return *this;
1882 }
1883
1884 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1885 {
1886     domdecOptions_ = options;
1887 }
1888
1889 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1890 {
1891     nstlist_ = nstlist;
1892 }
1893
1894 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters &params)
1895 {
1896     replicaExchangeParameters_ = params;
1897 }
1898
1899 Mdrunner Mdrunner::BuilderImplementation::build()
1900 {
1901     auto newRunner = Mdrunner(std::move(mdModules_));
1902
1903     newRunner.mdrunOptions          = mdrunOptions_;
1904     newRunner.pforce                = forceWarningThreshold_;
1905     newRunner.startingBehavior      = startingBehavior_;
1906     newRunner.domdecOptions         = domdecOptions_;
1907
1908     // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1909     newRunner.hw_opt          = hardwareOptions_;
1910
1911     // No invariant to check. This parameter exists to optionally override other behavior.
1912     newRunner.nstlist_cmdline = nstlist_;
1913
1914     newRunner.replExParams    = replicaExchangeParameters_;
1915
1916     newRunner.filenames = filenames_;
1917
1918     newRunner.communicator = communicator_;
1919
1920     // nullptr is a valid value for the multisim handle
1921     newRunner.ms = multiSimulation_;
1922
1923     // \todo Clarify ownership and lifetime management for gmx_output_env_t
1924     // \todo Update sanity checking when output environment has clearly specified invariants.
1925     // Initialization and default values for oenv are not well specified in the current version.
1926     if (outputEnvironment_)
1927     {
1928         newRunner.oenv  = outputEnvironment_;
1929     }
1930     else
1931     {
1932         GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1933     }
1934
1935     newRunner.logFileHandle = logFileHandle_;
1936
1937     if (nbpu_opt_)
1938     {
1939         newRunner.nbpu_opt    = nbpu_opt_;
1940     }
1941     else
1942     {
1943         GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1944     }
1945
1946     if (pme_opt_ && pme_fft_opt_)
1947     {
1948         newRunner.pme_opt     = pme_opt_;
1949         newRunner.pme_fft_opt = pme_fft_opt_;
1950     }
1951     else
1952     {
1953         GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1954     }
1955
1956     if (bonded_opt_)
1957     {
1958         newRunner.bonded_opt = bonded_opt_;
1959     }
1960     else
1961     {
1962         GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1963     }
1964
1965     if (update_opt_)
1966     {
1967         newRunner.update_opt = update_opt_;
1968     }
1969     else
1970     {
1971         GMX_THROW(gmx::APIError("MdrunnerBuilder::addUpdateTaskAssignment() is required before build()  "));
1972     }
1973
1974
1975     newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1976
1977     if (stopHandlerBuilder_)
1978     {
1979         newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1980     }
1981     else
1982     {
1983         newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1984     }
1985
1986     return newRunner;
1987 }
1988
1989 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1990 {
1991     nbpu_opt_ = nbpu_opt;
1992 }
1993
1994 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1995                                              const char* pme_fft_opt)
1996 {
1997     pme_opt_     = pme_opt;
1998     pme_fft_opt_ = pme_fft_opt;
1999 }
2000
2001 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2002 {
2003     bonded_opt_ = bonded_opt;
2004 }
2005
2006 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2007 {
2008     update_opt_ = update_opt;
2009 }
2010
2011 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2012 {
2013     hardwareOptions_ = hardwareOptions;
2014 }
2015
2016 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2017 {
2018     filenames_ = filenames;
2019 }
2020
2021 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2022 {
2023     outputEnvironment_ = outputEnvironment;
2024 }
2025
2026 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
2027 {
2028     logFileHandle_ = logFileHandle;
2029 }
2030
2031 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2032 {
2033     stopHandlerBuilder_ = std::move(builder);
2034 }
2035
2036 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules>           mdModules,
2037                                  compat::not_null<SimulationContext*> context) :
2038     impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context)}
2039 {
2040 }
2041
2042 MdrunnerBuilder::~MdrunnerBuilder() = default;
2043
2044 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions    &options,
2045                                                       real                   forceWarningThreshold,
2046                                                       const StartingBehavior startingBehavior)
2047 {
2048     impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2049     return *this;
2050 }
2051
2052 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
2053 {
2054     impl_->addDomdec(options);
2055     return *this;
2056 }
2057
2058 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
2059 {
2060     impl_->addVerletList(nstlist);
2061     return *this;
2062 }
2063
2064 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters &params)
2065 {
2066     impl_->addReplicaExchange(params);
2067     return *this;
2068 }
2069
2070 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2071 {
2072     impl_->addNonBonded(nbpu_opt);
2073     return *this;
2074 }
2075
2076 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
2077                                                     const char* pme_fft_opt)
2078 {
2079     // The builder method may become more general in the future, but in this version,
2080     // parameters for PME electrostatics are both required and the only parameters
2081     // available.
2082     if (pme_opt && pme_fft_opt)
2083     {
2084         impl_->addPME(pme_opt, pme_fft_opt);
2085     }
2086     else
2087     {
2088         GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2089     }
2090     return *this;
2091 }
2092
2093 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2094 {
2095     impl_->addBondedTaskAssignment(bonded_opt);
2096     return *this;
2097 }
2098
2099 MdrunnerBuilder &MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2100 {
2101     impl_->addUpdateTaskAssignment(update_opt);
2102     return *this;
2103 }
2104
2105 Mdrunner MdrunnerBuilder::build()
2106 {
2107     return impl_->build();
2108 }
2109
2110 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2111 {
2112     impl_->addHardwareOptions(hardwareOptions);
2113     return *this;
2114 }
2115
2116 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2117 {
2118     impl_->addFilenames(filenames);
2119     return *this;
2120 }
2121
2122 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2123 {
2124     impl_->addOutputEnvironment(outputEnvironment);
2125     return *this;
2126 }
2127
2128 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2129 {
2130     impl_->addLogFile(logFileHandle);
2131     return *this;
2132 }
2133
2134 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2135 {
2136     impl_->addStopHandlerBuilder(std::move(builder));
2137     return *this;
2138 }
2139
2140 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2141
2142 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;
2143
2144 } // namespace gmx