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, 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
58 #include "gromacs/commandline/filenm.h"
59 #include "gromacs/compat/make_unique.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/localatomsetmanager.h"
63 #include "gromacs/ewald/ewald-utils.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme-gpu-program.h"
66 #include "gromacs/fileio/checkpoint.h"
67 #include "gromacs/fileio/gmxfio.h"
68 #include "gromacs/fileio/oenv.h"
69 #include "gromacs/fileio/tpxio.h"
70 #include "gromacs/gmxlib/network.h"
71 #include "gromacs/gmxlib/nrnb.h"
72 #include "gromacs/gpu_utils/clfftinitializer.h"
73 #include "gromacs/gpu_utils/gpu_utils.h"
74 #include "gromacs/hardware/cpuinfo.h"
75 #include "gromacs/hardware/detecthardware.h"
76 #include "gromacs/hardware/printhardware.h"
77 #include "gromacs/listed-forces/disre.h"
78 #include "gromacs/listed-forces/manage-threading.h"
79 #include "gromacs/listed-forces/orires.h"
80 #include "gromacs/math/functions.h"
81 #include "gromacs/math/utilities.h"
82 #include "gromacs/math/vec.h"
83 #include "gromacs/mdlib/boxdeformation.h"
84 #include "gromacs/mdlib/calc_verletbuf.h"
85 #include "gromacs/mdlib/forcerec.h"
86 #include "gromacs/mdlib/gmx_omp_nthreads.h"
87 #include "gromacs/mdlib/makeconstraints.h"
88 #include "gromacs/mdlib/md_support.h"
89 #include "gromacs/mdlib/mdatoms.h"
90 #include "gromacs/mdlib/mdrun.h"
91 #include "gromacs/mdlib/membed.h"
92 #include "gromacs/mdlib/nb_verlet.h"
93 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
94 #include "gromacs/mdlib/nbnxn_search.h"
95 #include "gromacs/mdlib/nbnxn_tuning.h"
96 #include "gromacs/mdlib/ppforceworkload.h"
97 #include "gromacs/mdlib/qmmm.h"
98 #include "gromacs/mdlib/sighandler.h"
99 #include "gromacs/mdlib/sim_util.h"
100 #include "gromacs/mdlib/stophandler.h"
101 #include "gromacs/mdrun/legacymdrunoptions.h"
102 #include "gromacs/mdrun/logging.h"
103 #include "gromacs/mdrun/multisim.h"
104 #include "gromacs/mdrun/simulationcontext.h"
105 #include "gromacs/mdrunutility/mdmodules.h"
106 #include "gromacs/mdrunutility/threadaffinity.h"
107 #include "gromacs/mdtypes/commrec.h"
108 #include "gromacs/mdtypes/fcdata.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/md_enums.h"
111 #include "gromacs/mdtypes/observableshistory.h"
112 #include "gromacs/mdtypes/state.h"
113 #include "gromacs/pbcutil/pbc.h"
114 #include "gromacs/pulling/output.h"
115 #include "gromacs/pulling/pull.h"
116 #include "gromacs/pulling/pull_rotation.h"
117 #include "gromacs/restraint/manager.h"
118 #include "gromacs/restraint/restraintmdmodule.h"
119 #include "gromacs/restraint/restraintpotential.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/taskassignment/decidegpuusage.h"
122 #include "gromacs/taskassignment/resourcedivision.h"
123 #include "gromacs/taskassignment/taskassignment.h"
124 #include "gromacs/taskassignment/usergpuids.h"
125 #include "gromacs/timing/wallcycle.h"
126 #include "gromacs/topology/mtop_util.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basenetwork.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/exceptions.h"
131 #include "gromacs/utility/fatalerror.h"
132 #include "gromacs/utility/filestream.h"
133 #include "gromacs/utility/gmxassert.h"
134 #include "gromacs/utility/gmxmpi.h"
135 #include "gromacs/utility/logger.h"
136 #include "gromacs/utility/loggerbuilder.h"
137 #include "gromacs/utility/physicalnodecommunicator.h"
138 #include "gromacs/utility/pleasecite.h"
139 #include "gromacs/utility/programcontext.h"
140 #include "gromacs/utility/smalloc.h"
141 #include "gromacs/utility/stringutil.h"
143 #include "integrator.h"
144 #include "replicaexchange.h"
147 #include "corewrap.h"
153 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
155 * Used to ensure that the master thread does not modify mdrunner during copy
156 * on the spawned threads. */
157 static void threadMpiMdrunnerAccessBarrier()
160 MPI_Barrier(MPI_COMM_WORLD
);
164 Mdrunner
Mdrunner::cloneOnSpawnedThread() const
166 auto newRunner
= Mdrunner();
168 // All runners in the same process share a restraint manager resource because it is
169 // part of the interface to the client code, which is associated only with the
170 // original thread. Handles to the same resources can be obtained by copy.
172 newRunner
.restraintManager_
= compat::make_unique
<RestraintManager
>(*restraintManager_
);
175 // Copy original cr pointer before master thread can pass the thread barrier
176 newRunner
.cr
= reinitialize_commrec_for_this_thread(cr
);
178 // Copy members of master runner.
179 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
180 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
181 newRunner
.hw_opt
= hw_opt
;
182 newRunner
.filenames
= filenames
;
184 newRunner
.oenv
= oenv
;
185 newRunner
.mdrunOptions
= mdrunOptions
;
186 newRunner
.domdecOptions
= domdecOptions
;
187 newRunner
.nbpu_opt
= nbpu_opt
;
188 newRunner
.pme_opt
= pme_opt
;
189 newRunner
.pme_fft_opt
= pme_fft_opt
;
190 newRunner
.bonded_opt
= bonded_opt
;
191 newRunner
.nstlist_cmdline
= nstlist_cmdline
;
192 newRunner
.replExParams
= replExParams
;
193 newRunner
.pforce
= pforce
;
195 newRunner
.stopHandlerBuilder_
= compat::make_unique
<StopHandlerBuilder
>(*stopHandlerBuilder_
);
197 threadMpiMdrunnerAccessBarrier();
199 GMX_RELEASE_ASSERT(!MASTER(newRunner
.cr
), "cloneOnSpawnedThread should only be called on spawned threads");
204 /*! \brief The callback used for running on spawned threads.
206 * Obtains the pointer to the master mdrunner object from the one
207 * argument permitted to the thread-launch API call, copies it to make
208 * a new runner for this thread, reinitializes necessary data, and
209 * proceeds to the simulation. */
210 static void mdrunner_start_fn(const void *arg
)
214 auto masterMdrunner
= reinterpret_cast<const gmx::Mdrunner
*>(arg
);
215 /* copy the arg list to make sure that it's thread-local. This
216 doesn't copy pointed-to items, of course; fnm, cr and fplog
217 are reset in the call below, all others should be const. */
218 gmx::Mdrunner mdrunner
= masterMdrunner
->cloneOnSpawnedThread();
221 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
225 /*! \brief Start thread-MPI threads.
227 * Called by mdrunner() to start a specific number of threads
228 * (including the main thread) for thread-parallel runs. This in turn
229 * calls mdrunner() for each thread. All options are the same as for
231 t_commrec
*Mdrunner::spawnThreads(int numThreadsToLaunch
) const
234 /* first check whether we even need to start tMPI */
235 if (numThreadsToLaunch
< 2)
241 /* now spawn new threads that start mdrunner_start_fn(), while
242 the main thread returns, we set thread affinity later */
243 if (tMPI_Init_fn(TRUE
, numThreadsToLaunch
, TMPI_AFFINITY_NONE
,
244 mdrunner_start_fn
, static_cast<const void*>(this)) != TMPI_SUCCESS
)
246 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
249 threadMpiMdrunnerAccessBarrier();
251 GMX_UNUSED_VALUE(mdrunner_start_fn
);
254 return reinitialize_commrec_for_this_thread(cr
);
259 /*! \brief Initialize variables for Verlet scheme simulation */
260 static void prepare_verlet_scheme(FILE *fplog
,
264 const gmx_mtop_t
*mtop
,
266 bool makeGpuPairList
,
267 const gmx::CpuInfo
&cpuinfo
)
269 /* For NVE simulations, we will retain the initial list buffer */
270 if (EI_DYNAMICS(ir
->eI
) &&
271 ir
->verletbuf_tol
> 0 &&
272 !(EI_MD(ir
->eI
) && ir
->etc
== etcNO
))
274 /* Update the Verlet buffer size for the current run setup */
276 /* Here we assume SIMD-enabled kernels are being used. But as currently
277 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
278 * and 4x2 gives a larger buffer than 4x4, this is ok.
280 ListSetupType listType
= (makeGpuPairList
? ListSetupType::Gpu
: ListSetupType::CpuSimdWhenSupported
);
281 VerletbufListSetup listSetup
= verletbufGetSafeListSetup(listType
);
284 calc_verlet_buffer_size(mtop
, det(box
), ir
, ir
->nstlist
, ir
->nstlist
- 1, -1, &listSetup
, nullptr, &rlist_new
);
286 if (rlist_new
!= ir
->rlist
)
288 if (fplog
!= nullptr)
290 fprintf(fplog
, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
291 ir
->rlist
, rlist_new
,
292 listSetup
.cluster_size_i
, listSetup
.cluster_size_j
);
294 ir
->rlist
= rlist_new
;
298 if (nstlist_cmdline
> 0 && (!EI_DYNAMICS(ir
->eI
) || ir
->verletbuf_tol
<= 0))
300 gmx_fatal(FARGS
, "Can not set nstlist without %s",
301 !EI_DYNAMICS(ir
->eI
) ? "dynamics" : "verlet-buffer-tolerance");
304 if (EI_DYNAMICS(ir
->eI
))
306 /* Set or try nstlist values */
307 increaseNstlist(fplog
, cr
, ir
, nstlist_cmdline
, mtop
, box
, makeGpuPairList
, cpuinfo
);
311 /*! \brief Override the nslist value in inputrec
313 * with value passed on the command line (if any)
315 static void override_nsteps_cmdline(const gmx::MDLogger
&mdlog
,
316 int64_t nsteps_cmdline
,
321 /* override with anything else than the default -2 */
322 if (nsteps_cmdline
> -2)
324 char sbuf_steps
[STEPSTRSIZE
];
325 char sbuf_msg
[STRLEN
];
327 ir
->nsteps
= nsteps_cmdline
;
328 if (EI_DYNAMICS(ir
->eI
) && nsteps_cmdline
!= -1)
330 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
331 gmx_step_str(nsteps_cmdline
, sbuf_steps
),
332 fabs(nsteps_cmdline
*ir
->delta_t
));
336 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps",
337 gmx_step_str(nsteps_cmdline
, sbuf_steps
));
340 GMX_LOG(mdlog
.warning
).asParagraph().appendText(sbuf_msg
);
342 else if (nsteps_cmdline
< -2)
344 gmx_fatal(FARGS
, "Invalid nsteps value passed on the command line: %" PRId64
,
347 /* Do nothing if nsteps_cmdline == -2 */
353 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
355 * If not, and if a warning may be issued, logs a warning about
356 * falling back to CPU code. With thread-MPI, only the first
357 * call to this function should have \c issueWarning true. */
358 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger
&mdlog
,
359 const t_inputrec
*ir
,
362 if (ir
->opts
.ngener
- ir
->nwall
> 1)
364 /* The GPU code does not support more than one energy group.
365 * If the user requested GPUs explicitly, a fatal error is given later.
369 GMX_LOG(mdlog
.warning
).asParagraph()
370 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
371 "For better performance, run on the GPU without energy groups and then do "
372 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
379 //! Initializes the logger for mdrun.
380 static gmx::LoggerOwner
buildLogger(FILE *fplog
, const t_commrec
*cr
)
382 gmx::LoggerBuilder builder
;
383 if (fplog
!= nullptr)
385 builder
.addTargetFile(gmx::MDLogger::LogLevel::Info
, fplog
);
387 if (cr
== nullptr || SIMMASTER(cr
))
389 builder
.addTargetStream(gmx::MDLogger::LogLevel::Warning
,
390 &gmx::TextOutputFile::standardError());
392 return builder
.build();
395 //! Make a TaskTarget from an mdrun argument string.
396 static TaskTarget
findTaskTarget(const char *optionString
)
398 TaskTarget returnValue
= TaskTarget::Auto
;
400 if (strncmp(optionString
, "auto", 3) == 0)
402 returnValue
= TaskTarget::Auto
;
404 else if (strncmp(optionString
, "cpu", 3) == 0)
406 returnValue
= TaskTarget::Cpu
;
408 else if (strncmp(optionString
, "gpu", 3) == 0)
410 returnValue
= TaskTarget::Gpu
;
414 GMX_ASSERT(false, "Option string should have been checked for sanity already");
420 int Mdrunner::mdrunner()
424 t_forcerec
*fr
= nullptr;
425 t_fcdata
*fcd
= nullptr;
426 real ewaldcoeff_q
= 0;
427 real ewaldcoeff_lj
= 0;
428 int nChargePerturbed
= -1, nTypePerturbed
= 0;
429 gmx_wallcycle_t wcycle
;
430 gmx_walltime_accounting_t walltime_accounting
= nullptr;
432 int64_t reset_counters
;
433 int nthreads_pme
= 1;
434 gmx_membed_t
* membed
= nullptr;
435 gmx_hw_info_t
*hwinfo
= nullptr;
437 /* CAUTION: threads may be started later on in this function, so
438 cr doesn't reflect the final parallel state right now */
439 std::unique_ptr
<gmx::MDModules
> mdModules(new gmx::MDModules
);
440 t_inputrec inputrecInstance
;
441 t_inputrec
*inputrec
= &inputrecInstance
;
444 bool doMembed
= opt2bSet("-membed", filenames
.size(), filenames
.data());
445 bool doRerun
= mdrunOptions
.rerun
;
447 // Handle task-assignment related user options.
448 EmulateGpuNonbonded emulateGpuNonbonded
= (getenv("GMX_EMULATE_GPU") != nullptr ?
449 EmulateGpuNonbonded::Yes
: EmulateGpuNonbonded::No
);
450 std::vector
<int> gpuIdsAvailable
;
453 gpuIdsAvailable
= parseUserGpuIds(hw_opt
.gpuIdsAvailable
);
454 // TODO We could put the GPU IDs into a std::map to find
455 // duplicates, but for the small numbers of IDs involved, this
456 // code is simple and fast.
457 for (size_t i
= 0; i
!= gpuIdsAvailable
.size(); ++i
)
459 for (size_t j
= i
+1; j
!= gpuIdsAvailable
.size(); ++j
)
461 if (gpuIdsAvailable
[i
] == gpuIdsAvailable
[j
])
463 GMX_THROW(InvalidInputError(formatString("The string of available GPU device IDs '%s' may not contain duplicate device IDs", hw_opt
.gpuIdsAvailable
.c_str())));
468 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
470 std::vector
<int> userGpuTaskAssignment
;
473 userGpuTaskAssignment
= parseUserGpuIds(hw_opt
.userGpuTaskAssignment
);
475 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
476 auto nonbondedTarget
= findTaskTarget(nbpu_opt
);
477 auto pmeTarget
= findTaskTarget(pme_opt
);
478 auto pmeFftTarget
= findTaskTarget(pme_fft_opt
);
479 auto bondedTarget
= findTaskTarget(bonded_opt
);
480 PmeRunMode pmeRunMode
= PmeRunMode::None
;
482 // Here we assume that SIMMASTER(cr) does not change even after the
483 // threads are started.
485 FILE *fplog
= nullptr;
486 // If we are appending, we don't write log output because we need
487 // to check that the old log file matches what the checkpoint file
488 // expects. Otherwise, we should start to write log output now if
489 // there is a file ready for it.
490 if (logFileHandle
!= nullptr && !mdrunOptions
.continuationOptions
.appendFiles
)
492 fplog
= gmx_fio_getfp(logFileHandle
);
494 gmx::LoggerOwner
logOwner(buildLogger(fplog
, cr
));
495 gmx::MDLogger
mdlog(logOwner
.logger());
497 // TODO The thread-MPI master rank makes a working
498 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
499 // after the threads have been launched. This works because no use
500 // is made of that communicator until after the execution paths
501 // have rejoined. But it is likely that we can improve the way
502 // this is expressed, e.g. by expressly running detection only the
503 // master rank for thread-MPI, rather than relying on the mutex
504 // and reference count.
505 PhysicalNodeCommunicator
physicalNodeComm(MPI_COMM_WORLD
, gmx_physicalnode_id_hash());
506 hwinfo
= gmx_detect_hardware(mdlog
, physicalNodeComm
);
508 gmx_print_detected_hardware(fplog
, cr
, ms
, mdlog
, hwinfo
);
510 std::vector
<int> gpuIdsToUse
;
511 auto compatibleGpus
= getCompatibleGpus(hwinfo
->gpu_info
);
512 if (gpuIdsAvailable
.empty())
514 gpuIdsToUse
= compatibleGpus
;
518 for (const auto &availableGpuId
: gpuIdsAvailable
)
520 bool availableGpuIsCompatible
= false;
521 for (const auto &compatibleGpuId
: compatibleGpus
)
523 if (availableGpuId
== compatibleGpuId
)
525 availableGpuIsCompatible
= true;
529 if (!availableGpuIsCompatible
)
531 gmx_fatal(FARGS
, "You limited the set of compatible GPUs to a set that included ID #%d, but that ID is not for a compatible GPU. List only compatible GPUs.", availableGpuId
);
533 gpuIdsToUse
.push_back(availableGpuId
);
537 if (fplog
!= nullptr)
539 /* Print references after all software/hardware printing */
540 please_cite(fplog
, "Abraham2015");
541 please_cite(fplog
, "Pall2015");
542 please_cite(fplog
, "Pronk2013");
543 please_cite(fplog
, "Hess2008b");
544 please_cite(fplog
, "Spoel2005a");
545 please_cite(fplog
, "Lindahl2001a");
546 please_cite(fplog
, "Berendsen95a");
547 writeSourceDoi(fplog
);
550 std::unique_ptr
<t_state
> globalState
;
554 /* Only the master rank has the global state */
555 globalState
= compat::make_unique
<t_state
>();
557 /* Read (nearly) all data required for the simulation */
558 read_tpx_state(ftp2fn(efTPR
, filenames
.size(), filenames
.data()), inputrec
, globalState
.get(), &mtop
);
560 /* In rerun, set velocities to zero if present */
561 if (doRerun
&& ((globalState
->flags
& (1 << estV
)) != 0))
563 // rerun does not use velocities
564 GMX_LOG(mdlog
.info
).asParagraph().appendText(
565 "Rerun trajectory contains velocities. Rerun does only evaluate "
566 "potential energy and forces. The velocities will be ignored.");
567 for (int i
= 0; i
< globalState
->natoms
; i
++)
569 clear_rvec(globalState
->v
[i
]);
571 globalState
->flags
&= ~(1 << estV
);
574 if (inputrec
->cutoff_scheme
!= ecutsVERLET
)
576 if (nstlist_cmdline
> 0)
578 gmx_fatal(FARGS
, "Can not set nstlist with the group cut-off scheme");
581 if (!compatibleGpus
.empty())
583 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
584 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
585 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
590 /* Check and update the hardware options for internal consistency */
591 check_and_update_hw_opt_1(mdlog
, &hw_opt
, cr
, domdecOptions
.numPmeRanks
);
593 /* Early check for externally set process affinity. */
594 gmx_check_thread_affinity_set(mdlog
, cr
,
595 &hw_opt
, hwinfo
->nthreads_hw_avail
, FALSE
);
597 if (GMX_THREAD_MPI
&& SIMMASTER(cr
))
599 if (domdecOptions
.numPmeRanks
> 0 && hw_opt
.nthreads_tmpi
<= 0)
601 gmx_fatal(FARGS
, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
604 /* Since the master knows the cut-off scheme, update hw_opt for this.
605 * This is done later for normal MPI and also once more with tMPI
606 * for all tMPI ranks.
608 check_and_update_hw_opt_2(&hw_opt
, inputrec
->cutoff_scheme
);
610 bool useGpuForNonbonded
= false;
611 bool useGpuForPme
= false;
614 // If the user specified the number of ranks, then we must
615 // respect that, but in default mode, we need to allow for
616 // the number of GPUs to choose the number of ranks.
618 useGpuForNonbonded
= decideWhetherToUseGpusForNonbondedWithThreadMpi
619 (nonbondedTarget
, gpuIdsToUse
, userGpuTaskAssignment
, emulateGpuNonbonded
,
620 inputrec
->cutoff_scheme
== ecutsVERLET
,
621 gpuAccelerationOfNonbondedIsUseful(mdlog
, inputrec
, GMX_THREAD_MPI
),
622 hw_opt
.nthreads_tmpi
);
623 auto canUseGpuForPme
= pme_gpu_supports_build(*hwinfo
, nullptr) && pme_gpu_supports_input(*inputrec
, mtop
, nullptr);
624 useGpuForPme
= decideWhetherToUseGpusForPmeWithThreadMpi
625 (useGpuForNonbonded
, pmeTarget
, gpuIdsToUse
, userGpuTaskAssignment
,
626 canUseGpuForPme
, hw_opt
.nthreads_tmpi
, domdecOptions
.numPmeRanks
);
629 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
631 /* Determine how many thread-MPI ranks to start.
633 * TODO Over-writing the user-supplied value here does
634 * prevent any possible subsequent checks from working
636 hw_opt
.nthreads_tmpi
= get_nthreads_mpi(hwinfo
,
645 // Now start the threads for thread MPI.
646 cr
= spawnThreads(hw_opt
.nthreads_tmpi
);
647 /* The main thread continues here with a new cr. We don't deallocate
648 the old cr because other threads may still be reading it. */
649 // TODO Both master and spawned threads call dup_tfn and
650 // reinitialize_commrec_for_this_thread. Find a way to express
652 physicalNodeComm
= PhysicalNodeCommunicator(MPI_COMM_WORLD
, gmx_physicalnode_id_hash());
654 // END OF CAUTION: cr and physicalNodeComm are now reliable
658 /* now broadcast everything to the non-master nodes/threads: */
659 init_parallel(cr
, inputrec
, &mtop
);
662 // Now each rank knows the inputrec that SIMMASTER read and used,
663 // and (if applicable) cr->nnodes has been assigned the number of
664 // thread-MPI ranks that have been chosen. The ranks can now all
665 // run the task-deciding functions and will agree on the result
666 // without needing to communicate.
668 // TODO Should we do the communication in debug mode to support
669 // having an assertion?
671 // Note that these variables describe only their own node.
673 // Note that when bonded interactions run on a GPU they always run
674 // alongside a nonbonded task, so do not influence task assignment
675 // even though they affect the force calculation workload.
676 bool useGpuForNonbonded
= false;
677 bool useGpuForPme
= false;
678 bool useGpuForBonded
= false;
681 // It's possible that there are different numbers of GPUs on
682 // different nodes, which is the user's responsibilty to
683 // handle. If unsuitable, we will notice that during task
685 bool gpusWereDetected
= hwinfo
->ngpu_compatible_tot
> 0;
686 bool usingVerletScheme
= inputrec
->cutoff_scheme
== ecutsVERLET
;
687 useGpuForNonbonded
= decideWhetherToUseGpusForNonbonded(nonbondedTarget
, userGpuTaskAssignment
,
688 emulateGpuNonbonded
, usingVerletScheme
,
689 gpuAccelerationOfNonbondedIsUseful(mdlog
, inputrec
, !GMX_THREAD_MPI
),
691 auto canUseGpuForPme
= pme_gpu_supports_build(*hwinfo
, nullptr) && pme_gpu_supports_input(*inputrec
, mtop
, nullptr);
692 useGpuForPme
= decideWhetherToUseGpusForPme(useGpuForNonbonded
, pmeTarget
, userGpuTaskAssignment
,
693 canUseGpuForPme
, cr
->nnodes
, domdecOptions
.numPmeRanks
,
695 auto canUseGpuForBonded
= buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec
, mtop
, nullptr);
697 decideWhetherToUseGpusForBonded(useGpuForNonbonded
, useGpuForPme
, usingVerletScheme
,
698 bondedTarget
, canUseGpuForBonded
,
699 EVDW_PME(inputrec
->vdwtype
),
700 EEL_PME_EWALD(inputrec
->coulombtype
),
701 domdecOptions
.numPmeRanks
, gpusWereDetected
);
703 pmeRunMode
= (useGpuForPme
? PmeRunMode::GPU
: PmeRunMode::CPU
);
704 if (pmeRunMode
== PmeRunMode::GPU
)
706 if (pmeFftTarget
== TaskTarget::Cpu
)
708 pmeRunMode
= PmeRunMode::Mixed
;
711 else if (pmeFftTarget
== TaskTarget::Gpu
)
713 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.");
716 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
719 // TODO: hide restraint implementation details from Mdrunner.
720 // There is nothing unique about restraints at this point as far as the
721 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
722 // factory functions from the SimulationContext on which to call mdModules->add().
723 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
724 for (auto && restraint
: restraintManager_
->getRestraints())
726 auto module
= RestraintMDModule::create(restraint
,
728 mdModules
->add(std::move(module
));
731 // TODO: Error handling
732 mdModules
->assignOptionsToModules(*inputrec
->params
, nullptr);
734 if (fplog
!= nullptr)
736 pr_inputrec(fplog
, 0, "Input Parameters", inputrec
, FALSE
);
737 fprintf(fplog
, "\n");
742 /* now make sure the state is initialized and propagated */
743 set_state_entries(globalState
.get(), inputrec
);
746 /* NM and TPI parallelize over force/energy calculations, not atoms,
747 * so we need to initialize and broadcast the global state.
749 if (inputrec
->eI
== eiNM
|| inputrec
->eI
== eiTPI
)
753 globalState
= compat::make_unique
<t_state
>();
755 broadcastStateWithoutDynamics(cr
, globalState
.get());
758 /* A parallel command line option consistency check that we can
759 only do after any threads have started. */
760 if (!PAR(cr
) && (domdecOptions
.numCells
[XX
] > 1 ||
761 domdecOptions
.numCells
[YY
] > 1 ||
762 domdecOptions
.numCells
[ZZ
] > 1 ||
763 domdecOptions
.numPmeRanks
> 0))
766 "The -dd or -npme option request a parallel simulation, "
768 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv
));
771 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
773 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv
));
779 (EI_ENERGY_MINIMIZATION(inputrec
->eI
) || eiNM
== inputrec
->eI
))
781 gmx_fatal(FARGS
, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
784 if (can_use_allvsall(inputrec
, TRUE
, cr
, fplog
) && DOMAINDECOMP(cr
))
786 gmx_fatal(FARGS
, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
789 if (!(EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
)))
791 if (domdecOptions
.numPmeRanks
> 0)
793 gmx_fatal_collective(FARGS
, cr
->mpi_comm_mysim
, MASTER(cr
),
794 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
797 domdecOptions
.numPmeRanks
= 0;
800 if (useGpuForNonbonded
&& domdecOptions
.numPmeRanks
< 0)
802 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
803 * improve performance with many threads per GPU, since our OpenMP
804 * scaling is bad, but it's difficult to automate the setup.
806 domdecOptions
.numPmeRanks
= 0;
810 if (domdecOptions
.numPmeRanks
< 0)
812 domdecOptions
.numPmeRanks
= 0;
813 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
817 GMX_RELEASE_ASSERT(domdecOptions
.numPmeRanks
<= 1, "PME GPU decomposition is not supported");
824 fcRegisterSteps(inputrec
->nsteps
, inputrec
->init_step
);
828 /* NMR restraints must be initialized before load_checkpoint,
829 * since with time averaging the history is added to t_state.
830 * For proper consistency check we therefore need to extend
832 * So the PME-only nodes (if present) will also initialize
833 * the distance restraints.
837 /* This needs to be called before read_checkpoint to extend the state */
838 init_disres(fplog
, &mtop
, inputrec
, cr
, ms
, fcd
, globalState
.get(), replExParams
.exchangeInterval
> 0);
840 init_orires(fplog
, &mtop
, inputrec
, cr
, ms
, globalState
.get(), &(fcd
->orires
));
842 auto deform
= prepareBoxDeformation(globalState
->box
, cr
, *inputrec
);
844 ObservablesHistory observablesHistory
= {};
846 ContinuationOptions
&continuationOptions
= mdrunOptions
.continuationOptions
;
848 if (continuationOptions
.startedFromCheckpoint
)
850 /* Check if checkpoint file exists before doing continuation.
851 * This way we can use identical input options for the first and subsequent runs...
855 load_checkpoint(opt2fn_master("-cpi", filenames
.size(), filenames
.data(), cr
),
857 cr
, domdecOptions
.numCells
,
858 inputrec
, globalState
.get(),
859 &bReadEkin
, &observablesHistory
,
860 continuationOptions
.appendFiles
,
861 continuationOptions
.appendFilesOptionSet
,
862 mdrunOptions
.reproducible
);
866 continuationOptions
.haveReadEkin
= true;
869 if (continuationOptions
.appendFiles
&& logFileHandle
)
871 // Now we can start normal logging to the truncated log file.
872 fplog
= gmx_fio_getfp(logFileHandle
);
873 prepareLogAppending(cr
->nodeid
, cr
->nnodes
, fplog
);
874 logOwner
= buildLogger(fplog
, cr
);
875 mdlog
= logOwner
.logger();
879 if (mdrunOptions
.numStepsCommandline
> -2)
881 GMX_LOG(mdlog
.info
).asParagraph().
882 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
883 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
885 /* override nsteps with value set on the commamdline */
886 override_nsteps_cmdline(mdlog
, mdrunOptions
.numStepsCommandline
, inputrec
);
890 copy_mat(globalState
->box
, box
);
895 gmx_bcast(sizeof(box
), box
, cr
);
898 /* Update rlist and nstlist. */
899 if (inputrec
->cutoff_scheme
== ecutsVERLET
)
901 prepare_verlet_scheme(fplog
, cr
, inputrec
, nstlist_cmdline
, &mtop
, box
,
902 useGpuForNonbonded
|| (emulateGpuNonbonded
== EmulateGpuNonbonded::Yes
), *hwinfo
->cpuInfo
);
905 LocalAtomSetManager atomSets
;
907 if (PAR(cr
) && !(EI_TPI(inputrec
->eI
) ||
908 inputrec
->eI
== eiNM
))
910 cr
->dd
= init_domain_decomposition(mdlog
, cr
, domdecOptions
, mdrunOptions
,
912 box
, positionsFromStatePointer(globalState
.get()),
914 // Note that local state still does not exist yet.
918 /* PME, if used, is done on all nodes with 1D decomposition */
920 cr
->duty
= (DUTY_PP
| DUTY_PME
);
922 if (inputrec
->ePBC
== epbcSCREW
)
925 "pbc=screw is only implemented with domain decomposition");
931 /* After possible communicator splitting in make_dd_communicators.
932 * we can set up the intra/inter node communication.
934 gmx_setup_nodecomm(fplog
, cr
);
940 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
941 "This is simulation %d out of %d running as a composite GROMACS\n"
942 "multi-simulation job. Setup for this simulation:\n",
945 GMX_LOG(mdlog
.warning
).appendTextFormatted(
949 cr
->nnodes
== 1 ? "thread" : "threads"
951 cr
->nnodes
== 1 ? "process" : "processes"
957 /* Check and update hw_opt for the cut-off scheme */
958 check_and_update_hw_opt_2(&hw_opt
, inputrec
->cutoff_scheme
);
960 /* Check and update the number of OpenMP threads requested */
961 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt
, *hwinfo
, cr
, ms
, physicalNodeComm
.size_
,
964 gmx_omp_nthreads_init(mdlog
, cr
,
965 hwinfo
->nthreads_hw_avail
,
966 physicalNodeComm
.size_
,
968 hw_opt
.nthreads_omp_pme
,
969 !thisRankHasDuty(cr
, DUTY_PP
),
970 inputrec
->cutoff_scheme
== ecutsVERLET
);
972 // Enable FP exception detection for the Verlet scheme, but not in
973 // Release mode and not for compilers with known buggy FP
974 // exception support (clang with any optimization) or suspected
975 // buggy FP exception support (gcc 7.* with optimization).
976 #if !defined NDEBUG && \
977 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
978 && defined __OPTIMIZE__)
979 const bool bEnableFPE
= inputrec
->cutoff_scheme
== ecutsVERLET
;
981 const bool bEnableFPE
= false;
983 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
986 gmx_feenableexcept();
989 // Build a data structure that expresses which kinds of non-bonded
990 // task are handled by this rank.
992 // TODO Later, this might become a loop over all registered modules
993 // relevant to the mdp inputs, to find those that have such tasks.
995 // TODO This could move before init_domain_decomposition() as part
996 // of refactoring that separates the responsibility for duty
997 // assignment from setup for communication between tasks, and
998 // setup for tasks handled with a domain (ie including short-ranged
999 // tasks, bonded tasks, etc.).
1001 // Note that in general useGpuForNonbonded, etc. can have a value
1002 // that is inconsistent with the presence of actual GPUs on any
1003 // rank, and that is not known to be a problem until the
1004 // duty of the ranks on a node become known.
1006 // TODO Later we might need the concept of computeTasksOnThisRank,
1007 // from which we construct gpuTasksOnThisRank.
1009 // Currently the DD code assigns duty to ranks that can
1010 // include PP work that currently can be executed on a single
1011 // GPU, if present and compatible. This has to be coordinated
1012 // across PP ranks on a node, with possible multiple devices
1013 // or sharing devices on a node, either from the user
1014 // selection, or automatically.
1015 auto haveGpus
= !gpuIdsToUse
.empty();
1016 std::vector
<GpuTask
> gpuTasksOnThisRank
;
1017 if (thisRankHasDuty(cr
, DUTY_PP
))
1019 if (useGpuForNonbonded
)
1021 // Note that any bonded tasks on a GPU always accompany a
1025 gpuTasksOnThisRank
.push_back(GpuTask::Nonbonded
);
1027 else if (nonbondedTarget
== TaskTarget::Gpu
)
1029 gmx_fatal(FARGS
, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1031 else if (bondedTarget
== TaskTarget::Gpu
)
1033 gmx_fatal(FARGS
, "Cannot run bonded interactions on a GPU because there is none detected.");
1037 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1038 if (EEL_PME(inputrec
->coulombtype
) && (thisRankHasDuty(cr
, DUTY_PME
)))
1044 gpuTasksOnThisRank
.push_back(GpuTask::Pme
);
1046 else if (pmeTarget
== TaskTarget::Gpu
)
1048 gmx_fatal(FARGS
, "Cannot run PME on a GPU because there is none detected.");
1053 GpuTaskAssignment gpuTaskAssignment
;
1056 // Produce the task assignment for this rank.
1057 gpuTaskAssignment
= runTaskAssignment(gpuIdsToUse
, userGpuTaskAssignment
, *hwinfo
,
1058 mdlog
, cr
, ms
, physicalNodeComm
, gpuTasksOnThisRank
,
1059 useGpuForBonded
, pmeRunMode
);
1061 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1063 /* Prevent other ranks from continuing after an issue was found
1064 * and reported as a fatal error.
1066 * TODO This function implements a barrier so that MPI runtimes
1067 * can organize an orderly shutdown if one of the ranks has had to
1068 * issue a fatal error in various code already run. When we have
1069 * MPI-aware error handling and reporting, this should be
1074 MPI_Barrier(cr
->mpi_comm_mysim
);
1080 MPI_Barrier(ms
->mpi_comm_masters
);
1082 /* We need another barrier to prevent non-master ranks from contiuing
1083 * when an error occured in a different simulation.
1085 MPI_Barrier(cr
->mpi_comm_mysim
);
1089 /* Now that we know the setup is consistent, check for efficiency */
1090 check_resource_division_efficiency(hwinfo
, !gpuTaskAssignment
.empty(), mdrunOptions
.ntompOptionIsSet
,
1093 gmx_device_info_t
*nonbondedDeviceInfo
= nullptr;
1095 if (thisRankHasDuty(cr
, DUTY_PP
))
1097 // This works because only one task of each type is currently permitted.
1098 auto nbGpuTaskMapping
= std::find_if(gpuTaskAssignment
.begin(), gpuTaskAssignment
.end(),
1099 hasTaskType
<GpuTask::Nonbonded
>);
1100 if (nbGpuTaskMapping
!= gpuTaskAssignment
.end())
1102 int nonbondedDeviceId
= nbGpuTaskMapping
->deviceId_
;
1103 nonbondedDeviceInfo
= getDeviceInfo(hwinfo
->gpu_info
, nonbondedDeviceId
);
1104 init_gpu(nonbondedDeviceInfo
);
1106 if (DOMAINDECOMP(cr
))
1108 /* When we share GPUs over ranks, we need to know this for the DLB */
1109 dd_setup_dlb_resource_sharing(cr
, nonbondedDeviceId
);
1115 std::unique_ptr
<ClfftInitializer
> initializedClfftLibrary
;
1117 gmx_device_info_t
*pmeDeviceInfo
= nullptr;
1118 // Later, this program could contain kernels that might be later
1119 // re-used as auto-tuning progresses, or subsequent simulations
1121 PmeGpuProgramStorage pmeGpuProgram
;
1122 // This works because only one task of each type is currently permitted.
1123 auto pmeGpuTaskMapping
= std::find_if(gpuTaskAssignment
.begin(), gpuTaskAssignment
.end(), hasTaskType
<GpuTask::Pme
>);
1124 const bool thisRankHasPmeGpuTask
= (pmeGpuTaskMapping
!= gpuTaskAssignment
.end());
1125 if (thisRankHasPmeGpuTask
)
1127 pmeDeviceInfo
= getDeviceInfo(hwinfo
->gpu_info
, pmeGpuTaskMapping
->deviceId_
);
1128 init_gpu(pmeDeviceInfo
);
1129 pmeGpuProgram
= buildPmeGpuProgram(pmeDeviceInfo
);
1130 // TODO It would be nice to move this logic into the factory
1131 // function. See Redmine #2535.
1132 bool isMasterThread
= !GMX_THREAD_MPI
|| MASTER(cr
);
1133 if (pmeRunMode
== PmeRunMode::GPU
&& !initializedClfftLibrary
&& isMasterThread
)
1135 initializedClfftLibrary
= initializeClfftLibrary();
1139 /* getting number of PP/PME threads
1140 PME: env variable should be read only on one node to make sure it is
1141 identical everywhere;
1143 nthreads_pme
= gmx_omp_nthreads_get(emntPME
);
1145 int numThreadsOnThisRank
;
1146 /* threads on this MPI process or TMPI thread */
1147 if (thisRankHasDuty(cr
, DUTY_PP
))
1149 numThreadsOnThisRank
= gmx_omp_nthreads_get(emntNonbonded
);
1153 numThreadsOnThisRank
= nthreads_pme
;
1156 checkHardwareOversubscription(numThreadsOnThisRank
, cr
->nodeid
,
1157 *hwinfo
->hardwareTopology
,
1158 physicalNodeComm
, mdlog
);
1160 if (hw_opt
.thread_affinity
!= threadaffOFF
)
1162 /* Before setting affinity, check whether the affinity has changed
1163 * - which indicates that probably the OpenMP library has changed it
1164 * since we first checked).
1166 gmx_check_thread_affinity_set(mdlog
, cr
,
1167 &hw_opt
, hwinfo
->nthreads_hw_avail
, TRUE
);
1169 int numThreadsOnThisNode
, intraNodeThreadOffset
;
1170 analyzeThreadsOnThisNode(physicalNodeComm
, numThreadsOnThisRank
, &numThreadsOnThisNode
,
1171 &intraNodeThreadOffset
);
1173 /* Set the CPU affinity */
1174 gmx_set_thread_affinity(mdlog
, cr
, &hw_opt
, *hwinfo
->hardwareTopology
,
1175 numThreadsOnThisRank
, numThreadsOnThisNode
,
1176 intraNodeThreadOffset
, nullptr);
1179 if (mdrunOptions
.timingOptions
.resetStep
> -1)
1181 GMX_LOG(mdlog
.info
).asParagraph().
1182 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1184 wcycle
= wallcycle_init(fplog
, mdrunOptions
.timingOptions
.resetStep
, cr
);
1188 /* Master synchronizes its value of reset_counters with all nodes
1189 * including PME only nodes */
1190 reset_counters
= wcycle_get_reset_counters(wcycle
);
1191 gmx_bcast_sim(sizeof(reset_counters
), &reset_counters
, cr
);
1192 wcycle_set_reset_counters(wcycle
, reset_counters
);
1195 // Membrane embedding must be initialized before we call init_forcerec()
1200 fprintf(stderr
, "Initializing membed");
1202 /* Note that membed cannot work in parallel because mtop is
1203 * changed here. Fix this if we ever want to make it run with
1204 * multiple ranks. */
1205 membed
= init_membed(fplog
, filenames
.size(), filenames
.data(), &mtop
, inputrec
, globalState
.get(), cr
,
1207 .checkpointOptions
.period
);
1210 std::unique_ptr
<MDAtoms
> mdAtoms
;
1211 std::unique_ptr
<gmx_vsite_t
> vsite
;
1214 if (thisRankHasDuty(cr
, DUTY_PP
))
1216 /* Initiate forcerecord */
1218 fr
->forceProviders
= mdModules
->initForceProviders();
1219 init_forcerec(fplog
, mdlog
, fr
, fcd
,
1220 inputrec
, &mtop
, cr
, box
,
1221 opt2fn("-table", filenames
.size(), filenames
.data()),
1222 opt2fn("-tablep", filenames
.size(), filenames
.data()),
1223 opt2fns("-tableb", filenames
.size(), filenames
.data()),
1224 *hwinfo
, nonbondedDeviceInfo
,
1229 /* Initialize the mdAtoms structure.
1230 * mdAtoms is not filled with atom data,
1231 * as this can not be done now with domain decomposition.
1233 mdAtoms
= makeMDAtoms(fplog
, mtop
, *inputrec
, thisRankHasPmeGpuTask
);
1234 if (globalState
&& thisRankHasPmeGpuTask
)
1236 // The pinning of coordinates in the global state object works, because we only use
1237 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1238 // points to the global state object without DD.
1239 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1240 // which should also perform the pinning.
1241 changePinningPolicy(&globalState
->x
, pme_get_pinning_policy());
1244 /* Initialize the virtual site communication */
1245 vsite
= initVsite(mtop
, cr
);
1247 calc_shifts(box
, fr
->shift_vec
);
1249 /* With periodic molecules the charge groups should be whole at start up
1250 * and the virtual sites should not be far from their proper positions.
1252 if (!inputrec
->bContinuation
&& MASTER(cr
) &&
1253 !(inputrec
->ePBC
!= epbcNONE
&& inputrec
->bPeriodicMols
))
1255 /* Make molecules whole at start of run */
1256 if (fr
->ePBC
!= epbcNONE
)
1258 do_pbc_first_mtop(fplog
, inputrec
->ePBC
, box
, &mtop
, globalState
->x
.rvec_array());
1262 /* Correct initial vsite positions are required
1263 * for the initial distribution in the domain decomposition
1264 * and for the initial shell prediction.
1266 constructVsitesGlobal(mtop
, globalState
->x
);
1270 if (EEL_PME(fr
->ic
->eeltype
) || EVDW_PME(fr
->ic
->vdwtype
))
1272 ewaldcoeff_q
= fr
->ic
->ewaldcoeff_q
;
1273 ewaldcoeff_lj
= fr
->ic
->ewaldcoeff_lj
;
1278 /* This is a PME only node */
1280 GMX_ASSERT(globalState
== nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1282 ewaldcoeff_q
= calc_ewaldcoeff_q(inputrec
->rcoulomb
, inputrec
->ewald_rtol
);
1283 ewaldcoeff_lj
= calc_ewaldcoeff_lj(inputrec
->rvdw
, inputrec
->ewald_rtol_lj
);
1286 gmx_pme_t
*sepPmeData
= nullptr;
1287 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1288 GMX_ASSERT(thisRankHasDuty(cr
, DUTY_PP
) == (fr
!= nullptr), "Double-checking that only PME-only ranks have no forcerec");
1289 gmx_pme_t
* &pmedata
= fr
? fr
->pmedata
: sepPmeData
;
1291 /* Initiate PME if necessary,
1292 * either on all nodes or on dedicated PME nodes only. */
1293 if (EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
))
1295 if (mdAtoms
&& mdAtoms
->mdatoms())
1297 nChargePerturbed
= mdAtoms
->mdatoms()->nChargePerturbed
;
1298 if (EVDW_PME(inputrec
->vdwtype
))
1300 nTypePerturbed
= mdAtoms
->mdatoms()->nTypePerturbed
;
1303 if (cr
->npmenodes
> 0)
1305 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1306 gmx_bcast_sim(sizeof(nChargePerturbed
), &nChargePerturbed
, cr
);
1307 gmx_bcast_sim(sizeof(nTypePerturbed
), &nTypePerturbed
, cr
);
1310 if (thisRankHasDuty(cr
, DUTY_PME
))
1314 pmedata
= gmx_pme_init(cr
,
1315 getNumPmeDomains(cr
->dd
),
1317 mtop
.natoms
, nChargePerturbed
!= 0, nTypePerturbed
!= 0,
1318 mdrunOptions
.reproducible
,
1319 ewaldcoeff_q
, ewaldcoeff_lj
,
1321 pmeRunMode
, nullptr,
1322 pmeDeviceInfo
, pmeGpuProgram
.get(), mdlog
);
1324 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1329 if (EI_DYNAMICS(inputrec
->eI
))
1331 /* Turn on signal handling on all nodes */
1333 * (A user signal from the PME nodes (if any)
1334 * is communicated to the PP nodes.
1336 signal_handler_install();
1339 if (thisRankHasDuty(cr
, DUTY_PP
))
1341 /* Assumes uniform use of the number of OpenMP threads */
1342 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntDefault
));
1344 if (inputrec
->bPull
)
1346 /* Initialize pull code */
1347 inputrec
->pull_work
=
1348 init_pull(fplog
, inputrec
->pull
, inputrec
,
1349 &mtop
, cr
, &atomSets
, inputrec
->fepvals
->init_lambda
);
1350 if (inputrec
->pull
->bXOutAverage
|| inputrec
->pull
->bFOutAverage
)
1352 initPullHistory(inputrec
->pull_work
, &observablesHistory
);
1354 if (EI_DYNAMICS(inputrec
->eI
) && MASTER(cr
))
1356 init_pull_output_files(inputrec
->pull_work
,
1357 filenames
.size(), filenames
.data(), oenv
,
1358 continuationOptions
);
1362 std::unique_ptr
<EnforcedRotation
> enforcedRotation
;
1365 /* Initialize enforced rotation code */
1366 enforcedRotation
= init_rot(fplog
,
1378 if (inputrec
->eSwapCoords
!= eswapNO
)
1380 /* Initialize ion swapping code */
1381 init_swapcoords(fplog
, inputrec
, opt2fn_master("-swap", filenames
.size(), filenames
.data(), cr
),
1382 &mtop
, globalState
.get(), &observablesHistory
,
1383 cr
, &atomSets
, oenv
, mdrunOptions
);
1386 /* Let makeConstraints know whether we have essential dynamics constraints.
1387 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1389 bool doEssentialDynamics
= (opt2fn_null("-ei", filenames
.size(), filenames
.data()) != nullptr
1390 || observablesHistory
.edsamHistory
);
1391 auto constr
= makeConstraints(mtop
, *inputrec
, doEssentialDynamics
,
1392 fplog
, *mdAtoms
->mdatoms(),
1393 cr
, *ms
, nrnb
, wcycle
, fr
->bMolPBC
);
1395 if (DOMAINDECOMP(cr
))
1397 GMX_RELEASE_ASSERT(fr
, "fr was NULL while cr->duty was DUTY_PP");
1398 /* This call is not included in init_domain_decomposition mainly
1399 * because fr->cginfo_mb is set later.
1401 dd_init_bondeds(fplog
, cr
->dd
, &mtop
, vsite
.get(), inputrec
,
1402 domdecOptions
.checkBondedInteractions
,
1406 // TODO This is not the right place to manage the lifetime of
1407 // this data structure, but currently it's the easiest way to
1408 // make it work. Later, it should probably be made/updated
1409 // after the workload for the lifetime of a PP domain is
1411 PpForceWorkload ppForceWorkload
;
1413 GMX_ASSERT(stopHandlerBuilder_
, "Runner must provide StopHandlerBuilder to integrator.");
1414 /* Now do whatever the user wants us to do (how flexible...) */
1415 Integrator integrator
{
1416 fplog
, cr
, ms
, mdlog
, static_cast<int>(filenames
.size()), filenames
.data(),
1419 vsite
.get(), constr
.get(),
1420 enforcedRotation
? enforcedRotation
->getLegacyEnfrot() : nullptr,
1422 mdModules
->outputProvider(),
1426 &observablesHistory
,
1427 mdAtoms
.get(), nrnb
, wcycle
, fr
,
1431 walltime_accounting
,
1432 std::move(stopHandlerBuilder_
)
1434 integrator
.run(inputrec
->eI
, doRerun
);
1436 if (inputrec
->bPull
)
1438 finish_pull(inputrec
->pull_work
);
1444 GMX_RELEASE_ASSERT(pmedata
, "pmedata was NULL while cr->duty was not DUTY_PP");
1446 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntPME
));
1447 gmx_pmeonly(pmedata
, cr
, nrnb
, wcycle
, walltime_accounting
, inputrec
, pmeRunMode
);
1450 wallcycle_stop(wcycle
, ewcRUN
);
1452 /* Finish up, write some stuff
1453 * if rerunMD, don't write last frame again
1455 finish_run(fplog
, mdlog
, cr
,
1456 inputrec
, nrnb
, wcycle
, walltime_accounting
,
1457 fr
? fr
->nbv
: nullptr,
1459 EI_DYNAMICS(inputrec
->eI
) && !isMultiSim(ms
));
1464 gmx_pme_destroy(pmedata
);
1468 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1469 // before we destroy the GPU context(s) in free_gpu_resources().
1470 // Pinned buffers are associated with contexts in CUDA.
1471 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1472 mdAtoms
.reset(nullptr);
1473 globalState
.reset(nullptr);
1474 mdModules
.reset(nullptr); // destruct force providers here as they might also use the GPU
1476 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1477 free_gpu_resources(fr
, physicalNodeComm
);
1478 free_gpu(nonbondedDeviceInfo
);
1479 free_gpu(pmeDeviceInfo
);
1480 done_forcerec(fr
, mtop
.molblock
.size(), mtop
.groups
.grps
[egcENER
].nr
);
1485 free_membed(membed
);
1488 gmx_hardware_info_free();
1490 /* Does what it says */
1491 print_date_and_time(fplog
, cr
->nodeid
, "Finished mdrun", gmx_gettime());
1492 walltime_accounting_destroy(walltime_accounting
);
1495 // Ensure log file content is written
1498 gmx_fio_flush(logFileHandle
);
1501 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1502 * exceptions were enabled before function was called. */
1505 gmx_fedisableexcept();
1508 rc
= static_cast<int>(gmx_get_stop_condition());
1511 /* we need to join all threads. The sub-threads join when they
1512 exit this function, but the master thread needs to be told to
1514 if (PAR(cr
) && MASTER(cr
))
1524 Mdrunner::~Mdrunner()
1526 // Clean up of the Manager.
1527 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1528 // but okay as long as threads synchronize some time before adding or accessing
1529 // a new set of restraints.
1530 if (restraintManager_
)
1532 restraintManager_
->clear();
1533 GMX_ASSERT(restraintManager_
->countRestraints() == 0,
1534 "restraints added during runner life time should be cleared at runner destruction.");
1538 void Mdrunner::addPotential(std::shared_ptr
<gmx::IRestraintPotential
> puller
,
1541 GMX_ASSERT(restraintManager_
, "Mdrunner must have a restraint manager.");
1542 // Not sure if this should be logged through the md logger or something else,
1543 // but it is helpful to have some sort of INFO level message sent somewhere.
1544 // std::cout << "Registering restraint named " << name << std::endl;
1546 // When multiple restraints are used, it may be wasteful to register them separately.
1547 // Maybe instead register an entire Restraint Manager as a force provider.
1548 restraintManager_
->addToSpec(std::move(puller
),
1552 Mdrunner::Mdrunner(Mdrunner
&&) noexcept
= default;
1554 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1555 Mdrunner
&Mdrunner::operator=(Mdrunner
&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING
) = default;
1557 class Mdrunner::BuilderImplementation
1560 BuilderImplementation() = delete;
1561 explicit BuilderImplementation(SimulationContext
* context
);
1562 ~BuilderImplementation();
1564 BuilderImplementation
&setExtraMdrunOptions(const MdrunOptions
&options
,
1565 real forceWarningThreshold
);
1567 void addDomdec(const DomdecOptions
&options
);
1569 void addVerletList(int nstlist
);
1571 void addReplicaExchange(const ReplicaExchangeParameters
¶ms
);
1573 void addMultiSim(gmx_multisim_t
* multisim
);
1575 void addNonBonded(const char* nbpu_opt
);
1577 void addPME(const char* pme_opt_
, const char* pme_fft_opt_
);
1579 void addBondedTaskAssignment(const char* bonded_opt
);
1581 void addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
);
1583 void addFilenames(ArrayRef
<const t_filenm
> filenames
);
1585 void addOutputEnvironment(gmx_output_env_t
* outputEnvironment
);
1587 void addLogFile(t_fileio
*logFileHandle
);
1589 void addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
);
1594 // Default parameters copied from runner.h
1595 // \todo Clarify source(s) of default parameters.
1597 const char* nbpu_opt_
= nullptr;
1598 const char* pme_opt_
= nullptr;
1599 const char* pme_fft_opt_
= nullptr;
1600 const char *bonded_opt_
= nullptr;
1602 MdrunOptions mdrunOptions_
;
1604 DomdecOptions domdecOptions_
;
1606 ReplicaExchangeParameters replicaExchangeParameters_
;
1608 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1611 //! Non-owning multisim communicator handle.
1612 std::unique_ptr
<gmx_multisim_t
*> multisim_
= nullptr;
1614 //! Print a warning if any force is larger than this (in kJ/mol nm).
1615 real forceWarningThreshold_
= -1;
1617 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1620 * \todo Establish robust protocol to make sure resources remain valid.
1621 * SimulationContext will likely be separated into multiple layers for
1622 * different levels of access and different phases of execution. Ref
1623 * https://redmine.gromacs.org/issues/2375
1624 * https://redmine.gromacs.org/issues/2587
1626 SimulationContext
* context_
= nullptr;
1628 //! \brief Parallelism information.
1629 gmx_hw_opt_t hardwareOptions_
;
1631 //! filename options for simulation.
1632 ArrayRef
<const t_filenm
> filenames_
;
1634 /*! \brief Handle to output environment.
1636 * \todo gmx_output_env_t needs lifetime management.
1638 gmx_output_env_t
* outputEnvironment_
= nullptr;
1640 /*! \brief Non-owning handle to MD log file.
1642 * \todo Context should own output facilities for client.
1643 * \todo Improve log file handle management.
1645 * Code managing the FILE* relies on the ability to set it to
1646 * nullptr to check whether the filehandle is valid.
1648 t_fileio
* logFileHandle_
= nullptr;
1651 * \brief Builder for simulation stop signal handler.
1653 std::unique_ptr
<StopHandlerBuilder
> stopHandlerBuilder_
= nullptr;
1656 Mdrunner::BuilderImplementation::BuilderImplementation(SimulationContext
* context
) :
1659 GMX_ASSERT(context_
, "Bug found. It should not be possible to construct builder without a valid context.");
1662 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1664 Mdrunner::BuilderImplementation
&
1665 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions
&options
,
1666 real forceWarningThreshold
)
1668 mdrunOptions_
= options
;
1669 forceWarningThreshold_
= forceWarningThreshold
;
1673 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions
&options
)
1675 domdecOptions_
= options
;
1678 void Mdrunner::BuilderImplementation::addVerletList(int nstlist
)
1683 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1685 replicaExchangeParameters_
= params
;
1688 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t
* multisim
)
1690 multisim_
= compat::make_unique
<gmx_multisim_t
*>(multisim
);
1693 Mdrunner
Mdrunner::BuilderImplementation::build()
1695 auto newRunner
= Mdrunner();
1697 GMX_ASSERT(context_
, "Bug found. It should not be possible to call build() without a valid context.");
1699 newRunner
.mdrunOptions
= mdrunOptions_
;
1700 newRunner
.domdecOptions
= domdecOptions_
;
1702 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1703 newRunner
.hw_opt
= hardwareOptions_
;
1705 // No invariant to check. This parameter exists to optionally override other behavior.
1706 newRunner
.nstlist_cmdline
= nstlist_
;
1708 newRunner
.replExParams
= replicaExchangeParameters_
;
1710 newRunner
.filenames
= filenames_
;
1712 GMX_ASSERT(context_
->communicationRecord_
, "SimulationContext communications not initialized.");
1713 newRunner
.cr
= context_
->communicationRecord_
;
1717 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1718 newRunner
.ms
= *multisim_
;
1722 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1725 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1726 // \todo Update sanity checking when output environment has clearly specified invariants.
1727 // Initialization and default values for oenv are not well specified in the current version.
1728 if (outputEnvironment_
)
1730 newRunner
.oenv
= outputEnvironment_
;
1734 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1737 newRunner
.logFileHandle
= logFileHandle_
;
1741 newRunner
.nbpu_opt
= nbpu_opt_
;
1745 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1748 if (pme_opt_
&& pme_fft_opt_
)
1750 newRunner
.pme_opt
= pme_opt_
;
1751 newRunner
.pme_fft_opt
= pme_fft_opt_
;
1755 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1760 newRunner
.bonded_opt
= bonded_opt_
;
1764 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1767 newRunner
.restraintManager_
= compat::make_unique
<gmx::RestraintManager
>();
1769 if (stopHandlerBuilder_
)
1771 newRunner
.stopHandlerBuilder_
= std::move(stopHandlerBuilder_
);
1775 newRunner
.stopHandlerBuilder_
= compat::make_unique
<StopHandlerBuilder
>();
1781 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt
)
1783 nbpu_opt_
= nbpu_opt
;
1786 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt
,
1787 const char* pme_fft_opt
)
1790 pme_fft_opt_
= pme_fft_opt
;
1793 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt
)
1795 bonded_opt_
= bonded_opt
;
1798 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
1800 hardwareOptions_
= hardwareOptions
;
1803 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef
<const t_filenm
> filenames
)
1805 filenames_
= filenames
;
1808 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
1810 outputEnvironment_
= outputEnvironment
;
1813 void Mdrunner::BuilderImplementation::addLogFile(t_fileio
*logFileHandle
)
1815 logFileHandle_
= logFileHandle
;
1818 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
1820 stopHandlerBuilder_
= std::move(builder
);
1823 MdrunnerBuilder::MdrunnerBuilder(compat::not_null
<SimulationContext
*> context
) :
1824 impl_
{gmx::compat::make_unique
<Mdrunner::BuilderImplementation
>(context
)}
1828 MdrunnerBuilder::~MdrunnerBuilder() = default;
1830 MdrunnerBuilder
&MdrunnerBuilder::addSimulationMethod(const MdrunOptions
&options
,
1831 real forceWarningThreshold
)
1833 impl_
->setExtraMdrunOptions(options
, forceWarningThreshold
);
1837 MdrunnerBuilder
&MdrunnerBuilder::addDomainDecomposition(const DomdecOptions
&options
)
1839 impl_
->addDomdec(options
);
1843 MdrunnerBuilder
&MdrunnerBuilder::addNeighborList(int nstlist
)
1845 impl_
->addVerletList(nstlist
);
1849 MdrunnerBuilder
&MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1851 impl_
->addReplicaExchange(params
);
1855 MdrunnerBuilder
&MdrunnerBuilder::addMultiSim(gmx_multisim_t
* multisim
)
1857 impl_
->addMultiSim(multisim
);
1861 MdrunnerBuilder
&MdrunnerBuilder::addNonBonded(const char* nbpu_opt
)
1863 impl_
->addNonBonded(nbpu_opt
);
1867 MdrunnerBuilder
&MdrunnerBuilder::addElectrostatics(const char* pme_opt
,
1868 const char* pme_fft_opt
)
1870 // The builder method may become more general in the future, but in this version,
1871 // parameters for PME electrostatics are both required and the only parameters
1873 if (pme_opt
&& pme_fft_opt
)
1875 impl_
->addPME(pme_opt
, pme_fft_opt
);
1879 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
1884 MdrunnerBuilder
&MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt
)
1886 impl_
->addBondedTaskAssignment(bonded_opt
);
1890 Mdrunner
MdrunnerBuilder::build()
1892 return impl_
->build();
1895 MdrunnerBuilder
&MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
1897 impl_
->addHardwareOptions(hardwareOptions
);
1901 MdrunnerBuilder
&MdrunnerBuilder::addFilenames(ArrayRef
<const t_filenm
> filenames
)
1903 impl_
->addFilenames(filenames
);
1907 MdrunnerBuilder
&MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
1909 impl_
->addOutputEnvironment(outputEnvironment
);
1913 MdrunnerBuilder
&MdrunnerBuilder::addLogFile(t_fileio
*logFileHandle
)
1915 impl_
->addLogFile(logFileHandle
);
1919 MdrunnerBuilder
&MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
1921 impl_
->addStopHandlerBuilder(std::move(builder
));
1925 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder
&&) noexcept
= default;
1927 MdrunnerBuilder
&MdrunnerBuilder::operator=(MdrunnerBuilder
&&) noexcept
= default;