|
Gromacs
2026.0-dev-20251020-94eac3f
|
#include <gromacs/domdec/domainpaircomm.h>
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 cluster 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 | GridClusterRange |
| Struct for storing a cluster range for a grid column. 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 | getTargetZoneCorners (const gmx_domdec_t &dd, const matrix box, const DomainPairComm &domainPairComm) |
| Determines the corner for 2-body, corner_2b, and multi-body, corner_mb, communication distances. More... | |
| void | selectHaloAtoms (const gmx_domdec_t &dd, const Grid &grid, const real cutoffTwoBody, const real cutoffMultiBody, const ivec dimensionIsTriclinic, ArrayRef< const RVec > normal, const std::vector< bool > &isCellMissingLinks) |
| Determine which NBNxM grid clusters (and atoms) we need to send. More... | |
| FastVector< std::pair< int, int > > | makeColumnsSendBuffer () const |
| Creates and returns a buffer with column indices and cluster 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 | numAtomsPerCluster () const |
| Returns the number of atoms per cluster. | |
| ArrayRef< const GridClusterRange > | clusterRangesToSend () 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< RVec > | rvecBuffer () |
| Returns the buffer for sending coordinates and receiving forces. | |
| ArrayRef< const RVec > | rvecBuffer () const |
| Returns the buffer for sending coordinates and receiving forces. | |
| MPI_Comm | mpiCommAll () const |
| The MPI communication for halo exchange. | |
| gmx::DomainCommBackward::DomainCommBackward | ( | int | rank, |
| int | zone, | ||
| const IVec & | domainShift, | ||
| PbcType | pbcType, | ||
| bool | commOverPbc, | ||
| const IVec & | pbcCoordinateShift, | ||
| MPI_Comm | mpiCommAll | ||
| ) |
Constructor.
| [in] | rank | The MPI rank we communicate with |
| [in] | zone | The domain decomposition zone this pair of domains belongs to |
| [in] | domainShift | The shift >=0 in domain indices between the two domains |
| [in] | pbcType | The type of PBC |
| [in] | commOverPbc | Whether we are communicating over a periodic boundary |
| [in] | pbcCoordinateShift | The PBC coordinate shift, 0 or 1 for each dimension |
| [in] | mpiCommAll | MPI communicator for all PP ranks |
| void gmx::DomainCommBackward::getTargetZoneCorners | ( | const gmx_domdec_t & | dd, |
| const matrix | box, | ||
| const DomainPairComm & | domainPairComm | ||
| ) |
Determines the corner for 2-body, corner_2b, and multi-body, corner_mb, communication distances.
Communicates the computed corners between domain pairs
Note that the column bounding boxes are computed for the centers of update groups. Atoms in update groups can stick out by at most grid.dimensions().maxAtomGroupRadius.
| [in] | dd | The domain decomposition struct |
| [in] | box | The system unit cell |
| [in] | domainPairComm | Communication setup for a domain pair |
| void gmx::DomainCommBackward::selectHaloAtoms | ( | const gmx_domdec_t & | dd, |
| const Grid & | grid, | ||
| const real | cutoffTwoBody, | ||
| const real | cutoffMultiBody, | ||
| const ivec | dimensionIsTriclinic, | ||
| ArrayRef< const RVec > | normal, | ||
| const std::vector< bool > & | isCellMissingLinks | ||
| ) |
Determine which NBNxM grid clusters (and atoms) we need to send.
Should be called after calling getTargetZoneCorners().
| [in] | dd | The domain decomposition struct |
| [in] | grid | The local NBNxM pair-search grid |
| [in] | cutoffTwoBody | The cutoff for two-body interactions, includes 2 maxAtomGroupRadius |
| [in] | cutoffMultiBody | The cutoff for multi-body interactions includes 2 maxAtomGroupRadius |
| [in] | dimensionIsTriclinic | Tells whether the dimensions require triclinic distance checks |
| [in] | normal | The normal vectors to planes separating domains along the triclinic dimensions |
| [in] | isCellMissingLinks | Tells whether cells are missing bonded interactions, only used when filtering bonded communication |
1.8.5