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.
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
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/localatomsetmanager.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/ewald_utils.h"
65 #include "gromacs/ewald/pme.h"
66 #include "gromacs/ewald/pme_gpu_program.h"
67 #include "gromacs/fileio/checkpoint.h"
68 #include "gromacs/fileio/gmxfio.h"
69 #include "gromacs/fileio/oenv.h"
70 #include "gromacs/fileio/tpxio.h"
71 #include "gromacs/gmxlib/network.h"
72 #include "gromacs/gmxlib/nrnb.h"
73 #include "gromacs/gpu_utils/clfftinitializer.h"
74 #include "gromacs/gpu_utils/gpu_utils.h"
75 #include "gromacs/hardware/cpuinfo.h"
76 #include "gromacs/hardware/detecthardware.h"
77 #include "gromacs/hardware/printhardware.h"
78 #include "gromacs/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/ppforceworkload.h"
98 #include "gromacs/mdlib/qmmm.h"
99 #include "gromacs/mdlib/sighandler.h"
100 #include "gromacs/mdlib/stophandler.h"
101 #include "gromacs/mdrun/mdmodules.h"
102 #include "gromacs/mdrun/simulationcontext.h"
103 #include "gromacs/mdrunutility/handlerestart.h"
104 #include "gromacs/mdrunutility/logging.h"
105 #include "gromacs/mdrunutility/multisim.h"
106 #include "gromacs/mdrunutility/printtime.h"
107 #include "gromacs/mdrunutility/threadaffinity.h"
108 #include "gromacs/mdtypes/commrec.h"
109 #include "gromacs/mdtypes/enerdata.h"
110 #include "gromacs/mdtypes/fcdata.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/state.h"
116 #include "gromacs/nbnxm/gpu_data_mgmt.h"
117 #include "gromacs/nbnxm/nbnxm.h"
118 #include "gromacs/nbnxm/pairlist_tuning.h"
119 #include "gromacs/pbcutil/pbc.h"
120 #include "gromacs/pulling/output.h"
121 #include "gromacs/pulling/pull.h"
122 #include "gromacs/pulling/pull_rotation.h"
123 #include "gromacs/restraint/manager.h"
124 #include "gromacs/restraint/restraintmdmodule.h"
125 #include "gromacs/restraint/restraintpotential.h"
126 #include "gromacs/swap/swapcoords.h"
127 #include "gromacs/taskassignment/decidegpuusage.h"
128 #include "gromacs/taskassignment/resourcedivision.h"
129 #include "gromacs/taskassignment/taskassignment.h"
130 #include "gromacs/taskassignment/usergpuids.h"
131 #include "gromacs/timing/gpu_timing.h"
132 #include "gromacs/timing/wallcycle.h"
133 #include "gromacs/timing/wallcyclereporting.h"
134 #include "gromacs/topology/mtop_util.h"
135 #include "gromacs/trajectory/trajectoryframe.h"
136 #include "gromacs/utility/basenetwork.h"
137 #include "gromacs/utility/cstringutil.h"
138 #include "gromacs/utility/exceptions.h"
139 #include "gromacs/utility/fatalerror.h"
140 #include "gromacs/utility/filestream.h"
141 #include "gromacs/utility/gmxassert.h"
142 #include "gromacs/utility/gmxmpi.h"
143 #include "gromacs/utility/logger.h"
144 #include "gromacs/utility/loggerbuilder.h"
145 #include "gromacs/utility/physicalnodecommunicator.h"
146 #include "gromacs/utility/pleasecite.h"
147 #include "gromacs/utility/programcontext.h"
148 #include "gromacs/utility/smalloc.h"
149 #include "gromacs/utility/stringutil.h"
151 #include "isimulator.h"
152 #include "replicaexchange.h"
153 #include "simulatorbuilder.h"
156 #include "corewrap.h"
162 /*! \brief Log if development feature flags are encountered
164 * The use of dev features indicated by environment variables is logged
165 * in order to ensure that runs with such featrues enabled can be identified
166 * from their log and standard output.
168 * \param[in] mdlog Logger object.
170 static void reportDevelopmentFeatures(const gmx::MDLogger
&mdlog
)
172 const bool enableGpuBufOps
= (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
173 const bool useGpuUpdateConstrain
= (getenv("GMX_UPDATE_CONSTRAIN_GPU") != nullptr);
177 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
178 "NOTE: This run uses the 'GPU buffer ops' feature, enabled by the GMX_USE_GPU_BUFFER_OPS environment variable.");
181 if (useGpuUpdateConstrain
)
183 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
184 "NOTE: This run uses the 'GPU update/constraints' feature, enabled by the GMX_UPDATE_CONSTRAIN_GPU environment variable.");
188 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
190 * Used to ensure that the master thread does not modify mdrunner during copy
191 * on the spawned threads. */
192 static void threadMpiMdrunnerAccessBarrier()
195 MPI_Barrier(MPI_COMM_WORLD
);
199 Mdrunner
Mdrunner::cloneOnSpawnedThread() const
201 auto newRunner
= Mdrunner(std::make_unique
<MDModules
>());
203 // All runners in the same process share a restraint manager resource because it is
204 // part of the interface to the client code, which is associated only with the
205 // original thread. Handles to the same resources can be obtained by copy.
207 newRunner
.restraintManager_
= std::make_unique
<RestraintManager
>(*restraintManager_
);
210 // Copy original cr pointer before master thread can pass the thread barrier
211 newRunner
.cr
= reinitialize_commrec_for_this_thread(cr
);
213 // Copy members of master runner.
214 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
215 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
216 newRunner
.hw_opt
= hw_opt
;
217 newRunner
.filenames
= filenames
;
219 newRunner
.oenv
= oenv
;
220 newRunner
.mdrunOptions
= mdrunOptions
;
221 newRunner
.domdecOptions
= domdecOptions
;
222 newRunner
.nbpu_opt
= nbpu_opt
;
223 newRunner
.pme_opt
= pme_opt
;
224 newRunner
.pme_fft_opt
= pme_fft_opt
;
225 newRunner
.bonded_opt
= bonded_opt
;
226 newRunner
.nstlist_cmdline
= nstlist_cmdline
;
227 newRunner
.replExParams
= replExParams
;
228 newRunner
.pforce
= pforce
;
230 newRunner
.startingBehavior
= startingBehavior
;
231 newRunner
.stopHandlerBuilder_
= std::make_unique
<StopHandlerBuilder
>(*stopHandlerBuilder_
);
233 threadMpiMdrunnerAccessBarrier();
235 GMX_RELEASE_ASSERT(!MASTER(newRunner
.cr
), "cloneOnSpawnedThread should only be called on spawned threads");
240 /*! \brief The callback used for running on spawned threads.
242 * Obtains the pointer to the master mdrunner object from the one
243 * argument permitted to the thread-launch API call, copies it to make
244 * a new runner for this thread, reinitializes necessary data, and
245 * proceeds to the simulation. */
246 static void mdrunner_start_fn(const void *arg
)
250 auto masterMdrunner
= reinterpret_cast<const gmx::Mdrunner
*>(arg
);
251 /* copy the arg list to make sure that it's thread-local. This
252 doesn't copy pointed-to items, of course; fnm, cr and fplog
253 are reset in the call below, all others should be const. */
254 gmx::Mdrunner mdrunner
= masterMdrunner
->cloneOnSpawnedThread();
257 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
261 /*! \brief Start thread-MPI threads.
263 * Called by mdrunner() to start a specific number of threads
264 * (including the main thread) for thread-parallel runs. This in turn
265 * calls mdrunner() for each thread. All options are the same as for
267 t_commrec
*Mdrunner::spawnThreads(int numThreadsToLaunch
) const
270 /* first check whether we even need to start tMPI */
271 if (numThreadsToLaunch
< 2)
277 /* now spawn new threads that start mdrunner_start_fn(), while
278 the main thread returns. Thread affinity is handled later. */
279 if (tMPI_Init_fn(TRUE
, numThreadsToLaunch
, TMPI_AFFINITY_NONE
,
280 mdrunner_start_fn
, static_cast<const void*>(this)) != TMPI_SUCCESS
)
282 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
285 threadMpiMdrunnerAccessBarrier();
287 GMX_UNUSED_VALUE(mdrunner_start_fn
);
290 return reinitialize_commrec_for_this_thread(cr
);
295 /*! \brief Initialize variables for Verlet scheme simulation */
296 static void prepare_verlet_scheme(FILE *fplog
,
300 const gmx_mtop_t
*mtop
,
302 bool makeGpuPairList
,
303 const gmx::CpuInfo
&cpuinfo
)
305 /* For NVE simulations, we will retain the initial list buffer */
306 if (EI_DYNAMICS(ir
->eI
) &&
307 ir
->verletbuf_tol
> 0 &&
308 !(EI_MD(ir
->eI
) && ir
->etc
== etcNO
))
310 /* Update the Verlet buffer size for the current run setup */
312 /* Here we assume SIMD-enabled kernels are being used. But as currently
313 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
314 * and 4x2 gives a larger buffer than 4x4, this is ok.
316 ListSetupType listType
= (makeGpuPairList
? ListSetupType::Gpu
: ListSetupType::CpuSimdWhenSupported
);
317 VerletbufListSetup listSetup
= verletbufGetSafeListSetup(listType
);
319 const real rlist_new
=
320 calcVerletBufferSize(*mtop
, det(box
), *ir
, ir
->nstlist
, ir
->nstlist
- 1, -1, listSetup
);
322 if (rlist_new
!= ir
->rlist
)
324 if (fplog
!= nullptr)
326 fprintf(fplog
, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
327 ir
->rlist
, rlist_new
,
328 listSetup
.cluster_size_i
, listSetup
.cluster_size_j
);
330 ir
->rlist
= rlist_new
;
334 if (nstlist_cmdline
> 0 && (!EI_DYNAMICS(ir
->eI
) || ir
->verletbuf_tol
<= 0))
336 gmx_fatal(FARGS
, "Can not set nstlist without %s",
337 !EI_DYNAMICS(ir
->eI
) ? "dynamics" : "verlet-buffer-tolerance");
340 if (EI_DYNAMICS(ir
->eI
))
342 /* Set or try nstlist values */
343 increaseNstlist(fplog
, cr
, ir
, nstlist_cmdline
, mtop
, box
, makeGpuPairList
, cpuinfo
);
347 /*! \brief Override the nslist value in inputrec
349 * with value passed on the command line (if any)
351 static void override_nsteps_cmdline(const gmx::MDLogger
&mdlog
,
352 int64_t nsteps_cmdline
,
357 /* override with anything else than the default -2 */
358 if (nsteps_cmdline
> -2)
360 char sbuf_steps
[STEPSTRSIZE
];
361 char sbuf_msg
[STRLEN
];
363 ir
->nsteps
= nsteps_cmdline
;
364 if (EI_DYNAMICS(ir
->eI
) && nsteps_cmdline
!= -1)
366 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
367 gmx_step_str(nsteps_cmdline
, sbuf_steps
),
368 fabs(nsteps_cmdline
*ir
->delta_t
));
372 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps",
373 gmx_step_str(nsteps_cmdline
, sbuf_steps
));
376 GMX_LOG(mdlog
.warning
).asParagraph().appendText(sbuf_msg
);
378 else if (nsteps_cmdline
< -2)
380 gmx_fatal(FARGS
, "Invalid nsteps value passed on the command line: %" PRId64
,
383 /* Do nothing if nsteps_cmdline == -2 */
389 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
391 * If not, and if a warning may be issued, logs a warning about
392 * falling back to CPU code. With thread-MPI, only the first
393 * call to this function should have \c issueWarning true. */
394 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger
&mdlog
,
395 const t_inputrec
*ir
,
398 if (ir
->opts
.ngener
- ir
->nwall
> 1)
400 /* The GPU code does not support more than one energy group.
401 * If the user requested GPUs explicitly, a fatal error is given later.
405 GMX_LOG(mdlog
.warning
).asParagraph()
406 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
407 "For better performance, run on the GPU without energy groups and then do "
408 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
415 //! Initializes the logger for mdrun.
416 static gmx::LoggerOwner
buildLogger(FILE *fplog
, const t_commrec
*cr
)
418 gmx::LoggerBuilder builder
;
419 if (fplog
!= nullptr)
421 builder
.addTargetFile(gmx::MDLogger::LogLevel::Info
, fplog
);
423 if (cr
== nullptr || SIMMASTER(cr
))
425 builder
.addTargetStream(gmx::MDLogger::LogLevel::Warning
,
426 &gmx::TextOutputFile::standardError());
428 return builder
.build();
431 //! Make a TaskTarget from an mdrun argument string.
432 static TaskTarget
findTaskTarget(const char *optionString
)
434 TaskTarget returnValue
= TaskTarget::Auto
;
436 if (strncmp(optionString
, "auto", 3) == 0)
438 returnValue
= TaskTarget::Auto
;
440 else if (strncmp(optionString
, "cpu", 3) == 0)
442 returnValue
= TaskTarget::Cpu
;
444 else if (strncmp(optionString
, "gpu", 3) == 0)
446 returnValue
= TaskTarget::Gpu
;
450 GMX_ASSERT(false, "Option string should have been checked for sanity already");
456 //! Finish run, aggregate data to print performance info.
457 static void finish_run(FILE *fplog
,
458 const gmx::MDLogger
&mdlog
,
460 const t_inputrec
*inputrec
,
461 t_nrnb nrnb
[], gmx_wallcycle_t wcycle
,
462 gmx_walltime_accounting_t walltime_accounting
,
463 nonbonded_verlet_t
*nbv
,
464 const gmx_pme_t
*pme
,
468 double nbfs
= 0, mflop
= 0;
470 elapsed_time_over_all_ranks
,
471 elapsed_time_over_all_threads
,
472 elapsed_time_over_all_threads_over_all_ranks
;
473 /* Control whether it is valid to print a report. Only the
474 simulation master may print, but it should not do so if the run
475 terminated e.g. before a scheduled reset step. This is
476 complicated by the fact that PME ranks are unaware of the
477 reason why they were sent a pmerecvqxFINISH. To avoid
478 communication deadlocks, we always do the communication for the
479 report, even if we've decided not to write the report, because
480 how long it takes to finish the run is not important when we've
481 decided not to report on the simulation performance.
483 Further, we only report performance for dynamical integrators,
484 because those are the only ones for which we plan to
485 consider doing any optimizations. */
486 bool printReport
= EI_DYNAMICS(inputrec
->eI
) && SIMMASTER(cr
);
488 if (printReport
&& !walltime_accounting_get_valid_finish(walltime_accounting
))
490 GMX_LOG(mdlog
.warning
).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
495 std::unique_ptr
<t_nrnb
> nrnbTotalStorage
;
498 nrnbTotalStorage
= std::make_unique
<t_nrnb
>();
499 nrnb_tot
= nrnbTotalStorage
.get();
501 MPI_Allreduce(nrnb
->n
, nrnb_tot
->n
, eNRNB
, MPI_DOUBLE
, MPI_SUM
,
510 elapsed_time
= walltime_accounting_get_time_since_reset(walltime_accounting
);
511 elapsed_time_over_all_threads
= walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting
);
515 /* reduce elapsed_time over all MPI ranks in the current simulation */
516 MPI_Allreduce(&elapsed_time
,
517 &elapsed_time_over_all_ranks
,
518 1, MPI_DOUBLE
, MPI_SUM
,
520 elapsed_time_over_all_ranks
/= cr
->nnodes
;
521 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
522 * current simulation. */
523 MPI_Allreduce(&elapsed_time_over_all_threads
,
524 &elapsed_time_over_all_threads_over_all_ranks
,
525 1, MPI_DOUBLE
, MPI_SUM
,
531 elapsed_time_over_all_ranks
= elapsed_time
;
532 elapsed_time_over_all_threads_over_all_ranks
= elapsed_time_over_all_threads
;
537 print_flop(fplog
, nrnb_tot
, &nbfs
, &mflop
);
540 if (thisRankHasDuty(cr
, DUTY_PP
) && DOMAINDECOMP(cr
))
542 print_dd_statistics(cr
, inputrec
, fplog
);
545 /* TODO Move the responsibility for any scaling by thread counts
546 * to the code that handled the thread region, so that there's a
547 * mechanism to keep cycle counting working during the transition
548 * to task parallelism. */
549 int nthreads_pp
= gmx_omp_nthreads_get(emntNonbonded
);
550 int nthreads_pme
= gmx_omp_nthreads_get(emntPME
);
551 wallcycle_scale_by_num_threads(wcycle
, thisRankHasDuty(cr
, DUTY_PME
) && !thisRankHasDuty(cr
, DUTY_PP
), nthreads_pp
, nthreads_pme
);
552 auto cycle_sum(wallcycle_sum(cr
, wcycle
));
556 auto nbnxn_gpu_timings
= (nbv
!= nullptr && nbv
->useGpu()) ? Nbnxm::gpu_get_timings(nbv
->gpu_nbv
) : nullptr;
557 gmx_wallclock_gpu_pme_t pme_gpu_timings
= {};
559 if (pme_gpu_task_enabled(pme
))
561 pme_gpu_get_timings(pme
, &pme_gpu_timings
);
563 wallcycle_print(fplog
, mdlog
, cr
->nnodes
, cr
->npmenodes
, nthreads_pp
, nthreads_pme
,
564 elapsed_time_over_all_ranks
,
569 if (EI_DYNAMICS(inputrec
->eI
))
571 delta_t
= inputrec
->delta_t
;
576 print_perf(fplog
, elapsed_time_over_all_threads_over_all_ranks
,
577 elapsed_time_over_all_ranks
,
578 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting
),
579 delta_t
, nbfs
, mflop
);
583 print_perf(stderr
, elapsed_time_over_all_threads_over_all_ranks
,
584 elapsed_time_over_all_ranks
,
585 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting
),
586 delta_t
, nbfs
, mflop
);
591 int Mdrunner::mdrunner()
594 t_forcerec
*fr
= nullptr;
595 t_fcdata
*fcd
= nullptr;
596 real ewaldcoeff_q
= 0;
597 real ewaldcoeff_lj
= 0;
598 int nChargePerturbed
= -1, nTypePerturbed
= 0;
599 gmx_wallcycle_t wcycle
;
600 gmx_walltime_accounting_t walltime_accounting
= nullptr;
601 gmx_membed_t
* membed
= nullptr;
602 gmx_hw_info_t
*hwinfo
= nullptr;
604 /* CAUTION: threads may be started later on in this function, so
605 cr doesn't reflect the final parallel state right now */
606 t_inputrec inputrecInstance
;
607 t_inputrec
*inputrec
= &inputrecInstance
;
610 bool doMembed
= opt2bSet("-membed", filenames
.size(), filenames
.data());
611 bool doRerun
= mdrunOptions
.rerun
;
613 // Handle task-assignment related user options.
614 EmulateGpuNonbonded emulateGpuNonbonded
= (getenv("GMX_EMULATE_GPU") != nullptr ?
615 EmulateGpuNonbonded::Yes
: EmulateGpuNonbonded::No
);
617 std::vector
<int> userGpuTaskAssignment
;
620 userGpuTaskAssignment
= parseUserTaskAssignmentString(hw_opt
.userGpuTaskAssignment
);
622 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
623 auto nonbondedTarget
= findTaskTarget(nbpu_opt
);
624 auto pmeTarget
= findTaskTarget(pme_opt
);
625 auto pmeFftTarget
= findTaskTarget(pme_fft_opt
);
626 auto bondedTarget
= findTaskTarget(bonded_opt
);
627 PmeRunMode pmeRunMode
= PmeRunMode::None
;
629 // Here we assume that SIMMASTER(cr) does not change even after the
630 // threads are started.
632 FILE *fplog
= nullptr;
633 // If we are appending, we don't write log output because we need
634 // to check that the old log file matches what the checkpoint file
635 // expects. Otherwise, we should start to write log output now if
636 // there is a file ready for it.
637 if (logFileHandle
!= nullptr && startingBehavior
!= StartingBehavior::RestartWithAppending
)
639 fplog
= gmx_fio_getfp(logFileHandle
);
641 gmx::LoggerOwner
logOwner(buildLogger(fplog
, cr
));
642 gmx::MDLogger
mdlog(logOwner
.logger());
644 // report any development features that may be enabled by environment variables
645 reportDevelopmentFeatures(mdlog
);
647 // With thread-MPI, the communicator changes after threads are
648 // launched, so this is rebuilt for the master rank at that
649 // time. The non-master ranks are fine to keep the one made here.
650 PhysicalNodeCommunicator
physicalNodeComm(MPI_COMM_WORLD
, gmx_physicalnode_id_hash());
651 hwinfo
= gmx_detect_hardware(mdlog
, physicalNodeComm
);
653 gmx_print_detected_hardware(fplog
, isMasterSimMasterRank(ms
, MASTER(cr
)), mdlog
, hwinfo
);
655 std::vector
<int> gpuIdsToUse
= makeGpuIdsToUse(hwinfo
->gpu_info
, hw_opt
.gpuIdsAvailable
);
657 // Print citation requests after all software/hardware printing
658 pleaseCiteGromacs(fplog
);
660 std::unique_ptr
<t_state
> globalState
;
664 /* Only the master rank has the global state */
665 globalState
= std::make_unique
<t_state
>();
667 /* Read (nearly) all data required for the simulation */
668 read_tpx_state(ftp2fn(efTPR
, filenames
.size(), filenames
.data()), inputrec
, globalState
.get(), &mtop
);
671 /* Check and update the hardware options for internal consistency */
672 checkAndUpdateHardwareOptions(mdlog
, &hw_opt
, SIMMASTER(cr
), domdecOptions
.numPmeRanks
);
674 if (GMX_THREAD_MPI
&& SIMMASTER(cr
))
676 bool useGpuForNonbonded
= false;
677 bool useGpuForPme
= false;
680 // If the user specified the number of ranks, then we must
681 // respect that, but in default mode, we need to allow for
682 // the number of GPUs to choose the number of ranks.
683 auto canUseGpuForNonbonded
= buildSupportsNonbondedOnGpu(nullptr);
684 useGpuForNonbonded
= decideWhetherToUseGpusForNonbondedWithThreadMpi
685 (nonbondedTarget
, gpuIdsToUse
, userGpuTaskAssignment
, emulateGpuNonbonded
,
686 canUseGpuForNonbonded
,
687 gpuAccelerationOfNonbondedIsUseful(mdlog
, inputrec
, GMX_THREAD_MPI
),
688 hw_opt
.nthreads_tmpi
);
689 useGpuForPme
= decideWhetherToUseGpusForPmeWithThreadMpi
690 (useGpuForNonbonded
, pmeTarget
, gpuIdsToUse
, userGpuTaskAssignment
,
691 *hwinfo
, *inputrec
, mtop
, hw_opt
.nthreads_tmpi
, domdecOptions
.numPmeRanks
);
694 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
696 /* Determine how many thread-MPI ranks to start.
698 * TODO Over-writing the user-supplied value here does
699 * prevent any possible subsequent checks from working
701 hw_opt
.nthreads_tmpi
= get_nthreads_mpi(hwinfo
,
710 // Now start the threads for thread MPI.
711 cr
= spawnThreads(hw_opt
.nthreads_tmpi
);
712 /* The main thread continues here with a new cr. We don't deallocate
713 the old cr because other threads may still be reading it. */
714 // TODO Both master and spawned threads call dup_tfn and
715 // reinitialize_commrec_for_this_thread. Find a way to express
717 physicalNodeComm
= PhysicalNodeCommunicator(MPI_COMM_WORLD
, gmx_physicalnode_id_hash());
719 // END OF CAUTION: cr and physicalNodeComm are now reliable
723 /* now broadcast everything to the non-master nodes/threads: */
724 init_parallel(cr
, inputrec
, &mtop
);
727 // Now each rank knows the inputrec that SIMMASTER read and used,
728 // and (if applicable) cr->nnodes has been assigned the number of
729 // thread-MPI ranks that have been chosen. The ranks can now all
730 // run the task-deciding functions and will agree on the result
731 // without needing to communicate.
733 // TODO Should we do the communication in debug mode to support
734 // having an assertion?
736 // Note that these variables describe only their own node.
738 // Note that when bonded interactions run on a GPU they always run
739 // alongside a nonbonded task, so do not influence task assignment
740 // even though they affect the force calculation workload.
741 bool useGpuForNonbonded
= false;
742 bool useGpuForPme
= false;
743 bool useGpuForBonded
= false;
746 // It's possible that there are different numbers of GPUs on
747 // different nodes, which is the user's responsibilty to
748 // handle. If unsuitable, we will notice that during task
750 bool gpusWereDetected
= hwinfo
->ngpu_compatible_tot
> 0;
751 auto canUseGpuForNonbonded
= buildSupportsNonbondedOnGpu(nullptr);
752 useGpuForNonbonded
= decideWhetherToUseGpusForNonbonded(nonbondedTarget
, userGpuTaskAssignment
,
754 canUseGpuForNonbonded
,
755 gpuAccelerationOfNonbondedIsUseful(mdlog
, inputrec
, !GMX_THREAD_MPI
),
757 useGpuForPme
= decideWhetherToUseGpusForPme(useGpuForNonbonded
, pmeTarget
, userGpuTaskAssignment
,
758 *hwinfo
, *inputrec
, mtop
,
759 cr
->nnodes
, domdecOptions
.numPmeRanks
,
761 auto canUseGpuForBonded
= buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec
, mtop
, nullptr);
763 decideWhetherToUseGpusForBonded(useGpuForNonbonded
, useGpuForPme
,
764 bondedTarget
, canUseGpuForBonded
,
765 EVDW_PME(inputrec
->vdwtype
),
766 EEL_PME_EWALD(inputrec
->coulombtype
),
767 domdecOptions
.numPmeRanks
, gpusWereDetected
);
769 pmeRunMode
= (useGpuForPme
? PmeRunMode::GPU
: PmeRunMode::CPU
);
770 if (pmeRunMode
== PmeRunMode::GPU
)
772 if (pmeFftTarget
== TaskTarget::Cpu
)
774 pmeRunMode
= PmeRunMode::Mixed
;
777 else if (pmeFftTarget
== TaskTarget::Gpu
)
779 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.");
782 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
785 // TODO: hide restraint implementation details from Mdrunner.
786 // There is nothing unique about restraints at this point as far as the
787 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
788 // factory functions from the SimulationContext on which to call mdModules_->add().
789 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
790 for (auto && restraint
: restraintManager_
->getRestraints())
792 auto module
= RestraintMDModule::create(restraint
,
794 mdModules_
->add(std::move(module
));
797 // TODO: Error handling
798 mdModules_
->assignOptionsToModules(*inputrec
->params
, nullptr);
800 if (fplog
!= nullptr)
802 pr_inputrec(fplog
, 0, "Input Parameters", inputrec
, FALSE
);
803 fprintf(fplog
, "\n");
808 /* In rerun, set velocities to zero if present */
809 if (doRerun
&& ((globalState
->flags
& (1 << estV
)) != 0))
811 // rerun does not use velocities
812 GMX_LOG(mdlog
.info
).asParagraph().appendText(
813 "Rerun trajectory contains velocities. Rerun does only evaluate "
814 "potential energy and forces. The velocities will be ignored.");
815 for (int i
= 0; i
< globalState
->natoms
; i
++)
817 clear_rvec(globalState
->v
[i
]);
819 globalState
->flags
&= ~(1 << estV
);
822 /* now make sure the state is initialized and propagated */
823 set_state_entries(globalState
.get(), inputrec
);
826 /* NM and TPI parallelize over force/energy calculations, not atoms,
827 * so we need to initialize and broadcast the global state.
829 if (inputrec
->eI
== eiNM
|| inputrec
->eI
== eiTPI
)
833 globalState
= std::make_unique
<t_state
>();
835 broadcastStateWithoutDynamics(cr
, globalState
.get());
838 /* A parallel command line option consistency check that we can
839 only do after any threads have started. */
840 if (!PAR(cr
) && (domdecOptions
.numCells
[XX
] > 1 ||
841 domdecOptions
.numCells
[YY
] > 1 ||
842 domdecOptions
.numCells
[ZZ
] > 1 ||
843 domdecOptions
.numPmeRanks
> 0))
846 "The -dd or -npme option request a parallel simulation, "
848 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv
));
851 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
853 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv
));
859 (EI_ENERGY_MINIMIZATION(inputrec
->eI
) || eiNM
== inputrec
->eI
))
861 gmx_fatal(FARGS
, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
864 if (!(EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
)))
866 if (domdecOptions
.numPmeRanks
> 0)
868 gmx_fatal_collective(FARGS
, cr
->mpi_comm_mysim
, MASTER(cr
),
869 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
872 domdecOptions
.numPmeRanks
= 0;
875 if (useGpuForNonbonded
&& domdecOptions
.numPmeRanks
< 0)
877 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
878 * improve performance with many threads per GPU, since our OpenMP
879 * scaling is bad, but it's difficult to automate the setup.
881 domdecOptions
.numPmeRanks
= 0;
885 if (domdecOptions
.numPmeRanks
< 0)
887 domdecOptions
.numPmeRanks
= 0;
888 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
892 GMX_RELEASE_ASSERT(domdecOptions
.numPmeRanks
<= 1, "PME GPU decomposition is not supported");
899 fcRegisterSteps(inputrec
->nsteps
, inputrec
->init_step
);
903 /* NMR restraints must be initialized before load_checkpoint,
904 * since with time averaging the history is added to t_state.
905 * For proper consistency check we therefore need to extend
907 * So the PME-only nodes (if present) will also initialize
908 * the distance restraints.
912 /* This needs to be called before read_checkpoint to extend the state */
913 init_disres(fplog
, &mtop
, inputrec
, cr
, ms
, fcd
, globalState
.get(), replExParams
.exchangeInterval
> 0);
915 init_orires(fplog
, &mtop
, inputrec
, cr
, ms
, globalState
.get(), &(fcd
->orires
));
917 auto deform
= prepareBoxDeformation(globalState
->box
, cr
, *inputrec
);
919 ObservablesHistory observablesHistory
= {};
921 if (startingBehavior
!= StartingBehavior::NewSimulation
)
923 /* Check if checkpoint file exists before doing continuation.
924 * This way we can use identical input options for the first and subsequent runs...
926 if (mdrunOptions
.numStepsCommandline
> -2)
928 /* Temporarily set the number of steps to unmlimited to avoid
929 * triggering the nsteps check in load_checkpoint().
930 * This hack will go away soon when the -nsteps option is removed.
932 inputrec
->nsteps
= -1;
935 load_checkpoint(opt2fn_master("-cpi", filenames
.size(), filenames
.data(), cr
),
937 cr
, domdecOptions
.numCells
,
938 inputrec
, globalState
.get(),
940 mdrunOptions
.reproducible
);
942 if (startingBehavior
== StartingBehavior::RestartWithAppending
&& logFileHandle
)
944 // Now we can start normal logging to the truncated log file.
945 fplog
= gmx_fio_getfp(logFileHandle
);
946 prepareLogAppending(fplog
);
947 logOwner
= buildLogger(fplog
, cr
);
948 mdlog
= logOwner
.logger();
952 if (mdrunOptions
.numStepsCommandline
> -2)
954 GMX_LOG(mdlog
.info
).asParagraph().
955 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
956 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
958 /* override nsteps with value set on the commamdline */
959 override_nsteps_cmdline(mdlog
, mdrunOptions
.numStepsCommandline
, inputrec
);
963 copy_mat(globalState
->box
, box
);
968 gmx_bcast(sizeof(box
), box
, cr
);
971 if (inputrec
->cutoff_scheme
!= ecutsVERLET
)
973 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.");
975 /* Update rlist and nstlist. */
976 prepare_verlet_scheme(fplog
, cr
, inputrec
, nstlist_cmdline
, &mtop
, box
,
977 useGpuForNonbonded
|| (emulateGpuNonbonded
== EmulateGpuNonbonded::Yes
), *hwinfo
->cpuInfo
);
979 LocalAtomSetManager atomSets
;
981 if (PAR(cr
) && !(EI_TPI(inputrec
->eI
) ||
982 inputrec
->eI
== eiNM
))
984 cr
->dd
= init_domain_decomposition(mdlog
, cr
, domdecOptions
, mdrunOptions
,
986 box
, positionsFromStatePointer(globalState
.get()),
988 // Note that local state still does not exist yet.
992 /* PME, if used, is done on all nodes with 1D decomposition */
994 cr
->duty
= (DUTY_PP
| DUTY_PME
);
996 if (inputrec
->ePBC
== epbcSCREW
)
999 "pbc=screw is only implemented with domain decomposition");
1005 /* After possible communicator splitting in make_dd_communicators.
1006 * we can set up the intra/inter node communication.
1008 gmx_setup_nodecomm(fplog
, cr
);
1014 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1015 "This is simulation %d out of %d running as a composite GROMACS\n"
1016 "multi-simulation job. Setup for this simulation:\n",
1019 GMX_LOG(mdlog
.warning
).appendTextFormatted(
1020 "Using %d MPI %s\n",
1023 cr
->nnodes
== 1 ? "thread" : "threads"
1025 cr
->nnodes
== 1 ? "process" : "processes"
1031 // If mdrun -pin auto honors any affinity setting that already
1032 // exists. If so, it is nice to provide feedback about whether
1033 // that existing affinity setting was from OpenMP or something
1034 // else, so we run this code both before and after we initialize
1035 // the OpenMP support.
1036 gmx_check_thread_affinity_set(mdlog
,
1037 &hw_opt
, hwinfo
->nthreads_hw_avail
, FALSE
);
1038 /* Check and update the number of OpenMP threads requested */
1039 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt
, *hwinfo
, cr
, ms
, physicalNodeComm
.size_
,
1042 gmx_omp_nthreads_init(mdlog
, cr
,
1043 hwinfo
->nthreads_hw_avail
,
1044 physicalNodeComm
.size_
,
1045 hw_opt
.nthreads_omp
,
1046 hw_opt
.nthreads_omp_pme
,
1047 !thisRankHasDuty(cr
, DUTY_PP
));
1049 // Enable FP exception detection, but not in
1050 // Release mode and not for compilers with known buggy FP
1051 // exception support (clang with any optimization) or suspected
1052 // buggy FP exception support (gcc 7.* with optimization).
1053 #if !defined NDEBUG && \
1054 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1055 && defined __OPTIMIZE__)
1056 const bool bEnableFPE
= true;
1058 const bool bEnableFPE
= false;
1060 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1063 gmx_feenableexcept();
1066 // Build a data structure that expresses which kinds of non-bonded
1067 // task are handled by this rank.
1069 // TODO Later, this might become a loop over all registered modules
1070 // relevant to the mdp inputs, to find those that have such tasks.
1072 // TODO This could move before init_domain_decomposition() as part
1073 // of refactoring that separates the responsibility for duty
1074 // assignment from setup for communication between tasks, and
1075 // setup for tasks handled with a domain (ie including short-ranged
1076 // tasks, bonded tasks, etc.).
1078 // Note that in general useGpuForNonbonded, etc. can have a value
1079 // that is inconsistent with the presence of actual GPUs on any
1080 // rank, and that is not known to be a problem until the
1081 // duty of the ranks on a node become known.
1083 // TODO Later we might need the concept of computeTasksOnThisRank,
1084 // from which we construct gpuTasksOnThisRank.
1086 // Currently the DD code assigns duty to ranks that can
1087 // include PP work that currently can be executed on a single
1088 // GPU, if present and compatible. This has to be coordinated
1089 // across PP ranks on a node, with possible multiple devices
1090 // or sharing devices on a node, either from the user
1091 // selection, or automatically.
1092 auto haveGpus
= !gpuIdsToUse
.empty();
1093 std::vector
<GpuTask
> gpuTasksOnThisRank
;
1094 if (thisRankHasDuty(cr
, DUTY_PP
))
1096 if (useGpuForNonbonded
)
1098 // Note that any bonded tasks on a GPU always accompany a
1102 gpuTasksOnThisRank
.push_back(GpuTask::Nonbonded
);
1104 else if (nonbondedTarget
== TaskTarget::Gpu
)
1106 gmx_fatal(FARGS
, "Cannot run short-ranged nonbonded interactions on a GPU because no GPU is detected.");
1108 else if (bondedTarget
== TaskTarget::Gpu
)
1110 gmx_fatal(FARGS
, "Cannot run bonded interactions on a GPU because no GPU is detected.");
1114 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1115 if (EEL_PME(inputrec
->coulombtype
) && (thisRankHasDuty(cr
, DUTY_PME
)))
1121 gpuTasksOnThisRank
.push_back(GpuTask::Pme
);
1123 else if (pmeTarget
== TaskTarget::Gpu
)
1125 gmx_fatal(FARGS
, "Cannot run PME on a GPU because no GPU is detected.");
1130 GpuTaskAssignment gpuTaskAssignment
;
1133 // Produce the task assignment for this rank.
1134 gpuTaskAssignment
= runTaskAssignment(gpuIdsToUse
, userGpuTaskAssignment
, *hwinfo
,
1135 mdlog
, cr
, ms
, physicalNodeComm
, gpuTasksOnThisRank
,
1136 useGpuForBonded
, pmeRunMode
);
1138 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1140 /* Prevent other ranks from continuing after an issue was found
1141 * and reported as a fatal error.
1143 * TODO This function implements a barrier so that MPI runtimes
1144 * can organize an orderly shutdown if one of the ranks has had to
1145 * issue a fatal error in various code already run. When we have
1146 * MPI-aware error handling and reporting, this should be
1151 MPI_Barrier(cr
->mpi_comm_mysim
);
1157 MPI_Barrier(ms
->mpi_comm_masters
);
1159 /* We need another barrier to prevent non-master ranks from contiuing
1160 * when an error occured in a different simulation.
1162 MPI_Barrier(cr
->mpi_comm_mysim
);
1166 /* Now that we know the setup is consistent, check for efficiency */
1167 check_resource_division_efficiency(hwinfo
, !gpuTaskAssignment
.empty(), mdrunOptions
.ntompOptionIsSet
,
1170 gmx_device_info_t
*nonbondedDeviceInfo
= nullptr;
1172 if (thisRankHasDuty(cr
, DUTY_PP
))
1174 // This works because only one task of each type is currently permitted.
1175 auto nbGpuTaskMapping
= std::find_if(gpuTaskAssignment
.begin(), gpuTaskAssignment
.end(),
1176 hasTaskType
<GpuTask::Nonbonded
>);
1177 if (nbGpuTaskMapping
!= gpuTaskAssignment
.end())
1179 int nonbondedDeviceId
= nbGpuTaskMapping
->deviceId_
;
1180 nonbondedDeviceInfo
= getDeviceInfo(hwinfo
->gpu_info
, nonbondedDeviceId
);
1181 init_gpu(nonbondedDeviceInfo
);
1183 if (DOMAINDECOMP(cr
))
1185 /* When we share GPUs over ranks, we need to know this for the DLB */
1186 dd_setup_dlb_resource_sharing(cr
, nonbondedDeviceId
);
1192 std::unique_ptr
<ClfftInitializer
> initializedClfftLibrary
;
1194 gmx_device_info_t
*pmeDeviceInfo
= nullptr;
1195 // Later, this program could contain kernels that might be later
1196 // re-used as auto-tuning progresses, or subsequent simulations
1198 PmeGpuProgramStorage pmeGpuProgram
;
1199 // This works because only one task of each type is currently permitted.
1200 auto pmeGpuTaskMapping
= std::find_if(gpuTaskAssignment
.begin(), gpuTaskAssignment
.end(), hasTaskType
<GpuTask::Pme
>);
1201 const bool thisRankHasPmeGpuTask
= (pmeGpuTaskMapping
!= gpuTaskAssignment
.end());
1202 if (thisRankHasPmeGpuTask
)
1204 pmeDeviceInfo
= getDeviceInfo(hwinfo
->gpu_info
, pmeGpuTaskMapping
->deviceId_
);
1205 init_gpu(pmeDeviceInfo
);
1206 pmeGpuProgram
= buildPmeGpuProgram(pmeDeviceInfo
);
1207 // TODO It would be nice to move this logic into the factory
1208 // function. See Redmine #2535.
1209 bool isMasterThread
= !GMX_THREAD_MPI
|| MASTER(cr
);
1210 if (pmeRunMode
== PmeRunMode::GPU
&& !initializedClfftLibrary
&& isMasterThread
)
1212 initializedClfftLibrary
= initializeClfftLibrary();
1216 /* getting number of PP/PME threads on this MPI / tMPI rank.
1217 PME: env variable should be read only on one node to make sure it is
1218 identical everywhere;
1220 const int numThreadsOnThisRank
=
1221 thisRankHasDuty(cr
, DUTY_PP
) ? gmx_omp_nthreads_get(emntNonbonded
) : gmx_omp_nthreads_get(emntPME
);
1222 checkHardwareOversubscription(numThreadsOnThisRank
, cr
->nodeid
,
1223 *hwinfo
->hardwareTopology
,
1224 physicalNodeComm
, mdlog
);
1226 if (hw_opt
.threadAffinity
!= ThreadAffinity::Off
)
1228 /* Before setting affinity, check whether the affinity has changed
1229 * - which indicates that probably the OpenMP library has changed it
1230 * since we first checked).
1232 gmx_check_thread_affinity_set(mdlog
,
1233 &hw_opt
, hwinfo
->nthreads_hw_avail
, TRUE
);
1235 int numThreadsOnThisNode
, intraNodeThreadOffset
;
1236 analyzeThreadsOnThisNode(physicalNodeComm
, numThreadsOnThisRank
, &numThreadsOnThisNode
,
1237 &intraNodeThreadOffset
);
1239 /* Set the CPU affinity */
1240 gmx_set_thread_affinity(mdlog
, cr
, &hw_opt
, *hwinfo
->hardwareTopology
,
1241 numThreadsOnThisRank
, numThreadsOnThisNode
,
1242 intraNodeThreadOffset
, nullptr);
1245 if (mdrunOptions
.timingOptions
.resetStep
> -1)
1247 GMX_LOG(mdlog
.info
).asParagraph().
1248 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1250 wcycle
= wallcycle_init(fplog
, mdrunOptions
.timingOptions
.resetStep
, cr
);
1254 /* Master synchronizes its value of reset_counters with all nodes
1255 * including PME only nodes */
1256 int64_t reset_counters
= wcycle_get_reset_counters(wcycle
);
1257 gmx_bcast_sim(sizeof(reset_counters
), &reset_counters
, cr
);
1258 wcycle_set_reset_counters(wcycle
, reset_counters
);
1261 // Membrane embedding must be initialized before we call init_forcerec()
1266 fprintf(stderr
, "Initializing membed");
1268 /* Note that membed cannot work in parallel because mtop is
1269 * changed here. Fix this if we ever want to make it run with
1270 * multiple ranks. */
1271 membed
= init_membed(fplog
, filenames
.size(), filenames
.data(), &mtop
, inputrec
, globalState
.get(), cr
,
1273 .checkpointOptions
.period
);
1276 std::unique_ptr
<MDAtoms
> mdAtoms
;
1277 std::unique_ptr
<gmx_vsite_t
> vsite
;
1280 if (thisRankHasDuty(cr
, DUTY_PP
))
1282 /* Initiate forcerecord */
1283 fr
= new t_forcerec
;
1284 fr
->forceProviders
= mdModules_
->initForceProviders();
1285 init_forcerec(fplog
, mdlog
, fr
, fcd
,
1286 inputrec
, &mtop
, cr
, box
,
1287 opt2fn("-table", filenames
.size(), filenames
.data()),
1288 opt2fn("-tablep", filenames
.size(), filenames
.data()),
1289 opt2fns("-tableb", filenames
.size(), filenames
.data()),
1290 *hwinfo
, nonbondedDeviceInfo
,
1296 /* Initialize the mdAtoms structure.
1297 * mdAtoms is not filled with atom data,
1298 * as this can not be done now with domain decomposition.
1300 mdAtoms
= makeMDAtoms(fplog
, mtop
, *inputrec
, thisRankHasPmeGpuTask
);
1301 if (globalState
&& thisRankHasPmeGpuTask
)
1303 // The pinning of coordinates in the global state object works, because we only use
1304 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1305 // points to the global state object without DD.
1306 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1307 // which should also perform the pinning.
1308 changePinningPolicy(&globalState
->x
, pme_get_pinning_policy());
1311 /* Initialize the virtual site communication */
1312 vsite
= initVsite(mtop
, cr
);
1314 calc_shifts(box
, fr
->shift_vec
);
1316 /* With periodic molecules the charge groups should be whole at start up
1317 * and the virtual sites should not be far from their proper positions.
1319 if (!inputrec
->bContinuation
&& MASTER(cr
) &&
1320 !(inputrec
->ePBC
!= epbcNONE
&& inputrec
->bPeriodicMols
))
1322 /* Make molecules whole at start of run */
1323 if (fr
->ePBC
!= epbcNONE
)
1325 do_pbc_first_mtop(fplog
, inputrec
->ePBC
, box
, &mtop
, globalState
->x
.rvec_array());
1329 /* Correct initial vsite positions are required
1330 * for the initial distribution in the domain decomposition
1331 * and for the initial shell prediction.
1333 constructVsitesGlobal(mtop
, globalState
->x
);
1337 if (EEL_PME(fr
->ic
->eeltype
) || EVDW_PME(fr
->ic
->vdwtype
))
1339 ewaldcoeff_q
= fr
->ic
->ewaldcoeff_q
;
1340 ewaldcoeff_lj
= fr
->ic
->ewaldcoeff_lj
;
1345 /* This is a PME only node */
1347 GMX_ASSERT(globalState
== nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1349 ewaldcoeff_q
= calc_ewaldcoeff_q(inputrec
->rcoulomb
, inputrec
->ewald_rtol
);
1350 ewaldcoeff_lj
= calc_ewaldcoeff_lj(inputrec
->rvdw
, inputrec
->ewald_rtol_lj
);
1353 gmx_pme_t
*sepPmeData
= nullptr;
1354 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1355 GMX_ASSERT(thisRankHasDuty(cr
, DUTY_PP
) == (fr
!= nullptr), "Double-checking that only PME-only ranks have no forcerec");
1356 gmx_pme_t
* &pmedata
= fr
? fr
->pmedata
: sepPmeData
;
1358 /* Initiate PME if necessary,
1359 * either on all nodes or on dedicated PME nodes only. */
1360 if (EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
))
1362 if (mdAtoms
&& mdAtoms
->mdatoms())
1364 nChargePerturbed
= mdAtoms
->mdatoms()->nChargePerturbed
;
1365 if (EVDW_PME(inputrec
->vdwtype
))
1367 nTypePerturbed
= mdAtoms
->mdatoms()->nTypePerturbed
;
1370 if (cr
->npmenodes
> 0)
1372 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1373 gmx_bcast_sim(sizeof(nChargePerturbed
), &nChargePerturbed
, cr
);
1374 gmx_bcast_sim(sizeof(nTypePerturbed
), &nTypePerturbed
, cr
);
1377 if (thisRankHasDuty(cr
, DUTY_PME
))
1381 pmedata
= gmx_pme_init(cr
,
1382 getNumPmeDomains(cr
->dd
),
1384 nChargePerturbed
!= 0, nTypePerturbed
!= 0,
1385 mdrunOptions
.reproducible
,
1386 ewaldcoeff_q
, ewaldcoeff_lj
,
1387 gmx_omp_nthreads_get(emntPME
),
1388 pmeRunMode
, nullptr,
1389 pmeDeviceInfo
, pmeGpuProgram
.get(), mdlog
);
1391 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1396 if (EI_DYNAMICS(inputrec
->eI
))
1398 /* Turn on signal handling on all nodes */
1400 * (A user signal from the PME nodes (if any)
1401 * is communicated to the PP nodes.
1403 signal_handler_install();
1406 pull_t
*pull_work
= nullptr;
1407 if (thisRankHasDuty(cr
, DUTY_PP
))
1409 /* Assumes uniform use of the number of OpenMP threads */
1410 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntDefault
));
1412 if (inputrec
->bPull
)
1414 /* Initialize pull code */
1416 init_pull(fplog
, inputrec
->pull
, inputrec
,
1417 &mtop
, cr
, &atomSets
, inputrec
->fepvals
->init_lambda
);
1418 if (inputrec
->pull
->bXOutAverage
|| inputrec
->pull
->bFOutAverage
)
1420 initPullHistory(pull_work
, &observablesHistory
);
1422 if (EI_DYNAMICS(inputrec
->eI
) && MASTER(cr
))
1424 init_pull_output_files(pull_work
,
1425 filenames
.size(), filenames
.data(), oenv
,
1430 std::unique_ptr
<EnforcedRotation
> enforcedRotation
;
1433 /* Initialize enforced rotation code */
1434 enforcedRotation
= init_rot(fplog
,
1447 t_swap
*swap
= nullptr;
1448 if (inputrec
->eSwapCoords
!= eswapNO
)
1450 /* Initialize ion swapping code */
1451 swap
= init_swapcoords(fplog
, inputrec
,
1452 opt2fn_master("-swap", filenames
.size(), filenames
.data(), cr
),
1453 &mtop
, globalState
.get(), &observablesHistory
,
1454 cr
, &atomSets
, oenv
, mdrunOptions
,
1458 /* Let makeConstraints know whether we have essential dynamics constraints.
1459 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1461 bool doEssentialDynamics
= (opt2fn_null("-ei", filenames
.size(), filenames
.data()) != nullptr
1462 || observablesHistory
.edsamHistory
);
1463 auto constr
= makeConstraints(mtop
, *inputrec
, pull_work
, doEssentialDynamics
,
1464 fplog
, *mdAtoms
->mdatoms(),
1465 cr
, ms
, &nrnb
, wcycle
, fr
->bMolPBC
);
1467 /* Energy terms and groups */
1468 gmx_enerdata_t
enerd(mtop
.groups
.groups
[SimulationAtomGroupType::EnergyOutput
].size(), inputrec
->fepvals
->n_lambda
);
1470 /* Set up interactive MD (IMD) */
1471 auto imdSession
= makeImdSession(inputrec
, cr
, wcycle
, &enerd
, ms
, &mtop
, mdlog
,
1472 MASTER(cr
) ? globalState
->x
.rvec_array() : nullptr,
1473 filenames
.size(), filenames
.data(), oenv
, mdrunOptions
.imdOptions
,
1476 if (DOMAINDECOMP(cr
))
1478 GMX_RELEASE_ASSERT(fr
, "fr was NULL while cr->duty was DUTY_PP");
1479 /* This call is not included in init_domain_decomposition mainly
1480 * because fr->cginfo_mb is set later.
1482 dd_init_bondeds(fplog
, cr
->dd
, &mtop
, vsite
.get(), inputrec
,
1483 domdecOptions
.checkBondedInteractions
,
1487 // TODO This is not the right place to manage the lifetime of
1488 // this data structure, but currently it's the easiest way to
1489 // make it work. Later, it should probably be made/updated
1490 // after the workload for the lifetime of a PP domain is
1492 PpForceWorkload ppForceWorkload
;
1494 GMX_ASSERT(stopHandlerBuilder_
, "Runner must provide StopHandlerBuilder to simulator.");
1495 SimulatorBuilder simulatorBuilder
;
1497 // build and run simulator object based on user-input
1498 auto simulator
= simulatorBuilder
.build(
1499 fplog
, cr
, ms
, mdlog
, static_cast<int>(filenames
.size()), filenames
.data(),
1503 vsite
.get(), constr
.get(),
1504 enforcedRotation
? enforcedRotation
->getLegacyEnfrot() : nullptr,
1506 mdModules_
->outputProvider(),
1507 inputrec
, imdSession
.get(), pull_work
, swap
, &mtop
,
1510 &observablesHistory
,
1511 mdAtoms
.get(), &nrnb
, wcycle
, fr
,
1516 walltime_accounting
,
1517 std::move(stopHandlerBuilder_
),
1521 if (inputrec
->bPull
)
1523 finish_pull(pull_work
);
1525 finish_swapcoords(swap
);
1529 GMX_RELEASE_ASSERT(pmedata
, "pmedata was NULL while cr->duty was not DUTY_PP");
1531 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntPME
));
1532 gmx_pmeonly(pmedata
, cr
, &nrnb
, wcycle
, walltime_accounting
, inputrec
, pmeRunMode
);
1535 wallcycle_stop(wcycle
, ewcRUN
);
1537 /* Finish up, write some stuff
1538 * if rerunMD, don't write last frame again
1540 finish_run(fplog
, mdlog
, cr
,
1541 inputrec
, &nrnb
, wcycle
, walltime_accounting
,
1542 fr
? fr
->nbv
.get() : nullptr,
1544 EI_DYNAMICS(inputrec
->eI
) && !isMultiSim(ms
));
1546 // clean up cycle counter
1547 wallcycle_destroy(wcycle
);
1552 gmx_pme_destroy(pmedata
);
1556 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1557 // before we destroy the GPU context(s) in free_gpu_resources().
1558 // Pinned buffers are associated with contexts in CUDA.
1559 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1560 mdAtoms
.reset(nullptr);
1561 globalState
.reset(nullptr);
1562 mdModules_
.reset(nullptr); // destruct force providers here as they might also use the GPU
1564 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1565 free_gpu_resources(fr
, physicalNodeComm
, hwinfo
->gpu_info
);
1566 free_gpu(nonbondedDeviceInfo
);
1567 free_gpu(pmeDeviceInfo
);
1568 done_forcerec(fr
, mtop
.molblock
.size());
1573 free_membed(membed
);
1576 /* Does what it says */
1577 print_date_and_time(fplog
, cr
->nodeid
, "Finished mdrun", gmx_gettime());
1578 walltime_accounting_destroy(walltime_accounting
);
1580 // Ensure log file content is written
1583 gmx_fio_flush(logFileHandle
);
1586 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1587 * exceptions were enabled before function was called. */
1590 gmx_fedisableexcept();
1593 auto rc
= static_cast<int>(gmx_get_stop_condition());
1596 /* we need to join all threads. The sub-threads join when they
1597 exit this function, but the master thread needs to be told to
1599 if (PAR(cr
) && MASTER(cr
))
1603 //TODO free commrec in MPI simulations
1609 Mdrunner::~Mdrunner()
1611 // Clean up of the Manager.
1612 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1613 // but okay as long as threads synchronize some time before adding or accessing
1614 // a new set of restraints.
1615 if (restraintManager_
)
1617 restraintManager_
->clear();
1618 GMX_ASSERT(restraintManager_
->countRestraints() == 0,
1619 "restraints added during runner life time should be cleared at runner destruction.");
1623 void Mdrunner::addPotential(std::shared_ptr
<gmx::IRestraintPotential
> puller
,
1626 GMX_ASSERT(restraintManager_
, "Mdrunner must have a restraint manager.");
1627 // Not sure if this should be logged through the md logger or something else,
1628 // but it is helpful to have some sort of INFO level message sent somewhere.
1629 // std::cout << "Registering restraint named " << name << std::endl;
1631 // When multiple restraints are used, it may be wasteful to register them separately.
1632 // Maybe instead register an entire Restraint Manager as a force provider.
1633 restraintManager_
->addToSpec(std::move(puller
),
1637 Mdrunner::Mdrunner(std::unique_ptr
<MDModules
> mdModules
)
1638 : mdModules_(std::move(mdModules
))
1642 Mdrunner::Mdrunner(Mdrunner
&&) noexcept
= default;
1644 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1645 Mdrunner
&Mdrunner::operator=(Mdrunner
&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING
) = default;
1647 class Mdrunner::BuilderImplementation
1650 BuilderImplementation() = delete;
1651 BuilderImplementation(std::unique_ptr
<MDModules
> mdModules
,
1652 SimulationContext
* context
);
1653 ~BuilderImplementation();
1655 BuilderImplementation
&setExtraMdrunOptions(const MdrunOptions
&options
,
1656 real forceWarningThreshold
,
1657 StartingBehavior startingBehavior
);
1659 void addDomdec(const DomdecOptions
&options
);
1661 void addVerletList(int nstlist
);
1663 void addReplicaExchange(const ReplicaExchangeParameters
¶ms
);
1665 void addMultiSim(gmx_multisim_t
* multisim
);
1667 void addNonBonded(const char* nbpu_opt
);
1669 void addPME(const char* pme_opt_
, const char* pme_fft_opt_
);
1671 void addBondedTaskAssignment(const char* bonded_opt
);
1673 void addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
);
1675 void addFilenames(ArrayRef
<const t_filenm
> filenames
);
1677 void addOutputEnvironment(gmx_output_env_t
* outputEnvironment
);
1679 void addLogFile(t_fileio
*logFileHandle
);
1681 void addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
);
1687 // Default parameters copied from runner.h
1688 // \todo Clarify source(s) of default parameters.
1690 const char* nbpu_opt_
= nullptr;
1691 const char* pme_opt_
= nullptr;
1692 const char* pme_fft_opt_
= nullptr;
1693 const char *bonded_opt_
= nullptr;
1695 MdrunOptions mdrunOptions_
;
1697 DomdecOptions domdecOptions_
;
1699 ReplicaExchangeParameters replicaExchangeParameters_
;
1701 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1704 //! Non-owning multisim communicator handle.
1705 std::unique_ptr
<gmx_multisim_t
*> multisim_
= nullptr;
1707 //! Print a warning if any force is larger than this (in kJ/mol nm).
1708 real forceWarningThreshold_
= -1;
1710 //! Whether the simulation will start afresh, or restart with/without appending.
1711 StartingBehavior startingBehavior_
= StartingBehavior::NewSimulation
;
1713 //! The modules that comprise the functionality of mdrun.
1714 std::unique_ptr
<MDModules
> mdModules_
;
1716 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1719 * \todo Establish robust protocol to make sure resources remain valid.
1720 * SimulationContext will likely be separated into multiple layers for
1721 * different levels of access and different phases of execution. Ref
1722 * https://redmine.gromacs.org/issues/2375
1723 * https://redmine.gromacs.org/issues/2587
1725 SimulationContext
* context_
= nullptr;
1727 //! \brief Parallelism information.
1728 gmx_hw_opt_t hardwareOptions_
;
1730 //! filename options for simulation.
1731 ArrayRef
<const t_filenm
> filenames_
;
1733 /*! \brief Handle to output environment.
1735 * \todo gmx_output_env_t needs lifetime management.
1737 gmx_output_env_t
* outputEnvironment_
= nullptr;
1739 /*! \brief Non-owning handle to MD log file.
1741 * \todo Context should own output facilities for client.
1742 * \todo Improve log file handle management.
1744 * Code managing the FILE* relies on the ability to set it to
1745 * nullptr to check whether the filehandle is valid.
1747 t_fileio
* logFileHandle_
= nullptr;
1750 * \brief Builder for simulation stop signal handler.
1752 std::unique_ptr
<StopHandlerBuilder
> stopHandlerBuilder_
= nullptr;
1755 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr
<MDModules
> mdModules
,
1756 SimulationContext
* context
) :
1757 mdModules_(std::move(mdModules
)),
1760 GMX_ASSERT(context_
, "Bug found. It should not be possible to construct builder without a valid context.");
1763 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1765 Mdrunner::BuilderImplementation
&
1766 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions
&options
,
1767 const real forceWarningThreshold
,
1768 const StartingBehavior startingBehavior
)
1770 mdrunOptions_
= options
;
1771 forceWarningThreshold_
= forceWarningThreshold
;
1772 startingBehavior_
= startingBehavior
;
1776 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions
&options
)
1778 domdecOptions_
= options
;
1781 void Mdrunner::BuilderImplementation::addVerletList(int nstlist
)
1786 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1788 replicaExchangeParameters_
= params
;
1791 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t
* multisim
)
1793 multisim_
= std::make_unique
<gmx_multisim_t
*>(multisim
);
1796 Mdrunner
Mdrunner::BuilderImplementation::build()
1798 auto newRunner
= Mdrunner(std::move(mdModules_
));
1800 GMX_ASSERT(context_
, "Bug found. It should not be possible to call build() without a valid context.");
1802 newRunner
.mdrunOptions
= mdrunOptions_
;
1803 newRunner
.pforce
= forceWarningThreshold_
;
1804 newRunner
.startingBehavior
= startingBehavior_
;
1805 newRunner
.domdecOptions
= domdecOptions_
;
1807 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1808 newRunner
.hw_opt
= hardwareOptions_
;
1810 // No invariant to check. This parameter exists to optionally override other behavior.
1811 newRunner
.nstlist_cmdline
= nstlist_
;
1813 newRunner
.replExParams
= replicaExchangeParameters_
;
1815 newRunner
.filenames
= filenames_
;
1817 GMX_ASSERT(context_
->communicationRecord_
, "SimulationContext communications not initialized.");
1818 newRunner
.cr
= context_
->communicationRecord_
;
1822 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1823 newRunner
.ms
= *multisim_
;
1827 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1830 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1831 // \todo Update sanity checking when output environment has clearly specified invariants.
1832 // Initialization and default values for oenv are not well specified in the current version.
1833 if (outputEnvironment_
)
1835 newRunner
.oenv
= outputEnvironment_
;
1839 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1842 newRunner
.logFileHandle
= logFileHandle_
;
1846 newRunner
.nbpu_opt
= nbpu_opt_
;
1850 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1853 if (pme_opt_
&& pme_fft_opt_
)
1855 newRunner
.pme_opt
= pme_opt_
;
1856 newRunner
.pme_fft_opt
= pme_fft_opt_
;
1860 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1865 newRunner
.bonded_opt
= bonded_opt_
;
1869 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1872 newRunner
.restraintManager_
= std::make_unique
<gmx::RestraintManager
>();
1874 if (stopHandlerBuilder_
)
1876 newRunner
.stopHandlerBuilder_
= std::move(stopHandlerBuilder_
);
1880 newRunner
.stopHandlerBuilder_
= std::make_unique
<StopHandlerBuilder
>();
1886 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt
)
1888 nbpu_opt_
= nbpu_opt
;
1891 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt
,
1892 const char* pme_fft_opt
)
1895 pme_fft_opt_
= pme_fft_opt
;
1898 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt
)
1900 bonded_opt_
= bonded_opt
;
1903 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
1905 hardwareOptions_
= hardwareOptions
;
1908 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef
<const t_filenm
> filenames
)
1910 filenames_
= filenames
;
1913 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
1915 outputEnvironment_
= outputEnvironment
;
1918 void Mdrunner::BuilderImplementation::addLogFile(t_fileio
*logFileHandle
)
1920 logFileHandle_
= logFileHandle
;
1923 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
1925 stopHandlerBuilder_
= std::move(builder
);
1928 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr
<MDModules
> mdModules
,
1929 compat::not_null
<SimulationContext
*> context
) :
1930 impl_
{std::make_unique
<Mdrunner::BuilderImplementation
>(std::move(mdModules
), context
)}
1934 MdrunnerBuilder::~MdrunnerBuilder() = default;
1936 MdrunnerBuilder
&MdrunnerBuilder::addSimulationMethod(const MdrunOptions
&options
,
1937 real forceWarningThreshold
,
1938 const StartingBehavior startingBehavior
)
1940 impl_
->setExtraMdrunOptions(options
, forceWarningThreshold
, startingBehavior
);
1944 MdrunnerBuilder
&MdrunnerBuilder::addDomainDecomposition(const DomdecOptions
&options
)
1946 impl_
->addDomdec(options
);
1950 MdrunnerBuilder
&MdrunnerBuilder::addNeighborList(int nstlist
)
1952 impl_
->addVerletList(nstlist
);
1956 MdrunnerBuilder
&MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1958 impl_
->addReplicaExchange(params
);
1962 MdrunnerBuilder
&MdrunnerBuilder::addMultiSim(gmx_multisim_t
* multisim
)
1964 impl_
->addMultiSim(multisim
);
1968 MdrunnerBuilder
&MdrunnerBuilder::addNonBonded(const char* nbpu_opt
)
1970 impl_
->addNonBonded(nbpu_opt
);
1974 MdrunnerBuilder
&MdrunnerBuilder::addElectrostatics(const char* pme_opt
,
1975 const char* pme_fft_opt
)
1977 // The builder method may become more general in the future, but in this version,
1978 // parameters for PME electrostatics are both required and the only parameters
1980 if (pme_opt
&& pme_fft_opt
)
1982 impl_
->addPME(pme_opt
, pme_fft_opt
);
1986 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
1991 MdrunnerBuilder
&MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt
)
1993 impl_
->addBondedTaskAssignment(bonded_opt
);
1997 Mdrunner
MdrunnerBuilder::build()
1999 return impl_
->build();
2002 MdrunnerBuilder
&MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
2004 impl_
->addHardwareOptions(hardwareOptions
);
2008 MdrunnerBuilder
&MdrunnerBuilder::addFilenames(ArrayRef
<const t_filenm
> filenames
)
2010 impl_
->addFilenames(filenames
);
2014 MdrunnerBuilder
&MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
2016 impl_
->addOutputEnvironment(outputEnvironment
);
2020 MdrunnerBuilder
&MdrunnerBuilder::addLogFile(t_fileio
*logFileHandle
)
2022 impl_
->addLogFile(logFileHandle
);
2026 MdrunnerBuilder
&MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
2028 impl_
->addStopHandlerBuilder(std::move(builder
));
2032 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder
&&) noexcept
= default;
2034 MdrunnerBuilder
&MdrunnerBuilder::operator=(MdrunnerBuilder
&&) noexcept
= default;