48 #ifndef __IFPACK2_FASTILU_BASE_DECL_HPP__  
   49 #define __IFPACK2_FASTILU_BASE_DECL_HPP__  
   51 #include <Tpetra_RowMatrix.hpp> 
   52 #include <Tpetra_CrsMatrix.hpp> 
   53 #include <Kokkos_DefaultNode.hpp> 
   54 #include <KokkosSparse_CrsMatrix.hpp> 
   65 template<
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
   68                          Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
 
   76     typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> 
TRowMatrix;
 
   78     typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> 
TCrsMatrix;
 
   80     typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> 
TMultiVec;
 
   82     typedef KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, execution_space> 
KCrsMatrix;
 
   84     typedef Kokkos::View<LocalOrdinal *, execution_space> 
OrdinalArray;
 
  194     mutable double applyTime_;
 
  209       static Params getDefaults();
 
  222     virtual std::string 
getName() 
const = 0;
 
virtual std::string getName() const =0
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic") 
 
Mix-in interface for preconditioners that can change their matrix after construction. 
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
 
int getNumCompute() const 
Get the number of times compute() was called. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:231
 
virtual void checkLocalIC() const 
Verify and print debug information about the underlying IC preconditioner. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:281
 
double getComputeTime() const 
Get the time spent in the last compute() call. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:252
 
virtual int getNTrisol() const =0
Get the "triangular solve iterations" parameter. 
 
device_type::execution_space execution_space
Kokkos execution space. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:74
 
double getCopyTime() const 
Get the time spent deep copying local 3-array CRS out of the matrix. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:266
 
virtual int getSweeps() const =0
Get the "sweeps" parameter. 
 
double getInitializeTime() const 
Get the time spent in the last initialize() call. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:245
 
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const 
Get the range map of the matrix. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:85
 
double getApplyTime() const 
Get the time spent in the last apply() call. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:259
 
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TCrsMatrix
Tpetra CRS matrix. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:78
 
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const 
Get the domain map of the matrix. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:77
 
Teuchos::RCP< const TRowMatrix > getMatrix() const 
Get the current matrix. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:217
 
Kokkos::View< Scalar *, Kokkos::HostSpace > ScalarArrayHost
Array of Scalar on host. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:90
 
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:273
 
Kokkos::View< LocalOrdinal *, execution_space > OrdinalArray
Array of LocalOrdinal on device. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:84
 
int getNumApply() const 
Get the number of times apply() was called. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:238
 
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic) 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:66
 
void setParameters(const Teuchos::ParameterList &List)
Validate parameters, and set defaults when parameters are not provided. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:146
 
std::string description() const 
Return a brief description of the preconditioner, in YAML format. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:288
 
void compute()
Compute the preconditioner. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:184
 
Node::device_type device_type
Kokkos device type. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:72
 
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:107
 
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:61
 
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:209
 
Declaration of interface for preconditioners that can change their matrix after construction. 
 
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:177
 
KokkosSparse::CrsMatrix< Scalar, LocalOrdinal, execution_space > KCrsMatrix
Kokkos CRS matrix. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:82
 
void initialize()
Initialize the preconditioner. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:154
 
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TMultiVec
Tpetra multivector. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:80
 
virtual void applyLocalPrec(ScalarArray x, ScalarArray y) const =0
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only) ...
 
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TRowMatrix
Tpetra row matrix. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:76
 
Kokkos::View< LocalOrdinal *, Kokkos::HostSpace > OrdinalArrayHost
Array of LocalOrdinal on host. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:86
 
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition: Ifpack2_Details_FastILU_Base_def.hpp:311
 
Kokkos::View< Scalar *, execution_space > ScalarArray
Array of Scalar on device. 
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:88
 
virtual void initLocalPrec()=0
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
 
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:92
 
int getNumInitialize() const 
Get the number of times initialize() was called. 
Definition: Ifpack2_Details_FastILU_Base_def.hpp:224
 
virtual void computeLocalPrec()=0
Get values array from the matrix and then call compute() on the underlying preconditioner.