Gromacs  2025.0-dev-20241014-f673b97
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Enumerations | Functions
#include <cstddef>
#include <vector>
#include "gromacs/gpu_utils/devicebuffer_datatype.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for domdec.h:
+ This graph shows which files directly or indirectly include this file:

Description

This file declares functions for mdrun to call to manage the details of its domain decomposition.

Author
Berk Hess hess@.nosp@m.kth..nosp@m.se

Classes

class  gmx::ArrayRef< typename >
 STL-like interface to a C array of T (or part of a std container of T). More...
 
class  gmx::FixedCapacityVector< typename, size_t >
 Vector that behaves likes std::vector but has fixed capacity. More...
 
struct  NumPmeDomains
 Struct for passing around the number of PME domains. More...
 

Enumerations

enum  {
  ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU,
  ddCyclPME, ddCyclNr
}
 Cycle counter indices used internally in the domain decomposition.
 

Functions

int ddglatnr (const gmx_domdec_t *dd, int i)
 Returns the global topology atom number belonging to local atom index i. More...
 
void dd_store_state (const gmx_domdec_t &dd, t_state *state)
 Store the global cg indices of the home cgs in state,. More...
 
const gmx::DomdecZonesgetDomdecZones (const gmx_domdec_t &dd)
 Returns a const reference to the gmx::DomdecZones object.
 
int dd_numAtomsZones (const gmx_domdec_t &dd)
 Returns the range for atoms in zones.
 
int dd_numHomeAtoms (const gmx_domdec_t &dd)
 Returns the number of home atoms.
 
int dd_natoms_mdatoms (const gmx_domdec_t &dd)
 Returns the atom range in the local state for atoms that need to be present in mdatoms.
 
int dd_natoms_vsite (const gmx_domdec_t &dd)
 Returns the atom range in the local state for atoms involved in virtual sites.
 
void dd_get_constraint_range (const gmx_domdec_t &dd, int *at_start, int *at_end)
 Sets the atom range for atom in the local state for atoms received in constraints communication.
 
NumPmeDomains getNumPmeDomains (const gmx_domdec_t *dd)
 Returns the number of PME domains, can be called with dd=NULL.
 
std::vector< int > get_pme_ddranks (const t_commrec *cr, int pmenodeid)
 Returns the set of DD ranks that communicate with pme node cr->nodeid.
 
int dd_pme_maxshift_x (const gmx_domdec_t &dd)
 Returns the maximum shift for coordinate communication in PME, dim x.
 
int dd_pme_maxshift_y (const gmx_domdec_t &dd)
 Returns the maximum shift for coordinate communication in PME, dim y.
 
bool ddUsesUpdateGroups (const gmx_domdec_t &dd)
 Return whether update groups are used.
 
bool dd_moleculesAreAlwaysWhole (const gmx_domdec_t &dd)
 Returns whether molecules are always whole, i.e. not broken by PBC.
 
bool dd_bonded_molpbc (const gmx_domdec_t &dd, PbcType pbcType)
 Returns if we need to do pbc for calculating bonded interactions.
 
bool change_dd_cutoff (t_commrec *cr, const matrix box, gmx::ArrayRef< const gmx::RVec > x, real cutoffRequested, bool checkGpuDdLimitation)
 Change the DD non-bonded communication cut-off. More...
 
void dd_setup_dlb_resource_sharing (const t_commrec *cr, size_t uniqueDeviceId)
 Set up communication between PP ranks for averaging GPU wait times over domains. More...
 
void dd_cycles_add (const gmx_domdec_t *dd, float cycles, int ddCycl)
 Add the wallcycle count to the DD counter.
 
void dd_move_x (struct gmx_domdec_t *dd, const matrix box, gmx::ArrayRef< gmx::RVec > x, gmx_wallcycle *wcycle)
 Communicate the coordinates to the neighboring cells and do pbc.
 
void dd_move_f (struct gmx_domdec_t *dd, gmx::ForceWithShiftForces *forceWithShiftForces, gmx_wallcycle *wcycle)
 Sum the forces over the neighboring cells. More...
 
