Domi
Multi-dimensional, distributed data structures
 All Classes Files Functions Variables Typedefs Friends
List of all members
Domi::MDMap Class Reference

Multi-dimensional map. More...

#include <Domi_MDMap.hpp>

Public Member Functions

Constructors and destructor
 MDMap (const Teuchos::RCP< const MDComm > mdComm, const Teuchos::ArrayView< const dim_type > &dimensions, const Teuchos::ArrayView< const int > &commPad=Teuchos::ArrayView< const int >(), const Teuchos::ArrayView< const int > &bndryPad=Teuchos::ArrayView< const int >(), const Teuchos::ArrayView< const int > &replicatedBoundary=Teuchos::ArrayView< const int >(), const Layout layout=DEFAULT_ORDER)
 Main constructor. More...
 
 MDMap (Teuchos::ParameterList &plist)
 Constructor with ParameterList. More...
 
 MDMap (const Teuchos::RCP< const Teuchos::Comm< int > > teuchosComm, Teuchos::ParameterList &plist)
 Constructor with Teuchos::Comm and ParameterList. More...
 
 MDMap (const Teuchos::RCP< const MDComm > mdComm, Teuchos::ParameterList &plist)
 Constructor with MDComm and ParameterList. More...
 
 MDMap (const Teuchos::RCP< const MDComm > mdComm, const Teuchos::ArrayView< Slice > &myGlobalBounds, const Teuchos::ArrayView< padding_type > &padding=Teuchos::ArrayView< padding_type >(), const Teuchos::ArrayView< const int > &replicatedBoundary=Teuchos::ArrayView< const int >(), const Layout layout=DEFAULT_ORDER)
 Constructor with global bounds for this processor. More...
 
 MDMap (const MDMap &source)
 Copy constructor. More...
 
 MDMap (const MDMap &parent, int axis, dim_type index)
 Parent/single global ordinal sub-map constructor. More...
 
 MDMap (const MDMap &parent, int axis, const Slice &slice, int bndryPad=0)
 Parent/single slice sub-map constructor. More...
 
 MDMap (const MDMap &parent, const Teuchos::ArrayView< Slice > &slices, const Teuchos::ArrayView< int > &bndryPad=Teuchos::ArrayView< int >())
 Parent/array of slices sub-map constructor. More...
 
 ~MDMap ()
 Destructor.
 
MDMap standard operators
MDMapoperator= (const MDMap &source)
 Assignment operator. More...
 
MDComm accessor and pass-through methods
Teuchos::RCP< const
Teuchos::Comm< int > > 
getTeuchosComm () const
 Get the Teuchos communicator. More...
 
Teuchos::RCP< const MDCommgetMDComm () const
 Access the underlying MDComm object. More...
 
bool onSubcommunicator () const
 Query whether this processor is on the sub-communicator. More...
 
int numDims () const
 Get the number of dimensions. More...
 
Teuchos::Array< int > getCommDims () const
 Get the communicator sizes along each axis. More...
 
int getCommDim (int axis) const
 Get the communicator size along the given axis. More...
 
bool isPeriodic (int axis) const
 Return the periodic flag for the given axis. More...
 
int getCommIndex (int axis) const
 Get the axis rank of this processor. More...
 
int getLowerNeighbor (int axis) const
 Get the rank of the lower neighbor. More...
 
int getUpperNeighbor (int axis) const
 Get the rank of the upper neighbor. More...
 
Dimensions and looping bounds
Teuchos::Array< dim_type > getGlobalDims () const
 Get an array of the the global dimensions, including boundary padding.
 
dim_type getGlobalDim (int axis, bool withBndryPad=false) const
 Get the global dimension along the specified axis. More...
 
Slice getGlobalBounds (int axis, bool withBndryPad=false) const
 Get the bounds of the global problem. More...
 
Slice getGlobalRankBounds (int axis, bool withBndryPad=false) const
 Get the global loop bounds along the specified axis. More...
 
Teuchos::Array< dim_type > getLocalDims () const
 Get an array of the local dimensions, including padding.
 
dim_type getLocalDim (int axis, bool withPad=false) const
 Get the local dimension along the specified axis. More...
 
Teuchos::ArrayView< const SlicegetLocalBounds () const
 Get the local loop bounds along every axis. More...
 
