Teuchos - Trilinos Tools Package
Version of the Day
|
Abstract interface for distributed-memory communication. More...
#include <Teuchos_Comm.hpp>
Public Member Functions | |
virtual int | getTag () const =0 |
The current tag. More... | |
Destructor | |
virtual | ~Comm () |
Destructor, declared virtual for safety of derived classes. More... | |
Query functions | |
virtual int | getRank () const =0 |
Returns the rank of this process. More... | |
virtual int | getSize () const =0 |
Returns the number of processes that make up this communicator. More... | |
Collective Operations | |
virtual void | barrier () const =0 |
Pause every process in *this communicator until all the processes reach this point. More... | |
virtual void | broadcast (const int rootRank, const Ordinal bytes, char buffer[]) const =0 |
Broadcast values from the root process to the slave processes. More... | |
virtual void | gather (const Ordinal sendBytes, const char sendBuffer[], const Ordinal recvBytes, char recvBuffer[], const int root) const =0 |
Gather values from all processes to the root process. More... | |
virtual void | gatherAll (const Ordinal sendBytes, const char sendBuffer[], const Ordinal recvBytes, char recvBuffer[]) const =0 |
Gather values from each process to collect on all processes. More... | |
virtual void | reduceAll (const ValueTypeReductionOp< Ordinal, char > &reductOp, const Ordinal bytes, const char sendBuffer[], char globalReducts[]) const =0 |
Global reduction. More... | |
virtual void | scan (const ValueTypeReductionOp< Ordinal, char > &reductOp, const Ordinal bytes, const char sendBuffer[], char scanReducts[]) const =0 |
Scan reduction. More... | |
Blocking Point-to-Point Operations | |
virtual void | send (const Ordinal bytes, const char sendBuffer[], const int destRank) const =0 |
Possibly blocking send of data from this process to another process. More... | |
virtual void | send (const Ordinal bytes, const char sendBuffer[], const int destRank, const int tag) const =0 |
Variant of send() that takes a tag. More... | |
virtual void | ssend (const Ordinal bytes, const char sendBuffer[], const int destRank) const =0 |
Always blocking send of data from this process to another process. More... | |
virtual void | ssend (const Ordinal bytes, const char sendBuffer[], const int destRank, const int tag) const =0 |
Variant of ssend() that takes a message tag. More... | |
virtual int | receive (const int sourceRank, const Ordinal bytes, char recvBuffer[]) const =0 |
Blocking receive of data from this process to another process. More... | |
virtual void | readySend (const ArrayView< const char > &sendBuffer, const int destRank) const =0 |
Ready send of data from this process to another process. More... | |
virtual void | readySend (const Ordinal bytes, const char sendBuffer[], const int destRank, const int tag) const =0 |
Variant of readySend() that accepts a message tag. More... | |
Non-blocking Point-to-Point Operations | |
virtual RCP< CommRequest < Ordinal > > | isend (const ArrayView< const char > &sendBuffer, const int destRank) const =0 |
Non-blocking send. More... | |
virtual RCP< CommRequest < Ordinal > > | isend (const ArrayView< const char > &sendBuffer, const int destRank, const int tag) const =0 |
Variant of isend() that takes a tag. More... | |
virtual RCP< CommRequest < Ordinal > > | ireceive (const ArrayView< char > &recvBuffer, const int sourceRank) const =0 |
Non-blocking receive. More... | |
virtual RCP< CommRequest < Ordinal > > | ireceive (const ArrayView< char > &recvBuffer, const int sourceRank, const int tag) const =0 |
Variant of ireceive that takes a tag. More... | |
virtual void | waitAll (const ArrayView< RCP< CommRequest< Ordinal > > > &requests) const =0 |
Wait on a set of communication requests. More... | |
virtual void | waitAll (const ArrayView< RCP< CommRequest< Ordinal > > > &requests, const ArrayView< RCP< CommStatus< Ordinal > > > &statuses) const =0 |
Wait on communication requests, and return their statuses. More... | |
virtual RCP< CommStatus < Ordinal > > | wait (const Ptr< RCP< CommRequest< Ordinal > > > &request) const =0 |
Wait on a single communication request, and return its status. More... | |
Subcommunicator Operations | |
virtual RCP< Comm > | duplicate () const =0 |
Duplicate this communicator. More... | |
virtual RCP< Comm > | split (const int color, const int key) const =0 |
Split a communicator into subcommunicators based on color and key. More... | |
virtual RCP< Comm > | createSubcommunicator (const ArrayView< const int > &ranks) const =0 |
Create a subcommunicator containing the specified processes. More... | |
Public Member Functions inherited from Teuchos::Describable | |
virtual std::string | description () const |
Return a simple one-line description of this object. More... | |
virtual void | describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
Print the object with some verbosity level to a FancyOStream. More... | |
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
Version of describe() that takes an std::ostream instead of a FancyOStream. More... | |
virtual | ~Describable () |
Destructor (marked virtual for memory safety of derived classes). More... | |
Public Member Functions inherited from Teuchos::LabeledObject | |
LabeledObject () | |
Construct with an empty label. More... | |
virtual | ~LabeledObject () |
virtual void | setObjectLabel (const std::string &objectLabel) |
Set the object label (see LabeledObject). More... | |
virtual std::string | getObjectLabel () const |
Get the object label (see LabeledObject). More... | |
Related Functions | |
(Note that these are not member functions.) | |
enum | EReductionType { REDUCE_SUM, REDUCE_MIN, REDUCE_MAX, REDUCE_AND, REDUCE_BOR } |
Predefined reduction operations that Teuchos::Comm understands. More... | |
template<typename Ordinal > | |
int | rank (const Comm< Ordinal > &comm) |
Get the process rank. More... | |
template<typename Ordinal > | |
int | size (const Comm< Ordinal > &comm) |
Get the number of processes in the communicator. More... | |
template<typename Ordinal > | |
void | barrier (const Comm< Ordinal > &comm) |
Barrier. More... | |
template<typename Ordinal , typename Packet > | |
void | broadcast (const Comm< Ordinal > &comm, const int rootRank, const Ordinal count, Packet buffer[]) |
Broadcast array of objects that use value semantics. More... | |
template<typename Ordinal , typename Packet > | |
void | broadcast (const Comm< Ordinal > &comm, const int rootRank, const ArrayView< Packet > &buffer) |
Broadcast array of objects that use value semantics. More... | |
template<typename Ordinal , typename Packet > | |
void | broadcast (const Comm< Ordinal > &comm, const int rootRank, Packet *object) |
Broadcast single object that use value semantics. More... | |
template<typename Ordinal , typename Packet > | |
void | broadcast (const Comm< Ordinal > &comm, const int rootRank, const Ptr< Packet > &object) |
Broadcast single object that use value semantics. More... | |
template<typename Ordinal , typename Packet > | |
void | broadcast (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const int rootRank, const Ordinal count, Packet *const buffer[]) |
Broadcast array of objects that use reference semantics. More... | |
template<typename Ordinal , typename Packet > | |
void | broadcast (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const int rootRank, const ArrayView< const Ptr< Packet > > &buffer) |
Broadcast array of objects that use reference semantics. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | broadcast (const Comm< Ordinal > &comm, const Serializer &serializer, const int rootRank, const Ordinal count, Packet buffer[]) |
Broadcast array of objects that use value semantics using customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
void | gather (const Packet sendBuf[], const Ordinal sendCount, Packet recvBuf[], const Ordinal recvCount, const int root, const Comm< Ordinal > &comm) |
Gather values from each process to the root process. More... | |
template<typename Ordinal , typename Packet > | |
void | gatherv (const Packet sendBuf[], const Ordinal sendCount, Packet recvBuf[], const Ordinal recvCounts[], const Ordinal displs[], const int root, const Comm< Ordinal > &comm) |
Gather arrays of possibly different lengths from each process to the root process. More... | |
template<typename Ordinal , typename Packet > | |
void | gatherAll (const Comm< Ordinal > &comm, const Ordinal sendCount, const Packet sendBuffer[], const Ordinal recvCount, Packet recvBuffer[]) |
Gather array of objects that use value semantics from every process to every process. More... | |
template<typename Ordinal , typename Packet > | |
void | gatherAll (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const Ordinal sendCount, const Packet *const sendBuffer[], const Ordinal recvCount, Packet *const recvBuffer[]) |
Gather array of objects that use reference semantics from every process to every process. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | gatherAll (const Comm< Ordinal > &comm, const Serializer &serializer, const Ordinal sendCount, const Packet sendBuffer[], const Ordinal recvCount, Packet recvBuffer[]) |
Gather array of objects that use value semantics from every process to every process using customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
void | scatter (const Packet sendBuf[], const Ordinal sendCount, Packet recvBuf[], const Ordinal recvCount, const Ordinal root, const Comm< Ordinal > &comm) |
Wrapper for MPI_Scatter; scatter collective. More... | |
template<typename Ordinal , typename Packet > | |
void | reduce (const Packet sendBuf[], Packet recvBuf[], const Ordinal count, const EReductionType reductType, const Ordinal root, const Comm< Ordinal > &comm) |
Wrapper for MPI_Reduce; reduction to one process, using a built-in reduction operator selected by enum. More... | |
template<typename Ordinal , typename Packet > | |
void | reduceAll (const Comm< Ordinal > &comm, const ValueTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet sendBuffer[], Packet globalReducts[]) |
Wrapper for MPI_Allreduce that takes a custom reduction operator. More... | |
template<typename Ordinal , typename Packet > | |
void | reduceAll (const Comm< Ordinal > &comm, const EReductionType reductType, const Ordinal count, const Packet sendBuffer[], Packet globalReducts[]) |
Collective reduce all of array of objects using value semantics using a pre-defined reduction type. More... | |
template<typename Ordinal , typename Packet > | |
void | reduceAll (const Comm< Ordinal > &comm, const EReductionType reductType, const Packet &send, const Ptr< Packet > &globalReduct) |
Collective reduce all for single object using value semantics using a pre-defined reduction type. More... | |
template<typename Ordinal , typename Packet > | |
void | reduceAll (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const ReferenceTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet *const sendBuffer[], Packet *const globalReducts[]) |
Collective reduce all for array of objects using reference semantics. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | reduceAll (const Comm< Ordinal > &comm, const Serializer &serializer, const ValueTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet sendBuffer[], Packet globalReducts[]) |
Collective reduce all of array of objects using value semantics using a user-defined reduction operator and customized serializer. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | reduceAll (const Comm< Ordinal > &comm, const Serializer &serializer, const EReductionType reductType, const Ordinal count, const Packet sendBuffer[], Packet globalReducts[]) |
Collective reduce all of array of objects using value semantics using a pre-defined reduction type and customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
void | scan (const Comm< Ordinal > &comm, const ValueTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet sendBuffer[], Packet scanReducts[]) |
Scan/Reduce array of objects that use value semantics using a user-defined reduction operator. More... | |
template<typename Ordinal , typename Packet > | |
void | scan (const Comm< Ordinal > &comm, const EReductionType reductType, const Ordinal count, const Packet sendBuffer[], Packet scanReducts[]) |
Scan/Reduce array of objects using value semantics using a predefined reduction type. More... | |
template<typename Ordinal , typename Packet > | |
void | scan (const Comm< Ordinal > &comm, const EReductionType reductType, const Packet &send, const Ptr< Packet > &scanReduct) |
Scan/Reduce single object using value semantics using a predefined reduction type. More... | |
template<typename Ordinal , typename Packet > | |
void | scan (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const ReferenceTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet *const sendBuffer[], Packet *const scanReducts[]) |
Scan/Reduce array of objects that use reference semantics using a user-defined reduction operator. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | scan (const Comm< Ordinal > &comm, const Serializer &serializer, const ValueTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet sendBuffer[], Packet scanReducts[]) |
Scan/Reduce array of objects that use value semantics using a user-defined reduction operator and customized serializer. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | scan (const Comm< Ordinal > &comm, const Serializer &serializer, const EReductionType reductType, const Ordinal count, const Packet sendBuffer[], Packet scanReducts[]) |
Scan/Reduce array of objects using value semantics using a predefined reduction type and customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
void | send (const Comm< Ordinal > &comm, const Ordinal count, const Packet sendBuffer[], const int destRank) |
Send objects that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet > | |
void | ssend (const Comm< Ordinal > &comm, const Ordinal count, const Packet sendBuffer[], const int destRank) |
Synchronously send objects that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet > | |
void | send (const Comm< Ordinal > &comm, const Packet &send, const int destRank) |
Send a single object that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet > | |
void | ssend (const Comm< Ordinal > &comm, const Packet &send, const int destRank) |
Synchronously send a single object that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet > | |
void | send (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const Ordinal count, const Packet *const sendBuffer[], const int destRank) |
Send objects that use reference semantics to another process. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | send (const Comm< Ordinal > &comm, const Serializer &serializer, const Ordinal count, const Packet sendBuffer[], const int destRank) |
Send objects that use values semantics to another process using customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
int | receive (const Comm< Ordinal > &comm, const int sourceRank, const Ordinal count, Packet recvBuffer[]) |
Receive objects that use values semantics from another process. More... | |
template<typename Ordinal , typename Packet > | |
int | receive (const Comm< Ordinal > &comm, const int sourceRank, Packet *recv) |
Receive a single object that use values semantics from another process. More... | |
template<typename Ordinal , typename Packet > | |
int | receive (const Comm< Ordinal > &comm, const Serializer< Ordinal, Packet > &serializer, const int sourceRank, const Ordinal count, Packet *const recvBuffer[]) |
Receive objects that use reference semantics from another process. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
int | receive (const Comm< Ordinal > &comm, const Serializer &serializer, const int sourceRank, const Ordinal count, Packet recvBuffer[]) |
Receive objects that use values semantics from another process using customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
void | readySend (const Comm< Ordinal > &comm, const ArrayView< const Packet > &sendBuffer, const int destRank) |
Ready-Send an array of objects that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet > | |
void | readySend (const Comm< Ordinal > &comm, const Packet &send, const int destRank) |
Ready-Send a single object that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
void | readySend (const Comm< Ordinal > &comm, const Serializer &serializer, const ArrayView< const Packet > &sendBuffer, const int destRank) |
Ready-Send an array of objects that use values semantics to another process using customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
RCP< CommRequest< Ordinal > > | isend (const Comm< Ordinal > &comm, const ArrayRCP< const Packet > &sendBuffer, const int destRank) |
Send objects that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet > | |
RCP< CommRequest< Ordinal > > | isend (const Comm< Ordinal > &comm, const RCP< const Packet > &send, const int destRank) |
Send a single object that use values semantics to another process. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
RCP< CommRequest< Ordinal > > | isend (const Comm< Ordinal > &comm, const Serializer &serializer, const ArrayRCP< const Packet > &sendBuffer, const int destRank) |
Send objects that use values semantics to another process using customized serializer. More... | |
template<typename Ordinal , typename Packet > | |
RCP< CommRequest< Ordinal > > | ireceive (const Comm< Ordinal > &comm, const ArrayRCP< Packet > &recvBuffer, const int sourceRank) |
Receive one or more objects (that use values semantics) from another process. More... | |
template<typename Ordinal , typename Packet > | |
RCP< CommRequest< Ordinal > > | ireceive (const Comm< Ordinal > &comm, const RCP< Packet > &recv, const int sourceRank) |
Receive one object (that uses values semantics) from another process. More... | |
template<typename Ordinal , typename Packet , typename Serializer > | |
RCP< CommRequest< Ordinal > > | ireceive (const Comm< Ordinal > &comm, const Serializer &serializer, const ArrayRCP< Packet > &recvBuffer, const int sourceRank) |
Send objects that use values semantics to another process using customized serializer. More... | |
template<typename Ordinal > | |
void | waitAll (const Comm< Ordinal > &comm, const ArrayView< RCP< CommRequest< Ordinal > > > &requests) |
Wait for an array of Teuchos::CommRequest objects. More... | |
template<typename Ordinal > | |
void | waitAll (const Comm< Ordinal > &comm, const ArrayView< RCP< CommRequest< Ordinal > > > &requests, const ArrayView< RCP< CommStatus< Ordinal > > > &statuses) |
Wait on one or more communication requests, and return their statuses. More... | |
template<typename Ordinal > | |
RCP< CommStatus< Ordinal > > | wait (const Comm< Ordinal > &comm, const Ptr< RCP< CommRequest< Ordinal > > > &request) |
Wait on a single communication request, and return its status. More... | |
Related Functions inherited from Teuchos::Describable | |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
Describable output stream manipulator. More... | |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
Output stream operator for Describable manipulator. More... | |
Additional Inherited Members | |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default = VERB_DEFAULT |
Default value for the verbLevel argument of describe(). More... | |
Abstract interface for distributed-memory communication.
Ordinal | Type of indices used for communication. |
This class is Teuchos' interface to distributed-memory communication between one or more parallel processes. It presents an interface very much like that of MPI (the Message Passing Interface). Teuchos provides two implementations of Comm:
Comm is an abstract interface. You cannot create a Comm directly. You have to create one of the subclasses. The normal way to handle a Comm is to pass it around using RCP (a reference-counted "smart" pointer). For example:
Comm's communication methods that actually send or receive data accept that data as an array of char
. You should never call these methods directly. Instead, you should use the nonmember "helper" functions in Teuchos_CommHelpers.hpp. These methods are templated on the Packet
type, that is, the type of data you want to send or receive. See the example below.
You should consider an RCP<const Comm<int> >
as equivalent to the MPI_Comm opaque handle, except that the RCP also does reference counting to ensure memory safety when using the same communicator in different parts of the code. That is, copying the RCP does not create a new communicator; the following two codes do about the same thing, except with a different syntax (and reference counting in the second code).
Raw MPI_Comm handles:
Reference-counted pointers to Comm:
If you want to make a "new communicator" rather than just "copy the handle," you should call the duplicate() method. This has the same behavior as MPI_Comm_dup (which see).
The "reference counting" feature means that the subclass of Comm will take care of freeing the underlying MPI_Comm (and any other data structures it may use) by calling MPI_Comm_free if necessary, once the reference count of the RCP goes to zero.
Comm works whether or not you have build Trilinos with MPI support. If you want to make a "default" Comm that is the equivalent of MPI_COMM_WORLD, but you don't know if your Trilinos with MPI enabled, you may use GlobalMPISession to call MPI_Init if necessary, and DefaultComm to "get a default communicator." For example: int main (int argc, char* argv[]) { using Teuchos::Comm; using Teuchos::DefaultComm; using Teuchos::RCP;
// This replaces the call to MPI_Init. If you didn't // build with MPI, this doesn't call MPI functions. Teuchos::GlobalMPISesssion session (&argc, &argv, NULL); // comm is the equivalent of MPI_COMM_WORLD. RCP<const Comm<int> > comm = DefaultComm<int>::getComm ();
// ... use comm in your code as you would use MPI_COMM_WORLD ...
// We don't need to call MPI_Finalize, since the // destructor of GlobalMPISession does that for us. return EXIT_SUCCESS; } This code works whether or not you built Trilinos with MPI support. It is not necessary to use GlobalMPISession, but it's useful so you don't have to remember to call MPI_Finalize. If you don't want to use GlobalMPISession, you can still call DefaultComm<int>::getComm()
, though you must have called MPI_Init first if you build Trilinos with MPI support. Furthermore, if you know MPI is present, you don't need to use DefaultComm. You may simply pass MPI_COMM_WORLD directly to MpiComm, like this:
You may also pass an arbitrary MPI_Comm directly into MpiComm's constructor, though you are responsible for freeing it after use (via MPI_Comm_free) if necessary. You may automate the freeing of your MPI_Comm by using OpaqueWrapper (which see).
As we mentioned above, for communication of data with Comm, you you should use the nonmember "helper" functions in Teuchos_CommHelpers.hpp. These methods are templated on the Packet
type, that is, the type of data you want to send or receive. For example, suppose you have two processes (with ranks 0 and 1, respectively), and you want to send an array of 10 double
from Process 0 to Process 1. Both processes have defined RCP<const Comm<int> > comm
as above. Here is the code on Process 0:
Here is the code on Process 1:
Please refer to the documentation in Teuchos_CommHelpers.hpp for more details.
This interface is templated on the ordinal type but only deals with buffers of untyped data represented as arrays char
type. All reduction operations that are initiated by the concreate communicator object are performed by user-defined ReductOpBase
objects. It is the responsibility of the ReductOpBase
object to know what the currect data type is, to perform casts or serializations/unserializations to and from char[]
buffers, and to know how to reduce the objects correctly. It is strictly up to the client to correctly convert data types to char[]
arrays but there is a great deal of helper code to make this easy and safe.
Definition at line 85 of file Teuchos_Comm.hpp.
|
inlinevirtual |
Destructor, declared virtual for safety of derived classes.
Definition at line 323 of file Teuchos_Comm.hpp.
|
pure virtual |
The current tag.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Returns the rank of this process.
Postconditions:
0 <= return && return < this->getSize()
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Returns the number of processes that make up this communicator.
Postconditions:
return > 0
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Pause every process in *this
communicator until all the processes reach this point.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Broadcast values from the root process to the slave processes.
rootRank | [in] The rank of the root process. |
count | [in] The number of bytes in buffer[] . |
buffer | [in/out] Array (length bytes ) of packed data. Must be set on input on the root processes with rank root . On output, each processs, including the root process contains the data. |
Preconditions:
0 <= rootRank && rootRank < this->getSize()
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Gather values from all processes to the root process.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Gather values from each process to collect on all processes.
sendBytes | [in] Number of entires in sendBuffer[] on input. |
sendBuffer | [in] Array (length sendBytes ) of data being sent from each process. |
recvBytes | [in] Number of entires in recvBuffer[] which must be equal to sendBytes*this->getSize() . This field is just here for debug checking. |
recvBuffer | [out] Array (length recvBytes ) of all of the entires sent from each processes. Specifically, recvBuffer[sendBytes*j+i] , for j=0...this->getSize()-1 and i=0...sendBytes-1 , is the entry sendBuffer[i] from process with rank j . |
Preconditions:
recvBytes==sendBytes*this->getSize()
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Global reduction.
reductOp | [in] The user-defined reduction operation |
bytes | [in] The length of the buffers sendBuffer[] and globalReducts[] . |
sendBuffer | [in] Array (length bytes ) of the data contributed from each process. |
globalReducts | [out] Array (length bytes ) of the global reduction from each process. |
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Scan reduction.
reductOp | [in] The user-defined reduction operation |
bytes | [in] The length of the buffers sendBuffer[] and scanReducts[] . |
sendBuffer | [in] Array (length bytes ) of the data contributed from each process. |
scanReducts | [out] Array (length bytes ) of the reduction up to and including this process. |
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Possibly blocking send of data from this process to another process.
This routine does not return until you can reuse the send buffer. Whether this routine blocks depends on whether the MPI implementation buffers.
bytes | [in] The number of bytes of data being passed between processes. |
sendBuffer | [in] Array (length bytes ) of data being sent from this process. This buffer can be immediately destroyed or reused as soon as the function exits (that is why this function is "blocking"). |
destRank | [in] The rank of the process to receive the data. |
Preconditions:
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Variant of send() that takes a tag.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Always blocking send of data from this process to another process.
This routine blocks until the matching receive posts. After it returns, you are allowed to reuse the send buffer.
bytes | [in] The number of bytes of data being passed between processes. |
sendBuffer | [in] Array (length bytes ) of data being sent from this process. This buffer can be immediately destroyed or reused as soon as the function exits (that is why this function is "blocking"). |
destRank | [in] The rank of the process to receive the data. |
Preconditions:
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Variant of ssend() that takes a message tag.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Blocking receive of data from this process to another process.
sourceRank | [in] The rank of the process to receive the data from. If sourceRank < 0 then data will be received from any process. |
bytes | [in] The number of bytes of data being passed between processes. |
recvBuffer | [out] Array (length bytes ) of data being received from this process. This buffer can be immediately used to access the data as soon as the function exits (that is why this function is "blocking"). |
Preconditions:
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Ready send of data from this process to another process.
sendBuffer | [in] The data to be sent. |
destRank | [in] The rank of the process to receive the data. |
Preconditions:
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Variant of readySend() that accepts a message tag.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Non-blocking send.
sendBuffer | [in] The data buffer to be sent. |
destRank | [in] The rank of the process to receive the data. |
Preconditions:
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Variant of isend() that takes a tag.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Non-blocking receive.
recvBuffer | [out] The location for storing the received data. |
sourceRank | [in] The rank of the process to receive the data from. If sourceRank < 0 then data will be received from any process. |
Preconditions:
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Variant of ireceive that takes a tag.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Wait on a set of communication requests.
Preconditions:
requests.size() > 0
Postconditions:
is_null(request[i]))
for i=0...requests.size()-1
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Wait on communication requests, and return their statuses.
requests | [in/out] On input: the requests on which to wait. On output: all set to null. |
statuses | [out] The status results of waiting on the requests. |
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Wait on a single communication request, and return its status.
request | [in/out] On input: request is not null, and *request is either null (in which case this function does nothing and returns null) or an RCP of a valid CommRequest instance representing an outstanding communication request. On output: If the communication request completed successfully, we set *request to null, indicating that the request has completed. (This helps prevent common bugs like trying to complete the same request twice.) |
MPI_ANY_SOURCE
.)!is_null(request)
(that is, the Ptr is not null). is_null(*request)
(that is, the RCP is null).This function blocks until the communication operation associated with the CommRequest object has completed.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Duplicate this communicator.
Make a copy of this communicator with a duplicate communication space. Note that the returned communicator has the same properties (including process ranks, attributes and topologies) as this communicator, but is distinct from the original. "Distinct" means that if you send a message on the original communicator, you can't receive it on the new one, and vice versa. The new communicator represents a separate message space. This has the same semantics as MPI_Comm_dup. (In fact, the subclass MpiComm implements this using MPI_Comm_dup.)
Most users don't want to do this. The duplicate() method returns a new communicator. In MPI terms, it is a different MPI_Comm. If you want a shallow copy of the handle, you should pass the Comm<Ordinal>
around by const pointer, like this:
This behaves the same as the following "raw MPI" code:
The subclass of Comm ensures that the "raw" MPI handle is freed only after the last reference to it by a subclass instance disappears. (It does reference counting underneath.)
Please, please do not invoke the copy constructor or assignment operator of Comm. Of course it's not legal to do that anyway, because Comm is pure virtual. However, even if you could do it, you must never do this! For example, do not do this:
and do not do this:
This is bad because it ignores the subclass' data. Depending on the subclass of Comm that you are actually using, it may be appropriate to invoke the copy constructor or assignment operator of the specific subclass, but never those of Comm itself.
Users are not responsible for freeing the returned communicator. The destructor of the subclass of Comm handles that itself.
In an MPI implementation, the returned communicator is created using MPI_Comm_dup, with the resulting semantic implications.
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Split a communicator into subcommunicators based on color and key.
Partition this communicator into multiple disjoint groups, and return the communicator corresponding to the group to which this process belongs. There will be as many groups as there are globally many distinct values for the color
parameter. Within each subset of the partition, the ranks will be ordered according to the key value each process passed for the key
parameter. If multiple processes pass the same value for key
, then they will be ordered according to their rank in the original communicator. To return a valid communicator, this function requires a nonnegative value for color
. If color
is negative, this method will return a null communicator.
This method must be called as a collective on all processes in this communicator. That is, if this method is called at all, it must be called on all processes in the communicator.
Users are not responsible for freeing the returned communicator. The destructor of the subclass of Comm handles that itself.
In an MPI implementation, the returned communicator is created using MPI_Comm_split, with the resulting semantic implications.
color | [in] An integer representing the color for the local rank. In the MPI implementation, if this is negative, MPI_Comm_split gets MPI_UNDEFINED as the color. |
key | [in] A key value to order processes of the same color. In the MPI implementation, this is passed directly to MPI_Comm_split. |
Implemented in Teuchos::SerialComm< Ordinal >.
|
pure virtual |
Create a subcommunicator containing the specified processes.
Create and return a subcommunicator of this communicator. The subcommunicator contains the processes in this communicator with the given ranks, in which they are listed in the input vector. Processes whose ranks are not included in the input vector will be given a null communicator.
This method must be called as a collective on all processes in this communicator. That is, if this method is called at all, it must be called on all processes in the communicator.
Users are not responsible for freeing the returned communicator. The destructor of the subclass of Comm handles that itself.
In an MPI implementation, the subcommunicator is created using MPI_Comm_create, with the resulting semantic implications.
ranks | The ranks of the processes to include in the subcommunicator. |
Implemented in Teuchos::SerialComm< Ordinal >.
|
related |
Predefined reduction operations that Teuchos::Comm understands.
Teuchos' Comm class wraps MPI (the Message Passing Interface for distributed-memory parallel programming). If you do not have MPI, it imitates MPI's functionality, as if you were running with a single "parallel" process. This means that Teuchos must wrap a subset of MPI's functionality, so that it can build without MPI.
Comm provides wrappers for MPI_Reduce
, MPI_Allreduce
, and other collectives that take a reduction operator MPI_Op
. Teuchos wraps MPI_Op
in two different ways. The first way is this enum, which lets users pick from a set of common predefined MPI_Op
. The second way is through Teuchos' wrappers for custom MPI_Op
, namely ValueTypeReductionOp and ValueTypeReductionOp. Most users should find the reduction operators below sufficient.
Enumerator | |
---|---|
REDUCE_SUM |
Sum. |
REDUCE_MIN |
Min. |
REDUCE_MAX |
Max. |
REDUCE_AND |
Logical AND. |
REDUCE_BOR |
Bitwise OR. |
Definition at line 71 of file Teuchos_EReductionType.hpp.
|
related |
Get the process rank.
|
related |
Get the number of processes in the communicator.
|
related |
Barrier.
|
related |
Broadcast array of objects that use value semantics.
|
related |
Broadcast array of objects that use value semantics.
|
related |
Broadcast single object that use value semantics.
|
related |
Broadcast single object that use value semantics.
|
related |
Broadcast array of objects that use reference semantics.
|
related |
Broadcast array of objects that use reference semantics.
|
related |
Broadcast array of objects that use value semantics using customized serializer.
|
related |
Gather values from each process to the root process.
This wraps MPI_Gather in an MPI build, when Comm implements MpiComm.
|
related |
Gather arrays of possibly different lengths from each process to the root process.
This wraps MPI_Gatherv in an MPI build, when Comm implements MpiComm.
|
related |
Gather array of objects that use value semantics from every process to every process.
|
related |
Gather array of objects that use reference semantics from every process to every process.
|
related |
Gather array of objects that use value semantics from every process to every process using customized serializer.
|
related |
Wrapper for MPI_Scatter; scatter collective.
Ordinal | The template parameter of Comm. Should always be int unless you REALLY know what you are doing. |
Packet | The type of the thing this operation communicates. For this particular overload of scatter(), Packet must have "value semantics" and must not require a custom Serializer . Built-in types and structs thereof are generally OK here. |
If Comm is an MpiComm, then this wraps MPI_Scatter(). This function is a collective operation over the input communicator.
sendBuf | [in] Array of Packet things to scatter out. This may NOT alias recvBuf . This is ONLY read on the root process of the collective. |
sendCount | The number of Packet things for the root process to send to each process in the input communicator. |
recvBuf | [out] Array of recvCount*comm.getSize() Packet things to receive. This may NOT alias sendBuf . |
recvCount | [in] Number of Packet things to receive from each process. |
root | [in] Rank of the process from which to scatter (the root of the scatter operation). The rank is relative to the input communicator. |
comm | [in] The communicator over which to scatter. |
Definition at line 271 of file Teuchos_CommHelpers.hpp.
|
related |
Wrapper for MPI_Reduce; reduction to one process, using a built-in reduction operator selected by enum.
Ordinal | The template parameter of Comm. Should always be int unless you REALLY know what you are doing. |
Packet | The type of the thing this operation communicates. For this particular overload of reduce(), Packet must have "value semantics" and must not require a custom Serializer . Built-in types and structs thereof are generally OK here. |
If Comm is an MpiComm, then this wraps MPI_Reduce(). This function is a collective operation over the input communicator.
sendBuf | [in] Array of count Packet things to send. This may NOT alias recvBuf . |
recvBuf | [in] Array of count Packet reduction results. This may NOT alias sendBuf . |
count | [in] Number of entries in the sendBuf and recvBuf arrays. |
reductType | [in] Type of reduction operation. Valid values include REDUCE_SUM, REDUCE_MIN, REDUCE_MAX, and REDUCE_AND. See the documentation of the EReductionType enum for details. |
root | [in] Rank of the process on which to receive the reduction result. The rank is relative to the input communicator. |
comm | [in] The communicator over which to reduce. |
|
related |
Wrapper for MPI_Allreduce that takes a custom reduction operator.
Ordinal | The template parameter of Comm. Should always be int unless you REALLY know what you are doing. |
Packet | The type of the thing this operation communicates. For this particular overload of reduce(), Packet must have "value semantics" and must not require a custom Serializer . Built-in types and structs thereof are generally OK here. |
If Comm is an MpiComm, then this wraps MPI_Allreduce(). This function is a collective operation over the input communicator.
comm | [in] The communicator over which to reduce. |
reductOp | [in] The custom reduction operator. |
count | [in] Number of entries in the sendBuffer and globalReducts arrays. |
sendBuffer | [in] Array of count Packet things to send. This may NOT alias globalReducts . |
globalReducts | [in] Array of count Packet reduction results. This may NOT alias recvBuf . |
|
related |
Collective reduce all of array of objects using value semantics using a pre-defined reduction type.
|
related |
Collective reduce all for single object using value semantics using a pre-defined reduction type.
|
related |
Collective reduce all for array of objects using reference semantics.
|
related |
Collective reduce all of array of objects using value semantics using a user-defined reduction operator and customized serializer.
|
related |
Collective reduce all of array of objects using value semantics using a pre-defined reduction type and customized serializer.
|
related |
Scan/Reduce array of objects that use value semantics using a user-defined reduction operator.
|
related |
Scan/Reduce array of objects using value semantics using a predefined reduction type.
|
related |
Scan/Reduce single object using value semantics using a predefined reduction type.
|
related |
Scan/Reduce array of objects that use reference semantics using a user-defined reduction operator.
|
related |
Scan/Reduce array of objects that use value semantics using a user-defined reduction operator and customized serializer.
|
related |
Scan/Reduce array of objects using value semantics using a predefined reduction type and customized serializer.
|
related |
Send objects that use values semantics to another process.
|
related |
Synchronously send objects that use values semantics to another process.
|
related |
Send a single object that use values semantics to another process.
|
related |
Synchronously send a single object that use values semantics to another process.
|
related |
Send objects that use reference semantics to another process.
NOTE: Not implemented yet!
|
related |
Send objects that use values semantics to another process using customized serializer.
|
related |
Receive objects that use values semantics from another process.
|
related |
Receive a single object that use values semantics from another process.
|
related |
Receive objects that use reference semantics from another process.
|
related |
Receive objects that use values semantics from another process using customized serializer.
|
related |
Ready-Send an array of objects that use values semantics to another process.
|
related |
Ready-Send a single object that use values semantics to another process.
|
related |
Ready-Send an array of objects that use values semantics to another process using customized serializer.
|
related |
Send objects that use values semantics to another process.
|
related |
Send a single object that use values semantics to another process.
|
related |
Send objects that use values semantics to another process using customized serializer.
|
related |
Receive one or more objects (that use values semantics) from another process.
comm | [in] The communicator. |
recvBuffer | [out] The buffer into which to receive the data. |
sourceRank | [in] The rank of the sending process. A negative source rank means that this will accept an incoming message from any process on the given communicator. (This is the equivalent of MPI_ANY_SOURCE.) |
|
related |
Receive one object (that uses values semantics) from another process.
comm | [in] The communicator. |
recv | [out] The buffer into which to receive the object. |
sourceRank | [in] The rank of the sending process. A negative source rank means that this will accept an incoming message from any process on the given communicator. |
|
related |
Send objects that use values semantics to another process using customized serializer.
|
related |
Wait for an array of Teuchos::CommRequest objects.
Blocks until all communication operations associated with the CommRequest objects have completed.
|
related |
Wait on one or more communication requests, and return their statuses.
This function will block until all the communication operations complete. You are responsible for ordering messages in a way that avoids deadlock. An erroneous ordering of messages may cause this function to wait forever.
requests.size() == statuses.size()
requests.size() - 1
, requests[i]
is either null or it was the return value of a nonblocking communication request.requests.size() - 1
, requests[i]
is null. If requests[i]
were nonnull on input, then the corresponding nonblocking communication operation completed successfully.comm | [in] The communicator on which the communication requests were made. |
requests | [in/out] On input: the requests on which to wait. Null elements of the array are ignored. On output: all elements of the array are set to null. |
statuses | [out] CommStatus instances representing the results of completing the requests. You may query each for information about its corresponding communication operation. Any element of the array may be null, for example if its corresponding CommRequest was null on input. |
|
related |
Wait on a single communication request, and return its status.
This function will block until the communication operation completes. You are responsible for ordering messages in a way that avoids deadlock. An erroneous ordering of messages may cause this function to wait forever.
request | [in/out] A nonnull pointer to an RCP. On input, the RCP is either null (in which case this function does nothing) or it points to a valid CommRequest representing an outstanding communication request. On output, this function sets *request to null, which indicates that the request completed successfully. (This function will not return unless the request completes.) |
*request
was null on input.!is_null(request)
(that is, the Ptr is not null). is_null(*request)
(that is, the RCP is null).