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> >
102 typedef typename MatrixType::scalar_type
scalar_type;
106 typedef typename Kokkos::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
110 typedef typename MatrixType::local_ordinal_type LO;
114 typedef typename MatrixType::global_ordinal_type GO;
143 template <
class NewMatrixType>
friend class RBILUK;
145 typedef typename block_crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
146 typedef typename block_crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
147 typedef typename block_crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
153 typedef typename block_crs_matrix_type::local_matrix_device_type local_matrix_device_type;
154 typedef typename block_crs_matrix_type::local_matrix_host_type local_matrix_host_type;
155 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
156 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
157 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
158 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
159 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
160 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
161 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
162 <
typename lno_row_view_t::const_value_type,
typename lno_nonzero_view_t::const_value_type,
typename scalar_nonzero_view_t::value_type,
163 HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;
219 using RILUK<Tpetra::RowMatrix<
typename MatrixType::scalar_type,
220 typename MatrixType::local_ordinal_type,
221 typename MatrixType::global_ordinal_type,
290 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
291 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
303 const block_crs_matrix_type&
getLBlock ()
const;
306 const block_crs_matrix_type&
getDBlock ()
const;
309 const block_crs_matrix_type&
getUBlock ()
const;
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;
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;
326 void allocate_L_and_U_blocks();
327 void initAllValues ();
342 Kokkos::View<impl_scalar_type*, typename values_device_view_type::device_type> tmp_;
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:286
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:113
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_Experimental_RBILUK_def.hpp:103
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:160
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:103
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:962
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:117
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:109
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:146
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:543
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:108
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:132
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:126
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:1159
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:95