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