Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Attributes | List of all members
Ifpack2::OverlappingPartitioner< GraphType > Class Template Referenceabstract

Create overlapping partitions of a local graph. More...

#include <Ifpack2_OverlappingPartitioner_decl.hpp>

Inheritance diagram for Ifpack2::OverlappingPartitioner< GraphType >:
Inheritance graph
[legend]

Public Member Functions

 OverlappingPartitioner (const Teuchos::RCP< const row_graph_type > &graph)
 Constructor. More...
 
virtual ~OverlappingPartitioner ()
 Destructor. More...
 
int numLocalParts () const
 Number of computed local partitions. More...
 
int overlappingLevel () const
 The number of levels of overlap. More...
 
local_ordinal_type operator() (const local_ordinal_type MyRow) const
 Local index of the nonoverlapping partition of the given row. More...
 
local_ordinal_type operator() (const local_ordinal_type i, const local_ordinal_type j) const
 Local index of the overlapping partition of the j-th vertex in partition i. More...
 
size_t numRowsInPart (const local_ordinal_type Part) const
 the number of rows contained in the given partition. More...
 
void rowsInPart (const local_ordinal_type Part, Teuchos::ArrayRCP< local_ordinal_type > &List) const
 Fill List with the local indices of the rows in the (overlapping) partition Part. More...
 
virtual Teuchos::ArrayView
< const local_ordinal_type > 
nonOverlappingPartition () const
 A view of the local indices of the nonoverlapping partitions of each local row. More...
 
virtual void setParameters (Teuchos::ParameterList &List)
 Set all the parameters for the partitioner. More...
 
virtual void setPartitionParameters (Teuchos::ParameterList &List)=0
 Set all the parameters for the partitioner. More...
 
virtual void compute ()
 Computes the partitions. Returns 0 if successful. More...
 
virtual void computePartitions ()=0
 Computes the partitions. Returns 0 if successful. More...
 
virtual void computeOverlappingPartitions ()
 Computes the partitions. Returns 0 if successful. More...
 
virtual bool isComputed () const
 Returns true if partitions have been computed successfully. More...
 
virtual std::ostream & print (std::ostream &os) const
 Prints basic information on iostream. This function is used by operator<<. More...
 
Implementation of Teuchos::Describable
std::string description () const
 Return a simple one-line description of this object. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print the object with some verbosity level to an FancyOStream object. More...
 
- Public Member Functions inherited from Ifpack2::Partitioner< GraphType >
virtual ~Partitioner ()
 Destructor. More...
 

Protected Attributes

int NumLocalParts_
 Number of local subgraphs. More...
 
Teuchos::Array
< local_ordinal_type > 
Partition_
 Mapping from local row to partition number. More...
 
Teuchos::Array
< Teuchos::ArrayRCP
< local_ordinal_type > > 
Parts_
 Mapping from partition to all local rows it contains. More...
 
Teuchos::RCP< const
row_graph_type > 
Graph_
 The graph to be partitioned. More...
 
int OverlappingLevel_
 Level of overlap. More...
 
bool IsComputed_
 If true, the graph has been successfully partitioned. More...
 
bool verbose_
 If true, information are reported to stdout. More...
 
bool maintainSparsity_
 If true, only add row to partition (block) if doing so won't add new columns to the column map. More...
 

Detailed Description

template<class GraphType>
class Ifpack2::OverlappingPartitioner< GraphType >

Create overlapping partitions of a local graph.

Template Parameters
GraphTypeSpecialization of Tpetra::CrsGraph or Tpetra::RowGraph.

This class enables the extension of nonoverlapping partitions to an arbitrary amount of constant overlap. "Overlap" here refers to the overlap among the local parts, and not the overlap among the processes.

For partitions with varying overlap, use UserPartitioner, which inherits from OverlappingPartitioner.

Supported parameters are:

This class implements common functionality for derived classes like LinearPartitioner. Graphs given as input to any derived class must contain no singletons. Usually, this means that the graph is the graph of a Tpetra::RowMatrix that has been filtered using SingletonFilter.

Constructor & Destructor Documentation

template<class GraphType >
Ifpack2::OverlappingPartitioner< GraphType >::OverlappingPartitioner ( const Teuchos::RCP< const row_graph_type > &  graph)

Constructor.

template<class GraphType >
Ifpack2::OverlappingPartitioner< GraphType >::~OverlappingPartitioner ( )
virtual

Destructor.

Member Function Documentation

template<class GraphType >
int Ifpack2::OverlappingPartitioner< GraphType >::numLocalParts ( ) const
virtual

Number of computed local partitions.

This is a (signed) int because negative values are significant. We can't use local_ordinal_type here, because that type may be either signed or unsigned.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
int Ifpack2::OverlappingPartitioner< GraphType >::overlappingLevel ( ) const
virtual

The number of levels of overlap.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
GraphType::local_ordinal_type Ifpack2::OverlappingPartitioner< GraphType >::operator() ( const local_ordinal_type  MyRow) const
virtual

