13 #ifndef IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
14 #define IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
16 #include <Tpetra_BlockCrsMatrix.hpp>
18 #include <Ifpack2_RILUK.hpp>
22 namespace Experimental {
94 template <
class MatrixType>
96 typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> > {
100 typedef typename MatrixType::scalar_type
scalar_type;
104 typedef typename Kokkos::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
108 typedef typename MatrixType::local_ordinal_type LO;
112 typedef typename MatrixType::global_ordinal_type GO;
117 typedef typename node_type::execution_space execution_space;
144 block_crs_matrix_type;
146 template <
class NewMatrixType>
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;
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>
217 using RILUK<Tpetra::RowMatrix<
typename MatrixType::scalar_type,
218 typename MatrixType::local_ordinal_type,
219 typename MatrixType::global_ordinal_type,
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,
300 const block_crs_matrix_type&
getLBlock()
const;
303 const block_crs_matrix_type&
getDBlock()
const;
306 const block_crs_matrix_type&
getUBlock()
const;
309 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
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;
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>;
324 void allocate_L_and_U_blocks();
325 void initAllValues();
327 template <
typename X,
typename Y>
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_;
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_;
352 std::vector<Teuchos::RCP<block_crs_matrix_type> > L_block_v_;
355 std::vector<Teuchos::RCP<block_crs_matrix_type> > U_block_v_;
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:111
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:120
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:134
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: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:127
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