46 #ifndef IFPACK2_CRSRILUK_DECL_HPP 47 #define IFPACK2_CRSRILUK_DECL_HPP 51 #include "Tpetra_CrsMatrix_decl.hpp" 54 #include "Ifpack2_LocalSparseTriangularSolver_decl.hpp" 56 #include <type_traits> 242 template<
class MatrixType>
245 typename MatrixType::local_ordinal_type,
246 typename MatrixType::global_ordinal_type,
247 typename MatrixType::node_type>,
249 typename MatrixType::local_ordinal_type,
250 typename MatrixType::global_ordinal_type,
251 typename MatrixType::node_type> >
267 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
276 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
"Ifpack2::RILUK: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.");
284 template <
class NewMatrixType>
friend class RILUK;
289 RILUK (
const Teuchos::RCP<const row_matrix_type>& A_in);
298 RILUK (
const Teuchos::RCP<const crs_matrix_type>& A_in);
334 return isInitialized_;
343 return numInitialize_;
356 return initializeTime_;
398 setMatrix (
const Teuchos::RCP<const row_matrix_type>& A);
412 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
416 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
449 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
450 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
451 Teuchos::ETransp mode = Teuchos::NO_TRANS,
452 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one (),
453 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero ())
const;
479 multiply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
480 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
481 const Teuchos::ETransp mode = Teuchos::NO_TRANS)
const;
485 Teuchos::RCP<const row_matrix_type>
getMatrix ()
const;
503 TEUCHOS_TEST_FOR_EXCEPTION(
504 true, std::logic_error,
"Ifpack2::RILUK::SetOverlapMode: " 505 "RILUK no longer implements overlap on its own. " 506 "Use RILUK with AdditiveSchwarz if you want overlap.");
511 return getL ().getGlobalNumEntries () +
getU ().getGlobalNumEntries ();
515 Teuchos::RCP<Ifpack2::IlukGraph<Tpetra::CrsGraph<local_ordinal_type,global_ordinal_type,node_type> > >
getGraph ()
const {
523 const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>&
530 Teuchos::RCP<const crs_matrix_type>
getCrsMatrix ()
const;
533 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
534 typedef Teuchos::ScalarTraits<scalar_type> STS;
535 typedef Teuchos::ScalarTraits<magnitude_type> STM;
537 void allocateSolvers ();
538 void allocate_L_and_U ();
547 static Teuchos::RCP<const row_matrix_type>
548 makeLocalFilter (
const Teuchos::RCP<const row_matrix_type>& A);
551 typedef Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> vec_type;
554 Teuchos::RCP<const row_matrix_type>
A_;
566 Teuchos::RCP<crs_matrix_type>
L_;
568 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> >
L_solver_;
570 Teuchos::RCP<crs_matrix_type>
U_;
572 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> >
U_solver_;
574 Teuchos::RCP<vec_type>
D_;
585 mutable int numApply_;
587 double initializeTime_;
589 mutable double applyTime_;
604 template<
class MatrixType,
class NodeType>
605 struct setLocalSolveParams{
606 static Teuchos::RCP<Teuchos::ParameterList>
607 setParams (
const Teuchos::RCP<Teuchos::ParameterList>& param) {
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
Declaration and definition of IlukGraph.
Ifpack2::ScalingType enumerable type.
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:359
int getNumApply() const
Number of successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:350
Ifpack2::Preconditioner< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > ::scalar_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > ::local_ordinal_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > ::global_ordinal_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > ::node_type >::magnitude_type Teuchos::ScalarTraits< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > ::scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Preconditioner.hpp:111
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:261
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:566
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:396
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition: Ifpack2_RILUK_def.hpp:349
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_RILUK_decl.hpp:273
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
Definition: Ifpack2_RILUK_def.hpp:356
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_RILUK_def.hpp:107
void setParameters(const Teuchos::ParameterList ¶ms)
Definition: Ifpack2_RILUK_def.hpp:271
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:243
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition: Ifpack2_RILUK_decl.hpp:502
bool isComputed() const
Whether compute() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:337
double getApplyTime() const
Total time in seconds taken by all successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:363
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.
Definition: Ifpack2_RILUK_decl.hpp:276
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack2_RILUK_decl.hpp:493
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > getGraph() const
Return the Ifpack2::IlukGraph associated with this factored matrix.
Definition: Ifpack2_RILUK_decl.hpp:515
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:264
const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > & getD() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:160
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:258
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_RILUK_def.hpp:96
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_RILUK_def.hpp:187
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition: Ifpack2_RILUK_decl.hpp:572
int getNumCompute() const
Number of successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:346
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.
Definition: Ifpack2_RILUK_def.hpp:205
magnitude_type getRelaxValue() const
Get RILU(k) relaxation parameter.
Definition: Ifpack2_RILUK_decl.hpp:490
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > Graph_
The ILU(k) graph.
Definition: Ifpack2_RILUK_decl.hpp:559
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:107
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:653
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:174
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_RILUK_decl.hpp:267
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:97
Declaration of interface for preconditioners that can change their matrix after construction.
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition: Ifpack2_RILUK_decl.hpp:568
bool isInitialized() const
Whether initialize() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:333
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:554
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack2_RILUK_decl.hpp:510
Definition: Ifpack2_Container_decl.hpp:565
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition: Ifpack2_RILUK_decl.hpp:574
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.
Definition: Ifpack2_RILUK_def.hpp:835
std::string description() const
A one-line description of this object.
Definition: Ifpack2_RILUK_def.hpp:983
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition: Ifpack2_RILUK_decl.hpp:496
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:342
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:73
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition: Ifpack2_RILUK_decl.hpp:499
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object. ...
Definition: Ifpack2_RILUK_decl.hpp:355
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_t...
Definition: Ifpack2_RILUK_decl.hpp:563
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:570
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.
Definition: Ifpack2_RILUK_def.hpp:224
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:143
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:255