void reset_dd_statistics_counters (struct gmx_domdec_t *dd)
 Reset all the statistics and counters for total run counting.
 
void dd_move_f_vsites (const gmx_domdec_t &dd, gmx::ArrayRef< gmx::RVec > f, gmx::ArrayRef< gmx::RVec > fshift)
 Communicates the virtual site forces, reduces the shift forces when fshift != NULL.
 
void dd_clear_f_vsites (const gmx_domdec_t &dd, gmx::ArrayRef< gmx::RVec > f)
 Clears the forces for virtual sites.
 
void dd_move_x_constraints (struct gmx_domdec_t *dd, const matrix box, gmx::ArrayRef< gmx::RVec > x0, gmx::ArrayRef< gmx::RVec > x1, bool bX1IsCoord)
 Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1.
 
void dd_move_x_vsites (const gmx_domdec_t &dd, const matrix box, gmx::ArrayRef< gmx::RVec > x)
 Communicates the coordinates involved in virtual sites.
 
void dd_move_x_and_v_vsites (const gmx_domdec_t &dd, const matrix box, gmx::ArrayRef< gmx::RVec > x, gmx::ArrayRef< gmx::RVec > v)
 Communicates the positions and velocities involved in virtual sites.
 
gmx::ArrayRef< const int > dd_constraints_nlocalatoms (const gmx_domdec_t *dd)
 Returns the local atom count array for all constraints. More...
 
void dd_init_local_state (const gmx_domdec_t &dd, const t_state *state_global, t_state *local_state)
 Construct local state.
 
void constructGpuHaloExchange (const t_commrec &cr, const gmx::DeviceStreamManager &deviceStreamManager, gmx_wallcycle *wcycle)
 Construct the GPU halo exchange object(s). More...
 
void reinitGpuHaloExchange (const t_commrec &cr, DeviceBuffer< gmx::RVec > d_coordinatesBuffer, DeviceBuffer< gmx::RVec > d_forcesBuffer)
 (Re-) Initialization for GPU halo exchange More...
 
GpuEventSynchronizer * communicateGpuHaloCoordinates (const t_commrec &cr, const matrix box, GpuEventSynchronizer *dependencyEvent)
 GPU halo exchange of coordinates buffer. More...
 
void communicateGpuHaloForces (const t_commrec &cr, bool accumulateForces, gmx::FixedCapacityVector< GpuEventSynchronizer *, 2 > *dependencyEvents)
 Wait for copy of nonlocal part of coordinate array from GPU to CPU following coordinate halo exchange. More...
 
void putUpdateGroupAtomsInSamePeriodicImage (const gmx_domdec_t &dd, const gmx_mtop_t &mtop, const matrix box, gmx::ArrayRef< gmx::RVec > positions)
 Wraps the positions so that atoms from the same update group share the same periodic image wrt box. More...
 

Function Documentation

bool change_dd_cutoff ( t_commrec *  cr,
const matrix  box,
gmx::ArrayRef< const gmx::RVec x,
real  cutoffRequested,
bool  checkGpuDdLimitation 
)

Change the DD non-bonded communication cut-off.

This could fail when trying to increase the cut-off, then FALSE will be returned and the cut-off is not modified.

Parameters
[in]crCommunication recrod
[in]boxBox matrix, used for computing the dimensions of the system
[in]xPosition vector, used for computing the dimensions of the system
[in]cutoffRequestedThe requested atom to atom cut-off distance, usually the pair-list cutoff distance
[in]checkGpuDdLimitationWhether to check the GPU DD support limitation
GpuEventSynchronizer* communicateGpuHaloCoordinates ( const t_commrec &  cr,
const matrix  box,
GpuEventSynchronizer *  dependencyEvent 
)

GPU halo exchange of coordinates buffer.

Parameters
[in]crThe commrec object
[in]boxCoordinate box (from which shifts will be constructed)
[in]dependencyEventDependency event for this operation
Returns
Event recorded when this operation has been launched
void communicateGpuHaloForces ( const t_commrec &  cr,
bool  accumulateForces,
gmx::FixedCapacityVector< GpuEventSynchronizer *, 2 > *  dependencyEvents 
)

