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

MPI implementation of Epetra_Distributor. More...

#include <Epetra_MpiDistributor.h>

Inheritance diagram for Epetra_MpiDistributor:
Inheritance graph
[legend]

Private Member Functions

int CreateSendStructures_ (int my_proc, int nprocs, const int &NumExportIDs, const int *ExportPIDs)
 
int CreateRecvStructures_ (const int &NumRemoteIDs, const int *RemotePIDs)
 
int ComputeRecvs_ (int my_proc, int nprocs)
 
template<typename id_type >
int ComputeSends_ (int num_imports, const id_type *&import_ids, const int *&import_procs, int &num_exports, id_type *&export_ids, int *&export_procs, int my_proc)
 
int Resize_ (int *sizes)
 
int Sort_ints_ (int *vals, int *other, int nvals)
 
Epetra_MpiDistributoroperator= (const Epetra_MpiDistributor &src)
 
void CreateReverseDistributor ()
 

Private Attributes

int * lengths_to_
 
int * procs_to_
 
int * indices_to_
 
int size_indices_to_
 
int * lengths_from_
 
int * procs_from_
 
int * indices_from_
 
int size_indices_from_
 
bool resized_
 
int * sizes_
 
int * sizes_to_
 
int * starts_to_
 
int * starts_to_ptr_
 
int * indices_to_ptr_
 
int * sizes_from_
 
int * starts_from_
 
int * starts_from_ptr_
 
int * indices_from_ptr_
 
int nrecvs_
 
int nsends_
 
int nexports_
 
int self_msg_
 
int max_send_length_
 
int total_recv_length_
 
int tag_
 
const Epetra_MpiCommepComm_
 
const MPI_Comm comm_
 
MPI_Request * request_
 
MPI_Status * status_
 
bool no_delete_
 
char * send_array_
 
int send_array_size_
 
Epetra_MpiDistributorcomm_plan_reverse_
 
int lastRoundBytesSend_
 
int lastRoundBytesRecv_
 

Constructors/Destructor

 Epetra_MpiDistributor (const Epetra_MpiComm &Comm)
 Default constructor. More...
 
 Epetra_MpiDistributor (const Epetra_MpiDistributor &Distributor)
 Copy constructor. More...
 
Epetra_DistributorClone ()
 Clone method. More...
 
Epetra_DistributorReverseClone ()
 Create and extract the reverse version of the distributor. More...
 
virtual ~Epetra_MpiDistributor ()
 Destructor (declared virtual for memory safety). More...
 

Gather/Scatter Constructors

int CreateFromSends (const int &NumExportIDs, const int *ExportPIDs, bool Deterministic, int &NumRemoteIDs)
 Create a communication plan from send list. More...
 
int CreateFromRecvs (const int &NumRemoteIDs, const int *RemoteGIDs, const int *RemotePIDs, bool Deterministic, int &NumExportIDs, int *&ExportGIDs, int *&ExportPIDs)
 Create a communication plan from receive list. More...
 
int CreateFromRecvs (const int &NumRemoteIDs, const long long *RemoteGIDs, const int *RemotePIDs, bool Deterministic, int &NumExportIDs, long long *&ExportGIDs, int *&ExportPIDs)
 
int CreateFromSendsAndRecvs (const int &NumExportIDs, const int *ExportPIDs, const int &NumRemoteIDs, const int *RemoteGIDs, const int *RemotePIDs, bool Deterministic)
 Create a communication plan from send list and a recv list. More...
 
int CreateFromSendsAndRecvs (const int &NumExportIDs, const int *ExportPIDs, const int &NumRemoteIDs, const long long *RemoteGIDs, const int *RemotePIDs, bool Deterministic)
 

Execute Gather/Scatter Operations

int Do (char *export_objs, int obj_size, int &len_import_objs, char *&import_objs)
 Execute plan on buffer of export objects in a single step. More...
 
int DoReverse (char *export_objs, int obj_size, int &len_import_objs, char *&import_objs)
 Execute reverse of plan on buffer of export objects in a single step. More...
 
int DoPosts (char *export_objs, int obj_size, int &len_import_objs, char *&import_objs)
 Post buffer of export objects (can do other local work before executing Waits) More...
 
int DoWaits ()
 Wait on a set of posts. More...
 
int DoReversePosts (char *export_objs, int obj_size, int &len_import_objs, char *&import_objs)
 Do reverse post of buffer of export objects (can do other local work before executing Waits) More...
 
int DoReverseWaits ()
 Wait on a reverse set of posts. More...
 

Execute Gather/Scatter Operations (Non-constant size objects)

int Do (char *export_objs, int obj_size, int *&sizes, int &len_import_objs, char *&import_objs)
 Execute plan on buffer of export objects in a single step (object size may vary) More...
 
