12#ifndef __IFPACK2_FILDL_DEF_HPP__
13#define __IFPACK2_FILDL_DEF_HPP__
15#include "Ifpack2_Details_Fildl_decl.hpp"
17#include <Kokkos_Timer.hpp>
18#include <shylu_fastildl.hpp>
25template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
27Fildl(Teuchos::RCP<const TRowMatrix> A) :
28 FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
30template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
34 return localPrec_->getNFact();
37template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
41 return localPrec_->getSpTrsvType();
44template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
48 return localPrec_->getNTrisol();
51template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
55 localPrec_->checkIC();
58template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
62 auto nRows = this->mat_->getLocalNumRows();
63 auto&
p = this->params_;
64 localPrec_ = Teuchos::rcp(
new LocalFILDL(this->localRowPtrs_, this->localColInds_, this->localValues_,
nRows, (
p.sptrsv_algo != FastILU::SpTRSV::Fast),
65 p.nFact,
p.nTrisol,
p.level,
p.omega,
p.shift,
p.guessFlag ? 1 : 0,
p.blockSizeILU,
p.blockSize));
66 localPrec_->initialize();
67 this->initTime_ = localPrec_->getInitializeTime();
70template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
75 localPrec_->setValues(this->localValues_);
76 localPrec_->compute();
77 this->computeTime_ = localPrec_->getComputeTime();
80template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
84 localPrec_->apply(x, y);
86 this->applyTime_ += localPrec_->getApplyTime();
89template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
96#define IFPACK2_DETAILS_FILDL_INSTANT(S, L, G, N) \
97template class Ifpack2::Details::Fildl<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:41
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_Fildl_def.hpp:82
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition Ifpack2_Details_Fildl_def.hpp:72
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition Ifpack2_Details_Fildl_def.hpp:32
Fildl(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_Fildl_def.hpp:27
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition Ifpack2_Details_Fildl_def.hpp:91
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition Ifpack2_Details_Fildl_def.hpp:39
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition Ifpack2_Details_Fildl_def.hpp:46
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition Ifpack2_Details_Fildl_def.hpp:53
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition Ifpack2_Details_Fildl_def.hpp:60
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:77
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41