45 #ifndef __IFPACK2_FILU_DECL_HPP__
46 #define __IFPACK2_FILU_DECL_HPP__
48 #include <Ifpack2_Details_FastILU_Base.hpp>
51 template<
typename LocalOrdinal,
typename Scalar,
typename execution_space>
61 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
68 typedef FastILUPrec<LocalOrdinal, Scalar, typename Base::execution_space> LocalFILU;
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:86
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition: Ifpack2_Details_Filu_def.hpp:93
The Ifpack2 wrapper to the ILU preconditioner of ShyLU FastILU.
Definition: Ifpack2_Details_Filu_decl.hpp:62
void checkLocalILU() const
Verify and print debug info about the internal ILU preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:79
Filu(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_Filu_def.hpp:60
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:66
void applyLocalPrec(ScalarArray x, ScalarArray y) const
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only) ...
Definition: Ifpack2_Details_Filu_def.hpp:116
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TRowMatrix
Tpetra row matrix.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:76
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition: Ifpack2_Details_Filu_def.hpp:65
Kokkos::View< Scalar *, execution_space > ScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:88
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition: Ifpack2_Details_Filu_def.hpp:125
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition: Ifpack2_Details_Filu_def.hpp:72
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition: Ifpack2_Details_Filu_def.hpp:106