int DoReverse (char *export_objs, int obj_size, int *&sizes, int &len_import_objs, char *&import_objs)
 Execute reverse of plan on buffer of export objects in a single step (object size may vary) More...
 
int DoPosts (char *export_objs, int obj_size, int *&sizes, int &len_import_objs, char *&import_objs)
 Post buffer of export objects (can do other local work before executing Waits) More...
 
int DoReversePosts (char *export_objs, int obj_size, int *&sizes, int &len_import_objs, char *&import_objs)
 Do reverse post of buffer of export objects (can do other local work before executing Waits) More...
 

Attribute Accessor Methods

int NumReceives () const
 The number of procs from which we will receive data. More...
 
int NumSends () const
 The number of procs to which we will send data. More...
 
int MaxSendLength () const
 Maximum number of values that this proc is sending to another single proc. More...
 
int TotalReceiveLength () const
 Total number of values that this proc is receiving from other procs. More...
 
const int * ProcsFrom () const
 A list of procs sending values to this proc. More...
 
const int * ProcsTo () const
 A list of procs to which this proc is sending values. More...
 
const int * LengthsFrom () const
 Number of values we're receiving from each proc. More...
 
const int * LengthsTo () const
 Number of values we're sending to each proc. More...
 
void GetLastDoStatistics (int &bytes_sent, int &bytes_recvd) const
 Information on the last call to Do/DoReverse. More...
 

Print object to an output stream

void Print (std::ostream &os) const
 

Additional Inherited Members

- 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_Distributor
virtual ~Epetra_Distributor ()
 Epetra_Distributor Destructor. More...
 
- 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

MPI implementation of Epetra_Distributor.

This class is an MPI implementation of Epetra_Distributor. It encapsulates the general information and services needed for other Epetra classes to perform gather/scatter operations on a parallel computer. An Epetra_MpiDistributor instance is actually produced by calling a method in the Epetra_MpiComm class.

Definition at line 59 of file Epetra_MpiDistributor.h.

Constructor & Destructor Documentation

Epetra_MpiDistributor::Epetra_MpiDistributor ( const Epetra_MpiComm Comm)

Default constructor.

Definition at line 49 of file Epetra_MpiDistributor.cpp.

Epetra_MpiDistributor::Epetra_MpiDistributor ( const Epetra_MpiDistributor Distributor)

Copy constructor.

Definition at line 88 of file Epetra_MpiDistributor.cpp.

Epetra_MpiDistributor::~Epetra_MpiDistributor ( )
virtual

Destructor (declared virtual for memory safety).

Definition at line 189 of file Epetra_MpiDistributor.cpp.

Member Function Documentation

Epetra_Distributor* Epetra_MpiDistributor::Clone ( )
inlinevirtual

Clone method.

Implements Epetra_Distributor.

Definition at line 73 of file Epetra_MpiDistributor.h.

Epetra_Distributor * Epetra_MpiDistributor::ReverseClone ( )
virtual

Create and extract the reverse version of the distributor.

This is not a const method since a reverse distributor might need to be created. This works like Clone, returning a new object the user must deallocate.

Implements Epetra_Distributor.

Definition at line 1589 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateFromSends ( const int &  NumExportIDs,
const int *  ExportPIDs,
bool  Deterministic,
int &  NumRemoteIDs 
)
virtual

Create a communication plan from send list.

Given a list of process IDs to which to send the given number of data IDs, construct a communication plan for efficiently scattering data to these processes.

Returns
The number of data IDs being sent to me.
Parameters
NumExportIDs[in] Number of data IDs that need to be sent from the calling process.
ExportPIDs[in] List of process IDs that will get the exported data IDs.
Deterministic[in] Currently has no effect.
NumRemoteIDs[out] Number of data IDs the calling process will be receiving.

Implements Epetra_Distributor.

Definition at line 221 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateFromRecvs ( const int &  NumRemoteIDs,
const int *  RemoteGIDs,
const int *  RemotePIDs,
bool  Deterministic,
int &  NumExportIDs,
int *&  ExportGIDs,
int *&  ExportPIDs 
)
virtual

Create a communication plan from receive list.

Given a list of remote data IDs and corresponding process IDs from which to receive data, construct a communication plan for efficiently scattering data to these processes.

Returns
The number and list of data IDs being sent by me.
Parameters
NumRemoteIDs[in] Number of data IDs the calling process will be receiving.
RemoteGIDs[in] List of data IDs that the calling process wants to receive.
RemotePIDs[in] List of IDs of the processes that will send the remote data IDs to the calling process.
Deterministic[in] Currently has no effect.
NumExportIDs[out] Number of data IDs that need to be sent from the calling process.
ExportGIDs[out] List of data IDs that the calling process will send out.
ExportPIDs[out] List of IDs of the processes that will receive the data IDs sent by the calling process.
Note
This method allocates the output arrays using new. The caller is responsible for deallocating them after use.

Implements Epetra_Distributor.

Definition at line 252 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateFromRecvs ( const int &  NumRemoteIDs,
const long long *  RemoteGIDs,
const int *  RemotePIDs,
bool  Deterministic,
int &  NumExportIDs,
long long *&  ExportGIDs,
int *&  ExportPIDs 
)
virtual

Implements Epetra_Distributor.

Definition at line 282 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateFromSendsAndRecvs ( const int &  NumExportIDs,
const int *  ExportPIDs,
const int &  NumRemoteIDs,
const int *  RemoteGIDs,
const int *  RemotePIDs,
bool  Deterministic 
)

Create a communication plan from send list and a recv list.

Given a list of process IDs to which to send the given number of data IDs, and a list of remote data IDs and corresponding process IDs from which to receive data, construct a communication plan for efficiently scattering data to these processes.

Needless to say, knowing both of these at the same time is a pretty rare occurance. But if it happens, this routine will avoid a lot of communication that CreateFromSends or CreateFromRecvs would have to do.

Returns
zero if this worked.
Parameters
NumExportIDs[in] Number of data IDs that need to be sent from the calling process.
ExportPIDs[in] List of process IDs that will get the exported data IDs.
NumRemoteIDs[in] Number of data IDs the calling process will be receiving.
RemoteGIDs[in] List of data IDs that the calling process wants to receive.
RemotePIDs[in] List of IDs of the processes that will send the remote data IDs to the calling process.
Deterministic[in] Currently has no effect.

Definition at line 313 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateFromSendsAndRecvs ( const int &  NumExportIDs,
const int *  ExportPIDs,
const int &  NumRemoteIDs,
const long long *  RemoteGIDs,
const int *  RemotePIDs,
bool  Deterministic 
)

