Ifpack2 Templated Preconditioning Package
Version 1.0
|
Ifpack2::LinePartitioner: A class to define partitions into a set of lines. More...
#include <Ifpack2_LinePartitioner_decl.hpp>
Public Member Functions | |
LinePartitioner (const Teuchos::RCP< const row_graph_type > &graph) | |
Constructor. More... | |
virtual | ~LinePartitioner () |
Destructor. More... | |
void | setPartitionParameters (Teuchos::ParameterList &List) |
Set the partitioner's parameters (none for linear partitioning). More... | |
void | computePartitions () |
Compute the partitions. More... | |
Public Member Functions inherited from Ifpack2::OverlappingPartitioner< GraphType > | |
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 | compute () |
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... | |
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... | |
Additional Inherited Members | |
Protected Attributes inherited from Ifpack2::OverlappingPartitioner< GraphType > | |
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... | |
Ifpack2::LinePartitioner: A class to define partitions into a set of lines.
These "line" partitions could then be used in to do block Gauss-Seidel smoothing, for instance.
The current implementation uses a local line detection inspired by the work of Mavriplis for convection-diffusion (AIAA Journal, Vol 37(10), 1999).
Here we use coordinate information to pick "close" points if they are sufficiently far away from the "far" points. We also make sure the line can never double back on itself, so that the associated sub-matrix could (in theory) be handed off to a fast triangular solver. This implementation doesn't actual do that, however.
This implementation is derived from the related routine in ML.
Supported parameters:
"partitioner: line detection threshold"
: if \(||x_j - x_i||^2 < thresh * \max_k||x_k - x_i||^2\), then the points are close enough to line smooth <Scalar> "partitioner: coordinates"
: coordinates of local nodes < Teuchos::MultiVector<double> > "partitioner: PDE equations"
: number of equations per node <int> GraphType | Specialization of Tpetra::RowGraph or Tpetra::CrsGraph. |
Ifpack2::LinePartitioner< GraphType, Scalar >::LinePartitioner | ( | const Teuchos::RCP< const row_graph_type > & | graph | ) |
Constructor.
|
virtual |
Destructor.
|
virtual |
Set the partitioner's parameters (none for linear partitioning).
Implements Ifpack2::OverlappingPartitioner< GraphType >.
|
virtual |
Compute the partitions.
Implements Ifpack2::OverlappingPartitioner< GraphType >.