Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Enumerations | Functions
#include "gromacs/math/vectypes.h"
+ This graph shows which files directly or indirectly include this file:

Description

This file declares functions for (mostly) the domdec module to use MPI functionality.

Todo:
Wrap the raw dd_bcast in md.cpp into a higher-level function in the domdec module, then this file can be module-internal.
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...
 

Enumerations

enum  { dddirForward, dddirBackward }
 

Functions

template<typename T >
void ddSendrecv (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, gmx::ArrayRef< T > sendBuffer, gmx::ArrayRef< T > receiveBuffer)
 Move a view of T values in the communication region one cell along the domain decomposition. More...
 
template void ddSendrecv< int > (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, gmx::ArrayRef< int > sendBuffer, gmx::ArrayRef< int > receiveBuffer)
 Extern declaration for int specialization.
 
template void ddSendrecv< real > (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, gmx::ArrayRef< real > sendBuffer, gmx::ArrayRef< real > receiveBuffer)
 Extern declaration for real specialization.
 
template void ddSendrecv< gmx::RVec > (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, gmx::ArrayRef< gmx::RVec > sendBuffer, gmx::ArrayRef< gmx::RVec > receiveBuffer)
 Extern declaration for gmx::RVec specialization.
 
void dd_sendrecv2_rvec (const struct gmx_domdec_t *dd, int ddimind, rvec *buf_s_fw, int n_s_fw, rvec *buf_r_fw, int n_r_fw, rvec *buf_s_bw, int n_s_bw, rvec *buf_r_bw, int n_r_bw)
 Move revc's in the comm. region one cell along the domain decomposition. More...
 
void dd_bcast (const gmx_domdec_t *dd, int nbytes, void *data)
 Broadcasts nbytes from data on DDMAINRANK to all PP ranks.
 
void dd_scatter (const gmx_domdec_t *dd, int nbytes, const void *src, void *dest)
 Scatters nbytes from src on DDMAINRANK to all PP ranks, received in dest.
 
void dd_gather (const gmx_domdec_t *dd, int nbytes, const void *src, void *dest)
 Gathers nbytes from src on all PP ranks, received in dest on DDMAINRANK.
 
template<typename T >
void dd_scatterv (const gmx_domdec_t *dd, gmx::ArrayRef< const int > scounts, gmx::ArrayRef< const int > disps, const T *sbuf, int rcount, T *rbuf)
 Scatters scounts elements of type T from src on DDMAINRANK to all PP ranks, receiving rcount * sc_scattervSize elements in dest. More...
 
template void dd_scatterv (const gmx_domdec_t *dd, gmx::ArrayRef< const int > scounts, gmx::ArrayRef< const int > disps, const int *sbuf, int rcount, int *rbuf)
 Instantiation of dd_scatterv for type int.
 
template void dd_scatterv (const gmx_domdec_t *dd, gmx::ArrayRef< const int > scounts, gmx::ArrayRef< const int > disps, const gmx::RVec *sbuf, int rcount, gmx::RVec *rbuf)
 Instantiation of dd_scatterv for type real.
 
template<typename T >
void dd_gatherv (const gmx_domdec_t *dd, int scount, const T *sbuf, gmx::ArrayRef< const int > rcounts, gmx::ArrayRef< const int > disps, T *rbuf)
 Gathers rcount elements of type T from src on all PP ranks, received in scounts elements in dest on DDMAINRANK. More...
 
template void dd_gatherv (const gmx_domdec_t *dd, int scount, const int *sbuf, gmx::ArrayRef< const int > rcounts, gmx::ArrayRef< const int > disps, int *rbuf)
 Instantiation of dd_gatherv for type int.
 
template void dd_gatherv (const gmx_domdec_t *dd, int scount, const gmx::RVec *sbuf, gmx::ArrayRef< const int > rcounts, gmx::ArrayRef< const int > disps, gmx::RVec *rbuf)
 Instantiation of dd_gatherv for type real.
 

Function Documentation

template<typename T >
void dd_gatherv ( const gmx_domdec_t *  dd,
int  scount,
const T *  sbuf,
gmx::ArrayRef< const int >  rcounts,
gmx::ArrayRef< const int >  disps,
T *  rbuf 
)

Gathers rcount elements of type T from src on all PP ranks, received in scounts elements in dest on DDMAINRANK.

Template Parameters
TData type, can only be int or real

See man MPI_Gatherv for details of how to construct scounts and disps.

If scount==0, sbuf is allowed to be nullptr.

template<typename T >
void dd_scatterv ( const gmx_domdec_t *  dd,
gmx::ArrayRef< const int >  scounts,
gmx::ArrayRef< const int >  disps,
const T *  sbuf,
int  rcount,
T *  rbuf 
)

Scatters scounts elements of type T from src on DDMAINRANK to all PP ranks, receiving rcount * sc_scattervSize elements in dest.

Template Parameters
TData type, can only be int or real

See man MPI_Scatterv for details of how to construct scounts and disps. If rcount==0, rbuf is allowed to be nullptr.

void dd_sendrecv2_rvec ( const struct gmx_domdec_t *  dd,
int  ddimind,
rvec buf_s_fw,
int  n_s_fw,
rvec buf_r_fw,
int  n_r_fw,
rvec buf_s_bw,
int  n_s_bw,
rvec buf_r_bw,
int  n_r_bw 
)

Move revc's in the comm. region one cell along the domain decomposition.

Moves in dimension indexed by ddimind, simultaneously in the forward and backward directions.

template<typename T >
void ddSendrecv ( const gmx_domdec_t *  dd,
int  ddDimensionIndex,
int  direction,
gmx::ArrayRef< T >  sendBuffer,
gmx::ArrayRef< T >  receiveBuffer 
)

Move a view of T values in the communication region one cell along the domain decomposition.

Moves in the dimension indexed by ddDimensionIndex, either forward (direction=dddirFoward) or backward (direction=dddirBackward).