Slice getLocalBounds (int axis, bool withPad=false) const
 Get the local loop bounds along the specified axis. More...
 
Slice getLocalInteriorBounds (int axis) const
 Get the local interior loop bounds along the specified axis. More...
 
Storage order, communication and boundary padding
bool hasPadding () const
 Return true if there is any padding stored locally. More...
 
int getLowerPadSize (int axis) const
 Get the size of the lower padding along the given axis. More...
 
int getUpperPadSize (int axis) const
 Get the size of the upper padding along the given axis. More...
 
int getCommPadSize (int axis) const
 Get the communication padding size along the given axis. More...
 
int getLowerBndryPad (int axis) const
 Get the size of the lower boundary padding along the given axis. More...
 
int getUpperBndryPad (int axis) const
 Get the size of the upper boundary padding along the given axis. More...
 
Teuchos::Array< int > getBndryPadSizes () const
 Get the boundary padding sizes along each axis. More...
 
int getBndryPadSize (int axis) const
 Get the boundary padding size along the given axis. More...
 
bool isPad (const Teuchos::ArrayView< dim_type > &index) const
 
bool isCommPad (const Teuchos::ArrayView< dim_type > &index) const
 
bool isBndryPad (const Teuchos::ArrayView< dim_type > &index) const
 
bool isReplicatedBoundary (int axis) const
 Return whether the given axis supports a replicated boundary.
 
Layout getLayout () const
 Get the storage order.
 
Conversions to other Maps
Teuchos::ArrayView
< Teuchos::RCP< const
Domi::MDMap > > 
getAxisMaps () const
 Return an array of axis maps. More...
 
Teuchos::RCP< const Domi::MDMapgetAxisMap (int axis) const
 Return an axis map for the given axis. More...
 
Teuchos::RCP< const MDMapgetAugmentedMDMap (const dim_type leadingDim, const dim_type trailingDim=0) const
 Return an RCP to a new MDMap that is a simple augmentation of this MDMap. More...
 
Global/local index and ID conversion methods
Teuchos::Array< dim_type > getGlobalIndex (size_type globalID) const
 Convert a global ID to a global index. More...
 
Teuchos::Array< dim_type > getLocalIndex (size_type localID) const
 Convert a local ID to a local index. More...
 
size_type getGlobalID (size_type localID) const
 Convert a local ID to a global ID. More...
 
size_type getGlobalID (const Teuchos::ArrayView< dim_type > &globalIndex) const
 convert a global index to a global ID More...
 
size_type getLocalID (size_type globalID) const
 Convert a global ID to a local ID. More...
 
size_type getLocalID (const Teuchos::ArrayView< dim_type > &localIndex) const
 Convert a local index to a local ID. More...
 
Boolean comparison methods
bool isCompatible (const MDMap &mdMap) const
 True if two MDMaps are "compatible". More...
 
bool isSameAs (const MDMap &mdMap, const int verbose=0) const
 True if two MDMaps are "identical". More...
 
bool isContiguous () const
 True if there are no stride gaps on any processor. More...
 

Detailed Description

Multi-dimensional map.

The MDMap class is intended to perform the functions of Epetra or Tpetra maps, except that the data is multi-dimensional. Rather than taking a Teuchos::Comm object in its constructor, an MDMap object takes a Domi::MDComm object. An MDComm object is ordered in a structured way, with processors assigned to each axis, whose product of lengths is the total number of processors of its underlying communicator.

At a bare minimum, an MDMap constructor requires an MDComm and an array of dimensions indicating the size and shape of the data structure(s) the MDMap will describe. Each dimension will be decomposed along each axis, according to the number of processors along that axis in the MDComm, in as even a fashion as is possible.

Each axis may be flagged as periodic when constructing the MDComm. This attribute is transferred to the MDMap at construction.

There are two distinct buffer padding concepts employed by MDMaps, communication padding and boundary padding. Communication padding refers to local extensions of the data structure used to receive communication from neighboring processors. Boundary padding refers to points added to the exterior of the global domain that are useful, for example, in certain finite differencing algorithms. They are treated separately from the dimensions so that the decomposition will be performed relative to the dimensions only, and not influenced by the size of the boundary padding.

