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