Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_FastILU_Base_decl.hpp
Go to the documentation of this file.
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
14
15#ifndef __IFPACK2_FASTILU_BASE_DECL_HPP__
16#define __IFPACK2_FASTILU_BASE_DECL_HPP__
17
18#include <Tpetra_RowMatrix.hpp>
19#include <Tpetra_CrsMatrix.hpp>
20#include <Tpetra_KokkosCompat_DefaultNode.hpp>
21#include <KokkosSparse_CrsMatrix.hpp>
24#include <shylu_fastutil.hpp>
25
26#ifdef HAVE_IFPACK2_METIS
27#include "metis.h"
28#endif
29
30namespace Ifpack2 {
31namespace Details {
32
35template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36class FastILU_Base : public Ifpack2::Preconditioner<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
38 Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> {
39 public:
41 typedef typename Node::device_type device_type;
43 typedef typename device_type::execution_space execution_space;
45 typedef typename Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::impl_scalar_type ImplScalar;
47 typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TRowMatrix;
49 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TCrsMatrix;
51 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TMultiVec;
53 typedef KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, execution_space> KCrsMatrix;
55 typedef Kokkos::View<LocalOrdinal*, execution_space> OrdinalArray;
57 typedef typename Kokkos::View<LocalOrdinal*, execution_space>::host_mirror_type OrdinalArrayHost;
59 typedef Kokkos::View<ImplScalar*, execution_space> ImplScalarArray;
60 typedef Kokkos::View<Scalar*, execution_space> ScalarArray;
61 typedef Kokkos::View<const Scalar*, execution_space> ConstScalarArray;
62#ifdef HAVE_IFPACK2_METIS
63 typedef Kokkos::View<idx_t*, Kokkos::HostSpace> MetisArrayHost;
64#endif
65
67 FastILU_Base(Teuchos::RCP<const TRowMatrix> mat_);
68
70 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
71 getDomainMap() const;
72
74 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
75 getRangeMap() const;
76
78 void
79 apply(const TMultiVec& X,
80 TMultiVec& Y,
81 Teuchos::ETransp mode = Teuchos::NO_TRANS,
82 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
83 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const;
85
96 void setParameters(const Teuchos::ParameterList& List);
97
98 bool isBlockCrs() const;
99
101 void initialize();
102
104 bool isInitialized() const;
105
107 void compute();
108
110 bool isComputed() const;
111
113 Teuchos::RCP<const TRowMatrix> getMatrix() const;
114
116 int getNumInitialize() const;
117
119 int getNumCompute() const;
120
122 int getNumApply() const;
123
125 double getInitializeTime() const;
126
128 double getComputeTime() const;
129
131 double getApplyTime() const;
132
134 double getCopyTime() const;
135
137 virtual int getSweeps() const = 0;
138
140 virtual std::string getSpTrsvType() const = 0;
141
143 virtual int getNTrisol() const = 0;
144
146 virtual void checkLocalILU() const;
147
149 virtual void checkLocalIC() const;
150
152 std::string description() const;
153
157 void setMatrix(const Teuchos::RCP<const TRowMatrix>& A);
158
159 protected:
160 Teuchos::RCP<const TRowMatrix> mat_;
161 bool initFlag_;
162 bool computedFlag_;
163 int nInit_;
164 int nComputed_;
165 mutable int nApply_;
166 // store the local CRS components
167 ImplScalarArray localValues_; // set at beginning of compute()
168 OrdinalArray localRowPtrs_; // set in initialize()
169 OrdinalArray localColInds_; // set in initialize()
170 OrdinalArrayHost localRowPtrsHost_; // set in initialize() and used to get localValues_ in compute()
171
172 double initTime_;
173 double computeTime_;
174 mutable double applyTime_;
175 double crsCopyTime_; // total time spent deep copying values, rowptrs, colinds out of mat
176
177 // Store validated parameter values (instead of keeping a ParameterList around)
178 struct Params {
179 Params() {}
180 Params(const Teuchos::ParameterList& pL, std::string precType);
181 bool use_metis;
182 FastILU::SpTRSV sptrsv_algo;
183 int nFact;
184 int nTrisol;
185 int level;
186 int blkSize;
187 double omega;
188 double shift;
189 bool guessFlag;
190 int blockSizeILU;
191 int blockSize;
192 bool blockCrs;
193 int blockCrsSize;
194 bool fillBlocks;
195 static Params getDefaults();
196 };
197
198 Params params_;
199
200#ifdef HAVE_IFPACK2_METIS
201 MetisArrayHost metis_perm_;
202 MetisArrayHost metis_iperm_;
203#endif
204
206 // \pre !mat_.is_null()
207 virtual void initLocalPrec() = 0;
209 virtual void computeLocalPrec() = 0;
211 virtual void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const = 0;
213 virtual std::string getName() const = 0;
214};
215
216} // namespace Details
217} // namespace Ifpack2
218
219#endif
Declaration of interface for preconditioners that can change their matrix after construction.
Mix-in interface for preconditioners that can change their matrix after construction.
Definition Ifpack2_Details_CanChangeMatrix.hpp:60
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
double getInitializeTime() const
Get the time spent in the last initialize() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:292
virtual void checkLocalILU() const
Verify and print debug information about the underlying ILU preconditioner (only supported if this is...
Definition Ifpack2_Details_FastILU_Base_def.hpp:316
Node::device_type device_type
Kokkos device type.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:41
virtual std::string getSpTrsvType() const =0
Get the name of triangular solve algorithm.
void setMatrix(const Teuchos::RCP< const TRowMatrix > &A)
Definition Ifpack2_Details_FastILU_Base_def.hpp:351
virtual void computeLocalPrec()=0
Get values array from the matrix and then call compute() on the underlying preconditioner.
void compute()
Compute the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:238
double getComputeTime() const
Get the time spent in the last compute() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:298
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TCrsMatrix
Tpetra CRS matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:49
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TRowMatrix
Tpetra row matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:47
device_type::execution_space execution_space
Kokkos execution space.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:43
double getCopyTime() const
Get the time spent deep copying local 3-array CRS out of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:310
int getNumApply() const
Get the number of times apply() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:286
void initialize()
Initialize the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:121
KokkosSparse::CrsMatrix< Scalar, LocalOrdinal, execution_space > KCrsMatrix
Kokkos CRS matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:53
bool isInitialized() const
Whether initialize() has been called since the last time the matrix's structure was changed.
Definition Ifpack2_Details_FastILU_Base_def.hpp:232
virtual std::string getName() const =0
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic")
Kokkos::View< LocalOrdinal *, execution_space >::host_mirror_type OrdinalArrayHost
Array of LocalOrdinal on host.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:57
virtual int getSweeps() const =0
Get the "sweeps" parameter.
void setParameters(const Teuchos::ParameterList &List)
Validate parameters, and set defaults when parameters are not provided.
Definition Ifpack2_Details_FastILU_Base_def.hpp:108
void apply(const TMultiVec &X, TMultiVec &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Apply the preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:58
double getApplyTime() const
Get the time spent in the last apply() call.
Definition Ifpack2_Details_FastILU_Base_def.hpp:304
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Get the domain map of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:45
virtual void checkLocalIC() const
Verify and print debug information about the underlying IC preconditioner.
Definition Ifpack2_Details_FastILU_Base_def.hpp:323
virtual void initLocalPrec()=0
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
int getNumInitialize() const
Get the number of times initialize() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:274
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TMultiVec
Tpetra multivector.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:51
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Get the range map of the matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:52
Teuchos::RCP< const TRowMatrix > getMatrix() const
Get the current matrix.
Definition Ifpack2_Details_FastILU_Base_def.hpp:268
Kokkos::View< LocalOrdinal *, execution_space > OrdinalArray
Array of LocalOrdinal on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:55
std::string description() const
Return a brief description of the preconditioner, in YAML format.
Definition Ifpack2_Details_FastILU_Base_def.hpp:329
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:59
virtual int getNTrisol() const =0
Get the "triangular solve iterations" parameter.
int getNumCompute() const
Get the number of times compute() was called.
Definition Ifpack2_Details_FastILU_Base_def.hpp:280
bool isComputed() const
Whether compute() has been called since the last time the matrix's values or structure were changed.
Definition Ifpack2_Details_FastILU_Base_def.hpp:261
virtual void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const =0
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only)
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40