45 #ifndef __IFPACK2_FIC_DEF_HPP__
46 #define __IFPACK2_FIC_DEF_HPP__
48 #include "Ifpack2_Details_Fic_decl.hpp"
50 #include <Kokkos_Timer.hpp>
51 #include <shylu_fastic.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_->getSpTrsvType();
77 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
81 return localPrec_->getNTrisol();
84 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
88 localPrec_->checkIC();
91 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
95 auto nRows = this->mat_->getLocalNumRows();
96 auto& p = this->params_;
97 localPrec_ =
Teuchos::rcp(
new LocalFIC(this->localRowPtrs_, this->localColInds_, this->localValues_, nRows, (p.sptrsv_algo != FastILU::SpTRSV::Fast),
98 p.nFact, p.nTrisol, p.level, p.omega, p.shift, p.guessFlag ? 1 : 0, p.blockSizeILU, p.blockSize));
99 localPrec_->initialize();
100 this->initTime_ = localPrec_->getInitializeTime();
103 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
108 localPrec_->setValues(this->localValues_);
109 localPrec_->compute();
110 this->computeTime_ = localPrec_->getComputeTime();
113 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
117 localPrec_->apply(x, y);
119 this->applyTime_ += localPrec_->getApplyTime();
122 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
129 #define IFPACK2_DETAILS_FIC_INSTANT(S, L, G, N) \
130 template class Ifpack2::Details::Fic<S, L, G, N>;
Fic(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_Fic_def.hpp:60
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition: Ifpack2_Details_Fic_def.hpp:65
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition: Ifpack2_Details_Fic_def.hpp:124
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition: Ifpack2_Details_Fic_def.hpp:93
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:71
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition: Ifpack2_Details_Fic_def.hpp:105
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition: Ifpack2_Details_Fic_def.hpp:72
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:95
void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only) ...
Definition: Ifpack2_Details_Fic_def.hpp:115
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition: Ifpack2_Details_Fic_def.hpp:79
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition: Ifpack2_Details_Fic_def.hpp:86
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...