From a global perspective, the only padding points that are visible are boundary padding points. Communication padding is internal and essentially invisible from this perspective. From a local perspective, a local buffer can have padding on all sides. Some of this padding could be boundary padding, some could be communication padding, and some could be both.

The sizes of boundary and communication padding can be different along each axis, and the sizes of boundary and communication padding on any given axis can be different from each other.

When you want to loop over the elements of a data structure described by an MDMap, you can obtain the local looping bounds along a given axis with the getLocalBounds(int axis) method. This returns a Slice object, whose start() and stop() methods return the loop bounds (stop() is non-inclusive). This method also takes a boolean flag indicating whether the bounds should include the communication padding, which defaults to false.

There are two methods for obtaining global loop bounds, getGlobalBounds(int axis) and getGlobalRankBounds(int axis). The first returns the global bounds for the full data structure and the second returns the bounds on this processor.

Note that all indexes start at zero, whether global or local, whether unique 1D IDs or multidimensional indexes, whether boundary or communication padding is present. This allows negative indexes to represent reverse indexing from the end of a dimension.

Constructor & Destructor Documentation

Domi::MDMap::MDMap ( const Teuchos::RCP< const MDComm mdComm,
const Teuchos::ArrayView< const dim_type > &  dimensions,
const Teuchos::ArrayView< const int > &  commPad = Teuchos::ArrayView< const int >(),
const Teuchos::ArrayView< const int > &  bndryPad = Teuchos::ArrayView< const int >(),
const Teuchos::ArrayView< const int > &  replicatedBoundary = Teuchos::ArrayView< const int >(),
const Layout  layout = DEFAULT_ORDER 
)

Main constructor.

Parameters
mdComm[in] an RCP of an MDComm (multi-dimensional communicator), on which this MDMap will be built.
dimensions[in] the dimensions of the map along each axis. This array must be the same number of dimensions as the MDComm.
commPad[in] the number of indexes in the communication padding along each axis. If this array is less than the number of dimensions, unspecified communication padding will be set to zero.
bndryPad[in] the number of indexes in the boundary padding along each axis. If this array is less than the number of dimensions, unspecified unspecified boundary padding sizes will be set to zero.
replicatedBoundary[in] An array of ints which are simple flags denoting whether each axis contains replicated boundary points (RBPs). RBPs pertain only to periodic axes. As an example, consider a map that describes the decomposition of a longitude coordinate. If the global lower boundary of this domain is longitude = 0 degrees, and the global upper boundary of this domain is longitude = 360 degrees, then this axis has a replicated boundary. Note that communication of the boundary padding will skip these boundary points, and it is up to the user to ensure that their values are equal or compatible.
layout[in] the storage order of the map
Domi::MDMap::MDMap ( Teuchos::ParameterList &  plist)

Constructor with ParameterList.

Parameters
plist[in] ParameterList with construction information

This constructor uses the Teuchos::DefaultComm

Domi::MDMap::MDMap ( const Teuchos::RCP< const Teuchos::Comm< int > >  teuchosComm,
Teuchos::ParameterList &  plist 
)

Constructor with Teuchos::Comm and ParameterList.

Parameters
teuchosComm[in] The Teuchos Communicator. Note that an MDComm will be constructed from the information in plist.
plist[in] ParameterList with construction information
Domi::MDMap::MDMap ( const Teuchos::RCP< const MDComm mdComm,
Teuchos::ParameterList &  plist 
)

Constructor with MDComm and ParameterList.

Parameters
mdComm[in] an RCP of an MDComm (multi-dimensional communicator), on which this MDMap will be built.
plist[in] ParameterList with construction information
Domi::MDMap::MDMap ( const Teuchos::RCP< const MDComm mdComm,
const Teuchos::ArrayView< Slice > &  myGlobalBounds,
const Teuchos::ArrayView< padding_type > &  padding = Teuchos::ArrayView< padding_type >(),
const Teuchos::ArrayView< const int > &  replicatedBoundary = Teuchos::ArrayView< const int >(),
const Layout  layout = DEFAULT_ORDER 
)

Constructor with global bounds for this processor.

