15#ifndef __IFPACK2_FASTILU_BASE_DECL_HPP__
16#define __IFPACK2_FASTILU_BASE_DECL_HPP__
18#include <Tpetra_RowMatrix.hpp>
19#include <Tpetra_CrsMatrix.hpp>
20#include <Tpetra_KokkosCompat_DefaultNode.hpp>
21#include <KokkosSparse_CrsMatrix.hpp>
24#include <shylu_fastutil.hpp>
26#ifdef HAVE_IFPACK2_METIS
37template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40 Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
48 typedef typename Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::impl_scalar_type
ImplScalar;
50 typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
TRowMatrix;
52 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
TCrsMatrix;
54 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
TMultiVec;
56 typedef KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, execution_space>
KCrsMatrix;
58 typedef Kokkos::View<LocalOrdinal *, execution_space>
OrdinalArray;
60 typedef typename Kokkos::View<LocalOrdinal *, execution_space>::HostMirror
OrdinalArrayHost;
63 typedef Kokkos::View< Scalar *, execution_space> ScalarArray;
64 typedef Kokkos::View<const Scalar *, execution_space> ConstScalarArray;
65 #ifdef HAVE_IFPACK2_METIS
66 typedef Kokkos::View<idx_t*, Kokkos::HostSpace> MetisArrayHost;
73 Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
77 Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
84 Teuchos::ETransp mode = Teuchos::NO_TRANS,
85 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
86 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero())
const;
101 bool isBlockCrs()
const;
116 Teuchos::RCP<const TRowMatrix>
getMatrix()
const;
160 void setMatrix(
const Teuchos::RCP<const TRowMatrix>& A);
163 Teuchos::RCP<const TRowMatrix> mat_;
177 mutable double applyTime_;
184 Params(
const Teuchos::ParameterList& pL, std::string precType);
186 FastILU::SpTRSV sptrsv_algo;
199 static Params getDefaults();
204 #ifdef HAVE_IFPACK2_METIS
205 MetisArrayHost metis_perm_;
206 MetisArrayHost metis_iperm_;
Declaration of interface for preconditioners that can change their matrix after construction.
Mix-in interface for preconditioners that can change their matrix after construction.
Definition Ifpack2_Details_CanChangeMatrix.hpp:60
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type ImplScalar
Kokkos scalar type.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:48
double getInitializeTime() const
Get the time spent in the last initialize() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:319
virtual void checkLocalILU() const
Verify and print debug information about the underlying ILU preconditioner (only supported if this is...
Definition Ifpack2_Details_FastILU_Base_def.hpp:347
Node::device_type device_type
Kokkos device type.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:44
virtual std::string getSpTrsvType() const =0
Get the name of triangular solve algorithm.
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition Ifpack2_Details_FastILU_Base_def.hpp:388
virtual void computeLocalPrec()=0
Get values array from the matrix and then call compute() on the underlying preconditioner.
void compute()
Compute the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:258
double getComputeTime() const
Get the time spent in the last compute() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:326
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TCrsMatrix
Tpetra CRS matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:52
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TRowMatrix
Tpetra row matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:50
device_type::execution_space execution_space
Kokkos execution space.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:46
double getCopyTime() const
Get the time spent deep copying local 3-array CRS out of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:340
int getNumApply() const
Get the number of times apply() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:312
void initialize()
Initialize the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:135
KokkosSparse::CrsMatrix< Scalar, LocalOrdinal, execution_space > KCrsMatrix
Kokkos CRS matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:56
bool isInitialized() const
Whether initialize() has been called since the last time the matrix's structure was changed.
Definition Ifpack2_Details_FastILU_Base_def.hpp:251
virtual std::string getName() const =0
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic").
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_FastILU_Base_def.hpp:31
virtual int getSweeps() const =0
Get the "sweeps" parameter.
void setParameters(const Teuchos::ParameterList &List)
Validate parameters, and set defaults when parameters are not provided.
Definition Ifpack2_Details_FastILU_Base_def.hpp:120
Kokkos::View< LocalOrdinal *, execution_space > OrdinalArray
Array of LocalOrdinal on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:58
void apply(const TMultiVec &X, TMultiVec &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Apply the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:62
double getApplyTime() const
Get the time spent in the last apply() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:333
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Get the domain map of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:47
virtual void checkLocalIC() const
Verify and print debug information about the underlying IC preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:355
virtual void initLocalPrec()=0
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
int getNumInitialize() const
Get the number of times initialize() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:298
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TMultiVec
Tpetra multivector.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:54
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Get the range map of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:55
Teuchos::RCP< const TRowMatrix > getMatrix() const
Get the current matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:291
Kokkos::View< LocalOrdinal *, execution_space >::HostMirror OrdinalArrayHost
Array of LocalOrdinal on host.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:60
std::string description() const
Return a brief description of the preconditioner, in YAML format.
Definition Ifpack2_Details_FastILU_Base_def.hpp:362
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:62
virtual int getNTrisol() const =0
Get the "triangular solve iterations" parameter.
int getNumCompute() const
Get the number of times compute() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:305
bool isComputed() const
Whether compute() has been called since the last time the matrix's values or structure were changed.
Definition Ifpack2_Details_FastILU_Base_def.hpp:283
virtual void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const =0
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only).
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:75
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41