Split simulationWork.useGpuBufferOps into separate x and f flags
[alexxy/gromacs.git] / src / gromacs / mdrun / runner.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011-2019,2020,2021, 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/makebondedlinks.h"
66 #include "gromacs/domdec/partition.h"
67 #include "gromacs/domdec/reversetopology.h"
68 #include "gromacs/ewald/ewald_utils.h"
69 #include "gromacs/ewald/pme.h"
70 #include "gromacs/ewald/pme_gpu_program.h"
71 #include "gromacs/ewald/pme_only.h"
72 #include "gromacs/ewald/pme_pp_comm_gpu.h"
73 #include "gromacs/fileio/checkpoint.h"
74 #include "gromacs/fileio/gmxfio.h"
75 #include "gromacs/fileio/oenv.h"
76 #include "gromacs/fileio/tpxio.h"
77 #include "gromacs/gmxlib/network.h"
78 #include "gromacs/gmxlib/nrnb.h"
79 #include "gromacs/gpu_utils/device_stream_manager.h"
80 #include "gromacs/hardware/cpuinfo.h"
81 #include "gromacs/hardware/detecthardware.h"
82 #include "gromacs/hardware/device_management.h"
83 #include "gromacs/hardware/hardwaretopology.h"
84 #include "gromacs/hardware/printhardware.h"
85 #include "gromacs/imd/imd.h"
86 #include "gromacs/listed_forces/disre.h"
87 #include "gromacs/listed_forces/listed_forces_gpu.h"
88 #include "gromacs/listed_forces/listed_forces.h"
89 #include "gromacs/listed_forces/orires.h"
90 #include "gromacs/math/functions.h"
91 #include "gromacs/math/utilities.h"
92 #include "gromacs/math/vec.h"
93 #include "gromacs/mdlib/boxdeformation.h"
94 #include "gromacs/mdlib/broadcaststructs.h"
95 #include "gromacs/mdlib/calc_verletbuf.h"
96 #include "gromacs/mdlib/dispersioncorrection.h"
97 #include "gromacs/mdlib/enerdata_utils.h"
98 #include "gromacs/mdlib/force.h"
99 #include "gromacs/mdlib/forcerec.h"
100 #include "gromacs/mdlib/gmx_omp_nthreads.h"
101 #include "gromacs/mdlib/gpuforcereduction.h"
102 #include "gromacs/mdlib/makeconstraints.h"
103 #include "gromacs/mdlib/md_support.h"
104 #include "gromacs/mdlib/mdatoms.h"
105 #include "gromacs/mdlib/sighandler.h"
106 #include "gromacs/mdlib/stophandler.h"
107 #include "gromacs/mdlib/tgroup.h"
108 #include "gromacs/mdlib/updategroups.h"
109 #include "gromacs/mdlib/vsite.h"
110 #include "gromacs/mdrun/mdmodules.h"
111 #include "gromacs/mdrun/simulationcontext.h"
112 #include "gromacs/mdrun/simulationinput.h"
113 #include "gromacs/mdrun/simulationinputhandle.h"
114 #include "gromacs/mdrunutility/handlerestart.h"
115 #include "gromacs/mdrunutility/logging.h"
116 #include "gromacs/mdrunutility/multisim.h"
117 #include "gromacs/mdrunutility/printtime.h"
118 #include "gromacs/mdrunutility/threadaffinity.h"
119 #include "gromacs/mdtypes/checkpointdata.h"
120 #include "gromacs/mdtypes/commrec.h"
121 #include "gromacs/mdtypes/enerdata.h"
122 #include "gromacs/mdtypes/fcdata.h"
123 #include "gromacs/mdtypes/forcerec.h"
124 #include "gromacs/mdtypes/group.h"
125 #include "gromacs/mdtypes/inputrec.h"
126 #include "gromacs/mdtypes/interaction_const.h"
127 #include "gromacs/mdtypes/md_enums.h"
128 #include "gromacs/mdtypes/mdatom.h"
129 #include "gromacs/mdtypes/mdrunoptions.h"
130 #include "gromacs/mdtypes/observableshistory.h"
131 #include "gromacs/mdtypes/observablesreducer.h"
132 #include "gromacs/mdtypes/simulation_workload.h"
133 #include "gromacs/mdtypes/state.h"
134 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
135 #include "gromacs/modularsimulator/modularsimulator.h"
136 #include "gromacs/nbnxm/gpu_data_mgmt.h"
137 #include "gromacs/nbnxm/nbnxm.h"
138 #include "gromacs/nbnxm/pairlist_tuning.h"
139 #include "gromacs/pbcutil/pbc.h"
140 #include "gromacs/pulling/output.h"
141 #include "gromacs/pulling/pull.h"
142 #include "gromacs/pulling/pull_rotation.h"
143 #include "gromacs/restraint/manager.h"
144 #include "gromacs/restraint/restraintmdmodule.h"
145 #include "gromacs/restraint/restraintpotential.h"
146 #include "gromacs/swap/swapcoords.h"
147 #include "gromacs/taskassignment/decidegpuusage.h"
148 #include "gromacs/taskassignment/decidesimulationworkload.h"
149 #include "gromacs/taskassignment/resourcedivision.h"
150 #include "gromacs/taskassignment/taskassignment.h"
151 #include "gromacs/taskassignment/usergpuids.h"
152 #include "gromacs/timing/gpu_timing.h"
153 #include "gromacs/timing/wallcycle.h"
154 #include "gromacs/timing/wallcyclereporting.h"
155 #include "gromacs/topology/mtop_util.h"
156 #include "gromacs/trajectory/trajectoryframe.h"
157 #include "gromacs/utility/basenetwork.h"
158 #include "gromacs/utility/cstringutil.h"
159 #include "gromacs/utility/exceptions.h"
160 #include "gromacs/utility/fatalerror.h"
161 #include "gromacs/utility/filestream.h"
162 #include "gromacs/utility/gmxassert.h"
163 #include "gromacs/utility/gmxmpi.h"
164 #include "gromacs/utility/keyvaluetree.h"
165 #include "gromacs/utility/logger.h"
166 #include "gromacs/utility/loggerbuilder.h"
167 #include "gromacs/utility/mdmodulesnotifiers.h"
168 #include "gromacs/utility/physicalnodecommunicator.h"
169 #include "gromacs/utility/pleasecite.h"
170 #include "gromacs/utility/programcontext.h"
171 #include "gromacs/utility/smalloc.h"
172 #include "gromacs/utility/stringutil.h"
173 #include "gromacs/utility/mpiinfo.h"
174
175 #include "isimulator.h"
176 #include "membedholder.h"
177 #include "replicaexchange.h"
178 #include "simulatorbuilder.h"
179
180 namespace gmx
181 {
182
183
184 /*! \brief Manage any development feature flag variables encountered
185  *
186  * The use of dev features indicated by environment variables is
187  * logged in order to ensure that runs with such features enabled can
188  * be identified from their log and standard output. Any cross
189  * dependencies are also checked, and if unsatisfied, a fatal error
190  * issued.
191  *
192  * Note that some development features overrides are applied already here:
193  * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
194  *
195  * \param[in]  mdlog                Logger object.
196  * \param[in]  useGpuForNonbonded   True if the nonbonded task is offloaded in this run.
197  * \param[in]  pmeRunMode           The PME run mode for this run
198  * \returns                         The object populated with development feature flags.
199  */
200 static DevelopmentFeatureFlags manageDevelopmentFeatures(const gmx::MDLogger& mdlog,
201                                                          const bool           useGpuForNonbonded,
202                                                          const PmeRunMode     pmeRunMode)
203 {
204     DevelopmentFeatureFlags devFlags;
205
206     devFlags.enableGpuBufferOps = (GMX_GPU_CUDA || GMX_GPU_SYCL) && useGpuForNonbonded
207                                   && (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
208     devFlags.enableGpuHaloExchange = GMX_MPI && GMX_GPU_CUDA && getenv("GMX_GPU_DD_COMMS") != nullptr;
209     devFlags.forceGpuUpdateDefault = (getenv("GMX_FORCE_UPDATE_DEFAULT_GPU") != nullptr) || GMX_FAHCORE;
210     devFlags.enableGpuPmePPComm = GMX_MPI && GMX_GPU_CUDA && getenv("GMX_GPU_PME_PP_COMMS") != nullptr;
211
212     // Direct GPU comm path is being used with CUDA_AWARE_MPI
213     // make sure underlying MPI implementation is CUDA-aware
214     if (!GMX_THREAD_MPI && (devFlags.enableGpuPmePPComm || devFlags.enableGpuHaloExchange))
215     {
216         const bool haveDetectedCudaAwareMpi =
217                 (checkMpiCudaAwareSupport() == CudaAwareMpiStatus::Supported);
218         const bool forceCudaAwareMpi = (getenv("GMX_FORCE_CUDA_AWARE_MPI") != nullptr);
219
220         if (!haveDetectedCudaAwareMpi && forceCudaAwareMpi)
221         {
222             // CUDA-aware support not detected in MPI library but, user has forced it's use
223             GMX_LOG(mdlog.warning)
224                     .asParagraph()
225                     .appendTextFormatted(
226                             "This run has forced use of 'CUDA-aware MPI'. "
227                             "But, GROMACS cannot determine if underlying MPI "
228                             "is CUDA-aware. GROMACS recommends use of latest openMPI version "
229                             "for CUDA-aware support. "
230                             "If you observe failures at runtime, try unsetting "
231                             "GMX_FORCE_CUDA_AWARE_MPI environment variable.");
232         }
233
234         if (haveDetectedCudaAwareMpi || forceCudaAwareMpi)
235         {
236             devFlags.usingCudaAwareMpi = true;
237             GMX_LOG(mdlog.warning)
238                     .asParagraph()
239                     .appendTextFormatted(
240                             "Using CUDA-aware MPI for 'GPU halo exchange' or 'GPU PME-PP "
241                             "communications' feature.");
242         }
243         else
244         {
245             if (devFlags.enableGpuHaloExchange)
246             {
247                 GMX_LOG(mdlog.warning)
248                         .asParagraph()
249                         .appendTextFormatted(
250                                 "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
251                                 "halo exchange' feature will not be enabled as GROMACS couldn't "
252                                 "detect CUDA_aware support in underlying MPI implementation.");
253                 devFlags.enableGpuHaloExchange = false;
254             }
255             if (devFlags.enableGpuPmePPComm)
256             {
257                 GMX_LOG(mdlog.warning)
258                         .asParagraph()
259                         .appendText(
260                                 "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
261                                 "'GPU PME-PP communications' feature will not be enabled as "
262                                 "GROMACS couldn't "
263                                 "detect CUDA_aware support in underlying MPI implementation.");
264                 devFlags.enableGpuPmePPComm = false;
265             }
266
267             GMX_LOG(mdlog.warning)
268                     .asParagraph()
269                     .appendTextFormatted(
270                             "GROMACS recommends use of latest OpenMPI version for CUDA-aware "
271                             "support. "
272                             "If you are certain about CUDA-aware support in your MPI library, "
273                             "you can force it's use by setting environment variable "
274                             " GMX_FORCE_CUDA_AWARE_MPI.");
275         }
276     }
277
278     if (devFlags.enableGpuBufferOps)
279     {
280         GMX_LOG(mdlog.warning)
281                 .asParagraph()
282                 .appendTextFormatted(
283                         "This run uses the 'GPU buffer ops' feature, enabled by the "
284                         "GMX_USE_GPU_BUFFER_OPS environment variable.");
285     }
286
287     if (devFlags.forceGpuUpdateDefault)
288     {
289         GMX_LOG(mdlog.warning)
290                 .asParagraph()
291                 .appendTextFormatted(
292                         "This run will default to '-update gpu' as requested by the "
293                         "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable. GPU update with domain "
294                         "decomposition lacks substantial testing and should be used with caution.");
295     }
296
297     if (devFlags.enableGpuHaloExchange)
298     {
299         if (useGpuForNonbonded)
300         {
301             if (!devFlags.enableGpuBufferOps)
302             {
303                 GMX_LOG(mdlog.warning)
304                         .asParagraph()
305                         .appendTextFormatted(
306                                 "Enabling GPU buffer operations required by GMX_GPU_DD_COMMS "
307                                 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
308                 devFlags.enableGpuBufferOps = true;
309             }
310             GMX_LOG(mdlog.warning)
311                     .asParagraph()
312                     .appendTextFormatted(
313                             "This run has requested the 'GPU halo exchange' feature, enabled by "
314                             "the "
315                             "GMX_GPU_DD_COMMS environment variable.");
316         }
317         else
318         {
319             GMX_LOG(mdlog.warning)
320                     .asParagraph()
321                     .appendTextFormatted(
322                             "GMX_GPU_DD_COMMS environment variable detected, but the 'GPU "
323                             "halo exchange' feature will not be enabled as nonbonded interactions "
324                             "are not offloaded.");
325             devFlags.enableGpuHaloExchange = false;
326         }
327     }
328
329     if (devFlags.enableGpuPmePPComm)
330     {
331         if (pmeRunMode == PmeRunMode::GPU)
332         {
333             if (!devFlags.enableGpuBufferOps)
334             {
335                 GMX_LOG(mdlog.warning)
336                         .asParagraph()
337                         .appendTextFormatted(
338                                 "Enabling GPU buffer operations required by GMX_GPU_PME_PP_COMMS "
339                                 "(equivalent with GMX_USE_GPU_BUFFER_OPS=1).");
340                 devFlags.enableGpuBufferOps = true;
341             }
342             GMX_LOG(mdlog.warning)
343                     .asParagraph()
344                     .appendTextFormatted(
345                             "This run uses the 'GPU PME-PP communications' feature, enabled "
346                             "by the GMX_GPU_PME_PP_COMMS environment variable.");
347         }
348         else
349         {
350             std::string clarification;
351             if (pmeRunMode == PmeRunMode::Mixed)
352             {
353                 clarification =
354                         "PME FFT and gather are not offloaded to the GPU (PME is running in mixed "
355                         "mode).";
356             }
357             else
358             {
359                 clarification = "PME is not offloaded to the GPU.";
360             }
361             GMX_LOG(mdlog.warning)
362                     .asParagraph()
363                     .appendText(
364                             "GMX_GPU_PME_PP_COMMS environment variable detected, but the "
365                             "'GPU PME-PP communications' feature was not enabled as "
366                             + clarification);
367             devFlags.enableGpuPmePPComm = false;
368         }
369     }
370
371     return devFlags;
372 }
373
374 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
375  *
376  * Used to ensure that the master thread does not modify mdrunner during copy
377  * on the spawned threads. */
378 static void threadMpiMdrunnerAccessBarrier()
379 {
380 #if GMX_THREAD_MPI
381     MPI_Barrier(MPI_COMM_WORLD);
382 #endif
383 }
384
385 Mdrunner Mdrunner::cloneOnSpawnedThread() const
386 {
387     auto newRunner = Mdrunner(std::make_unique<MDModules>());
388
389     // All runners in the same process share a restraint manager resource because it is
390     // part of the interface to the client code, which is associated only with the
391     // original thread. Handles to the same resources can be obtained by copy.
392     {
393         newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
394     }
395
396     // Copy members of master runner.
397     // \todo Replace with builder when Simulation context and/or runner phases are better defined.
398     // Ref https://gitlab.com/gromacs/gromacs/-/issues/2587 and https://gitlab.com/gromacs/gromacs/-/issues/2375
399     newRunner.hw_opt    = hw_opt;
400     newRunner.filenames = filenames;
401
402     newRunner.hwinfo_         = hwinfo_;
403     newRunner.oenv            = oenv;
404     newRunner.mdrunOptions    = mdrunOptions;
405     newRunner.domdecOptions   = domdecOptions;
406     newRunner.nbpu_opt        = nbpu_opt;
407     newRunner.pme_opt         = pme_opt;
408     newRunner.pme_fft_opt     = pme_fft_opt;
409     newRunner.bonded_opt      = bonded_opt;
410     newRunner.update_opt      = update_opt;
411     newRunner.nstlist_cmdline = nstlist_cmdline;
412     newRunner.replExParams    = replExParams;
413     newRunner.pforce          = pforce;
414     // Give the spawned thread the newly created valid communicator
415     // for the simulation.
416     newRunner.libraryWorldCommunicator = MPI_COMM_WORLD;
417     newRunner.simulationCommunicator   = MPI_COMM_WORLD;
418     newRunner.ms                       = ms;
419     newRunner.startingBehavior         = startingBehavior;
420     newRunner.stopHandlerBuilder_      = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
421     newRunner.inputHolder_             = inputHolder_;
422
423     threadMpiMdrunnerAccessBarrier();
424
425     return newRunner;
426 }
427
428 /*! \brief The callback used for running on spawned threads.
429  *
430  * Obtains the pointer to the master mdrunner object from the one
431  * argument permitted to the thread-launch API call, copies it to make
432  * a new runner for this thread, reinitializes necessary data, and
433  * proceeds to the simulation. */
434 static void mdrunner_start_fn(const void* arg)
435 {
436     try
437     {
438         const auto* masterMdrunner = reinterpret_cast<const gmx::Mdrunner*>(arg);
439         /* copy the arg list to make sure that it's thread-local. This
440            doesn't copy pointed-to items, of course; fnm, cr and fplog
441            are reset in the call below, all others should be const. */
442         gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
443         mdrunner.mdrunner();
444     }
445     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
446 }
447
448
449 void Mdrunner::spawnThreads(int numThreadsToLaunch)
450 {
451 #if GMX_THREAD_MPI
452     /* now spawn new threads that start mdrunner_start_fn(), while
453        the main thread returns. Thread affinity is handled later. */
454     if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE, mdrunner_start_fn, static_cast<const void*>(this))
455         != TMPI_SUCCESS)
456     {
457         GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
458     }
459
460     // Give the master thread the newly created valid communicator for
461     // the simulation.
462     libraryWorldCommunicator = MPI_COMM_WORLD;
463     simulationCommunicator   = MPI_COMM_WORLD;
464     threadMpiMdrunnerAccessBarrier();
465 #else
466     GMX_UNUSED_VALUE(numThreadsToLaunch);
467     GMX_UNUSED_VALUE(mdrunner_start_fn);
468 #endif
469 }
470
471 } // namespace gmx
472
473 /*! \brief Initialize variables for Verlet scheme simulation */
474 static void prepare_verlet_scheme(FILE*               fplog,
475                                   t_commrec*          cr,
476                                   t_inputrec*         ir,
477                                   int                 nstlist_cmdline,
478                                   const gmx_mtop_t&   mtop,
479                                   const matrix        box,
480                                   bool                makeGpuPairList,
481                                   const gmx::CpuInfo& cpuinfo)
482 {
483     // We checked the cut-offs in grompp, but double-check here.
484     // We have PME+LJcutoff kernels for rcoulomb>rvdw.
485     if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == VanDerWaalsType::Cut)
486     {
487         GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
488                            "With Verlet lists and PME we should have rcoulomb>=rvdw");
489     }
490     else
491     {
492         GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw,
493                            "With Verlet lists and no PME rcoulomb and rvdw should be identical");
494     }
495     /* For NVE simulations, we will retain the initial list buffer */
496     if (EI_DYNAMICS(ir->eI) && ir->verletbuf_tol > 0
497         && !(EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No))
498     {
499         /* Update the Verlet buffer size for the current run setup */
500
501         /* Here we assume SIMD-enabled kernels are being used. But as currently
502          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
503          * and 4x2 gives a larger buffer than 4x4, this is ok.
504          */
505         ListSetupType listType =
506                 (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
507         VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
508
509         const real rlist_new =
510                 calcVerletBufferSize(mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
511
512         if (rlist_new != ir->rlist)
513         {
514             if (fplog != nullptr)
515             {
516                 fprintf(fplog,
517                         "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
518                         ir->rlist,
519                         rlist_new,
520                         listSetup.cluster_size_i,
521                         listSetup.cluster_size_j);
522             }
523             ir->rlist = rlist_new;
524         }
525     }
526
527     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
528     {
529         gmx_fatal(FARGS,
530                   "Can not set nstlist without %s",
531                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
532     }
533
534     if (EI_DYNAMICS(ir->eI))
535     {
536         /* Set or try nstlist values */
537         increaseNstlist(fplog, cr, ir, nstlist_cmdline, &mtop, box, makeGpuPairList, cpuinfo);
538     }
539 }
540
541 /*! \brief Override the nslist value in inputrec
542  *
543  * with value passed on the command line (if any)
544  */
545 static void override_nsteps_cmdline(const gmx::MDLogger& mdlog, int64_t nsteps_cmdline, t_inputrec* ir)
546 {
547     assert(ir);
548
549     /* override with anything else than the default -2 */
550     if (nsteps_cmdline > -2)
551     {
552         char sbuf_steps[STEPSTRSIZE];
553         char sbuf_msg[STRLEN];
554
555         ir->nsteps = nsteps_cmdline;
556         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
557         {
558             sprintf(sbuf_msg,
559                     "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
560                     gmx_step_str(nsteps_cmdline, sbuf_steps),
561                     fabs(nsteps_cmdline * ir->delta_t));
562         }
563         else
564         {
565             sprintf(sbuf_msg,
566                     "Overriding nsteps with value passed on the command line: %s steps",
567                     gmx_step_str(nsteps_cmdline, sbuf_steps));
568         }
569
570         GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
571     }
572     else if (nsteps_cmdline < -2)
573     {
574         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline);
575     }
576     /* Do nothing if nsteps_cmdline == -2 */
577 }
578
579 namespace gmx
580 {
581
582 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
583  *
584  * If not, and if a warning may be issued, logs a warning about
585  * falling back to CPU code. With thread-MPI, only the first
586  * call to this function should have \c issueWarning true. */
587 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger& mdlog, const t_inputrec& ir, bool issueWarning)
588 {
589     bool        gpuIsUseful = true;
590     std::string warning;
591
592     if (ir.opts.ngener - ir.nwall > 1)
593     {
594         /* The GPU code does not support more than one energy group.
595          * If the user requested GPUs explicitly, a fatal error is given later.
596          */
597         gpuIsUseful = false;
598         warning =
599                 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
600                 "For better performance, run on the GPU without energy groups and then do "
601                 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
602     }
603
604     if (EI_TPI(ir.eI))
605     {
606         gpuIsUseful = false;
607         warning     = "TPI is not implemented for GPUs.";
608     }
609
610     if (!gpuIsUseful && issueWarning)
611     {
612         GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
613     }
614
615     return gpuIsUseful;
616 }
617
618 //! Initializes the logger for mdrun.
619 static gmx::LoggerOwner buildLogger(FILE* fplog, const bool isSimulationMasterRank)
620 {
621     gmx::LoggerBuilder builder;
622     if (fplog != nullptr)
623     {
624         builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
625     }
626     if (isSimulationMasterRank)
627     {
628         builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
629     }
630     return builder.build();
631 }
632
633 //! Make a TaskTarget from an mdrun argument string.
634 static TaskTarget findTaskTarget(const char* optionString)
635 {
636     TaskTarget returnValue = TaskTarget::Auto;
637
638     if (strncmp(optionString, "auto", 3) == 0)
639     {
640         returnValue = TaskTarget::Auto;
641     }
642     else if (strncmp(optionString, "cpu", 3) == 0)
643     {
644         returnValue = TaskTarget::Cpu;
645     }
646     else if (strncmp(optionString, "gpu", 3) == 0)
647     {
648         returnValue = TaskTarget::Gpu;
649     }
650     else
651     {
652         GMX_ASSERT(false, "Option string should have been checked for sanity already");
653     }
654
655     return returnValue;
656 }
657
658 //! Finish run, aggregate data to print performance info.
659 static void finish_run(FILE*                     fplog,
660                        const gmx::MDLogger&      mdlog,
661                        const t_commrec*          cr,
662                        const t_inputrec&         inputrec,
663                        t_nrnb                    nrnb[],
664                        gmx_wallcycle*            wcycle,
665                        gmx_walltime_accounting_t walltime_accounting,
666                        nonbonded_verlet_t*       nbv,
667                        const gmx_pme_t*          pme,
668                        gmx_bool                  bWriteStat)
669 {
670     double delta_t = 0;
671     double nbfs = 0, mflop = 0;
672     double elapsed_time, elapsed_time_over_all_ranks, elapsed_time_over_all_threads,
673             elapsed_time_over_all_threads_over_all_ranks;
674     /* Control whether it is valid to print a report. Only the
675        simulation master may print, but it should not do so if the run
676        terminated e.g. before a scheduled reset step. This is
677        complicated by the fact that PME ranks are unaware of the
678        reason why they were sent a pmerecvqxFINISH. To avoid
679        communication deadlocks, we always do the communication for the
680        report, even if we've decided not to write the report, because
681        how long it takes to finish the run is not important when we've
682        decided not to report on the simulation performance.
683
684        Further, we only report performance for dynamical integrators,
685        because those are the only ones for which we plan to
686        consider doing any optimizations. */
687     bool printReport = EI_DYNAMICS(inputrec.eI) && SIMMASTER(cr);
688
689     if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
690     {
691         GMX_LOG(mdlog.warning)
692                 .asParagraph()
693                 .appendText("Simulation ended prematurely, no performance report will be written.");
694         printReport = false;
695     }
696
697     t_nrnb*                 nrnb_tot;
698     std::unique_ptr<t_nrnb> nrnbTotalStorage;
699     if (cr->nnodes > 1)
700     {
701         nrnbTotalStorage = std::make_unique<t_nrnb>();
702         nrnb_tot         = nrnbTotalStorage.get();
703 #if GMX_MPI
704         MPI_Allreduce(nrnb->n.data(), nrnb_tot->n.data(), eNRNB, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
705 #endif
706     }
707     else
708     {
709         nrnb_tot = nrnb;
710     }
711
712     elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
713     elapsed_time_over_all_threads =
714             walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
715     if (cr->nnodes > 1)
716     {
717 #if GMX_MPI
718         /* reduce elapsed_time over all MPI ranks in the current simulation */
719         MPI_Allreduce(&elapsed_time, &elapsed_time_over_all_ranks, 1, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mysim);
720         elapsed_time_over_all_ranks /= cr->nnodes;
721         /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
722          * current simulation. */
723         MPI_Allreduce(&elapsed_time_over_all_threads,
724                       &elapsed_time_over_all_threads_over_all_ranks,
725                       1,
726                       MPI_DOUBLE,
727                       MPI_SUM,
728                       cr->mpi_comm_mysim);
729 #endif
730     }
731     else
732     {
733         elapsed_time_over_all_ranks                  = elapsed_time;
734         elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
735     }
736
737     if (printReport)
738     {
739         print_flop(fplog, nrnb_tot, &nbfs, &mflop);
740     }
741
742     if (thisRankHasDuty(cr, DUTY_PP) && haveDDAtomOrdering(*cr))
743     {
744         print_dd_statistics(cr, inputrec, fplog);
745     }
746
747     /* TODO Move the responsibility for any scaling by thread counts
748      * to the code that handled the thread region, so that there's a
749      * mechanism to keep cycle counting working during the transition
750      * to task parallelism. */
751     int nthreads_pp  = gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded);
752     int nthreads_pme = gmx_omp_nthreads_get(ModuleMultiThread::Pme);
753     wallcycle_scale_by_num_threads(
754             wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
755     auto cycle_sum(wallcycle_sum(cr, wcycle));
756
757     if (printReport)
758     {
759         auto* nbnxn_gpu_timings =
760                 (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
761         gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
762
763         if (pme_gpu_task_enabled(pme))
764         {
765             pme_gpu_get_timings(pme, &pme_gpu_timings);
766         }
767         wallcycle_print(fplog,
768                         mdlog,
769                         cr->nnodes,
770                         cr->npmenodes,
771                         nthreads_pp,
772                         nthreads_pme,
773                         elapsed_time_over_all_ranks,
774                         wcycle,
775                         cycle_sum,
776                         nbnxn_gpu_timings,
777                         &pme_gpu_timings);
778
779         if (EI_DYNAMICS(inputrec.eI))
780         {
781             delta_t = inputrec.delta_t;
782         }
783
784         if (fplog)
785         {
786             print_perf(fplog,
787                        elapsed_time_over_all_threads_over_all_ranks,
788                        elapsed_time_over_all_ranks,
789                        walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
790                        delta_t,
791                        nbfs,
792                        mflop);
793         }
794         if (bWriteStat)
795         {
796             print_perf(stderr,
797                        elapsed_time_over_all_threads_over_all_ranks,
798                        elapsed_time_over_all_ranks,
799                        walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
800                        delta_t,
801                        nbfs,
802                        mflop);
803         }
804     }
805 }
806
807 int Mdrunner::mdrunner()
808 {
809     matrix                      box;
810     std::unique_ptr<t_forcerec> fr;
811     real                        ewaldcoeff_q     = 0;
812     real                        ewaldcoeff_lj    = 0;
813     int                         nChargePerturbed = -1, nTypePerturbed = 0;
814     gmx_walltime_accounting_t   walltime_accounting = nullptr;
815     MembedHolder                membedHolder(filenames.size(), filenames.data());
816
817     /* CAUTION: threads may be started later on in this function, so
818        cr doesn't reflect the final parallel state right now */
819     gmx_mtop_t mtop;
820
821     /* TODO: inputrec should tell us whether we use an algorithm, not a file option */
822     const bool doEssentialDynamics = opt2bSet("-ei", filenames.size(), filenames.data());
823     const bool doRerun             = mdrunOptions.rerun;
824
825     // Handle task-assignment related user options.
826     EmulateGpuNonbonded emulateGpuNonbonded =
827             (getenv("GMX_EMULATE_GPU") != nullptr ? EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
828
829     std::vector<int> userGpuTaskAssignment;
830     try
831     {
832         userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
833     }
834     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
835     auto nonbondedTarget = findTaskTarget(nbpu_opt);
836     auto pmeTarget       = findTaskTarget(pme_opt);
837     auto pmeFftTarget    = findTaskTarget(pme_fft_opt);
838     auto bondedTarget    = findTaskTarget(bonded_opt);
839     auto updateTarget    = findTaskTarget(update_opt);
840
841     FILE* fplog = nullptr;
842     // If we are appending, we don't write log output because we need
843     // to check that the old log file matches what the checkpoint file
844     // expects. Otherwise, we should start to write log output now if
845     // there is a file ready for it.
846     if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
847     {
848         fplog = gmx_fio_getfp(logFileHandle);
849     }
850     const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, simulationCommunicator);
851     gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
852     gmx::MDLogger    mdlog(logOwner.logger());
853
854     gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo_);
855
856     std::vector<int> availableDevices =
857             makeListOfAvailableDevices(hwinfo_->deviceInfoList, hw_opt.devicesSelectedByUser);
858     const int numAvailableDevices = gmx::ssize(availableDevices);
859
860     // Print citation requests after all software/hardware printing
861     pleaseCiteGromacs(fplog);
862
863     // Note: legacy program logic relies on checking whether these pointers are assigned.
864     // Objects may or may not be allocated later.
865     std::unique_ptr<t_inputrec> inputrec;
866     std::unique_ptr<t_state>    globalState;
867
868     auto partialDeserializedTpr = std::make_unique<PartialDeserializedTprFile>();
869
870     if (isSimulationMasterRank)
871     {
872         // Allocate objects to be initialized by later function calls.
873         /* Only the master rank has the global state */
874         globalState = std::make_unique<t_state>();
875         inputrec    = std::make_unique<t_inputrec>();
876
877         /* Read (nearly) all data required for the simulation
878          * and keep the partly serialized tpr contents to send to other ranks later
879          */
880         applyGlobalSimulationState(
881                 *inputHolder_.get(), partialDeserializedTpr.get(), globalState.get(), inputrec.get(), &mtop);
882     }
883
884     /* Check and update the hardware options for internal consistency */
885     checkAndUpdateHardwareOptions(
886             mdlog, &hw_opt, isSimulationMasterRank, domdecOptions.numPmeRanks, inputrec.get());
887
888     if (GMX_THREAD_MPI && isSimulationMasterRank)
889     {
890         bool useGpuForNonbonded = false;
891         bool useGpuForPme       = false;
892         try
893         {
894             GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
895
896             // If the user specified the number of ranks, then we must
897             // respect that, but in default mode, we need to allow for
898             // the number of GPUs to choose the number of ranks.
899             auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
900             useGpuForNonbonded         = decideWhetherToUseGpusForNonbondedWithThreadMpi(
901                     nonbondedTarget,
902                     numAvailableDevices > 0,
903                     userGpuTaskAssignment,
904                     emulateGpuNonbonded,
905                     canUseGpuForNonbonded,
906                     gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
907                     hw_opt.nthreads_tmpi);
908             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(useGpuForNonbonded,
909                                                                      pmeTarget,
910                                                                      pmeFftTarget,
911                                                                      numAvailableDevices,
912                                                                      userGpuTaskAssignment,
913                                                                      *hwinfo_,
914                                                                      *inputrec,
915                                                                      hw_opt.nthreads_tmpi,
916                                                                      domdecOptions.numPmeRanks);
917         }
918         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
919
920         /* Determine how many thread-MPI ranks to start.
921          *
922          * TODO Over-writing the user-supplied value here does
923          * prevent any possible subsequent checks from working
924          * correctly. */
925         hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo_,
926                                                 &hw_opt,
927                                                 numAvailableDevices,
928                                                 useGpuForNonbonded,
929                                                 useGpuForPme,
930                                                 inputrec.get(),
931                                                 mtop,
932                                                 mdlog,
933                                                 membedHolder.doMembed());
934
935         // Now start the threads for thread MPI.
936         spawnThreads(hw_opt.nthreads_tmpi);
937         // The spawned threads enter mdrunner() and execution of
938         // master and spawned threads joins at the end of this block.
939     }
940
941     GMX_RELEASE_ASSERT(!GMX_MPI || ms || simulationCommunicator != MPI_COMM_NULL,
942                        "Must have valid communicator unless running a multi-simulation");
943     CommrecHandle crHandle = init_commrec(simulationCommunicator);
944     t_commrec*    cr       = crHandle.get();
945     GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
946
947     PhysicalNodeCommunicator physicalNodeComm(libraryWorldCommunicator, gmx_physicalnode_id_hash());
948
949     // If we detected the topology on this system, double-check that it makes sense
950     if (hwinfo_->hardwareTopology->isThisSystem())
951     {
952         hardwareTopologyDoubleCheckDetection(mdlog, *hwinfo_->hardwareTopology);
953     }
954
955     if (PAR(cr))
956     {
957         /* now broadcast everything to the non-master nodes/threads: */
958         if (!isSimulationMasterRank)
959         {
960             // Until now, only the master rank has a non-null pointer.
961             // On non-master ranks, allocate the object that will receive data in the following call.
962             inputrec = std::make_unique<t_inputrec>();
963         }
964         init_parallel(cr->mpiDefaultCommunicator,
965                       MASTER(cr),
966                       inputrec.get(),
967                       &mtop,
968                       partialDeserializedTpr.get());
969     }
970     GMX_RELEASE_ASSERT(inputrec != nullptr, "All ranks should have a valid inputrec now");
971     partialDeserializedTpr.reset(nullptr);
972
973     // Now the number of ranks is known to all ranks, and each knows
974     // the inputrec read by the master rank. The ranks can now all run
975     // the task-deciding functions and will agree on the result
976     // without needing to communicate.
977     // The LBFGS minimizer, test-particle insertion, normal modes and shell dynamics don't support DD
978     const bool useDomainDecomposition =
979             !(inputrec->eI == IntegrationAlgorithm::LBFGS || EI_TPI(inputrec->eI)
980               || inputrec->eI == IntegrationAlgorithm::NM
981               || gmx_mtop_particletype_count(mtop)[ParticleType::Shell] > 0);
982
983     // Note that these variables describe only their own node.
984     //
985     // Note that when bonded interactions run on a GPU they always run
986     // alongside a nonbonded task, so do not influence task assignment
987     // even though they affect the force calculation workload.
988     bool useGpuForNonbonded = false;
989     bool useGpuForPme       = false;
990     bool useGpuForBonded    = false;
991     bool useGpuForUpdate    = false;
992     bool gpusWereDetected   = hwinfo_->ngpu_compatible_tot > 0;
993     try
994     {
995         // It's possible that there are different numbers of GPUs on
996         // different nodes, which is the user's responsibility to
997         // handle. If unsuitable, we will notice that during task
998         // assignment.
999         auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
1000         useGpuForNonbonded         = decideWhetherToUseGpusForNonbonded(
1001                 nonbondedTarget,
1002                 userGpuTaskAssignment,
1003                 emulateGpuNonbonded,
1004                 canUseGpuForNonbonded,
1005                 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI),
1006                 gpusWereDetected);
1007         useGpuForPme    = decideWhetherToUseGpusForPme(useGpuForNonbonded,
1008                                                     pmeTarget,
1009                                                     pmeFftTarget,
1010                                                     userGpuTaskAssignment,
1011                                                     *hwinfo_,
1012                                                     *inputrec,
1013                                                     cr->sizeOfDefaultCommunicator,
1014                                                     domdecOptions.numPmeRanks,
1015                                                     gpusWereDetected);
1016         useGpuForBonded = decideWhetherToUseGpusForBonded(
1017                 useGpuForNonbonded, useGpuForPme, bondedTarget, *inputrec, mtop, domdecOptions.numPmeRanks, gpusWereDetected);
1018     }
1019     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1020
1021     const PmeRunMode pmeRunMode = determinePmeRunMode(useGpuForPme, pmeFftTarget, *inputrec);
1022
1023     // Initialize development feature flags that enabled by environment variable
1024     // and report those features that are enabled.
1025     const DevelopmentFeatureFlags devFlags =
1026             manageDevelopmentFeatures(mdlog, useGpuForNonbonded, pmeRunMode);
1027
1028     const bool useModularSimulator = checkUseModularSimulator(false,
1029                                                               inputrec.get(),
1030                                                               doRerun,
1031                                                               mtop,
1032                                                               ms,
1033                                                               replExParams,
1034                                                               nullptr,
1035                                                               doEssentialDynamics,
1036                                                               membedHolder.doMembed());
1037
1038     ObservablesReducerBuilder observablesReducerBuilder;
1039
1040     // Build restraints.
1041     // TODO: hide restraint implementation details from Mdrunner.
1042     // There is nothing unique about restraints at this point as far as the
1043     // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
1044     // factory functions from the SimulationContext on which to call mdModules_->add().
1045     // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
1046     for (auto&& restraint : restraintManager_->getRestraints())
1047     {
1048         auto module = RestraintMDModule::create(restraint, restraint->sites());
1049         mdModules_->add(std::move(module));
1050     }
1051
1052     // TODO: Error handling
1053     mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
1054     // now that the MDModules know their options, they know which callbacks to sign up to
1055     mdModules_->subscribeToSimulationSetupNotifications();
1056     const auto& setupNotifier = mdModules_->notifiers().simulationSetupNotifier_;
1057
1058     // Notify MdModules of existing logger
1059     setupNotifier.notify(mdlog);
1060
1061     // Notify MdModules of internal parameters, saved into KVT
1062     if (inputrec->internalParameters != nullptr)
1063     {
1064         setupNotifier.notify(*inputrec->internalParameters);
1065     }
1066
1067     // Let MdModules know the .tpr filename
1068     {
1069         gmx::MdRunInputFilename mdRunInputFilename = { ftp2fn(efTPR, filenames.size(), filenames.data()) };
1070         setupNotifier.notify(mdRunInputFilename);
1071     }
1072
1073     if (fplog != nullptr)
1074     {
1075         pr_inputrec(fplog, 0, "Input Parameters", inputrec.get(), FALSE);
1076         fprintf(fplog, "\n");
1077     }
1078
1079     if (SIMMASTER(cr))
1080     {
1081         /* In rerun, set velocities to zero if present */
1082         if (doRerun && ((globalState->flags & enumValueToBitMask(StateEntry::V)) != 0))
1083         {
1084             // rerun does not use velocities
1085             GMX_LOG(mdlog.info)
1086                     .asParagraph()
1087                     .appendText(
1088                             "Rerun trajectory contains velocities. Rerun does only evaluate "
1089                             "potential energy and forces. The velocities will be ignored.");
1090             for (int i = 0; i < globalState->natoms; i++)
1091             {
1092                 clear_rvec(globalState->v[i]);
1093             }
1094             globalState->flags &= ~enumValueToBitMask(StateEntry::V);
1095         }
1096
1097         /* now make sure the state is initialized and propagated */
1098         set_state_entries(globalState.get(), inputrec.get(), useModularSimulator);
1099     }
1100
1101     /* NM and TPI parallelize over force/energy calculations, not atoms,
1102      * so we need to initialize and broadcast the global state.
1103      */
1104     if (inputrec->eI == IntegrationAlgorithm::NM || inputrec->eI == IntegrationAlgorithm::TPI)
1105     {
1106         if (!MASTER(cr))
1107         {
1108             globalState = std::make_unique<t_state>();
1109         }
1110         broadcastStateWithoutDynamics(
1111                 cr->mpiDefaultCommunicator, haveDDAtomOrdering(*cr), PAR(cr), globalState.get());
1112     }
1113
1114     /* A parallel command line option consistency check that we can
1115        only do after any threads have started. */
1116     if (!PAR(cr)
1117         && (domdecOptions.numCells[XX] > 1 || domdecOptions.numCells[YY] > 1
1118             || domdecOptions.numCells[ZZ] > 1 || domdecOptions.numPmeRanks > 0))
1119     {
1120         gmx_fatal(FARGS,
1121                   "The -dd or -npme option request a parallel simulation, "
1122 #if !GMX_MPI
1123                   "but %s was compiled without threads or MPI enabled",
1124                   output_env_get_program_display_name(oenv));
1125 #elif GMX_THREAD_MPI
1126                   "but the number of MPI-threads (option -ntmpi) is not set or is 1");
1127 #else
1128                   "but %s was not started through mpirun/mpiexec or only one rank was requested "
1129                   "through mpirun/mpiexec",
1130                   output_env_get_program_display_name(oenv));
1131 #endif
1132     }
1133
1134     if (doRerun && (EI_ENERGY_MINIMIZATION(inputrec->eI) || IntegrationAlgorithm::NM == inputrec->eI))
1135     {
1136         gmx_fatal(FARGS,
1137                   "The .mdp file specified an energy mininization or normal mode algorithm, and "
1138                   "these are not compatible with mdrun -rerun");
1139     }
1140
1141     /* NMR restraints must be initialized before load_checkpoint,
1142      * since with time averaging the history is added to t_state.
1143      * For proper consistency check we therefore need to extend
1144      * t_state here.
1145      * So the PME-only nodes (if present) will also initialize
1146      * the distance restraints.
1147      */
1148
1149     /* This needs to be called before read_checkpoint to extend the state */
1150     t_disresdata* disresdata;
1151     snew(disresdata, 1);
1152     init_disres(fplog,
1153                 mtop,
1154                 inputrec.get(),
1155                 DisResRunMode::MDRun,
1156                 MASTER(cr) ? DDRole::Master : DDRole::Agent,
1157                 PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1158                 cr->mpi_comm_mysim,
1159                 ms,
1160                 disresdata,
1161                 globalState.get(),
1162                 replExParams.exchangeInterval > 0);
1163
1164     if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0 && isSimulationMasterRank)
1165     {
1166         extendStateWithOriresHistory(mtop, *inputrec, globalState.get());
1167     }
1168
1169     auto deform = prepareBoxDeformation(globalState != nullptr ? globalState->box : box,
1170                                         MASTER(cr) ? DDRole::Master : DDRole::Agent,
1171                                         PAR(cr) ? NumRanks::Multiple : NumRanks::Single,
1172                                         cr->mpi_comm_mygroup,
1173                                         *inputrec);
1174
1175 #if GMX_FAHCORE
1176     /* We have to remember the generation's first step before reading checkpoint.
1177        This way, we can report to the F@H core both the generation's first step
1178        and the restored first step, thus making it able to distinguish between
1179        an interruption/resume and start of the n-th generation simulation.
1180        Having this information, the F@H core can correctly calculate and report
1181        the progress.
1182      */
1183     int gen_first_step = 0;
1184     if (MASTER(cr))
1185     {
1186         gen_first_step = inputrec->init_step;
1187     }
1188 #endif
1189
1190     ObservablesHistory observablesHistory = {};
1191
1192     auto modularSimulatorCheckpointData = std::make_unique<ReadCheckpointDataHolder>();
1193     if (startingBehavior != StartingBehavior::NewSimulation)
1194     {
1195         /* Check if checkpoint file exists before doing continuation.
1196          * This way we can use identical input options for the first and subsequent runs...
1197          */
1198         if (mdrunOptions.numStepsCommandline > -2)
1199         {
1200             /* Temporarily set the number of steps to unlimited to avoid
1201              * triggering the nsteps check in load_checkpoint().
1202              * This hack will go away soon when the -nsteps option is removed.
1203              */
1204             inputrec->nsteps = -1;
1205         }
1206
1207         // Finish applying initial simulation state information from external sources on all ranks.
1208         // Reconcile checkpoint file data with Mdrunner state established up to this point.
1209         applyLocalState(*inputHolder_.get(),
1210                         logFileHandle,
1211                         cr,
1212                         domdecOptions.numCells,
1213                         inputrec.get(),
1214                         globalState.get(),
1215                         &observablesHistory,
1216                         mdrunOptions.reproducible,
1217                         mdModules_->notifiers(),
1218                         modularSimulatorCheckpointData.get(),
1219                         useModularSimulator);
1220         // TODO: (#3652) Synchronize filesystem state, SimulationInput contents, and program
1221         //  invariants
1222         //  on all code paths.
1223         // Write checkpoint or provide hook to update SimulationInput.
1224         // If there was a checkpoint file, SimulationInput contains more information
1225         // than if there wasn't. At this point, we have synchronized the in-memory
1226         // state with the filesystem state only for restarted simulations. We should
1227         // be calling applyLocalState unconditionally and expect that the completeness
1228         // of SimulationInput is not dependent on its creation method.
1229
1230         if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
1231         {
1232             // Now we can start normal logging to the truncated log file.
1233             fplog = gmx_fio_getfp(logFileHandle);
1234             prepareLogAppending(fplog);
1235             logOwner = buildLogger(fplog, MASTER(cr));
1236             mdlog    = logOwner.logger();
1237         }
1238     }
1239
1240 #if GMX_FAHCORE
1241     if (MASTER(cr))
1242     {
1243         fcRegisterSteps(inputrec->nsteps + inputrec->init_step, gen_first_step);
1244     }
1245 #endif
1246
1247     if (mdrunOptions.numStepsCommandline > -2)
1248     {
1249         GMX_LOG(mdlog.info)
1250                 .asParagraph()
1251                 .appendText(
1252                         "The -nsteps functionality is deprecated, and may be removed in a future "
1253                         "version. "
1254                         "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp "
1255                         "file field.");
1256     }
1257     /* override nsteps with value set on the commandline */
1258     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec.get());
1259
1260     if (isSimulationMasterRank)
1261     {
1262         copy_mat(globalState->box, box);
1263     }
1264
1265     if (PAR(cr))
1266     {
1267         gmx_bcast(sizeof(box), box, cr->mpiDefaultCommunicator);
1268     }
1269
1270     if (inputrec->cutoff_scheme != CutoffScheme::Verlet)
1271     {
1272         gmx_fatal(FARGS,
1273                   "This group-scheme .tpr file can no longer be run by mdrun. Please update to the "
1274                   "Verlet scheme, or use an earlier version of GROMACS if necessary.");
1275     }
1276     /* Update rlist and nstlist. */
1277     /* Note: prepare_verlet_scheme is calling increaseNstlist(...), which (while attempting to
1278      * increase rlist) tries to check if the newly chosen value fits with the DD scheme. As this is
1279      * run before any DD scheme is set up, this check is never executed. See #3334 for more details.
1280      */
1281     prepare_verlet_scheme(fplog,
1282                           cr,
1283                           inputrec.get(),
1284                           nstlist_cmdline,
1285                           mtop,
1286                           box,
1287                           useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes),
1288                           *hwinfo_->cpuInfo);
1289
1290     // We need to decide on update groups early, as this affects
1291     // inter-domain communication distances.
1292     auto       updateGroupingsPerMoleculeType = makeUpdateGroupingsPerMoleculeType(mtop);
1293     const real maxUpdateGroupRadius           = computeMaxUpdateGroupRadius(
1294             mtop, updateGroupingsPerMoleculeType, maxReferenceTemperature(*inputrec));
1295     const real   cutoffMargin = std::sqrt(max_cutoff2(inputrec->pbcType, box)) - inputrec->rlist;
1296     UpdateGroups updateGroups = makeUpdateGroups(mdlog,
1297                                                  std::move(updateGroupingsPerMoleculeType),
1298                                                  maxUpdateGroupRadius,
1299                                                  useDomainDecomposition,
1300                                                  systemHasConstraintsOrVsites(mtop),
1301                                                  cutoffMargin);
1302
1303     try
1304     {
1305         const bool haveFrozenAtoms = inputrecFrozenAtoms(inputrec.get());
1306
1307         useGpuForUpdate = decideWhetherToUseGpuForUpdate(useDomainDecomposition,
1308                                                          updateGroups.useUpdateGroups(),
1309                                                          pmeRunMode,
1310                                                          domdecOptions.numPmeRanks > 0,
1311                                                          useGpuForNonbonded,
1312                                                          updateTarget,
1313                                                          gpusWereDetected,
1314                                                          *inputrec,
1315                                                          mtop,
1316                                                          doEssentialDynamics,
1317                                                          gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
1318                                                          haveFrozenAtoms,
1319                                                          doRerun,
1320                                                          devFlags,
1321                                                          mdlog);
1322     }
1323     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1324
1325     bool useGpuDirectHalo = false;
1326
1327     if (useGpuForNonbonded)
1328     {
1329         // cr->npmenodes is not yet initialized.
1330         // domdecOptions.numPmeRanks == -1 results in 0 separate PME ranks when useGpuForNonbonded is true.
1331         // Todo: remove this assumption later once auto mode has support for separate PME rank
1332         const int numPmeRanks = domdecOptions.numPmeRanks > 0 ? domdecOptions.numPmeRanks : 0;
1333         bool      havePPDomainDecomposition = (cr->sizeOfDefaultCommunicator - numPmeRanks) > 1;
1334         useGpuDirectHalo                    = decideWhetherToUseGpuForHalo(devFlags,
1335                                                         havePPDomainDecomposition,
1336                                                         useGpuForNonbonded,
1337                                                         useModularSimulator,
1338                                                         doRerun,
1339                                                         EI_ENERGY_MINIMIZATION(inputrec->eI));
1340     }
1341
1342     // This builder is necessary while we have multi-part construction
1343     // of DD. Before DD is constructed, we use the existence of
1344     // the builder object to indicate that further construction of DD
1345     // is needed.
1346     std::unique_ptr<DomainDecompositionBuilder> ddBuilder;
1347     if (useDomainDecomposition)
1348     {
1349         // P2P GPU comm + GPU update leads to case in which we enqueue async work for multiple
1350         // timesteps. DLB needs to be disabled in that case
1351         const bool directGpuCommUsedWithGpuUpdate = GMX_THREAD_MPI && useGpuDirectHalo && useGpuForUpdate;
1352         ddBuilder                                 = std::make_unique<DomainDecompositionBuilder>(
1353                 mdlog,
1354                 cr,
1355                 domdecOptions,
1356                 mdrunOptions,
1357                 mtop,
1358                 *inputrec,
1359                 mdModules_->notifiers(),
1360                 box,
1361                 updateGroups.updateGroupingPerMoleculeType(),
1362                 updateGroups.useUpdateGroups(),
1363                 updateGroups.maxUpdateGroupRadius(),
1364                 positionsFromStatePointer(globalState.get()),
1365                 useGpuForNonbonded,
1366                 useGpuForPme,
1367                 directGpuCommUsedWithGpuUpdate);
1368     }
1369     else
1370     {
1371         /* PME, if used, is done on all nodes with 1D decomposition */
1372         cr->nnodes     = cr->sizeOfDefaultCommunicator;
1373         cr->sim_nodeid = cr->rankInDefaultCommunicator;
1374         cr->nodeid     = cr->rankInDefaultCommunicator;
1375         cr->npmenodes  = 0;
1376         cr->duty       = (DUTY_PP | DUTY_PME);
1377
1378         if (inputrec->pbcType == PbcType::Screw)
1379         {
1380             gmx_fatal(FARGS, "pbc=screw is only implemented with domain decomposition");
1381         }
1382     }
1383
1384     // Produce the task assignment for this rank - done after DD is constructed
1385     GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
1386             availableDevices,
1387             userGpuTaskAssignment,
1388             *hwinfo_,
1389             simulationCommunicator,
1390             physicalNodeComm,
1391             nonbondedTarget,
1392             pmeTarget,
1393             bondedTarget,
1394             updateTarget,
1395             useGpuForNonbonded,
1396             useGpuForPme,
1397             thisRankHasDuty(cr, DUTY_PP),
1398             // TODO cr->duty & DUTY_PME should imply that a PME
1399             // algorithm is active, but currently does not.
1400             EEL_PME(inputrec->coulombtype) && thisRankHasDuty(cr, DUTY_PME));
1401
1402     // Get the device handles for the modules, nullptr when no task is assigned.
1403     int                deviceId   = -1;
1404     DeviceInformation* deviceInfo = gpuTaskAssignments.initDevice(&deviceId);
1405
1406     // TODO Currently this is always built, yet DD partition code
1407     // checks if it is built before using it. Probably it should
1408     // become an MDModule that is made only when another module
1409     // requires it (e.g. pull, CompEl, density fitting), so that we
1410     // don't update the local atom sets unilaterally every step.
1411     LocalAtomSetManager atomSets;
1412
1413     // Local state and topology are declared (and perhaps constructed)
1414     // now, because DD needs them for the LocalTopologyChecker, but
1415     // they do not contain valid data until after the first DD
1416     // partition.
1417     std::unique_ptr<t_state> localStateInstance;
1418     t_state*                 localState;
1419     gmx_localtop_t           localTopology(mtop.ffparams);
1420
1421     if (ddBuilder)
1422     {
1423         localStateInstance = std::make_unique<t_state>();
1424         localState         = localStateInstance.get();
1425         // TODO Pass the GPU streams to ddBuilder to use in buffer
1426         // transfers (e.g. halo exchange)
1427         cr->dd = ddBuilder->build(&atomSets, localTopology, *localState, &observablesReducerBuilder);
1428         // The builder's job is done, so destruct it
1429         ddBuilder.reset(nullptr);
1430         // Note that local state still does not exist yet.
1431     }
1432     else
1433     {
1434         // Without DD, the local state is merely an alias to the global state,
1435         // so we don't need to allocate anything.
1436         localState = globalState.get();
1437     }
1438
1439     // Ensure that all atoms within the same update group are in the
1440     // same periodic image. Otherwise, a simulation that did not use
1441     // update groups (e.g. a single-rank simulation) cannot always be
1442     // correctly restarted in a way that does use update groups
1443     // (e.g. a multi-rank simulation).
1444     if (isSimulationMasterRank)
1445     {
1446         const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
1447         if (useUpdateGroups)
1448         {
1449             putUpdateGroupAtomsInSamePeriodicImage(*cr->dd, mtop, globalState->box, globalState->x);
1450         }
1451     }
1452
1453     const bool printHostName = (cr->nnodes > 1);
1454     gpuTaskAssignments.reportGpuUsage(mdlog, printHostName, useGpuForBonded, pmeRunMode, useGpuForUpdate);
1455
1456     const bool disableNonbondedCalculation = (getenv("GMX_NO_NONBONDED") != nullptr);
1457     if (disableNonbondedCalculation)
1458     {
1459         /* turn off non-bonded calculations */
1460         GMX_LOG(mdlog.warning)
1461                 .asParagraph()
1462                 .appendText(
1463                         "Found environment variable GMX_NO_NONBONDED.\n"
1464                         "Disabling nonbonded calculations.");
1465     }
1466
1467     MdrunScheduleWorkload runScheduleWork;
1468
1469     // Also populates the simulation constant workload description.
1470     // Note: currently the default duty is DUTY_PP | DUTY_PME for all simulations, including those without PME,
1471     // so this boolean is sufficient on all ranks to determine whether separate PME ranks are used,
1472     // but this will no longer be the case if cr->duty is changed for !EEL_PME(fr->ic->eeltype).
1473     const bool haveSeparatePmeRank = (!thisRankHasDuty(cr, DUTY_PP) || !thisRankHasDuty(cr, DUTY_PME));
1474     runScheduleWork.simulationWork = createSimulationWorkload(*inputrec,
1475                                                               disableNonbondedCalculation,
1476                                                               devFlags,
1477                                                               havePPDomainDecomposition(cr),
1478                                                               haveSeparatePmeRank,
1479                                                               useGpuForNonbonded,
1480                                                               pmeRunMode,
1481                                                               useGpuForBonded,
1482                                                               useGpuForUpdate,
1483                                                               useGpuDirectHalo);
1484
1485     std::unique_ptr<DeviceStreamManager> deviceStreamManager = nullptr;
1486
1487     if (deviceInfo != nullptr)
1488     {
1489         if (runScheduleWork.simulationWork.havePpDomainDecomposition && thisRankHasDuty(cr, DUTY_PP))
1490         {
1491             dd_setup_dlb_resource_sharing(cr, deviceId);
1492         }
1493         const bool useGpuTiming = decideGpuTimingsUsage();
1494         deviceStreamManager     = std::make_unique<DeviceStreamManager>(
1495                 *deviceInfo, havePPDomainDecomposition(cr), runScheduleWork.simulationWork, useGpuTiming);
1496     }
1497
1498     // If the user chose a task assignment, give them some hints
1499     // where appropriate.
1500     if (!userGpuTaskAssignment.empty())
1501     {
1502         gpuTaskAssignments.logPerformanceHints(mdlog, numAvailableDevices);
1503     }
1504
1505     if (PAR(cr))
1506     {
1507         /* After possible communicator splitting in make_dd_communicators.
1508          * we can set up the intra/inter node communication.
1509          */
1510         gmx_setup_nodecomm(fplog, cr);
1511     }
1512
1513 #if GMX_MPI
1514     if (isMultiSim(ms))
1515     {
1516         GMX_LOG(mdlog.warning)
1517                 .asParagraph()
1518                 .appendTextFormatted(
1519                         "This is simulation %d out of %d running as a composite GROMACS\n"
1520                         "multi-simulation job. Setup for this simulation:\n",
1521                         ms->simulationIndex_,
1522                         ms->numSimulations_);
1523     }
1524     GMX_LOG(mdlog.warning)
1525             .appendTextFormatted("Using %d MPI %s\n",
1526                                  cr->nnodes,
1527 #    if GMX_THREAD_MPI
1528                                  cr->nnodes == 1 ? "thread" : "threads"
1529 #    else
1530                                  cr->nnodes == 1 ? "process" : "processes"
1531 #    endif
1532             );
1533     fflush(stderr);
1534 #endif
1535
1536     // If mdrun -pin auto honors any affinity setting that already
1537     // exists. If so, it is nice to provide feedback about whether
1538     // that existing affinity setting was from OpenMP or something
1539     // else, so we run this code both before and after we initialize
1540     // the OpenMP support.
1541     gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, FALSE);
1542     /* Check and update the number of OpenMP threads requested */
1543     checkAndUpdateRequestedNumOpenmpThreads(
1544             &hw_opt, *hwinfo_, cr, ms, physicalNodeComm.size_, pmeRunMode, mtop, *inputrec);
1545
1546     gmx_omp_nthreads_init(mdlog,
1547                           cr,
1548                           hwinfo_->nthreads_hw_avail,
1549                           physicalNodeComm.size_,
1550                           hw_opt.nthreads_omp,
1551                           hw_opt.nthreads_omp_pme,
1552                           !thisRankHasDuty(cr, DUTY_PP));
1553
1554     const bool bEnableFPE = gmxShouldEnableFPExceptions();
1555     // FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1556     if (bEnableFPE)
1557     {
1558         gmx_feenableexcept();
1559     }
1560
1561     /* Now that we know the setup is consistent, check for efficiency */
1562     check_resource_division_efficiency(
1563             hwinfo_, gpuTaskAssignments.thisRankHasAnyGpuTask(), mdrunOptions.ntompOptionIsSet, cr, mdlog);
1564
1565     /* getting number of PP/PME threads on this MPI / tMPI rank.
1566        PME: env variable should be read only on one node to make sure it is
1567        identical everywhere;
1568      */
1569     const int numThreadsOnThisRank = thisRankHasDuty(cr, DUTY_PP)
1570                                              ? gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded)
1571                                              : gmx_omp_nthreads_get(ModuleMultiThread::Pme);
1572     checkHardwareOversubscription(
1573             numThreadsOnThisRank, cr->nodeid, *hwinfo_->hardwareTopology, physicalNodeComm, mdlog);
1574
1575     // Enable Peer access between GPUs where available
1576     // Only for DD, only master PP rank needs to perform setup, and only if thread MPI plus
1577     // any of the GPU communication features are active.
1578     if (haveDDAtomOrdering(*cr) && MASTER(cr) && thisRankHasDuty(cr, DUTY_PP) && GMX_THREAD_MPI
1579         && (runScheduleWork.simulationWork.useGpuHaloExchange
1580             || runScheduleWork.simulationWork.useGpuPmePpCommunication))
1581     {
1582         setupGpuDevicePeerAccess(gpuTaskAssignments.deviceIdsAssigned(), mdlog);
1583     }
1584
1585     if (hw_opt.threadAffinity != ThreadAffinity::Off)
1586     {
1587         /* Before setting affinity, check whether the affinity has changed
1588          * - which indicates that probably the OpenMP library has changed it
1589          * since we first checked).
1590          */
1591         gmx_check_thread_affinity_set(mdlog, &hw_opt, hwinfo_->nthreads_hw_avail, TRUE);
1592
1593         int numThreadsOnThisNode, intraNodeThreadOffset;
1594         analyzeThreadsOnThisNode(
1595                 physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode, &intraNodeThreadOffset);
1596
1597         /* Set the CPU affinity */
1598         gmx_set_thread_affinity(mdlog,
1599                                 cr,
1600                                 &hw_opt,
1601                                 *hwinfo_->hardwareTopology,
1602                                 numThreadsOnThisRank,
1603                                 numThreadsOnThisNode,
1604                                 intraNodeThreadOffset,
1605                                 nullptr);
1606     }
1607
1608     if (mdrunOptions.timingOptions.resetStep > -1)
1609     {
1610         GMX_LOG(mdlog.info)
1611                 .asParagraph()
1612                 .appendText(
1613                         "The -resetstep functionality is deprecated, and may be removed in a "
1614                         "future version.");
1615     }
1616     std::unique_ptr<gmx_wallcycle> wcycle =
1617             wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1618
1619     if (PAR(cr))
1620     {
1621         /* Master synchronizes its value of reset_counters with all nodes
1622          * including PME only nodes */
1623         int64_t reset_counters = wcycle_get_reset_counters(wcycle.get());
1624         gmx_bcast(sizeof(reset_counters), &reset_counters, cr->mpi_comm_mysim);
1625         wcycle_set_reset_counters(wcycle.get(), reset_counters);
1626     }
1627
1628     // Membrane embedding must be initialized before we call init_forcerec()
1629     membedHolder.initializeMembed(fplog,
1630                                   filenames.size(),
1631                                   filenames.data(),
1632                                   &mtop,
1633                                   inputrec.get(),
1634                                   globalState.get(),
1635                                   cr,
1636                                   &mdrunOptions.checkpointOptions.period);
1637
1638     const bool               thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
1639     std::unique_ptr<MDAtoms> mdAtoms;
1640     std::unique_ptr<VirtualSitesHandler> vsite;
1641
1642     t_nrnb nrnb;
1643     if (thisRankHasDuty(cr, DUTY_PP))
1644     {
1645         setupNotifier.notify(*cr);
1646         setupNotifier.notify(&atomSets);
1647         setupNotifier.notify(mtop);
1648         setupNotifier.notify(inputrec->pbcType);
1649         setupNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
1650         /* Initiate forcerecord */
1651         fr                 = std::make_unique<t_forcerec>();
1652         fr->forceProviders = mdModules_->initForceProviders();
1653         init_forcerec(fplog,
1654                       mdlog,
1655                       runScheduleWork.simulationWork,
1656                       fr.get(),
1657                       *inputrec,
1658                       mtop,
1659                       cr,
1660                       box,
1661                       opt2fn("-table", filenames.size(), filenames.data()),
1662                       opt2fn("-tablep", filenames.size(), filenames.data()),
1663                       opt2fns("-tableb", filenames.size(), filenames.data()),
1664                       pforce);
1665         // Dirty hack, for fixing disres and orires should be made mdmodules
1666         fr->fcdata->disres = disresdata;
1667         if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1668         {
1669             fr->fcdata->orires = std::make_unique<t_oriresdata>(
1670                     fplog, mtop, *inputrec, ms, globalState.get(), &atomSets);
1671         }
1672
1673         // Save a handle to device stream manager to use elsewhere in the code
1674         // TODO: Forcerec is not a correct place to store it.
1675         fr->deviceStreamManager = deviceStreamManager.get();
1676
1677         if (runScheduleWork.simulationWork.useGpuPmePpCommunication && !thisRankHasDuty(cr, DUTY_PME))
1678         {
1679             GMX_RELEASE_ASSERT(
1680                     deviceStreamManager != nullptr,
1681                     "GPU device stream manager should be valid in order to use PME-PP direct "
1682                     "communications.");
1683             GMX_RELEASE_ASSERT(
1684                     deviceStreamManager->streamIsValid(DeviceStreamType::PmePpTransfer),
1685                     "GPU PP-PME stream should be valid in order to use GPU PME-PP direct "
1686                     "communications.");
1687             fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(
1688                     cr->mpi_comm_mysim,
1689                     cr->dd->pme_nodeid,
1690                     &cr->dd->pmeForceReceiveBuffer,
1691                     deviceStreamManager->context(),
1692                     deviceStreamManager->stream(DeviceStreamType::PmePpTransfer));
1693         }
1694
1695         fr->nbv = Nbnxm::init_nb_verlet(mdlog,
1696                                         *inputrec,
1697                                         *fr,
1698                                         cr,
1699                                         *hwinfo_,
1700                                         runScheduleWork.simulationWork.useGpuNonbonded,
1701                                         deviceStreamManager.get(),
1702                                         mtop,
1703                                         box,
1704                                         wcycle.get());
1705         // TODO: Move the logic below to a GPU bonded builder
1706         if (runScheduleWork.simulationWork.useGpuBonded)
1707         {
1708             GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
1709                                "GPU device stream manager should be valid in order to use GPU "
1710                                "version of bonded forces.");
1711             fr->listedForcesGpu = std::make_unique<ListedForcesGpu>(
1712                     mtop.ffparams,
1713                     fr->ic->epsfac * fr->fudgeQQ,
1714                     deviceStreamManager->context(),
1715                     deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)),
1716                     wcycle.get());
1717         }
1718         fr->longRangeNonbondeds = std::make_unique<CpuPpLongRangeNonbondeds>(fr->n_tpi,
1719                                                                              fr->ic->ewaldcoeff_q,
1720                                                                              fr->ic->epsilon_r,
1721                                                                              fr->qsum,
1722                                                                              fr->ic->eeltype,
1723                                                                              fr->ic->vdwtype,
1724                                                                              *inputrec,
1725                                                                              &nrnb,
1726                                                                              wcycle.get(),
1727                                                                              fplog);
1728
1729         /* Initialize the mdAtoms structure.
1730          * mdAtoms is not filled with atom data,
1731          * as this can not be done now with domain decomposition.
1732          */
1733         mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1734         if (globalState && thisRankHasPmeGpuTask)
1735         {
1736             // The pinning of coordinates in the global state object works, because we only use
1737             // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1738             // points to the global state object without DD.
1739             // FIXME: MD and EM separately set up the local state - this should happen in the same
1740             // function, which should also perform the pinning.
1741             changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1742         }
1743
1744         /* Initialize the virtual site communication */
1745         vsite = makeVirtualSitesHandler(
1746                 mtop, cr, fr->pbcType, updateGroups.updateGroupingPerMoleculeType());
1747
1748         calc_shifts(box, fr->shift_vec);
1749
1750         /* With periodic molecules the charge groups should be whole at start up
1751          * and the virtual sites should not be far from their proper positions.
1752          */
1753         if (!inputrec->bContinuation && MASTER(cr)
1754             && !(inputrec->pbcType != PbcType::No && inputrec->bPeriodicMols))
1755         {
1756             /* Make molecules whole at start of run */
1757             if (fr->pbcType != PbcType::No)
1758             {
1759                 do_pbc_first_mtop(fplog, inputrec->pbcType, box, &mtop, globalState->x.rvec_array());
1760             }
1761             if (vsite)
1762             {
1763                 /* Correct initial vsite positions are required
1764                  * for the initial distribution in the domain decomposition
1765                  * and for the initial shell prediction.
1766                  */
1767                 constructVirtualSitesGlobal(mtop, globalState->x);
1768             }
1769         }
1770         // Make the DD reverse topology, now that any vsites that are present are available
1771         if (haveDDAtomOrdering(*cr))
1772         {
1773             dd_make_reverse_top(fplog, cr->dd, mtop, vsite.get(), *inputrec, domdecOptions.ddBondedChecking);
1774         }
1775
1776         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1777         {
1778             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1779             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1780         }
1781     }
1782     else
1783     {
1784         /* This is a PME only node */
1785
1786         GMX_ASSERT(globalState == nullptr,
1787                    "We don't need the state on a PME only rank and expect it to be unitialized");
1788
1789         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1790         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1791     }
1792
1793     gmx_pme_t* sepPmeData = nullptr;
1794     // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1795     GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr),
1796                "Double-checking that only PME-only ranks have no forcerec");
1797     gmx_pme_t*& pmedata = fr ? fr->pmedata : sepPmeData;
1798
1799     // TODO should live in ewald module once its testing is improved
1800     //
1801     // Later, this program could contain kernels that might be later
1802     // re-used as auto-tuning progresses, or subsequent simulations
1803     // are invoked.
1804     PmeGpuProgramStorage pmeGpuProgram;
1805     if (thisRankHasPmeGpuTask)
1806     {
1807         GMX_RELEASE_ASSERT(
1808                 (deviceStreamManager != nullptr),
1809                 "GPU device stream manager should be initialized in order to use GPU for PME.");
1810         GMX_RELEASE_ASSERT((deviceInfo != nullptr),
1811                            "GPU device should be initialized in order to use GPU for PME.");
1812         pmeGpuProgram = buildPmeGpuProgram(deviceStreamManager->context());
1813     }
1814
1815     /* Initiate PME if necessary,
1816      * either on all nodes or on dedicated PME nodes only. */
1817     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1818     {
1819         if (mdAtoms && mdAtoms->mdatoms())
1820         {
1821             nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1822             if (EVDW_PME(inputrec->vdwtype))
1823             {
1824                 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1825             }
1826         }
1827         if (cr->npmenodes > 0)
1828         {
1829             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1830             gmx_bcast(sizeof(nChargePerturbed), &nChargePerturbed, cr->mpi_comm_mysim);
1831             gmx_bcast(sizeof(nTypePerturbed), &nTypePerturbed, cr->mpi_comm_mysim);
1832         }
1833
1834         if (thisRankHasDuty(cr, DUTY_PME))
1835         {
1836             try
1837             {
1838                 // TODO: This should be in the builder.
1839                 GMX_RELEASE_ASSERT(!runScheduleWork.simulationWork.useGpuPme
1840                                            || (deviceStreamManager != nullptr),
1841                                    "Device stream manager should be valid in order to use GPU "
1842                                    "version of PME.");
1843                 GMX_RELEASE_ASSERT(
1844                         !runScheduleWork.simulationWork.useGpuPme
1845                                 || deviceStreamManager->streamIsValid(DeviceStreamType::Pme),
1846                         "GPU PME stream should be valid in order to use GPU version of PME.");
1847
1848                 const DeviceContext* deviceContext = runScheduleWork.simulationWork.useGpuPme
1849                                                              ? &deviceStreamManager->context()
1850                                                              : nullptr;
1851                 const DeviceStream*  pmeStream =
1852                         runScheduleWork.simulationWork.useGpuPme
1853                                  ? &deviceStreamManager->stream(DeviceStreamType::Pme)
1854                                  : nullptr;
1855
1856                 pmedata = gmx_pme_init(cr,
1857                                        getNumPmeDomains(cr->dd),
1858                                        inputrec.get(),
1859                                        nChargePerturbed != 0,
1860                                        nTypePerturbed != 0,
1861                                        mdrunOptions.reproducible,
1862                                        ewaldcoeff_q,
1863                                        ewaldcoeff_lj,
1864                                        gmx_omp_nthreads_get(ModuleMultiThread::Pme),
1865                                        pmeRunMode,
1866                                        nullptr,
1867                                        deviceContext,
1868                                        pmeStream,
1869                                        pmeGpuProgram.get(),
1870                                        mdlog);
1871             }
1872             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1873         }
1874     }
1875
1876
1877     if (EI_DYNAMICS(inputrec->eI))
1878     {
1879         /* Turn on signal handling on all nodes */
1880         /*
1881          * (A user signal from the PME nodes (if any)
1882          * is communicated to the PP nodes.
1883          */
1884         signal_handler_install();
1885     }
1886
1887     pull_t* pull_work = nullptr;
1888     if (thisRankHasDuty(cr, DUTY_PP))
1889     {
1890         /* Assumes uniform use of the number of OpenMP threads */
1891         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(ModuleMultiThread::Default));
1892
1893         if (inputrec->bPull)
1894         {
1895             /* Initialize pull code */
1896             pull_work = init_pull(fplog,
1897                                   inputrec->pull.get(),
1898                                   inputrec.get(),
1899                                   mtop,
1900                                   cr,
1901                                   &atomSets,
1902                                   inputrec->fepvals->init_lambda);
1903             if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1904             {
1905                 initPullHistory(pull_work, &observablesHistory);
1906             }
1907             if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1908             {
1909                 init_pull_output_files(pull_work, filenames.size(), filenames.data(), oenv, startingBehavior);
1910             }
1911         }
1912
1913         std::unique_ptr<EnforcedRotation> enforcedRotation;
1914         if (inputrec->bRot)
1915         {
1916             /* Initialize enforced rotation code */
1917             enforcedRotation = init_rot(fplog,
1918                                         inputrec.get(),
1919                                         filenames.size(),
1920                                         filenames.data(),
1921                                         cr,
1922                                         &atomSets,
1923                                         globalState.get(),
1924                                         mtop,
1925                                         oenv,
1926                                         mdrunOptions,
1927                                         startingBehavior);
1928         }
1929
1930         t_swap* swap = nullptr;
1931         if (inputrec->eSwapCoords != SwapType::No)
1932         {
1933             /* Initialize ion swapping code */
1934             swap = init_swapcoords(fplog,
1935                                    inputrec.get(),
1936                                    opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1937                                    mtop,
1938                                    globalState.get(),
1939                                    &observablesHistory,
1940                                    cr,
1941                                    &atomSets,
1942                                    oenv,
1943                                    mdrunOptions,
1944                                    startingBehavior);
1945         }
1946
1947         /* Let makeConstraints know whether we have essential dynamics constraints. */
1948         auto constr = makeConstraints(mtop,
1949                                       *inputrec,
1950                                       pull_work,
1951                                       doEssentialDynamics,
1952                                       fplog,
1953                                       cr,
1954                                       updateGroups.useUpdateGroups(),
1955                                       ms,
1956                                       &nrnb,
1957                                       wcycle.get(),
1958                                       fr->bMolPBC,
1959                                       &observablesReducerBuilder);
1960
1961         /* Energy terms and groups */
1962         gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1963                              inputrec->fepvals->n_lambda);
1964
1965         // cos acceleration is only supported by md, but older tpr
1966         // files might still combine it with other integrators
1967         GMX_RELEASE_ASSERT(inputrec->cos_accel == 0.0 || inputrec->eI == IntegrationAlgorithm::MD,
1968                            "cos_acceleration is only supported by integrator=md");
1969
1970         /* Kinetic energy data */
1971         gmx_ekindata_t ekind(inputrec->opts.ngtc,
1972                              inputrec->cos_accel,
1973                              gmx_omp_nthreads_get(ModuleMultiThread::Update));
1974
1975         /* Set up interactive MD (IMD) */
1976         auto imdSession = makeImdSession(inputrec.get(),
1977                                          cr,
1978                                          wcycle.get(),
1979                                          &enerd,
1980                                          ms,
1981                                          mtop,
1982                                          mdlog,
1983                                          MASTER(cr) ? globalState->x : gmx::ArrayRef<gmx::RVec>(),
1984                                          filenames.size(),
1985                                          filenames.data(),
1986                                          oenv,
1987                                          mdrunOptions.imdOptions,
1988                                          startingBehavior);
1989
1990         if (haveDDAtomOrdering(*cr))
1991         {
1992             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1993             /* This call is not included in init_domain_decomposition
1994              * because fr->atomInfoForEachMoleculeBlock is set later.
1995              */
1996             makeBondedLinks(cr->dd, mtop, fr->atomInfoForEachMoleculeBlock);
1997         }
1998
1999         if (runScheduleWork.simulationWork.useGpuFBufferOps)
2000         {
2001             fr->gpuForceReduction[gmx::AtomLocality::Local] = std::make_unique<gmx::GpuForceReduction>(
2002                     deviceStreamManager->context(),
2003                     deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedLocal),
2004                     wcycle.get());
2005             fr->gpuForceReduction[gmx::AtomLocality::NonLocal] = std::make_unique<gmx::GpuForceReduction>(
2006                     deviceStreamManager->context(),
2007                     deviceStreamManager->stream(gmx::DeviceStreamType::NonBondedNonLocal),
2008                     wcycle.get());
2009         }
2010
2011         std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
2012         if (gpusWereDetected
2013             && ((runScheduleWork.simulationWork.useGpuPme && thisRankHasDuty(cr, DUTY_PME))
2014                 || runScheduleWork.simulationWork.useGpuXBufferOps))
2015         {
2016             GpuApiCallBehavior transferKind =
2017                     (inputrec->eI == IntegrationAlgorithm::MD && !doRerun && !useModularSimulator)
2018                             ? GpuApiCallBehavior::Async
2019                             : GpuApiCallBehavior::Sync;
2020             GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
2021                                "GPU device stream manager should be initialized to use GPU.");
2022             stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
2023                     *deviceStreamManager, transferKind, pme_gpu_get_block_size(fr->pmedata), wcycle.get());
2024             fr->stateGpu = stateGpu.get();
2025         }
2026
2027         GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
2028         SimulatorBuilder simulatorBuilder;
2029
2030         simulatorBuilder.add(SimulatorStateData(
2031                 globalState.get(), localState, &observablesHistory, &enerd, &ekind));
2032         simulatorBuilder.add(std::move(membedHolder));
2033         simulatorBuilder.add(std::move(stopHandlerBuilder_));
2034         simulatorBuilder.add(SimulatorConfig(mdrunOptions, startingBehavior, &runScheduleWork));
2035
2036
2037         simulatorBuilder.add(SimulatorEnv(fplog, cr, ms, mdlog, oenv, &observablesReducerBuilder));
2038         simulatorBuilder.add(Profiling(&nrnb, walltime_accounting, wcycle.get()));
2039         simulatorBuilder.add(ConstraintsParam(
2040                 constr.get(), enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr, vsite.get()));
2041         // TODO: Separate `fr` to a separate add, and make the `build` handle the coupling sensibly.
2042         simulatorBuilder.add(LegacyInput(
2043                 static_cast<int>(filenames.size()), filenames.data(), inputrec.get(), fr.get()));
2044         simulatorBuilder.add(ReplicaExchangeParameters(replExParams));
2045         simulatorBuilder.add(InteractiveMD(imdSession.get()));
2046         simulatorBuilder.add(SimulatorModules(mdModules_->outputProvider(), mdModules_->notifiers()));
2047         simulatorBuilder.add(CenterOfMassPulling(pull_work));
2048         // Todo move to an MDModule
2049         simulatorBuilder.add(IonSwapping(swap));
2050         simulatorBuilder.add(TopologyData(mtop, &localTopology, mdAtoms.get()));
2051         simulatorBuilder.add(BoxDeformationHandle(deform.get()));
2052         simulatorBuilder.add(std::move(modularSimulatorCheckpointData));
2053
2054         // build and run simulator object based on user-input
2055         auto simulator = simulatorBuilder.build(useModularSimulator);
2056         simulator->run();
2057
2058         if (fr->pmePpCommGpu)
2059         {
2060             // destroy object since it is no longer required. (This needs to be done while the GPU context still exists.)
2061             fr->pmePpCommGpu.reset();
2062         }
2063
2064         if (inputrec->bPull)
2065         {
2066             finish_pull(pull_work);
2067         }
2068         finish_swapcoords(swap);
2069     }
2070     else
2071     {
2072         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
2073         /* do PME only */
2074         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(ModuleMultiThread::Pme));
2075         gmx_pmeonly(pmedata,
2076                     cr,
2077                     &nrnb,
2078                     wcycle.get(),
2079                     walltime_accounting,
2080                     inputrec.get(),
2081                     pmeRunMode,
2082                     runScheduleWork.simulationWork.useGpuPmePpCommunication,
2083                     deviceStreamManager.get());
2084     }
2085
2086     wallcycle_stop(wcycle.get(), WallCycleCounter::Run);
2087
2088     /* Finish up, write some stuff
2089      * if rerunMD, don't write last frame again
2090      */
2091     finish_run(fplog,
2092                mdlog,
2093                cr,
2094                *inputrec,
2095                &nrnb,
2096                wcycle.get(),
2097                walltime_accounting,
2098                fr ? fr->nbv.get() : nullptr,
2099                pmedata,
2100                EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
2101
2102
2103     deviceStreamManager.reset(nullptr);
2104     // Free PME data
2105     if (pmedata)
2106     {
2107         gmx_pme_destroy(pmedata);
2108         pmedata = nullptr;
2109     }
2110
2111     // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
2112     // before we destroy the GPU context(s)
2113     // Pinned buffers are associated with contexts in CUDA.
2114     // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
2115     mdAtoms.reset(nullptr);
2116     globalState.reset(nullptr);
2117     localStateInstance.reset(nullptr);
2118     mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
2119     fr.reset(nullptr);         // destruct forcerec before gpu
2120     // TODO convert to C++ so we can get rid of these frees
2121     sfree(disresdata);
2122
2123     if (!hwinfo_->deviceInfoList.empty())
2124     {
2125         /* stop the GPU profiler (only CUDA) */
2126         stopGpuProfiler();
2127     }
2128
2129     /* With tMPI we need to wait for all ranks to finish deallocation before
2130      * destroying the CUDA context as some tMPI ranks may be sharing
2131      * GPU and context.
2132      *
2133      * This is not a concern in OpenCL where we use one context per rank.
2134      *
2135      * Note: it is safe to not call the barrier on the ranks which do not use GPU,
2136      * but it is easier and more futureproof to call it on the whole node.
2137      *
2138      * Note that this function needs to be called even if GPUs are not used
2139      * in this run because the PME ranks have no knowledge of whether GPUs
2140      * are used or not, but all ranks need to enter the barrier below.
2141      * \todo Remove this physical node barrier after making sure
2142      * that it's not needed anymore (with a shared GPU run).
2143      */
2144     if (GMX_THREAD_MPI)
2145     {
2146         physicalNodeComm.barrier();
2147     }
2148
2149     if (!devFlags.usingCudaAwareMpi)
2150     {
2151         // Don't reset GPU in case of CUDA-AWARE MPI
2152         // UCX creates CUDA buffers which are cleaned-up as part of MPI_Finalize()
2153         // resetting the device before MPI_Finalize() results in crashes inside UCX
2154         releaseDevice(deviceInfo);
2155     }
2156
2157     /* Does what it says */
2158     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
2159     walltime_accounting_destroy(walltime_accounting);
2160
2161     // Ensure log file content is written
2162     if (logFileHandle)
2163     {
2164         gmx_fio_flush(logFileHandle);
2165     }
2166
2167     /* Reset FPEs (important for unit tests) by disabling them. Assumes no
2168      * exceptions were enabled before function was called. */
2169     if (bEnableFPE)
2170     {
2171         gmx_fedisableexcept();
2172     }
2173
2174     auto rc = static_cast<int>(gmx_get_stop_condition());
2175
2176 #if GMX_THREAD_MPI
2177     /* we need to join all threads. The sub-threads join when they
2178        exit this function, but the master thread needs to be told to
2179        wait for that. */
2180     if (MASTER(cr))
2181     {
2182         tMPI_Finalize();
2183     }
2184 #endif
2185     return rc;
2186 } // namespace gmx
2187
2188 Mdrunner::~Mdrunner()
2189 {
2190     // Clean up of the Manager.
2191     // This will end up getting called on every thread-MPI rank, which is unnecessary,
2192     // but okay as long as threads synchronize some time before adding or accessing
2193     // a new set of restraints.
2194     if (restraintManager_)
2195     {
2196         restraintManager_->clear();
2197         GMX_ASSERT(restraintManager_->countRestraints() == 0,
2198                    "restraints added during runner life time should be cleared at runner "
2199                    "destruction.");
2200     }
2201 };
2202
2203 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller, const std::string& name)
2204 {
2205     GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
2206     // Not sure if this should be logged through the md logger or something else,
2207     // but it is helpful to have some sort of INFO level message sent somewhere.
2208     //    std::cout << "Registering restraint named " << name << std::endl;
2209
2210     // When multiple restraints are used, it may be wasteful to register them separately.
2211     // Maybe instead register an entire Restraint Manager as a force provider.
2212     restraintManager_->addToSpec(std::move(puller), name);
2213 }
2214
2215 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules) : mdModules_(std::move(mdModules)) {}
2216
2217 Mdrunner::Mdrunner(Mdrunner&&) noexcept = default;
2218
2219 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265 in CentOS 7
2220 Mdrunner& Mdrunner::operator=(Mdrunner&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
2221
2222 class Mdrunner::BuilderImplementation
2223 {
2224 public:
2225     BuilderImplementation() = delete;
2226     BuilderImplementation(std::unique_ptr<MDModules> mdModules, compat::not_null<SimulationContext*> context);
2227     ~BuilderImplementation();
2228
2229     BuilderImplementation& setExtraMdrunOptions(const MdrunOptions& options,
2230                                                 real                forceWarningThreshold,
2231                                                 StartingBehavior    startingBehavior);
2232
2233     void addHardwareDetectionResult(const gmx_hw_info_t* hwinfo);
2234
2235     void addDomdec(const DomdecOptions& options);
2236
2237     void addInput(SimulationInputHandle inputHolder);
2238
2239     void addVerletList(int nstlist);
2240
2241     void addReplicaExchange(const ReplicaExchangeParameters& params);
2242
2243     void addNonBonded(const char* nbpu_opt);
2244
2245     void addPME(const char* pme_opt_, const char* pme_fft_opt_);
2246
2247     void addBondedTaskAssignment(const char* bonded_opt);
2248
2249     void addUpdateTaskAssignment(const char* update_opt);
2250
2251     void addHardwareOptions(const gmx_hw_opt_t& hardwareOptions);
2252
2253     void addFilenames(ArrayRef<const t_filenm> filenames);
2254
2255     void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
2256
2257     void addLogFile(t_fileio* logFileHandle);
2258
2259     void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
2260
2261     Mdrunner build();
2262
2263 private:
2264     // Default parameters copied from runner.h
2265     // \todo Clarify source(s) of default parameters.
2266
2267     const char* nbpu_opt_    = nullptr;
2268     const char* pme_opt_     = nullptr;
2269     const char* pme_fft_opt_ = nullptr;
2270     const char* bonded_opt_  = nullptr;
2271     const char* update_opt_  = nullptr;
2272
2273     MdrunOptions mdrunOptions_;
2274
2275     DomdecOptions domdecOptions_;
2276
2277     ReplicaExchangeParameters replicaExchangeParameters_;
2278
2279     //! Command-line override for the duration of a neighbor list with the Verlet scheme.
2280     int nstlist_ = 0;
2281
2282     //! World communicator, used for hardware detection and task assignment
2283     MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
2284
2285     //! Multisim communicator handle.
2286     gmx_multisim_t* multiSimulation_;
2287
2288     //! mdrun communicator
2289     MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
2290
2291     //! Print a warning if any force is larger than this (in kJ/mol nm).
2292     real forceWarningThreshold_ = -1;
2293
2294     //! Whether the simulation will start afresh, or restart with/without appending.
2295     StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
2296
2297     //! The modules that comprise the functionality of mdrun.
2298     std::unique_ptr<MDModules> mdModules_;
2299
2300     //! Detected hardware.
2301     const gmx_hw_info_t* hwinfo_ = nullptr;
2302
2303     //! \brief Parallelism information.
2304     gmx_hw_opt_t hardwareOptions_;
2305
2306     //! filename options for simulation.
2307     ArrayRef<const t_filenm> filenames_;
2308
2309     /*! \brief Handle to output environment.
2310      *
2311      * \todo gmx_output_env_t needs lifetime management.
2312      */
2313     gmx_output_env_t* outputEnvironment_ = nullptr;
2314
2315     /*! \brief Non-owning handle to MD log file.
2316      *
2317      * \todo Context should own output facilities for client.
2318      * \todo Improve log file handle management.
2319      * \internal
2320      * Code managing the FILE* relies on the ability to set it to
2321      * nullptr to check whether the filehandle is valid.
2322      */
2323     t_fileio* logFileHandle_ = nullptr;
2324
2325     /*!
2326      * \brief Builder for simulation stop signal handler.
2327      */
2328     std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
2329
2330     /*!
2331      * \brief Sources for initial simulation state.
2332      *
2333      * See issue #3652 for near-term refinements to the SimulationInput interface.
2334      *
2335      * See issue #3379 for broader discussion on API aspects of simulation inputs and outputs.
2336      */
2337     SimulationInputHandle inputHolder_;
2338 };
2339
2340 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules,
2341                                                        compat::not_null<SimulationContext*> context) :
2342     mdModules_(std::move(mdModules))
2343 {
2344     libraryWorldCommunicator_ = context->libraryWorldCommunicator_;
2345     simulationCommunicator_   = context->simulationCommunicator_;
2346     multiSimulation_          = context->multiSimulation_.get();
2347 }
2348
2349 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
2350
2351 Mdrunner::BuilderImplementation&
2352 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions&    options,
2353                                                       const real             forceWarningThreshold,
2354                                                       const StartingBehavior startingBehavior)
2355 {
2356     mdrunOptions_          = options;
2357     forceWarningThreshold_ = forceWarningThreshold;
2358     startingBehavior_      = startingBehavior;
2359     return *this;
2360 }
2361
2362 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions& options)
2363 {
2364     domdecOptions_ = options;
2365 }
2366
2367 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
2368 {
2369     nstlist_ = nstlist;
2370 }
2371
2372 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters& params)
2373 {
2374     replicaExchangeParameters_ = params;
2375 }
2376
2377 Mdrunner Mdrunner::BuilderImplementation::build()
2378 {
2379     auto newRunner = Mdrunner(std::move(mdModules_));
2380
2381     newRunner.mdrunOptions     = mdrunOptions_;
2382     newRunner.pforce           = forceWarningThreshold_;
2383     newRunner.startingBehavior = startingBehavior_;
2384     newRunner.domdecOptions    = domdecOptions_;
2385
2386     // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
2387     newRunner.hw_opt = hardwareOptions_;
2388
2389     // No invariant to check. This parameter exists to optionally override other behavior.
2390     newRunner.nstlist_cmdline = nstlist_;
2391
2392     newRunner.replExParams = replicaExchangeParameters_;
2393
2394     newRunner.filenames = filenames_;
2395
2396     newRunner.libraryWorldCommunicator = libraryWorldCommunicator_;
2397
2398     newRunner.simulationCommunicator = simulationCommunicator_;
2399
2400     // nullptr is a valid value for the multisim handle
2401     newRunner.ms = multiSimulation_;
2402
2403     if (hwinfo_)
2404     {
2405         newRunner.hwinfo_ = hwinfo_;
2406     }
2407     else
2408     {
2409         GMX_THROW(gmx::APIError(
2410                 "MdrunnerBuilder::addHardwareDetectionResult() is required before build()"));
2411     }
2412
2413     if (inputHolder_)
2414     {
2415         newRunner.inputHolder_ = std::move(inputHolder_);
2416     }
2417     else
2418     {
2419         GMX_THROW(gmx::APIError("MdrunnerBuilder::addInput() is required before build()."));
2420     }
2421
2422     // \todo Clarify ownership and lifetime management for gmx_output_env_t
2423     // \todo Update sanity checking when output environment has clearly specified invariants.
2424     // Initialization and default values for oenv are not well specified in the current version.
2425     if (outputEnvironment_)
2426     {
2427         newRunner.oenv = outputEnvironment_;
2428     }
2429     else
2430     {
2431         GMX_THROW(gmx::APIError(
2432                 "MdrunnerBuilder::addOutputEnvironment() is required before build()"));
2433     }
2434
2435     newRunner.logFileHandle = logFileHandle_;
2436
2437     if (nbpu_opt_)
2438     {
2439         newRunner.nbpu_opt = nbpu_opt_;
2440     }
2441     else
2442     {
2443         GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
2444     }
2445
2446     if (pme_opt_ && pme_fft_opt_)
2447     {
2448         newRunner.pme_opt     = pme_opt_;
2449         newRunner.pme_fft_opt = pme_fft_opt_;
2450     }
2451     else
2452     {
2453         GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
2454     }
2455
2456     if (bonded_opt_)
2457     {
2458         newRunner.bonded_opt = bonded_opt_;
2459     }
2460     else
2461     {
2462         GMX_THROW(gmx::APIError(
2463                 "MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
2464     }
2465
2466     if (update_opt_)
2467     {
2468         newRunner.update_opt = update_opt_;
2469     }
2470     else
2471     {
2472         GMX_THROW(gmx::APIError(
2473                 "MdrunnerBuilder::addUpdateTaskAssignment() is required before build()  "));
2474     }
2475
2476
2477     newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
2478
2479     if (stopHandlerBuilder_)
2480     {
2481         newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
2482     }
2483     else
2484     {
2485         newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
2486     }
2487
2488     return newRunner;
2489 }
2490
2491 void Mdrunner::BuilderImplementation::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2492 {
2493     hwinfo_ = hwinfo;
2494 }
2495
2496 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
2497 {
2498     nbpu_opt_ = nbpu_opt;
2499 }
2500
2501 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt, const char* pme_fft_opt)
2502 {
2503     pme_opt_     = pme_opt;
2504     pme_fft_opt_ = pme_fft_opt;
2505 }
2506
2507 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
2508 {
2509     bonded_opt_ = bonded_opt;
2510 }
2511
2512 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
2513 {
2514     update_opt_ = update_opt;
2515 }
2516
2517 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2518 {
2519     hardwareOptions_ = hardwareOptions;
2520 }
2521
2522 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2523 {
2524     filenames_ = filenames;
2525 }
2526
2527 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2528 {
2529     outputEnvironment_ = outputEnvironment;
2530 }
2531
2532 void Mdrunner::BuilderImplementation::addLogFile(t_fileio* logFileHandle)
2533 {
2534     logFileHandle_ = logFileHandle;
2535 }
2536
2537 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2538 {
2539     stopHandlerBuilder_ = std::move(builder);
2540 }
2541
2542 void Mdrunner::BuilderImplementation::addInput(SimulationInputHandle inputHolder)
2543 {
2544     inputHolder_ = std::move(inputHolder);
2545 }
2546
2547 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules>           mdModules,
2548                                  compat::not_null<SimulationContext*> context) :
2549     impl_{ std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules), context) }
2550 {
2551 }
2552
2553 MdrunnerBuilder::~MdrunnerBuilder() = default;
2554
2555 MdrunnerBuilder& MdrunnerBuilder::addHardwareDetectionResult(const gmx_hw_info_t* hwinfo)
2556 {
2557     impl_->addHardwareDetectionResult(hwinfo);
2558     return *this;
2559 }
2560
2561 MdrunnerBuilder& MdrunnerBuilder::addSimulationMethod(const MdrunOptions&    options,
2562                                                       real                   forceWarningThreshold,
2563                                                       const StartingBehavior startingBehavior)
2564 {
2565     impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2566     return *this;
2567 }
2568
2569 MdrunnerBuilder& MdrunnerBuilder::addDomainDecomposition(const DomdecOptions& options)
2570 {
2571     impl_->addDomdec(options);
2572     return *this;
2573 }
2574
2575 MdrunnerBuilder& MdrunnerBuilder::addNeighborList(int nstlist)
2576 {
2577     impl_->addVerletList(nstlist);
2578     return *this;
2579 }
2580
2581 MdrunnerBuilder& MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters& params)
2582 {
2583     impl_->addReplicaExchange(params);
2584     return *this;
2585 }
2586
2587 MdrunnerBuilder& MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2588 {
2589     impl_->addNonBonded(nbpu_opt);
2590     return *this;
2591 }
2592
2593 MdrunnerBuilder& MdrunnerBuilder::addElectrostatics(const char* pme_opt, const char* pme_fft_opt)
2594 {
2595     // The builder method may become more general in the future, but in this version,
2596     // parameters for PME electrostatics are both required and the only parameters
2597     // available.
2598     if (pme_opt && pme_fft_opt)
2599     {
2600         impl_->addPME(pme_opt, pme_fft_opt);
2601     }
2602     else
2603     {
2604         GMX_THROW(
2605                 gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2606     }
2607     return *this;
2608 }
2609
2610 MdrunnerBuilder& MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2611 {
2612     impl_->addBondedTaskAssignment(bonded_opt);
2613     return *this;
2614 }
2615
2616 MdrunnerBuilder& MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2617 {
2618     impl_->addUpdateTaskAssignment(update_opt);
2619     return *this;
2620 }
2621
2622 Mdrunner MdrunnerBuilder::build()
2623 {
2624     return impl_->build();
2625 }
2626
2627 MdrunnerBuilder& MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t& hardwareOptions)
2628 {
2629     impl_->addHardwareOptions(hardwareOptions);
2630     return *this;
2631 }
2632
2633 MdrunnerBuilder& MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2634 {
2635     impl_->addFilenames(filenames);
2636     return *this;
2637 }
2638
2639 MdrunnerBuilder& MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2640 {
2641     impl_->addOutputEnvironment(outputEnvironment);
2642     return *this;
2643 }
2644
2645 MdrunnerBuilder& MdrunnerBuilder::addLogFile(t_fileio* logFileHandle)
2646 {
2647     impl_->addLogFile(logFileHandle);
2648     return *this;
2649 }
2650
2651 MdrunnerBuilder& MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2652 {
2653     impl_->addStopHandlerBuilder(std::move(builder));
2654     return *this;
2655 }
2656
2657 MdrunnerBuilder& MdrunnerBuilder::addInput(SimulationInputHandle input)
2658 {
2659     impl_->addInput(std::move(input));
2660     return *this;
2661 }
2662
2663 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder&&) noexcept = default;
2664
2665 MdrunnerBuilder& MdrunnerBuilder::operator=(MdrunnerBuilder&&) noexcept = default;
2666
2667 } // namespace gmx