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