12#ifndef __IFPACK2_FIC_DEF_HPP__
13#define __IFPACK2_FIC_DEF_HPP__
15#include "Ifpack2_Details_Fic_decl.hpp"
17#include <Kokkos_Timer.hpp>
18#include <shylu_fastic.hpp>
23template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
25 Fic(Teuchos::RCP<const TRowMatrix> A)
26 :
FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
28template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
31 return localPrec_->getNFact();
34template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
37 return localPrec_->getSpTrsvType();
40template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
43 return localPrec_->getNTrisol();
46template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
49 localPrec_->checkIC();
52template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
55 auto nRows = this->mat_->getLocalNumRows();
56 auto&
p = this->params_;
57 localPrec_ = Teuchos::rcp(
new LocalFIC(this->localRowPtrs_, this->localColInds_, this->localValues_,
nRows, (
p.sptrsv_algo != FastILU::SpTRSV::Fast),
58 p.nFact,
p.nTrisol,
p.level,
p.omega,
p.shift,
p.guessFlag ? 1 : 0,
p.blockSizeILU,
p.blockSize));
59 localPrec_->initialize();
60 this->initTime_ = localPrec_->getInitializeTime();
63template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
67 localPrec_->setValues(this->localValues_);
68 localPrec_->compute();
69 this->computeTime_ = localPrec_->getComputeTime();
72template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
75 localPrec_->apply(x, y);
77 this->applyTime_ += localPrec_->getApplyTime();
80template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
86#define IFPACK2_DETAILS_FIC_INSTANT(S, L, G, N) \
87 template class Ifpack2::Details::Fic<S, L, G, N>;
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition Ifpack2_Details_FastILU_Base_decl.hpp:38
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition Ifpack2_Details_Fic_def.hpp:42
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition Ifpack2_Details_Fic_def.hpp:48
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition Ifpack2_Details_Fic_def.hpp:30
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition Ifpack2_Details_Fic_def.hpp:65
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition Ifpack2_Details_Fic_def.hpp:36
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:74
Fic(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_Fic_def.hpp:25
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition Ifpack2_Details_Fic_def.hpp:54
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition Ifpack2_Details_Fic_def.hpp:82
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40