Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_Fildl_decl.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_DECL_HPP__
13#define __IFPACK2_FILDL_DECL_HPP__
14
15#include <Ifpack2_Details_FastILU_Base.hpp>
16
17// forward-declare the local preconditioner type
18template <typename LocalOrdinal, typename Scalar, typename execution_space>
19class FastILDLPrec;
20
21namespace Ifpack2 {
22namespace Details {
23
26template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
27class Fildl : public FastILU_Base<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
28 public:
30 typedef typename Base::TRowMatrix TRowMatrix;
31 typedef typename Base::ImplScalar ImplScalar;
32 typedef typename Base::ImplScalarArray ImplScalarArray;
34
36 Fildl(Teuchos::RCP<const TRowMatrix> mat_);
37
39 int getSweeps() const;
40
42 std::string getSpTrsvType() const;
43
45 int getNTrisol() const;
46
48 void checkLocalIC() const;
49
50 protected:
51 mutable Teuchos::RCP<LocalFILDL> localPrec_;
52
53 void initLocalPrec();
54 // compute() takes A's local values
55 void computeLocalPrec();
56 void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const;
57 std::string getName() const;
58};
59
60} // namespace Details
61} // namespace Ifpack2
62
63#endif
The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)
Definition Ifpack2_Details_FastILU_Base_decl.hpp:38
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type ImplScalar
Kokkos scalar type.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:45
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TRowMatrix
Tpetra row matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:47
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:59
The Ifpack2 wrapper to the ILDL preconditioner of ShyLU FastILU.
Definition Ifpack2_Details_Fildl_decl.hpp:27
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
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition Ifpack2_Details_Fildl_def.hpp:65
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition Ifpack2_Details_Fildl_def.hpp:30
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Definition Ifpack2_Details_Fildl_def.hpp:82
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition Ifpack2_Details_Fildl_def.hpp:36
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition Ifpack2_Details_Fildl_def.hpp:42
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition Ifpack2_Details_Fildl_def.hpp:48
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition Ifpack2_Details_Fildl_def.hpp:54
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