Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Enumerations | Functions
#include "gromacs/math/vectypes.h"
#include "gromacs/utility/arrayref.h"
+ Include dependency graph for domdec_network.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

Enumerations

enum  { dddirForward, dddirBackward }
 

Functions

template<typename T >
void ddSendrecv (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, T *sendBuffer, int numElementsToSend, T *receiveBuffer, int numElementsToReceive)
 Move 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, int *buf_s, int n_s, int *buf_r, int n_r)
 Extern declaration for int specialization.
 
template void ddSendrecv< real > (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, real *buf_s, int n_s, real *buf_r, int n_r)
 Extern declaration for real specialization.
 
template void ddSendrecv< rvec > (const gmx_domdec_t *dd, int ddDimensionIndex, int direction, rvec *buf_s, int n_s, rvec *buf_r, int n_r)
 Extern declaration for rvec specialization.
 
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 DDMASTERRANK to all PP ranks.
 
void dd_bcastc (const gmx_domdec_t *dd, int nbytes, void *src, void *dest)
 Copies nbytes from src to dest on DDMASTERRANK and then broadcasts to dest on all PP ranks.
 
void dd_scatter (const gmx_domdec_t *dd, int nbytes, const void *src, void *dest)
 Scatters nbytes from src on DDMASTERRANK 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 DDMASTERRANK.
 
void dd_scatterv (const gmx_domdec_t *dd, int *scounts, int *disps, const void *sbuf, int rcount, void *rbuf)
 Scatters scounts bytes from src on DDMASTERRANK to all PP ranks, receiving rcount bytes in dest. More...
 
void dd_gatherv (const gmx_domdec_t *dd, int scount, const void *sbuf, int *rcounts, int *disps, void *rbuf)
 Gathers rcount bytes from src on all PP ranks, received in scounts bytes in dest on DDMASTERRANK. More...
 

Function Documentation

void dd_gatherv ( const gmx_domdec_t *  dd,
int  scount,
const void *  sbuf,
int *  rcounts,
int *  disps,
void *  rbuf 
)

Gathers rcount bytes from src on all PP ranks, received in scounts bytes in dest on DDMASTERRANK.

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

If scount==0, sbuf is allowed to be NULL

void dd_scatterv ( const gmx_domdec_t *  dd,
int *  scounts,
int *  disps,
const void *  sbuf,
int  rcount,
void *  rbuf 
)

Scatters scounts bytes from src on DDMASTERRANK to all PP ranks, receiving rcount bytes in dest.

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

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,
T *  sendBuffer,
int  numElementsToSend,
T *  receiveBuffer,
int  numElementsToReceive 
)

Move 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).

Todo:
This function template is deprecated, new calls should be made to the version taking ArrayRef parameters and this function template removed when unused.
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).