Parameters
mdComm[in] an RCP of an MDComm (multi-dimensional communicator), on which the MDMap will be built
myGlobalBounds[in] an array of Slices, one for each axis, that represent the global indexes of the bounds on this processor, excluding padding
padding[in] an array of padding_type (a 2-tuple of integers) specifying the local padding along each axis. Since this is local, the MDMap constructor can determine from context whether each pad refers to communication padding or boundary padding.
replicatedBoundary[in] An array of ints which are simple flags denoting whether each axis contains replicated boundary points (RBPs). RBPs pertain only to periodic axes. As an example, consider a map that describes the decomposition of a longitude coordinate. If the global lower boundary of this domain is longitude = 0 degrees, and the global upper boundary of this domain is longitude = 360 degrees, then this axis has a replicated boundary. Note that communication of the boundary padding will skip these boundary points, and it is up to the user to ensure that their values are equal or compatible.
layout[in] the storage order of the map
Domi::MDMap::MDMap ( const MDMap source)

Copy constructor.

Parameters
source[in] the source MDMap to be copied
Domi::MDMap::MDMap ( const MDMap parent,
int  axis,
dim_type  index 
)

Parent/single global ordinal sub-map constructor.

Parameters
parent[in] an MDMap, from which this sub-map will be derived.
axis[in] the axis to which the slice argument applies
index[in] a global ordinal that defines the sub-map.

This constructor will return an MDMap that is one dimension less than the dimensions of the parent MDMap (unless the parent MDMap is already one dimension).

Domi::MDMap::MDMap ( const MDMap parent,
int  axis,
const Slice slice,
int  bndryPad = 0 
)

Parent/single slice sub-map constructor.

Parameters
parent[in] an MDMap, from which this sub-map will be derived.
axis[in] the axis to which the slice argument applies
slice[in] a Slice of global axis indexes that defines the sub-map. This slice must not include indexes from the boundary padding along each axis.
bndryPad[in] The size of the boundary padding along the altered axis of the new sub-map. This may include indexes from the boundary padding of the parent MDMap, but it does not have to.

This constructor will return an MDMap that is the same number of dimensions as the parent MDMap, but with a smaller dimension along the given axis (unless the given Slice represents the entire axis).

Domi::MDMap::MDMap ( const MDMap parent,
const Teuchos::ArrayView< Slice > &  slices,
const Teuchos::ArrayView< int > &  bndryPad = Teuchos::ArrayView< int >() 
)

Parent/array of slices sub-map constructor.

Parameters
parent[in] an MDMap, from which this sub-map will be derived.
slices[in] an array of Slices of global axis indexes that defines the sub-map. These slices must not include indexes from the boundary padding along each axis.
bndryPad[in] The boundary padding of the new sub-map. These may include indexes from the boundary padding of the parent MDMap, but they do not have to.

Member Function Documentation

Teuchos::RCP< const MDMap > Domi::MDMap::getAugmentedMDMap ( const dim_type  leadingDim,
const dim_type  trailingDim = 0 
) const

Return an RCP to a new MDMap that is a simple augmentation of this MDMap.

Parameters
leadingDim[in] If this value is greater than 0, then add a dimension to the new MDMap that comes before the existing dimensions and has this many degrees of freedom.
trailingDim[in] If this value is greater than 0, then add a dimension to the new MDMap that comes after the existing dimensions and has this many degrees of freedom.

Note that the new dimension(s) will not be distributed; i.e. the commDim for the new leading dimension (if requested) will be 1 and the commDim for the trailing dimension (if requested) will be 1.

Teuchos::RCP< const Domi::MDMap > Domi::MDMap::getAxisMap ( int  axis) const

Return an axis map for the given axis.

An axis map is a 1D map along a given axis.

Teuchos::ArrayView< Teuchos::RCP< const Domi::MDMap > > Domi::MDMap::getAxisMaps ( ) const

Return an array of axis maps.

An axis map is a 1D map along a given axis.

int Domi::MDMap::getBndryPadSize ( int  axis) const

Get the boundary padding size along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This returns the value of the boundary padding along the given axis at the time of construction, regardless of whether a sub-map reduced these values.

Teuchos::Array< int > Domi::MDMap::getBndryPadSizes ( ) const

Get the boundary padding sizes along each axis.

This returns an array of the boundary padding along each axis at the time of construction, regardless of whether a subsequent sub-map reduced these values.

