Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_Details_Fildl_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_FILDL_DEF_HPP__
13 #define __IFPACK2_FILDL_DEF_HPP__
14 
15 #include "Ifpack2_Details_Fildl_decl.hpp"
17 #include <Kokkos_Timer.hpp>
18 #include <shylu_fastildl.hpp>
19 
20 namespace Ifpack2 {
21 namespace Details {
22 
23 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
26  : FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A) {}
27 
28 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
30  getSweeps() const {
31  return localPrec_->getNFact();
32 }
33 
34 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
36  getSpTrsvType() const {
37  return localPrec_->getSpTrsvType();
38 }
39 
40 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
42  getNTrisol() const {
43  return localPrec_->getNTrisol();
44 }
45 
46 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
48  checkLocalIC() const {
49  localPrec_->checkIC();
50 }
51 
52 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
55  auto nRows = this->mat_->getLocalNumRows();
56  auto& p = this->params_;
57  localPrec_ = Teuchos::rcp(new LocalFILDL(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();
61 }
62 
63 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
66  // update values in local prec (until compute(), values aren't needed)
67  localPrec_->setValues(this->localValues_);
68  localPrec_->compute();
69  this->computeTime_ = localPrec_->getComputeTime();
70 }
71 
72 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
75  localPrec_->apply(x, y);
76  // since this may be applied to multiple vectors, add to applyTime_ instead of setting it
77  this->applyTime_ += localPrec_->getApplyTime();
78 }
79 
80 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
82  getName() const {
83  return "Fildl";
84 }
85 
86 #define IFPACK2_DETAILS_FILDL_INSTANT(S, L, G, N) \
87  template class Ifpack2::Details::Fildl<S, L, G, N>;
88 
89 } // namespace Details
90 } // namespace Ifpack2
91 
92 #endif
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_-&gt;ini...
Definition: Ifpack2_Details_Fildl_def.hpp:54
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition: Ifpack2_Details_Fildl_def.hpp:36
std::string getName() const
Get the name of the underlying preconditioner (&quot;Filu&quot;, &quot;Fildl&quot; or &quot;Fic&quot;)
Definition: Ifpack2_Details_Fildl_def.hpp:82
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:36
Fildl(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition: Ifpack2_Details_Fildl_def.hpp:25
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition: Ifpack2_Details_Fildl_def.hpp:65
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition: Ifpack2_Details_Fildl_def.hpp:42
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition: Ifpack2_Details_Fildl_def.hpp:30
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition: Ifpack2_Details_FastILU_Base_decl.hpp:59
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:74
Provides functions for retrieving local CRS arrays (row pointers, column indices, and values) from Tp...
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition: Ifpack2_Details_Fildl_def.hpp:48