Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_MDF_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
12
13#ifndef IFPACK2_MDF_DECL_HPP
14#define IFPACK2_MDF_DECL_HPP
15
18#include "Tpetra_CrsMatrix_decl.hpp"
20#include "Ifpack2_IlukGraph.hpp"
21#include "Ifpack2_LocalSparseTriangularSolver_decl.hpp"
22#include "KokkosSparse_mdf.hpp"
23
24#include <type_traits>
25
26namespace Teuchos {
27class ParameterList; // forward declaration
28}
29namespace Ifpack2 {
30
49template <class MatrixType>
50class MDF : virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
51 typename MatrixType::local_ordinal_type,
52 typename MatrixType::global_ordinal_type,
53 typename MatrixType::node_type>,
54 virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
55 typename MatrixType::local_ordinal_type,
56 typename MatrixType::global_ordinal_type,
57 typename MatrixType::node_type> > {
58 public:
60 typedef typename MatrixType::scalar_type scalar_type;
61
63 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
64
66 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
67
69 typedef typename MatrixType::node_type node_type;
70
72 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
73
75 typedef typename node_type::device_type device_type;
76
78 typedef typename node_type::execution_space execution_space;
79
81 typedef Tpetra::RowMatrix<scalar_type,
86
87 static_assert(std::is_same<MatrixType, row_matrix_type>::value, "Ifpack2::MDF: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.");
88
90 typedef Tpetra::CrsMatrix<scalar_type,
95
97 typedef typename crs_matrix_type::impl_scalar_type impl_scalar_type;
98
99 template <class NewMatrixType>
100 friend class MDF;
101
102 typedef typename crs_matrix_type::global_inds_host_view_type global_inds_host_view_type;
103 typedef typename crs_matrix_type::local_inds_host_view_type local_inds_host_view_type;
104 typedef typename crs_matrix_type::values_host_view_type values_host_view_type;
105
106 typedef typename crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
107 typedef typename crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
108 typedef typename crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
109
111
113
114 typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
115 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
116 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
117 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
118 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
119 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
120 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
121 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle<typename lno_row_view_t::const_value_type, typename lno_nonzero_view_t::const_value_type, typename scalar_nonzero_view_t::value_type,
122 HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace>
123 kk_handle_type;
124
128 MDF(const Teuchos::RCP<const row_matrix_type>& A_in);
129
137 MDF(const Teuchos::RCP<const crs_matrix_type>& A_in);
138
139 private:
142 MDF(const MDF<MatrixType>& src);
143
144 public:
146 virtual ~MDF() = default;
147
154 void setParameters(const Teuchos::ParameterList& params);
155
157 void initialize();
158
167 void compute();
168
170 bool isInitialized() const {
171 return isInitialized_;
172 }
174 bool isComputed() const {
175 return isComputed_;
176 }
177
179 int getNumInitialize() const {
180 return numInitialize_;
181 }
183 int getNumCompute() const {
184 return numCompute_;
185 }
187 int getNumApply() const {
188 return numApply_;
189 }
190
192 double getInitializeTime() const {
193 return initializeTime_;
194 }
196 double getComputeTime() const {
197 return computeTime_;
198 }
200 double getApplyTime() const {
201 return applyTime_;
202 }
203
205 size_t getNodeSmootherComplexity() const;
206
208
209
232 virtual void
233 setMatrix(const Teuchos::RCP<const row_matrix_type>& A);
234
236
238
240 std::string description() const;
241
243
245
247 Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
248 getDomainMap() const;
249
251 Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
252 getRangeMap() const;
253
283 void
284 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
285 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
286 Teuchos::ETransp mode = Teuchos::NO_TRANS,
287 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
288 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
290
291 private:
292 // Split off to a different impl call so that nested apply calls don't mess up apply counts/timers
293 void apply_impl(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
294 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
295 Teuchos::ETransp mode = Teuchos::NO_TRANS,
296 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
297 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
298
320 void
321 multiply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
322 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
323 const Teuchos::ETransp mode = Teuchos::NO_TRANS) const;
324
325 public:
326 using MDF_handle_device_type = KokkosSparse::Experimental::MDF_handle<local_matrix_device_type>;
327 using permutations_type = Teuchos::ArrayRCP<local_ordinal_type>;
328
330 Teuchos::RCP<const row_matrix_type> getMatrix() const;
331
333 int getLevelOfFill() const { return LevelOfFill_; }
334
336 Tpetra::CombineMode getOverlapMode() {
337 TEUCHOS_TEST_FOR_EXCEPTION(
338 true, std::logic_error,
339 "Ifpack2::MDF::SetOverlapMode: "
340 "MDF no longer implements overlap on its own. "
341 "Use MDF with AdditiveSchwarz if you want overlap.");
342 }
343
345 Tpetra::global_size_t getGlobalNumEntries() const {
346 return getL().getGlobalNumEntries() + getU().getGlobalNumEntries();
347 }
348
350 const crs_matrix_type& getL() const;
351
353 const crs_matrix_type& getU() const;
354
356 permutations_type& getPermutations() const;
357
359 permutations_type& getReversePermutations() const;
360
362 Teuchos::RCP<const crs_matrix_type> getCrsMatrix() const;
363
364 private:
365 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
366 typedef Teuchos::ScalarTraits<scalar_type> STS;
367 typedef Teuchos::ScalarTraits<magnitude_type> STM;
368
369 void allocateSolvers();
370 void allocatePermutations(bool force = false);
371 static void checkOrderingConsistency(const row_matrix_type& A);
372 // void initAllValues (const row_matrix_type& A);
373
379 static Teuchos::RCP<const row_matrix_type>
380 makeLocalFilter(const Teuchos::RCP<const row_matrix_type>& A);
381
382 protected:
383 typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vec_type;
384
386 Teuchos::RCP<const row_matrix_type> A_;
387
388 // The MDF handle
389 Teuchos::RCP<MDF_handle_device_type> MDF_handle_;
390
394 Teuchos::RCP<const row_matrix_type> A_local_;
395 lno_row_view_t A_local_rowmap_;
396 lno_nonzero_view_t A_local_entries_;
397 scalar_nonzero_view_t A_local_values_;
398
400 Teuchos::RCP<crs_matrix_type> L_;
402 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> > L_solver_;
404 Teuchos::RCP<crs_matrix_type> U_;
406 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> > U_solver_;
407
409 permutations_type permutations_;
410
412 permutations_type reversePermutations_;
413
414 int Verbosity_;
415
416 int LevelOfFill_;
417 double Overalloc_;
418
419 bool isAllocated_;
420 bool isInitialized_;
421 bool isComputed_;
422
423 int numInitialize_;
424 int numCompute_;
425 mutable int numApply_;
426
427 double initializeTime_;
428 double computeTime_;
429 mutable double applyTime_;
430};
431
432} // namespace Ifpack2
433
434#endif /* IFPACK2_MDF_DECL_HPP */
Declaration of interface for preconditioners that can change their matrix after construction.
Declaration and definition of IlukGraph.
Ifpack2::ScalingType enumerable type.
Mix-in interface for preconditioners that can change their matrix after construction.
Definition Ifpack2_Details_CanChangeMatrix.hpp:60
MDF (incomplete LU factorization with minimum discarded fill reordering) of a Tpetra sparse matrix.
Definition Ifpack2_MDF_decl.hpp:57
void setParameters(const Teuchos::ParameterList &params)
Definition Ifpack2_MDF_def.hpp:318
double getApplyTime() const
Total time in seconds taken by all successful apply() calls for this object.
Definition Ifpack2_MDF_decl.hpp:200
const crs_matrix_type & getL() const
Return the L factor of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:213
node_type::execution_space execution_space
The Kokkos execution space of the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:78
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the (inverse of the) incomplete factorization to X, resulting in Y.
Definition Ifpack2_MDF_def.hpp:626
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_MDF_decl.hpp:72
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition Ifpack2_MDF_decl.hpp:333
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition Ifpack2_MDF_decl.hpp:402
bool isComputed() const
Whether compute() has been called on this object.
Definition Ifpack2_MDF_decl.hpp:174
permutations_type & getReversePermutations() const
Return the reverse permutations of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:238
permutations_type reversePermutations_
The reverse permuations from MDF factorization.
Definition Ifpack2_MDF_decl.hpp:412
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition Ifpack2_MDF_decl.hpp:336
const crs_matrix_type & getU() const
Return the U factor of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:251
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_MDF_def.hpp:504
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition Ifpack2_MDF_decl.hpp:404
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:63
bool isInitialized() const
Whether initialize() has been called on this object.
Definition Ifpack2_MDF_decl.hpp:170
Teuchos::RCP< const row_matrix_type > A_local_
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_t...
Definition Ifpack2_MDF_decl.hpp:394
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_MDF_def.hpp:410
int getNumCompute() const
Number of successful compute() calls for this object.
Definition Ifpack2_MDF_decl.hpp:183
std::string description() const
A one-line description of this object.
Definition Ifpack2_MDF_def.hpp:696
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition Ifpack2_MDF_decl.hpp:196
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition Ifpack2_MDF_def.hpp:300
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition Ifpack2_MDF_decl.hpp:345
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U.
Definition Ifpack2_MDF_decl.hpp:94
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition Ifpack2_MDF_def.hpp:280
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition Ifpack2_MDF_decl.hpp:85
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
Definition Ifpack2_MDF_def.hpp:374
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:69
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:66
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object.
Definition Ifpack2_MDF_decl.hpp:192
permutations_type permutations_
The computed permuations from MDF factorization.
Definition Ifpack2_MDF_decl.hpp:409
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition Ifpack2_MDF_decl.hpp:179
crs_matrix_type::impl_scalar_type impl_scalar_type
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector)
Definition Ifpack2_MDF_decl.hpp:97
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition Ifpack2_MDF_decl.hpp:406
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:60
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition Ifpack2_MDF_def.hpp:368
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_MDF_def.hpp:179
node_type::device_type device_type
The Kokkos device type of the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:75
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_MDF_def.hpp:263
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition Ifpack2_MDF_decl.hpp:400
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition Ifpack2_MDF_decl.hpp:386
permutations_type & getPermutations() const
Return the permutations of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:226
virtual ~MDF()=default
Destructor (declared virtual for memory safety).
int getNumApply() const
Number of successful apply() calls for this object.
Definition Ifpack2_MDF_decl.hpp:187
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40