int Domi::MDMap::getCommDim ( int  axis) const
inline

Get the communicator size along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This method will throw a Domi::SubcommunicatorError if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

Teuchos::Array< int > Domi::MDMap::getCommDims ( ) const
inline

Get the communicator sizes along each axis.

This method will throw a Domi::SubcommunicatorError if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

int Domi::MDMap::getCommIndex ( int  axis) const
inline

Get the axis rank of this processor.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This method will throw a Domi::SubcommunicatorError if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

int Domi::MDMap::getCommPadSize ( int  axis) const

Get the communication padding size along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This returns the value of the communication padding along the given axis at the time of construction, regardless of the processor's position relative to the domain boundary.

Slice Domi::MDMap::getGlobalBounds ( int  axis,
bool  withBndryPad = false 
) const

Get the bounds of the global problem.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
withBndryPad[in] specify whether the bounds should include boundary padding or not
dim_type Domi::MDMap::getGlobalDim ( int  axis,
bool  withBndryPad = false 
) const

Get the global dimension along the specified axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
withBndryPad[in] specify whether the dimension should include boundary padding or not
size_type Domi::MDMap::getGlobalID ( size_type  localID) const

Convert a local ID to a global ID.

Parameters
localID[in] a unique 1D local identifier
size_type Domi::MDMap::getGlobalID ( const Teuchos::ArrayView< dim_type > &  globalIndex) const

convert a global index to a global ID

Parameters
globalIndex[in] a global index
Teuchos::Array< dim_type > Domi::MDMap::getGlobalIndex ( size_type  globalID) const

Convert a global ID to a global index.

Parameters
globalID[in] a unique 1D global identifier
Slice Domi::MDMap::getGlobalRankBounds ( int  axis,
bool  withBndryPad = false 
) const

Get the global loop bounds along the specified axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
withBndryPad[in] specify whether the dimension should include boundary padding or not

The loop bounds are returned in the form of a Slice, in which the start() method returns the loop begin value, and the stop() method returns the non-inclusive end value.

Teuchos::ArrayView< const Slice > Domi::MDMap::getLocalBounds ( ) const

Get the local loop bounds along every axis.

The loop bounds are returned in the form of a Slice, in which the start() method returns the loop begin value, and the stop() method returns the non-inclusive end value. For this method, padding is included in the bounds.

Slice Domi::MDMap::getLocalBounds ( int  axis,
bool  withPad = false 
) const

Get the local loop bounds along the specified axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
withPad[in] specify whether the dimension should include padding or not

The loop bounds are returned in the form of a Slice, in which the start() method returns the loop begin value, and the stop() method returns the non-inclusive end value.

dim_type Domi::MDMap::getLocalDim ( int  axis,
bool  withPad = false 
) const

Get the local dimension along the specified axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
withPad[in] specify whether the dimension should include padding or not
size_type Domi::MDMap::getLocalID ( size_type  globalID) const

Convert a global ID to a local ID.

Parameters
globalID[in] a unique 1D global identifier

This method can throw a Domi::RangeError if the global ID is not on the current processor.

size_type Domi::MDMap::getLocalID ( const Teuchos::ArrayView< dim_type > &  localIndex) const

Convert a local index to a local ID.

Parameters
localIndex[in] a local index
Teuchos::Array< dim_type > Domi::MDMap::getLocalIndex ( size_type  localID) const

Convert a local ID to a local index.

Parameters
localID[in] a unique 1D local identifier
Slice Domi::MDMap::getLocalInteriorBounds ( int  axis) const

Get the local interior loop bounds along the specified axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

Local interior loop bounds are the same as local loop bounds without padding, except that for non-periodic axes the global end points of the given axis are excluded. For periodic axes, the local interior loop bounds are exactly the same as local loop bounds without padding.

The loop bounds are returned in the form of a Slice, in which the start() method returns the loop begin value, and the stop() method returns the non-inclusive end value.

int Domi::MDMap::getLowerBndryPad ( int  axis) const

Get the size of the lower boundary padding along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
int Domi::MDMap::getLowerNeighbor ( int  axis) const
inline

Get the rank of the lower neighbor.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This method will throw a Domi::SubcommunicatorError if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

If the periodic flag for the given axis is off, and the axis rank of the calling processor is zero, then this method returns -1.

