Ifpack2 Templated Preconditioning Package
Version 1.0
|
ILUT (incomplete LU factorization with threshold) of a Tpetra sparse matrix. More...
#include <Ifpack2_ILUT_decl.hpp>
Public Types | |
Typedefs | |
typedef MatrixType::scalar_type | scalar_type |
The type of the entries of the input MatrixType. More... | |
typedef MatrixType::local_ordinal_type | local_ordinal_type |
The type of local indices in the input MatrixType. More... | |
typedef MatrixType::global_ordinal_type | global_ordinal_type |
The type of global indices in the input MatrixType. More... | |
typedef MatrixType::node_type | node_type |
The Node type used by the input MatrixType. More... | |
typedef Teuchos::ScalarTraits < scalar_type >::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
typedef Tpetra::RowMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | row_matrix_type |
Type of the Tpetra::RowMatrix specialization that this class uses. More... | |
typedef Tpetra::CrsMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | crs_matrix_type |
Type of the Tpetra::CrsMatrix specialization that this class uses for the L and U factors. More... | |
typedef Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > | map_type |
Type of the Tpetra::Map specialization that this class uses. More... | |
Public Types inherited from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > | |
typedef Teuchos::ScalarTraits < MatrixType::scalar_type > ::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
Public Member Functions | |
Constructors and Destructors | |
ILUT (const Teuchos::RCP< const row_matrix_type > &A) | |
Constructor. More... | |
virtual | ~ILUT ()=default |
Destructor. More... | |
Methods for setting up and computing the incomplete factorization | |
void | setParameters (const Teuchos::ParameterList ¶ms) |
Set preconditioner parameters. More... | |
void | initialize () |
Clear any previously computed factors. More... | |
bool | isInitialized () const |
Returns true if the preconditioner has been successfully initialized. More... | |
void | compute () |
Compute factors L and U using the specified diagonal perturbation thresholds and relaxation parameters. More... | |
bool | isComputed () const |
If compute() is completed, this query returns true, otherwise it returns false. More... | |
Implementation of Ifpack2::Details::CanChangeMatrix | |
virtual void | setMatrix (const Teuchos::RCP< const row_matrix_type > &A) |
Change the matrix to be preconditioned. More... | |
Implementation of Tpetra::Operator | |
void | apply (const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const |
Apply the ILUT preconditioner to X, resulting in Y. More... | |
Teuchos::RCP< const map_type > | getDomainMap () const |
Tpetra::Map representing the domain of this operator. More... | |
Teuchos::RCP< const map_type > | getRangeMap () const |
Tpetra::Map representing the range of this operator. More... | |
bool | hasTransposeApply () const |
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable). More... | |
Mathematical functions | |
Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const |
Returns the input matrix's communicator. More... | |
Teuchos::RCP< const row_matrix_type > | getMatrix () const |
Returns a reference to the matrix to be preconditioned. More... | |
Teuchos::RCP< const crs_matrix_type > | getL () const |
Returns a reference to the L factor. More... | |
Teuchos::RCP< const crs_matrix_type > | getU () const |
Returns a reference to the U factor. More... | |
int | getNumInitialize () const |
Returns the number of calls to Initialize(). More... | |
int | getNumCompute () const |
Returns the number of calls to Compute(). More... | |
int | getNumApply () const |
Returns the number of calls to apply(). More... | |
double | getInitializeTime () const |
Returns the time spent in Initialize(). More... | |
double | getComputeTime () const |
Returns the time spent in Compute(). More... | |
double | getApplyTime () const |
Returns the time spent in apply(). More... | |
size_t | getNodeSmootherComplexity () const |
Get a rough estimate of cost per iteration. More... | |
double | getLevelOfFill () const |
The level of fill. More... | |
magnitude_type | getAbsoluteThreshold () const |
Get absolute threshold value. More... | |
magnitude_type | getRelativeThreshold () const |
Get relative threshold value. More... | |
magnitude_type | getRelaxValue () const |
Get the relax value. More... | |
magnitude_type | getDropTolerance () const |
Gets the dropping tolerance. More... | |
global_size_t | getGlobalNumEntries () const |
Returns the number of nonzero entries in the global graph. More... | |
size_t | getNodeNumEntries () const |
Returns the number of nonzero entries in the local graph. 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::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > | |
virtual | ~Preconditioner () |
Destructor. More... | |
Public Member Functions inherited from Ifpack2::Details::CanChangeMatrix< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > | |
virtual void | setMatrix (const Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > &A)=0 |
Set the new matrix. More... | |
virtual | ~CanChangeMatrix () |
Destructor. More... | |
ILUT (incomplete LU factorization with threshold) of a Tpetra sparse matrix.
A | specialization of Tpetra::RowMatrix. |
This class computes a sparse ILUT (incomplete LU) factorization with specified fill and drop tolerance, of the local part of a given sparse matrix represented as a Tpetra::RowMatrix or Tpetra::CrsMatrix. The "local part" is the square diagonal block of the matrix owned by the calling process. Thus, if the input matrix is distributed over multiple MPI processes, this preconditioner is equivalent to nonoverlapping additive Schwarz domain decomposition over the MPI processes, with ILUT as the subdomain solver on each process.
typedef MatrixType::scalar_type Ifpack2::ILUT< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::ILUT< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::ILUT< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type Ifpack2::ILUT< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::ILUT< MatrixType >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::ILUT< MatrixType >::row_matrix_type |
Type of the Tpetra::RowMatrix specialization that this class uses.
typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::ILUT< MatrixType >::crs_matrix_type |
Type of the Tpetra::CrsMatrix specialization that this class uses for the L and U factors.
typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::ILUT< MatrixType >::map_type |
Type of the Tpetra::Map specialization that this class uses.
|
explicit |
Constructor.
A | [in] The sparse matrix to factor, as a Tpetra::RowMatrix. (Tpetra::CrsMatrix inherits from this, so you may use a Tpetra::CrsMatrix here instead.) |
The factorization will not modify the input matrix. It stores the L and U factors in the incomplete factorization separately.
|
virtualdefault |
Destructor.
|
virtual |
Set preconditioner parameters.
ILUT implements the following parameters:
double
) magnitude_type
) magnitude_type
) magnitude_type
) magnitude_type
) "fact: drop tolerance" is the magnitude threshold for dropping entries. It corresponds to the \(\tau\) parameter in Saad's original description of ILUT. "fact: ilut level-of-fill" controls the number of entries to keep in the strict upper triangle of the current row, and in the strict lower triangle of the current row. It does not correspond to the \(p\) parameter in Saad's original description. This parameter represents a maximum fill fraction. In this implementation, the L and U factors always contains nonzeros corresponding to the original sparsity pattern of A, so this value should be >= 1.0. Letting \(fill = \frac{(level-of-fill - 1)*nnz(A)}{2*N}\), each row of the computed L and U factors contains at most \(fill\) nonzero elements in addition to those from the sparsity pattern of A. ILUT always keeps the diagonal entry in the current row, regardless of the drop tolerance or fill level.
The absolute and relative threshold parameters affect how this code modifies the diagonal entry of the output factor. These parameters are not part of the original ILUT algorithm, but we include them for consistency with other Ifpack2 preconditioners.
The "fact: relax value" parameter currently has no effect.
|
virtual |
Clear any previously computed factors.
You may call this before calling compute(). The compute() method will call this automatically if it has not yet been called. If you call this after calling compute(), you must recompute the factorization (by calling compute() again) before you may call apply().
|
inlinevirtual |
Returns true
if the preconditioner has been successfully initialized.
|
virtual |
Compute factors L and U using the specified diagonal perturbation thresholds and relaxation parameters.
This function computes the ILUT factors L and U using the current:
|
inlinevirtual |
If compute() is completed, this query returns true, otherwise it returns false.
|
virtual |
Change the matrix to be preconditioned.
A | [in] The new matrix. |
! isInitialized ()
! isComputed ()
Calling this method resets the preconditioner's state. After calling this method with a nonnull input, you must first call initialize() and compute() (in that order) before you may call apply().
You may call this method with a null input. If A is null, then you may not call initialize() or compute() until you first call this method again with a nonnull input. This method invalidates any previous factorization whether or not A is null, so calling setMatrix() with a null input is one way to clear the preconditioner's state (and free any memory that it may be using).
The new matrix A need not necessarily have the same Maps or even the same communicator as the original matrix.
|
virtual |
Apply the ILUT preconditioner to X, resulting in Y.
X | [in] Input multivector; "right-hand side" of the solve. |
Y | [out] Output multivector; result of the solve. |
|
virtual |
Tpetra::Map representing the domain of this operator.
|
virtual |
Tpetra::Map representing the range of this operator.
bool Ifpack2::ILUT< MatrixType >::hasTransposeApply | ( | ) | const |
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable).
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::ILUT< MatrixType >::getComm | ( | ) | const |
Returns the input matrix's communicator.
|
virtual |
Returns a reference to the matrix to be preconditioned.
|
inline |
Returns a reference to the L factor.
|
inline |
Returns a reference to the U factor.
|
virtual |
Returns the number of calls to Initialize().
|
virtual |
Returns the number of calls to Compute().
|
virtual |
Returns the number of calls to apply().
|
virtual |
Returns the time spent in Initialize().
|
virtual |
Returns the time spent in Compute().
|
virtual |
Returns the time spent in apply().
size_t Ifpack2::ILUT< MatrixType >::getNodeSmootherComplexity | ( | ) | const |
Get a rough estimate of cost per iteration.
|
inline |
|
inline |
Get absolute threshold value.
|
inline |
Get relative threshold value.
|
inline |
Get the relax value.
|
inline |
Gets the dropping tolerance.
global_size_t Ifpack2::ILUT< MatrixType >::getGlobalNumEntries | ( | ) | const |
Returns the number of nonzero entries in the global graph.
size_t Ifpack2::ILUT< MatrixType >::getNodeNumEntries | ( | ) | const |
Returns the number of nonzero entries in the local graph.
std::string Ifpack2::ILUT< MatrixType >::description | ( | ) | const |
Return a simple one-line description of this object.
void Ifpack2::ILUT< MatrixType >::describe | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
) | const |
Print the object with some verbosity level to an FancyOStream object.