Ifpack2 Templated Preconditioning Package
Version 1.0
|
MDF (incomplete LU factorization with minimum discarded fill reordering) of a Tpetra sparse matrix. More...
#include <Ifpack2_MDF_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... | |
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... | |
const crs_matrix_type & | getL () const |
Return the L factor of the MDF factorization. More... | |
const crs_matrix_type & | getU () const |
Return the U factor of the MDF factorization. More... | |
permutations_type & | getPermutations () const |
Return the permutations of the MDF factorization. More... | |
permutations_type & | getReversePermutations () const |
Return the reverse permutations of the MDF 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... | |
Protected Attributes | |
Teuchos::RCP< const row_matrix_type > | A_ |
The (original) input matrix for which to compute ILU(k). 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... | |
permutations_type | permutations_ |
The computed permuations from MDF factorization. More... | |
permutations_type | reversePermutations_ |
The reverse permuations from MDF factorization. More... | |
Implementation of Kokkos Kernels MDF. | |
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 |
MDF (const Teuchos::RCP< const row_matrix_type > &A_in) | |
Constructor that takes a Tpetra::RowMatrix. More... | |
MDF (const Teuchos::RCP< const crs_matrix_type > &A_in) | |
Constructor that takes a Tpetra::CrsMatrix. More... | |
virtual | ~MDF ()=default |
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... | |
MDF (incomplete LU factorization with minimum discarded fill reordering) of a Tpetra sparse matrix.
A | specialization of Tpetra::RowMatrix. |
This class computes a sparse MDF (incomplete LU) factorization with a reordering that minimizes the discarded fill 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 MDF as the subdomain solver on each process.
typedef MatrixType::scalar_type Ifpack2::MDF< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::MDF< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::MDF< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type Ifpack2::MDF< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::MDF< MatrixType >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
typedef node_type::device_type Ifpack2::MDF< MatrixType >::device_type |
The Kokkos device type of the input MatrixType.
typedef node_type::execution_space Ifpack2::MDF< 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::MDF< 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::MDF< MatrixType >::crs_matrix_type |
Tpetra::CrsMatrix specialization used by this class for representing L and U.
typedef crs_matrix_type::impl_scalar_type Ifpack2::MDF< MatrixType >::impl_scalar_type |
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector)
Ifpack2::MDF< MatrixType >::MDF | ( | const Teuchos::RCP< const row_matrix_type > & | A_in | ) |
Constructor that takes a Tpetra::RowMatrix.
A_in | [in] The input matrix. |
Ifpack2::MDF< MatrixType >::MDF | ( | 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. |
|
virtualdefault |
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::MDF< 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::MDF< 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 level of fill (the "k" in ILU(k)).
|
inline |
Get overlap mode type.
|
inline |
Returns the number of nonzero entries in the global graph.
const MDF< MatrixType >::crs_matrix_type & Ifpack2::MDF< MatrixType >::getL | ( | ) | const |
Return the L factor of the MDF factorization.
const MDF< MatrixType >::crs_matrix_type & Ifpack2::MDF< MatrixType >::getU | ( | ) | const |
Return the U factor of the MDF factorization.
MDF< MatrixType >::permutations_type & Ifpack2::MDF< MatrixType >::getPermutations | ( | ) | const |
Return the permutations of the MDF factorization.
MDF< MatrixType >::permutations_type & Ifpack2::MDF< MatrixType >::getReversePermutations | ( | ) | const |
Return the reverse permutations of the MDF factorization.
Teuchos::RCP< const typename MDF< MatrixType >::crs_matrix_type > Ifpack2::MDF< MatrixType >::getCrsMatrix | ( | ) | const |
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
|
protected |
The (original) input matrix for which to compute ILU(k).
|
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 computed permuations from MDF factorization.
|
protected |
The reverse permuations from MDF factorization.