MADNESS
0.10.1
|
This class wraps/extends the MPI interface for World
.
More...
#include <worldmpi.h>
Public Member Functions | |
WorldMpiInterface (const SafeMPI::Intracomm &comm) | |
Constructs an interface in the specified SafeMPI communicator. More... | |
~WorldMpiInterface ()=default | |
template<typename T > | |
std::enable_if<!std::is_pointer< T >::value, void >::type | Bcast (T &buffer, int root) const |
MPI broadcast a datum. More... | |
template<typename T > | |
void | Bcast (T *buffer, int count, int root) const |
MPI broadcast an array of count elements. More... | |
SafeMPI::Intracomm & | comm () |
Returns the associated SafeMPI communicator. More... | |
template<typename T > | |
std::enable_if<!std::is_pointer< T >::value, SafeMPI::Request >::type | Irecv (T &buf, int source, int tag=SafeMPI::DEFAULT_SEND_RECV_TAG) const |
Async receive datum from process source with default tag=1 . More... | |
template<typename T > | |
SafeMPI::Request | Irecv (T *buf, int count, int source, int tag=SafeMPI::DEFAULT_SEND_RECV_TAG) const |
Async receive data of up to count elements from process source . More... | |
template<typename T > | |
std::enable_if<!std::is_pointer< T >::value, SafeMPI::Request >::type | Isend (const T &datum, int dest, int tag=SafeMPI::DEFAULT_SEND_RECV_TAG) const |
Isend one element. More... | |
int | nproc () const |
Access the total number of processes. More... | |
int | rank () const |
Access the rank of this process. More... | |
template<typename T > | |
std::enable_if<!std::is_pointer< T >::value, void >::type | Recv (T &buf, int src, int tag=SafeMPI::DEFAULT_SEND_RECV_TAG) const |
Receive datum from process src . More... | |
template<typename T > | |
void | Recv (T *buf, long lenbuf, int src, int tag) const |
Receive data of up to lenbuf elements from process src . More... | |
template<typename T > | |
void | Recv (T *buf, long lenbuf, int src, int tag, SafeMPI::Status &status) const |
Receive data of up to lenbuf elements from process dest with status . More... | |
template<typename T > | |
std::enable_if<!std::is_pointer< T >::value, void >::type | Send (const T &datum, int dest, int tag=SafeMPI::DEFAULT_SEND_RECV_TAG) const |
Send element to process dest with default tag=1001 . More... | |
template<class T > | |
void | Send (const T *buf, long lenbuf, int dest, int tag=SafeMPI::DEFAULT_SEND_RECV_TAG) const |
Send array of lenbuf elements to process dest . More... | |
int | size () const |
Access the total number of processes. More... | |
Public Member Functions inherited from SafeMPI::Intracomm | |
Intracomm (const Intracomm &other) | |
Intracomm (const MPI_Comm &comm, bool take_ownership_of_comm=true) | |
Intracomm (const WorldInitObject &) | |
~Intracomm () | |
void | Abort (int code=1) const |
void | Allreduce (const void *sendbuf, void *recvbuf, const int count, const MPI_Datatype datatype, const MPI_Op op) const |
void | Barrier () const |
void | Bcast (void *buf, size_t count, const MPI_Datatype datatype, const int root) const |
void | binary_tree_info (int root, int &parent, int &child0, int &child1) |
Construct info about a binary tree with given root. More... | |
Intracomm | Clone () const |
Intracomm | Create (Group group) const |
bool | Get_attr (int key, void *value) const |
Group | Get_group () const |
MPI_Comm & | Get_mpi_comm () const |
int | Get_rank () const |
int | Get_size () const |
Request | Irecv (void *buf, const int count, const MPI_Datatype datatype, const int src, const int tag) const |
Request | Isend (const void *buf, const int count, const MPI_Datatype datatype, const int dest, const int tag) const |
Request | Issend (const void *buf, const int count, const MPI_Datatype datatype, const int dest, const int tag) const |
bool | operator== (const Intracomm &other) const |
void | Recv (void *buf, const int count, const MPI_Datatype datatype, const int source, const int tag) const |
void | Recv (void *buf, const int count, const MPI_Datatype datatype, const int source, const int tag, MPI_Status &status) const |
void | Reduce (const void *sendbuf, void *recvbuf, const int count, const MPI_Datatype datatype, const MPI_Op op, const int root) const |
void | Send (const void *buf, const int count, const MPI_Datatype datatype, int dest, int tag) const |
Intracomm | Split (int Color, int Key=0) const |
Intracomm | Split_type (int Type, int Key=0) const |
int | unique_reserved_tag () |
Returns a unique tag reserved for long-term use (0<tag<1000) More... | |
int | unique_tag () |
Returns a unique tag for temporary use (1023<tag<4095) More... | |
Private Member Functions | |
WorldMpiInterface (const WorldMpiInterface &)=delete | |
WorldMpiInterface & | operator= (const WorldMpiInterface &)=delete |
Private Member Functions inherited from madness::detail::WorldMpiRuntime | |
WorldMpiRuntime () | |
Constructor. More... | |
~WorldMpiRuntime () | |
Destructor. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from SafeMPI::Intracomm | |
static int | unique_tag_period () |
Static Public Attributes inherited from SafeMPI::Intracomm | |
static const int | SHARED_SPLIT_TYPE = MPI_COMM_TYPE_SHARED |
static const int | UNDEFINED_COLOR = MPI_UNDEFINED |
static const int | UNDEFINED_SPLIT_TYPE = MPI_UNDEFINED |
This class wraps/extends the MPI interface for World
.
|
privatedelete |
|
inline |
|
default |
|
inline |
MPI broadcast a datum.
T | The type of data to broadcast. |
[in,out] | buffer | The datum to send (if this is the root process); otherwise, the received datum. |
[in] | root | The process that is sending the data to other processes. |
References SafeMPI::Intracomm::Bcast(), MPI_BYTE, and T().
|
inline |
MPI broadcast an array of count
elements.
T | The type of data to broadcast. |
[in,out] | buffer | The data to send (if this is the root process); otherwise, a buffer to receive the data. |
[in] | count | The number of data elements being broadcast. |
[in] | root | The process that is sending the data to other processes. |
References SafeMPI::Intracomm::Bcast(), MPI_BYTE, and T().
Referenced by madness::World::World(), and madness::Cloud::replicate().
|
inline |
Returns the associated SafeMPI
communicator.
SafeMPI
communicator. Referenced by madness::WorldAmInterface::WorldAmInterface(), madness::MacroTaskQ::create_worlds(), madness::World::find_instance(), madness::archive::ArchiveStoreImpl< ParallelOutputArchive< VectorOutputArchive >, WorldContainer< keyT, valueT > >::store(), and test_multi_world().
|
inline |
Async receive datum from process source
with default tag=1
.
T | The type of data to receive. |
[out] | buf | Where to put the received datum. |
[in] | source | The source process. |
[in] | tag | The MPI tag. |
References SafeMPI::Intracomm::Irecv(), MPI_BYTE, source(), and T().
|
inline |
Async receive data of up to count
elements from process source
.
T | The type of data to receive. |
[out] | buf | Where to put the received data. |
[in] | count | The number of data elements to receive. |
[in] | source | The source process. |
[in] | tag | The MPI tag. |
References SafeMPI::Intracomm::Irecv(), MPI_BYTE, source(), and T().
Referenced by madness::WorldGopInterface::broadcast(), madness::WorldGopInterface::concat0(), madness::SystolicMatrixAlgorithm< T >::cycle(), madness::WorldGopInterface::fence_impl(), and madness::WorldGopInterface::reduce().
|
inline |
Isend one element.
T | The type of data to send. |
[in] | datum | The element to send. |
[in] | dest | The destination process. |
[in] | tag | The MPI tag. |
References SafeMPI::Intracomm::Isend(), MPI_BYTE, and T().
Referenced by madness::WorldGopInterface::broadcast(), madness::WorldGopInterface::concat0(), madness::WorldGopInterface::fence_impl(), and madness::WorldGopInterface::reduce().
|
inline |
Access the total number of processes.
References SafeMPI::Intracomm::Get_size().
Referenced by madness::World::nproc().
|
privatedelete |
|
inline |
Access the rank of this process.
References SafeMPI::Intracomm::Get_rank().
Referenced by madness::World::rank(), and madness::redirectio().
|
inline |
Receive datum from process src
.
T | The type of data to receive. |
[out] | buf | The received datum. |
[in] | src | The source process. |
[in] | tag | The MPI tag. |
References MPI_BYTE, SafeMPI::Intracomm::Recv(), and T().
|
inline |
Receive data of up to lenbuf
elements from process src
.
T | The type of data to receive. |
[out] | buf | Where to put the received data. |
[in] | lenbuf | The maximum number of data elements to receive. |
[in] | src | The source process. |
[in] | tag | The MPI tag. |
References MPI_BYTE, SafeMPI::Intracomm::Recv(), and T().
Referenced by madness::WorldGopInterface::concat0(), madness::SystolicMatrixAlgorithm< T >::cycle(), madness::archive::MPIRawInputArchive::load(), madness::archive::MPIInputArchive::load(), and madness::archive::ArchiveStoreImpl< ParallelOutputArchive< localarchiveT >, WorldContainer< keyT, valueT > >::store().
|
inline |
Receive data of up to lenbuf
elements from process dest
with status
.
T | The type of data to receive. |
[out] | buf | Where to put the received data. |
[in] | lenbuf | The maximum number of data elements to receive. |
[in] | src | The source process. |
[in] | tag | The MPI tag. |
[in] | status | The status. |
References MPI_BYTE, SafeMPI::Intracomm::Recv(), status, and T().
|
inline |
Send element to process dest
with default tag=1001
.
T | The type of data to send. |
[in] | datum | The data element to send. |
[in] | dest | The destination process. |
[in] | tag | The MPI tag. |
References MPI_BYTE, SafeMPI::Intracomm::Send(), and T().
|
inline |
Send array of lenbuf
elements to process dest
.
T | The type of data to send. |
[in] | buf | Pointer to the data. |
[in] | lenbuf | The number of data elements to send. |
[in] | dest | The destination process. |
[in] | tag | The MPI tag. |
References MPI_BYTE, SafeMPI::Intracomm::Send(), and T().
Referenced by madness::WorldGopInterface::concat0(), madness::SystolicMatrixAlgorithm< T >::cycle(), madness::archive::MPIOutputArchive::flush(), madness::archive::ArchiveStoreImpl< ParallelOutputArchive< localarchiveT >, WorldContainer< keyT, valueT > >::store(), and madness::archive::MPIRawOutputArchive::store().
|
inline |
Access the total number of processes.
References SafeMPI::Intracomm::Get_size().
Referenced by madness::World::size().