Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Epetra_MpiComm Class Reference

Epetra_MpiComm: The Epetra MPI Communication Class. More...

#include <Epetra_MpiComm.h>

Inheritance diagram for Epetra_MpiComm:
Inheritance graph
[legend]

Public Member Functions

Epetra_MpiCommoperator= (const Epetra_MpiComm &Comm)
 Assignment Operator. More...
 
- Public Member Functions inherited from Epetra_Object
 Epetra_Object (int TracebackModeIn=-1, bool set_label=true)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const char *const Label, int TracebackModeIn=-1)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const Epetra_Object &Object)
 Epetra_Object Copy Constructor. More...
 
virtual ~Epetra_Object ()
 Epetra_Object Destructor. More...
 
virtual void SetLabel (const char *const Label)
 Epetra_Object Label definition using char *. More...
 
virtual const char * Label () const
 Epetra_Object Label access funtion. More...
 
virtual int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method. More...
 
- Public Member Functions inherited from Epetra_Comm
virtual ~Epetra_Comm ()
 Epetra_Comm Destructor. More...
 

Private Member Functions

int CheckInput (double *ptr, int count) const
 
int CheckInput (int *ptr, int count) const
 
int CheckInput (long *ptr, int count) const
 
int CheckInput (char *ptr, int count) const
 
int CheckInput (long long *ptr, int count) const
 
void CleanupData ()
 

Private Attributes

Epetra_MpiCommDataMpiCommData_
 

Constructor/Destructor Methods

 Epetra_MpiComm (MPI_Comm comm)
 Epetra_MpiComm MPI Constructor. More...
 
 Epetra_MpiComm (const Epetra_MpiComm &Comm)
 Epetra_MpiComm Copy Constructor. More...
 
Epetra_CommClone () const
 Clone method. More...
 
virtual ~Epetra_MpiComm ()
 Epetra_MpiComm Destructor. More...
 

Barrier Methods

void Barrier () const
 Epetra_MpiComm Barrier function. More...
 

Broadcast Methods

int Broadcast (double *MyVals, int Count, int Root) const
 Epetra_MpiComm Broadcast function. More...
 
int Broadcast (int *MyVals, int Count, int Root) const
 Epetra_MpiComm Broadcast function. More...
 
int Broadcast (long *MyVals, int Count, int Root) const
 Epetra_MpiComm Broadcast function. More...
 
int Broadcast (long long *MyVals, int Count, int Root) const
 Epetra_MpiComm Broadcast function. More...
 
int Broadcast (char *MyVals, int Count, int Root) const
 Epetra_MpiComm Broadcast function. More...
 

Gather Methods

int GatherAll (double *MyVals, double *AllVals, int Count) const
 Epetra_MpiComm All Gather function. More...
 
int GatherAll (int *MyVals, int *AllVals, int Count) const
 Epetra_MpiComm All Gather function. More...
 
int GatherAll (long *MyVals, long *AllVals, int Count) const
 Epetra_MpiComm All Gather function. More...
 
int GatherAll (long long *MyVals, long long *AllVals, int Count) const
 Epetra_MpiComm All Gather function. More...
 

Sum Methods

int SumAll (double *PartialSums, double *GlobalSums, int Count) const
 Epetra_MpiComm Global Sum function. More...
 
int SumAll (int *PartialSums, int *GlobalSums, int Count) const
 Epetra_MpiComm Global Sum function. More...
 
int SumAll (long *PartialSums, long *GlobalSums, int Count) const
 Epetra_MpiComm Global Sum function. More...
 
int SumAll (long long *PartialSums, long long *GlobalSums, int Count) const
 Epetra_MpiComm Global Sum function. More...
 

Max/Min Methods

int MaxAll (double *PartialMaxs, double *GlobalMaxs, int Count) const
 Epetra_MpiComm Global Max function. More...
 
int MaxAll (int *PartialMaxs, int *GlobalMaxs, int Count) const
 Epetra_MpiComm Global Max function. More...
 
int MaxAll (long *PartialMaxs, long *GlobalMaxs, int Count) const
 Epetra_MpiComm Global Max function. More...
 
