Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Experimental_RBILUK_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
12
13#ifndef IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
14#define IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
15
16#include <Tpetra_BlockCrsMatrix.hpp>
17
18#include <Ifpack2_RILUK.hpp>
19
20namespace Ifpack2 {
21
22namespace Experimental {
23
94template <class MatrixType>
95class RBILUK : virtual public Ifpack2::RILUK<Tpetra::RowMatrix<typename MatrixType::scalar_type,
96 typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> > {
97 public:
99
100
101 typedef typename MatrixType::scalar_type scalar_type;
102
103 // typedef typename MatrixType::impl_scalar_type impl_scalar_type;
104 typedef typename KokkosKernels::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
105
107 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
108 typedef typename MatrixType::local_ordinal_type LO;
109
111 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
112 typedef typename MatrixType::global_ordinal_type GO;
113
115 typedef typename MatrixType::node_type node_type;
116
117 typedef typename node_type::execution_space execution_space;
118
120 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
121
123 typedef Tpetra::RowMatrix<scalar_type,
126 node_type>
128
130 typedef Tpetra::CrsMatrix<scalar_type,
133 node_type>
135
136 using crs_graph_type = Tpetra::CrsGraph<local_ordinal_type,
138 node_type>;
139
140 typedef Tpetra::BlockCrsMatrix<scalar_type,
143 node_type>
144 block_crs_matrix_type;
145
146 template <class NewMatrixType>
147 friend class RBILUK;
148
149 typedef typename block_crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
150 typedef typename block_crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
151 typedef typename block_crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
152
154
156
157 typedef typename block_crs_matrix_type::local_matrix_device_type local_matrix_device_type;
158 typedef typename block_crs_matrix_type::local_matrix_host_type local_matrix_host_type;
159 typedef typename crs_matrix_type::local_matrix_device_type local_crs_matrix_device_type;
160 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
161 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
162 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
163 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
164 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
165 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
166 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,
167 HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace>
168 kk_handle_type;
169
171
173
174
177 RBILUK(const Teuchos::RCP<const row_matrix_type>& A_in);
178
182 RBILUK(const Teuchos::RCP<const block_crs_matrix_type>& A_in);
183
184 private:
187 RBILUK(const RBILUK<MatrixType>& src);
188
189 public:
191 virtual ~RBILUK();
193
195 void initialize();
196
205 void compute();
206
208
209
210 // Declare that we intend to overload RILUK::setMatrix, not hide it.
211 // This avoids build warnings that the method below "hides
212 // overloaded virtual function" (e.g., Clang 3.5).
213 //
214 // NOTE: If the base class of this class changes, e.g., if its
215 // template parameter changes, then be sure to change the code below
216 // to refer to the proper base class.
217 using RILUK<Tpetra::RowMatrix<typename MatrixType::scalar_type,
218 typename MatrixType::local_ordinal_type,
219 typename MatrixType::global_ordinal_type,
220 typename MatrixType::node_type> >::setMatrix;
221
244 void
245 setMatrix(const Teuchos::RCP<const block_crs_matrix_type>& A);
246
248
250
252 std::string description() const;
253
255
257
287 void
288 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
289 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
290 Teuchos::ETransp mode = Teuchos::NO_TRANS,
291 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
292 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
294
295 public:
297 Teuchos::RCP<const block_crs_matrix_type> getBlockMatrix() const;
298
300 const block_crs_matrix_type& getLBlock() const;
301
303 const block_crs_matrix_type& getDBlock() const;
304
306 const block_crs_matrix_type& getUBlock() const;
307
308 private:
309 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
310 typedef Teuchos::ScalarTraits<impl_scalar_type> STS;
311 typedef Teuchos::ScalarTraits<magnitude_type> STM;
312 typedef typename block_crs_matrix_type::little_block_type little_block_type;
313 typedef typename block_crs_matrix_type::little_block_host_type little_block_host_type;
314 typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
315 typedef typename block_crs_matrix_type::little_host_vec_type little_host_vec_type;
316 typedef typename block_crs_matrix_type::const_host_little_vec_type const_host_little_vec_type;
317
318 using local_inds_host_view_type = typename block_crs_matrix_type::local_inds_host_view_type;
319 using values_host_view_type = typename block_crs_matrix_type::values_host_view_type;
320 using local_inds_device_view_type = typename block_crs_matrix_type::local_inds_device_view_type;
321 using values_device_view_type = typename block_crs_matrix_type::values_device_view_type;
322 using view1d = Kokkos::View<impl_scalar_type*, typename values_device_view_type::device_type>;
323
324 void allocate_L_and_U_blocks();
325 void initAllValues();
326
327 template <typename X, typename Y>
328 void stream_apply_impl(const X& X_views, Y& Y_views, const bool use_temp_x, const bool use_temp_y, const std::vector<Teuchos::RCP<block_crs_matrix_type> >& LU_block_v, const std::vector<Teuchos::RCP<kk_handle_type> >& LU_sptrsv_handle_v, const LO numVecs) const;
329
330 Teuchos::RCP<kk_handle_type> KernelHandle_block_;
331 Teuchos::RCP<kk_handle_type> L_Sptrsv_KernelHandle_;
332 Teuchos::RCP<kk_handle_type> U_Sptrsv_KernelHandle_;
333
334 std::vector<Teuchos::RCP<kk_handle_type> > KernelHandle_block_v_;
335 std::vector<Teuchos::RCP<kk_handle_type> > L_Sptrsv_KernelHandle_v_;
336 std::vector<Teuchos::RCP<kk_handle_type> > U_Sptrsv_KernelHandle_v_;
337
339 local_ordinal_type blockSize_;
340
341 Teuchos::RCP<const block_crs_matrix_type> A_local_bcrs_;
342 Teuchos::RCP<crs_matrix_type> A_local_crs_nc_;
343
345 std::vector<local_matrix_device_type> A_block_local_diagblks_v_;
346 std::vector<lno_row_view_t> A_block_local_diagblks_rowmap_v_;
347 std::vector<lno_nonzero_view_t> A_block_local_diagblks_entries_v_;
348 std::vector<scalar_nonzero_view_t> A_block_local_diagblks_values_v_;
349
351 Teuchos::RCP<block_crs_matrix_type> L_block_;
352 std::vector<Teuchos::RCP<block_crs_matrix_type> > L_block_v_;
354 Teuchos::RCP<block_crs_matrix_type> U_block_;
355 std::vector<Teuchos::RCP<block_crs_matrix_type> > U_block_v_;
357 Teuchos::RCP<block_crs_matrix_type> D_block_;
358
360 Teuchos::RCP<block_crs_matrix_type> D_block_inverse_;
361
362 view1d tmp_;
363};
364
365} // namespace Experimental
366
367} // namespace Ifpack2
368
369#endif /* IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP */
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:96
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:115
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition Ifpack2_Experimental_RBILUK_def.hpp:208
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:120
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:107
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_Experimental_RBILUK_def.hpp:162
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_Experimental_RBILUK_def.hpp:1200
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_Experimental_RBILUK_def.hpp:701
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:127
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:101
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:111
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition Ifpack2_Experimental_RBILUK_def.hpp:159
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition Ifpack2_Experimental_RBILUK_def.hpp:195
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_Experimental_RBILUK_def.hpp:344
std::string description() const
A one-line description of this object.
Definition Ifpack2_Experimental_RBILUK_def.hpp:1383
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition Ifpack2_Experimental_RBILUK_def.hpp:182
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_Experimental_RBILUK_decl.hpp:134
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition Ifpack2_RILUK_decl.hpp:218
Ifpack2 features that are experimental. Use at your own risk.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40