If the periodic flag for the given axis is on, and the axis rank of the calling processor is the highest axis rank processor along this axis, then the returned lower neighbor will be zero.

int Domi::MDMap::getLowerPadSize ( int  axis) const

Get the size of the lower padding along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

Note that the returned padding can be either communication padding or boundary padding as appropriate.

Teuchos::RCP< const MDComm > Domi::MDMap::getMDComm ( ) const
inline

Access the underlying MDComm object.

Return a reference counted pointer to the Domi::MDComm object upon which this MDMap is built.

Teuchos::RCP< const Teuchos::Comm< int > > Domi::MDMap::getTeuchosComm ( ) const
inline

Get the Teuchos communicator.

Note that if the communicator is not a full communicator, i.e. a sub-communicator, that the underlying Comm pointer may be NULL, depending on this processor's rank.

int Domi::MDMap::getUpperBndryPad ( int  axis) const

Get the size of the upper boundary padding along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)
int Domi::MDMap::getUpperNeighbor ( int  axis) const
inline

Get the rank of the upper neighbor.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This method will throw a Domi::SubcommunicatorError if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

If the periodic flag for the given axis is off, and the axis rank of the calling processor is the highest axis rank processor along this axis, then this method returns -1.

If the periodic flag for the given axis is on, and the axis rank of the calling processor is zero, then the returned lower neighbor will be the highest axis rank processor along this axis.

int Domi::MDMap::getUpperPadSize ( int  axis) const

Get the size of the upper padding along the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

Note that the returned padding can be either communication padding or boundary padding as appropriate.

bool Domi::MDMap::hasPadding ( ) const

Return true if there is any padding stored locally.

Note that it is not as simple as whether there were communication padding specified in the constructor. Suppose non-zero communication padding was specified, but the problem is on one processor and no boundary padding was specified: this method will return false. Similarly, if no communication padding was specified but boundary padding was, then boundary processors will have padding and this method will return true.

bool Domi::MDMap::isCompatible ( const MDMap mdMap) const

True if two MDMaps are "compatible".

Parameters
mdMap[in] MDMap to compare against

Two MDMaps are considered "compatible" if all of the following are true:

  1. The commDims of their underlying MDComms are identical.
  2. Their global dimensions are identical.
  3. Their local dimensions, not including padding, are identical.
bool Domi::MDMap::isContiguous ( ) const

True if there are no stride gaps on any processor.

An MDMap constructed from a communicator and dimensions will always be contiguous. An MDMap that is a slice of a parent MDMap will generally be non-contiguous, with some exceptions. There are cases where some local data is contiguous and some is not, but this method returns True only if all processes' local data is contiguous.

bool Domi::MDMap::isPeriodic ( int  axis) const
inline

Return the periodic flag for the given axis.

Parameters
axis[in] the index of the axis (from zero to the number of dimensions - 1)

This method will throw a Domi::SubcommunicatorError if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

bool Domi::MDMap::isSameAs ( const MDMap mdMap,
const int  verbose = 0 
) const

True if two MDMaps are "identical".

Parameters
mdMap[in] MDMap to compare against
verbose[in] set to one to print why MDMaps are not the same

Two MDMaps are considered "identical" if all of the following are true:

  1. They are compatible.
  2. Their underlying communicators are congruent (have the same number of processes, in the same order: this corresponds to the MPI_IDENT or MPI_CONGRUENT return values of MPI_Comm_compare).
  3. Both are distributed or not distributed.
  4. Their global bounds identical on each process.
  5. Their local dimensions, including padding, are identical.
int Domi::MDMap::numDims ( ) const
inline

Get the number of dimensions.

This method will return 0 if the communicator is a sub-communicator and this processor does not belong to the sub-communicator.

bool Domi::MDMap::onSubcommunicator ( ) const
inline

Query whether this processor is on the sub-communicator.

Sub-communicators are formed when a parent MDComm is sliced by using the (parent,ordinal) or (parent,slices) constructor. For a full communicator, this method will always return true.

MDMap& Domi::MDMap::operator= ( const MDMap source)

Assignment operator.

Parameters
source[in] source MDComm

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

Generated on Thu Mar 28 2024 09:28:09 for Domi by doxygen 1.8.5