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

Ifpack2::Partitioner: More...

#include <Ifpack2_Partitioner.hpp>

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

Public Member Functions

virtual ~Partitioner ()
 Destructor. More...
 
virtual int numLocalParts () const =0
 Number of computed local partitions. More...
 
virtual int overlappingLevel () const =0
 The level of overlap. More...
 
virtual LocalOrdinal operator() (LocalOrdinal MyRow) const =0
 The local (nonoverlapping) partition index of the specified local row. More...
 
virtual LocalOrdinal operator() (LocalOrdinal i, LocalOrdinal j) const =0
 The local overlapping partition index of the j-th node in partition i. More...
 
virtual size_t numRowsInPart (const LocalOrdinal Part) const =0
 The number of rows contained in the specified partition. More...
 
virtual void rowsInPart (const LocalOrdinal Part, Teuchos::ArrayRCP< LocalOrdinal > &List) const =0
 Copy into List the rows in the (overlapping) partition Part. More...
 
virtual Teuchos::ArrayView
< const LocalOrdinal > 
nonOverlappingPartition () const =0
 The nonoverlapping partition indices of each local row. More...
 
virtual void setParameters (Teuchos::ParameterList &List)=0
 Set all the parameters for the partitioner. More...
 
virtual void compute ()=0
 Compute the partitions. More...
 
virtual bool isComputed () const =0
 Return true if partitions have been computed successfully. More...
 
virtual std::ostream & print (std::ostream &os) const =0
 Print basic information about the partitioning object. More...
 

Detailed Description

template<class GraphType>
class Ifpack2::Partitioner< GraphType >

Ifpack2::Partitioner:

A class to decompose local graphs.

Summary

Most Ifpack2 users will not need to create or use a Partitioner instance, or write a Partitioner subclass. However, those implementing their own block relaxation algorithms may need to interact with Partitioner or write their own subclass thereof. Ifpack2's main application of this class is in BlockRelaxation. Partitioner defines the diagonal blocks of the matrix that BlockRelaxation uses. BlockRelaxation creates a Partitioner subclass instance internally.

Partitions

A Partitioner instance can partition a local graph by rows. A local graph is one for which, on every process, the column Map contains no entries not in the domain Map on that process. You may use LocalFilter on the graph of a matrix to make a local graph; it excludes entries in the column Map not in the domain Map. This class assumes that the graph is local.

The partitions created by Partitioner implementations are nonoverlapping in the graph sense. This means that each row (or, more appropriately, vertex) of the graph belongs to at most one partition. Furthermore, these nonoverlapping partitions are local: partitions do not cross process boundaries.

operator () (LocalOrdinal i) returns the local partition index corresponding to local row i of the graph. The above implies that the local partition index is unique.

Overlapping decomposition

The OverlappingPartitioner subclass extends the nonoverlapping partitions by the required amount of overlap, considering local vertices only. That is, this overlap does not modify the overlap among the processes. (The mathematical definition of "partition" does not allow overlap, but we accept this as a useful extension.)

Subclasses of Partitioner

Partitioner is just an interface; it does not implement any fucntionality. You cannot create an instance of Partitioner; you must instantiate a concrete implementation thereof. Concrete implementations include:

The constructor takes a Tpetra::RowGraph instance, which is the graph to partition.

Example code

The following example code shows how to use LinearPartitioner, a subclass of Partitioner. Note that most Ifpack2 users will not need to create or interact with Partitioner instances. We only show this example for the sake of developers who might need to implement or use Preconditioner subclasses, and want an example of a Partitioner "in action."

#include "Ifpack2_LinearPartitioner.hpp"
#include "Tpetra_CrsMatrix.hpp"
// ...
// The matrix A must be fill complete.
void example (Tpetra::CrsMatrix<double, int>& A) {
using std::cout;
using std::endl;
typedef Tpetra::CrsGraph<int> graph_type;
typedef Ifpack2::LinearPartitioner<graph_type> partitioner_type;
// Create the partitioner.
partitioner_type partitioner (A.getGraph ());
// Set up the partitioner's parameters.
// We want 16 local partitions,
// and an overlap of 0 among the local partitions.
params.set ("partitioner: local parts", 16);
params.set ("partitioner: overlap", 0);
partitioner.setParameters (params);
// Partition the graph. If the structure of the
// graph changes, you must call compute() again,
// but you need not call setParameters() again.
partitioner.compute ();
// Get the number of partitions created on the calling process.
const int numParts = partitioner.numLocalParts ();
// Get the number of rows in each partition.
for (int i = 0; i < numParts; ++i) {
cout << "Number of rows in partition " << i << ": "
<< partitioner.numRowsInPart (i) << endl;
}
// For nonoverlapping partitions only, operator()(i)
// returns the partition index for each local row.
const size_t numLocalRows = A.getNodeNumRows ();
for (size_t i = 0; i < numLocalRows; ++i) {
cout << "Partition[" << i <<"] = " << partitioner(i) << endl;
}
}

When overlapping partitions are created, users can get the local indices of the rows in each partition:

const int numLocalParts = partitioner.numLocalParts ();
for (int i = 0; i < numLocalParts; ++i) {
cout << "Local rows of Partition " << i << ": [";
for (size_t j = 0; j < partitioner.numRowsInPart (i); ++j) {
cout << partitioner(i,j) << " ";
}
cout << "]";
}

Constructor & Destructor Documentation

template<class GraphType>
virtual Ifpack2::Partitioner< GraphType >::~Partitioner ( )
inlinevirtual

Destructor.

Member Function Documentation

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

Number of computed local partitions.

See Ifpack2_OverlappingPartitioner_decl.hpp for explanation of why this is an int instead of LocalOrdinal.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

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

The level of overlap.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

template<class GraphType>
virtual LocalOrdinal Ifpack2::Partitioner< GraphType >::operator() ( LocalOrdinal  MyRow) const
pure virtual

The local (nonoverlapping) partition index of the specified local row.

Parameters
MyRow[in] Local index of the row.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

template<class GraphType>
virtual LocalOrdinal Ifpack2::Partitioner< GraphType >::operator() ( LocalOrdinal  i,
LocalOrdinal  j 
) const
pure virtual

The local overlapping partition index of the j-th node in partition i.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

template<class GraphType>
virtual size_t Ifpack2::Partitioner< GraphType >::numRowsInPart ( const LocalOrdinal  Part) const
pure virtual

The number of rows contained in the specified partition.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

template<class GraphType>
virtual void Ifpack2::Partitioner< GraphType >::rowsInPart ( const LocalOrdinal  Part,
Teuchos::ArrayRCP< LocalOrdinal > &  List 
) const
pure virtual

Copy into List the rows in the (overlapping) partition Part.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

template<class GraphType>
virtual Teuchos::ArrayView<const LocalOrdinal> Ifpack2::Partitioner< GraphType >::nonOverlappingPartition ( ) const
pure virtual

The nonoverlapping partition indices of each local row.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

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

Set all the parameters for the partitioner.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

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

Compute the partitions.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

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

Return true if partitions have been computed successfully.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.

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

Print basic information about the partitioning object.

Implemented in Ifpack2::OverlappingPartitioner< GraphType >.


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