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