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