2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements the MD runner routine calling all integrators.
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_mdrun
59 #include "gromacs/commandline/filenm.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/localatomsetmanager.h"
63 #include "gromacs/ewald/ewald_utils.h"
64 #include "gromacs/ewald/pme.h"
65 #include "gromacs/ewald/pme_gpu_program.h"
66 #include "gromacs/fileio/checkpoint.h"
67 #include "gromacs/fileio/gmxfio.h"
68 #include "gromacs/fileio/oenv.h"
69 #include "gromacs/fileio/tpxio.h"
70 #include "gromacs/gmxlib/network.h"
71 #include "gromacs/gmxlib/nrnb.h"
72 #include "gromacs/gpu_utils/clfftinitializer.h"
73 #include "gromacs/gpu_utils/gpu_utils.h"
74 #include "gromacs/hardware/cpuinfo.h"
75 #include "gromacs/hardware/detecthardware.h"
76 #include "gromacs/hardware/printhardware.h"
77 #include "gromacs/listed_forces/disre.h"
78 #include "gromacs/listed_forces/gpubonded.h"
79 #include "gromacs/listed_forces/orires.h"
80 #include "gromacs/math/functions.h"
81 #include "gromacs/math/utilities.h"
82 #include "gromacs/math/vec.h"
83 #include "gromacs/mdlib/boxdeformation.h"
84 #include "gromacs/mdlib/calc_verletbuf.h"
85 #include "gromacs/mdlib/forcerec.h"
86 #include "gromacs/mdlib/gmx_omp_nthreads.h"
87 #include "gromacs/mdlib/makeconstraints.h"
88 #include "gromacs/mdlib/md_support.h"
89 #include "gromacs/mdlib/mdatoms.h"
90 #include "gromacs/mdlib/mdrun.h"
91 #include "gromacs/mdlib/membed.h"
92 #include "gromacs/mdlib/ppforceworkload.h"
93 #include "gromacs/mdlib/qmmm.h"
94 #include "gromacs/mdlib/sighandler.h"
95 #include "gromacs/mdlib/sim_util.h"
96 #include "gromacs/mdlib/stophandler.h"
97 #include "gromacs/mdrun/legacymdrunoptions.h"
98 #include "gromacs/mdrun/logging.h"
99 #include "gromacs/mdrun/multisim.h"
100 #include "gromacs/mdrun/simulationcontext.h"
101 #include "gromacs/mdrunutility/mdmodules.h"
102 #include "gromacs/mdrunutility/threadaffinity.h"
103 #include "gromacs/mdtypes/commrec.h"
104 #include "gromacs/mdtypes/fcdata.h"
105 #include "gromacs/mdtypes/inputrec.h"
106 #include "gromacs/mdtypes/md_enums.h"
107 #include "gromacs/mdtypes/observableshistory.h"
108 #include "gromacs/mdtypes/state.h"
109 #include "gromacs/nbnxm/nbnxm.h"
110 #include "gromacs/nbnxm/pairlist_tuning.h"
111 #include "gromacs/pbcutil/pbc.h"
112 #include "gromacs/pulling/output.h"
113 #include "gromacs/pulling/pull.h"
114 #include "gromacs/pulling/pull_rotation.h"
115 #include "gromacs/restraint/manager.h"
116 #include "gromacs/restraint/restraintmdmodule.h"
117 #include "gromacs/restraint/restraintpotential.h"
118 #include "gromacs/swap/swapcoords.h"
119 #include "gromacs/taskassignment/decidegpuusage.h"
120 #include "gromacs/taskassignment/resourcedivision.h"
121 #include "gromacs/taskassignment/taskassignment.h"
122 #include "gromacs/taskassignment/usergpuids.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/topology/mtop_util.h"
125 #include "gromacs/trajectory/trajectoryframe.h"
126 #include "gromacs/utility/basenetwork.h"
127 #include "gromacs/utility/cstringutil.h"
128 #include "gromacs/utility/exceptions.h"
129 #include "gromacs/utility/fatalerror.h"
130 #include "gromacs/utility/filestream.h"
131 #include "gromacs/utility/gmxassert.h"
132 #include "gromacs/utility/gmxmpi.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/loggerbuilder.h"
135 #include "gromacs/utility/physicalnodecommunicator.h"
136 #include "gromacs/utility/pleasecite.h"
137 #include "gromacs/utility/programcontext.h"
138 #include "gromacs/utility/smalloc.h"
139 #include "gromacs/utility/stringutil.h"
141 #include "integrator.h"
142 #include "replicaexchange.h"
145 #include "corewrap.h"
151 /*! \brief Barrier for safe simultaneous thread access to mdrunner data
153 * Used to ensure that the master thread does not modify mdrunner during copy
154 * on the spawned threads. */
155 static void threadMpiMdrunnerAccessBarrier()
158 MPI_Barrier(MPI_COMM_WORLD
);
162 Mdrunner
Mdrunner::cloneOnSpawnedThread() const
164 auto newRunner
= Mdrunner();
166 // All runners in the same process share a restraint manager resource because it is
167 // part of the interface to the client code, which is associated only with the
168 // original thread. Handles to the same resources can be obtained by copy.
170 newRunner
.restraintManager_
= std::make_unique
<RestraintManager
>(*restraintManager_
);
173 // Copy original cr pointer before master thread can pass the thread barrier
174 newRunner
.cr
= reinitialize_commrec_for_this_thread(cr
);
176 // Copy members of master runner.
177 // \todo Replace with builder when Simulation context and/or runner phases are better defined.
178 // Ref https://redmine.gromacs.org/issues/2587 and https://redmine.gromacs.org/issues/2375
179 newRunner
.hw_opt
= hw_opt
;
180 newRunner
.filenames
= filenames
;
182 newRunner
.oenv
= oenv
;
183 newRunner
.mdrunOptions
= mdrunOptions
;
184 newRunner
.domdecOptions
= domdecOptions
;
185 newRunner
.nbpu_opt
= nbpu_opt
;
186 newRunner
.pme_opt
= pme_opt
;
187 newRunner
.pme_fft_opt
= pme_fft_opt
;
188 newRunner
.bonded_opt
= bonded_opt
;
189 newRunner
.nstlist_cmdline
= nstlist_cmdline
;
190 newRunner
.replExParams
= replExParams
;
191 newRunner
.pforce
= pforce
;
193 newRunner
.stopHandlerBuilder_
= std::make_unique
<StopHandlerBuilder
>(*stopHandlerBuilder_
);
195 threadMpiMdrunnerAccessBarrier();
197 GMX_RELEASE_ASSERT(!MASTER(newRunner
.cr
), "cloneOnSpawnedThread should only be called on spawned threads");
202 /*! \brief The callback used for running on spawned threads.
204 * Obtains the pointer to the master mdrunner object from the one
205 * argument permitted to the thread-launch API call, copies it to make
206 * a new runner for this thread, reinitializes necessary data, and
207 * proceeds to the simulation. */
208 static void mdrunner_start_fn(const void *arg
)
212 auto masterMdrunner
= reinterpret_cast<const gmx::Mdrunner
*>(arg
);
213 /* copy the arg list to make sure that it's thread-local. This
214 doesn't copy pointed-to items, of course; fnm, cr and fplog
215 are reset in the call below, all others should be const. */
216 gmx::Mdrunner mdrunner
= masterMdrunner
->cloneOnSpawnedThread();
219 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
223 /*! \brief Start thread-MPI threads.
225 * Called by mdrunner() to start a specific number of threads
226 * (including the main thread) for thread-parallel runs. This in turn
227 * calls mdrunner() for each thread. All options are the same as for
229 t_commrec
*Mdrunner::spawnThreads(int numThreadsToLaunch
) const
232 /* first check whether we even need to start tMPI */
233 if (numThreadsToLaunch
< 2)
239 /* now spawn new threads that start mdrunner_start_fn(), while
240 the main thread returns, we set thread affinity later */
241 if (tMPI_Init_fn(TRUE
, numThreadsToLaunch
, TMPI_AFFINITY_NONE
,
242 mdrunner_start_fn
, static_cast<const void*>(this)) != TMPI_SUCCESS
)
244 GMX_THROW(gmx::InternalError("Failed to spawn thread-MPI threads"));
247 threadMpiMdrunnerAccessBarrier();
249 GMX_UNUSED_VALUE(mdrunner_start_fn
);
252 return reinitialize_commrec_for_this_thread(cr
);
257 /*! \brief Initialize variables for Verlet scheme simulation */
258 static void prepare_verlet_scheme(FILE *fplog
,
262 const gmx_mtop_t
*mtop
,
264 bool makeGpuPairList
,
265 const gmx::CpuInfo
&cpuinfo
)
267 /* For NVE simulations, we will retain the initial list buffer */
268 if (EI_DYNAMICS(ir
->eI
) &&
269 ir
->verletbuf_tol
> 0 &&
270 !(EI_MD(ir
->eI
) && ir
->etc
== etcNO
))
272 /* Update the Verlet buffer size for the current run setup */
274 /* Here we assume SIMD-enabled kernels are being used. But as currently
275 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
276 * and 4x2 gives a larger buffer than 4x4, this is ok.
278 ListSetupType listType
= (makeGpuPairList
? ListSetupType::Gpu
: ListSetupType::CpuSimdWhenSupported
);
279 VerletbufListSetup listSetup
= verletbufGetSafeListSetup(listType
);
282 calc_verlet_buffer_size(mtop
, det(box
), ir
, ir
->nstlist
, ir
->nstlist
- 1, -1, &listSetup
, nullptr, &rlist_new
);
284 if (rlist_new
!= ir
->rlist
)
286 if (fplog
!= nullptr)
288 fprintf(fplog
, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
289 ir
->rlist
, rlist_new
,
290 listSetup
.cluster_size_i
, listSetup
.cluster_size_j
);
292 ir
->rlist
= rlist_new
;
296 if (nstlist_cmdline
> 0 && (!EI_DYNAMICS(ir
->eI
) || ir
->verletbuf_tol
<= 0))
298 gmx_fatal(FARGS
, "Can not set nstlist without %s",
299 !EI_DYNAMICS(ir
->eI
) ? "dynamics" : "verlet-buffer-tolerance");
302 if (EI_DYNAMICS(ir
->eI
))
304 /* Set or try nstlist values */
305 increaseNstlist(fplog
, cr
, ir
, nstlist_cmdline
, mtop
, box
, makeGpuPairList
, cpuinfo
);
309 /*! \brief Override the nslist value in inputrec
311 * with value passed on the command line (if any)
313 static void override_nsteps_cmdline(const gmx::MDLogger
&mdlog
,
314 int64_t nsteps_cmdline
,
319 /* override with anything else than the default -2 */
320 if (nsteps_cmdline
> -2)
322 char sbuf_steps
[STEPSTRSIZE
];
323 char sbuf_msg
[STRLEN
];
325 ir
->nsteps
= nsteps_cmdline
;
326 if (EI_DYNAMICS(ir
->eI
) && nsteps_cmdline
!= -1)
328 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
329 gmx_step_str(nsteps_cmdline
, sbuf_steps
),
330 fabs(nsteps_cmdline
*ir
->delta_t
));
334 sprintf(sbuf_msg
, "Overriding nsteps with value passed on the command line: %s steps",
335 gmx_step_str(nsteps_cmdline
, sbuf_steps
));
338 GMX_LOG(mdlog
.warning
).asParagraph().appendText(sbuf_msg
);
340 else if (nsteps_cmdline
< -2)
342 gmx_fatal(FARGS
, "Invalid nsteps value passed on the command line: %" PRId64
,
345 /* Do nothing if nsteps_cmdline == -2 */
351 /*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
353 * If not, and if a warning may be issued, logs a warning about
354 * falling back to CPU code. With thread-MPI, only the first
355 * call to this function should have \c issueWarning true. */
356 static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger
&mdlog
,
357 const t_inputrec
*ir
,
360 if (ir
->opts
.ngener
- ir
->nwall
> 1)
362 /* The GPU code does not support more than one energy group.
363 * If the user requested GPUs explicitly, a fatal error is given later.
367 GMX_LOG(mdlog
.warning
).asParagraph()
368 .appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
369 "For better performance, run on the GPU without energy groups and then do "
370 "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
377 //! Initializes the logger for mdrun.
378 static gmx::LoggerOwner
buildLogger(FILE *fplog
, const t_commrec
*cr
)
380 gmx::LoggerBuilder builder
;
381 if (fplog
!= nullptr)
383 builder
.addTargetFile(gmx::MDLogger::LogLevel::Info
, fplog
);
385 if (cr
== nullptr || SIMMASTER(cr
))
387 builder
.addTargetStream(gmx::MDLogger::LogLevel::Warning
,
388 &gmx::TextOutputFile::standardError());
390 return builder
.build();
393 //! Make a TaskTarget from an mdrun argument string.
394 static TaskTarget
findTaskTarget(const char *optionString
)
396 TaskTarget returnValue
= TaskTarget::Auto
;
398 if (strncmp(optionString
, "auto", 3) == 0)
400 returnValue
= TaskTarget::Auto
;
402 else if (strncmp(optionString
, "cpu", 3) == 0)
404 returnValue
= TaskTarget::Cpu
;
406 else if (strncmp(optionString
, "gpu", 3) == 0)
408 returnValue
= TaskTarget::Gpu
;
412 GMX_ASSERT(false, "Option string should have been checked for sanity already");
418 int Mdrunner::mdrunner()
422 t_forcerec
*fr
= nullptr;
423 t_fcdata
*fcd
= nullptr;
424 real ewaldcoeff_q
= 0;
425 real ewaldcoeff_lj
= 0;
426 int nChargePerturbed
= -1, nTypePerturbed
= 0;
427 gmx_wallcycle_t wcycle
;
428 gmx_walltime_accounting_t walltime_accounting
= nullptr;
429 gmx_membed_t
* membed
= nullptr;
430 gmx_hw_info_t
*hwinfo
= nullptr;
432 /* CAUTION: threads may be started later on in this function, so
433 cr doesn't reflect the final parallel state right now */
434 std::unique_ptr
<gmx::MDModules
> mdModules(new gmx::MDModules
);
435 t_inputrec inputrecInstance
;
436 t_inputrec
*inputrec
= &inputrecInstance
;
439 bool doMembed
= opt2bSet("-membed", filenames
.size(), filenames
.data());
440 bool doRerun
= mdrunOptions
.rerun
;
442 // Handle task-assignment related user options.
443 EmulateGpuNonbonded emulateGpuNonbonded
= (getenv("GMX_EMULATE_GPU") != nullptr ?
444 EmulateGpuNonbonded::Yes
: EmulateGpuNonbonded::No
);
445 std::vector
<int> gpuIdsAvailable
;
448 gpuIdsAvailable
= parseUserGpuIdString(hw_opt
.gpuIdsAvailable
);
450 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
452 std::vector
<int> userGpuTaskAssignment
;
455 userGpuTaskAssignment
= parseUserTaskAssignmentString(hw_opt
.userGpuTaskAssignment
);
457 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
458 auto nonbondedTarget
= findTaskTarget(nbpu_opt
);
459 auto pmeTarget
= findTaskTarget(pme_opt
);
460 auto pmeFftTarget
= findTaskTarget(pme_fft_opt
);
461 auto bondedTarget
= findTaskTarget(bonded_opt
);
462 PmeRunMode pmeRunMode
= PmeRunMode::None
;
464 // Here we assume that SIMMASTER(cr) does not change even after the
465 // threads are started.
467 FILE *fplog
= nullptr;
468 // If we are appending, we don't write log output because we need
469 // to check that the old log file matches what the checkpoint file
470 // expects. Otherwise, we should start to write log output now if
471 // there is a file ready for it.
472 if (logFileHandle
!= nullptr && !mdrunOptions
.continuationOptions
.appendFiles
)
474 fplog
= gmx_fio_getfp(logFileHandle
);
476 gmx::LoggerOwner
logOwner(buildLogger(fplog
, cr
));
477 gmx::MDLogger
mdlog(logOwner
.logger());
479 // TODO The thread-MPI master rank makes a working
480 // PhysicalNodeCommunicator here, but it gets rebuilt by all ranks
481 // after the threads have been launched. This works because no use
482 // is made of that communicator until after the execution paths
483 // have rejoined. But it is likely that we can improve the way
484 // this is expressed, e.g. by expressly running detection only the
485 // master rank for thread-MPI, rather than relying on the mutex
486 // and reference count.
487 PhysicalNodeCommunicator
physicalNodeComm(MPI_COMM_WORLD
, gmx_physicalnode_id_hash());
488 hwinfo
= gmx_detect_hardware(mdlog
, physicalNodeComm
);
490 gmx_print_detected_hardware(fplog
, cr
, ms
, mdlog
, hwinfo
);
492 std::vector
<int> gpuIdsToUse
;
493 auto compatibleGpus
= getCompatibleGpus(hwinfo
->gpu_info
);
494 if (gpuIdsAvailable
.empty())
496 gpuIdsToUse
= compatibleGpus
;
500 for (const auto &availableGpuId
: gpuIdsAvailable
)
502 bool availableGpuIsCompatible
= false;
503 for (const auto &compatibleGpuId
: compatibleGpus
)
505 if (availableGpuId
== compatibleGpuId
)
507 availableGpuIsCompatible
= true;
511 if (!availableGpuIsCompatible
)
513 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
);
515 gpuIdsToUse
.push_back(availableGpuId
);
519 if (fplog
!= nullptr)
521 /* Print references after all software/hardware printing */
522 please_cite(fplog
, "Abraham2015");
523 please_cite(fplog
, "Pall2015");
524 please_cite(fplog
, "Pronk2013");
525 please_cite(fplog
, "Hess2008b");
526 please_cite(fplog
, "Spoel2005a");
527 please_cite(fplog
, "Lindahl2001a");
528 please_cite(fplog
, "Berendsen95a");
529 writeSourceDoi(fplog
);
532 std::unique_ptr
<t_state
> globalState
;
536 /* Only the master rank has the global state */
537 globalState
= std::make_unique
<t_state
>();
539 /* Read (nearly) all data required for the simulation */
540 read_tpx_state(ftp2fn(efTPR
, filenames
.size(), filenames
.data()), inputrec
, globalState
.get(), &mtop
);
542 /* In rerun, set velocities to zero if present */
543 if (doRerun
&& ((globalState
->flags
& (1 << estV
)) != 0))
545 // rerun does not use velocities
546 GMX_LOG(mdlog
.info
).asParagraph().appendText(
547 "Rerun trajectory contains velocities. Rerun does only evaluate "
548 "potential energy and forces. The velocities will be ignored.");
549 for (int i
= 0; i
< globalState
->natoms
; i
++)
551 clear_rvec(globalState
->v
[i
]);
553 globalState
->flags
&= ~(1 << estV
);
556 if (inputrec
->cutoff_scheme
!= ecutsVERLET
)
558 if (nstlist_cmdline
> 0)
560 gmx_fatal(FARGS
, "Can not set nstlist with the group cut-off scheme");
563 if (!compatibleGpus
.empty())
565 GMX_LOG(mdlog
.warning
).asParagraph().appendText(
566 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
567 " To use a GPU, set the mdp option: cutoff-scheme = Verlet");
572 /* Check and update the hardware options for internal consistency */
573 check_and_update_hw_opt_1(mdlog
, &hw_opt
, cr
, domdecOptions
.numPmeRanks
);
575 /* Early check for externally set process affinity. */
576 gmx_check_thread_affinity_set(mdlog
, cr
,
577 &hw_opt
, hwinfo
->nthreads_hw_avail
, FALSE
);
579 if (GMX_THREAD_MPI
&& SIMMASTER(cr
))
581 if (domdecOptions
.numPmeRanks
> 0 && hw_opt
.nthreads_tmpi
<= 0)
583 gmx_fatal(FARGS
, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
586 /* Since the master knows the cut-off scheme, update hw_opt for this.
587 * This is done later for normal MPI and also once more with tMPI
588 * for all tMPI ranks.
590 check_and_update_hw_opt_2(&hw_opt
, inputrec
->cutoff_scheme
);
592 bool useGpuForNonbonded
= false;
593 bool useGpuForPme
= false;
596 // If the user specified the number of ranks, then we must
597 // respect that, but in default mode, we need to allow for
598 // the number of GPUs to choose the number of ranks.
599 auto canUseGpuForNonbonded
= buildSupportsNonbondedOnGpu(nullptr);
600 useGpuForNonbonded
= decideWhetherToUseGpusForNonbondedWithThreadMpi
601 (nonbondedTarget
, gpuIdsToUse
, userGpuTaskAssignment
, emulateGpuNonbonded
,
602 canUseGpuForNonbonded
,
603 inputrec
->cutoff_scheme
== ecutsVERLET
,
604 gpuAccelerationOfNonbondedIsUseful(mdlog
, inputrec
, GMX_THREAD_MPI
),
605 hw_opt
.nthreads_tmpi
);
606 auto canUseGpuForPme
= pme_gpu_supports_build(*hwinfo
, nullptr) && pme_gpu_supports_input(*inputrec
, mtop
, nullptr);
607 useGpuForPme
= decideWhetherToUseGpusForPmeWithThreadMpi
608 (useGpuForNonbonded
, pmeTarget
, gpuIdsToUse
, userGpuTaskAssignment
,
609 canUseGpuForPme
, hw_opt
.nthreads_tmpi
, domdecOptions
.numPmeRanks
);
612 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
614 /* Determine how many thread-MPI ranks to start.
616 * TODO Over-writing the user-supplied value here does
617 * prevent any possible subsequent checks from working
619 hw_opt
.nthreads_tmpi
= get_nthreads_mpi(hwinfo
,
628 // Now start the threads for thread MPI.
629 cr
= spawnThreads(hw_opt
.nthreads_tmpi
);
630 /* The main thread continues here with a new cr. We don't deallocate
631 the old cr because other threads may still be reading it. */
632 // TODO Both master and spawned threads call dup_tfn and
633 // reinitialize_commrec_for_this_thread. Find a way to express
635 physicalNodeComm
= PhysicalNodeCommunicator(MPI_COMM_WORLD
, gmx_physicalnode_id_hash());
637 // END OF CAUTION: cr and physicalNodeComm are now reliable
641 /* now broadcast everything to the non-master nodes/threads: */
642 init_parallel(cr
, inputrec
, &mtop
);
645 // Now each rank knows the inputrec that SIMMASTER read and used,
646 // and (if applicable) cr->nnodes has been assigned the number of
647 // thread-MPI ranks that have been chosen. The ranks can now all
648 // run the task-deciding functions and will agree on the result
649 // without needing to communicate.
651 // TODO Should we do the communication in debug mode to support
652 // having an assertion?
654 // Note that these variables describe only their own node.
656 // Note that when bonded interactions run on a GPU they always run
657 // alongside a nonbonded task, so do not influence task assignment
658 // even though they affect the force calculation workload.
659 bool useGpuForNonbonded
= false;
660 bool useGpuForPme
= false;
661 bool useGpuForBonded
= false;
664 // It's possible that there are different numbers of GPUs on
665 // different nodes, which is the user's responsibilty to
666 // handle. If unsuitable, we will notice that during task
668 bool gpusWereDetected
= hwinfo
->ngpu_compatible_tot
> 0;
669 bool usingVerletScheme
= inputrec
->cutoff_scheme
== ecutsVERLET
;
670 auto canUseGpuForNonbonded
= buildSupportsNonbondedOnGpu(nullptr);
671 useGpuForNonbonded
= decideWhetherToUseGpusForNonbonded(nonbondedTarget
, userGpuTaskAssignment
,
673 canUseGpuForNonbonded
,
675 gpuAccelerationOfNonbondedIsUseful(mdlog
, inputrec
, !GMX_THREAD_MPI
),
677 auto canUseGpuForPme
= pme_gpu_supports_build(*hwinfo
, nullptr) && pme_gpu_supports_input(*inputrec
, mtop
, nullptr);
678 useGpuForPme
= decideWhetherToUseGpusForPme(useGpuForNonbonded
, pmeTarget
, userGpuTaskAssignment
,
679 canUseGpuForPme
, cr
->nnodes
, domdecOptions
.numPmeRanks
,
681 auto canUseGpuForBonded
= buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec
, mtop
, nullptr);
683 decideWhetherToUseGpusForBonded(useGpuForNonbonded
, useGpuForPme
, usingVerletScheme
,
684 bondedTarget
, canUseGpuForBonded
,
685 EVDW_PME(inputrec
->vdwtype
),
686 EEL_PME_EWALD(inputrec
->coulombtype
),
687 domdecOptions
.numPmeRanks
, gpusWereDetected
);
689 pmeRunMode
= (useGpuForPme
? PmeRunMode::GPU
: PmeRunMode::CPU
);
690 if (pmeRunMode
== PmeRunMode::GPU
)
692 if (pmeFftTarget
== TaskTarget::Cpu
)
694 pmeRunMode
= PmeRunMode::Mixed
;
697 else if (pmeFftTarget
== TaskTarget::Gpu
)
699 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.");
702 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
705 // TODO: hide restraint implementation details from Mdrunner.
706 // There is nothing unique about restraints at this point as far as the
707 // Mdrunner is concerned. The Mdrunner should just be getting a sequence of
708 // factory functions from the SimulationContext on which to call mdModules->add().
709 // TODO: capture all restraints into a single RestraintModule, passed to the runner builder.
710 for (auto && restraint
: restraintManager_
->getRestraints())
712 auto module
= RestraintMDModule::create(restraint
,
714 mdModules
->add(std::move(module
));
717 // TODO: Error handling
718 mdModules
->assignOptionsToModules(*inputrec
->params
, nullptr);
720 if (fplog
!= nullptr)
722 pr_inputrec(fplog
, 0, "Input Parameters", inputrec
, FALSE
);
723 fprintf(fplog
, "\n");
728 /* now make sure the state is initialized and propagated */
729 set_state_entries(globalState
.get(), inputrec
);
732 /* NM and TPI parallelize over force/energy calculations, not atoms,
733 * so we need to initialize and broadcast the global state.
735 if (inputrec
->eI
== eiNM
|| inputrec
->eI
== eiTPI
)
739 globalState
= std::make_unique
<t_state
>();
741 broadcastStateWithoutDynamics(cr
, globalState
.get());
744 /* A parallel command line option consistency check that we can
745 only do after any threads have started. */
746 if (!PAR(cr
) && (domdecOptions
.numCells
[XX
] > 1 ||
747 domdecOptions
.numCells
[YY
] > 1 ||
748 domdecOptions
.numCells
[ZZ
] > 1 ||
749 domdecOptions
.numPmeRanks
> 0))
752 "The -dd or -npme option request a parallel simulation, "
754 "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv
));
757 "but the number of MPI-threads (option -ntmpi) is not set or is 1");
759 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv
));
765 (EI_ENERGY_MINIMIZATION(inputrec
->eI
) || eiNM
== inputrec
->eI
))
767 gmx_fatal(FARGS
, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
770 if (!(EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
)))
772 if (domdecOptions
.numPmeRanks
> 0)
774 gmx_fatal_collective(FARGS
, cr
->mpi_comm_mysim
, MASTER(cr
),
775 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
778 domdecOptions
.numPmeRanks
= 0;
781 if (useGpuForNonbonded
&& domdecOptions
.numPmeRanks
< 0)
783 /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
784 * improve performance with many threads per GPU, since our OpenMP
785 * scaling is bad, but it's difficult to automate the setup.
787 domdecOptions
.numPmeRanks
= 0;
791 if (domdecOptions
.numPmeRanks
< 0)
793 domdecOptions
.numPmeRanks
= 0;
794 // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
798 GMX_RELEASE_ASSERT(domdecOptions
.numPmeRanks
<= 1, "PME GPU decomposition is not supported");
805 fcRegisterSteps(inputrec
->nsteps
, inputrec
->init_step
);
809 /* NMR restraints must be initialized before load_checkpoint,
810 * since with time averaging the history is added to t_state.
811 * For proper consistency check we therefore need to extend
813 * So the PME-only nodes (if present) will also initialize
814 * the distance restraints.
818 /* This needs to be called before read_checkpoint to extend the state */
819 init_disres(fplog
, &mtop
, inputrec
, cr
, ms
, fcd
, globalState
.get(), replExParams
.exchangeInterval
> 0);
821 init_orires(fplog
, &mtop
, inputrec
, cr
, ms
, globalState
.get(), &(fcd
->orires
));
823 auto deform
= prepareBoxDeformation(globalState
->box
, cr
, *inputrec
);
825 ObservablesHistory observablesHistory
= {};
827 ContinuationOptions
&continuationOptions
= mdrunOptions
.continuationOptions
;
829 if (continuationOptions
.startedFromCheckpoint
)
831 /* Check if checkpoint file exists before doing continuation.
832 * This way we can use identical input options for the first and subsequent runs...
836 load_checkpoint(opt2fn_master("-cpi", filenames
.size(), filenames
.data(), cr
),
838 cr
, domdecOptions
.numCells
,
839 inputrec
, globalState
.get(),
840 &bReadEkin
, &observablesHistory
,
841 continuationOptions
.appendFiles
,
842 continuationOptions
.appendFilesOptionSet
,
843 mdrunOptions
.reproducible
);
847 continuationOptions
.haveReadEkin
= true;
850 if (continuationOptions
.appendFiles
&& logFileHandle
)
852 // Now we can start normal logging to the truncated log file.
853 fplog
= gmx_fio_getfp(logFileHandle
);
854 prepareLogAppending(fplog
);
855 logOwner
= buildLogger(fplog
, cr
);
856 mdlog
= logOwner
.logger();
860 if (mdrunOptions
.numStepsCommandline
> -2)
862 GMX_LOG(mdlog
.info
).asParagraph().
863 appendText("The -nsteps functionality is deprecated, and may be removed in a future version. "
864 "Consider using gmx convert-tpr -nsteps or changing the appropriate .mdp file field.");
866 /* override nsteps with value set on the commamdline */
867 override_nsteps_cmdline(mdlog
, mdrunOptions
.numStepsCommandline
, inputrec
);
871 copy_mat(globalState
->box
, box
);
876 gmx_bcast(sizeof(box
), box
, cr
);
879 /* Update rlist and nstlist. */
880 if (inputrec
->cutoff_scheme
== ecutsVERLET
)
882 prepare_verlet_scheme(fplog
, cr
, inputrec
, nstlist_cmdline
, &mtop
, box
,
883 useGpuForNonbonded
|| (emulateGpuNonbonded
== EmulateGpuNonbonded::Yes
), *hwinfo
->cpuInfo
);
886 LocalAtomSetManager atomSets
;
888 if (PAR(cr
) && !(EI_TPI(inputrec
->eI
) ||
889 inputrec
->eI
== eiNM
))
891 cr
->dd
= init_domain_decomposition(mdlog
, cr
, domdecOptions
, mdrunOptions
,
893 box
, positionsFromStatePointer(globalState
.get()),
895 // Note that local state still does not exist yet.
899 /* PME, if used, is done on all nodes with 1D decomposition */
901 cr
->duty
= (DUTY_PP
| DUTY_PME
);
903 if (inputrec
->ePBC
== epbcSCREW
)
906 "pbc=screw is only implemented with domain decomposition");
912 /* After possible communicator splitting in make_dd_communicators.
913 * we can set up the intra/inter node communication.
915 gmx_setup_nodecomm(fplog
, cr
);
921 GMX_LOG(mdlog
.warning
).asParagraph().appendTextFormatted(
922 "This is simulation %d out of %d running as a composite GROMACS\n"
923 "multi-simulation job. Setup for this simulation:\n",
926 GMX_LOG(mdlog
.warning
).appendTextFormatted(
930 cr
->nnodes
== 1 ? "thread" : "threads"
932 cr
->nnodes
== 1 ? "process" : "processes"
938 /* Check and update hw_opt for the cut-off scheme */
939 check_and_update_hw_opt_2(&hw_opt
, inputrec
->cutoff_scheme
);
941 /* Check and update the number of OpenMP threads requested */
942 checkAndUpdateRequestedNumOpenmpThreads(&hw_opt
, *hwinfo
, cr
, ms
, physicalNodeComm
.size_
,
945 gmx_omp_nthreads_init(mdlog
, cr
,
946 hwinfo
->nthreads_hw_avail
,
947 physicalNodeComm
.size_
,
949 hw_opt
.nthreads_omp_pme
,
950 !thisRankHasDuty(cr
, DUTY_PP
),
951 inputrec
->cutoff_scheme
== ecutsVERLET
);
953 // Enable FP exception detection for the Verlet scheme, but not in
954 // Release mode and not for compilers with known buggy FP
955 // exception support (clang with any optimization) or suspected
956 // buggy FP exception support (gcc 7.* with optimization).
957 #if !defined NDEBUG && \
958 !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
959 && defined __OPTIMIZE__)
960 const bool bEnableFPE
= inputrec
->cutoff_scheme
== ecutsVERLET
;
962 const bool bEnableFPE
= false;
964 //FIXME - reconcile with gmx_feenableexcept() call from CommandLineModuleManager::run()
967 gmx_feenableexcept();
970 // Build a data structure that expresses which kinds of non-bonded
971 // task are handled by this rank.
973 // TODO Later, this might become a loop over all registered modules
974 // relevant to the mdp inputs, to find those that have such tasks.
976 // TODO This could move before init_domain_decomposition() as part
977 // of refactoring that separates the responsibility for duty
978 // assignment from setup for communication between tasks, and
979 // setup for tasks handled with a domain (ie including short-ranged
980 // tasks, bonded tasks, etc.).
982 // Note that in general useGpuForNonbonded, etc. can have a value
983 // that is inconsistent with the presence of actual GPUs on any
984 // rank, and that is not known to be a problem until the
985 // duty of the ranks on a node become known.
987 // TODO Later we might need the concept of computeTasksOnThisRank,
988 // from which we construct gpuTasksOnThisRank.
990 // Currently the DD code assigns duty to ranks that can
991 // include PP work that currently can be executed on a single
992 // GPU, if present and compatible. This has to be coordinated
993 // across PP ranks on a node, with possible multiple devices
994 // or sharing devices on a node, either from the user
995 // selection, or automatically.
996 auto haveGpus
= !gpuIdsToUse
.empty();
997 std::vector
<GpuTask
> gpuTasksOnThisRank
;
998 if (thisRankHasDuty(cr
, DUTY_PP
))
1000 if (useGpuForNonbonded
)
1002 // Note that any bonded tasks on a GPU always accompany a
1006 gpuTasksOnThisRank
.push_back(GpuTask::Nonbonded
);
1008 else if (nonbondedTarget
== TaskTarget::Gpu
)
1010 gmx_fatal(FARGS
, "Cannot run short-ranged nonbonded interactions on a GPU because there is none detected.");
1012 else if (bondedTarget
== TaskTarget::Gpu
)
1014 gmx_fatal(FARGS
, "Cannot run bonded interactions on a GPU because there is none detected.");
1018 // TODO cr->duty & DUTY_PME should imply that a PME algorithm is active, but currently does not.
1019 if (EEL_PME(inputrec
->coulombtype
) && (thisRankHasDuty(cr
, DUTY_PME
)))
1025 gpuTasksOnThisRank
.push_back(GpuTask::Pme
);
1027 else if (pmeTarget
== TaskTarget::Gpu
)
1029 gmx_fatal(FARGS
, "Cannot run PME on a GPU because there is none detected.");
1034 GpuTaskAssignment gpuTaskAssignment
;
1037 // Produce the task assignment for this rank.
1038 gpuTaskAssignment
= runTaskAssignment(gpuIdsToUse
, userGpuTaskAssignment
, *hwinfo
,
1039 mdlog
, cr
, ms
, physicalNodeComm
, gpuTasksOnThisRank
,
1040 useGpuForBonded
, pmeRunMode
);
1042 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1044 /* Prevent other ranks from continuing after an issue was found
1045 * and reported as a fatal error.
1047 * TODO This function implements a barrier so that MPI runtimes
1048 * can organize an orderly shutdown if one of the ranks has had to
1049 * issue a fatal error in various code already run. When we have
1050 * MPI-aware error handling and reporting, this should be
1055 MPI_Barrier(cr
->mpi_comm_mysim
);
1061 MPI_Barrier(ms
->mpi_comm_masters
);
1063 /* We need another barrier to prevent non-master ranks from contiuing
1064 * when an error occured in a different simulation.
1066 MPI_Barrier(cr
->mpi_comm_mysim
);
1070 /* Now that we know the setup is consistent, check for efficiency */
1071 check_resource_division_efficiency(hwinfo
, !gpuTaskAssignment
.empty(), mdrunOptions
.ntompOptionIsSet
,
1074 gmx_device_info_t
*nonbondedDeviceInfo
= nullptr;
1076 if (thisRankHasDuty(cr
, DUTY_PP
))
1078 // This works because only one task of each type is currently permitted.
1079 auto nbGpuTaskMapping
= std::find_if(gpuTaskAssignment
.begin(), gpuTaskAssignment
.end(),
1080 hasTaskType
<GpuTask::Nonbonded
>);
1081 if (nbGpuTaskMapping
!= gpuTaskAssignment
.end())
1083 int nonbondedDeviceId
= nbGpuTaskMapping
->deviceId_
;
1084 nonbondedDeviceInfo
= getDeviceInfo(hwinfo
->gpu_info
, nonbondedDeviceId
);
1085 init_gpu(nonbondedDeviceInfo
);
1087 if (DOMAINDECOMP(cr
))
1089 /* When we share GPUs over ranks, we need to know this for the DLB */
1090 dd_setup_dlb_resource_sharing(cr
, nonbondedDeviceId
);
1096 std::unique_ptr
<ClfftInitializer
> initializedClfftLibrary
;
1098 gmx_device_info_t
*pmeDeviceInfo
= nullptr;
1099 // Later, this program could contain kernels that might be later
1100 // re-used as auto-tuning progresses, or subsequent simulations
1102 PmeGpuProgramStorage pmeGpuProgram
;
1103 // This works because only one task of each type is currently permitted.
1104 auto pmeGpuTaskMapping
= std::find_if(gpuTaskAssignment
.begin(), gpuTaskAssignment
.end(), hasTaskType
<GpuTask::Pme
>);
1105 const bool thisRankHasPmeGpuTask
= (pmeGpuTaskMapping
!= gpuTaskAssignment
.end());
1106 if (thisRankHasPmeGpuTask
)
1108 pmeDeviceInfo
= getDeviceInfo(hwinfo
->gpu_info
, pmeGpuTaskMapping
->deviceId_
);
1109 init_gpu(pmeDeviceInfo
);
1110 pmeGpuProgram
= buildPmeGpuProgram(pmeDeviceInfo
);
1111 // TODO It would be nice to move this logic into the factory
1112 // function. See Redmine #2535.
1113 bool isMasterThread
= !GMX_THREAD_MPI
|| MASTER(cr
);
1114 if (pmeRunMode
== PmeRunMode::GPU
&& !initializedClfftLibrary
&& isMasterThread
)
1116 initializedClfftLibrary
= initializeClfftLibrary();
1120 /* getting number of PP/PME threads on this MPI / tMPI rank.
1121 PME: env variable should be read only on one node to make sure it is
1122 identical everywhere;
1124 const int numThreadsOnThisRank
=
1125 thisRankHasDuty(cr
, DUTY_PP
) ? gmx_omp_nthreads_get(emntNonbonded
) : gmx_omp_nthreads_get(emntPME
);
1126 checkHardwareOversubscription(numThreadsOnThisRank
, cr
->nodeid
,
1127 *hwinfo
->hardwareTopology
,
1128 physicalNodeComm
, mdlog
);
1130 if (hw_opt
.thread_affinity
!= threadaffOFF
)
1132 /* Before setting affinity, check whether the affinity has changed
1133 * - which indicates that probably the OpenMP library has changed it
1134 * since we first checked).
1136 gmx_check_thread_affinity_set(mdlog
, cr
,
1137 &hw_opt
, hwinfo
->nthreads_hw_avail
, TRUE
);
1139 int numThreadsOnThisNode
, intraNodeThreadOffset
;
1140 analyzeThreadsOnThisNode(physicalNodeComm
, numThreadsOnThisRank
, &numThreadsOnThisNode
,
1141 &intraNodeThreadOffset
);
1143 /* Set the CPU affinity */
1144 gmx_set_thread_affinity(mdlog
, cr
, &hw_opt
, *hwinfo
->hardwareTopology
,
1145 numThreadsOnThisRank
, numThreadsOnThisNode
,
1146 intraNodeThreadOffset
, nullptr);
1149 if (mdrunOptions
.timingOptions
.resetStep
> -1)
1151 GMX_LOG(mdlog
.info
).asParagraph().
1152 appendText("The -resetstep functionality is deprecated, and may be removed in a future version.");
1154 wcycle
= wallcycle_init(fplog
, mdrunOptions
.timingOptions
.resetStep
, cr
);
1158 /* Master synchronizes its value of reset_counters with all nodes
1159 * including PME only nodes */
1160 int64_t reset_counters
= wcycle_get_reset_counters(wcycle
);
1161 gmx_bcast_sim(sizeof(reset_counters
), &reset_counters
, cr
);
1162 wcycle_set_reset_counters(wcycle
, reset_counters
);
1165 // Membrane embedding must be initialized before we call init_forcerec()
1170 fprintf(stderr
, "Initializing membed");
1172 /* Note that membed cannot work in parallel because mtop is
1173 * changed here. Fix this if we ever want to make it run with
1174 * multiple ranks. */
1175 membed
= init_membed(fplog
, filenames
.size(), filenames
.data(), &mtop
, inputrec
, globalState
.get(), cr
,
1177 .checkpointOptions
.period
);
1180 std::unique_ptr
<MDAtoms
> mdAtoms
;
1181 std::unique_ptr
<gmx_vsite_t
> vsite
;
1184 if (thisRankHasDuty(cr
, DUTY_PP
))
1186 /* Initiate forcerecord */
1188 fr
->forceProviders
= mdModules
->initForceProviders();
1189 init_forcerec(fplog
, mdlog
, fr
, fcd
,
1190 inputrec
, &mtop
, cr
, box
,
1191 opt2fn("-table", filenames
.size(), filenames
.data()),
1192 opt2fn("-tablep", filenames
.size(), filenames
.data()),
1193 opt2fns("-tableb", filenames
.size(), filenames
.data()),
1194 *hwinfo
, nonbondedDeviceInfo
,
1199 /* Initialize the mdAtoms structure.
1200 * mdAtoms is not filled with atom data,
1201 * as this can not be done now with domain decomposition.
1203 mdAtoms
= makeMDAtoms(fplog
, mtop
, *inputrec
, thisRankHasPmeGpuTask
);
1204 if (globalState
&& thisRankHasPmeGpuTask
)
1206 // The pinning of coordinates in the global state object works, because we only use
1207 // PME on GPU without DD or on a separate PME rank, and because the local state pointer
1208 // points to the global state object without DD.
1209 // FIXME: MD and EM separately set up the local state - this should happen in the same function,
1210 // which should also perform the pinning.
1211 changePinningPolicy(&globalState
->x
, pme_get_pinning_policy());
1214 /* Initialize the virtual site communication */
1215 vsite
= initVsite(mtop
, cr
);
1217 calc_shifts(box
, fr
->shift_vec
);
1219 /* With periodic molecules the charge groups should be whole at start up
1220 * and the virtual sites should not be far from their proper positions.
1222 if (!inputrec
->bContinuation
&& MASTER(cr
) &&
1223 !(inputrec
->ePBC
!= epbcNONE
&& inputrec
->bPeriodicMols
))
1225 /* Make molecules whole at start of run */
1226 if (fr
->ePBC
!= epbcNONE
)
1228 do_pbc_first_mtop(fplog
, inputrec
->ePBC
, box
, &mtop
, globalState
->x
.rvec_array());
1232 /* Correct initial vsite positions are required
1233 * for the initial distribution in the domain decomposition
1234 * and for the initial shell prediction.
1236 constructVsitesGlobal(mtop
, globalState
->x
);
1240 if (EEL_PME(fr
->ic
->eeltype
) || EVDW_PME(fr
->ic
->vdwtype
))
1242 ewaldcoeff_q
= fr
->ic
->ewaldcoeff_q
;
1243 ewaldcoeff_lj
= fr
->ic
->ewaldcoeff_lj
;
1248 /* This is a PME only node */
1250 GMX_ASSERT(globalState
== nullptr, "We don't need the state on a PME only rank and expect it to be unitialized");
1252 ewaldcoeff_q
= calc_ewaldcoeff_q(inputrec
->rcoulomb
, inputrec
->ewald_rtol
);
1253 ewaldcoeff_lj
= calc_ewaldcoeff_lj(inputrec
->rvdw
, inputrec
->ewald_rtol_lj
);
1256 gmx_pme_t
*sepPmeData
= nullptr;
1257 // This reference hides the fact that PME data is owned by runner on PME-only ranks and by forcerec on other ranks
1258 GMX_ASSERT(thisRankHasDuty(cr
, DUTY_PP
) == (fr
!= nullptr), "Double-checking that only PME-only ranks have no forcerec");
1259 gmx_pme_t
* &pmedata
= fr
? fr
->pmedata
: sepPmeData
;
1261 /* Initiate PME if necessary,
1262 * either on all nodes or on dedicated PME nodes only. */
1263 if (EEL_PME(inputrec
->coulombtype
) || EVDW_PME(inputrec
->vdwtype
))
1265 if (mdAtoms
&& mdAtoms
->mdatoms())
1267 nChargePerturbed
= mdAtoms
->mdatoms()->nChargePerturbed
;
1268 if (EVDW_PME(inputrec
->vdwtype
))
1270 nTypePerturbed
= mdAtoms
->mdatoms()->nTypePerturbed
;
1273 if (cr
->npmenodes
> 0)
1275 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1276 gmx_bcast_sim(sizeof(nChargePerturbed
), &nChargePerturbed
, cr
);
1277 gmx_bcast_sim(sizeof(nTypePerturbed
), &nTypePerturbed
, cr
);
1280 if (thisRankHasDuty(cr
, DUTY_PME
))
1284 pmedata
= gmx_pme_init(cr
,
1285 getNumPmeDomains(cr
->dd
),
1287 mtop
.natoms
, nChargePerturbed
!= 0, nTypePerturbed
!= 0,
1288 mdrunOptions
.reproducible
,
1289 ewaldcoeff_q
, ewaldcoeff_lj
,
1290 gmx_omp_nthreads_get(emntPME
),
1291 pmeRunMode
, nullptr,
1292 pmeDeviceInfo
, pmeGpuProgram
.get(), mdlog
);
1294 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1299 if (EI_DYNAMICS(inputrec
->eI
))
1301 /* Turn on signal handling on all nodes */
1303 * (A user signal from the PME nodes (if any)
1304 * is communicated to the PP nodes.
1306 signal_handler_install();
1309 if (thisRankHasDuty(cr
, DUTY_PP
))
1311 /* Assumes uniform use of the number of OpenMP threads */
1312 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntDefault
));
1314 if (inputrec
->bPull
)
1316 /* Initialize pull code */
1317 inputrec
->pull_work
=
1318 init_pull(fplog
, inputrec
->pull
, inputrec
,
1319 &mtop
, cr
, &atomSets
, inputrec
->fepvals
->init_lambda
);
1320 if (inputrec
->pull
->bXOutAverage
|| inputrec
->pull
->bFOutAverage
)
1322 initPullHistory(inputrec
->pull_work
, &observablesHistory
);
1324 if (EI_DYNAMICS(inputrec
->eI
) && MASTER(cr
))
1326 init_pull_output_files(inputrec
->pull_work
,
1327 filenames
.size(), filenames
.data(), oenv
,
1328 continuationOptions
);
1332 std::unique_ptr
<EnforcedRotation
> enforcedRotation
;
1335 /* Initialize enforced rotation code */
1336 enforcedRotation
= init_rot(fplog
,
1348 if (inputrec
->eSwapCoords
!= eswapNO
)
1350 /* Initialize ion swapping code */
1351 init_swapcoords(fplog
, inputrec
, opt2fn_master("-swap", filenames
.size(), filenames
.data(), cr
),
1352 &mtop
, globalState
.get(), &observablesHistory
,
1353 cr
, &atomSets
, oenv
, mdrunOptions
);
1356 /* Let makeConstraints know whether we have essential dynamics constraints.
1357 * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1359 bool doEssentialDynamics
= (opt2fn_null("-ei", filenames
.size(), filenames
.data()) != nullptr
1360 || observablesHistory
.edsamHistory
);
1361 auto constr
= makeConstraints(mtop
, *inputrec
, doEssentialDynamics
,
1362 fplog
, *mdAtoms
->mdatoms(),
1363 cr
, ms
, nrnb
, wcycle
, fr
->bMolPBC
);
1365 if (DOMAINDECOMP(cr
))
1367 GMX_RELEASE_ASSERT(fr
, "fr was NULL while cr->duty was DUTY_PP");
1368 /* This call is not included in init_domain_decomposition mainly
1369 * because fr->cginfo_mb is set later.
1371 dd_init_bondeds(fplog
, cr
->dd
, &mtop
, vsite
.get(), inputrec
,
1372 domdecOptions
.checkBondedInteractions
,
1376 // TODO This is not the right place to manage the lifetime of
1377 // this data structure, but currently it's the easiest way to
1378 // make it work. Later, it should probably be made/updated
1379 // after the workload for the lifetime of a PP domain is
1381 PpForceWorkload ppForceWorkload
;
1383 GMX_ASSERT(stopHandlerBuilder_
, "Runner must provide StopHandlerBuilder to integrator.");
1384 /* Now do whatever the user wants us to do (how flexible...) */
1385 Integrator integrator
{
1386 fplog
, cr
, ms
, mdlog
, static_cast<int>(filenames
.size()), filenames
.data(),
1389 vsite
.get(), constr
.get(),
1390 enforcedRotation
? enforcedRotation
->getLegacyEnfrot() : nullptr,
1392 mdModules
->outputProvider(),
1396 &observablesHistory
,
1397 mdAtoms
.get(), nrnb
, wcycle
, fr
,
1401 walltime_accounting
,
1402 std::move(stopHandlerBuilder_
)
1404 integrator
.run(inputrec
->eI
, doRerun
);
1406 if (inputrec
->bPull
)
1408 finish_pull(inputrec
->pull_work
);
1414 GMX_RELEASE_ASSERT(pmedata
, "pmedata was NULL while cr->duty was not DUTY_PP");
1416 walltime_accounting
= walltime_accounting_init(gmx_omp_nthreads_get(emntPME
));
1417 gmx_pmeonly(pmedata
, cr
, nrnb
, wcycle
, walltime_accounting
, inputrec
, pmeRunMode
);
1420 wallcycle_stop(wcycle
, ewcRUN
);
1422 /* Finish up, write some stuff
1423 * if rerunMD, don't write last frame again
1425 finish_run(fplog
, mdlog
, cr
,
1426 inputrec
, nrnb
, wcycle
, walltime_accounting
,
1427 fr
? fr
->nbv
: nullptr,
1429 EI_DYNAMICS(inputrec
->eI
) && !isMultiSim(ms
));
1431 // clean up cycle counter
1432 wallcycle_destroy(wcycle
);
1437 gmx_pme_destroy(pmedata
);
1441 // FIXME: this is only here to manually unpin mdAtoms->chargeA_ and state->x,
1442 // before we destroy the GPU context(s) in free_gpu_resources().
1443 // Pinned buffers are associated with contexts in CUDA.
1444 // As soon as we destroy GPU contexts after mdrunner() exits, these lines should go.
1445 mdAtoms
.reset(nullptr);
1446 globalState
.reset(nullptr);
1447 mdModules
.reset(nullptr); // destruct force providers here as they might also use the GPU
1449 /* Free GPU memory and set a physical node tMPI barrier (which should eventually go away) */
1450 free_gpu_resources(fr
, physicalNodeComm
);
1451 free_gpu(nonbondedDeviceInfo
);
1452 free_gpu(pmeDeviceInfo
);
1453 done_forcerec(fr
, mtop
.molblock
.size(), mtop
.groups
.grps
[egcENER
].nr
);
1458 free_membed(membed
);
1461 gmx_hardware_info_free();
1463 /* Does what it says */
1464 print_date_and_time(fplog
, cr
->nodeid
, "Finished mdrun", gmx_gettime());
1465 walltime_accounting_destroy(walltime_accounting
);
1468 // Ensure log file content is written
1471 gmx_fio_flush(logFileHandle
);
1474 /* Reset FPEs (important for unit tests) by disabling them. Assumes no
1475 * exceptions were enabled before function was called. */
1478 gmx_fedisableexcept();
1481 auto rc
= static_cast<int>(gmx_get_stop_condition());
1484 /* we need to join all threads. The sub-threads join when they
1485 exit this function, but the master thread needs to be told to
1487 if (PAR(cr
) && MASTER(cr
))
1491 //TODO free commrec in MPI simulations
1497 Mdrunner::~Mdrunner()
1499 // Clean up of the Manager.
1500 // This will end up getting called on every thread-MPI rank, which is unnecessary,
1501 // but okay as long as threads synchronize some time before adding or accessing
1502 // a new set of restraints.
1503 if (restraintManager_
)
1505 restraintManager_
->clear();
1506 GMX_ASSERT(restraintManager_
->countRestraints() == 0,
1507 "restraints added during runner life time should be cleared at runner destruction.");
1511 void Mdrunner::addPotential(std::shared_ptr
<gmx::IRestraintPotential
> puller
,
1514 GMX_ASSERT(restraintManager_
, "Mdrunner must have a restraint manager.");
1515 // Not sure if this should be logged through the md logger or something else,
1516 // but it is helpful to have some sort of INFO level message sent somewhere.
1517 // std::cout << "Registering restraint named " << name << std::endl;
1519 // When multiple restraints are used, it may be wasteful to register them separately.
1520 // Maybe instead register an entire Restraint Manager as a force provider.
1521 restraintManager_
->addToSpec(std::move(puller
),
1525 Mdrunner::Mdrunner(Mdrunner
&&) noexcept
= default;
1527 //NOLINTNEXTLINE(performance-noexcept-move-constructor) working around GCC bug 58265
1528 Mdrunner
&Mdrunner::operator=(Mdrunner
&& /*handle*/) noexcept(BUGFREE_NOEXCEPT_STRING
) = default;
1530 class Mdrunner::BuilderImplementation
1533 BuilderImplementation() = delete;
1534 explicit BuilderImplementation(SimulationContext
* context
);
1535 ~BuilderImplementation();
1537 BuilderImplementation
&setExtraMdrunOptions(const MdrunOptions
&options
,
1538 real forceWarningThreshold
);
1540 void addDomdec(const DomdecOptions
&options
);
1542 void addVerletList(int nstlist
);
1544 void addReplicaExchange(const ReplicaExchangeParameters
¶ms
);
1546 void addMultiSim(gmx_multisim_t
* multisim
);
1548 void addNonBonded(const char* nbpu_opt
);
1550 void addPME(const char* pme_opt_
, const char* pme_fft_opt_
);
1552 void addBondedTaskAssignment(const char* bonded_opt
);
1554 void addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
);
1556 void addFilenames(ArrayRef
<const t_filenm
> filenames
);
1558 void addOutputEnvironment(gmx_output_env_t
* outputEnvironment
);
1560 void addLogFile(t_fileio
*logFileHandle
);
1562 void addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
);
1567 // Default parameters copied from runner.h
1568 // \todo Clarify source(s) of default parameters.
1570 const char* nbpu_opt_
= nullptr;
1571 const char* pme_opt_
= nullptr;
1572 const char* pme_fft_opt_
= nullptr;
1573 const char *bonded_opt_
= nullptr;
1575 MdrunOptions mdrunOptions_
;
1577 DomdecOptions domdecOptions_
;
1579 ReplicaExchangeParameters replicaExchangeParameters_
;
1581 //! Command-line override for the duration of a neighbor list with the Verlet scheme.
1584 //! Non-owning multisim communicator handle.
1585 std::unique_ptr
<gmx_multisim_t
*> multisim_
= nullptr;
1587 //! Print a warning if any force is larger than this (in kJ/mol nm).
1588 real forceWarningThreshold_
= -1;
1590 /*! \brief Non-owning pointer to SimulationContext (owned and managed by client)
1593 * \todo Establish robust protocol to make sure resources remain valid.
1594 * SimulationContext will likely be separated into multiple layers for
1595 * different levels of access and different phases of execution. Ref
1596 * https://redmine.gromacs.org/issues/2375
1597 * https://redmine.gromacs.org/issues/2587
1599 SimulationContext
* context_
= nullptr;
1601 //! \brief Parallelism information.
1602 gmx_hw_opt_t hardwareOptions_
;
1604 //! filename options for simulation.
1605 ArrayRef
<const t_filenm
> filenames_
;
1607 /*! \brief Handle to output environment.
1609 * \todo gmx_output_env_t needs lifetime management.
1611 gmx_output_env_t
* outputEnvironment_
= nullptr;
1613 /*! \brief Non-owning handle to MD log file.
1615 * \todo Context should own output facilities for client.
1616 * \todo Improve log file handle management.
1618 * Code managing the FILE* relies on the ability to set it to
1619 * nullptr to check whether the filehandle is valid.
1621 t_fileio
* logFileHandle_
= nullptr;
1624 * \brief Builder for simulation stop signal handler.
1626 std::unique_ptr
<StopHandlerBuilder
> stopHandlerBuilder_
= nullptr;
1629 Mdrunner::BuilderImplementation::BuilderImplementation(SimulationContext
* context
) :
1632 GMX_ASSERT(context_
, "Bug found. It should not be possible to construct builder without a valid context.");
1635 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
1637 Mdrunner::BuilderImplementation
&
1638 Mdrunner::BuilderImplementation::setExtraMdrunOptions(const MdrunOptions
&options
,
1639 real forceWarningThreshold
)
1641 mdrunOptions_
= options
;
1642 forceWarningThreshold_
= forceWarningThreshold
;
1646 void Mdrunner::BuilderImplementation::addDomdec(const DomdecOptions
&options
)
1648 domdecOptions_
= options
;
1651 void Mdrunner::BuilderImplementation::addVerletList(int nstlist
)
1656 void Mdrunner::BuilderImplementation::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1658 replicaExchangeParameters_
= params
;
1661 void Mdrunner::BuilderImplementation::addMultiSim(gmx_multisim_t
* multisim
)
1663 multisim_
= std::make_unique
<gmx_multisim_t
*>(multisim
);
1666 Mdrunner
Mdrunner::BuilderImplementation::build()
1668 auto newRunner
= Mdrunner();
1670 GMX_ASSERT(context_
, "Bug found. It should not be possible to call build() without a valid context.");
1672 newRunner
.mdrunOptions
= mdrunOptions_
;
1673 newRunner
.domdecOptions
= domdecOptions_
;
1675 // \todo determine an invariant to check or confirm that all gmx_hw_opt_t objects are valid
1676 newRunner
.hw_opt
= hardwareOptions_
;
1678 // No invariant to check. This parameter exists to optionally override other behavior.
1679 newRunner
.nstlist_cmdline
= nstlist_
;
1681 newRunner
.replExParams
= replicaExchangeParameters_
;
1683 newRunner
.filenames
= filenames_
;
1685 GMX_ASSERT(context_
->communicationRecord_
, "SimulationContext communications not initialized.");
1686 newRunner
.cr
= context_
->communicationRecord_
;
1690 // nullptr is a valid value for the multisim handle, so we don't check the pointed-to pointer.
1691 newRunner
.ms
= *multisim_
;
1695 GMX_THROW(gmx::APIError("MdrunnerBuilder::addMultiSim() is required before build()"));
1698 // \todo Clarify ownership and lifetime management for gmx_output_env_t
1699 // \todo Update sanity checking when output environment has clearly specified invariants.
1700 // Initialization and default values for oenv are not well specified in the current version.
1701 if (outputEnvironment_
)
1703 newRunner
.oenv
= outputEnvironment_
;
1707 GMX_THROW(gmx::APIError("MdrunnerBuilder::addOutputEnvironment() is required before build()"));
1710 newRunner
.logFileHandle
= logFileHandle_
;
1714 newRunner
.nbpu_opt
= nbpu_opt_
;
1718 GMX_THROW(gmx::APIError("MdrunnerBuilder::addNonBonded() is required before build()"));
1721 if (pme_opt_
&& pme_fft_opt_
)
1723 newRunner
.pme_opt
= pme_opt_
;
1724 newRunner
.pme_fft_opt
= pme_fft_opt_
;
1728 GMX_THROW(gmx::APIError("MdrunnerBuilder::addElectrostatics() is required before build()"));
1733 newRunner
.bonded_opt
= bonded_opt_
;
1737 GMX_THROW(gmx::APIError("MdrunnerBuilder::addBondedTaskAssignment() is required before build()"));
1740 newRunner
.restraintManager_
= std::make_unique
<gmx::RestraintManager
>();
1742 if (stopHandlerBuilder_
)
1744 newRunner
.stopHandlerBuilder_
= std::move(stopHandlerBuilder_
);
1748 newRunner
.stopHandlerBuilder_
= std::make_unique
<StopHandlerBuilder
>();
1754 void Mdrunner::BuilderImplementation::addNonBonded(const char* nbpu_opt
)
1756 nbpu_opt_
= nbpu_opt
;
1759 void Mdrunner::BuilderImplementation::addPME(const char* pme_opt
,
1760 const char* pme_fft_opt
)
1763 pme_fft_opt_
= pme_fft_opt
;
1766 void Mdrunner::BuilderImplementation::addBondedTaskAssignment(const char* bonded_opt
)
1768 bonded_opt_
= bonded_opt
;
1771 void Mdrunner::BuilderImplementation::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
1773 hardwareOptions_
= hardwareOptions
;
1776 void Mdrunner::BuilderImplementation::addFilenames(ArrayRef
<const t_filenm
> filenames
)
1778 filenames_
= filenames
;
1781 void Mdrunner::BuilderImplementation::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
1783 outputEnvironment_
= outputEnvironment
;
1786 void Mdrunner::BuilderImplementation::addLogFile(t_fileio
*logFileHandle
)
1788 logFileHandle_
= logFileHandle
;
1791 void Mdrunner::BuilderImplementation::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
1793 stopHandlerBuilder_
= std::move(builder
);
1796 MdrunnerBuilder::MdrunnerBuilder(compat::not_null
<SimulationContext
*> context
) :
1797 impl_
{std::make_unique
<Mdrunner::BuilderImplementation
>(context
)}
1801 MdrunnerBuilder::~MdrunnerBuilder() = default;
1803 MdrunnerBuilder
&MdrunnerBuilder::addSimulationMethod(const MdrunOptions
&options
,
1804 real forceWarningThreshold
)
1806 impl_
->setExtraMdrunOptions(options
, forceWarningThreshold
);
1810 MdrunnerBuilder
&MdrunnerBuilder::addDomainDecomposition(const DomdecOptions
&options
)
1812 impl_
->addDomdec(options
);
1816 MdrunnerBuilder
&MdrunnerBuilder::addNeighborList(int nstlist
)
1818 impl_
->addVerletList(nstlist
);
1822 MdrunnerBuilder
&MdrunnerBuilder::addReplicaExchange(const ReplicaExchangeParameters
¶ms
)
1824 impl_
->addReplicaExchange(params
);
1828 MdrunnerBuilder
&MdrunnerBuilder::addMultiSim(gmx_multisim_t
* multisim
)
1830 impl_
->addMultiSim(multisim
);
1834 MdrunnerBuilder
&MdrunnerBuilder::addNonBonded(const char* nbpu_opt
)
1836 impl_
->addNonBonded(nbpu_opt
);
1840 MdrunnerBuilder
&MdrunnerBuilder::addElectrostatics(const char* pme_opt
,
1841 const char* pme_fft_opt
)
1843 // The builder method may become more general in the future, but in this version,
1844 // parameters for PME electrostatics are both required and the only parameters
1846 if (pme_opt
&& pme_fft_opt
)
1848 impl_
->addPME(pme_opt
, pme_fft_opt
);
1852 GMX_THROW(gmx::InvalidInputError("addElectrostatics() arguments must be non-null pointers."));
1857 MdrunnerBuilder
&MdrunnerBuilder::addBondedTaskAssignment(const char* bonded_opt
)
1859 impl_
->addBondedTaskAssignment(bonded_opt
);
1863 Mdrunner
MdrunnerBuilder::build()
1865 return impl_
->build();
1868 MdrunnerBuilder
&MdrunnerBuilder::addHardwareOptions(const gmx_hw_opt_t
&hardwareOptions
)
1870 impl_
->addHardwareOptions(hardwareOptions
);
1874 MdrunnerBuilder
&MdrunnerBuilder::addFilenames(ArrayRef
<const t_filenm
> filenames
)
1876 impl_
->addFilenames(filenames
);
1880 MdrunnerBuilder
&MdrunnerBuilder::addOutputEnvironment(gmx_output_env_t
* outputEnvironment
)
1882 impl_
->addOutputEnvironment(outputEnvironment
);
1886 MdrunnerBuilder
&MdrunnerBuilder::addLogFile(t_fileio
*logFileHandle
)
1888 impl_
->addLogFile(logFileHandle
);
1892 MdrunnerBuilder
&MdrunnerBuilder::addStopHandlerBuilder(std::unique_ptr
<StopHandlerBuilder
> builder
)
1894 impl_
->addStopHandlerBuilder(std::move(builder
));
1898 MdrunnerBuilder::MdrunnerBuilder(MdrunnerBuilder
&&) noexcept
= default;
1900 MdrunnerBuilder
&MdrunnerBuilder::operator=(MdrunnerBuilder
&&) noexcept
= default;