Wait for copy of nonlocal part of coordinate array from GPU to CPU following coordinate halo exchange.

Parameters
[in]crThe commrec object
[in]accumulateForcesTrue if forces should accumulate, otherwise they are set
[in]dependencyEventsDependency events for this operation
void constructGpuHaloExchange ( const t_commrec &  cr,
const gmx::DeviceStreamManager deviceStreamManager,
gmx_wallcycle *  wcycle 
)

Construct the GPU halo exchange object(s).

Parameters
[in]crThe commrec object.
[in]deviceStreamManagerManager of the GPU context and streams.
[in]wcycleThe wallclock counter.
gmx::ArrayRef<const int> dd_constraints_nlocalatoms ( const gmx_domdec_t *  dd)

Returns the local atom count array for all constraints.

The local atom count for a constraint, possible values 2/1/0, is needed to avoid not/double-counting contributions linked to the Lagrange multiplier, such as the virial and free-energy derivatives.

Note
When dd = nullptr, an empty reference is returned.
void dd_move_f ( struct gmx_domdec_t *  dd,
gmx::ForceWithShiftForces forceWithShiftForces,
gmx_wallcycle *  wcycle 
)

Sum the forces over the neighboring cells.

When fshift!=NULL the shift forces are updated to obtain the correct virial from the single sum including f.

void dd_setup_dlb_resource_sharing ( const t_commrec *  cr,
size_t  uniqueDeviceId 
)

Set up communication between PP ranks for averaging GPU wait times over domains.

When domains (PP MPI ranks) share a GPU, the individual GPU wait times are meaningless, as it depends on the order in which tasks on the same GPU finish. Therefore there wait times need to be averaged over the ranks sharing the same GPU. This function sets up the communication for that.

When there is no PP work on this rank or only one rank, no sharing is set up. It's the caller's responsibility to call this method on a PP rank only when that rank is using a GPU.

It is not necessary that hashForDevice is truly the hash of a real device UUID, merely that its value will be the same on each rank when a device is shared between ranks, and otherwise different. An index into the array of devices visible to all ranks of the same node is sufficient if all such ranks see all devices. However, that index is unreliable when ranks on the same node see different subsets of the devices on the node because of different local environments. In the case where each rank sees only one device, DLB will behave as if all ranks on the node share the same device, which will not be optimal.

void dd_store_state ( const gmx_domdec_t &  dd,
t_state state 
)

Store the global cg indices of the home cgs in state,.

This means it can be reset, even after a new DD partitioning.

int ddglatnr ( const gmx_domdec_t *  dd,
int  i 
)

Returns the global topology atom number belonging to local atom index i.

This function is intended for writing ASCII output and returns atom numbers starting at 1. When dd=NULL returns i+1.

void putUpdateGroupAtomsInSamePeriodicImage ( const gmx_domdec_t &  dd,
const gmx_mtop_t &  mtop,
const matrix  box,
gmx::ArrayRef< gmx::RVec positions 
)

Wraps the positions so that atoms from the same update group share the same periodic image wrt box.

When DD and update groups are in use, the simulation main rank should call this to ensure that e.g. when restarting a simulation that did not use update groups that the coordinates satisfy the new requirements.

This function can probably be removed when even single-rank simulations use domain decomposition, because then the choice of whether update groups are used is probably going to be the same regardless of the rank count.

Parameters
[in]ddThe DD manager
[in]mtopThe system topology
[in]boxThe global system box
[in]positionsThe global system positions
void reinitGpuHaloExchange ( const t_commrec &  cr,
DeviceBuffer< gmx::RVec d_coordinatesBuffer,
DeviceBuffer< gmx::RVec d_forcesBuffer 
)

(Re-) Initialization for GPU halo exchange

Parameters
[in]crThe commrec object
[in]d_coordinatesBufferpointer to coordinates buffer in GPU memory
[in]d_forcesBufferpointer to forces buffer in GPU memory