Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 
20 namespace Ifpack2 {
21 
22 namespace Experimental {
23 
94 template <class MatrixType>
95 class 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  typedef typename MatrixType::scalar_type scalar_type;
102 
103  // typedef typename MatrixType::impl_scalar_type impl_scalar_type;
104 #if KOKKOS_VERSION >= 40799
105  typedef typename KokkosKernels::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
106 #else
107  typedef typename Kokkos::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
108 #endif
109 
111  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
112  typedef typename MatrixType::local_ordinal_type LO;
113 
115  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
116  typedef typename MatrixType::global_ordinal_type GO;
117 
119  typedef typename MatrixType::node_type node_type;
120 
121  typedef typename node_type::execution_space execution_space;
122 
125 
127  typedef Tpetra::RowMatrix<scalar_type,
130  node_type>
132 
134  typedef Tpetra::CrsMatrix<scalar_type,
137  node_type>
139 
140  using crs_graph_type = Tpetra::CrsGraph<local_ordinal_type,
142  node_type>;
143 
144  typedef Tpetra::BlockCrsMatrix<scalar_type,
147  node_type>
148  block_crs_matrix_type;
149 
150  template <class NewMatrixType>
151  friend class RBILUK;
152 
153  typedef typename block_crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
154  typedef typename block_crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
155  typedef typename block_crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
156 
158 
160 
161  typedef typename block_crs_matrix_type::local_matrix_device_type local_matrix_device_type;
162  typedef typename block_crs_matrix_type::local_matrix_host_type local_matrix_host_type;
163  typedef typename crs_matrix_type::local_matrix_device_type local_crs_matrix_device_type;
164  typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
165  typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
166  typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
167  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
168  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
169  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
170  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,
171  HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace>
172  kk_handle_type;
173 
175 
177 
182 
187 
188  private:
191  RBILUK(const RBILUK<MatrixType>& src);
192 
193  public:
195  virtual ~RBILUK();
197 
199  void initialize();
200 
209  void compute();
210 
212 
213 
214  // Declare that we intend to overload RILUK::setMatrix, not hide it.
215  // This avoids build warnings that the method below "hides
216  // overloaded virtual function" (e.g., Clang 3.5).
217  //
218  // NOTE: If the base class of this class changes, e.g., if its
219  // template parameter changes, then be sure to change the code below
220  // to refer to the proper base class.
221  using RILUK<Tpetra::RowMatrix<typename MatrixType::scalar_type,
222  typename MatrixType::local_ordinal_type,
223  typename MatrixType::global_ordinal_type,
224  typename MatrixType::node_type> >::setMatrix;
225 
248  void
250 
252 
254 
256  std::string description() const;
257 
259 
261 
291  void
292  apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
293  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
298 
299  public:
302 
304  const block_crs_matrix_type& getLBlock() const;
305 
307  const block_crs_matrix_type& getDBlock() const;
308 
310  const block_crs_matrix_type& getUBlock() const;
311 
312  private:
313  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
316  typedef typename block_crs_matrix_type::little_block_type little_block_type;
317  typedef typename block_crs_matrix_type::little_block_host_type little_block_host_type;
318  typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
319  typedef typename block_crs_matrix_type::little_host_vec_type little_host_vec_type;
320  typedef typename block_crs_matrix_type::const_host_little_vec_type const_host_little_vec_type;
321 
322  using local_inds_host_view_type = typename block_crs_matrix_type::local_inds_host_view_type;
323  using values_host_view_type = typename block_crs_matrix_type::values_host_view_type;
324  using local_inds_device_view_type = typename block_crs_matrix_type::local_inds_device_view_type;
325  using values_device_view_type = typename block_crs_matrix_type::values_device_view_type;
326  using view1d = Kokkos::View<impl_scalar_type*, typename values_device_view_type::device_type>;
327 
328  void allocate_L_and_U_blocks();
329  void initAllValues();
330 
331  template <typename X, typename Y>
332  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;
333 
334  Teuchos::RCP<kk_handle_type> KernelHandle_block_;
335  Teuchos::RCP<kk_handle_type> L_Sptrsv_KernelHandle_;
336  Teuchos::RCP<kk_handle_type> U_Sptrsv_KernelHandle_;
337 
338  std::vector<Teuchos::RCP<kk_handle_type> > KernelHandle_block_v_;
339  std::vector<Teuchos::RCP<kk_handle_type> > L_Sptrsv_KernelHandle_v_;
340  std::vector<Teuchos::RCP<kk_handle_type> > U_Sptrsv_KernelHandle_v_;
341 
343  local_ordinal_type blockSize_;
344 
346  Teuchos::RCP<crs_matrix_type> A_local_crs_nc_;
347 
349  std::vector<local_matrix_device_type> A_block_local_diagblks_v_;
350  std::vector<lno_row_view_t> A_block_local_diagblks_rowmap_v_;
351  std::vector<lno_nonzero_view_t> A_block_local_diagblks_entries_v_;
352  std::vector<scalar_nonzero_view_t> A_block_local_diagblks_values_v_;
353 
356  std::vector<Teuchos::RCP<block_crs_matrix_type> > L_block_v_;
359  std::vector<Teuchos::RCP<block_crs_matrix_type> > U_block_v_;
362 
364  Teuchos::RCP<block_crs_matrix_type> D_block_inverse_;
365 
366  view1d tmp_;
367 };
368 
369 } // namespace Experimental
370 
371 } // namespace Ifpack2
372 
373 #endif /* IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP */
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:344
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:115
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_Experimental_RBILUK_def.hpp:159
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:208
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:101
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
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:124
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:213
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:119
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:138
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:111
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:195
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:701
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:162
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:182
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:131
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:1383
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:95