int MaxAll (long long *PartialMaxs, long long *GlobalMaxs, int Count) const
 Epetra_MpiComm Global Max function. More...
 
int MinAll (double *PartialMins, double *GlobalMins, int Count) const
 Epetra_MpiComm Global Min function. More...
 
int MinAll (int *PartialMins, int *GlobalMins, int Count) const
 Epetra_MpiComm Global Min function. More...
 
int MinAll (long *PartialMins, long *GlobalMins, int Count) const
 Epetra_MpiComm Global Min function. More...
 
int MinAll (long long *PartialMins, long long *GlobalMins, int Count) const
 Epetra_MpiComm Global Min function. More...
 

Parallel Prefix Methods

int ScanSum (double *MyVals, double *ScanSums, int Count) const
 Epetra_MpiComm Scan Sum function. More...
 
int ScanSum (int *MyVals, int *ScanSums, int Count) const
 Epetra_MpiComm Scan Sum function. More...
 
int ScanSum (long *MyVals, long *ScanSums, int Count) const
 Epetra_MpiComm Scan Sum function. More...
 
int ScanSum (long long *MyVals, long long *ScanSums, int Count) const
 Epetra_MpiComm Scan Sum function. More...
 

Attribute Accessor Methods

MPI_Comm Comm () const
 Extract MPI Communicator from a Epetra_MpiComm object. More...
 
int MyPID () const
 Return my process ID. More...
 
int NumProc () const
 Returns total number of processes. More...
 

Gather/Scatter and Directory Constructors

Epetra_DistributorCreateDistributor () const
 Create a distributor object. More...
 
Epetra_DirectoryCreateDirectory (const Epetra_BlockMap &Map) const
 Create a directory object for the given Epetra_BlockMap. More...
 

MPI-specific Methods

int GetMpiTag () const
 Acquire an MPI tag from the Epetra range of 24050-24099, increment tag. More...
 
