12#ifndef __IFPACK2_FILU_DEF_HPP__
13#define __IFPACK2_FILU_DEF_HPP__
15#include "Ifpack2_Details_Filu_decl.hpp"
17#include "Ifpack2_Details_getCrsMatrix.hpp"
18#include <Kokkos_Timer.hpp>
19#include <shylu_fastilu.hpp>
26template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
28Filu(Teuchos::RCP<const TRowMatrix> A) :
29 FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
31template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
35 return localPrec_->getNFact();
38template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
42 return localPrec_->getSpTrsvType();
45template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
49 return localPrec_->getNTrisol();
52template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
56 localPrec_->checkILU();
59template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
63 localPrec_->checkIC();
66template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
70 auto nRows = this->mat_->getLocalNumRows();
71 auto& p = this->params_;
72 auto matCrs = Ifpack2::Details::getCrsMatrix(this->mat_);
74 if (p.blockCrsSize > 1 && !BlockCrsEnabled) {
75 throw std::runtime_error(
"Must use prec type FAST_ILU_B if you want blockCrs support");
78 bool skipSortMatrix = !matCrs.is_null() && matCrs->getCrsGraph()->isSorted() &&
81 Teuchos::rcp(
new LocalFILU(skipSortMatrix, this->localRowPtrs_, this->localColInds_, this->localValues_, nRows, p.sptrsv_algo,
82 p.nFact, p.nTrisol, p.level, p.omega, p.shift, p.guessFlag ? 1 : 0, p.blockSizeILU, p.blockSize,
85 #ifdef HAVE_IFPACK2_METIS
87 localPrec_->setMetisPerm(this->metis_perm_, this->metis_iperm_);
91 localPrec_->initialize();
92 this->initTime_ = localPrec_->getInitializeTime();
95template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
100 localPrec_->setValues(this->localValues_);
101 localPrec_->compute();
102 this->computeTime_ = localPrec_->getComputeTime();
105template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
109 localPrec_->apply(x, y);
112 this->applyTime_ += localPrec_->getApplyTime();
115template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node,
bool BlockCrsEnabled>
122#define IFPACK2_DETAILS_FILU_INSTANT(S, L, G, N) \
123template class Ifpack2::Details::Filu<S, L, G, N,false>; \
124template class Ifpack2::Details::Filu<S, L, G, N,true>;
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_FastILU_Base_def.hpp:31
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic").
Definition Ifpack2_Details_Filu_def.hpp:117
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition Ifpack2_Details_Filu_def.hpp:61
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition Ifpack2_Details_Filu_def.hpp:47
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition Ifpack2_Details_Filu_def.hpp:68
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition Ifpack2_Details_Filu_def.hpp:97
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition Ifpack2_Details_Filu_def.hpp:33
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition Ifpack2_Details_Filu_def.hpp:40
void checkLocalILU() const
Verify and print debug info about the internal ILU preconditioner.
Definition Ifpack2_Details_Filu_def.hpp:54
Filu(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_Filu_def.hpp:28
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_Filu_def.hpp:107
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41