Local index of the nonoverlapping partition of the given row.

Parameters
MyRow[in] Local index of the row.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
GraphType::local_ordinal_type Ifpack2::OverlappingPartitioner< GraphType >::operator() ( const local_ordinal_type  i,
const local_ordinal_type  j 
) const
virtual

Local index of the overlapping partition of the j-th vertex in partition i.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
size_t Ifpack2::OverlappingPartitioner< GraphType >::numRowsInPart ( const local_ordinal_type  Part) const
virtual

the number of rows contained in the given partition.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
void Ifpack2::OverlappingPartitioner< GraphType >::rowsInPart ( const local_ordinal_type  Part,
Teuchos::ArrayRCP< local_ordinal_type > &  List 
) const
virtual

Fill List with the local indices of the rows in the (overlapping) partition Part.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
Teuchos::ArrayView< const typename GraphType::local_ordinal_type > Ifpack2::OverlappingPartitioner< GraphType >::nonOverlappingPartition ( ) const
virtual

A view of the local indices of the nonoverlapping partitions of each local row.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
void Ifpack2::OverlappingPartitioner< GraphType >::setParameters ( Teuchos::ParameterList List)
virtual

Set all the parameters for the partitioner.

The supported parameters are:

  • "partitioner: overlap" (int, default = 0).
  • "partitioner: local parts" (int, default = 1).
  • "partitioner: print level" (bool, default = false).

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
virtual void Ifpack2::OverlappingPartitioner< GraphType >::setPartitionParameters ( Teuchos::ParameterList List)
pure virtual

Set all the parameters for the partitioner.

This function is used by derived classes to set their own parameters. These classes should not derive SetParameters(), so that common parameters can be set just once.

Implemented in Ifpack2::Details::UserPartitioner< GraphType >, Ifpack2::LinePartitioner< GraphType, Scalar >, and Ifpack2::LinearPartitioner< GraphType >.

template<class GraphType >
void Ifpack2::OverlappingPartitioner< GraphType >::compute ( )
virtual

Computes the partitions. Returns 0 if successful.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
virtual void Ifpack2::OverlappingPartitioner< GraphType >::computePartitions ( )
pure virtual
template<class GraphType >
void Ifpack2::OverlappingPartitioner< GraphType >::computeOverlappingPartitions ( )
virtual

Computes the partitions. Returns 0 if successful.

template<class GraphType >
bool Ifpack2::OverlappingPartitioner< GraphType >::isComputed ( ) const
virtual

Returns true if partitions have been computed successfully.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
std::ostream & Ifpack2::OverlappingPartitioner< GraphType >::print ( std::ostream &  os) const
virtual

Prints basic information on iostream. This function is used by operator<<.

Implements Ifpack2::Partitioner< GraphType >.

template<class GraphType >
std::string Ifpack2::OverlappingPartitioner< GraphType >::description ( ) const
virtual

Return a simple one-line description of this object.

Reimplemented from Teuchos::Describable.

template<class GraphType >
void Ifpack2::OverlappingPartitioner< GraphType >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
virtual

Print the object with some verbosity level to an FancyOStream object.

Reimplemented from Teuchos::Describable.

Member Data Documentation

template<class GraphType >
int Ifpack2::OverlappingPartitioner< GraphType >::NumLocalParts_
protected

Number of local subgraphs.

This is a (signed) int because negative values are significant. We can't use local_ordinal_type here, because that type may be either signed or unsigned.

template<class GraphType >
Teuchos::Array<local_ordinal_type> Ifpack2::OverlappingPartitioner< GraphType >::Partition_
protected

Mapping from local row to partition number.

Partition_[i] contains the local index of the nonoverlapping partition to which local row i belongs. If the application has explicitly defined Parts_, then Partition_ is unused.

template<class GraphType >
Teuchos::Array<Teuchos::ArrayRCP<local_ordinal_type> > Ifpack2::OverlappingPartitioner< GraphType >::Parts_
protected

Mapping from partition to all local rows it contains.

Parts_[i][j] is the local index of the j-th row contained in the (overlapping) partition i.

template<class GraphType >
Teuchos::RCP<const row_graph_type> Ifpack2::OverlappingPartitioner< GraphType >::Graph_
protected

The graph to be partitioned.

template<class GraphType >
int Ifpack2::OverlappingPartitioner< GraphType >::OverlappingLevel_
protected

Level of overlap.

template<class GraphType >
bool Ifpack2::OverlappingPartitioner< GraphType >::IsComputed_
protected

If true, the graph has been successfully partitioned.

template<class GraphType >
bool Ifpack2::OverlappingPartitioner< GraphType >::verbose_
protected

If true, information are reported to stdout.

template<class GraphType >
bool Ifpack2::OverlappingPartitioner< GraphType >::maintainSparsity_
protected

If true, only add row to partition (block) if doing so won't add new columns to the column map.

Default is false. This option should decrease computational cost of the apply, but convergence will likely degrade.


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