Remove charge group data from t_forcerec
[gromacs.git] / src / gromacs / mdrun / runner.cpp
blob690bfd748d286a2d6b6dc936bc0cc07888198cca
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /*! \internal \file
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
44 #include "gmxpre.h"
46 #include "runner.h"
48 #include "config.h"
50 #include <cassert>
51 #include <cinttypes>
52 #include <csignal>
53 #include <cstdlib>
54 #include <cstring>
56 #include <algorithm>
57 #include <memory>
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/gpuhaloexchange.h"
63 #include "gromacs/domdec/localatomsetmanager.h"
64 #include "gromacs/domdec/partition.h"
65 #include "gromacs/ewald/ewald_utils.h"
66 #include "gromacs/ewald/pme.h"
67 #include "gromacs/ewald/pme_gpu_program.h"
68 #include "gromacs/fileio/checkpoint.h"
69 #include "gromacs/fileio/gmxfio.h"
70 #include "gromacs/fileio/oenv.h"
71 #include "gromacs/fileio/tpxio.h"
72 #include "gromacs/gmxlib/network.h"
73 #include "gromacs/gmxlib/nrnb.h"
74 #include "gromacs/gpu_utils/gpu_utils.h"
75 #include "gromacs/hardware/cpuinfo.h"
76 #include "gromacs/hardware/detecthardware.h"
77 #include "gromacs/hardware/printhardware.h"
78 #include "gromacs/imd/imd.h"
79 #include "gromacs/listed_forces/disre.h"
80 #include "gromacs/listed_forces/gpubonded.h"
81 #include "gromacs/listed_forces/orires.h"
82 #include "gromacs/math/functions.h"
83 #include "gromacs/math/utilities.h"
84 #include "gromacs/math/vec.h"
85 #include "gromacs/mdlib/boxdeformation.h"
86 #include "gromacs/mdlib/broadcaststructs.h"
87 #include "gromacs/mdlib/calc_verletbuf.h"
88 #include "gromacs/mdlib/dispersioncorrection.h"
89 #include "gromacs/mdlib/enerdata_utils.h"
90 #include "gromacs/mdlib/force.h"
91 #include "gromacs/mdlib/forcerec.h"
92 #include "gromacs/mdlib/gmx_omp_nthreads.h"
93 #include "gromacs/mdlib/makeconstraints.h"
94 #include "gromacs/mdlib/md_support.h"
95 #include "gromacs/mdlib/mdatoms.h"
96 #include "gromacs/mdlib/membed.h"
97 #include "gromacs/mdlib/qmmm.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/stophandler.h"
100 #include "gromacs/mdrun/mdmodules.h"
101 #include "gromacs/mdrun/simulationcontext.h"
102 #include "gromacs/mdrunutility/handlerestart.h"
103 #include "gromacs/mdrunutility/logging.h"
104 #include "gromacs/mdrunutility/multisim.h"
105 #include "gromacs/mdrunutility/printtime.h"
106 #include "gromacs/mdrunutility/threadaffinity.h"
107 #include "gromacs/mdtypes/commrec.h"
108 #include "gromacs/mdtypes/enerdata.h"
109 #include "gromacs/mdtypes/fcdata.h"
110 #include "gromacs/mdtypes/group.h"
111 #include "gromacs/mdtypes/inputrec.h"
112 #include "gromacs/mdtypes/md_enums.h"
113 #include "gromacs/mdtypes/mdrunoptions.h"
114 #include "gromacs/mdtypes/observableshistory.h"
115 #include "gromacs/mdtypes/simulation_workload.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/nbnxm/gpu_data_mgmt.h"
118 #include "gromacs/nbnxm/nbnxm.h"
119 #include "gromacs/nbnxm/pairlist_tuning.h"
120 #include "gromacs/pbcutil/pbc.h"
121 #include "gromacs/pulling/output.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/pulling/pull_rotation.h"
124 #include "gromacs/restraint/manager.h"
125 #include "gromacs/restraint/restraintmdmodule.h"
126 #include "gromacs/restraint/restraintpotential.h"
127 #include "gromacs/swap/swapcoords.h"
128 #include "gromacs/taskassignment/decidegpuusage.h"
129 #include "gromacs/taskassignment/resourcedivision.h"
130 #include "gromacs/taskassignment/taskassignment.h"
131 #include "gromacs/taskassignment/usergpuids.h"
132 #include "gromacs/timing/gpu_timing.h"
133 #include "gromacs/timing/wallcycle.h"
134 #include "gromacs/timing/wallcyclereporting.h"
135 #include "gromacs/topology/mtop_util.h"
136 #include "gromacs/trajectory/trajectoryframe.h"
137 #include "gromacs/utility/basenetwork.h"
138 #include "gromacs/utility/cstringutil.h"
139 #include "gromacs/utility/exceptions.h"
140 #include "gromacs/utility/fatalerror.h"
141 #include "gromacs/utility/filestream.h"
142 #include "gromacs/utility/gmxassert.h"
143 #include "gromacs/utility/gmxmpi.h"
144 #include "gromacs/utility/keyvaluetree.h"
145 #include "gromacs/utility/logger.h"
146 #include "gromacs/utility/loggerbuilder.h"
147 #include "gromacs/utility/mdmodulenotification.h"
148 #include "gromacs/utility/physicalnodecommunicator.h"
149 #include "gromacs/utility/pleasecite.h"
150 #include "gromacs/utility/programcontext.h"
151 #include "gromacs/utility/smalloc.h"
152 #include "gromacs/utility/stringutil.h"
154 #include "isimulator.h"
155 #include "replicaexchange.h"
156 #include "simulatorbuilder.h"
158 #if GMX_FAHCORE
159 #include "corewrap.h"
160 #endif
162 namespace gmx
165 /*! \brief environment variable to enable GPU P2P communication */
166 static const bool c_enableGpuHaloExchange = (getenv("GMX_GPU_DD_COMMS") != nullptr)
167 && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA);
169 /*! \brief Manage any development feature flag variables encountered
171 * The use of dev features indicated by environment variables is
172 * logged in order to ensure that runs with such featrues enabled can
173 * be identified from their log and standard output. Any cross
174 * dependencies are also checked, and if unsatisified, a fatal error
175 * issued.
177 * \param[in] mdlog Logger object.
179 static void manageDevelopmentFeatures(const gmx::MDLogger &mdlog)
181 const bool enableGpuBufOps = (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
182 const bool useGpuUpdateConstrain = (getenv("GMX_UPDATE_CONSTRAIN_GPU") != nullptr);
183 const bool enableGpuHaloExchange = (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA));
185 if (enableGpuBufOps)
187 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
188 "NOTE: This run uses the 'GPU buffer ops' feature, enabled by the GMX_USE_GPU_BUFFER_OPS environment variable.");
191 if (enableGpuHaloExchange)
193 if (!enableGpuBufOps)
195 gmx_fatal(FARGS, "Cannot enable GPU halo exchange without GPU buffer operations, set GMX_USE_GPU_BUFFER_OPS=1\n");
197 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
198 "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the GMX_GPU_DD_COMMS environment variable.");
201 if (useGpuUpdateConstrain)
203 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
204 "NOTE: This run uses the 'GPU update/constraints' feature, enabled by the GMX_UPDATE_CONSTRAIN_GPU environment variable.");
209 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
211 * Used to ensure that the master thread does not modify mdrunner during copy
212 * on the spawned threads. */
213 static void threadMpiMdrunnerAccessBarrier()
215 #if GMX_THREAD_MPI
216 MPI_Barrier(MPI_COMM_WORLD);
217 #endif
220 Mdrunner Mdrunner::cloneOnSpawnedThread() const
222 auto newRunner = Mdrunner(std::make_unique<MDModules>());
224 // All runners in the same process share a restraint manager resource because it is
225 // part of the interface to the client code, which is associated only with the
226 // original thread. Handles to the same resources can be obtained by copy.
228 newRunner.restraintManager_ = std::make_unique<RestraintManager>(*restraintManager_);
231 // Copy original cr pointer before master thread can pass the thread barrier
232 newRunner.cr = reinitialize_commrec_for_this_thread(cr);
234 // Copy members of master runner.
235 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
236 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
237 newRunner.hw_opt = hw_opt;
238 newRunner.filenames = filenames;
240 newRunner.oenv = oenv;
241 newRunner.mdrunOptions = mdrunOptions;
242 newRunner.domdecOptions = domdecOptions;
243 newRunner.nbpu_opt = nbpu_opt;
244 newRunner.pme_opt = pme_opt;
245 newRunner.pme_fft_opt = pme_fft_opt;
246 newRunner.bonded_opt = bonded_opt;
247 newRunner.update_opt = update_opt;
248 newRunner.nstlist_cmdline = nstlist_cmdline;
249 newRunner.replExParams = replExParams;
250 newRunner.pforce = pforce;
251 newRunner.ms = ms;
252 newRunner.startingBehavior = startingBehavior;
253 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
255 threadMpiMdrunnerAccessBarrier();
257 GMX_RELEASE_ASSERT(!MASTER(newRunner.cr), "cloneOnSpawnedThread should only be called on spawned threads");
259 return newRunner;
262 /*! \brief The callback used for running on spawned threads.
264 * Obtains the pointer to the master mdrunner object from the one
265 * argument permitted to the thread-launch API call, copies it to make
266 * a new runner for this thread, reinitializes necessary data, and
267 * proceeds to the simulation. */
268 static void mdrunner_start_fn(const void *arg)
272 auto masterMdrunner = reinterpret_cast<const gmx::Mdrunner *>(arg);
273 /* copy the arg list to make sure that it's thread-local. This
274 doesn't copy pointed-to items, of course; fnm, cr and fplog
275 are reset in the call below, all others should be const. */
276 gmx::Mdrunner mdrunner = masterMdrunner->cloneOnSpawnedThread();
277 mdrunner.mdrunner();
279 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
283 /*! \brief Start thread-MPI threads.
285 * Called by mdrunner() to start a specific number of threads
286 * (including the main thread) for thread-parallel runs. This in turn
287 * calls mdrunner() for each thread. All options are the same as for
288 * mdrunner(). */
289 t_commrec *Mdrunner::spawnThreads(int numThreadsToLaunch) const
292 /* first check whether we even need to start tMPI */
293 if (numThreadsToLaunch < 2)
295 return cr;
298 #if GMX_THREAD_MPI
299 /* now spawn new threads that start mdrunner_start_fn(), while
300 the main thread returns. Thread affinity is handled later. */
301 if (tMPI_Init_fn(TRUE, numThreadsToLaunch, TMPI_AFFINITY_NONE,
302 mdrunner_start_fn, static_cast<const void*>(this)) != TMPI_SUCCESS)
304 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
307 threadMpiMdrunnerAccessBarrier();
308 #else
309 GMX_UNUSED_VALUE(mdrunner_start_fn);
310 #endif
312 return reinitialize_commrec_for_this_thread(cr);
315 } // namespace gmx
317 /*! \brief Initialize variables for Verlet scheme simulation */
318 static void prepare_verlet_scheme(FILE *fplog,
319 t_commrec *cr,
320 t_inputrec *ir,
321 int nstlist_cmdline,
322 const gmx_mtop_t *mtop,
323 const matrix box,
324 bool makeGpuPairList,
325 const gmx::CpuInfo &cpuinfo)
327 /* For NVE simulations, we will retain the initial list buffer */
328 if (EI_DYNAMICS(ir->eI) &&
329 ir->verletbuf_tol > 0 &&
330 !(EI_MD(ir->eI) && ir->etc == etcNO))
332 /* Update the Verlet buffer size for the current run setup */
334 /* Here we assume SIMD-enabled kernels are being used. But as currently
335 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
336 * and 4x2 gives a larger buffer than 4x4, this is ok.
338 ListSetupType listType = (makeGpuPairList ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
339 VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
341 const real rlist_new =
342 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
344 if (rlist_new != ir->rlist)
346 if (fplog != nullptr)
348 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
349 ir->rlist, rlist_new,
350 listSetup.cluster_size_i, listSetup.cluster_size_j);
352 ir->rlist = rlist_new;
356 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
358 gmx_fatal(FARGS, "Can not set nstlist without %s",
359 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
362 if (EI_DYNAMICS(ir->eI))
364 /* Set or try nstlist values */
365 increaseNstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, makeGpuPairList, cpuinfo);
369 /*! \brief Override the nslist value in inputrec
371 * with value passed on the command line (if any)
373 static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
374 int64_t nsteps_cmdline,
375 t_inputrec *ir)
377 assert(ir);
379 /* override with anything else than the default -2 */
380 if (nsteps_cmdline > -2)
382 char sbuf_steps[STEPSTRSIZE];
383 char sbuf_msg[STRLEN];
385 ir->nsteps = nsteps_cmdline;
386 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
388 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
389 gmx_step_str(nsteps_cmdline, sbuf_steps),
390 fabs(nsteps_cmdline*ir->delta_t));
392 else
394 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
395 gmx_step_str(nsteps_cmdline, sbuf_steps));
398 GMX_LOG(mdlog.warning).asParagraph().appendText(sbuf_msg);
400 else if (nsteps_cmdline < -2)
402 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64,
403 nsteps_cmdline);
405 /* Do nothing if nsteps_cmdline == -2 */
408 namespace gmx
411 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
413 * If not, and if a warning may be issued, logs a warning about
414 * falling back to CPU code. With thread-MPI, only the first
415 * call to this function should have \c issueWarning true. */
416 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger &mdlog,
417 const t_inputrec &ir,
418 bool issueWarning)
420 bool gpuIsUseful = true;
421 std::string warning;
423 if (ir.opts.ngener - ir.nwall > 1)
425 /* The GPU code does not support more than one energy group.
426 * If the user requested GPUs explicitly, a fatal error is given later.
428 gpuIsUseful = false;
429 warning =
430 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
431 "For better performance, run on the GPU without energy groups and then do "
432 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
435 if (EI_TPI(ir.eI))
437 gpuIsUseful = false;
438 warning = "TPI is not implemented for GPUs.";
441 if (!gpuIsUseful && issueWarning)
443 GMX_LOG(mdlog.warning).asParagraph().appendText(warning);
446 return gpuIsUseful;
449 //! Initializes the logger for mdrun.
450 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
452 gmx::LoggerBuilder builder;
453 if (fplog != nullptr)
455 builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
457 if (cr == nullptr || SIMMASTER(cr))
459 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
460 &gmx::TextOutputFile::standardError());
462 return builder.build();
465 //! Make a TaskTarget from an mdrun argument string.
466 static TaskTarget findTaskTarget(const char *optionString)
468 TaskTarget returnValue = TaskTarget::Auto;
470 if (strncmp(optionString, "auto", 3) == 0)
472 returnValue = TaskTarget::Auto;
474 else if (strncmp(optionString, "cpu", 3) == 0)
476 returnValue = TaskTarget::Cpu;
478 else if (strncmp(optionString, "gpu", 3) == 0)
480 returnValue = TaskTarget::Gpu;
482 else
484 GMX_ASSERT(false, "Option string should have been checked for sanity already");
487 return returnValue;
490 //! Finish run, aggregate data to print performance info.
491 static void finish_run(FILE *fplog,
492 const gmx::MDLogger &mdlog,
493 const t_commrec *cr,
494 const t_inputrec *inputrec,
495 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
496 gmx_walltime_accounting_t walltime_accounting,
497 nonbonded_verlet_t *nbv,
498 const gmx_pme_t *pme,
499 gmx_bool bWriteStat)
501 double delta_t = 0;
502 double nbfs = 0, mflop = 0;
503 double elapsed_time,
504 elapsed_time_over_all_ranks,
505 elapsed_time_over_all_threads,
506 elapsed_time_over_all_threads_over_all_ranks;
507 /* Control whether it is valid to print a report. Only the
508 simulation master may print, but it should not do so if the run
509 terminated e.g. before a scheduled reset step. This is
510 complicated by the fact that PME ranks are unaware of the
511 reason why they were sent a pmerecvqxFINISH. To avoid
512 communication deadlocks, we always do the communication for the
513 report, even if we've decided not to write the report, because
514 how long it takes to finish the run is not important when we've
515 decided not to report on the simulation performance.
517 Further, we only report performance for dynamical integrators,
518 because those are the only ones for which we plan to
519 consider doing any optimizations. */
520 bool printReport = EI_DYNAMICS(inputrec->eI) && SIMMASTER(cr);
522 if (printReport && !walltime_accounting_get_valid_finish(walltime_accounting))
524 GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
525 printReport = false;
528 t_nrnb *nrnb_tot;
529 std::unique_ptr<t_nrnb> nrnbTotalStorage;
530 if (cr->nnodes > 1)
532 nrnbTotalStorage = std::make_unique<t_nrnb>();
533 nrnb_tot = nrnbTotalStorage.get();
534 #if GMX_MPI
535 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
536 cr->mpi_comm_mysim);
537 #endif
539 else
541 nrnb_tot = nrnb;
544 elapsed_time = walltime_accounting_get_time_since_reset(walltime_accounting);
545 elapsed_time_over_all_threads = walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting);
546 if (cr->nnodes > 1)
548 #if GMX_MPI
549 /* reduce elapsed_time over all MPI ranks in the current simulation */
550 MPI_Allreduce(&elapsed_time,
551 &elapsed_time_over_all_ranks,
552 1, MPI_DOUBLE, MPI_SUM,
553 cr->mpi_comm_mysim);
554 elapsed_time_over_all_ranks /= cr->nnodes;
555 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
556 * current simulation. */
557 MPI_Allreduce(&elapsed_time_over_all_threads,
558 &elapsed_time_over_all_threads_over_all_ranks,
559 1, MPI_DOUBLE, MPI_SUM,
560 cr->mpi_comm_mysim);
561 #endif
563 else
565 elapsed_time_over_all_ranks = elapsed_time;
566 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
569 if (printReport)
571 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
574 if (thisRankHasDuty(cr, DUTY_PP) && DOMAINDECOMP(cr))
576 print_dd_statistics(cr, inputrec, fplog);
579 /* TODO Move the responsibility for any scaling by thread counts
580 * to the code that handled the thread region, so that there's a
581 * mechanism to keep cycle counting working during the transition
582 * to task parallelism. */
583 int nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
584 int nthreads_pme = gmx_omp_nthreads_get(emntPME);
585 wallcycle_scale_by_num_threads(wcycle, thisRankHasDuty(cr, DUTY_PME) && !thisRankHasDuty(cr, DUTY_PP), nthreads_pp, nthreads_pme);
586 auto cycle_sum(wallcycle_sum(cr, wcycle));
588 if (printReport)
590 auto nbnxn_gpu_timings = (nbv != nullptr && nbv->useGpu()) ? Nbnxm::gpu_get_timings(nbv->gpu_nbv) : nullptr;
591 gmx_wallclock_gpu_pme_t pme_gpu_timings = {};
593 if (pme_gpu_task_enabled(pme))
595 pme_gpu_get_timings(pme, &pme_gpu_timings);
597 wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
598 elapsed_time_over_all_ranks,
599 wcycle, cycle_sum,
600 nbnxn_gpu_timings,
601 &pme_gpu_timings);
603 if (EI_DYNAMICS(inputrec->eI))
605 delta_t = inputrec->delta_t;
608 if (fplog)
610 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
611 elapsed_time_over_all_ranks,
612 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
613 delta_t, nbfs, mflop);
615 if (bWriteStat)
617 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
618 elapsed_time_over_all_ranks,
619 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting),
620 delta_t, nbfs, mflop);
625 int Mdrunner::mdrunner()
627 matrix box;
628 t_forcerec *fr = nullptr;
629 t_fcdata *fcd = nullptr;
630 real ewaldcoeff_q = 0;
631 real ewaldcoeff_lj = 0;
632 int nChargePerturbed = -1, nTypePerturbed = 0;
633 gmx_wallcycle_t wcycle;
634 gmx_walltime_accounting_t walltime_accounting = nullptr;
635 gmx_membed_t * membed = nullptr;
636 gmx_hw_info_t *hwinfo = nullptr;
638 /* CAUTION: threads may be started later on in this function, so
639 cr doesn't reflect the final parallel state right now */
640 gmx_mtop_t mtop;
642 bool doMembed = opt2bSet("-membed", filenames.size(), filenames.data());
643 bool doRerun = mdrunOptions.rerun;
645 // Handle task-assignment related user options.
646 EmulateGpuNonbonded emulateGpuNonbonded = (getenv("GMX_EMULATE_GPU") != nullptr ?
647 EmulateGpuNonbonded::Yes : EmulateGpuNonbonded::No);
649 std::vector<int> userGpuTaskAssignment;
652 userGpuTaskAssignment = parseUserTaskAssignmentString(hw_opt.userGpuTaskAssignment);
654 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
655 auto nonbondedTarget = findTaskTarget(nbpu_opt);
656 auto pmeTarget = findTaskTarget(pme_opt);
657 auto pmeFftTarget = findTaskTarget(pme_fft_opt);
658 auto bondedTarget = findTaskTarget(bonded_opt);
659 auto updateTarget = findTaskTarget(update_opt);
660 PmeRunMode pmeRunMode = PmeRunMode::None;
662 // Here we assume that SIMMASTER(cr) does not change even after the
663 // threads are started.
665 FILE *fplog = nullptr;
666 // If we are appending, we don't write log output because we need
667 // to check that the old log file matches what the checkpoint file
668 // expects. Otherwise, we should start to write log output now if
669 // there is a file ready for it.
670 if (logFileHandle != nullptr && startingBehavior != StartingBehavior::RestartWithAppending)
672 fplog = gmx_fio_getfp(logFileHandle);
674 gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
675 gmx::MDLogger mdlog(logOwner.logger());
677 // report any development features that may be enabled by environment variables
678 manageDevelopmentFeatures(mdlog);
680 // With thread-MPI, the communicator changes after threads are
681 // launched, so this is rebuilt for the master rank at that
682 // time. The non-master ranks are fine to keep the one made here.
683 PhysicalNodeCommunicator physicalNodeComm(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
684 hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
686 gmx_print_detected_hardware(fplog, isMasterSimMasterRank(ms, MASTER(cr)), mdlog, hwinfo);
688 std::vector<int> gpuIdsToUse = makeGpuIdsToUse(hwinfo->gpu_info, hw_opt.gpuIdsAvailable);
690 // Print citation requests after all software/hardware printing
691 pleaseCiteGromacs(fplog);
693 // TODO Replace this by unique_ptr once t_inputrec is C++
694 t_inputrec inputrecInstance;
695 t_inputrec *inputrec = nullptr;
696 std::unique_ptr<t_state> globalState;
698 if (SIMMASTER(cr))
700 /* Only the master rank has the global state */
701 globalState = std::make_unique<t_state>();
703 /* Read (nearly) all data required for the simulation */
704 read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()),
705 &inputrecInstance, globalState.get(), &mtop);
706 inputrec = &inputrecInstance;
709 /* Check and update the hardware options for internal consistency */
710 checkAndUpdateHardwareOptions(mdlog, &hw_opt, SIMMASTER(cr), domdecOptions.numPmeRanks, inputrec);
712 if (GMX_THREAD_MPI && SIMMASTER(cr))
714 bool useGpuForNonbonded = false;
715 bool useGpuForPme = false;
718 GMX_RELEASE_ASSERT(inputrec != nullptr, "Keep the compiler happy");
720 // If the user specified the number of ranks, then we must
721 // respect that, but in default mode, we need to allow for
722 // the number of GPUs to choose the number of ranks.
723 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
724 useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
725 (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
726 canUseGpuForNonbonded,
727 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
728 hw_opt.nthreads_tmpi);
729 useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
730 (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
731 *hwinfo, *inputrec, mtop, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
734 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
736 /* Determine how many thread-MPI ranks to start.
738 * TODO Over-writing the user-supplied value here does
739 * prevent any possible subsequent checks from working
740 * correctly. */
741 hw_opt.nthreads_tmpi = get_nthreads_mpi(hwinfo,
742 &hw_opt,
743 gpuIdsToUse,
744 useGpuForNonbonded,
745 useGpuForPme,
746 inputrec, &mtop,
747 mdlog,
748 doMembed);
750 // Now start the threads for thread MPI.
751 cr = spawnThreads(hw_opt.nthreads_tmpi);
752 /* The main thread continues here with a new cr. We don't deallocate
753 the old cr because other threads may still be reading it. */
754 // TODO Both master and spawned threads call dup_tfn and
755 // reinitialize_commrec_for_this_thread. Find a way to express
756 // this better.
757 physicalNodeComm = PhysicalNodeCommunicator(MPI_COMM_WORLD, gmx_physicalnode_id_hash());
759 // END OF CAUTION: cr and physicalNodeComm are now reliable
761 if (PAR(cr))
763 /* now broadcast everything to the non-master nodes/threads: */
764 if (!SIMMASTER(cr))
766 inputrec = &inputrecInstance;
768 init_parallel(cr, inputrec, &mtop);
770 GMX_RELEASE_ASSERT(inputrec != nullptr, "All range should have a valid inputrec now");
772 // Now each rank knows the inputrec that SIMMASTER read and used,
773 // and (if applicable) cr->nnodes has been assigned the number of
774 // thread-MPI ranks that have been chosen. The ranks can now all
775 // run the task-deciding functions and will agree on the result
776 // without needing to communicate.
778 // TODO Should we do the communication in debug mode to support
779 // having an assertion?
781 // Note that these variables describe only their own node.
783 // Note that when bonded interactions run on a GPU they always run
784 // alongside a nonbonded task, so do not influence task assignment
785 // even though they affect the force calculation workload.
786 bool useGpuForNonbonded = false;
787 bool useGpuForPme = false;
788 bool useGpuForBonded = false;
789 bool useGpuForUpdate = false;
790 bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
793 // It's possible that there are different numbers of GPUs on
794 // different nodes, which is the user's responsibilty to
795 // handle. If unsuitable, we will notice that during task
796 // assignment.
797 auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
798 useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
799 emulateGpuNonbonded,
800 canUseGpuForNonbonded,
801 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI),
802 gpusWereDetected);
803 useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
804 *hwinfo, *inputrec, mtop,
805 cr->nnodes, domdecOptions.numPmeRanks,
806 gpusWereDetected);
807 auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
808 useGpuForBonded =
809 decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme,
810 bondedTarget, canUseGpuForBonded,
811 EVDW_PME(inputrec->vdwtype),
812 EEL_PME_EWALD(inputrec->coulombtype),
813 domdecOptions.numPmeRanks, gpusWereDetected);
815 pmeRunMode = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
816 if (pmeRunMode == PmeRunMode::GPU)
818 if (pmeFftTarget == TaskTarget::Cpu)
820 pmeRunMode = PmeRunMode::Mixed;
823 else if (pmeFftTarget == TaskTarget::Gpu)
825 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.");
828 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
830 // Build restraints.
831 // TODO: hide restraint implementation details from Mdrunner.
832 // There is nothing unique about restraints at this point as far as the
833 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
834 // factory functions from the SimulationContext on which to call mdModules_->add().
835 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
836 for (auto && restraint : restraintManager_->getRestraints())
838 auto module = RestraintMDModule::create(restraint,
839 restraint->sites());
840 mdModules_->add(std::move(module));
843 // TODO: Error handling
844 mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
845 const auto &mdModulesNotifier = mdModules_->notifier().notifier_;
847 if (inputrec->internalParameters != nullptr)
849 mdModulesNotifier.notify(*inputrec->internalParameters);
852 if (fplog != nullptr)
854 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
855 fprintf(fplog, "\n");
858 if (SIMMASTER(cr))
860 /* In rerun, set velocities to zero if present */
861 if (doRerun && ((globalState->flags & (1 << estV)) != 0))
863 // rerun does not use velocities
864 GMX_LOG(mdlog.info).asParagraph().appendText(
865 "Rerun trajectory contains velocities. Rerun does only evaluate "
866 "potential energy and forces. The velocities will be ignored.");
867 for (int i = 0; i < globalState->natoms; i++)
869 clear_rvec(globalState->v[i]);
871 globalState->flags &= ~(1 << estV);
874 /* now make sure the state is initialized and propagated */
875 set_state_entries(globalState.get(), inputrec);
878 /* NM and TPI parallelize over force/energy calculations, not atoms,
879 * so we need to initialize and broadcast the global state.
881 if (inputrec->eI == eiNM || inputrec->eI == eiTPI)
883 if (!MASTER(cr))
885 globalState = std::make_unique<t_state>();
887 broadcastStateWithoutDynamics(cr, globalState.get());
890 /* A parallel command line option consistency check that we can
891 only do after any threads have started. */
892 if (!PAR(cr) && (domdecOptions.numCells[XX] > 1 ||
893 domdecOptions.numCells[YY] > 1 ||
894 domdecOptions.numCells[ZZ] > 1 ||
895 domdecOptions.numPmeRanks > 0))
897 gmx_fatal(FARGS,
898 "The -dd or -npme option request a parallel simulation, "
899 #if !GMX_MPI
900 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
901 #else
902 #if GMX_THREAD_MPI
903 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
904 #else
905 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
906 #endif
907 #endif
910 if (doRerun &&
911 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
913 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
916 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
918 if (domdecOptions.numPmeRanks > 0)
920 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
921 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
924 domdecOptions.numPmeRanks = 0;
927 if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
929 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
930 * improve performance with many threads per GPU, since our OpenMP
931 * scaling is bad, but it's difficult to automate the setup.
933 domdecOptions.numPmeRanks = 0;
935 if (useGpuForPme)
937 if (domdecOptions.numPmeRanks < 0)
939 domdecOptions.numPmeRanks = 0;
940 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
942 else
944 GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1, "PME GPU decomposition is not supported");
948 #if GMX_FAHCORE
949 if (MASTER(cr))
951 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
953 #endif
955 /* NMR restraints must be initialized before load_checkpoint,
956 * since with time averaging the history is added to t_state.
957 * For proper consistency check we therefore need to extend
958 * t_state here.
959 * So the PME-only nodes (if present) will also initialize
960 * the distance restraints.
962 snew(fcd, 1);
964 /* This needs to be called before read_checkpoint to extend the state */
965 init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
967 init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
969 auto deform = prepareBoxDeformation(globalState->box, cr, *inputrec);
971 ObservablesHistory observablesHistory = {};
973 if (startingBehavior != StartingBehavior::NewSimulation)
975 /* Check if checkpoint file exists before doing continuation.
976 * This way we can use identical input options for the first and subsequent runs...
978 if (mdrunOptions.numStepsCommandline > -2)
980 /* Temporarily set the number of steps to unmlimited to avoid
981 * triggering the nsteps check in load_checkpoint().
982 * This hack will go away soon when the -nsteps option is removed.
984 inputrec->nsteps = -1;
987 load_checkpoint(opt2fn_master("-cpi", filenames.size(), filenames.data(), cr),
988 logFileHandle,
989 cr, domdecOptions.numCells,
990 inputrec, globalState.get(),
991 &observablesHistory,
992 mdrunOptions.reproducible, mdModules_->notifier());
994 if (startingBehavior == StartingBehavior::RestartWithAppending && logFileHandle)
996 // Now we can start normal logging to the truncated log file.
997 fplog = gmx_fio_getfp(logFileHandle);
998 prepareLogAppending(fplog);
999 logOwner = buildLogger(fplog, cr);
1000 mdlog = logOwner.logger();
1004 if (mdrunOptions.numStepsCommandline > -2)
1006 GMX_LOG(mdlog.info).asParagraph().
1007 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
1008 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
1010 /* override nsteps with value set on the commamdline */
1011 override_nsteps_cmdline(mdlog, mdrunOptions.numStepsCommandline, inputrec);
1013 if (SIMMASTER(cr))
1015 copy_mat(globalState->box, box);
1018 if (PAR(cr))
1020 gmx_bcast(sizeof(box), box, cr);
1023 if (inputrec->cutoff_scheme != ecutsVERLET)
1025 gmx_fatal(FARGS, "This group-scheme .tpr file can no longer be run by mdrun. Please update to the Verlet scheme, or use an earlier version of GROMACS if necessary.");
1027 /* Update rlist and nstlist. */
1028 prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
1029 useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
1031 LocalAtomSetManager atomSets;
1032 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1033 inputrec->eI == eiNM))
1035 cr->dd = init_domain_decomposition(mdlog, cr, domdecOptions, mdrunOptions,
1036 &mtop, inputrec,
1037 box, positionsFromStatePointer(globalState.get()),
1038 &atomSets);
1039 // Note that local state still does not exist yet.
1041 else
1043 /* PME, if used, is done on all nodes with 1D decomposition */
1044 cr->npmenodes = 0;
1045 cr->duty = (DUTY_PP | DUTY_PME);
1047 if (inputrec->ePBC == epbcSCREW)
1049 gmx_fatal(FARGS,
1050 "pbc=screw is only implemented with domain decomposition");
1054 if (PAR(cr))
1056 /* After possible communicator splitting in make_dd_communicators.
1057 * we can set up the intra/inter node communication.
1059 gmx_setup_nodecomm(fplog, cr);
1062 #if GMX_MPI
1063 if (isMultiSim(ms))
1065 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1066 "This is simulation %d out of %d running as a composite GROMACS\n"
1067 "multi-simulation job. Setup for this simulation:\n",
1068 ms->sim, ms->nsim);
1070 GMX_LOG(mdlog.warning).appendTextFormatted(
1071 "Using %d MPI %s\n",
1072 cr->nnodes,
1073 #if GMX_THREAD_MPI
1074 cr->nnodes == 1 ? "thread" : "threads"
1075 #else
1076 cr->nnodes == 1 ? "process" : "processes"
1077 #endif
1079 fflush(stderr);
1080 #endif
1082 // If mdrun -pin auto honors any affinity setting that already
1083 // exists. If so, it is nice to provide feedback about whether
1084 // that existing affinity setting was from OpenMP or something
1085 // else, so we run this code both before and after we initialize
1086 // the OpenMP support.
1087 gmx_check_thread_affinity_set(mdlog,
1088 &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1089 /* Check and update the number of OpenMP threads requested */
1090 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
1091 pmeRunMode, mtop, *inputrec);
1093 gmx_omp_nthreads_init(mdlog, cr,
1094 hwinfo->nthreads_hw_avail,
1095 physicalNodeComm.size_,
1096 hw_opt.nthreads_omp,
1097 hw_opt.nthreads_omp_pme,
1098 !thisRankHasDuty(cr, DUTY_PP));
1100 // Enable FP exception detection, but not in
1101 // Release mode and not for compilers with known buggy FP
1102 // exception support (clang with any optimization) or suspected
1103 // buggy FP exception support (gcc 7.* with optimization).
1104 #if !defined NDEBUG && \
1105 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1106 && defined __OPTIMIZE__)
1107 const bool bEnableFPE = true;
1108 #else
1109 const bool bEnableFPE = false;
1110 #endif
1111 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1112 if (bEnableFPE)
1114 gmx_feenableexcept();
1117 // Build a data structure that expresses which kinds of non-bonded
1118 // task are handled by this rank.
1120 // TODO Later, this might become a loop over all registered modules
1121 // relevant to the mdp inputs, to find those that have such tasks.
1123 // TODO This could move before init_domain_decomposition() as part
1124 // of refactoring that separates the responsibility for duty
1125 // assignment from setup for communication between tasks, and
1126 // setup for tasks handled with a domain (ie including short-ranged
1127 // tasks, bonded tasks, etc.).
1129 // Note that in general useGpuForNonbonded, etc. can have a value
1130 // that is inconsistent with the presence of actual GPUs on any
1131 // rank, and that is not known to be a problem until the
1132 // duty of the ranks on a node become known.
1134 // TODO Later we might need the concept of computeTasksOnThisRank,
1135 // from which we construct gpuTasksOnThisRank.
1137 // Currently the DD code assigns duty to ranks that can
1138 // include PP work that currently can be executed on a single
1139 // GPU, if present and compatible. This has to be coordinated
1140 // across PP ranks on a node, with possible multiple devices
1141 // or sharing devices on a node, either from the user
1142 // selection, or automatically.
1143 auto haveGpus = !gpuIdsToUse.empty();
1144 std::vector<GpuTask> gpuTasksOnThisRank;
1145 if (thisRankHasDuty(cr, DUTY_PP))
1147 if (useGpuForNonbonded)
1149 // Note that any bonded tasks on a GPU always accompany a
1150 // non-bonded task.
1151 if (haveGpus)
1153 gpuTasksOnThisRank.push_back(GpuTask::Nonbonded);
1155 else if (nonbondedTarget == TaskTarget::Gpu)
1157 gmx_fatal(FARGS, "Cannot run short-ranged nonbonded interactions on a GPU because no GPU is detected.");
1159 else if (bondedTarget == TaskTarget::Gpu)
1161 gmx_fatal(FARGS, "Cannot run bonded interactions on a GPU because no GPU is detected.");
1165 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1166 if (EEL_PME(inputrec->coulombtype) && (thisRankHasDuty(cr, DUTY_PME)))
1168 if (useGpuForPme)
1170 if (haveGpus)
1172 gpuTasksOnThisRank.push_back(GpuTask::Pme);
1174 else if (pmeTarget == TaskTarget::Gpu)
1176 gmx_fatal(FARGS, "Cannot run PME on a GPU because no GPU is detected.");
1181 GpuTaskAssignment gpuTaskAssignment;
1184 // Produce the task assignment for this rank.
1185 gpuTaskAssignment = runTaskAssignment(gpuIdsToUse, userGpuTaskAssignment, *hwinfo,
1186 mdlog, cr, ms, physicalNodeComm, gpuTasksOnThisRank,
1187 useGpuForBonded, pmeRunMode);
1189 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1191 /* Prevent other ranks from continuing after an issue was found
1192 * and reported as a fatal error.
1194 * TODO This function implements a barrier so that MPI runtimes
1195 * can organize an orderly shutdown if one of the ranks has had to
1196 * issue a fatal error in various code already run. When we have
1197 * MPI-aware error handling and reporting, this should be
1198 * improved. */
1199 #if GMX_MPI
1200 if (PAR(cr))
1202 MPI_Barrier(cr->mpi_comm_mysim);
1204 if (isMultiSim(ms))
1206 if (SIMMASTER(cr))
1208 MPI_Barrier(ms->mpi_comm_masters);
1210 /* We need another barrier to prevent non-master ranks from contiuing
1211 * when an error occured in a different simulation.
1213 MPI_Barrier(cr->mpi_comm_mysim);
1215 #endif
1217 /* Now that we know the setup is consistent, check for efficiency */
1218 check_resource_division_efficiency(hwinfo, !gpuTaskAssignment.empty(), mdrunOptions.ntompOptionIsSet,
1219 cr, mdlog);
1221 gmx_device_info_t *nonbondedDeviceInfo = nullptr;
1223 if (thisRankHasDuty(cr, DUTY_PP))
1225 // This works because only one task of each type is currently permitted.
1226 auto nbGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(),
1227 hasTaskType<GpuTask::Nonbonded>);
1228 if (nbGpuTaskMapping != gpuTaskAssignment.end())
1230 int nonbondedDeviceId = nbGpuTaskMapping->deviceId_;
1231 nonbondedDeviceInfo = getDeviceInfo(hwinfo->gpu_info, nonbondedDeviceId);
1232 init_gpu(nonbondedDeviceInfo);
1234 if (DOMAINDECOMP(cr))
1236 /* When we share GPUs over ranks, we need to know this for the DLB */
1237 dd_setup_dlb_resource_sharing(cr, nonbondedDeviceId);
1243 gmx_device_info_t *pmeDeviceInfo = nullptr;
1244 // Later, this program could contain kernels that might be later
1245 // re-used as auto-tuning progresses, or subsequent simulations
1246 // are invoked.
1247 PmeGpuProgramStorage pmeGpuProgram;
1248 // This works because only one task of each type is currently permitted.
1249 auto pmeGpuTaskMapping = std::find_if(gpuTaskAssignment.begin(), gpuTaskAssignment.end(), hasTaskType<GpuTask::Pme>);
1250 const bool thisRankHasPmeGpuTask = (pmeGpuTaskMapping != gpuTaskAssignment.end());
1251 if (thisRankHasPmeGpuTask)
1253 pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
1254 init_gpu(pmeDeviceInfo);
1255 pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
1258 /* getting number of PP/PME threads on this MPI / tMPI rank.
1259 PME: env variable should be read only on one node to make sure it is
1260 identical everywhere;
1262 const int numThreadsOnThisRank =
1263 thisRankHasDuty(cr, DUTY_PP) ? gmx_omp_nthreads_get(emntNonbonded) : gmx_omp_nthreads_get(emntPME);
1264 checkHardwareOversubscription(numThreadsOnThisRank, cr->nodeid,
1265 *hwinfo->hardwareTopology,
1266 physicalNodeComm, mdlog);
1268 if (hw_opt.threadAffinity != ThreadAffinity::Off)
1270 /* Before setting affinity, check whether the affinity has changed
1271 * - which indicates that probably the OpenMP library has changed it
1272 * since we first checked).
1274 gmx_check_thread_affinity_set(mdlog,
1275 &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1277 int numThreadsOnThisNode, intraNodeThreadOffset;
1278 analyzeThreadsOnThisNode(physicalNodeComm, numThreadsOnThisRank, &numThreadsOnThisNode,
1279 &intraNodeThreadOffset);
1281 /* Set the CPU affinity */
1282 gmx_set_thread_affinity(mdlog, cr, &hw_opt, *hwinfo->hardwareTopology,
1283 numThreadsOnThisRank, numThreadsOnThisNode,
1284 intraNodeThreadOffset, nullptr);
1287 if (mdrunOptions.timingOptions.resetStep > -1)
1289 GMX_LOG(mdlog.info).asParagraph().
1290 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1292 wcycle = wallcycle_init(fplog, mdrunOptions.timingOptions.resetStep, cr);
1294 if (PAR(cr))
1296 /* Master synchronizes its value of reset_counters with all nodes
1297 * including PME only nodes */
1298 int64_t reset_counters = wcycle_get_reset_counters(wcycle);
1299 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1300 wcycle_set_reset_counters(wcycle, reset_counters);
1303 // Membrane embedding must be initialized before we call init_forcerec()
1304 if (doMembed)
1306 if (MASTER(cr))
1308 fprintf(stderr, "Initializing membed");
1310 /* Note that membed cannot work in parallel because mtop is
1311 * changed here. Fix this if we ever want to make it run with
1312 * multiple ranks. */
1313 membed = init_membed(fplog, filenames.size(), filenames.data(), &mtop, inputrec, globalState.get(), cr,
1314 &mdrunOptions
1315 .checkpointOptions.period);
1318 std::unique_ptr<MDAtoms> mdAtoms;
1319 std::unique_ptr<gmx_vsite_t> vsite;
1321 t_nrnb nrnb;
1322 if (thisRankHasDuty(cr, DUTY_PP))
1324 mdModulesNotifier.notify(*cr);
1325 mdModulesNotifier.notify(&atomSets);
1326 mdModulesNotifier.notify(PeriodicBoundaryConditionType {inputrec->ePBC});
1327 /* Initiate forcerecord */
1328 fr = new t_forcerec;
1329 fr->forceProviders = mdModules_->initForceProviders();
1330 init_forcerec(fplog, mdlog, fr, fcd,
1331 inputrec, &mtop, cr, box,
1332 opt2fn("-table", filenames.size(), filenames.data()),
1333 opt2fn("-tablep", filenames.size(), filenames.data()),
1334 opt2fns("-tableb", filenames.size(), filenames.data()),
1335 *hwinfo, nonbondedDeviceInfo,
1336 useGpuForBonded,
1337 pforce,
1338 wcycle);
1340 // TODO Move this to happen during domain decomposition setup,
1341 // once stream and event handling works well with that.
1342 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
1343 if (havePPDomainDecomposition(cr) && c_enableGpuHaloExchange && useGpuForNonbonded)
1345 void *streamLocal = Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
1346 void *streamNonLocal =
1347 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
1348 void *coordinatesOnDeviceEvent = fr->nbv->get_x_on_device_event();
1349 cr->dd->gpuHaloExchange = std::make_unique<GpuHaloExchange>(cr->dd, cr->mpi_comm_mysim, streamLocal,
1350 streamNonLocal, coordinatesOnDeviceEvent);
1353 /* Initialize the mdAtoms structure.
1354 * mdAtoms is not filled with atom data,
1355 * as this can not be done now with domain decomposition.
1357 mdAtoms = makeMDAtoms(fplog, mtop, *inputrec, thisRankHasPmeGpuTask);
1358 if (globalState && thisRankHasPmeGpuTask)
1360 // The pinning of coordinates in the global state object works, because we only use
1361 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1362 // points to the global state object without DD.
1363 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1364 // which should also perform the pinning.
1365 changePinningPolicy(&globalState->x, pme_get_pinning_policy());
1368 /* Initialize the virtual site communication */
1369 vsite = initVsite(mtop, cr);
1371 calc_shifts(box, fr->shift_vec);
1373 /* With periodic molecules the charge groups should be whole at start up
1374 * and the virtual sites should not be far from their proper positions.
1376 if (!inputrec->bContinuation && MASTER(cr) &&
1377 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1379 /* Make molecules whole at start of run */
1380 if (fr->ePBC != epbcNONE)
1382 do_pbc_first_mtop(fplog, inputrec->ePBC, box, &mtop, globalState->x.rvec_array());
1384 if (vsite)
1386 /* Correct initial vsite positions are required
1387 * for the initial distribution in the domain decomposition
1388 * and for the initial shell prediction.
1390 constructVsitesGlobal(mtop, globalState->x);
1394 if (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
1396 ewaldcoeff_q = fr->ic->ewaldcoeff_q;
1397 ewaldcoeff_lj = fr->ic->ewaldcoeff_lj;
1400 else
1402 /* This is a PME only node */
1404 GMX_ASSERT(globalState == nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1406 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1407 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1410 gmx_pme_t *sepPmeData = nullptr;
1411 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1412 GMX_ASSERT(thisRankHasDuty(cr, DUTY_PP) == (fr != nullptr), "Double-checking that only PME-only ranks have no forcerec");
1413 gmx_pme_t * &pmedata = fr ? fr->pmedata : sepPmeData;
1415 /* Initiate PME if necessary,
1416 * either on all nodes or on dedicated PME nodes only. */
1417 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1419 if (mdAtoms && mdAtoms->mdatoms())
1421 nChargePerturbed = mdAtoms->mdatoms()->nChargePerturbed;
1422 if (EVDW_PME(inputrec->vdwtype))
1424 nTypePerturbed = mdAtoms->mdatoms()->nTypePerturbed;
1427 if (cr->npmenodes > 0)
1429 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1430 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1431 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1434 if (thisRankHasDuty(cr, DUTY_PME))
1438 pmedata = gmx_pme_init(cr,
1439 getNumPmeDomains(cr->dd),
1440 inputrec,
1441 nChargePerturbed != 0, nTypePerturbed != 0,
1442 mdrunOptions.reproducible,
1443 ewaldcoeff_q, ewaldcoeff_lj,
1444 gmx_omp_nthreads_get(emntPME),
1445 pmeRunMode, nullptr,
1446 pmeDeviceInfo, pmeGpuProgram.get(), mdlog);
1448 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1453 if (EI_DYNAMICS(inputrec->eI))
1455 /* Turn on signal handling on all nodes */
1457 * (A user signal from the PME nodes (if any)
1458 * is communicated to the PP nodes.
1460 signal_handler_install();
1463 pull_t *pull_work = nullptr;
1464 if (thisRankHasDuty(cr, DUTY_PP))
1466 /* Assumes uniform use of the number of OpenMP threads */
1467 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1469 if (inputrec->bPull)
1471 /* Initialize pull code */
1472 pull_work =
1473 init_pull(fplog, inputrec->pull, inputrec,
1474 &mtop, cr, &atomSets, inputrec->fepvals->init_lambda);
1475 if (inputrec->pull->bXOutAverage || inputrec->pull->bFOutAverage)
1477 initPullHistory(pull_work, &observablesHistory);
1479 if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
1481 init_pull_output_files(pull_work,
1482 filenames.size(), filenames.data(), oenv,
1483 startingBehavior);
1487 std::unique_ptr<EnforcedRotation> enforcedRotation;
1488 if (inputrec->bRot)
1490 /* Initialize enforced rotation code */
1491 enforcedRotation = init_rot(fplog,
1492 inputrec,
1493 filenames.size(),
1494 filenames.data(),
1496 &atomSets,
1497 globalState.get(),
1498 &mtop,
1499 oenv,
1500 mdrunOptions,
1501 startingBehavior);
1504 t_swap *swap = nullptr;
1505 if (inputrec->eSwapCoords != eswapNO)
1507 /* Initialize ion swapping code */
1508 swap = init_swapcoords(fplog, inputrec,
1509 opt2fn_master("-swap", filenames.size(), filenames.data(), cr),
1510 &mtop, globalState.get(), &observablesHistory,
1511 cr, &atomSets, oenv, mdrunOptions,
1512 startingBehavior);
1515 /* Let makeConstraints know whether we have essential dynamics constraints.
1516 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1518 bool doEssentialDynamics = (opt2fn_null("-ei", filenames.size(), filenames.data()) != nullptr
1519 || observablesHistory.edsamHistory);
1520 auto constr = makeConstraints(mtop, *inputrec, pull_work, doEssentialDynamics,
1521 fplog, *mdAtoms->mdatoms(),
1522 cr, ms, &nrnb, wcycle, fr->bMolPBC);
1524 /* Energy terms and groups */
1525 gmx_enerdata_t enerd(mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(), inputrec->fepvals->n_lambda);
1527 /* Kinetic energy data */
1528 gmx_ekindata_t ekind;
1529 init_ekindata(fplog, &mtop, &(inputrec->opts), &ekind, inputrec->cos_accel);
1531 /* Set up interactive MD (IMD) */
1532 auto imdSession = makeImdSession(inputrec, cr, wcycle, &enerd, ms, &mtop, mdlog,
1533 MASTER(cr) ? globalState->x.rvec_array() : nullptr,
1534 filenames.size(), filenames.data(), oenv, mdrunOptions.imdOptions,
1535 startingBehavior);
1537 if (DOMAINDECOMP(cr))
1539 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1540 /* This call is not included in init_domain_decomposition mainly
1541 * because fr->cginfo_mb is set later.
1543 dd_init_bondeds(fplog, cr->dd, &mtop, vsite.get(), inputrec,
1544 domdecOptions.checkBondedInteractions,
1545 fr->cginfo_mb);
1548 if (updateTarget == TaskTarget::Gpu)
1550 if (SIMMASTER(cr))
1552 gmx_fatal(FARGS, "It is currently not possible to redirect the calculation "
1553 "of update and constraints to the GPU!");
1557 // Before we start the actual simulator, try if we can run the update task on the GPU.
1558 useGpuForUpdate = decideWhetherToUseGpuForUpdate(DOMAINDECOMP(cr),
1559 useGpuForNonbonded,
1560 updateTarget,
1561 gpusWereDetected,
1562 *inputrec,
1563 *mdAtoms,
1564 doEssentialDynamics,
1565 fcd->orires.nr != 0,
1566 fcd->disres.nsystems != 0);
1568 // TODO This is not the right place to manage the lifetime of
1569 // this data structure, but currently it's the easiest way to
1570 // make it work.
1571 MdrunScheduleWorkload runScheduleWork;
1573 GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
1574 SimulatorBuilder simulatorBuilder;
1576 // build and run simulator object based on user-input
1577 const bool inputIsCompatibleWithModularSimulator = ModularSimulator::isInputCompatible(
1578 false,
1579 inputrec, doRerun, vsite.get(), ms, replExParams,
1580 fcd, static_cast<int>(filenames.size()), filenames.data(),
1581 &observablesHistory, membed);
1582 auto simulator = simulatorBuilder.build(
1583 inputIsCompatibleWithModularSimulator,
1584 fplog, cr, ms, mdlog, static_cast<int>(filenames.size()), filenames.data(),
1585 oenv,
1586 mdrunOptions,
1587 startingBehavior,
1588 vsite.get(), constr.get(),
1589 enforcedRotation ? enforcedRotation->getLegacyEnfrot() : nullptr,
1590 deform.get(),
1591 mdModules_->outputProvider(),
1592 mdModules_->notifier(),
1593 inputrec, imdSession.get(), pull_work, swap, &mtop,
1594 fcd,
1595 globalState.get(),
1596 &observablesHistory,
1597 mdAtoms.get(), &nrnb, wcycle, fr,
1598 &enerd,
1599 &ekind,
1600 &runScheduleWork,
1601 replExParams,
1602 membed,
1603 walltime_accounting,
1604 std::move(stopHandlerBuilder_),
1605 doRerun,
1606 useGpuForUpdate);
1607 simulator->run();
1609 if (inputrec->bPull)
1611 finish_pull(pull_work);
1613 finish_swapcoords(swap);
1615 else
1617 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1618 /* do PME only */
1619 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1620 gmx_pmeonly(pmedata, cr, &nrnb, wcycle, walltime_accounting, inputrec, pmeRunMode);
1623 wallcycle_stop(wcycle, ewcRUN);
1625 /* Finish up, write some stuff
1626 * if rerunMD, don't write last frame again
1628 finish_run(fplog, mdlog, cr,
1629 inputrec, &nrnb, wcycle, walltime_accounting,
1630 fr ? fr->nbv.get() : nullptr,
1631 pmedata,
1632 EI_DYNAMICS(inputrec->eI) && !isMultiSim(ms));
1634 // clean up cycle counter
1635 wallcycle_destroy(wcycle);
1637 // Free PME data
1638 if (pmedata)
1640 gmx_pme_destroy(pmedata);
1641 pmedata = nullptr;
1644 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1645 // before we destroy the GPU context(s) in free_gpu_resources().
1646 // Pinned buffers are associated with contexts in CUDA.
1647 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1648 mdAtoms.reset(nullptr);
1649 globalState.reset(nullptr);
1650 mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
1652 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1653 free_gpu_resources(fr, physicalNodeComm, hwinfo->gpu_info);
1654 free_gpu(nonbondedDeviceInfo);
1655 free_gpu(pmeDeviceInfo);
1656 done_forcerec(fr, mtop.molblock.size());
1657 sfree(fcd);
1659 if (doMembed)
1661 free_membed(membed);
1664 /* Does what it says */
1665 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1666 walltime_accounting_destroy(walltime_accounting);
1668 // Ensure log file content is written
1669 if (logFileHandle)
1671 gmx_fio_flush(logFileHandle);
1674 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1675 * exceptions were enabled before function was called. */
1676 if (bEnableFPE)
1678 gmx_fedisableexcept();
1681 auto rc = static_cast<int>(gmx_get_stop_condition());
1683 #if GMX_THREAD_MPI
1684 /* we need to join all threads. The sub-threads join when they
1685 exit this function, but the master thread needs to be told to
1686 wait for that. */
1687 if (PAR(cr) && MASTER(cr))
1689 tMPI_Finalize();
1691 #endif
1692 return rc;
1695 Mdrunner::~Mdrunner()
1697 // Clean up of the Manager.
1698 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1699 // but okay as long as threads synchronize some time before adding or accessing
1700 // a new set of restraints.
1701 if (restraintManager_)
1703 restraintManager_->clear();
1704 GMX_ASSERT(restraintManager_->countRestraints() == 0,
1705 "restraints added during runner life time should be cleared at runner destruction.");
1709 void Mdrunner::addPotential(std::shared_ptr<gmx::IRestraintPotential> puller,
1710 const std::string &name)
1712 GMX_ASSERT(restraintManager_, "Mdrunner must have a restraint manager.");
1713 // Not sure if this should be logged through the md logger or something else,
1714 // but it is helpful to have some sort of INFO level message sent somewhere.
1715 // std::cout << "Registering restraint named " << name << std::endl;
1717 // When multiple restraints are used, it may be wasteful to register them separately.
1718 // Maybe instead register an entire Restraint Manager as a force provider.
1719 restraintManager_->addToSpec(std::move(puller),
1720 name);
1723 Mdrunner::Mdrunner(std::unique_ptr<MDModules> mdModules)
1724 : mdModules_(std::move(mdModules))
1728 Mdrunner::Mdrunner(Mdrunner &&) noexcept = default;
1730 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1731 Mdrunner &Mdrunner::operator=(Mdrunner && /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING) = default;
1733 class Mdrunner::BuilderImplementation
1735 public:
1736 BuilderImplementation() = delete;
1737 BuilderImplementation(std::unique_ptr<MDModules> mdModules);
1738 ~BuilderImplementation();
1740 BuilderImplementation &setExtraMdrunOptions(const MdrunOptions &options,
1741 real forceWarningThreshold,
1742 StartingBehavior startingBehavior);
1744 void addDomdec(const DomdecOptions &options);
1746 void addVerletList(int nstlist);
1748 void addReplicaExchange(const ReplicaExchangeParameters &params);
1750 void addMultiSim(gmx_multisim_t* multisim);
1752 void addCommunicationRecord(t_commrec *cr);
1754 void addNonBonded(const char* nbpu_opt);
1756 void addPME(const char* pme_opt_, const char* pme_fft_opt_);
1758 void addBondedTaskAssignment(const char* bonded_opt);
1760 void addUpdateTaskAssignment(const char* update_opt);
1762 void addHardwareOptions(const gmx_hw_opt_t &hardwareOptions);
1764 void addFilenames(ArrayRef <const t_filenm> filenames);
1766 void addOutputEnvironment(gmx_output_env_t* outputEnvironment);
1768 void addLogFile(t_fileio *logFileHandle);
1770 void addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder);
1772 Mdrunner build();
1774 private:
1776 // Default parameters copied from runner.h
1777 // \todo Clarify source(s) of default parameters.
1779 const char* nbpu_opt_ = nullptr;
1780 const char* pme_opt_ = nullptr;
1781 const char* pme_fft_opt_ = nullptr;
1782 const char *bonded_opt_ = nullptr;
1783 const char *update_opt_ = nullptr;
1785 MdrunOptions mdrunOptions_;
1787 DomdecOptions domdecOptions_;
1789 ReplicaExchangeParameters replicaExchangeParameters_;
1791 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1792 int nstlist_ = 0;
1794 //! Multisim communicator handle.
1795 std::unique_ptr<gmx_multisim_t*> multisim_ = nullptr;
1797 //! Non-owning communication record (only used when building command-line mdrun)
1798 t_commrec *communicationRecord_ = nullptr;
1800 //! Print a warning if any force is larger than this (in kJ/mol nm).
1801 real forceWarningThreshold_ = -1;
1803 //! Whether the simulation will start afresh, or restart with/without appending.
1804 StartingBehavior startingBehavior_ = StartingBehavior::NewSimulation;
1806 //! The modules that comprise the functionality of mdrun.
1807 std::unique_ptr<MDModules> mdModules_;
1809 //! \brief Parallelism information.
1810 gmx_hw_opt_t hardwareOptions_;
1812 //! filename options for simulation.
1813 ArrayRef<const t_filenm> filenames_;
1815 /*! \brief Handle to output environment.
1817 * \todo gmx_output_env_t needs lifetime management.
1819 gmx_output_env_t* outputEnvironment_ = nullptr;
1821 /*! \brief Non-owning handle to MD log file.
1823 * \todo Context should own output facilities for client.
1824 * \todo Improve log file handle management.
1825 * \internal
1826 * Code managing the FILE* relies on the ability to set it to
1827 * nullptr to check whether the filehandle is valid.
1829 t_fileio* logFileHandle_ = nullptr;
1832 * \brief Builder for simulation stop signal handler.
1834 std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_ = nullptr;
1837 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules> mdModules) :
1838 mdModules_(std::move(mdModules))
1842 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1844 Mdrunner::BuilderImplementation &
1845 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions &options,
1846 const real forceWarningThreshold,
1847 const StartingBehavior startingBehavior)
1849 mdrunOptions_ = options;
1850 forceWarningThreshold_ = forceWarningThreshold;
1851 startingBehavior_ = startingBehavior;
1852 return *this;
1855 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions &options)
1857 domdecOptions_ = options;
1860 void Mdrunner::BuilderImplementation::addVerletList(int nstlist)
1862 nstlist_ = nstlist;
1865 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters &params)
1867 replicaExchangeParameters_ = params;
1870 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t* multisim)
1872 multisim_ = std::make_unique<gmx_multisim_t*>(multisim);
1875 void Mdrunner::BuilderImplementation::addCommunicationRecord(t_commrec *cr)
1877 communicationRecord_ = cr;
1880 Mdrunner Mdrunner::BuilderImplementation::build()
1882 auto newRunner = Mdrunner(std::move(mdModules_));
1884 newRunner.mdrunOptions = mdrunOptions_;
1885 newRunner.pforce = forceWarningThreshold_;
1886 newRunner.startingBehavior = startingBehavior_;
1887 newRunner.domdecOptions = domdecOptions_;
1889 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1890 newRunner.hw_opt = hardwareOptions_;
1892 // No invariant to check. This parameter exists to optionally override other behavior.
1893 newRunner.nstlist_cmdline = nstlist_;
1895 newRunner.replExParams = replicaExchangeParameters_;
1897 newRunner.filenames = filenames_;
1899 GMX_ASSERT(communicationRecord_, "Bug found. It should not be possible to call build() without a valid communicationRecord_.");
1900 newRunner.cr = communicationRecord_;
1902 if (multisim_)
1904 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1905 newRunner.ms = *multisim_;
1907 else
1909 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1912 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1913 // \todo Update sanity checking when output environment has clearly specified invariants.
1914 // Initialization and default values for oenv are not well specified in the current version.
1915 if (outputEnvironment_)
1917 newRunner.oenv = outputEnvironment_;
1919 else
1921 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1924 newRunner.logFileHandle = logFileHandle_;
1926 if (nbpu_opt_)
1928 newRunner.nbpu_opt = nbpu_opt_;
1930 else
1932 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1935 if (pme_opt_ && pme_fft_opt_)
1937 newRunner.pme_opt = pme_opt_;
1938 newRunner.pme_fft_opt = pme_fft_opt_;
1940 else
1942 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1945 if (bonded_opt_)
1947 newRunner.bonded_opt = bonded_opt_;
1949 else
1951 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1954 if (update_opt_)
1956 newRunner.update_opt = update_opt_;
1958 else
1960 GMX_THROW(gmx::APIError("MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1964 newRunner.restraintManager_ = std::make_unique<gmx::RestraintManager>();
1966 if (stopHandlerBuilder_)
1968 newRunner.stopHandlerBuilder_ = std::move(stopHandlerBuilder_);
1970 else
1972 newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>();
1975 return newRunner;
1978 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt)
1980 nbpu_opt_ = nbpu_opt;
1983 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt,
1984 const char* pme_fft_opt)
1986 pme_opt_ = pme_opt;
1987 pme_fft_opt_ = pme_fft_opt;
1990 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt)
1992 bonded_opt_ = bonded_opt;
1995 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt)
1997 update_opt_ = update_opt;
2000 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2002 hardwareOptions_ = hardwareOptions;
2005 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef<const t_filenm> filenames)
2007 filenames_ = filenames;
2010 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2012 outputEnvironment_ = outputEnvironment;
2015 void Mdrunner::BuilderImplementation::addLogFile(t_fileio *logFileHandle)
2017 logFileHandle_ = logFileHandle;
2020 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2022 stopHandlerBuilder_ = std::move(builder);
2025 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr<MDModules> mdModules) :
2026 impl_ {std::make_unique<Mdrunner::BuilderImplementation>(std::move(mdModules))}
2030 MdrunnerBuilder::~MdrunnerBuilder() = default;
2032 MdrunnerBuilder &MdrunnerBuilder::addSimulationMethod(const MdrunOptions &options,
2033 real forceWarningThreshold,
2034 const StartingBehavior startingBehavior)
2036 impl_->setExtraMdrunOptions(options, forceWarningThreshold, startingBehavior);
2037 return *this;
2040 MdrunnerBuilder &MdrunnerBuilder::addDomainDecomposition(const DomdecOptions &options)
2042 impl_->addDomdec(options);
2043 return *this;
2046 MdrunnerBuilder &MdrunnerBuilder::addNeighborList(int nstlist)
2048 impl_->addVerletList(nstlist);
2049 return *this;
2052 MdrunnerBuilder &MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters &params)
2054 impl_->addReplicaExchange(params);
2055 return *this;
2058 MdrunnerBuilder &MdrunnerBuilder::addMultiSim(gmx_multisim_t* multisim)
2060 impl_->addMultiSim(multisim);
2061 return *this;
2064 MdrunnerBuilder &MdrunnerBuilder::addCommunicationRecord(t_commrec *cr)
2066 impl_->addCommunicationRecord(cr);
2067 return *this;
2070 MdrunnerBuilder &MdrunnerBuilder::addNonBonded(const char* nbpu_opt)
2072 impl_->addNonBonded(nbpu_opt);
2073 return *this;
2076 MdrunnerBuilder &MdrunnerBuilder::addElectrostatics(const char* pme_opt,
2077 const char* pme_fft_opt)
2079 // The builder method may become more general in the future, but in this version,
2080 // parameters for PME electrostatics are both required and the only parameters
2081 // available.
2082 if (pme_opt && pme_fft_opt)
2084 impl_->addPME(pme_opt, pme_fft_opt);
2086 else
2088 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2090 return *this;
2093 MdrunnerBuilder &MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt)
2095 impl_->addBondedTaskAssignment(bonded_opt);
2096 return *this;
2099 MdrunnerBuilder &MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt)
2101 impl_->addUpdateTaskAssignment(update_opt);
2102 return *this;
2105 Mdrunner MdrunnerBuilder::build()
2107 return impl_->build();
2110 MdrunnerBuilder &MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t &hardwareOptions)
2112 impl_->addHardwareOptions(hardwareOptions);
2113 return *this;
2116 MdrunnerBuilder &MdrunnerBuilder::addFilenames(ArrayRef<const t_filenm> filenames)
2118 impl_->addFilenames(filenames);
2119 return *this;
2122 MdrunnerBuilder &MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t* outputEnvironment)
2124 impl_->addOutputEnvironment(outputEnvironment);
2125 return *this;
2128 MdrunnerBuilder &MdrunnerBuilder::addLogFile(t_fileio *logFileHandle)
2130 impl_->addLogFile(logFileHandle);
2131 return *this;
2134 MdrunnerBuilder &MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr<StopHandlerBuilder> builder)
2136 impl_->addStopHandlerBuilder(std::move(builder));
2137 return *this;
2140 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder &&) noexcept = default;
2142 MdrunnerBuilder &MdrunnerBuilder::operator=(MdrunnerBuilder &&) noexcept = default;
2144 } // namespace gmx