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  typedef typename Kokkos::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 
119 
121  typedef Tpetra::RowMatrix<scalar_type,
124  node_type>
126 
128  typedef Tpetra::CrsMatrix<scalar_type,
131  node_type>
133 
134  using crs_graph_type = Tpetra::CrsGraph<local_ordinal_type,
136  node_type>;
137 
138  typedef Tpetra::BlockCrsMatrix<scalar_type,
141  node_type>
142  block_crs_matrix_type;
143 
144  template <class NewMatrixType>
145  friend class RBILUK;
146 
147  typedef typename block_crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
148  typedef typename block_crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
149  typedef typename block_crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
150 
152 
154 
155  typedef typename block_crs_matrix_type::local_matrix_device_type local_matrix_device_type;
156  typedef typename block_crs_matrix_type::local_matrix_host_type local_matrix_host_type;
157  typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
158  typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
159  typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
160  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
161  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
162  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
163  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,
164  HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace>
165  kk_handle_type;
166  // typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
167  // <typename lno_row_view_t::non_const_value_type, typename lno_nonzero_view_t::non_const_value_type, typename scalar_nonzero_view_t::value_type,
168  // HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;//test
169  Teuchos::RCP<kk_handle_type> KernelHandle_;
170  Teuchos::RCP<kk_handle_type> L_Sptrsv_KernelHandle_;
171  Teuchos::RCP<kk_handle_type> U_Sptrsv_KernelHandle_;
172 
174 
176 
181 
186 
187  private:
190  RBILUK(const RBILUK<MatrixType>& src);
191 
192  public:
194  virtual ~RBILUK();
196 
198  void initialize();
199 
208  void compute();
209 
211 
212 
213  // Declare that we intend to overload RILUK::setMatrix, not hide it.
214  // This avoids build warnings that the method below "hides
215  // overloaded virtual function" (e.g., Clang 3.5).
216  //
217  // NOTE: If the base class of this class changes, e.g., if its
218  // template parameter changes, then be sure to change the code below
219  // to refer to the proper base class.
220  using RILUK<Tpetra::RowMatrix<typename MatrixType::scalar_type,
221  typename MatrixType::local_ordinal_type,
222  typename MatrixType::global_ordinal_type,
223  typename MatrixType::node_type> >::setMatrix;
224 
247  void
249 
251 
253 
255  std::string description() const;
256 
258 
260 
290  void
291  apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
292  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
297 
298  public:
301 
303  const block_crs_matrix_type& getLBlock() const;
304 
306  const block_crs_matrix_type& getDBlock() const;
307 
309  const block_crs_matrix_type& getUBlock() const;
310 
311  private:
312  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
315  typedef typename block_crs_matrix_type::little_block_type little_block_type;
316  typedef typename block_crs_matrix_type::little_block_host_type little_block_host_type;
317  typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
318  typedef typename block_crs_matrix_type::little_host_vec_type little_host_vec_type;
319  typedef typename block_crs_matrix_type::const_host_little_vec_type const_host_little_vec_type;
320 
321  using local_inds_host_view_type = typename block_crs_matrix_type::local_inds_host_view_type;
322  using values_host_view_type = typename block_crs_matrix_type::values_host_view_type;
323  using local_inds_device_view_type = typename block_crs_matrix_type::local_inds_device_view_type;
324  using values_device_view_type = typename block_crs_matrix_type::values_device_view_type;
325 
326  void allocate_L_and_U_blocks();
327  void initAllValues();
328 
330  local_ordinal_type blockSize_;
331 
338 
340  Teuchos::RCP<block_crs_matrix_type> D_block_inverse_;
341 
342  Kokkos::View<impl_scalar_type*, typename values_device_view_type::device_type> tmp_;
343 };
344 
345 } // namespace Experimental
346 
347 } // namespace Ifpack2
348 
349 #endif /* IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP */
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:257
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:91
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:140
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:913
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:118
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:115
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:132
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:107
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:127
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:508
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:94
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:114
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:125
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:1101
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:95