MPI_Comm GetMpiComm () const
 Get the MPI Communicator (identical to Comm() method; used when we know we are MPI. More...
 

Print object to an output stream

void Print (std::ostream &os) const
 Print method that implements Epetra_Object virtual Print method. More...
 
void PrintInfo (std::ostream &os) const
 Print method that implements Epetra_Comm virtual PrintInfo method. More...
 

Expert Users and Developers Only

int ReferenceCount () const
 Returns the reference count of MpiCommData. More...
 
const Epetra_MpiCommDataDataPtr () const
 Returns a pointer to the MpiCommData instance this MpiComm uses. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from Epetra_Object
static void SetTracebackMode (int TracebackModeValue)
 Set the value of the Epetra_Object error traceback report mode. More...
 
static int GetTracebackMode ()
 Get the value of the Epetra_Object error report mode. More...
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting. More...
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 
- Protected Member Functions inherited from Epetra_Object
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 

Detailed Description

Epetra_MpiComm: The Epetra MPI Communication Class.

The Epetra_MpiComm class is an implementation of Epetra_Comm that encapsulates the general information and services needed for other Epetra classes to run on a parallel computer using MPI.

Definition at line 64 of file Epetra_MpiComm.h.

Constructor & Destructor Documentation

Epetra_MpiComm::Epetra_MpiComm ( MPI_Comm  comm)

Epetra_MpiComm MPI Constructor.

Creates a Epetra_MpiComm instance for use with MPI. If no specialized MPI communicator is needed, this constuctor can be called with the argument MPI_COMM_WORLD.

Definition at line 46 of file Epetra_MpiComm.cpp.

Epetra_MpiComm::Epetra_MpiComm ( const Epetra_MpiComm Comm)

Epetra_MpiComm Copy Constructor.

Makes an exact copy of an existing Epetra_MpiComm instance.

Definition at line 53 of file Epetra_MpiComm.cpp.

Epetra_MpiComm::~Epetra_MpiComm ( )
virtual

Epetra_MpiComm Destructor.

Completely deletes a Epetra_MpiComm object.

Warning
Note: All objects that depend on a Epetra_MpiComm instance should be destroyed prior to calling this function.

Definition at line 247 of file Epetra_MpiComm.cpp.

Member Function Documentation

Epetra_Comm* Epetra_MpiComm::Clone ( ) const
inlinevirtual

Clone method.

Implements Epetra_Comm.

Definition at line 84 of file Epetra_MpiComm.h.

void Epetra_MpiComm::Barrier ( ) const
virtual

Epetra_MpiComm Barrier function.

Causes each processor in the communicator to wait until all processors have arrived.

Implements Epetra_Comm.

Definition at line 61 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::Broadcast ( double *  MyVals,
int  Count,
int  Root 
) const
virtual

Epetra_MpiComm Broadcast function.

Takes list of input values from the root processor and sends to all other processors.

Parameters
ValuesInOut On entry, the root processor contains the list of values. On exit, all processors will have the same list of values. Note that values must be allocated on all processor before the broadcast.
CountIn On entry, contains the length of the list of Values.
RootIn On entry, contains the processor from which all processors will receive a copy of Values.

Implements Epetra_Comm.

Definition at line 65 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::Broadcast ( int *  MyVals,
int  Count,
int  Root 
) const
virtual

Epetra_MpiComm Broadcast function.

Take list of input values from the root processor and sends to all other processors.

Parameters
ValuesInOut On entry, the root processor contains the list of values. On exit, all processors will have the same list of values. Note that values must be allocated on all processor before the broadcast.
CountIn On entry, contains the length of the list of Values.
RootIn On entry, contains the processor from which all processors will receive a copy of Values.

Implements Epetra_Comm.

Definition at line 71 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::Broadcast ( long *  MyVals,
int  Count,
int  Root 
) const
virtual

Epetra_MpiComm Broadcast function.

Take list of input values from the root processor and sends to all other processors.

Parameters
ValuesInOut On entry, the root processor contains the list of values. On exit, all processors will have the same list of values. Note that values must be allocated on all processor before the broadcast.
CountIn On entry, contains the length of the list of Values.
RootIn On entry, contains the processor from which all processors will receive a copy of Values.

Implements Epetra_Comm.

Definition at line 77 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::Broadcast ( long long *  MyVals,
int  Count,
int  Root 
) const
virtual

Epetra_MpiComm Broadcast function.

Take list of input values from the root processor and sends to all other processors.

Parameters
ValuesInOut On entry, the root processor contains the list of values. On exit, all processors will have the same list of values. Note that values must be allocated on all processor before the broadcast.
CountIn On entry, contains the length of the list of Values.
RootIn On entry, contains the processor from which all processors will receive a copy of Values.

Implements Epetra_Comm.

Definition at line 83 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::Broadcast ( char *  MyVals,
int  Count,
int  Root 
) const
virtual

Epetra_MpiComm Broadcast function.

Take list of input values from the root processor and sends to all other processors.

Parameters
ValuesInOut On entry, the root processor contains the list of values. On exit, all processors will have the same list of values. Note that values must be allocated on all processor before the broadcast.
CountIn On entry, contains the length of the list of Values.
RootIn On entry, contains the processor from which all processors will receive a copy of Values.

Implements Epetra_Comm.

Definition at line 89 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::GatherAll ( double *  MyVals,
double *  AllVals,
int  Count 
) const
virtual

Epetra_MpiComm All Gather function.

Take list of input values from all processors in the communicator and creates an ordered contiguous list of those values on each processor.

Parameters
MyValsIn On entry, contains the list of values, to be sent to all processors.
AllValsOut On exit, contains the list of values from all processors. Must by of size NumProc*Count.
CountIn On entry, contains the length of the list of MyVals.

Implements Epetra_Comm.

Definition at line 95 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::GatherAll ( int *  MyVals,
int *  AllVals,
int  Count 
) const
virtual

Epetra_MpiComm All Gather function.

Take list of input values from all processors in the communicator and creates an ordered contiguous list of those values on each processor.

Parameters
MyValsIn On entry, contains the list of values, to be sent to all processors.
AllValsOut On exit, contains the list of values from all processors. Must by of size NumProc*Count.
CountIn On entry, contains the length of the list of MyVals.

Implements Epetra_Comm.

Definition at line 102 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::GatherAll ( long *  MyVals,
long *  AllVals,
int  Count 
) const
virtual

Epetra_MpiComm All Gather function.

Take list of input values from all processors in the communicator and creates an ordered contiguous list of those values on each processor.

Parameters
MyValsIn On entry, contains the list of values, to be sent to all processors.
AllValsOut On exit, contains the list of values from all processors. Must by of size NumProc*Count.
CountIn On entry, contains the length of the list of MyVals.

Implements Epetra_Comm.

Definition at line 109 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::GatherAll ( long long *  MyVals,
long long *  AllVals,
int  Count 
) const
virtual

Epetra_MpiComm All Gather function.

Take list of input values from all processors in the communicator and creates an ordered contiguous list of those values on each processor.

Parameters
MyValsIn On entry, contains the list of values, to be sent to all processors.
AllValsOut On exit, contains the list of values from all processors. Must by of size NumProc*Count.
CountIn On entry, contains the length of the list of MyVals.

Implements Epetra_Comm.

Definition at line 116 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::SumAll ( double *  PartialSums,
double *  GlobalSums,
int  Count 
) const
virtual

Epetra_MpiComm Global Sum function.

Take list of input values from all processors in the communicator, computes the sum and returns the sum to all processors.

Parameters
PartialSumsIn On entry, contains the list of values, usually partial sums computed locally, to be summed across all processors.
GlobalSumsOut On exit, contains the list of values summed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 123 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::SumAll ( int *  PartialSums,
int *  GlobalSums,
int  Count 
) const
virtual

Epetra_MpiComm Global Sum function.

Take list of input values from all processors in the communicator, computes the sum and returns the sum to all processors.

Parameters
PartialSumsIn On entry, contains the list of values, usually partial sums computed locally, to be summed across all processors.
GlobalSumsOut On exit, contains the list of values summed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 130 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::SumAll ( long *  PartialSums,
long *  GlobalSums,
int  Count 
) const
virtual

Epetra_MpiComm Global Sum function.

Take list of input values from all processors in the communicator, computes the sum and returns the sum to all processors.

Parameters
PartialSumsIn On entry, contains the list of values, usually partial sums computed locally, to be summed across all processors.
GlobalSumsOut On exit, contains the list of values summed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 137 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::SumAll ( long long *  PartialSums,
long long *  GlobalSums,
int  Count 
) const
virtual

Epetra_MpiComm Global Sum function.

Take list of input values from all processors in the communicator, computes the sum and returns the sum to all processors.

Parameters
PartialSumsIn On entry, contains the list of values, usually partial sums computed locally, to be summed across all processors.
GlobalSumsOut On exit, contains the list of values summed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 144 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MaxAll ( double *  PartialMaxs,
double *  GlobalMaxs,
int  Count 
) const
virtual

Epetra_MpiComm Global Max function.

Take list of input values from all processors in the communicator, computes the max and returns the max to all processors.

Parameters
PartialMaxsIn On entry, contains the list of values, usually partial maxs computed locally; using these Partial Maxs, the max across all processors will be computed.
GlobalMaxsOut On exit, contains the list of maxs computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 151 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MaxAll ( int *  PartialMaxs,
int *  GlobalMaxs,
int  Count 
) const
virtual

Epetra_MpiComm Global Max function.

Take list of input values from all processors in the communicator, computes the max and returns the max to all processors.

Parameters
PartialMaxsIn On entry, contains the list of values, usually partial maxs computed locally; using these Partial Maxs, the max across all processors will be computed.
GlobalMaxsOut On exit, contains the list of maxs computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 158 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MaxAll ( long *  PartialMaxs,
long *  GlobalMaxs,
int  Count 
) const
virtual

Epetra_MpiComm Global Max function.

Take list of input values from all processors in the communicator, computes the max and returns the max to all processors.

Parameters
PartialMaxsIn On entry, contains the list of values, usually partial maxs computed locally; using these Partial Maxs, the max across all processors will be computed.
GlobalMaxsOut On exit, contains the list of maxs computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 165 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MaxAll ( long long *  PartialMaxs,
long long *  GlobalMaxs,
int  Count 
) const
virtual

Epetra_MpiComm Global Max function.

Take list of input values from all processors in the communicator, computes the max and returns the max to all processors.

Parameters
PartialMaxsIn On entry, contains the list of values, usually partial maxs computed locally; using these Partial Maxs, the max across all processors will be computed.
GlobalMaxsOut On exit, contains the list of maxs computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 172 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MinAll ( double *  PartialMins,
double *  GlobalMins,
int  Count 
) const
virtual

Epetra_MpiComm Global Min function.

Take list of input values from all processors in the communicator, computes the min and returns the min to all processors.

Parameters
PartialMinsIn On entry, contains the list of values, usually partial mins computed locally; using these Partial Mins, the min across all processors will be computed.
GlobalMinsOut On exit, contains the list of mins computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 179 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MinAll ( int *  PartialMins,
int *  GlobalMins,
int  Count 
) const
virtual

Epetra_MpiComm Global Min function.

Take list of input values from all processors in the communicator, computes the min and returns the min to all processors.

Parameters
PartialMinsIn On entry, contains the list of values, usually partial mins computed locally; using these Partial Mins, the min across all processors will be computed.
GlobalMinsOut On exit, contains the list of mins computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 186 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MinAll ( long *  PartialMins,
long *  GlobalMins,
int  Count 
) const
virtual

Epetra_MpiComm Global Min function.

Take list of input values from all processors in the communicator, computes the min and returns the min to all processors.

Parameters
PartialMinsIn On entry, contains the list of values, usually partial mins computed locally; using these Partial Mins, the min across all processors will be computed.
GlobalMinsOut On exit, contains the list of mins computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 193 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::MinAll ( long long *  PartialMins,
long long *  GlobalMins,
int  Count 
) const
virtual

Epetra_MpiComm Global Min function.

Take list of input values from all processors in the communicator, computes the min and returns the min to all processors.

Parameters
PartialMinsIn On entry, contains the list of values, usually partial mins computed locally; using these Partial Mins, the min across all processors will be computed.
GlobalMinsOut On exit, contains the list of mins computed across all processors.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 200 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::ScanSum ( double *  MyVals,
double *  ScanSums,
int  Count 
) const
virtual

Epetra_MpiComm Scan Sum function.

Take list of input values from all processors in the communicator, computes the scan sum and returns it to all processors such that processor i contains the sum of values from processor 0 up to and including processor i.

Parameters
MyValsIn On entry, contains the list of values to be summed across all processors.
ScanSumsOut On exit, contains the list of values summed across processors 0 through i.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 207 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::ScanSum ( int *  MyVals,
int *  ScanSums,
int  Count 
) const
virtual

Epetra_MpiComm Scan Sum function.

Take list of input values from all processors in the communicator, computes the scan sum and returns it to all processors such that processor i contains the sum of values from processor 0 up to and including processor i.

Parameters
MyValsIn On entry, contains the list of values to be summed across all processors.
ScanSumsOut On exit, contains the list of values summed across processors 0 through i.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 214 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::ScanSum ( long *  MyVals,
long *  ScanSums,
int  Count 
) const
virtual

Epetra_MpiComm Scan Sum function.

Take list of input values from all processors in the communicator, computes the scan sum and returns it to all processors such that processor i contains the sum of values from processor 0 up to and including processor i.

Parameters
MyValsIn On entry, contains the list of values to be summed across all processors.
ScanSumsOut On exit, contains the list of values summed across processors 0 through i.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 221 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::ScanSum ( long long *  MyVals,
long long *  ScanSums,
int  Count 
) const
virtual

Epetra_MpiComm Scan Sum function.

Take list of input values from all processors in the communicator, computes the scan sum and returns it to all processors such that processor i contains the sum of values from processor 0 up to and including processor i.

Parameters
MyValsIn On entry, contains the list of values to be summed across all processors.
ScanSumsOut On exit, contains the list of values summed across processors 0 through i.
CountIn On entry, contains the length of the list of values.

Implements Epetra_Comm.

Definition at line 228 of file Epetra_MpiComm.cpp.

MPI_Comm Epetra_MpiComm::Comm ( ) const
inline

Extract MPI Communicator from a Epetra_MpiComm object.

Definition at line 457 of file Epetra_MpiComm.h.

int Epetra_MpiComm::MyPID ( ) const
inlinevirtual

Return my process ID.

In MPI mode returns the rank of the calling process. In serial mode returns 0.

Implements Epetra_Comm.

Definition at line 463 of file Epetra_MpiComm.h.

int Epetra_MpiComm::NumProc ( ) const
inlinevirtual

Returns total number of processes.

In MPI mode returns the size of the MPI communicator. In serial mode returns 1.

Implements Epetra_Comm.

Definition at line 469 of file Epetra_MpiComm.h.

Epetra_Distributor * Epetra_MpiComm::CreateDistributor ( ) const
virtual

Create a distributor object.

Implements Epetra_Comm.

Definition at line 235 of file Epetra_MpiComm.cpp.

Epetra_Directory * Epetra_MpiComm::CreateDirectory ( const Epetra_BlockMap Map) const
virtual

Create a directory object for the given Epetra_BlockMap.

Implements Epetra_Comm.

Definition at line 241 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::GetMpiTag ( ) const
inline

Acquire an MPI tag from the Epetra range of 24050-24099, increment tag.

Definition at line 487 of file Epetra_MpiComm.h.

MPI_Comm Epetra_MpiComm::GetMpiComm ( ) const
inline

Get the MPI Communicator (identical to Comm() method; used when we know we are MPI.

Definition at line 490 of file Epetra_MpiComm.h.

void Epetra_MpiComm::Print ( std::ostream &  os) const
inlinevirtual

Print method that implements Epetra_Object virtual Print method.

Reimplemented from Epetra_Object.

Definition at line 495 of file Epetra_MpiComm.h.

void Epetra_MpiComm::PrintInfo ( std::ostream &  os) const
inlinevirtual

Print method that implements Epetra_Comm virtual PrintInfo method.

Implements Epetra_Comm.

Definition at line 510 of file Epetra_MpiComm.h.

int Epetra_MpiComm::ReferenceCount ( ) const
inline

Returns the reference count of MpiCommData.

(Intended for testing purposes.)

Definition at line 518 of file Epetra_MpiComm.h.

const Epetra_MpiCommData* Epetra_MpiComm::DataPtr ( ) const
inline

Returns a pointer to the MpiCommData instance this MpiComm uses.

(Intended for developer use only for testing purposes.)

Definition at line 522 of file Epetra_MpiComm.h.

Epetra_MpiComm & Epetra_MpiComm::operator= ( const Epetra_MpiComm Comm)

Assignment Operator.

Definition at line 261 of file Epetra_MpiComm.cpp.

int Epetra_MpiComm::CheckInput ( double *  ptr,
int  count 
) const
inlineprivate

Definition at line 531 of file Epetra_MpiComm.h.

int Epetra_MpiComm::CheckInput ( int *  ptr,
int  count 
) const
inlineprivate

Definition at line 532 of file Epetra_MpiComm.h.

int Epetra_MpiComm::CheckInput ( long *  ptr,
int  count 
) const
inlineprivate

Definition at line 533 of file Epetra_MpiComm.h.

int Epetra_MpiComm::CheckInput ( char *  ptr,
int  count 
) const
inlineprivate

Definition at line 534 of file Epetra_MpiComm.h.

int Epetra_MpiComm::CheckInput ( long long *  ptr,
int  count 
) const
inlineprivate

Definition at line 535 of file Epetra_MpiComm.h.

void Epetra_MpiComm::CleanupData ( )
private

Definition at line 251 of file Epetra_MpiComm.cpp.

Member Data Documentation

Epetra_MpiCommData* Epetra_MpiComm::MpiCommData_
private

Definition at line 538 of file Epetra_MpiComm.h.


The documentation for this class was generated from the following files: