Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
Epetra_BlockMap Class Reference

Epetra_BlockMap: A class for partitioning block element vectors and matrices. More...

#include <Epetra_BlockMap.h>

Inheritance diagram for Epetra_BlockMap:
Inheritance graph
[legend]
Collaboration diagram for Epetra_BlockMap:
Collaboration graph
[legend]

Public Member Functions

template<>
bool GlobalIndicesIsType () const
 
template<>
bool GlobalIndicesIsType () const
 
Constructors/destructors
 Epetra_BlockMap (int NumGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
 Epetra_BlockMap constructor for a Epetra-defined uniform linear distribution of constant size elements. More...
 
 Epetra_BlockMap (long long NumGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (long long NumGlobalElements, int ElementSize, long long IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (int NumGlobalElements, int NumMyElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
 Epetra_BlockMap constructor for a user-defined linear distribution of constant size elements. More...
 
 Epetra_BlockMap (long long NumGlobalElements, int NumMyElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (long long NumGlobalElements, int NumMyElements, int ElementSize, long long IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (int NumGlobalElements, int NumMyElements, const int *MyGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
 Epetra_BlockMap constructor for a user-defined arbitrary distribution of constant size elements. More...
 
 Epetra_BlockMap (long long NumGlobalElements, int NumMyElements, const long long *MyGlobalElements, int ElementSize, int IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (long long NumGlobalElements, int NumMyElements, const long long *MyGlobalElements, int ElementSize, long long IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (int NumGlobalElements, int NumMyElements, const int *MyGlobalElements, const int *ElementSizeList, int IndexBase, const Epetra_Comm &Comm)
 Epetra_BlockMap constructor for a user-defined arbitrary distribution of variable size elements. More...
 
 Epetra_BlockMap (long long NumGlobalElements, int NumMyElements, const long long *MyGlobalElements, const int *ElementSizeList, int IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (long long NumGlobalElements, int NumMyElements, const long long *MyGlobalElements, const int *ElementSizeList, long long IndexBase, const Epetra_Comm &Comm)
 
 Epetra_BlockMap (long long NumGlobal_Elements, int NumMy_Elements, const long long *myGlobalElements, int ElementSize, int indexBase, const Epetra_Comm &comm, bool UserIsDistributedGlobal, long long UserMinAllGID, long long UserMaxAllGID)
 Epetra_BlockMap constructor for a user-defined arbitrary distribution of constant size elements, where the user provides all the globals. More...
 
 Epetra_BlockMap (long long NumGlobal_Elements, int NumMy_Elements, const long long *myGlobalElements, int ElementSize, long long indexBase, const Epetra_Comm &comm, bool UserIsDistributedGlobal, long long UserMinAllGID, long long UserMaxAllGID)
 
 Epetra_BlockMap (int NumGlobal_Elements, int NumMy_Elements, const int *myGlobalElements, int ElementSize, int indexBase, const Epetra_Comm &comm, bool UserIsDistributedGlobal, int UserMinAllGID, int UserMaxAllGID)
 
 Epetra_BlockMap (const Epetra_BlockMap &map)
 Epetra_BlockMap copy constructor.
 
virtual ~Epetra_BlockMap (void)
 Epetra_BlockMap destructor.
 
Local/Global ID accessor methods
int RemoteIDList (int NumIDs, const int *GIDList, int *PIDList, int *LIDList) const
 Returns the processor IDs and corresponding local index value for a given list of global indices. More...
 
int RemoteIDList (int NumIDs, const long long *GIDList, int *PIDList, int *LIDList) const
 
int RemoteIDList (int NumIDs, const int *GIDList, int *PIDList, int *LIDList, int *SizeList) const
 Returns the processor IDs, corresponding local index value, and element size for a given list of global indices. More...
 
int RemoteIDList (int NumIDs, const long long *GIDList, int *PIDList, int *LIDList, int *SizeList) const
 
int LID (int GID) const
 Returns local ID of global ID, return -1 if not found on this processor.
 
int LID (long long GID) const
 
int GID (int LID) const
 Returns global ID of local ID, return IndexBase-1 if not found on this processor.
 
long long GID64 (int LID) const
 
int FindLocalElementID (int PointID, int &ElementID, int &ElementOffset) const
 Returns the LID of the element that contains the given local PointID, and the Offset of the point in that element.
 
bool MyGID (int GID_in) const
 Returns true if the GID passed in belongs to the calling processor in this map, otherwise returns false.
 
bool MyGID (long long GID_in) const
 
bool MyLID (int lid) const
 Returns true if the LID passed in belongs to the calling processor in this map, otherwise returns false.
 
int MinAllGID () const
 Returns the minimum global ID across the entire map.
 
long long MinAllGID64 () const
 
int MaxAllGID () const
 Returns the maximum global ID across the entire map.
 
long long MaxAllGID64 () const
 
int MinMyGID () const
 Returns the minimum global ID owned by this processor.
 
long long MinMyGID64 () const
 
int MaxMyGID () const
 Returns the maximum global ID owned by this processor.
 
long long MaxMyGID64 () const
 
int MinLID () const
 The minimum local index value on the calling processor.
 
int MaxLID () const
 The maximum local index value on the calling processor.
 
Size and dimension accessor functions
int NumGlobalElements () const
 Number of elements across all processors.
 
long long NumGlobalElements64 () const
 
int NumMyElements () const
 Number of elements on the calling processor.
 
int MyGlobalElements (int *MyGlobalElementList) const
 Puts list of global elements on this processor into the user-provided array.
 
int MyGlobalElements (long long *MyGlobalElementList) const
 
int MyGlobalElementsPtr (int *&MyGlobalElementList) const
 
int MyGlobalElementsPtr (long long *&MyGlobalElementList) const
 
int ElementSize () const
 Returns the size of elements in the map; only valid if map has constant element size.
 
int ElementSize (int LID) const
 Size of element for specified LID.
 
int FirstPointInElement (int LID) const
 Returns the requested entry in the FirstPointInElementList; see FirstPointInElementList() for details. More...
 
int IndexBase () const
 Index base for this map.
 
long long IndexBase64 () const
 
int NumGlobalPoints () const
 Number of global points for this map; equals the sum of all element sizes across all processors.
 
long long NumGlobalPoints64 () const
 
int NumMyPoints () const
 Number of local points for this map; equals the sum of all element sizes on the calling processor.
 
int MinMyElementSize () const
 Minimum element size on the calling processor.
 
int MaxMyElementSize () const
 Maximum element size on the calling processor.
 
int MinElementSize () const
 Minimum element size across all processors.
 
int MaxElementSize () const
 Maximum element size across all processors.
 
Miscellaneous boolean tests
bool UniqueGIDs () const
 Returns true if map GIDs are 1-to-1. More...
 
bool GlobalIndicesInt () const
 Returns true if map create with int NumGlobalElements.
 
bool GlobalIndicesLongLong () const
 Returns true if map create with long long NumGlobalElements.
 
template<typename int_type >
bool GlobalIndicesIsType () const
 
bool GlobalIndicesTypeValid () const
 
bool GlobalIndicesTypeMatch (const Epetra_BlockMap &other) const
 
bool ConstantElementSize () const
 Returns true if map has constant element size.
 
bool SameBlockMapDataAs (const Epetra_BlockMap &Map) const
 
bool SameAs (const Epetra_BlockMap &Map) const
 Returns true if this and Map are identical maps.
 
bool PointSameAs (const Epetra_BlockMap &Map) const
 Returns true if this and Map have identical point-wise structure. More...
 
bool LinearMap () const
 Returns true if the global ID space is contiguously divided (but not necessarily uniformly) across all processors.
 
bool DistributedGlobal () const
 Returns true if map is defined across more than one processor.
 
Array accessor functions
int * MyGlobalElements () const
 Pointer to internal array containing list of global IDs assigned to the calling processor.
 
long long * MyGlobalElements64 () const
 
void MyGlobalElements (const int *&IntGIDs, const long long *&LLGIDs) const
 
void MyGlobalElements (int *&IntGIDs, long long *&LLGIDs)
 
int * FirstPointInElementList () const
 Pointer to internal array containing a mapping between the local elements and the first local point number in each element. More...
 
int * ElementSizeList () const
 List of the element sizes corresponding to the array MyGlobalElements().
 
int * PointToElementList () const
 For each local point, indicates the local element ID that the point belongs to.
 
int ElementSizeList (int *ElementSizeList) const
 Same as ElementSizeList() except it fills the user array that is passed in.
 
int FirstPointInElementList (int *FirstPointInElementList) const
 Same as FirstPointInElementList() except it fills the user array that is passed in.
 
int PointToElementList (int *PointToElementList) const
 Same as PointToElementList() except it fills the user array that is passed in.
 
Miscellaneous
virtual void Print (std::ostream &os) const
 Print object to an output stream.
 
const Epetra_CommComm () const
 Access function for Epetra_Comm communicator.
 
bool IsOneToOne () const
 
Epetra_BlockMapoperator= (const Epetra_BlockMap &map)
 Assignment Operator.
 
Expert Users and Developers Only
int ReferenceCount () const
 Returns the reference count of BlockMapData. More...
 
const Epetra_BlockMapDataDataPtr () const
 Returns a pointer to the BlockMapData instance this BlockMap uses. More...
 
Epetra_BlockMapRemoveEmptyProcesses () const
 Return a new BlockMap with processes with zero elements removed. More...
 
Epetra_BlockMapReplaceCommWithSubset (const Epetra_Comm *Comm) const
 Replace this BlockMap's communicator with a subset communicator. 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 int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method.
 
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...
 

Protected Member Functions

void CleanupData ()
 
- 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
 

Protected Attributes

Epetra_BlockMapDataBlockMapData_
 

Friends

class Epetra_Directory
 
class Epetra_LocalMap
 

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.
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting.
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 

Detailed Description

Epetra_BlockMap: A class for partitioning block element vectors and matrices.

It is often the case that multiple matrix and vector objects have an identical distribution

of elements on a parallel machine. The Epetra_BlockMap class keeps information that describes this distribution for matrices and vectors that have block elements. The definition of an element can vary depending on the situation. For vectors (and multi-vectors), an element is a span of one or more contiguous entries. For matrices, it is a span of one or more matrix rows. More generally, an element in the BlockMap class is an ordered list of points. (NOTE: Points do not have global ID's.) Two additional definitions useful in understanding the BlockMap class follow:

This class has a variety of constructors that can be separated into two categories:

  1. Fixed element size constructors: All map elements have an identical size. This corresponds to a block partitioning of matrices and vectors where the element size is the same for all elements. A common example is multiple degrees of freedom per mesh node in finite element computations where the number of degrees of freedom is the same for all nodes.
  2. Variable element size constructor: Map element sizes may vary and are individually defined via a list of element sizes. This is the most general case and corresponds to a variable block partitioning of the matrices and vectors. A common example is multiple degrees of freedom per mesh node in finite element computations where the number of degrees of freedom varies. This happens, for example, if regions have differing material types or there are chemical reactions in the simulation.

Epetra_BlockMap allows the storage and retrieval of the following information. Depending on the constructor that is used, some of the information is defined by the user and some is determined by the constructor. Once an Epetra_BlockMap is constructed any of the following can be obtained by calling a query function that has the same name as the attribute, e.g. to get the value of NumGlobalElements, you can call a function NumGlobalElements(). For attributes that are lists, the query functions return the list values in a user allocated array.

In addition to the information above that is passed in to or created by the Epetra_BlockMap constructor, the following attributes are computed and available via query to the user using the same scheme as above, e.g., use NumGlobalPoints() to get the value of NumGlobalPoints.

The following functions allow boolean tests for certain properties.

Warning
A Epetra_Comm object is required for all Epetra_BlockMap constructors.

error handling

Most methods in Epetra_BlockMap return an integer error code. If the error code is 0, then no error occurred. If > 0 then a warning error occurred. If < 0 then a fatal error occurred.

Epetra_BlockMap constructors will throw an exception of an error occurrs. These exceptions will alway be negative integer values as follows:

  1. -1 NumGlobalElements < -1. Should be >= -1 (Should be >= 0 for first BlockMap constructor).
  2. -2 NumMyElements < 0. Should be >= 0.
  3. -3 ElementSize <= 0. Should be > 0.
  4. -4 Invalid NumGlobalElements. Should equal sum of MyGlobalElements, or set to -1 to compute automatically.
  5. -5 Minimum global element index is less than index base.
  6. -99 Internal Epetra_BlockMap error. Contact developer.

For robust code, Epetra_BlockMap constructor calls should be caught using the try {...} catch {...} mechanism. For example:

  try {

    Epetra_BlockMap * map = new Epetra_BlockMap(NumGlobalElements, ElementSize, IndexBase, Comm);
  }
  catch (int Error) {
    if (Error==-1) { // handle error }
    if (Error==-2) ...
Note
{ In the current implementation, Epetra_BlockMap is the base class for: }

Constructor & Destructor Documentation

Epetra_BlockMap::Epetra_BlockMap ( int  NumGlobalElements,
int  ElementSize,
int  IndexBase,
const Epetra_Comm Comm 
)

Epetra_BlockMap constructor for a Epetra-defined uniform linear distribution of constant size elements.

Creates a map that distributes NumGlobalElements elements evenly across all processors in the Epetra_Comm communicator. If NumGlobalElements does not divide exactly into the number of processors, the first processors in the communicator get one extra element until the remainder is gone.

The elements are defined to have a constant fixed size specified by ElementSize.

Parameters
InNumGlobalElements - Number of elements to distribute.
InElementSize - Number of points or vector entries per element.
InIndexBase - Minimum index value used for arrays that use this map. Typically 0 for C/C++ and 1 for Fortran.
InComm - Epetra_Comm communicator containing information on the number of processors.
Returns
Pointer to a Epetra_BlockMap object.
Epetra_BlockMap::Epetra_BlockMap ( int  NumGlobalElements,
int  NumMyElements,
int  ElementSize,
int  IndexBase,
const Epetra_Comm Comm 
)

Epetra_BlockMap constructor for a user-defined linear distribution of constant size elements.

Creates a map that puts NumMyElements on the calling processor. If NumGlobalElements=-1, the number of global elements will be the computed sum of NumMyElements across all processors in the Epetra_Comm communicator.

The elements are defined to have a constant fixed size specified by ElementSize.

Parameters
InNumGlobalElements - Number of elements to distribute. Must be either -1 or equal to the computed sum of NumMyElements across all processors in the Epetra_Comm communicator.
InNumMyElements - Number of elements owned by the calling processor.
InElementSize - Number of points or vector entries per element.
InIndexBase - Minimum index value used for arrays that use this map. Typically 0 for C/C++ and 1 for Fortran.
InComm - Epetra_Comm communicator containing information on the number of processors.
Returns
Pointer to a Epetra_BlockMap object.
Epetra_BlockMap::Epetra_BlockMap ( int  NumGlobalElements,
int  NumMyElements,
const int *  MyGlobalElements,
int  ElementSize,
int  IndexBase,
const Epetra_Comm Comm 
)

Epetra_BlockMap constructor for a user-defined arbitrary distribution of constant size elements.

Creates a map that puts NumMyElements on the calling processor. The indices of the elements are determined from the list MyGlobalElements. If NumGlobalElements=-1, the number of global elements will be the computed sum of NumMyElements across all processors in the Epetra_Comm communicator.

The elements are defined to have a constant fixed size specified by ElementSize.

Parameters
InNumGlobalElements - Number of elements to distribute. Must be either -1 or equal to the computed sum of NumMyElements across all processors in the Epetra_Comm communicator.
InNumMyElements - Number of elements owned by the calling processor.
InMyGlobalElements - Integer array of length NumMyElements. The ith entry contains the global index value of the ith element on this processor. Index values are not required to be contiguous on a processor, or to be within the range of 0 to NumGlobalElements. As long as the index values are consistently defined and used, any set of NumGlobalElements distinct integer values is acceptable.
InElementSize - Number of points or vector entries per element.
InIndexBase - Minimum index value used for arrays that use this map. Typically 0 for C/C++ and 1 for Fortran.
InComm - Epetra_Comm communicator containing information on the number of processors.
Returns
Pointer to a Epetra_BlockMap object.
Epetra_BlockMap::Epetra_BlockMap ( int  NumGlobalElements,
int  NumMyElements,
const int *  MyGlobalElements,
const int *  ElementSizeList,
int  IndexBase,
const Epetra_Comm Comm 
)

Epetra_BlockMap constructor for a user-defined arbitrary distribution of variable size elements.

Creates a map that puts NumMyElements on the calling processor. If NumGlobalElements=-1, the number of global elements will be the computed sum of NumMyElements across all processors in the Epetra_Comm communicator.

The elements are defined to have a variable size defined by ElementSizeList.

Parameters
InNumGlobalElements - Number of elements to distribute. Must be either -1 or equal to the computed sum of NumMyElements across all processors in the Epetra_Comm communicator.
InNumMyElements - Number of elements owned by the calling processor.
InMyGlobalElements - Integer array of length NumMyElements. The ith entry contains the global index value of the ith element on this processor. Index values are not required to be contiguous on a processor, or to be within the range of 0 to NumGlobalElements. As long as the index values are consistently defined and used, any set of NumGlobalElements distinct integer values is acceptable.
InElementSizeList - A list of the element sizes for elements owned by the calling processor. The ith entry contains the element size of the ith element on this processor.
InIndexBase - Minimum index value used for arrays that use this map. Typically 0 for C/C++ and 1 for Fortran.
InComm - Epetra_Comm communicator containing information on the number of processors.
Returns
Pointer to a Epetra_BlockMap object.
Epetra_BlockMap::Epetra_BlockMap ( long long  NumGlobal_Elements,
int  NumMy_Elements,
const long long *  myGlobalElements,
int  ElementSize,
int  indexBase,
const Epetra_Comm comm,
bool  UserIsDistributedGlobal,
long long  UserMinAllGID,
long long  UserMaxAllGID 
)

Epetra_BlockMap constructor for a user-defined arbitrary distribution of constant size elements, where the user provides all the globals.

Warning
This method is intended for expert developer use only, and should never be called by user code.

Member Function Documentation

const Epetra_BlockMapData* Epetra_BlockMap::DataPtr ( ) const
inline

Returns a pointer to the BlockMapData instance this BlockMap uses.

(Intended for developer use only for testing purposes.)

int Epetra_BlockMap::FirstPointInElement ( int  LID) const

Returns the requested entry in the FirstPointInElementList; see FirstPointInElementList() for details.

This function provides similar functionality to FirstPointInElementList(), but for simple maps may avoid the explicit construction of the FirstPointInElementList array. Returns -1 if LID is out-of-range.

int* Epetra_BlockMap::FirstPointInElementList ( ) const

Pointer to internal array containing a mapping between the local elements and the first local point number in each element.

This array is a scan sum of the ElementSizeList such that the ith entry in FirstPointInElementList is the sum of the first i-1 entries of ElementSizeList().

bool Epetra_BlockMap::PointSameAs ( const Epetra_BlockMap Map) const

Returns true if this and Map have identical point-wise structure.

If both maps have the same number of global points and the same point distribution across processors then this method returns true.

int Epetra_BlockMap::ReferenceCount ( ) const
inline

Returns the reference count of BlockMapData.

(Intended for testing purposes.)

int Epetra_BlockMap::RemoteIDList ( int  NumIDs,
const int *  GIDList,
int *  PIDList,
int *  LIDList 
) const
inline

Returns the processor IDs and corresponding local index value for a given list of global indices.

For each element (GID) of a given list of global element numbers (stored in GIDList) of length NumIDs, this function returns (in PIDList) the ID (rank) of the processor that owns the GID for this map and returns the local index (in LIDList) of the GID on that processor.

If a GID is present on more than one processor, the lowest rank processor ID is used, as is the LID for that processor. If a GID is not present on any processor, the corresponding PID will return as -1.

int Epetra_BlockMap::RemoteIDList ( int  NumIDs,
const int *  GIDList,
int *  PIDList,
int *  LIDList,
int *  SizeList 
) const

Returns the processor IDs, corresponding local index value, and element size for a given list of global indices.

For each element (GID) of a given a list of global element numbers (stored in GIDList) of length NumIDs, this function returns (in PIDList) the with processor that owns the GID for this map and returns the local index (in LIDList) of the GID on that processor. Finally it returns the element sizes in SizeList.

Epetra_BlockMap* Epetra_BlockMap::RemoveEmptyProcesses ( ) const

Return a new BlockMap with processes with zero elements removed.

Warning
This method is only for expert users. Understanding how to use this method correctly requires some familiarity with semantics of MPI communicators.
We make no promises of backwards compatibility for this method. It may go away or change at any time.

This method first computes a new communicator, which contains only those processes in this Map's communicator (the "original communicator") that have a nonzero number of elements in this BlockMap (the "original BlockMap"). It then returns a new BlockMap distributed over the new communicator. The new BlockMap represents the same distribution as the original BlockMap, except that processes containing zero elements are not included in the new BlockMap or its communicator. On processes not included in the new BlockMap or communicator, this method returns NULL.

The returned BlockMap always has a distinct communicator from this BlockMap's original communicator. The new communicator contains a subset of processes from the original communicator. Even if the number of processes in the new communicator equals the number of processes in the original communicator, the new communicator is distinct. (In an MPI implementation, the new communicator is created using MPI_Comm_split.)

This method must be called collectively on the original communicator. It leaves the original Map and communicator unchanged.

This method was intended for applications such as algebraic multigrid or other multilevel preconditioners. Construction of each level of the multilevel preconditioner typically requires constructing sparse matrices, which in turn requires all-reduces over all participating processes at that level. Matrix sizes at successively coarser levels shrink geometrically. At the coarsest levels, some processes might be left with zero rows of the matrix, or the multigrid implementation might "rebalance" (redistribute the matrix) and intentionally leave some processes with zero rows. Removing processes with zero rows makes the all-reduces and other communication operations cheaper.

Epetra_BlockMap* Epetra_BlockMap::ReplaceCommWithSubset ( const Epetra_Comm Comm) const

Replace this BlockMap's communicator with a subset communicator.

Warning
This method is only for expert users. Understanding how to use this method correctly requires some familiarity with semantics of MPI communicators.
We make no promises of backwards compatibility for this method. It may go away or change at any time.
Precondition
The input communicator's processes are a subset of this BlockMap's current communicator's processes.
On processes which are not included in the input communicator, the input communicator is null.

This method must be called collectively on the original communicator. It leaves the original BlockMap and communicator unchanged.

Note
This method differs from removeEmptyProcesses(), in that it does not assume that excluded processes have zero entries. For example, one might wish to remove empty processes from the row BlockMap of a CrsGraph using removeEmptyProcesses(), and then apply the resulting subset communicator to the column, domain, and range Maps of the same graph. For the latter three Maps, one would in general use this method instead of removeEmptyProcesses(), giving the new row BlockMap's communicator to this method.
bool Epetra_BlockMap::SameBlockMapDataAs ( const Epetra_BlockMap Map) const

Returns true if maps share same block map data underneath. This is a very cheap check whether two maps are identical.

bool Epetra_BlockMap::UniqueGIDs ( ) const
inline

Returns true if map GIDs are 1-to-1.

Certain operations involving Epetra_BlockMap and Epetra_Map objects are well-defined only if the map GIDs are uniquely present in the map. In other words, if a GID occurs in the map, it occurs only once on a single processor and nowhere else. This boolean test returns true if this property is true, otherwise it returns false.


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