Gromacs  2026.0-dev-20250116-fa3fd9d
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Classes | Public Member Functions
gmx::DomainCommBackward Class Reference

#include <gromacs/domdec/domainpaircomm.h>

Description

Setup for selecting halo atoms to be sent and sending coordinates to another domain.

Also used for receiving forces.

This object is used for communicating with a domain that resides in backward direction along the domain decomposition grid. This object holds the NBNxM grid cell indices that needs to be sent and does the packing of coordinates into a buffer in this object. It also does the reverse indexing for the halo force reduction of forces received.

Classes

struct  ColumnInfo
 Struct for collecting information on grid columns. More...
 

Public Member Functions

 DomainCommBackward (int rank, int zone, const IVec &domainShift, PbcType pbcType, bool commOverPbc, const IVec &pbcCoordinateShift, MPI_Comm mpiCommAll)
 Constructor. More...
 
void clear ()
 Clears this communication, set no columns and atoms to send.
 
void selectHaloAtoms (const gmx_domdec_t &dd, const Grid &grid, const real cutoffSquaredNonbonded, const real cutoffSquaredBonded, const matrix box, const ivec dimensionIsTriclinic, ArrayRef< const RVec > normal, const std::vector< bool > &isCellMissingLinks)
 Determine which NBNxM grid cells (and atoms) we need to send. More...
 
FastVector< std::pair< int, int > > makeColumnsSendBuffer () const
 Creates and returns a buffer with column indices and cell counts to be sent.
 
void packCoordinateSendBuffer (const matrix box, ArrayRef< const RVec > x, ArrayRef< RVec > sendBuffer) const
 Copies the coordinates to commnicate to the send buffer.
 
void accumulateReceivedForces (ArrayRef< RVec > forces, ArrayRef< RVec > shiftForces) const
 Accumulates the received forces to the force buffer and shift force buffer, when not empty.
 
int rank () const
 Returns the rank we communicate with.
 
int zone () const
 Return the zone this domain pair resides in.
 
int domainShift (int dim) const
 Returns the shift in domains in Cartesian coordinates between the communicating domains.
 
bool shiftMultipleDomains () const
 Returns whether the domain shift is more than one domain along at least one dimension.
 
int numAtomsPerCell () const
 Returns the number of atoms per cell.
 
ArrayRef< const ColumnInfocolumnsToSend () const
 Returns the list of all columns to send with column information.
 
int numAtoms () const
 Returns the number of atoms to send.
 
ArrayRef< const int > globalAtomIndices () const
 Returns the buffer with global atom indices to communicate.
 
ArrayRef< RVecrvecBuffer ()
 Returns the buffer for sending coordinates and receiving forces.
 
ArrayRef< const RVecrvecBuffer () const
 Returns the buffer for sending coordinates and receiving forces.
 
MPI_Comm mpiCommAll () const
 The MPI communication for halo exchange.
 

Constructor & Destructor Documentation

gmx::DomainCommBackward::DomainCommBackward ( int  rank,
int  zone,
const IVec domainShift,
PbcType  pbcType,
bool  commOverPbc,
const IVec pbcCoordinateShift,
MPI_Comm  mpiCommAll 
)

Constructor.

Parameters
[in]rankThe MPI rank we communicate with
[in]zoneThe domain decomposition zone this pair of domains belongs to
[in]domainShiftThe shift >=0 in domain indices between the two domains
[in]pbcTypeThe type of PBC
[in]commOverPbcWhether we are communicating over a periodic boundary
[in]pbcCoordinateShiftThe PBC coordinate shift, 0 or 1 for each dimension
[in]mpiCommAllMPI communicator for all PP ranks

Member Function Documentation

void gmx::DomainCommBackward::selectHaloAtoms ( const gmx_domdec_t &  dd,
const Grid grid,
const real  cutoffSquaredNonbonded,
const real  cutoffSquaredBonded,
const matrix  box,
const ivec  dimensionIsTriclinic,
ArrayRef< const RVec normal,
const std::vector< bool > &  isCellMissingLinks 
)

Determine which NBNxM grid cells (and atoms) we need to send.

Parameters
[in]ddThe domain decomposition struct
[in]gridThe local NBNxM pair-search grid
[in]cutoffSquaredNonbondedThe cutoff^2 for non-bonded interactions
[in]cutoffSquaredBondedThe cutoff^2 for bonded interactions
[in]boxThe box
[in]dimensionIsTriclinicTells whether the dimensions require triclinic distance checks
[in]normalThe normal vectors to planes separating domains along the triclinic dimensions
[in]isCellMissingLinksTells whether cells are missing bonded interactions, only used when filtering bonded communication

The documentation for this class was generated from the following files: