Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_Filu_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
11
12#ifndef __IFPACK2_FILU_DEF_HPP__
13#define __IFPACK2_FILU_DEF_HPP__
14
15#include "Ifpack2_Details_Filu_decl.hpp"
17#include "Ifpack2_Details_getCrsMatrix.hpp"
18#include <Kokkos_Timer.hpp>
19#include <shylu_fastilu.hpp>
20
21namespace Ifpack2
22{
23namespace Details
24{
25
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) {}
30
31template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
33getSweeps() const
34{
35 return localPrec_->getNFact();
36}
37
38template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
40getSpTrsvType() const
41{
42 return localPrec_->getSpTrsvType();
43}
44
45template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
47getNTrisol() const
48{
49 return localPrec_->getNTrisol();
50}
51
52template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
58
59template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
65
66template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
69{
70 auto nRows = this->mat_->getLocalNumRows();
71 auto& p = this->params_;
72 auto matCrs = Ifpack2::Details::getCrsMatrix(this->mat_);
73
74 if (p.blockCrsSize > 1 && !BlockCrsEnabled) {
75 throw std::runtime_error("Must use prec type FAST_ILU_B if you want blockCrs support");
76 }
77
78 bool skipSortMatrix = !matCrs.is_null() && matCrs->getCrsGraph()->isSorted() &&
79 !p.use_metis;
80 localPrec_ =
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,
83 p.blockCrsSize));
84
85 #ifdef HAVE_IFPACK2_METIS
86 if (p.use_metis) {
87 localPrec_->setMetisPerm(this->metis_perm_, this->metis_iperm_);
88 }
89 #endif
90
91 localPrec_->initialize();
92 this->initTime_ = localPrec_->getInitializeTime();
93}
94
95template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
98{
99 //update values in local prec (until compute(), values aren't needed)
100 localPrec_->setValues(this->localValues_);
101 localPrec_->compute();
102 this->computeTime_ = localPrec_->getComputeTime();
103}
104
105template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
107applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const
108{
109 localPrec_->apply(x, y);
110
111 //since this may be applied to multiple vectors, add to applyTime_ instead of setting it
112 this->applyTime_ += localPrec_->getApplyTime();
113}
114
115template<typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node, bool BlockCrsEnabled>
117getName() const
118{
119 return "Filu";
120}
121
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>;
125
126} //namespace Details
127} //namespace Ifpack2
128
129#endif
130
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