Ifpack2 Templated Preconditioning Package
Version 1.0
|
ILU(k) factorization of a given Tpetra::RowMatrix. More...
#include <Ifpack2_RILUK_decl.hpp>
Public Types | |
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 node_type::device_type | device_type |
The Kokkos device type of the input MatrixType. More... | |
typedef node_type::execution_space | execution_space |
The Kokkos execution space of the input MatrixType. More... | |
typedef Tpetra::RowMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | row_matrix_type |
Tpetra::RowMatrix specialization used by this class. More... | |
typedef Tpetra::CrsMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | crs_matrix_type |
Tpetra::CrsMatrix specialization used by this class for representing L and U. More... | |
typedef crs_matrix_type::impl_scalar_type | impl_scalar_type |
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector) 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 | |
Teuchos::RCP< const row_matrix_type > | getMatrix () const |
Get the input matrix. More... | |
magnitude_type | getRelaxValue () const |
Get RILU(k) relaxation parameter. More... | |
magnitude_type | getAbsoluteThreshold () const |
Get absolute threshold value. More... | |
magnitude_type | getRelativeThreshold () const |
Get relative threshold value. More... | |
int | getLevelOfFill () const |
Get level of fill (the "k" in ILU(k)). More... | |
Tpetra::CombineMode | getOverlapMode () |
Get overlap mode type. More... | |
Tpetra::global_size_t | getGlobalNumEntries () const |
Returns the number of nonzero entries in the global graph. More... | |
Teuchos::RCP < Ifpack2::IlukGraph < Tpetra::CrsGraph < local_ordinal_type, global_ordinal_type, node_type > , kk_handle_type > > | getGraph () const |
Return the Ifpack2::IlukGraph associated with this factored matrix. More... | |
const crs_matrix_type & | getL () const |
Return the L factor of the ILU factorization. More... | |
const Tpetra::Vector < scalar_type, local_ordinal_type, global_ordinal_type, node_type > & | getD () const |
Return the diagonal entries of the ILU factorization. More... | |
const crs_matrix_type & | getU () const |
Return the U factor of the ILU factorization. More... | |
Teuchos::RCP< const crs_matrix_type > | getCrsMatrix () const |
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws. 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 Teuchos::Describable interface | |
std::string | description () const |
A one-line description of this object. More... | |
Implementation of Tpetra::Operator | |
Teuchos::RCP< const Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > > | getDomainMap () const |
Returns the Tpetra::Map object associated with the domain of this operator. More... | |
Teuchos::RCP< const Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > > | getRangeMap () const |
Returns the Tpetra::Map object associated with the range of this operator. More... | |
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 (inverse of the) incomplete factorization to X, resulting in Y. 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... | |
virtual void | setZeroStartingSolution (bool zeroStartingSolution) |
Set this preconditioner's parameters. 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... | |
Static Public Member Functions | |
static Teuchos::RCP< const row_matrix_type > | makeLocalFilter (const Teuchos::RCP< const row_matrix_type > &A) |
Return A, wrapped in a LocalFilter, if necessary. More... | |
Protected Attributes | |
Teuchos::RCP< const row_matrix_type > | A_ |
The (original) input matrix for which to compute ILU(k). More... | |
Teuchos::RCP< iluk_graph_type > | Graph_ |
The ILU(k) graph. More... | |
Teuchos::RCP< const row_matrix_type > | A_local_ |
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_type that initialize() constructs temporarily. More... | |
Teuchos::RCP< crs_matrix_type > | L_ |
The L (lower triangular) factor of ILU(k). More... | |
Teuchos::RCP < LocalSparseTriangularSolver < row_matrix_type > > | L_solver_ |
Sparse triangular solver for L. More... | |
Teuchos::RCP< crs_matrix_type > | U_ |
The U (upper triangular) factor of ILU(k). More... | |
Teuchos::RCP < LocalSparseTriangularSolver < row_matrix_type > > | U_solver_ |
Sparse triangular solver for U. More... | |
Teuchos::RCP< vec_type > | D_ |
The diagonal entries of the ILU(k) factorization. More... | |
bool | isKokkosKernelsSpiluk_ |
Optional KokkosKernels implementation. More... | |
Implementation of Kokkos Kernels ILU(k). | |
typedef crs_matrix_type::local_matrix_device_type | local_matrix_device_type |
typedef local_matrix_device_type::StaticCrsGraphType::row_map_type | lno_row_view_t |
typedef local_matrix_device_type::StaticCrsGraphType::entries_type | lno_nonzero_view_t |
typedef local_matrix_device_type::values_type | scalar_nonzero_view_t |
typedef local_matrix_device_type::StaticCrsGraphType::device_type::memory_space | TemporaryMemorySpace |
typedef local_matrix_device_type::StaticCrsGraphType::device_type::memory_space | PersistentMemorySpace |
typedef local_matrix_device_type::StaticCrsGraphType::device_type::execution_space | HandleExecSpace |
typedef KokkosKernels::Experimental::KokkosKernelsHandle < typename lno_row_view_t::const_value_type, typename lno_nonzero_view_t::const_value_type, typename scalar_nonzero_view_t::value_type, HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace > | kk_handle_type |
typedef Ifpack2::IlukGraph < Tpetra::CrsGraph < local_ordinal_type, global_ordinal_type, node_type > , kk_handle_type > | iluk_graph_type |
RILUK (const Teuchos::RCP< const row_matrix_type > &A_in) | |
Constructor that takes a Tpetra::RowMatrix. More... | |
RILUK (const Teuchos::RCP< const crs_matrix_type > &A_in) | |
Constructor that takes a Tpetra::CrsMatrix. More... | |
virtual | ~RILUK () |
Destructor (declared virtual for memory safety). More... | |
void | setParameters (const Teuchos::ParameterList ¶ms) |
void | initialize () |
Initialize by computing the symbolic incomplete factorization. More... | |
void | compute () |
Compute the (numeric) incomplete factorization. More... | |
bool | isInitialized () const |
Whether initialize() has been called on this object. More... | |
bool | isComputed () const |
Whether compute() has been called on this object. More... | |
int | getNumInitialize () const |
Number of successful initialize() calls for this object. More... | |
int | getNumCompute () const |
Number of successful compute() calls for this object. More... | |
int | getNumApply () const |
Number of successful apply() calls for this object. More... | |
double | getInitializeTime () const |
Total time in seconds taken by all successful initialize() calls for this object. More... | |
double | getComputeTime () const |
Total time in seconds taken by all successful compute() calls for this object. More... | |
double | getApplyTime () const |
Total time in seconds taken by all successful apply() calls for this object. More... | |
size_t | getNodeSmootherComplexity () const |
Get a rough estimate of cost per iteration. More... | |
ILU(k) factorization of a given Tpetra::RowMatrix.
MatrixType | A specialization of Tpetra::RowMatrix. |
This class implements a "relaxed" incomplete ILU (ILU) factorization with level k fill. It is based upon the ILU algorithms outlined in Yousef Saad's "Iterative Methods for Sparse Linear Systems", 2nd edition, Chapter 10.
For a complete list of valid parameters, see the documentation of setParameters().
The computed factorization is a function of several parameters:
The graph structure (sparsity pattern) of the matrix: All fill is derived from the original matrix nonzero structure. Level zero fill is defined as the original matrix pattern (nonzero structure), even if the matrix value at an entry is stored as a zero. (Thus it is possible to add entries to the ILU factors by adding zero entries to the original matrix.)
Level of fill: Starting with the original matrix pattern as level fill of zero, the next level of fill is determined by analyzing the graph of the previous level and determining nonzero fill that is a result of combining entries that were from previous level only (not the current level). This rule limits fill to entries that are direct decendents from the previous level graph. Fill for level k is determined by applying this rule recursively. For sufficiently large values of k, the fill would eventually be complete and an exact LU factorization would be computed.
Fraction of relaxation: Ifpack2::RILUK computes the ILU factorization row-by-row. As entries at a given row are computed, some number of them will be dropped because they do match the prescribed sparsity pattern. The relaxation factor determines how these dropped values will be handled. If the factor is zero, then these extra entries will by dropped. This is a classical ILU approach. If the RelaxValue is 1, then the sum of the extra entries will be added to the diagonal. This is a classical Modified ILU (MILU) approach. If RelaxValue is between 0 and 1, then the factor times the sum of extra entries will be added to the diagonal.
For most situations, the relaxation factor should be set to zero. For certain kinds of problems, e.g., reservoir modeling, there is a conservation principle involved such that any operator should obey a zero row-sum property. MILU was designed for these cases and you should set the relaxation factor to 1. For other situations, setting RelaxValue to some nonzero value may improve the stability of factorization, and can be used if the computed ILU factors are poorly conditioned.
Diagonal perturbation: Prior to computing the factorization, it is possible to modify the diagonal entries of the matrix for which the factorization will be computing. If the absolute and relative perturbation values are zero and one, respectively, the factorization will be compute for the original user matrix A. Otherwise, the factorization will computed for a matrix that differs from the original user matrix in the diagonal values only. Below we discuss the details of diagonal perturbations.
Note that the factorization is calculated based upon local ordering. This means that the ordering of the GIDs in the row map is ignored. Initial entries in \(L\), the strictly lower triangular part of A, and \(U\), the strictly upper triangular part of A, are given by
\(L(i,j) = A(i,j)\) if \(j < i\), for local IDs \(i\) and \(j\), even if GID \((j)\) \(>\) GID \((i)\),
and
\(U(i,j) = A(i,j)\) if \(i < j\), for local IDs \(i\) and \(j\), even if GID \((j)\) \(<\) GID \((i)\).
In particular, if the row map GIDs are not in ascending order on processor, then the incomplete factors will be different than those produced by ILU(k) using global IDs. If the row map GIDs are in ascending order, then the factors produced based on LID and GID ordering are the same.
For ill-conditioned matrices, we often have difficulty computing usable incomplete factorizations. The most common source of problems is that the factorization may encounter a small or zero pivot. In that case, the factorization may fail. Even if the factorization succeeds, the factors may be so poorly conditioned that use of them in the iterative phase produces meaningless results. Before we can fix this problem, we must be able to detect it. To this end, we use a simple but effective condition number estimate for \((LU)^{-1}\).
The condition number of a matrix \(B\), called \(cond_p(B)\), is defined as \(cond_p(B) = \|B\|_p\|B^{-1}\|_p\) in some appropriate norm \(p\). \(cond_p(B)\) gives some indication of how many accurate floating point digits can be expected from operations involving the matrix and its inverse. A condition number approaching the accuracy of a given floating point number system, about 15 decimal digits in IEEE double precision, means that any results involving \(B\) or \(B^{-1}\) may be meaningless.
The \(\infty\)-norm of a vector \(y\) is defined as the maximum of the absolute values of the vector entries, and the \(\infty\)-norm of a matrix C is defined as \(\|C\|_\infty = \max_{\|y\|_\infty = 1} \|Cy\|_\infty\). A crude lower bound for the \(cond_\infty(C)\) is \(\|C^{-1}e\|_\infty\) where \(e = (1, 1, \ldots, 1)^T\). It is a lower bound because \(cond_\infty(C) = \|C\|_\infty\|C^{-1}\|_\infty \ge \|C^{-1}\|_\infty \ge |C^{-1}e\|_\infty\).
For our purposes, we want to estimate \(cond_\infty(LU)\), where \(L\) and \(U\) are our incomplete factors. Edmond in his Ph.D. thesis demonstrates that \(\|(LU)^{-1}e\|_\infty\) provides an effective estimate for \(cond_\infty(LU)\). Furthermore, since finding \(z\) such that \(LUz = y\) is a basic kernel for applying the preconditioner, computing this estimate of \(cond_\infty(LU)\) is performed by setting \(y = e\), calling the solve kernel to compute \(z\) and then computing \(\|z\|_\infty\).
If we detect using the above method that our factorization is too ill-conditioned, we can improve the conditioning by perturbing the matrix diagonal and restarting the factorization using this more diagonally dominant matrix. In order to apply perturbation, prior to starting the factorization, we compute a diagonal perturbation of our matrix \(A\) and perform the factorization on this perturbed matrix. The overhead cost of perturbing the diagonal is minimal since the first step in computing the incomplete factors is to copy the matrix \(A\) into the memory space for the incomplete factors. We simply compute the perturbed diagonal at this point.
The actual perturbation values we use are the diagonal values \((d_1, d_2, \ldots, d_n)\) with \(d_i = sgn(d_i)\alpha + d_i\rho\), \(i=1, 2, \ldots, n\), where \(n\) is the matrix dimension and \(sgn(d_i)\) returns the sign of the diagonal entry. This has the effect of forcing the diagonal values to have minimal magnitude of \(\alpha\) and to increase each by an amount proportional to \(\rho\), and still keep the sign of the original diagonal entry.
Every Ifpack2 preconditioner has the following phases of computation:
RILUK constructs the symbolic incomplete factorization (that is, the structure of the incomplete factors) in the initialize() phase. It computes the numerical incomplete factorization (that is, it fills in the factors' entries with their correct values) in the compute() phase. The apply() phase applies the incomplete factorization to a given multivector using two triangular solves.
Each RILUK object keeps track of both the time required for various operations, and the number of times those operations have been applied for that object. The operations tracked include:
The getNum*
methods return the number of times that operation was called. The get*Time
methods return the total number of seconds spent in all invocations of that operation. For example, getApplyTime() returns the number of seconds spent in all apply() calls. For an average time per apply() call, divide by getNumApply(), the total number of calls to apply().
typedef MatrixType::scalar_type Ifpack2::RILUK< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::RILUK< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::RILUK< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type Ifpack2::RILUK< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::RILUK< MatrixType >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
typedef node_type::device_type Ifpack2::RILUK< MatrixType >::device_type |
The Kokkos device type of the input MatrixType.
typedef node_type::execution_space Ifpack2::RILUK< MatrixType >::execution_space |
The Kokkos execution space of the input MatrixType.
typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::RILUK< MatrixType >::row_matrix_type |
Tpetra::RowMatrix specialization used by this class.
typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::RILUK< MatrixType >::crs_matrix_type |
Tpetra::CrsMatrix specialization used by this class for representing L and U.
typedef crs_matrix_type::impl_scalar_type Ifpack2::RILUK< MatrixType >::impl_scalar_type |
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector)
Ifpack2::RILUK< MatrixType >::RILUK | ( | const Teuchos::RCP< const row_matrix_type > & | A_in | ) |
Constructor that takes a Tpetra::RowMatrix.
A_in | [in] The input matrix. |
Ifpack2::RILUK< MatrixType >::RILUK | ( | const Teuchos::RCP< const crs_matrix_type > & | A_in | ) |
Constructor that takes a Tpetra::CrsMatrix.
This constructor exists to avoid "ambiguous constructor" warnings. It does the same thing as the constructor that takes a Tpetra::RowMatrix.
A_in | [in] The input matrix. |
|
virtual |
Destructor (declared virtual for memory safety).
|
virtual |
Set parameters for the incomplete factorization.
This preconditioner supports the following parameters:
|
virtual |
Initialize by computing the symbolic incomplete factorization.
|
virtual |
Compute the (numeric) incomplete factorization.
This function computes the RILU(k) factors L and U using the current:
initialize() must be called first, before this method may be called.
|
inlinevirtual |
Whether initialize() has been called on this object.
|
inlinevirtual |
Whether compute() has been called on this object.
|
inlinevirtual |
Number of successful initialize() calls for this object.
|
inlinevirtual |
Number of successful compute() calls for this object.
|
inlinevirtual |
Number of successful apply() calls for this object.
|
inlinevirtual |
Total time in seconds taken by all successful initialize() calls for this object.
|
inlinevirtual |
Total time in seconds taken by all successful compute() calls for this object.
|
inlinevirtual |
Total time in seconds taken by all successful apply() calls for this object.
size_t Ifpack2::RILUK< MatrixType >::getNodeSmootherComplexity | ( | ) | const |
Get a rough estimate of cost per iteration.
|
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.
std::string Ifpack2::RILUK< MatrixType >::description | ( | ) | const |
A one-line description of this object.
|
virtual |
Returns the Tpetra::Map object associated with the domain of this operator.
|
virtual |
Returns the Tpetra::Map object associated with the range of this operator.
|
virtual |
Apply the (inverse of the) incomplete factorization to X, resulting in Y.
For an incomplete factorization \(A \approx LDU\), this method computes the following, depending on the value of mode:
Y = beta*Y + alpha*(U \ (D \ (L \ X)))
Y = beta*Y + alpha*(L^T \ (D^T \ (U^T \ X)))
Y = beta*Y + alpha*(L^* \ (D^* \ (U^* \ X)))
, where the asterisk indicates the conjugate transpose. If alpha is zero, then the result of applying the operator to a vector is ignored. This matters because zero times NaN (not a number) is NaN, not zero. Analogously, if beta is zero, then any values in Y on input are ignored.
X | [in] The input multivector. |
Y | [in/out] The output multivector. |
mode | [in] If Teuchos::TRANS resp. Teuchos::CONJ_TRANS, apply the transpose resp. conjugate transpose of the incomplete factorization. Otherwise, don't apply the tranpose. |
alpha | [in] Scaling factor for the result of applying the preconditioner. |
beta | [in] Scaling factor for the initial value of Y. |
|
virtual |
|
inline |
Get RILU(k) relaxation parameter.
|
inline |
Get absolute threshold value.
|
inline |
Get relative threshold value.
|
inline |
Get level of fill (the "k" in ILU(k)).
|
inline |
Get overlap mode type.
|
inline |
Returns the number of nonzero entries in the global graph.
|
inline |
Return the Ifpack2::IlukGraph associated with this factored matrix.
const RILUK< MatrixType >::crs_matrix_type & Ifpack2::RILUK< MatrixType >::getL | ( | ) | const |
Return the L factor of the ILU factorization.
const Tpetra::Vector< typename RILUK< MatrixType >::scalar_type, typename RILUK< MatrixType >::local_ordinal_type, typename RILUK< MatrixType >::global_ordinal_type, typename RILUK< MatrixType >::node_type > & Ifpack2::RILUK< MatrixType >::getD | ( | ) | const |
Return the diagonal entries of the ILU factorization.
const RILUK< MatrixType >::crs_matrix_type & Ifpack2::RILUK< MatrixType >::getU | ( | ) | const |
Return the U factor of the ILU factorization.
Teuchos::RCP< const typename RILUK< MatrixType >::crs_matrix_type > Ifpack2::RILUK< MatrixType >::getCrsMatrix | ( | ) | const |
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
|
static |
Return A, wrapped in a LocalFilter, if necessary.
"If necessary" means that if A is already a LocalFilter, or if its communicator only has one process, then we don't need to wrap it, so we just return A.
|
protected |
The (original) input matrix for which to compute ILU(k).
|
protected |
The ILU(k) graph.
|
protected |
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_type that initialize() constructs temporarily.
|
protected |
The L (lower triangular) factor of ILU(k).
|
protected |
Sparse triangular solver for L.
|
protected |
The U (upper triangular) factor of ILU(k).
|
protected |
Sparse triangular solver for U.
|
protected |
The diagonal entries of the ILU(k) factorization.
|
protected |
Optional KokkosKernels implementation.