Definition at line 338 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::Do ( char *  export_objs,
int  obj_size,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Execute plan on buffer of export objects in a single step.

Implements Epetra_Distributor.

Definition at line 743 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoReverse ( char *  export_objs,
int  obj_size,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Execute reverse of plan on buffer of export objects in a single step.

Implements Epetra_Distributor.

Definition at line 761 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoPosts ( char *  export_objs,
int  obj_size,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Post buffer of export objects (can do other local work before executing Waits)

Implements Epetra_Distributor.

Definition at line 781 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoWaits ( )
virtual

Wait on a set of posts.

Implements Epetra_Distributor.

Definition at line 941 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoReversePosts ( char *  export_objs,
int  obj_size,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Do reverse post of buffer of export objects (can do other local work before executing Waits)

Implements Epetra_Distributor.

Definition at line 949 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoReverseWaits ( )
virtual

Wait on a reverse set of posts.

Implements Epetra_Distributor.

Definition at line 967 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::Do ( char *  export_objs,
int  obj_size,
int *&  sizes,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Execute plan on buffer of export objects in a single step (object size may vary)

Implements Epetra_Distributor.

Definition at line 1166 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoReverse ( char *  export_objs,
int  obj_size,
int *&  sizes,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Execute reverse of plan on buffer of export objects in a single step (object size may vary)

Implements Epetra_Distributor.

Definition at line 1180 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoPosts ( char *  export_objs,
int  obj_size,
int *&  sizes,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Post buffer of export objects (can do other local work before executing Waits)

Implements Epetra_Distributor.

Definition at line 1194 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::DoReversePosts ( char *  export_objs,
int  obj_size,
int *&  sizes,
int &  len_import_objs,
char *&  import_objs 
)
virtual

Do reverse post of buffer of export objects (can do other local work before executing Waits)

Implements Epetra_Distributor.

Definition at line 1330 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::NumReceives ( ) const
inline

The number of procs from which we will receive data.

Definition at line 268 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::NumSends ( ) const
inline

The number of procs to which we will send data.

Definition at line 271 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::MaxSendLength ( ) const
inline

Maximum number of values that this proc is sending to another single proc.

Definition at line 274 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::TotalReceiveLength ( ) const
inline

Total number of values that this proc is receiving from other procs.

Definition at line 277 of file Epetra_MpiDistributor.h.

const int* Epetra_MpiDistributor::ProcsFrom ( ) const
inline

A list of procs sending values to this proc.

Definition at line 280 of file Epetra_MpiDistributor.h.

const int* Epetra_MpiDistributor::ProcsTo ( ) const
inline

A list of procs to which this proc is sending values.

Definition at line 283 of file Epetra_MpiDistributor.h.

const int* Epetra_MpiDistributor::LengthsFrom ( ) const
inline

Number of values we're receiving from each proc.

We will receive LengthsFrom[i] values from proc ProcsFrom[i].

Definition at line 287 of file Epetra_MpiDistributor.h.

const int* Epetra_MpiDistributor::LengthsTo ( ) const
inline

Number of values we're sending to each proc.

We will send LengthsTo[i] values to procs ProcsTo[i].

Definition at line 291 of file Epetra_MpiDistributor.h.

void Epetra_MpiDistributor::GetLastDoStatistics ( int &  bytes_sent,
int &  bytes_recvd 
) const
inline

Information on the last call to Do/DoReverse.

Returns the amount of data sent & recv'd on this processor in bytes

Definition at line 296 of file Epetra_MpiDistributor.h.

void Epetra_MpiDistributor::Print ( std::ostream &  os) const
virtual

Implements Epetra_Distributor.

Definition at line 1349 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateSendStructures_ ( int  my_proc,
int  nprocs,
const int &  NumExportIDs,
const int *  ExportPIDs 
)
private

Definition at line 369 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::CreateRecvStructures_ ( const int &  NumRemoteIDs,
const int *  RemotePIDs 
)
private

Definition at line 504 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::ComputeRecvs_ ( int  my_proc,
int  nprocs 
)
private

Definition at line 552 of file Epetra_MpiDistributor.cpp.

template<typename id_type >
int Epetra_MpiDistributor::ComputeSends_ ( int  num_imports,
const id_type *&  import_ids,
const int *&  import_procs,
int &  num_exports,
id_type *&  export_ids,
int *&  export_procs,
int  my_proc 
)
private

Definition at line 678 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::Resize_ ( int *  sizes)
private

Definition at line 980 of file Epetra_MpiDistributor.cpp.

int Epetra_MpiDistributor::Sort_ints_ ( int *  vals,
int *  other,
int  nvals 
)
private

Definition at line 1476 of file Epetra_MpiDistributor.cpp.

Epetra_MpiDistributor & Epetra_MpiDistributor::operator= ( const Epetra_MpiDistributor src)
private

Definition at line 1531 of file Epetra_MpiDistributor.cpp.

void Epetra_MpiDistributor::CreateReverseDistributor ( )
private

Definition at line 1544 of file Epetra_MpiDistributor.cpp.

Member Data Documentation

int* Epetra_MpiDistributor::lengths_to_
private

Definition at line 340 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::procs_to_
private

Definition at line 341 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::indices_to_
private

Definition at line 342 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::size_indices_to_
private

Definition at line 343 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::lengths_from_
private

Definition at line 345 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::procs_from_
private

Definition at line 346 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::indices_from_
private

Definition at line 347 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::size_indices_from_
private

Definition at line 348 of file Epetra_MpiDistributor.h.

bool Epetra_MpiDistributor::resized_
private

Definition at line 350 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::sizes_
private

Definition at line 351 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::sizes_to_
private

Definition at line 353 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::starts_to_
private

Definition at line 354 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::starts_to_ptr_
private

Definition at line 355 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::indices_to_ptr_
private

Definition at line 356 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::sizes_from_
private

Definition at line 358 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::starts_from_
private

Definition at line 359 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::starts_from_ptr_
private

Definition at line 360 of file Epetra_MpiDistributor.h.

int* Epetra_MpiDistributor::indices_from_ptr_
private

Definition at line 361 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::nrecvs_
private

Definition at line 363 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::nsends_
private

Definition at line 364 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::nexports_
private

Definition at line 365 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::self_msg_
private

Definition at line 367 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::max_send_length_
private

Definition at line 369 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::total_recv_length_
private

Definition at line 370 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::tag_
private

Definition at line 372 of file Epetra_MpiDistributor.h.

const Epetra_MpiComm* Epetra_MpiDistributor::epComm_
private

Definition at line 374 of file Epetra_MpiDistributor.h.

const MPI_Comm Epetra_MpiDistributor::comm_
private

Definition at line 375 of file Epetra_MpiDistributor.h.

MPI_Request* Epetra_MpiDistributor::request_
private

Definition at line 377 of file Epetra_MpiDistributor.h.

MPI_Status* Epetra_MpiDistributor::status_
private

Definition at line 378 of file Epetra_MpiDistributor.h.

bool Epetra_MpiDistributor::no_delete_
private

Definition at line 380 of file Epetra_MpiDistributor.h.

char* Epetra_MpiDistributor::send_array_
private

Definition at line 382 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::send_array_size_
private

Definition at line 383 of file Epetra_MpiDistributor.h.

Epetra_MpiDistributor* Epetra_MpiDistributor::comm_plan_reverse_
private

Definition at line 385 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::lastRoundBytesSend_
private

Definition at line 387 of file Epetra_MpiDistributor.h.

int Epetra_MpiDistributor::lastRoundBytesRecv_
private

Definition at line 388 of file Epetra_MpiDistributor.h.


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