45 #ifndef __IFPACK2_FILDL_DEF_HPP__
46 #define __IFPACK2_FILDL_DEF_HPP__
48 #include "Ifpack2_Details_Fildl_decl.hpp"
50 #include <impl/Kokkos_Timer.hpp>
51 #include <shylu_fastildl.hpp>
58 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
61 FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
63 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
67 return localPrec_->getNFact();
70 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
74 return localPrec_->getNTrisol();
77 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
81 localPrec_->checkIC();
84 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
88 auto nRows = this->mat_->getNodeNumRows();
89 auto& p = this->params_;
90 localPrec_ =
Teuchos::rcp(
new LocalFILDL(this->localRowPtrs_, this->localColInds_, this->localValues_, nRows,
91 p.nFact, p.nTrisol, p.level, p.omega,
92 p.shift, p.guessFlag ? 1 : 0, p.blockSize));
93 localPrec_->initialize();
94 this->initTime_ = localPrec_->getInitializeTime();
97 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
102 localPrec_->setValues(this->localValues_);
103 localPrec_->compute();
104 this->computeTime_ = localPrec_->getComputeTime();
107 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
111 localPrec_->apply(x, y);
113 this->applyTime_ += localPrec_->getApplyTime();
116 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
123 #define IFPACK2_DETAILS_FILDL_INSTANT(S, L, G, N) \
124 template class Ifpack2::Details::Fildl<S, L, G, N>;
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition: Ifpack2_Details_Fildl_def.hpp:86
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition: Ifpack2_Details_Fildl_def.hpp:118
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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_Fildl_def.hpp:109
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:66
Fildl(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_Fildl_def.hpp:60
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition: Ifpack2_Details_Fildl_def.hpp:99
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition: Ifpack2_Details_Fildl_def.hpp:72
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition: Ifpack2_Details_Fildl_def.hpp:65
Kokkos::View< Scalar *, execution_space > ScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:88
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition: Ifpack2_Details_Fildl_def.hpp:79