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/builder.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
64 #include "gromacs/domdec/localatomsetmanager.h"
65 #include "gromacs/domdec/partition.h"
66 #include "gromacs/ewald/ewald_utils.h"
67 #include "gromacs/ewald/pme.h"
68 #include "gromacs/ewald/pme_gpu_program.h"
69 #include "gromacs/fileio/checkpoint.h"
70 #include "gromacs/fileio/gmxfio.h"
71 #include "gromacs/fileio/oenv.h"
72 #include "gromacs/fileio/tpxio.h"
73 #include "gromacs/gmxlib/network.h"
74 #include "gromacs/gmxlib/nrnb.h"
75 #include "gromacs/gpu_utils/gpu_utils.h"
76 #include "gromacs/hardware/cpuinfo.h"
77 #include "gromacs/hardware/detecthardware.h"
78 #include "gromacs/hardware/printhardware.h"
79 #include "gromacs/imd/imd.h"
80 #include "gromacs/listed_forces/disre.h"
81 #include "gromacs/listed_forces/gpubonded.h"
82 #include "gromacs/listed_forces/orires.h"
83 #include "gromacs/math/functions.h"
84 #include "gromacs/math/utilities.h"
85 #include "gromacs/math/vec.h"
86 #include "gromacs/mdlib/boxdeformation.h"
87 #include "gromacs/mdlib/broadcaststructs.h"
88 #include "gromacs/mdlib/calc_verletbuf.h"
89 #include "gromacs/mdlib/dispersioncorrection.h"
90 #include "gromacs/mdlib/enerdata_utils.h"
91 #include "gromacs/mdlib/force.h"
92 #include "gromacs/mdlib/forcerec.h"
93 #include "gromacs/mdlib/gmx_omp_nthreads.h"
94 #include "gromacs/mdlib/makeconstraints.h"
95 #include "gromacs/mdlib/md_support.h"
96 #include "gromacs/mdlib/mdatoms.h"
97 #include "gromacs/mdlib/membed.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/group.h"
112 #include "gromacs/mdtypes/inputrec.h"
113 #include "gromacs/mdtypes/md_enums.h"
114 #include "gromacs/mdtypes/mdrunoptions.h"
115 #include "gromacs/mdtypes/observableshistory.h"
116 #include "gromacs/mdtypes/simulation_workload.h"
117 #include "gromacs/mdtypes/state.h"
118 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
119 #include "gromacs/nbnxm/gpu_data_mgmt.h"
120 #include "gromacs/nbnxm/nbnxm.h"
121 #include "gromacs/nbnxm/pairlist_tuning.h"
122 #include "gromacs/pbcutil/pbc.h"
123 #include "gromacs/pulling/output.h"
124 #include "gromacs/pulling/pull.h"
125 #include "gromacs/pulling/pull_rotation.h"
126 #include "gromacs/restraint/manager.h"
127 #include "gromacs/restraint/restraintmdmodule.h"
128 #include "gromacs/restraint/restraintpotential.h"
129 #include "gromacs/swap/swapcoords.h"
130 #include "gromacs/taskassignment/decidegpuusage.h"
131 #include "gromacs/taskassignment/decidesimulationworkload.h"
132 #include "gromacs/taskassignment/resourcedivision.h"
133 #include "gromacs/taskassignment/taskassignment.h"
134 #include "gromacs/taskassignment/usergpuids.h"
135 #include "gromacs/timing/gpu_timing.h"
136 #include "gromacs/timing/wallcycle.h"
137 #include "gromacs/timing/wallcyclereporting.h"
138 #include "gromacs/topology/mtop_util.h"
139 #include "gromacs/trajectory/trajectoryframe.h"
140 #include "gromacs/utility/basenetwork.h"
141 #include "gromacs/utility/cstringutil.h"
142 #include "gromacs/utility/exceptions.h"
143 #include "gromacs/utility/fatalerror.h"
144 #include "gromacs/utility/filestream.h"
145 #include "gromacs/utility/gmxassert.h"
146 #include "gromacs/utility/gmxmpi.h"
147 #include "gromacs/utility/keyvaluetree.h"
148 #include "gromacs/utility/logger.h"
149 #include "gromacs/utility/loggerbuilder.h"
150 #include "gromacs/utility/mdmodulenotification.h"
151 #include "gromacs/utility/physicalnodecommunicator.h"
152 #include "gromacs/utility/pleasecite.h"
153 #include "gromacs/utility/programcontext.h"
154 #include "gromacs/utility/smalloc.h"
155 #include "gromacs/utility/stringutil.h"
157 #include "isimulator.h"
158 #include "replicaexchange.h"
159 #include "simulatorbuilder.h"
162 #include "corewrap.h"
168 /*! \brief Structure that holds boolean flags corresponding to the development
169 * features present enabled through environemnt variables.
172 struct DevelopmentFeatureFlags
174 ///! True if the Buffer ops development feature is enabled
175 // TODO: when the trigger of the buffer ops offload is fully automated this should go away
176 bool enableGpuBufferOps
= false;
177 ///! True if the update-constraints development feature is enabled
178 // TODO This needs to be reomved when the code gets cleaned up of GMX_UPDATE_CONSTRAIN_GPU
179 bool useGpuUpdateConstrain
= false;
180 ///! True if the GPU halo exchange development feature is enabled
181 bool enableGpuHaloExchange
= false;
182 ///! True if the PME PP direct commuinication GPU development feature is enabled
183 bool enableGpuPmePPComm
= false;
186 /*! \brief Manage any development feature flag variables encountered
188 * The use of dev features indicated by environment variables is
189 * logged in order to ensure that runs with such features enabled can
190 * be identified from their log and standard output. Any cross
191 * dependencies are also checked, and if unsatisfied, a fatal error
194 * Note that some development features overrides are applied already here:
195 * the GPU communication flags are set to false in non-tMPI and non-CUDA builds.
197 * \param[in] mdlog Logger object.
198 * \returns The object populated with development feature flags.
200 static DevelopmentFeatureFlags
manageDevelopmentFeatures(const gmx::MDLogger
&mdlog
)
202 DevelopmentFeatureFlags devFlags
;
204 devFlags
.enableGpuBufferOps
= (getenv("GMX_USE_GPU_BUFFER_OPS") != nullptr);
205 devFlags
.useGpuUpdateConstrain
= (getenv("GMX_UPDATE_CONSTRAIN_GPU") != nullptr);
206 devFlags
.enableGpuHaloExchange
= (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI
&& (GMX_GPU
== GMX_GPU_CUDA
));
207 devFlags
.enableGpuPmePPComm
= (getenv("GMX_GPU_DD_COMMS") != nullptr && GMX_THREAD_MPI
&& (GMX_GPU
== GMX_GPU_CUDA
));
209 if (devFlags
.enableGpuBufferOps
)
211 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
212 "NOTE: This run uses the 'GPU buffer ops' feature, enabled by the GMX_USE_GPU_BUFFER_OPS environment variable.");
215 if (devFlags
.enableGpuHaloExchange
)
217 if (!devFlags
.enableGpuBufferOps
)
219 gmx_fatal(FARGS
, "Cannot enable GPU halo exchange without GPU buffer operations, set GMX_USE_GPU_BUFFER_OPS=1\n");
223 if (devFlags
.useGpuUpdateConstrain
)
225 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
226 "NOTE: This run uses the 'GPU update/constraints' feature, enabled by the GMX_UPDATE_CONSTRAIN_GPU environment variable.");
232 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
234 * Used to ensure that the master thread does not modify mdrunner during copy
235 * on the spawned threads. */
236 static void threadMpiMdrunnerAccessBarrier()
239 MPI_Barrier(MPI_COMM_WORLD
);
243 Mdrunner
Mdrunner::cloneOnSpawnedThread() const
245 auto newRunner
= Mdrunner(std::make_unique
<MDModules
>());
247 // All runners in the same process share a restraint manager resource because it is
248 // part of the interface to the client code, which is associated only with the
249 // original thread. Handles to the same resources can be obtained by copy.
251 newRunner
.restraintManager_
= std::make_unique
<RestraintManager
>(*restraintManager_
);
254 // Copy members of master runner.
255 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
256 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
257 newRunner
.hw_opt
= hw_opt
;
258 newRunner
.filenames
= filenames
;
260 newRunner
.oenv
= oenv
;
261 newRunner
.mdrunOptions
= mdrunOptions
;
262 newRunner
.domdecOptions
= domdecOptions
;
263 newRunner
.nbpu_opt
= nbpu_opt
;
264 newRunner
.pme_opt
= pme_opt
;
265 newRunner
.pme_fft_opt
= pme_fft_opt
;
266 newRunner
.bonded_opt
= bonded_opt
;
267 newRunner
.update_opt
= update_opt
;
268 newRunner
.nstlist_cmdline
= nstlist_cmdline
;
269 newRunner
.replExParams
= replExParams
;
270 newRunner
.pforce
= pforce
;
271 // Give the spawned thread the newly created valid communicator
272 // for the simulation.
273 newRunner
.communicator
= MPI_COMM_WORLD
;
275 newRunner
.startingBehavior
= startingBehavior
;
276 newRunner
.stopHandlerBuilder_
= std::make_unique
<StopHandlerBuilder
>(*stopHandlerBuilder_
);
278 threadMpiMdrunnerAccessBarrier();
283 /*! \brief The callback used for running on spawned threads.
285 * Obtains the pointer to the master mdrunner object from the one
286 * argument permitted to the thread-launch API call, copies it to make
287 * a new runner for this thread, reinitializes necessary data, and
288 * proceeds to the simulation. */
289 static void mdrunner_start_fn(const void *arg
)
293 auto masterMdrunner
= reinterpret_cast<const gmx::Mdrunner
*>(arg
);
294 /* copy the arg list to make sure that it's thread-local. This
295 doesn't copy pointed-to items, of course; fnm, cr and fplog
296 are reset in the call below, all others should be const. */
297 gmx::Mdrunner mdrunner
= masterMdrunner
->cloneOnSpawnedThread();
300 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
304 void Mdrunner::spawnThreads(int numThreadsToLaunch
)
307 /* now spawn new threads that start mdrunner_start_fn(), while
308 the main thread returns. Thread affinity is handled later. */
309 if (tMPI_Init_fn(TRUE
, numThreadsToLaunch
, TMPI_AFFINITY_NONE
,
310 mdrunner_start_fn
, static_cast<const void*>(this)) != TMPI_SUCCESS
)
312 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
315 // Give the master thread the newly created valid communicator for
317 communicator
= MPI_COMM_WORLD
;
318 threadMpiMdrunnerAccessBarrier();
320 GMX_UNUSED_VALUE(numThreadsToLaunch
);
321 GMX_UNUSED_VALUE(mdrunner_start_fn
);
327 /*! \brief Initialize variables for Verlet scheme simulation */
328 static void prepare_verlet_scheme(FILE *fplog
,
332 const gmx_mtop_t
*mtop
,
334 bool makeGpuPairList
,
335 const gmx::CpuInfo
&cpuinfo
)
337 /* For NVE simulations, we will retain the initial list buffer */
338 if (EI_DYNAMICS(ir
->eI
) &&
339 ir
->verletbuf_tol
> 0 &&
340 !(EI_MD(ir
->eI
) && ir
->etc
== etcNO
))
342 /* Update the Verlet buffer size for the current run setup */
344 /* Here we assume SIMD-enabled kernels are being used. But as currently
345 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
346 * and 4x2 gives a larger buffer than 4x4, this is ok.
348 ListSetupType listType
= (makeGpuPairList
? ListSetupType::Gpu
: ListSetupType::CpuSimdWhenSupported
);
349 VerletbufListSetup listSetup
= verletbufGetSafeListSetup(listType
);
351 const real rlist_new
=
352 calcVerletBufferSize(*mtop
, det(box
), *ir
, ir
->nstlist
, ir
->nstlist
- 1, -1, listSetup
);
354 if (rlist_new
!= ir
->rlist
)
356 if (fplog
!= nullptr)
358 fprintf(fplog
, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
359 ir
->rlist
, rlist_new
,
360 listSetup
.cluster_size_i
, listSetup
.cluster_size_j
);
362 ir
->rlist
= rlist_new
;
366 if (nstlist_cmdline
> 0 && (!EI_DYNAMICS(ir
->eI
) || ir
->verletbuf_tol
<= 0))
368 gmx_fatal(FARGS
, "Can not set nstlist without %s",
369 !EI_DYNAMICS(ir
->eI
) ? "dynamics" : "verlet-buffer-tolerance");
372 if (EI_DYNAMICS(ir
->eI
))
374 /* Set or try nstlist values */
375 increaseNstlist(fplog
, cr
, ir
, nstlist_cmdline
, mtop
, box
, makeGpuPairList
, cpuinfo
);
379 /*! \brief Override the nslist value in inputrec
381 * with value passed on the command line (if any)
383 static void override_nsteps_cmdline(const gmx::MDLogger
&mdlog
,
384 int64_t nsteps_cmdline
,
389 /* override with anything else than the default -2 */
390 if (nsteps_cmdline
> -2)
392 char sbuf_steps
[STEPSTRSIZE
];
393 char sbuf_msg
[STRLEN
];
395 ir
->nsteps
= nsteps_cmdline
;
396 if (EI_DYNAMICS(ir
->eI
) && nsteps_cmdline
!= -1)
398 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
399 gmx_step_str(nsteps_cmdline
, sbuf_steps
),
400 fabs(nsteps_cmdline
*ir
->delta_t
));
404 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps",
405 gmx_step_str(nsteps_cmdline
, sbuf_steps
));
408 GMX_LOG(mdlog
.warning
).asParagraph().appendText(sbuf_msg
);
410 else if (nsteps_cmdline
< -2)
412 gmx_fatal(FARGS
, "Invalid nsteps value passed on the command line: %" PRId64
,
415 /* Do nothing if nsteps_cmdline == -2 */
421 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
423 * If not, and if a warning may be issued, logs a warning about
424 * falling back to CPU code. With thread-MPI, only the first
425 * call to this function should have \c issueWarning true. */
426 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger
&mdlog
,
427 const t_inputrec
&ir
,
430 bool gpuIsUseful
= true;
433 if (ir
.opts
.ngener
- ir
.nwall
> 1)
435 /* The GPU code does not support more than one energy group.
436 * If the user requested GPUs explicitly, a fatal error is given later.
440 "Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
441 "For better performance, run on the GPU without energy groups and then do "
442 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.";
448 warning
= "TPI is not implemented for GPUs.";
451 if (!gpuIsUseful
&& issueWarning
)
453 GMX_LOG(mdlog
.warning
).asParagraph().appendText(warning
);
459 //! Initializes the logger for mdrun.
460 static gmx::LoggerOwner
buildLogger(FILE *fplog
,
461 const bool isSimulationMasterRank
)
463 gmx::LoggerBuilder builder
;
464 if (fplog
!= nullptr)
466 builder
.addTargetFile(gmx::MDLogger::LogLevel::Info
, fplog
);
468 if (isSimulationMasterRank
)
470 builder
.addTargetStream(gmx::MDLogger::LogLevel::Warning
,
471 &gmx::TextOutputFile::standardError());
473 return builder
.build();
476 //! Make a TaskTarget from an mdrun argument string.
477 static TaskTarget
findTaskTarget(const char *optionString
)
479 TaskTarget returnValue
= TaskTarget::Auto
;
481 if (strncmp(optionString
, "auto", 3) == 0)
483 returnValue
= TaskTarget::Auto
;
485 else if (strncmp(optionString
, "cpu", 3) == 0)
487 returnValue
= TaskTarget::Cpu
;
489 else if (strncmp(optionString
, "gpu", 3) == 0)
491 returnValue
= TaskTarget::Gpu
;
495 GMX_ASSERT(false, "Option string should have been checked for sanity already");
501 //! Finish run, aggregate data to print performance info.
502 static void finish_run(FILE *fplog
,
503 const gmx::MDLogger
&mdlog
,
505 const t_inputrec
*inputrec
,
506 t_nrnb nrnb
[], gmx_wallcycle_t wcycle
,
507 gmx_walltime_accounting_t walltime_accounting
,
508 nonbonded_verlet_t
*nbv
,
509 const gmx_pme_t
*pme
,
513 double nbfs
= 0, mflop
= 0;
515 elapsed_time_over_all_ranks
,
516 elapsed_time_over_all_threads
,
517 elapsed_time_over_all_threads_over_all_ranks
;
518 /* Control whether it is valid to print a report. Only the
519 simulation master may print, but it should not do so if the run
520 terminated e.g. before a scheduled reset step. This is
521 complicated by the fact that PME ranks are unaware of the
522 reason why they were sent a pmerecvqxFINISH. To avoid
523 communication deadlocks, we always do the communication for the
524 report, even if we've decided not to write the report, because
525 how long it takes to finish the run is not important when we've
526 decided not to report on the simulation performance.
528 Further, we only report performance for dynamical integrators,
529 because those are the only ones for which we plan to
530 consider doing any optimizations. */
531 bool printReport
= EI_DYNAMICS(inputrec
->eI
) && SIMMASTER(cr
);
533 if (printReport
&& !walltime_accounting_get_valid_finish(walltime_accounting
))
535 GMX_LOG(mdlog
.warning
).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
540 std::unique_ptr
<t_nrnb
> nrnbTotalStorage
;
543 nrnbTotalStorage
= std::make_unique
<t_nrnb
>();
544 nrnb_tot
= nrnbTotalStorage
.get();
546 MPI_Allreduce(nrnb
->n
, nrnb_tot
->n
, eNRNB
, MPI_DOUBLE
, MPI_SUM
,
555 elapsed_time
= walltime_accounting_get_time_since_reset(walltime_accounting
);
556 elapsed_time_over_all_threads
= walltime_accounting_get_time_since_reset_over_all_threads(walltime_accounting
);
560 /* reduce elapsed_time over all MPI ranks in the current simulation */
561 MPI_Allreduce(&elapsed_time
,
562 &elapsed_time_over_all_ranks
,
563 1, MPI_DOUBLE
, MPI_SUM
,
565 elapsed_time_over_all_ranks
/= cr
->nnodes
;
566 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
567 * current simulation. */
568 MPI_Allreduce(&elapsed_time_over_all_threads
,
569 &elapsed_time_over_all_threads_over_all_ranks
,
570 1, MPI_DOUBLE
, MPI_SUM
,
576 elapsed_time_over_all_ranks
= elapsed_time
;
577 elapsed_time_over_all_threads_over_all_ranks
= elapsed_time_over_all_threads
;
582 print_flop(fplog
, nrnb_tot
, &nbfs
, &mflop
);
585 if (thisRankHasDuty(cr
, DUTY_PP
) && DOMAINDECOMP(cr
))
587 print_dd_statistics(cr
, inputrec
, fplog
);
590 /* TODO Move the responsibility for any scaling by thread counts
591 * to the code that handled the thread region, so that there's a
592 * mechanism to keep cycle counting working during the transition
593 * to task parallelism. */
594 int nthreads_pp
= gmx_omp_nthreads_get(emntNonbonded
);
595 int nthreads_pme
= gmx_omp_nthreads_get(emntPME
);
596 wallcycle_scale_by_num_threads(wcycle
, thisRankHasDuty(cr
, DUTY_PME
) && !thisRankHasDuty(cr
, DUTY_PP
), nthreads_pp
, nthreads_pme
);
597 auto cycle_sum(wallcycle_sum(cr
, wcycle
));
601 auto nbnxn_gpu_timings
= (nbv
!= nullptr && nbv
->useGpu()) ? Nbnxm::gpu_get_timings(nbv
->gpu_nbv
) : nullptr;
602 gmx_wallclock_gpu_pme_t pme_gpu_timings
= {};
604 if (pme_gpu_task_enabled(pme
))
606 pme_gpu_get_timings(pme
, &pme_gpu_timings
);
608 wallcycle_print(fplog
, mdlog
, cr
->nnodes
, cr
->npmenodes
, nthreads_pp
, nthreads_pme
,
609 elapsed_time_over_all_ranks
,
614 if (EI_DYNAMICS(inputrec
->eI
))
616 delta_t
= inputrec
->delta_t
;
621 print_perf(fplog
, elapsed_time_over_all_threads_over_all_ranks
,
622 elapsed_time_over_all_ranks
,
623 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting
),
624 delta_t
, nbfs
, mflop
);
628 print_perf(stderr
, elapsed_time_over_all_threads_over_all_ranks
,
629 elapsed_time_over_all_ranks
,
630 walltime_accounting_get_nsteps_done_since_reset(walltime_accounting
),
631 delta_t
, nbfs
, mflop
);
636 int Mdrunner::mdrunner()
639 t_forcerec
*fr
= nullptr;
640 t_fcdata
*fcd
= nullptr;
641 real ewaldcoeff_q
= 0;
642 real ewaldcoeff_lj
= 0;
643 int nChargePerturbed
= -1, nTypePerturbed
= 0;
644 gmx_wallcycle_t wcycle
;
645 gmx_walltime_accounting_t walltime_accounting
= nullptr;
646 gmx_membed_t
* membed
= nullptr;
647 gmx_hw_info_t
*hwinfo
= nullptr;
649 /* CAUTION: threads may be started later on in this function, so
650 cr doesn't reflect the final parallel state right now */
653 bool doMembed
= opt2bSet("-membed", filenames
.size(), filenames
.data());
654 bool doRerun
= mdrunOptions
.rerun
;
656 // Handle task-assignment related user options.
657 EmulateGpuNonbonded emulateGpuNonbonded
= (getenv("GMX_EMULATE_GPU") != nullptr ?
658 EmulateGpuNonbonded::Yes
: EmulateGpuNonbonded::No
);
660 std::vector
<int> userGpuTaskAssignment
;
663 userGpuTaskAssignment
= parseUserTaskAssignmentString(hw_opt
.userGpuTaskAssignment
);
665 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
666 auto nonbondedTarget
= findTaskTarget(nbpu_opt
);
667 auto pmeTarget
= findTaskTarget(pme_opt
);
668 auto pmeFftTarget
= findTaskTarget(pme_fft_opt
);
669 auto bondedTarget
= findTaskTarget(bonded_opt
);
670 auto updateTarget
= findTaskTarget(update_opt
);
671 PmeRunMode pmeRunMode
= PmeRunMode::None
;
673 FILE *fplog
= nullptr;
674 // If we are appending, we don't write log output because we need
675 // to check that the old log file matches what the checkpoint file
676 // expects. Otherwise, we should start to write log output now if
677 // there is a file ready for it.
678 if (logFileHandle
!= nullptr && startingBehavior
!= StartingBehavior::RestartWithAppending
)
680 fplog
= gmx_fio_getfp(logFileHandle
);
682 const bool isSimulationMasterRank
= findIsSimulationMasterRank(ms
, communicator
);
683 gmx::LoggerOwner
logOwner(buildLogger(fplog
, isSimulationMasterRank
));
684 gmx::MDLogger
mdlog(logOwner
.logger());
686 // report any development features that may be enabled by environment variables
687 const DevelopmentFeatureFlags devFlags
= manageDevelopmentFeatures(mdlog
);
689 // TODO The thread-MPI master rank makes a working
690 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
691 // after the threads have been launched. This works because no use
692 // is made of that communicator until after the execution paths
693 // have rejoined. But it is likely that we can improve the way
694 // this is expressed, e.g. by expressly running detection only the
695 // master rank for thread-MPI, rather than relying on the mutex
696 // and reference count.
697 PhysicalNodeCommunicator
physicalNodeComm(communicator
, gmx_physicalnode_id_hash());
698 hwinfo
= gmx_detect_hardware(mdlog
, physicalNodeComm
);
700 gmx_print_detected_hardware(fplog
, isSimulationMasterRank
&& isMasterSim(ms
), mdlog
, hwinfo
);
702 std::vector
<int> gpuIdsToUse
= makeGpuIdsToUse(hwinfo
->gpu_info
, hw_opt
.gpuIdsAvailable
);
704 // Print citation requests after all software/hardware printing
705 pleaseCiteGromacs(fplog
);
707 // TODO Replace this by unique_ptr once t_inputrec is C++
708 t_inputrec inputrecInstance
;
709 t_inputrec
*inputrec
= nullptr;
710 std::unique_ptr
<t_state
> globalState
;
712 auto partialDeserializedTpr
= std::make_unique
<PartialDeserializedTprFile
>();
714 if (isSimulationMasterRank
)
716 /* Only the master rank has the global state */
717 globalState
= std::make_unique
<t_state
>();
719 /* Read (nearly) all data required for the simulation
720 * and keep the partly serialized tpr contents to send to other ranks later
722 *partialDeserializedTpr
= read_tpx_state(ftp2fn(efTPR
, filenames
.size(), filenames
.data()), &inputrecInstance
, globalState
.get(), &mtop
);
723 inputrec
= &inputrecInstance
;
726 /* Check and update the hardware options for internal consistency */
727 checkAndUpdateHardwareOptions(mdlog
, &hw_opt
, isSimulationMasterRank
, domdecOptions
.numPmeRanks
, inputrec
);
729 if (GMX_THREAD_MPI
&& isSimulationMasterRank
)
731 bool useGpuForNonbonded
= false;
732 bool useGpuForPme
= false;
735 GMX_RELEASE_ASSERT(inputrec
!= nullptr, "Keep the compiler happy");
737 // If the user specified the number of ranks, then we must
738 // respect that, but in default mode, we need to allow for
739 // the number of GPUs to choose the number of ranks.
740 auto canUseGpuForNonbonded
= buildSupportsNonbondedOnGpu(nullptr);
741 useGpuForNonbonded
= decideWhetherToUseGpusForNonbondedWithThreadMpi
742 (nonbondedTarget
, gpuIdsToUse
, userGpuTaskAssignment
, emulateGpuNonbonded
,
743 canUseGpuForNonbonded
,
744 gpuAccelerationOfNonbondedIsUseful(mdlog
, *inputrec
, GMX_THREAD_MPI
),
745 hw_opt
.nthreads_tmpi
);
746 useGpuForPme
= decideWhetherToUseGpusForPmeWithThreadMpi
747 (useGpuForNonbonded
, pmeTarget
, gpuIdsToUse
, userGpuTaskAssignment
,
748 *hwinfo
, *inputrec
, mtop
, hw_opt
.nthreads_tmpi
, domdecOptions
.numPmeRanks
);
751 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
753 /* Determine how many thread-MPI ranks to start.
755 * TODO Over-writing the user-supplied value here does
756 * prevent any possible subsequent checks from working
758 hw_opt
.nthreads_tmpi
= get_nthreads_mpi(hwinfo
,
767 // Now start the threads for thread MPI.
768 spawnThreads(hw_opt
.nthreads_tmpi
);
769 // The spawned threads enter mdrunner() and execution of
770 // master and spawned threads joins at the end of this block.
771 physicalNodeComm
= PhysicalNodeCommunicator(communicator
, gmx_physicalnode_id_hash());
774 GMX_RELEASE_ASSERT(communicator
== MPI_COMM_WORLD
, "Must have valid world communicator");
775 CommrecHandle crHandle
= init_commrec(communicator
, ms
);
776 t_commrec
*cr
= crHandle
.get();
777 GMX_RELEASE_ASSERT(cr
!= nullptr, "Must have valid commrec");
781 /* now broadcast everything to the non-master nodes/threads: */
782 if (!isSimulationMasterRank
)
784 inputrec
= &inputrecInstance
;
786 init_parallel(cr
, inputrec
, &mtop
, partialDeserializedTpr
.get());
788 GMX_RELEASE_ASSERT(inputrec
!= nullptr, "All ranks should have a valid inputrec now");
789 partialDeserializedTpr
.reset(nullptr);
791 // Now the number of ranks is known to all ranks, and each knows
792 // the inputrec read by the master rank. The ranks can now all run
793 // the task-deciding functions and will agree on the result
794 // without needing to communicate.
796 // TODO Should we do the communication in debug mode to support
797 // having an assertion?
799 // Note that these variables describe only their own node.
801 // Note that when bonded interactions run on a GPU they always run
802 // alongside a nonbonded task, so do not influence task assignment
803 // even though they affect the force calculation workload.
804 bool useGpuForNonbonded
= false;
805 bool useGpuForPme
= false;
806 bool useGpuForBonded
= false;
807 bool useGpuForUpdate
= false;
808 bool gpusWereDetected
= hwinfo
->ngpu_compatible_tot
> 0;
811 // It's possible that there are different numbers of GPUs on
812 // different nodes, which is the user's responsibility to
813 // handle. If unsuitable, we will notice that during task
815 auto canUseGpuForNonbonded
= buildSupportsNonbondedOnGpu(nullptr);
816 useGpuForNonbonded
= decideWhetherToUseGpusForNonbonded(nonbondedTarget
, userGpuTaskAssignment
,
818 canUseGpuForNonbonded
,
819 gpuAccelerationOfNonbondedIsUseful(mdlog
, *inputrec
, !GMX_THREAD_MPI
),
821 useGpuForPme
= decideWhetherToUseGpusForPme(useGpuForNonbonded
, pmeTarget
, userGpuTaskAssignment
,
822 *hwinfo
, *inputrec
, mtop
,
823 cr
->nnodes
, domdecOptions
.numPmeRanks
,
825 auto canUseGpuForBonded
= buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec
, mtop
, nullptr);
827 decideWhetherToUseGpusForBonded(useGpuForNonbonded
, useGpuForPme
,
828 bondedTarget
, canUseGpuForBonded
,
829 EVDW_PME(inputrec
->vdwtype
),
830 EEL_PME_EWALD(inputrec
->coulombtype
),
831 domdecOptions
.numPmeRanks
, gpusWereDetected
);
833 pmeRunMode
= (useGpuForPme
? PmeRunMode::GPU
: PmeRunMode::CPU
);
834 if (pmeRunMode
== PmeRunMode::GPU
)
836 if (pmeFftTarget
== TaskTarget::Cpu
)
838 pmeRunMode
= PmeRunMode::Mixed
;
841 else if (pmeFftTarget
== TaskTarget::Gpu
)
843 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.");
846 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
849 // TODO: hide restraint implementation details from Mdrunner.
850 // There is nothing unique about restraints at this point as far as the
851 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
852 // factory functions from the SimulationContext on which to call mdModules_->add().
853 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
854 for (auto && restraint
: restraintManager_
->getRestraints())
856 auto module
= RestraintMDModule::create(restraint
,
858 mdModules_
->add(std::move(module
));
861 // TODO: Error handling
862 mdModules_
->assignOptionsToModules(*inputrec
->params
, nullptr);
863 const auto &mdModulesNotifier
= mdModules_
->notifier().notifier_
;
865 if (inputrec
->internalParameters
!= nullptr)
867 mdModulesNotifier
.notify(*inputrec
->internalParameters
);
870 if (fplog
!= nullptr)
872 pr_inputrec(fplog
, 0, "Input Parameters", inputrec
, FALSE
);
873 fprintf(fplog
, "\n");
878 /* In rerun, set velocities to zero if present */
879 if (doRerun
&& ((globalState
->flags
& (1 << estV
)) != 0))
881 // rerun does not use velocities
882 GMX_LOG(mdlog
.info
).asParagraph().appendText(
883 "Rerun trajectory contains velocities. Rerun does only evaluate "
884 "potential energy and forces. The velocities will be ignored.");
885 for (int i
= 0; i
< globalState
->natoms
; i
++)
887 clear_rvec(globalState
->v
[i
]);
889 globalState
->flags
&= ~(1 << estV
);
892 /* now make sure the state is initialized and propagated */
893 set_state_entries(globalState
.get(), inputrec
);
896 /* NM and TPI parallelize over force/energy calculations, not atoms,
897 * so we need to initialize and broadcast the global state.
899 if (inputrec
->eI
== eiNM
|| inputrec
->eI
== eiTPI
)
903 globalState
= std::make_unique
<t_state
>();
905 broadcastStateWithoutDynamics(cr
, globalState
.get());
908 /* A parallel command line option consistency check that we can
909 only do after any threads have started. */
910 if (!PAR(cr
) && (domdecOptions
.numCells
[XX
] > 1 ||
911 domdecOptions
.numCells
[YY
] > 1 ||
912 domdecOptions
.numCells
[ZZ
] > 1 ||
913 domdecOptions
.numPmeRanks
> 0))
916 "The -dd or -npme option request a parallel simulation, "
918 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv
));
920 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
922 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv
));
927 (EI_ENERGY_MINIMIZATION(inputrec
->eI
) || eiNM
== inputrec
->eI
))
929 gmx_fatal(FARGS
, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
932 if (!(EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
)))
934 if (domdecOptions
.numPmeRanks
> 0)
936 gmx_fatal_collective(FARGS
, cr
->mpi_comm_mysim
, MASTER(cr
),
937 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
940 domdecOptions
.numPmeRanks
= 0;
943 if (useGpuForNonbonded
&& domdecOptions
.numPmeRanks
< 0)
945 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
946 * improve performance with many threads per GPU, since our OpenMP
947 * scaling is bad, but it's difficult to automate the setup.
949 domdecOptions
.numPmeRanks
= 0;
953 if (domdecOptions
.numPmeRanks
< 0)
955 domdecOptions
.numPmeRanks
= 0;
956 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
960 GMX_RELEASE_ASSERT(domdecOptions
.numPmeRanks
<= 1, "PME GPU decomposition is not supported");
967 fcRegisterSteps(inputrec
->nsteps
, inputrec
->init_step
);
971 /* NMR restraints must be initialized before load_checkpoint,
972 * since with time averaging the history is added to t_state.
973 * For proper consistency check we therefore need to extend
975 * So the PME-only nodes (if present) will also initialize
976 * the distance restraints.
980 /* This needs to be called before read_checkpoint to extend the state */
981 init_disres(fplog
, &mtop
, inputrec
, cr
, ms
, fcd
, globalState
.get(), replExParams
.exchangeInterval
> 0);
983 init_orires(fplog
, &mtop
, inputrec
, cr
, ms
, globalState
.get(), &(fcd
->orires
));
985 auto deform
= prepareBoxDeformation(globalState
->box
, cr
, *inputrec
);
987 ObservablesHistory observablesHistory
= {};
989 if (startingBehavior
!= StartingBehavior::NewSimulation
)
991 /* Check if checkpoint file exists before doing continuation.
992 * This way we can use identical input options for the first and subsequent runs...
994 if (mdrunOptions
.numStepsCommandline
> -2)
996 /* Temporarily set the number of steps to unlimited to avoid
997 * triggering the nsteps check in load_checkpoint().
998 * This hack will go away soon when the -nsteps option is removed.
1000 inputrec
->nsteps
= -1;
1003 load_checkpoint(opt2fn_master("-cpi", filenames
.size(), filenames
.data(), cr
),
1005 cr
, domdecOptions
.numCells
,
1006 inputrec
, globalState
.get(),
1007 &observablesHistory
,
1008 mdrunOptions
.reproducible
, mdModules_
->notifier());
1010 if (startingBehavior
== StartingBehavior::RestartWithAppending
&& logFileHandle
)
1012 // Now we can start normal logging to the truncated log file.
1013 fplog
= gmx_fio_getfp(logFileHandle
);
1014 prepareLogAppending(fplog
);
1015 logOwner
= buildLogger(fplog
, MASTER(cr
));
1016 mdlog
= logOwner
.logger();
1020 if (mdrunOptions
.numStepsCommandline
> -2)
1022 GMX_LOG(mdlog
.info
).asParagraph().
1023 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
1024 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
1026 /* override nsteps with value set on the commandline */
1027 override_nsteps_cmdline(mdlog
, mdrunOptions
.numStepsCommandline
, inputrec
);
1031 copy_mat(globalState
->box
, box
);
1036 gmx_bcast(sizeof(box
), box
, cr
);
1039 if (inputrec
->cutoff_scheme
!= ecutsVERLET
)
1041 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.");
1043 /* Update rlist and nstlist. */
1044 prepare_verlet_scheme(fplog
, cr
, inputrec
, nstlist_cmdline
, &mtop
, box
,
1045 useGpuForNonbonded
|| (emulateGpuNonbonded
== EmulateGpuNonbonded::Yes
), *hwinfo
->cpuInfo
);
1047 const bool prefer1DAnd1PulseDD
= (devFlags
.enableGpuHaloExchange
&& useGpuForNonbonded
);
1048 // This builder is necessary while we have multi-part construction
1049 // of DD. Before DD is constructed, we use the existence of
1050 // the builder object to indicate that further construction of DD
1052 std::unique_ptr
<DomainDecompositionBuilder
> ddBuilder
;
1053 if (PAR(cr
) && !(EI_TPI(inputrec
->eI
) ||
1054 inputrec
->eI
== eiNM
))
1056 ddBuilder
= std::make_unique
<DomainDecompositionBuilder
>
1057 (mdlog
, cr
, domdecOptions
, mdrunOptions
,
1058 prefer1DAnd1PulseDD
, mtop
, *inputrec
,
1059 box
, positionsFromStatePointer(globalState
.get()));
1063 /* PME, if used, is done on all nodes with 1D decomposition */
1065 cr
->duty
= (DUTY_PP
| DUTY_PME
);
1067 if (inputrec
->ePBC
== epbcSCREW
)
1070 "pbc=screw is only implemented with domain decomposition");
1074 // Produce the task assignment for this rank.
1075 GpuTaskAssignmentsBuilder gpuTaskAssignmentsBuilder
;
1076 GpuTaskAssignments gpuTaskAssignments
=
1077 gpuTaskAssignmentsBuilder
.build(gpuIdsToUse
,
1078 userGpuTaskAssignment
,
1089 thisRankHasDuty(cr
, DUTY_PP
),
1090 // TODO cr->duty & DUTY_PME should imply that a PME
1091 // algorithm is active, but currently does not.
1092 EEL_PME(inputrec
->coulombtype
) &&
1093 thisRankHasDuty(cr
, DUTY_PME
));
1095 const bool printHostName
= (cr
->nnodes
> 1);
1096 gpuTaskAssignments
.reportGpuUsage(mdlog
, printHostName
, useGpuForBonded
, pmeRunMode
);
1098 // If the user chose a task assignment, give them some hints
1099 // where appropriate.
1100 if (!userGpuTaskAssignment
.empty())
1102 gpuTaskAssignments
.logPerformanceHints(mdlog
,
1103 ssize(gpuIdsToUse
));
1106 // Get the device handles for the modules, nullptr when no task is assigned.
1107 gmx_device_info_t
*nonbondedDeviceInfo
= gpuTaskAssignments
.initNonbondedDevice(cr
);
1108 gmx_device_info_t
*pmeDeviceInfo
= gpuTaskAssignments
.initPmeDevice();
1110 // TODO Initialize GPU streams here.
1112 // TODO Currently this is always built, yet DD partition code
1113 // checks if it is built before using it. Probably it should
1114 // become an MDModule that is made only when another module
1115 // requires it (e.g. pull, CompEl, density fitting), so that we
1116 // don't update the local atom sets unilaterally every step.
1117 LocalAtomSetManager atomSets
;
1120 // TODO Pass the GPU streams to ddBuilder to use in buffer
1121 // transfers (e.g. halo exchange)
1122 cr
->dd
= ddBuilder
->build(&atomSets
);
1123 // The builder's job is done, so destruct it
1124 ddBuilder
.reset(nullptr);
1125 // Note that local state still does not exist yet.
1130 /* After possible communicator splitting in make_dd_communicators.
1131 * we can set up the intra/inter node communication.
1133 gmx_setup_nodecomm(fplog
, cr
);
1139 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1140 "This is simulation %d out of %d running as a composite GROMACS\n"
1141 "multi-simulation job. Setup for this simulation:\n",
1144 GMX_LOG(mdlog
.warning
).appendTextFormatted(
1145 "Using %d MPI %s\n",
1148 cr
->nnodes
== 1 ? "thread" : "threads"
1150 cr
->nnodes
== 1 ? "process" : "processes"
1156 // If mdrun -pin auto honors any affinity setting that already
1157 // exists. If so, it is nice to provide feedback about whether
1158 // that existing affinity setting was from OpenMP or something
1159 // else, so we run this code both before and after we initialize
1160 // the OpenMP support.
1161 gmx_check_thread_affinity_set(mdlog
,
1162 &hw_opt
, hwinfo
->nthreads_hw_avail
, FALSE
);
1163 /* Check and update the number of OpenMP threads requested */
1164 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt
, *hwinfo
, cr
, ms
, physicalNodeComm
.size_
,
1165 pmeRunMode
, mtop
, *inputrec
);
1167 gmx_omp_nthreads_init(mdlog
, cr
,
1168 hwinfo
->nthreads_hw_avail
,
1169 physicalNodeComm
.size_
,
1170 hw_opt
.nthreads_omp
,
1171 hw_opt
.nthreads_omp_pme
,
1172 !thisRankHasDuty(cr
, DUTY_PP
));
1174 // Enable FP exception detection, but not in
1175 // Release mode and not for compilers with known buggy FP
1176 // exception support (clang with any optimization) or suspected
1177 // buggy FP exception support (gcc 7.* with optimization).
1178 #if !defined NDEBUG && \
1179 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
1180 && defined __OPTIMIZE__)
1181 const bool bEnableFPE
= true;
1183 const bool bEnableFPE
= false;
1185 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
1188 gmx_feenableexcept();
1191 /* Now that we know the setup is consistent, check for efficiency */
1192 check_resource_division_efficiency(hwinfo
,
1193 gpuTaskAssignments
.thisRankHasAnyGpuTask(),
1194 mdrunOptions
.ntompOptionIsSet
,
1198 /* getting number of PP/PME threads on this MPI / tMPI rank.
1199 PME: env variable should be read only on one node to make sure it is
1200 identical everywhere;
1202 const int numThreadsOnThisRank
=
1203 thisRankHasDuty(cr
, DUTY_PP
) ? gmx_omp_nthreads_get(emntNonbonded
) : gmx_omp_nthreads_get(emntPME
);
1204 checkHardwareOversubscription(numThreadsOnThisRank
, cr
->nodeid
,
1205 *hwinfo
->hardwareTopology
,
1206 physicalNodeComm
, mdlog
);
1208 if (hw_opt
.threadAffinity
!= ThreadAffinity::Off
)
1210 /* Before setting affinity, check whether the affinity has changed
1211 * - which indicates that probably the OpenMP library has changed it
1212 * since we first checked).
1214 gmx_check_thread_affinity_set(mdlog
,
1215 &hw_opt
, hwinfo
->nthreads_hw_avail
, TRUE
);
1217 int numThreadsOnThisNode
, intraNodeThreadOffset
;
1218 analyzeThreadsOnThisNode(physicalNodeComm
, numThreadsOnThisRank
, &numThreadsOnThisNode
,
1219 &intraNodeThreadOffset
);
1221 /* Set the CPU affinity */
1222 gmx_set_thread_affinity(mdlog
, cr
, &hw_opt
, *hwinfo
->hardwareTopology
,
1223 numThreadsOnThisRank
, numThreadsOnThisNode
,
1224 intraNodeThreadOffset
, nullptr);
1227 if (mdrunOptions
.timingOptions
.resetStep
> -1)
1229 GMX_LOG(mdlog
.info
).asParagraph().
1230 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1232 wcycle
= wallcycle_init(fplog
, mdrunOptions
.timingOptions
.resetStep
, cr
);
1236 /* Master synchronizes its value of reset_counters with all nodes
1237 * including PME only nodes */
1238 int64_t reset_counters
= wcycle_get_reset_counters(wcycle
);
1239 gmx_bcast_sim(sizeof(reset_counters
), &reset_counters
, cr
);
1240 wcycle_set_reset_counters(wcycle
, reset_counters
);
1243 // Membrane embedding must be initialized before we call init_forcerec()
1248 fprintf(stderr
, "Initializing membed");
1250 /* Note that membed cannot work in parallel because mtop is
1251 * changed here. Fix this if we ever want to make it run with
1252 * multiple ranks. */
1253 membed
= init_membed(fplog
, filenames
.size(), filenames
.data(), &mtop
, inputrec
, globalState
.get(), cr
,
1255 .checkpointOptions
.period
);
1258 const bool thisRankHasPmeGpuTask
= gpuTaskAssignments
.thisRankHasPmeGpuTask();
1259 std::unique_ptr
<MDAtoms
> mdAtoms
;
1260 std::unique_ptr
<gmx_vsite_t
> vsite
;
1263 if (thisRankHasDuty(cr
, DUTY_PP
))
1265 mdModulesNotifier
.notify(*cr
);
1266 mdModulesNotifier
.notify(&atomSets
);
1267 mdModulesNotifier
.notify(PeriodicBoundaryConditionType
{inputrec
->ePBC
});
1268 /* Initiate forcerecord */
1269 fr
= new t_forcerec
;
1270 fr
->forceProviders
= mdModules_
->initForceProviders();
1271 init_forcerec(fplog
, mdlog
, fr
, fcd
,
1272 inputrec
, &mtop
, cr
, box
,
1273 opt2fn("-table", filenames
.size(), filenames
.data()),
1274 opt2fn("-tablep", filenames
.size(), filenames
.data()),
1275 opt2fns("-tableb", filenames
.size(), filenames
.data()),
1276 *hwinfo
, nonbondedDeviceInfo
,
1281 // TODO Move this to happen during domain decomposition setup,
1282 // once stream and event handling works well with that.
1283 // TODO remove need to pass local stream into GPU halo exchange - Redmine #3093
1284 if (havePPDomainDecomposition(cr
) && prefer1DAnd1PulseDD
&& is1DAnd1PulseDD(*cr
->dd
))
1286 GMX_RELEASE_ASSERT(devFlags
.enableGpuBufferOps
, "Must use GMX_GPU_BUFFER_OPS=1 to use GMX_GPU_DD_COMMS=1");
1287 void *streamLocal
= Nbnxm::gpu_get_command_stream(fr
->nbv
->gpu_nbv
, Nbnxm::InteractionLocality::Local
);
1288 void *streamNonLocal
= Nbnxm::gpu_get_command_stream(fr
->nbv
->gpu_nbv
, Nbnxm::InteractionLocality::NonLocal
);
1289 void *coordinatesOnDeviceEvent
= fr
->nbv
->get_x_on_device_event();
1290 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
1291 "NOTE: This run uses the 'GPU halo exchange' feature, enabled by the GMX_GPU_DD_COMMS environment variable.");
1292 cr
->dd
->gpuHaloExchange
= std::make_unique
<GpuHaloExchange
>(cr
->dd
, cr
->mpi_comm_mysim
, streamLocal
,
1293 streamNonLocal
, coordinatesOnDeviceEvent
);
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 // TODO should live in ewald module once its testing is improved
1360 // Later, this program could contain kernels that might be later
1361 // re-used as auto-tuning progresses, or subsequent simulations
1363 PmeGpuProgramStorage pmeGpuProgram
;
1364 if (thisRankHasPmeGpuTask
)
1366 pmeGpuProgram
= buildPmeGpuProgram(pmeDeviceInfo
);
1369 /* Initiate PME if necessary,
1370 * either on all nodes or on dedicated PME nodes only. */
1371 if (EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
))
1373 if (mdAtoms
&& mdAtoms
->mdatoms())
1375 nChargePerturbed
= mdAtoms
->mdatoms()->nChargePerturbed
;
1376 if (EVDW_PME(inputrec
->vdwtype
))
1378 nTypePerturbed
= mdAtoms
->mdatoms()->nTypePerturbed
;
1381 if (cr
->npmenodes
> 0)
1383 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1384 gmx_bcast_sim(sizeof(nChargePerturbed
), &nChargePerturbed
, cr
);
1385 gmx_bcast_sim(sizeof(nTypePerturbed
), &nTypePerturbed
, cr
);
1388 if (thisRankHasDuty(cr
, DUTY_PME
))
1392 pmedata
= gmx_pme_init(cr
,
1393 getNumPmeDomains(cr
->dd
),
1395 nChargePerturbed
!= 0, nTypePerturbed
!= 0,
1396 mdrunOptions
.reproducible
,
1397 ewaldcoeff_q
, ewaldcoeff_lj
,
1398 gmx_omp_nthreads_get(emntPME
),
1399 pmeRunMode
, nullptr,
1400 pmeDeviceInfo
, pmeGpuProgram
.get(), mdlog
);
1402 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1407 if (EI_DYNAMICS(inputrec
->eI
))
1409 /* Turn on signal handling on all nodes */
1411 * (A user signal from the PME nodes (if any)
1412 * is communicated to the PP nodes.
1414 signal_handler_install();
1417 pull_t
*pull_work
= nullptr;
1418 if (thisRankHasDuty(cr
, DUTY_PP
))
1420 /* Assumes uniform use of the number of OpenMP threads */
1421 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntDefault
));
1423 if (inputrec
->bPull
)
1425 /* Initialize pull code */
1427 init_pull(fplog
, inputrec
->pull
, inputrec
,
1428 &mtop
, cr
, &atomSets
, inputrec
->fepvals
->init_lambda
);
1429 if (inputrec
->pull
->bXOutAverage
|| inputrec
->pull
->bFOutAverage
)
1431 initPullHistory(pull_work
, &observablesHistory
);
1433 if (EI_DYNAMICS(inputrec
->eI
) && MASTER(cr
))
1435 init_pull_output_files(pull_work
,
1436 filenames
.size(), filenames
.data(), oenv
,
1441 std::unique_ptr
<EnforcedRotation
> enforcedRotation
;
1444 /* Initialize enforced rotation code */
1445 enforcedRotation
= init_rot(fplog
,
1458 t_swap
*swap
= nullptr;
1459 if (inputrec
->eSwapCoords
!= eswapNO
)
1461 /* Initialize ion swapping code */
1462 swap
= init_swapcoords(fplog
, inputrec
,
1463 opt2fn_master("-swap", filenames
.size(), filenames
.data(), cr
),
1464 &mtop
, globalState
.get(), &observablesHistory
,
1465 cr
, &atomSets
, oenv
, mdrunOptions
,
1469 /* Let makeConstraints know whether we have essential dynamics constraints.
1470 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1472 bool doEssentialDynamics
= (opt2fn_null("-ei", filenames
.size(), filenames
.data()) != nullptr
1473 || observablesHistory
.edsamHistory
);
1474 auto constr
= makeConstraints(mtop
, *inputrec
, pull_work
, doEssentialDynamics
,
1475 fplog
, *mdAtoms
->mdatoms(),
1476 cr
, ms
, &nrnb
, wcycle
, fr
->bMolPBC
);
1478 /* Energy terms and groups */
1479 gmx_enerdata_t
enerd(mtop
.groups
.groups
[SimulationAtomGroupType::EnergyOutput
].size(), inputrec
->fepvals
->n_lambda
);
1481 /* Kinetic energy data */
1482 gmx_ekindata_t ekind
;
1483 init_ekindata(fplog
, &mtop
, &(inputrec
->opts
), &ekind
, inputrec
->cos_accel
);
1485 /* Set up interactive MD (IMD) */
1486 auto imdSession
= makeImdSession(inputrec
, cr
, wcycle
, &enerd
, ms
, &mtop
, mdlog
,
1487 MASTER(cr
) ? globalState
->x
.rvec_array() : nullptr,
1488 filenames
.size(), filenames
.data(), oenv
, mdrunOptions
.imdOptions
,
1491 if (DOMAINDECOMP(cr
))
1493 GMX_RELEASE_ASSERT(fr
, "fr was NULL while cr->duty was DUTY_PP");
1494 /* This call is not included in init_domain_decomposition mainly
1495 * because fr->cginfo_mb is set later.
1497 dd_init_bondeds(fplog
, cr
->dd
, &mtop
, vsite
.get(), inputrec
,
1498 domdecOptions
.checkBondedInteractions
,
1502 if (updateTarget
== TaskTarget::Gpu
)
1506 gmx_fatal(FARGS
, "It is currently not possible to redirect the calculation "
1507 "of update and constraints to the GPU!");
1511 // Before we start the actual simulator, try if we can run the update task on the GPU.
1512 useGpuForUpdate
= decideWhetherToUseGpuForUpdate(DOMAINDECOMP(cr
),
1515 devFlags
.enableGpuBufferOps
,
1520 doEssentialDynamics
,
1521 fcd
->orires
.nr
!= 0,
1522 fcd
->disres
.nsystems
!= 0);
1524 const bool inputIsCompatibleWithModularSimulator
= ModularSimulator::isInputCompatible(
1526 inputrec
, doRerun
, vsite
.get(), ms
, replExParams
,
1527 fcd
, static_cast<int>(filenames
.size()), filenames
.data(),
1528 &observablesHistory
, membed
);
1530 const bool useModularSimulator
= inputIsCompatibleWithModularSimulator
&& !(getenv("GMX_DISABLE_MODULAR_SIMULATOR") != nullptr);
1532 std::unique_ptr
<gmx::StatePropagatorDataGpu
> stateGpu
;
1533 if (gpusWereDetected
&& ((useGpuForPme
&& thisRankHasDuty(cr
, DUTY_PME
)) || devFlags
.enableGpuBufferOps
))
1535 const void *pmeStream
= pme_gpu_get_device_stream(fr
->pmedata
);
1536 const void *localStream
= fr
->nbv
->gpu_nbv
!= nullptr ? Nbnxm::gpu_get_command_stream(fr
->nbv
->gpu_nbv
, Nbnxm::InteractionLocality::Local
) : nullptr;
1537 const void *nonLocalStream
= fr
->nbv
->gpu_nbv
!= nullptr ? Nbnxm::gpu_get_command_stream(fr
->nbv
->gpu_nbv
, Nbnxm::InteractionLocality::NonLocal
) : nullptr;
1538 const void *deviceContext
= pme_gpu_get_device_context(fr
->pmedata
);
1539 const int paddingSize
= pme_gpu_get_padding_size(fr
->pmedata
);
1540 GpuApiCallBehavior transferKind
= (inputrec
->eI
== eiMD
&& !doRerun
&& !useModularSimulator
) ? GpuApiCallBehavior::Async
: GpuApiCallBehavior::Sync
;
1542 stateGpu
= std::make_unique
<gmx::StatePropagatorDataGpu
>(pmeStream
,
1548 fr
->stateGpu
= stateGpu
.get();
1551 // TODO This is not the right place to manage the lifetime of
1552 // this data structure, but currently it's the easiest way to
1554 MdrunScheduleWorkload runScheduleWork
;
1555 // Also populates the simulation constant workload description.
1556 runScheduleWork
.simulationWork
= createSimulationWorkload(useGpuForNonbonded
,
1558 (pmeRunMode
== PmeRunMode::GPU
),
1561 devFlags
.enableGpuBufferOps
,
1562 devFlags
.enableGpuHaloExchange
,
1563 devFlags
.enableGpuPmePPComm
);
1566 GMX_ASSERT(stopHandlerBuilder_
, "Runner must provide StopHandlerBuilder to simulator.");
1567 SimulatorBuilder simulatorBuilder
;
1569 // build and run simulator object based on user-input
1570 auto simulator
= simulatorBuilder
.build(
1571 inputIsCompatibleWithModularSimulator
,
1572 fplog
, cr
, ms
, mdlog
, static_cast<int>(filenames
.size()), filenames
.data(),
1576 vsite
.get(), constr
.get(),
1577 enforcedRotation
? enforcedRotation
->getLegacyEnfrot() : nullptr,
1579 mdModules_
->outputProvider(),
1580 mdModules_
->notifier(),
1581 inputrec
, imdSession
.get(), pull_work
, swap
, &mtop
,
1584 &observablesHistory
,
1585 mdAtoms
.get(), &nrnb
, wcycle
, fr
,
1591 walltime_accounting
,
1592 std::move(stopHandlerBuilder_
),
1596 if (inputrec
->bPull
)
1598 finish_pull(pull_work
);
1600 finish_swapcoords(swap
);
1604 GMX_RELEASE_ASSERT(pmedata
, "pmedata was NULL while cr->duty was not DUTY_PP");
1606 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntPME
));
1607 gmx_pmeonly(pmedata
, cr
, &nrnb
, wcycle
, walltime_accounting
, inputrec
, pmeRunMode
);
1610 wallcycle_stop(wcycle
, ewcRUN
);
1612 /* Finish up, write some stuff
1613 * if rerunMD, don't write last frame again
1615 finish_run(fplog
, mdlog
, cr
,
1616 inputrec
, &nrnb
, wcycle
, walltime_accounting
,
1617 fr
? fr
->nbv
.get() : nullptr,
1619 EI_DYNAMICS(inputrec
->eI
) && !isMultiSim(ms
));
1621 // clean up cycle counter
1622 wallcycle_destroy(wcycle
);
1627 gmx_pme_destroy(pmedata
);
1631 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1632 // before we destroy the GPU context(s) in free_gpu_resources().
1633 // Pinned buffers are associated with contexts in CUDA.
1634 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1635 mdAtoms
.reset(nullptr);
1636 globalState
.reset(nullptr);
1637 mdModules_
.reset(nullptr); // destruct force providers here as they might also use the GPU
1639 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1640 free_gpu_resources(fr
, physicalNodeComm
, hwinfo
->gpu_info
);
1641 free_gpu(nonbondedDeviceInfo
);
1642 free_gpu(pmeDeviceInfo
);
1643 done_forcerec(fr
, mtop
.molblock
.size());
1648 free_membed(membed
);
1651 /* Does what it says */
1652 print_date_and_time(fplog
, cr
->nodeid
, "Finished mdrun", gmx_gettime());
1653 walltime_accounting_destroy(walltime_accounting
);
1655 // Ensure log file content is written
1658 gmx_fio_flush(logFileHandle
);
1661 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1662 * exceptions were enabled before function was called. */
1665 gmx_fedisableexcept();
1668 auto rc
= static_cast<int>(gmx_get_stop_condition());
1671 /* we need to join all threads. The sub-threads join when they
1672 exit this function, but the master thread needs to be told to
1674 if (PAR(cr
) && MASTER(cr
))
1682 Mdrunner::~Mdrunner()
1684 // Clean up of the Manager.
1685 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1686 // but okay as long as threads synchronize some time before adding or accessing
1687 // a new set of restraints.
1688 if (restraintManager_
)
1690 restraintManager_
->clear();
1691 GMX_ASSERT(restraintManager_
->countRestraints() == 0,
1692 "restraints added during runner life time should be cleared at runner destruction.");
1696 void Mdrunner::addPotential(std::shared_ptr
<gmx::IRestraintPotential
> puller
,
1697 const std::string
&name
)
1699 GMX_ASSERT(restraintManager_
, "Mdrunner must have a restraint manager.");
1700 // Not sure if this should be logged through the md logger or something else,
1701 // but it is helpful to have some sort of INFO level message sent somewhere.
1702 // std::cout << "Registering restraint named " << name << std::endl;
1704 // When multiple restraints are used, it may be wasteful to register them separately.
1705 // Maybe instead register an entire Restraint Manager as a force provider.
1706 restraintManager_
->addToSpec(std::move(puller
),
1710 Mdrunner::Mdrunner(std::unique_ptr
<MDModules
> mdModules
)
1711 : mdModules_(std::move(mdModules
))
1715 Mdrunner::Mdrunner(Mdrunner
&&) noexcept
= default;
1717 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1718 Mdrunner
&Mdrunner::operator=(Mdrunner
&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING
) = default;
1720 class Mdrunner::BuilderImplementation
1723 BuilderImplementation() = delete;
1724 BuilderImplementation(std::unique_ptr
<MDModules
> mdModules
,
1725 compat::not_null
<SimulationContext
*> context
);
1726 ~BuilderImplementation();
1728 BuilderImplementation
&setExtraMdrunOptions(const MdrunOptions
&options
,
1729 real forceWarningThreshold
,
1730 StartingBehavior startingBehavior
);
1732 void addDomdec(const DomdecOptions
&options
);
1734 void addVerletList(int nstlist
);
1736 void addReplicaExchange(const ReplicaExchangeParameters
¶ms
);
1738 void addNonBonded(const char* nbpu_opt
);
1740 void addPME(const char* pme_opt_
, const char* pme_fft_opt_
);
1742 void addBondedTaskAssignment(const char* bonded_opt
);
1744 void addUpdateTaskAssignment(const char* update_opt
);
1746 void addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
);
1748 void addFilenames(ArrayRef
<const t_filenm
> filenames
);
1750 void addOutputEnvironment(gmx_output_env_t
* outputEnvironment
);
1752 void addLogFile(t_fileio
*logFileHandle
);
1754 void addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
);
1760 // Default parameters copied from runner.h
1761 // \todo Clarify source(s) of default parameters.
1763 const char* nbpu_opt_
= nullptr;
1764 const char* pme_opt_
= nullptr;
1765 const char* pme_fft_opt_
= nullptr;
1766 const char *bonded_opt_
= nullptr;
1767 const char *update_opt_
= nullptr;
1769 MdrunOptions mdrunOptions_
;
1771 DomdecOptions domdecOptions_
;
1773 ReplicaExchangeParameters replicaExchangeParameters_
;
1775 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1778 //! Multisim communicator handle.
1779 gmx_multisim_t
*multiSimulation_
;
1781 //! mdrun communicator
1782 MPI_Comm communicator_
= MPI_COMM_NULL
;
1784 //! Print a warning if any force is larger than this (in kJ/mol nm).
1785 real forceWarningThreshold_
= -1;
1787 //! Whether the simulation will start afresh, or restart with/without appending.
1788 StartingBehavior startingBehavior_
= StartingBehavior::NewSimulation
;
1790 //! The modules that comprise the functionality of mdrun.
1791 std::unique_ptr
<MDModules
> mdModules_
;
1793 //! \brief Parallelism information.
1794 gmx_hw_opt_t hardwareOptions_
;
1796 //! filename options for simulation.
1797 ArrayRef
<const t_filenm
> filenames_
;
1799 /*! \brief Handle to output environment.
1801 * \todo gmx_output_env_t needs lifetime management.
1803 gmx_output_env_t
* outputEnvironment_
= nullptr;
1805 /*! \brief Non-owning handle to MD log file.
1807 * \todo Context should own output facilities for client.
1808 * \todo Improve log file handle management.
1810 * Code managing the FILE* relies on the ability to set it to
1811 * nullptr to check whether the filehandle is valid.
1813 t_fileio
* logFileHandle_
= nullptr;
1816 * \brief Builder for simulation stop signal handler.
1818 std::unique_ptr
<StopHandlerBuilder
> stopHandlerBuilder_
= nullptr;
1821 Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr
<MDModules
> mdModules
,
1822 compat::not_null
<SimulationContext
*> context
) :
1823 mdModules_(std::move(mdModules
))
1825 communicator_
= context
->communicator_
;
1826 multiSimulation_
= context
->multiSimulation_
.get();
1829 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1831 Mdrunner::BuilderImplementation
&
1832 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions
&options
,
1833 const real forceWarningThreshold
,
1834 const StartingBehavior startingBehavior
)
1836 mdrunOptions_
= options
;
1837 forceWarningThreshold_
= forceWarningThreshold
;
1838 startingBehavior_
= startingBehavior
;
1842 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions
&options
)
1844 domdecOptions_
= options
;
1847 void Mdrunner::BuilderImplementation::addVerletList(int nstlist
)
1852 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1854 replicaExchangeParameters_
= params
;
1857 Mdrunner
Mdrunner::BuilderImplementation::build()
1859 auto newRunner
= Mdrunner(std::move(mdModules_
));
1861 newRunner
.mdrunOptions
= mdrunOptions_
;
1862 newRunner
.pforce
= forceWarningThreshold_
;
1863 newRunner
.startingBehavior
= startingBehavior_
;
1864 newRunner
.domdecOptions
= domdecOptions_
;
1866 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1867 newRunner
.hw_opt
= hardwareOptions_
;
1869 // No invariant to check. This parameter exists to optionally override other behavior.
1870 newRunner
.nstlist_cmdline
= nstlist_
;
1872 newRunner
.replExParams
= replicaExchangeParameters_
;
1874 newRunner
.filenames
= filenames_
;
1876 newRunner
.communicator
= communicator_
;
1878 // nullptr is a valid value for the multisim handle
1879 newRunner
.ms
= multiSimulation_
;
1881 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1882 // \todo Update sanity checking when output environment has clearly specified invariants.
1883 // Initialization and default values for oenv are not well specified in the current version.
1884 if (outputEnvironment_
)
1886 newRunner
.oenv
= outputEnvironment_
;
1890 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1893 newRunner
.logFileHandle
= logFileHandle_
;
1897 newRunner
.nbpu_opt
= nbpu_opt_
;
1901 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1904 if (pme_opt_
&& pme_fft_opt_
)
1906 newRunner
.pme_opt
= pme_opt_
;
1907 newRunner
.pme_fft_opt
= pme_fft_opt_
;
1911 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1916 newRunner
.bonded_opt
= bonded_opt_
;
1920 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1925 newRunner
.update_opt
= update_opt_
;
1929 GMX_THROW(gmx::APIError("MdrunnerBuilder::addUpdateTaskAssignment() is required before build() "));
1933 newRunner
.restraintManager_
= std::make_unique
<gmx::RestraintManager
>();
1935 if (stopHandlerBuilder_
)
1937 newRunner
.stopHandlerBuilder_
= std::move(stopHandlerBuilder_
);
1941 newRunner
.stopHandlerBuilder_
= std::make_unique
<StopHandlerBuilder
>();
1947 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt
)
1949 nbpu_opt_
= nbpu_opt
;
1952 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt
,
1953 const char* pme_fft_opt
)
1956 pme_fft_opt_
= pme_fft_opt
;
1959 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt
)
1961 bonded_opt_
= bonded_opt
;
1964 void Mdrunner::BuilderImplementation::addUpdateTaskAssignment(const char* update_opt
)
1966 update_opt_
= update_opt
;
1969 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
1971 hardwareOptions_
= hardwareOptions
;
1974 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef
<const t_filenm
> filenames
)
1976 filenames_
= filenames
;
1979 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
1981 outputEnvironment_
= outputEnvironment
;
1984 void Mdrunner::BuilderImplementation::addLogFile(t_fileio
*logFileHandle
)
1986 logFileHandle_
= logFileHandle
;
1989 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
1991 stopHandlerBuilder_
= std::move(builder
);
1994 MdrunnerBuilder::MdrunnerBuilder(std::unique_ptr
<MDModules
> mdModules
,
1995 compat::not_null
<SimulationContext
*> context
) :
1996 impl_
{std::make_unique
<Mdrunner::BuilderImplementation
>(std::move(mdModules
), context
)}
2000 MdrunnerBuilder::~MdrunnerBuilder() = default;
2002 MdrunnerBuilder
&MdrunnerBuilder::addSimulationMethod(const MdrunOptions
&options
,
2003 real forceWarningThreshold
,
2004 const StartingBehavior startingBehavior
)
2006 impl_
->setExtraMdrunOptions(options
, forceWarningThreshold
, startingBehavior
);
2010 MdrunnerBuilder
&MdrunnerBuilder::addDomainDecomposition(const DomdecOptions
&options
)
2012 impl_
->addDomdec(options
);
2016 MdrunnerBuilder
&MdrunnerBuilder::addNeighborList(int nstlist
)
2018 impl_
->addVerletList(nstlist
);
2022 MdrunnerBuilder
&MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
2024 impl_
->addReplicaExchange(params
);
2028 MdrunnerBuilder
&MdrunnerBuilder::addNonBonded(const char* nbpu_opt
)
2030 impl_
->addNonBonded(nbpu_opt
);
2034 MdrunnerBuilder
&MdrunnerBuilder::addElectrostatics(const char* pme_opt
,
2035 const char* pme_fft_opt
)
2037 // The builder method may become more general in the future, but in this version,
2038 // parameters for PME electrostatics are both required and the only parameters
2040 if (pme_opt
&& pme_fft_opt
)
2042 impl_
->addPME(pme_opt
, pme_fft_opt
);
2046 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
2051 MdrunnerBuilder
&MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt
)
2053 impl_
->addBondedTaskAssignment(bonded_opt
);
2057 MdrunnerBuilder
&MdrunnerBuilder::addUpdateTaskAssignment(const char* update_opt
)
2059 impl_
->addUpdateTaskAssignment(update_opt
);
2063 Mdrunner
MdrunnerBuilder::build()
2065 return impl_
->build();
2068 MdrunnerBuilder
&MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
2070 impl_
->addHardwareOptions(hardwareOptions
);
2074 MdrunnerBuilder
&MdrunnerBuilder::addFilenames(ArrayRef
<const t_filenm
> filenames
)
2076 impl_
->addFilenames(filenames
);
2080 MdrunnerBuilder
&MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
2082 impl_
->addOutputEnvironment(outputEnvironment
);
2086 MdrunnerBuilder
&MdrunnerBuilder::addLogFile(t_fileio
*logFileHandle
)
2088 impl_
->addLogFile(logFileHandle
);
2092 MdrunnerBuilder
&MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
2094 impl_
->addStopHandlerBuilder(std::move(builder
));
2098 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder
&&) noexcept
= default;
2100 MdrunnerBuilder
&MdrunnerBuilder::operator=(MdrunnerBuilder
&&) noexcept
= default;