Fix segmentation fault with cosine pulling
[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
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/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_ = compat::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_ = compat::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 = compat::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             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
625                     (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
626                     *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
627
628         }
629         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
630
631         /* Determine how many thread-MPI ranks to start.
632          *
633          * TODO Over-writing the user-supplied value here does
634          * prevent any possible subsequent checks from working
635          * correctly. */
636         hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
637                                                 &hw_opt,
638                                                 gpuIdsToUse,
639                                                 useGpuForNonbonded,
640                                                 useGpuForPme,
641                                                 inputrec, &mtop,
642                                                 mdlog,
643                                                 doMembed);
644
645         // Now start the threads for thread MPI.
646         cr = spawnThreads(hw_opt.nthreads_tmpi);
647         /* The main thread continues here with a new cr. We don't deallocate
648            the old cr because other threads may still be reading it. */
649         // TODO Both master and spawned threads call dup_tfn and
650         // reinitialize_commrec_for_this_thread. Find a way to express
651         // this better.
652         physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
653     }
654     // END OF CAUTION: cr and physicalNodeComm are now reliable
655
656     if (PAR(cr))
657     {
658         /* now broadcast everything to the non-master nodes/threads: */
659         init_parallel(cr, inputrec, &mtop);
660     }
661
662     // Now each rank knows the inputrec that SIMMASTER read and used,
663     // and (if applicable) cr->nnodes has been assigned the number of
664     // thread-MPI ranks that have been chosen. The ranks can now all
665     // run the task-deciding functions and will agree on the result
666     // without needing to communicate.
667     //
668     // TODO Should we do the communication in debug mode to support
669     // having an assertion?
670     //
671     // Note that these variables describe only their own node.
672     //
673     // Note that when bonded interactions run on a GPU they always run
674     // alongside a nonbonded task, so do not influence task assignment
675     // even though they affect the force calculation workload.
676     bool useGpuForNonbonded = false;
677     bool useGpuForPme       = false;
678     bool useGpuForBonded    = false;
679     try
680     {
681         // It's possible that there are different numbers of GPUs on
682         // different nodes, which is the user's responsibilty to
683         // handle. If unsuitable, we will notice that during task
684         // assignment.
685         bool gpusWereDetected      = hwinfo->ngpu_compatible_tot > 0;
686         bool usingVerletScheme     = inputrec->cutoff_scheme == ecutsVERLET;
687         auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
688         useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
689                                                                 emulateGpuNonbonded,
690                                                                 canUseGpuForNonbonded,
691                                                                 usingVerletScheme,
692                                                                 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
693                                                                 gpusWereDetected);
694         useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
695                                                     *hwinfo, *inputrec, mtop,
696                                                     cr->nnodes, domdecOptions.numPmeRanks,
697                                                     gpusWereDetected);
698         auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
699         useGpuForBonded =
700             decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
701                                             bondedTarget, canUseGpuForBonded,
702                                             EVDW_PME(inputrec->vdwtype),
703                                             EEL_PME_EWALD(inputrec->coulombtype),
704                                             domdecOptions.numPmeRanks, gpusWereDetected);
705
706         pmeRunMode   = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
707         if (pmeRunMode == PmeRunMode::GPU)
708         {
709             if (pmeFftTarget == TaskTarget::Cpu)
710             {
711                 pmeRunMode = PmeRunMode::Mixed;
712             }
713         }
714         else if (pmeFftTarget == TaskTarget::Gpu)
715         {
716             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.");
717         }
718     }
719     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
720
721     // Build restraints.
722     // TODO: hide restraint implementation details from Mdrunner.
723     // There is nothing unique about restraints at this point as far as the
724     // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
725     // factory functions from the SimulationContext on which to call mdModules->add().
726     // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
727     for (auto && restraint : restraintManager_->getRestraints())
728     {
729         auto module = RestraintMDModule::create(restraint,
730                                                 restraint->sites());
731         mdModules->add(std::move(module));
732     }
733
734     // TODO: Error handling
735     mdModules->assignOptionsToModules(*inputrec->params, nullptr);
736
737     if (fplog != nullptr)
738     {
739         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
740         fprintf(fplog, "\n");
741     }
742
743     if (SIMMASTER(cr))
744     {
745         /* now make sure the state is initialized and propagated */
746         set_state_entries(globalState.get(), inputrec);
747     }
748
749     /* NM and TPI parallelize over force/energy calculations, not atoms,
750      * so we need to initialize and broadcast the global state.
751      */
752     if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
753     {
754         if (!MASTER(cr))
755         {
756             globalState = compat::make_unique<t_state>();
757         }
758         broadcastStateWithoutDynamics(cr, globalState.get());
759     }
760
761     /* A parallel command line option consistency check that we can
762        only do after any threads have started. */
763     if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
764                      domdecOptions.numCells[YY] > 1 ||
765                      domdecOptions.numCells[ZZ] > 1 ||
766                      domdecOptions.numPmeRanks > 0))
767     {
768         gmx_fatal(FARGS,
769                   "The -dd or -npme option request a parallel simulation, "
770 #if !GMX_MPI
771                   "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
772 #else
773 #if GMX_THREAD_MPI
774                   "but the number of MPI-threads (option -ntmpi) is not set or is 1");
775 #else
776                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
777 #endif
778 #endif
779     }
780
781     if (doRerun &&
782         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
783     {
784         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
785     }
786
787     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
788     {
789         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
790     }
791
792     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
793     {
794         if (domdecOptions.numPmeRanks > 0)
795         {
796             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
797                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
798         }
799
800         domdecOptions.numPmeRanks = 0;
801     }
802
803     if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
804     {
805         /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
806          * improve performance with many threads per GPU, since our OpenMP
807          * scaling is bad, but it's difficult to automate the setup.
808          */
809         domdecOptions.numPmeRanks = 0;
810     }
811     if (useGpuForPme)
812     {
813         if (domdecOptions.numPmeRanks < 0)
814         {
815             domdecOptions.numPmeRanks = 0;
816             // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
817         }
818         else
819         {
820             GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
821         }
822     }
823
824 #if GMX_FAHCORE
825     if (MASTER(cr))
826     {
827         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
828     }
829 #endif
830
831     /* NMR restraints must be initialized before load_checkpoint,
832      * since with time averaging the history is added to t_state.
833      * For proper consistency check we therefore need to extend
834      * t_state here.
835      * So the PME-only nodes (if present) will also initialize
836      * the distance restraints.
837      */
838     snew(fcd, 1);
839
840     /* This needs to be called before read_checkpoint to extend the state */
841     init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
842
843     init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
844
845     auto                 deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
846
847     ObservablesHistory   observablesHistory = {};
848
849     ContinuationOptions &continuationOptions = mdrunOptions.continuationOptions;
850
851     if (continuationOptions.startedFromCheckpoint)
852     {
853         /* Check if checkpoint file exists before doing continuation.
854          * This way we can use identical input options for the first and subsequent runs...
855          */
856         gmx_bool bReadEkin;
857
858         if (mdrunOptions.numStepsCommandline > -2)
859         {
860             /* Temporarily set the number of steps to unmlimited to avoid
861              * triggering the nsteps check in load_checkpoint().
862              * This hack will go away soon when the -nsteps option is removed.
863              */
864             inputrec->nsteps = -1;
865         }
866
867         load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
868                         logFileHandle,
869                         cr, domdecOptions.numCells,
870                         inputrec, globalState.get(),
871                         &bReadEkin, &observablesHistory,
872                         continuationOptions.appendFiles,
873                         continuationOptions.appendFilesOptionSet,
874                         mdrunOptions.reproducible);
875
876         if (bReadEkin)
877         {
878             continuationOptions.haveReadEkin = true;
879         }
880
881         if (continuationOptions.appendFiles && logFileHandle)
882         {
883             // Now we can start normal logging to the truncated log file.
884             fplog    = gmx_fio_getfp(logFileHandle);
885             prepareLogAppending(fplog);
886             logOwner = buildLogger(fplog, cr);
887             mdlog    = logOwner.logger();
888         }
889     }
890
891     if (mdrunOptions.numStepsCommandline > -2)
892     {
893         GMX_LOG(mdlog.info).asParagraph().
894             appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
895                        "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
896     }
897     /* override nsteps with value set on the commamdline */
898     override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
899
900     if (SIMMASTER(cr))
901     {
902         copy_mat(globalState->box, box);
903     }
904
905     if (PAR(cr))
906     {
907         gmx_bcast(sizeof(box), box, cr);
908     }
909
910     /* Update rlist and nstlist. */
911     if (inputrec->cutoff_scheme == ecutsVERLET)
912     {
913         prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
914                               useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
915     }
916
917     LocalAtomSetManager atomSets;
918
919     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
920                      inputrec->eI == eiNM))
921     {
922         cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
923                                            &mtop, inputrec,
924                                            box, positionsFromStatePointer(globalState.get()),
925                                            &atomSets);
926         // Note that local state still does not exist yet.
927     }
928     else
929     {
930         /* PME, if used, is done on all nodes with 1D decomposition */
931         cr->npmenodes = 0;
932         cr->duty      = (DUTY_PP | DUTY_PME);
933
934         if (inputrec->ePBC == epbcSCREW)
935         {
936             gmx_fatal(FARGS,
937                       "pbc=screw is only implemented with domain decomposition");
938         }
939     }
940
941     if (PAR(cr))
942     {
943         /* After possible communicator splitting in make_dd_communicators.
944          * we can set up the intra/inter node communication.
945          */
946         gmx_setup_nodecomm(fplog, cr);
947     }
948
949 #if GMX_MPI
950     if (isMultiSim(ms))
951     {
952         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
953                 "This is simulation %d out of %d running as a composite GROMACS\n"
954                 "multi-simulation job. Setup for this simulation:\n",
955                 ms->sim, ms->nsim);
956     }
957     GMX_LOG(mdlog.warning).appendTextFormatted(
958             "Using %d MPI %s\n",
959             cr->nnodes,
960 #if GMX_THREAD_MPI
961             cr->nnodes == 1 ? "thread" : "threads"
962 #else
963             cr->nnodes == 1 ? "process" : "processes"
964 #endif
965             );
966     fflush(stderr);
967 #endif
968
969     /* Check and update hw_opt for the cut-off scheme */
970     check_and_update_hw_opt_2(&hw_opt, inputrec->cutoff_scheme);
971
972     /* Check and update the number of OpenMP threads requested */
973     checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
974                                             pmeRunMode, mtop);
975
976     gmx_omp_nthreads_init(mdlog, cr,
977                           hwinfo->nthreads_hw_avail,
978                           physicalNodeComm.size_,
979                           hw_opt.nthreads_omp,
980                           hw_opt.nthreads_omp_pme,
981                           !thisRankHasDuty(cr, DUTY_PP),
982                           inputrec->cutoff_scheme == ecutsVERLET);
983
984     // Enable FP exception detection for the Verlet scheme, but not in
985     // Release mode and not for compilers with known buggy FP
986     // exception support (clang with any optimization) or suspected
987     // buggy FP exception support (gcc 7.* with optimization).
988 #if !defined NDEBUG && \
989     !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
990     && defined __OPTIMIZE__)
991     const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
992 #else
993     const bool bEnableFPE = false;
994 #endif
995     //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
996     if (bEnableFPE)
997     {
998         gmx_feenableexcept();
999     }
1000
1001     // Build a data structure that expresses which kinds of non-bonded
1002     // task are handled by this rank.
1003     //
1004     // TODO Later, this might become a loop over all registered modules
1005     // relevant to the mdp inputs, to find those that have such tasks.
1006     //
1007     // TODO This could move before init_domain_decomposition() as part
1008     // of refactoring that separates the responsibility for duty
1009     // assignment from setup for communication between tasks, and
1010     // setup for tasks handled with a domain (ie including short-ranged
1011     // tasks, bonded tasks, etc.).
1012     //
1013     // Note that in general useGpuForNonbonded, etc. can have a value
1014     // that is inconsistent with the presence of actual GPUs on any
1015     // rank, and that is not known to be a problem until the
1016     // duty of the ranks on a node become known.
1017     //
1018     // TODO Later we might need the concept of computeTasksOnThisRank,
1019     // from which we construct gpuTasksOnThisRank.
1020     //
1021     // Currently the DD code assigns duty to ranks that can
1022     // include PP work that currently can be executed on a single
1023     // GPU, if present and compatible.  This has to be coordinated
1024     // across PP ranks on a node, with possible multiple devices
1025     // or sharing devices on a node, either from the user
1026     // selection, or automatically.
1027     auto                 haveGpus = !gpuIdsToUse.empty();
1028     std::vector<GpuTask> gpuTasksOnThisRank;
1029     if (thisRankHasDuty(cr, DUTY_PP))
1030     {
1031         if (useGpuForNonbonded)
1032         {
1033             // Note that any bonded tasks on a GPU always accompany a
1034             // non-bonded task.
1035             if (haveGpus)
1036             {
1037                 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1038             }
1039             else if (nonbondedTarget == TaskTarget::Gpu)
1040             {
1041                 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1042             }
1043             else if (bondedTarget == TaskTarget::Gpu)
1044             {
1045                 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because there is none detected.");
1046             }
1047         }
1048     }
1049     // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1050     if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1051     {
1052         if (useGpuForPme)
1053         {
1054             if (haveGpus)
1055             {
1056                 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1057             }
1058             else if (pmeTarget == TaskTarget::Gpu)
1059             {
1060                 gmx_fatal(FARGS, "Cannot run PME on a GPU because there is none detected.");
1061             }
1062         }
1063     }
1064
1065     GpuTaskAssignment gpuTaskAssignment;
1066     try
1067     {
1068         // Produce the task assignment for this rank.
1069         gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1070                                               mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1071                                               useGpuForBonded, pmeRunMode);
1072     }
1073     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1074
1075     /* Prevent other ranks from continuing after an issue was found
1076      * and reported as a fatal error.
1077      *
1078      * TODO This function implements a barrier so that MPI runtimes
1079      * can organize an orderly shutdown if one of the ranks has had to
1080      * issue a fatal error in various code already run. When we have
1081      * MPI-aware error handling and reporting, this should be
1082      * improved. */
1083 #if GMX_MPI
1084     if (PAR(cr))
1085     {
1086         MPI_Barrier(cr->mpi_comm_mysim);
1087     }
1088     if (isMultiSim(ms))
1089     {
1090         if (SIMMASTER(cr))
1091         {
1092             MPI_Barrier(ms->mpi_comm_masters);
1093         }
1094         /* We need another barrier to prevent non-master ranks from contiuing
1095          * when an error occured in a different simulation.
1096          */
1097         MPI_Barrier(cr->mpi_comm_mysim);
1098     }
1099 #endif
1100
1101     /* Now that we know the setup is consistent, check for efficiency */
1102     check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1103                                        cr, mdlog);
1104
1105     gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1106
1107     if (thisRankHasDuty(cr, DUTY_PP))
1108     {
1109         // This works because only one task of each type is currently permitted.
1110         auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1111                                              hasTaskType<GpuTask::Nonbonded>);
1112         if (nbGpuTaskMapping != gpuTaskAssignment.end())
1113         {
1114             int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1115             nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1116             init_gpu(nonbondedDeviceInfo);
1117
1118             if (DOMAINDECOMP(cr))
1119             {
1120                 /* When we share GPUs over ranks, we need to know this for the DLB */
1121                 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1122             }
1123
1124         }
1125     }
1126
1127     std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
1128
1129     gmx_device_info_t                *pmeDeviceInfo = nullptr;
1130     // Later, this program could contain kernels that might be later
1131     // re-used as auto-tuning progresses, or subsequent simulations
1132     // are invoked.
1133     PmeGpuProgramStorage pmeGpuProgram;
1134     // This works because only one task of each type is currently permitted.
1135     auto                 pmeGpuTaskMapping     = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1136     const bool           thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1137     if (thisRankHasPmeGpuTask)
1138     {
1139         pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1140         init_gpu(pmeDeviceInfo);
1141         pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1142         // TODO It would be nice to move this logic into the factory
1143         // function. See Redmine #2535.
1144         bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
1145         if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
1146         {
1147             initializedClfftLibrary = initializeClfftLibrary();
1148         }
1149     }
1150
1151     /* getting number of PP/PME threads
1152        PME: env variable should be read only on one node to make sure it is
1153        identical everywhere;
1154      */
1155     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1156
1157     int numThreadsOnThisRank;
1158     /* threads on this MPI process or TMPI thread */
1159     if (thisRankHasDuty(cr, DUTY_PP))
1160     {
1161         numThreadsOnThisRank = gmx_omp_nthreads_get(emntNonbonded);
1162     }
1163     else
1164     {
1165         numThreadsOnThisRank = nthreads_pme;
1166     }
1167
1168     checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1169                                   *hwinfo->hardwareTopology,
1170                                   physicalNodeComm, mdlog);
1171
1172     if (hw_opt.thread_affinity != threadaffOFF)
1173     {
1174         /* Before setting affinity, check whether the affinity has changed
1175          * - which indicates that probably the OpenMP library has changed it
1176          * since we first checked).
1177          */
1178         gmx_check_thread_affinity_set(mdlog, cr,
1179                                       &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1180
1181         int numThreadsOnThisNode, intraNodeThreadOffset;
1182         analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1183                                  &intraNodeThreadOffset);
1184
1185         /* Set the CPU affinity */
1186         gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1187                                 numThreadsOnThisRank, numThreadsOnThisNode,
1188                                 intraNodeThreadOffset, nullptr);
1189     }
1190
1191     if (mdrunOptions.timingOptions.resetStep > -1)
1192     {
1193         GMX_LOG(mdlog.info).asParagraph().
1194             appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1195     }
1196     wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1197
1198     if (PAR(cr))
1199     {
1200         /* Master synchronizes its value of reset_counters with all nodes
1201          * including PME only nodes */
1202         reset_counters = wcycle_get_reset_counters(wcycle);
1203         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1204         wcycle_set_reset_counters(wcycle, reset_counters);
1205     }
1206
1207     // Membrane embedding must be initialized before we call init_forcerec()
1208     if (doMembed)
1209     {
1210         if (MASTER(cr))
1211         {
1212             fprintf(stderr, "Initializing membed");
1213         }
1214         /* Note that membed cannot work in parallel because mtop is
1215          * changed here. Fix this if we ever want to make it run with
1216          * multiple ranks. */
1217         membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1218                              &mdrunOptions
1219                                  .checkpointOptions.period);
1220     }
1221
1222     std::unique_ptr<MDAtoms>     mdAtoms;
1223     std::unique_ptr<gmx_vsite_t> vsite;
1224
1225     snew(nrnb, 1);
1226     if (thisRankHasDuty(cr, DUTY_PP))
1227     {
1228         /* Initiate forcerecord */
1229         fr                 = mk_forcerec();
1230         fr->forceProviders = mdModules->initForceProviders();
1231         init_forcerec(fplog, mdlog, fr, fcd,
1232                       inputrec, &mtop, cr, box,
1233                       opt2fn("-table", filenames.size(), filenames.data()),
1234                       opt2fn("-tablep", filenames.size(), filenames.data()),
1235                       opt2fns("-tableb", filenames.size(), filenames.data()),
1236                       *hwinfo, nonbondedDeviceInfo,
1237                       useGpuForBonded,
1238                       FALSE,
1239                       pforce);
1240
1241         /* Initialize the mdAtoms structure.
1242          * mdAtoms is not filled with atom data,
1243          * as this can not be done now with domain decomposition.
1244          */
1245         mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1246         if (globalState && thisRankHasPmeGpuTask)
1247         {
1248             // The pinning of coordinates in the global state object works, because we only use
1249             // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1250             // points to the global state object without DD.
1251             // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1252             // which should also perform the pinning.
1253             changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1254         }
1255
1256         /* Initialize the virtual site communication */
1257         vsite = initVsite(mtop, cr);
1258
1259         calc_shifts(box, fr->shift_vec);
1260
1261         /* With periodic molecules the charge groups should be whole at start up
1262          * and the virtual sites should not be far from their proper positions.
1263          */
1264         if (!inputrec->bContinuation && MASTER(cr) &&
1265             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1266         {
1267             /* Make molecules whole at start of run */
1268             if (fr->ePBC != epbcNONE)
1269             {
1270                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1271             }
1272             if (vsite)
1273             {
1274                 /* Correct initial vsite positions are required
1275                  * for the initial distribution in the domain decomposition
1276                  * and for the initial shell prediction.
1277                  */
1278                 constructVsitesGlobal(mtop, globalState->x);
1279             }
1280         }
1281
1282         if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1283         {
1284             ewaldcoeff_q  = fr->ic->ewaldcoeff_q;
1285             ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1286         }
1287     }
1288     else
1289     {
1290         /* This is a PME only node */
1291
1292         GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1293
1294         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1295         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1296     }
1297
1298     gmx_pme_t *sepPmeData = nullptr;
1299     // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1300     GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1301     gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1302
1303     /* Initiate PME if necessary,
1304      * either on all nodes or on dedicated PME nodes only. */
1305     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1306     {
1307         if (mdAtoms && mdAtoms->mdatoms())
1308         {
1309             nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1310             if (EVDW_PME(inputrec->vdwtype))
1311             {
1312                 nTypePerturbed   = mdAtoms->mdatoms()->nTypePerturbed;
1313             }
1314         }
1315         if (cr->npmenodes > 0)
1316         {
1317             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1318             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1319             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1320         }
1321
1322         if (thisRankHasDuty(cr, DUTY_PME))
1323         {
1324             try
1325             {
1326                 pmedata = gmx_pme_init(cr,
1327                                        getNumPmeDomains(cr->dd),
1328                                        inputrec,
1329                                        mtop.natoms, nChargePerturbed != 0, nTypePerturbed != 0,
1330                                        mdrunOptions.reproducible,
1331                                        ewaldcoeff_q, ewaldcoeff_lj,
1332                                        nthreads_pme,
1333                                        pmeRunMode, nullptr,
1334                                        pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1335             }
1336             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1337         }
1338     }
1339
1340
1341     if (EI_DYNAMICS(inputrec->eI))
1342     {
1343         /* Turn on signal handling on all nodes */
1344         /*
1345          * (A user signal from the PME nodes (if any)
1346          * is communicated to the PP nodes.
1347          */
1348         signal_handler_install();
1349     }
1350
1351     if (thisRankHasDuty(cr, DUTY_PP))
1352     {
1353         /* Assumes uniform use of the number of OpenMP threads */
1354         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1355
1356         if (inputrec->bPull)
1357         {
1358             /* Initialize pull code */
1359             inputrec->pull_work =
1360                 init_pull(fplog, inputrec->pull, inputrec,
1361                           &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1362             if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1363             {
1364                 initPullHistory(inputrec->pull_work, &observablesHistory);
1365             }
1366             if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1367             {
1368                 init_pull_output_files(inputrec->pull_work,
1369                                        filenames.size(), filenames.data(), oenv,
1370                                        continuationOptions);
1371             }
1372         }
1373
1374         std::unique_ptr<EnforcedRotation> enforcedRotation;
1375         if (inputrec->bRot)
1376         {
1377             /* Initialize enforced rotation code */
1378             enforcedRotation = init_rot(fplog,
1379                                         inputrec,
1380                                         filenames.size(),
1381                                         filenames.data(),
1382                                         cr,
1383                                         &atomSets,
1384                                         globalState.get(),
1385                                         &mtop,
1386                                         oenv,
1387                                         mdrunOptions);
1388         }
1389
1390         if (inputrec->eSwapCoords != eswapNO)
1391         {
1392             /* Initialize ion swapping code */
1393             init_swapcoords(fplog, inputrec, opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1394                             &mtop, globalState.get(), &observablesHistory,
1395                             cr, &atomSets, oenv, mdrunOptions);
1396         }
1397
1398         /* Let makeConstraints know whether we have essential dynamics constraints.
1399          * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1400          */
1401         bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1402                                     || observablesHistory.edsamHistory);
1403         auto constr              = makeConstraints(mtop, *inputrec, doEssentialDynamics,
1404                                                    fplog, *mdAtoms->mdatoms(),
1405                                                    cr, ms, nrnb, wcycle, fr->bMolPBC);
1406
1407         if (DOMAINDECOMP(cr))
1408         {
1409             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1410             /* This call is not included in init_domain_decomposition mainly
1411              * because fr->cginfo_mb is set later.
1412              */
1413             dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1414                             domdecOptions.checkBondedInteractions,
1415                             fr->cginfo_mb);
1416         }
1417
1418         // TODO This is not the right place to manage the lifetime of
1419         // this data structure, but currently it's the easiest way to
1420         // make it work. Later, it should probably be made/updated
1421         // after the workload for the lifetime of a PP domain is
1422         // understood.
1423         PpForceWorkload ppForceWorkload;
1424
1425         GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to integrator.");
1426         /* Now do whatever the user wants us to do (how flexible...) */
1427         Integrator integrator {
1428             fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1429             oenv,
1430             mdrunOptions,
1431             vsite.get(), constr.get(),
1432             enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1433             deform.get(),
1434             mdModules->outputProvider(),
1435             inputrec, &mtop,
1436             fcd,
1437             globalState.get(),
1438             &observablesHistory,
1439             mdAtoms.get(), nrnb, wcycle, fr,
1440             &ppForceWorkload,
1441             replExParams,
1442             membed,
1443             walltime_accounting,
1444             std::move(stopHandlerBuilder_)
1445         };
1446         integrator.run(inputrec->eI, doRerun);
1447
1448         if (inputrec->bPull)
1449         {
1450             finish_pull(inputrec->pull_work);
1451         }
1452
1453     }
1454     else
1455     {
1456         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1457         /* do PME only */
1458         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1459         gmx_pmeonly(pmedata, cr, nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1460     }
1461
1462     wallcycle_stop(wcycle, ewcRUN);
1463
1464     /* Finish up, write some stuff
1465      * if rerunMD, don't write last frame again
1466      */
1467     finish_run(fplog, mdlog, cr,
1468                inputrec, nrnb, wcycle, walltime_accounting,
1469                fr ? fr->nbv : nullptr,
1470                pmedata,
1471                EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1472
1473     // Free PME data
1474     if (pmedata)
1475     {
1476         gmx_pme_destroy(pmedata);
1477         pmedata = nullptr;
1478     }
1479
1480     // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1481     // before we destroy the GPU context(s) in free_gpu_resources().
1482     // Pinned buffers are associated with contexts in CUDA.
1483     // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1484     mdAtoms.reset(nullptr);
1485     globalState.reset(nullptr);
1486     mdModules.reset(nullptr);   // destruct force providers here as they might also use the GPU
1487
1488     /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1489     free_gpu_resources(fr, physicalNodeComm);
1490     free_gpu(nonbondedDeviceInfo);
1491     free_gpu(pmeDeviceInfo);
1492     done_forcerec(fr, mtop.molblock.size(), mtop.groups.grps[egcENER].nr);
1493     sfree(fcd);
1494
1495     if (doMembed)
1496     {
1497         free_membed(membed);
1498     }
1499
1500     gmx_hardware_info_free();
1501
1502     /* Does what it says */
1503     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1504     walltime_accounting_destroy(walltime_accounting);
1505     sfree(nrnb);
1506
1507     // Ensure log file content is written
1508     if (logFileHandle)
1509     {
1510         gmx_fio_flush(logFileHandle);
1511     }
1512
1513     /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1514      * exceptions were enabled before function was called. */
1515     if (bEnableFPE)
1516     {
1517         gmx_fedisableexcept();
1518     }
1519
1520     rc = static_cast<int>(gmx_get_stop_condition());
1521
1522 #if GMX_THREAD_MPI
1523     /* we need to join all threads. The sub-threads join when they
1524        exit this function, but the master thread needs to be told to
1525        wait for that. */
1526     if (PAR(cr) && MASTER(cr))
1527     {
1528         done_commrec(cr);
1529         tMPI_Finalize();
1530     }
1531 #endif
1532
1533     return rc;
1534 }
1535
1536 Mdrunner::~Mdrunner()
1537 {
1538     // Clean up of the Manager.
1539     // This will end up getting called on every thread-MPI rank, which is unnecessary,
1540     // but okay as long as threads synchronize some time before adding or accessing
1541     // a new set of restraints.
1542     if (restraintManager_)
1543     {
1544         restraintManager_->clear();
1545         GMX_ASSERT(restraintManager_->countRestraints() == 0,
1546                    "restraints added during runner life time should be cleared at runner destruction.");
1547     }
1548 };
1549
1550 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1551                             std::string                               name)
1552 {
1553     GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1554     // Not sure if this should be logged through the md logger or something else,
1555     // but it is helpful to have some sort of INFO level message sent somewhere.
1556     //    std::cout << "Registering restraint named " << name << std::endl;
1557
1558     // When multiple restraints are used, it may be wasteful to register them separately.
1559     // Maybe instead register an entire Restraint Manager as a force provider.
1560     restraintManager_->addToSpec(std::move(puller),
1561                                  std::move(name));
1562 }
1563
1564 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1565
1566 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1567 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1568
1569 class Mdrunner::BuilderImplementation
1570 {
1571     public:
1572         BuilderImplementation() = delete;
1573         explicit BuilderImplementation(SimulationContext* context);
1574         ~BuilderImplementation();
1575
1576         BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1577                                                     real                forceWarningThreshold);
1578
1579         void addDomdec(const DomdecOptions &options);
1580
1581         void addVerletList(int nstlist);
1582
1583         void addReplicaExchange(const ReplicaExchangeParameters &params);
1584
1585         void addMultiSim(gmx_multisim_t* multisim);
1586
1587         void addNonBonded(const char* nbpu_opt);
1588
1589         void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1590
1591         void addBondedTaskAssignment(const char* bonded_opt);
1592
1593         void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1594
1595         void addFilenames(ArrayRef <const t_filenm> filenames);
1596
1597         void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1598
1599         void addLogFile(t_fileio *logFileHandle);
1600
1601         void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1602
1603         Mdrunner build();
1604
1605     private:
1606         // Default parameters copied from runner.h
1607         // \todo Clarify source(s) of default parameters.
1608
1609         const char* nbpu_opt_    = nullptr;
1610         const char* pme_opt_     = nullptr;
1611         const char* pme_fft_opt_ = nullptr;
1612         const char *bonded_opt_  = nullptr;
1613
1614         MdrunOptions                          mdrunOptions_;
1615
1616         DomdecOptions                         domdecOptions_;
1617
1618         ReplicaExchangeParameters             replicaExchangeParameters_;
1619
1620         //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1621         int         nstlist_ = 0;
1622
1623         //! Non-owning multisim communicator handle.
1624         std::unique_ptr<gmx_multisim_t*>      multisim_ = nullptr;
1625
1626         //! Print a warning if any force is larger than this (in kJ/mol nm).
1627         real forceWarningThreshold_ = -1;
1628
1629         /*! \brief  Non-owning pointer to SimulationContext (owned and managed by client)
1630          *
1631          * \internal
1632          * \todo Establish robust protocol to make sure resources remain valid.
1633          * SimulationContext will likely be separated into multiple layers for
1634          * different levels of access and different phases of execution. Ref
1635          * https://redmine.gromacs.org/issues/2375
1636          * https://redmine.gromacs.org/issues/2587
1637          */
1638         SimulationContext* context_ = nullptr;
1639
1640         //! \brief Parallelism information.
1641         gmx_hw_opt_t hardwareOptions_;
1642
1643         //! filename options for simulation.
1644         ArrayRef<const t_filenm> filenames_;
1645
1646         /*! \brief Handle to output environment.
1647          *
1648          * \todo gmx_output_env_t needs lifetime management.
1649          */
1650         gmx_output_env_t*    outputEnvironment_ = nullptr;
1651
1652         /*! \brief Non-owning handle to MD log file.
1653          *
1654          * \todo Context should own output facilities for client.
1655          * \todo Improve log file handle management.
1656          * \internal
1657          * Code managing the FILE* relies on the ability to set it to
1658          * nullptr to check whether the filehandle is valid.
1659          */
1660         t_fileio* logFileHandle_ = nullptr;
1661
1662         /*!
1663          * \brief Builder for simulation stop signal handler.
1664          */
1665         std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1666 };
1667
1668 Mdrunner::BuilderImplementation::BuilderImplementation(SimulationContext* context) :
1669     context_(context)
1670 {
1671     GMX_ASSERT(context_, "Bug found. It should not be possible to construct builder without a valid context.");
1672 }
1673
1674 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1675
1676 Mdrunner::BuilderImplementation &
1677 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1678                                                       real                forceWarningThreshold)
1679 {
1680     mdrunOptions_          = options;
1681     forceWarningThreshold_ = forceWarningThreshold;
1682     return *this;
1683 }
1684
1685 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1686 {
1687     domdecOptions_ = options;
1688 }
1689
1690 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1691 {
1692     nstlist_ = nstlist;
1693 }
1694
1695 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters &params)
1696 {
1697     replicaExchangeParameters_ = params;
1698 }
1699
1700 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1701 {
1702     multisim_ = compat::make_unique<gmx_multisim_t*>(multisim);
1703 }
1704
1705 Mdrunner Mdrunner::BuilderImplementation::build()
1706 {
1707     auto newRunner = Mdrunner();
1708
1709     GMX_ASSERT(context_, "Bug found. It should not be possible to call build() without a valid context.");
1710
1711     newRunner.mdrunOptions    = mdrunOptions_;
1712     newRunner.domdecOptions   = domdecOptions_;
1713
1714     // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1715     newRunner.hw_opt          = hardwareOptions_;
1716
1717     // No invariant to check. This parameter exists to optionally override other behavior.
1718     newRunner.nstlist_cmdline = nstlist_;
1719
1720     newRunner.replExParams    = replicaExchangeParameters_;
1721
1722     newRunner.filenames = filenames_;
1723
1724     GMX_ASSERT(context_->communicationRecord_, "SimulationContext communications not initialized.");
1725     newRunner.cr = context_->communicationRecord_;
1726
1727     if (multisim_)
1728     {
1729         // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1730         newRunner.ms = *multisim_;
1731     }
1732     else
1733     {
1734         GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1735     }
1736
1737     // \todo Clarify ownership and lifetime management for gmx_output_env_t
1738     // \todo Update sanity checking when output environment has clearly specified invariants.
1739     // Initialization and default values for oenv are not well specified in the current version.
1740     if (outputEnvironment_)
1741     {
1742         newRunner.oenv  = outputEnvironment_;
1743     }
1744     else
1745     {
1746         GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1747     }
1748
1749     newRunner.logFileHandle = logFileHandle_;
1750
1751     if (nbpu_opt_)
1752     {
1753         newRunner.nbpu_opt    = nbpu_opt_;
1754     }
1755     else
1756     {
1757         GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1758     }
1759
1760     if (pme_opt_ && pme_fft_opt_)
1761     {
1762         newRunner.pme_opt     = pme_opt_;
1763         newRunner.pme_fft_opt = pme_fft_opt_;
1764     }
1765     else
1766     {
1767         GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1768     }
1769
1770     if (bonded_opt_)
1771     {
1772         newRunner.bonded_opt = bonded_opt_;
1773     }
1774     else
1775     {
1776         GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1777     }
1778
1779     newRunner.restraintManager_ = compat::make_unique<gmx::RestraintManager>();
1780
1781     if (stopHandlerBuilder_)
1782     {
1783         newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1784     }
1785     else
1786     {
1787         newRunner.stopHandlerBuilder_ = compat::make_unique<StopHandlerBuilder>();
1788     }
1789
1790     return newRunner;
1791 }
1792
1793 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1794 {
1795     nbpu_opt_ = nbpu_opt;
1796 }
1797
1798 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1799                                              const char* pme_fft_opt)
1800 {
1801     pme_opt_     = pme_opt;
1802     pme_fft_opt_ = pme_fft_opt;
1803 }
1804
1805 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1806 {
1807     bonded_opt_ = bonded_opt;
1808 }
1809
1810 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1811 {
1812     hardwareOptions_ = hardwareOptions;
1813 }
1814
1815 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
1816 {
1817     filenames_ = filenames;
1818 }
1819
1820 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1821 {
1822     outputEnvironment_ = outputEnvironment;
1823 }
1824
1825 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
1826 {
1827     logFileHandle_ = logFileHandle;
1828 }
1829
1830 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1831 {
1832     stopHandlerBuilder_ = std::move(builder);
1833 }
1834
1835 MdrunnerBuilder::MdrunnerBuilder(compat::not_null<SimulationContext*> context) :
1836     impl_ {gmx::compat::make_unique<Mdrunner::BuilderImplementation>(context)}
1837 {
1838 }
1839
1840 MdrunnerBuilder::~MdrunnerBuilder() = default;
1841
1842 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
1843                                                       real                forceWarningThreshold)
1844 {
1845     impl_->setExtraMdrunOptions(options, forceWarningThreshold);
1846     return *this;
1847 }
1848
1849 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
1850 {
1851     impl_->addDomdec(options);
1852     return *this;
1853 }
1854
1855 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
1856 {
1857     impl_->addVerletList(nstlist);
1858     return *this;
1859 }
1860
1861 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters &params)
1862 {
1863     impl_->addReplicaExchange(params);
1864     return *this;
1865 }
1866
1867 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
1868 {
1869     impl_->addMultiSim(multisim);
1870     return *this;
1871 }
1872
1873 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
1874 {
1875     impl_->addNonBonded(nbpu_opt);
1876     return *this;
1877 }
1878
1879 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
1880                                                     const char* pme_fft_opt)
1881 {
1882     // The builder method may become more general in the future, but in this version,
1883     // parameters for PME electrostatics are both required and the only parameters
1884     // available.
1885     if (pme_opt && pme_fft_opt)
1886     {
1887         impl_->addPME(pme_opt, pme_fft_opt);
1888     }
1889     else
1890     {
1891         GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
1892     }
1893     return *this;
1894 }
1895
1896 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
1897 {
1898     impl_->addBondedTaskAssignment(bonded_opt);
1899     return *this;
1900 }
1901
1902 Mdrunner MdrunnerBuilder::build()
1903 {
1904     return impl_->build();
1905 }
1906
1907 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
1908 {
1909     impl_->addHardwareOptions(hardwareOptions);
1910     return *this;
1911 }
1912
1913 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
1914 {
1915     impl_->addFilenames(filenames);
1916     return *this;
1917 }
1918
1919 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
1920 {
1921     impl_->addOutputEnvironment(outputEnvironment);
1922     return *this;
1923 }
1924
1925 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
1926 {
1927     impl_->addLogFile(logFileHandle);
1928     return *this;
1929 }
1930
1931 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
1932 {
1933     impl_->addStopHandlerBuilder(std::move(builder));
1934     return *this;
1935 }
1936
1937 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
1938
1939 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;
1940
1941 } // namespace gmx