46 #ifndef IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
47 #define IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
49 #include <Tpetra_BlockCrsMatrix.hpp>
51 #include <Ifpack2_RILUK.hpp>
55 namespace Experimental {
127 template<
class MatrixType>
129 typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> >
135 typedef typename MatrixType::scalar_type
scalar_type;
139 typedef typename Kokkos::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
143 typedef typename MatrixType::local_ordinal_type LO;
147 typedef typename MatrixType::global_ordinal_type GO;
176 template <
class NewMatrixType>
friend class RBILUK;
178 typedef typename block_crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
179 typedef typename block_crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
180 typedef typename block_crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
186 typedef typename block_crs_matrix_type::local_matrix_device_type local_matrix_device_type;
187 typedef typename block_crs_matrix_type::local_matrix_host_type local_matrix_host_type;
188 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
189 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
190 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
191 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
192 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
193 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
194 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
195 <
typename lno_row_view_t::const_value_type,
typename lno_nonzero_view_t::const_value_type,
typename scalar_nonzero_view_t::value_type,
196 HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;
250 using RILUK<Tpetra::RowMatrix<
typename MatrixType::scalar_type,
251 typename MatrixType::local_ordinal_type,
252 typename MatrixType::global_ordinal_type,
321 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
322 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
334 const block_crs_matrix_type&
getLBlock ()
const;
337 const block_crs_matrix_type&
getDBlock ()
const;
340 const block_crs_matrix_type&
getUBlock ()
const;
343 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
346 typedef typename block_crs_matrix_type::little_block_type little_block_type;
347 typedef typename block_crs_matrix_type::little_block_host_type little_block_host_type;
348 typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
349 typedef typename block_crs_matrix_type::little_host_vec_type little_host_vec_type;
350 typedef typename block_crs_matrix_type::const_host_little_vec_type const_host_little_vec_type;
352 using local_inds_host_view_type =
typename block_crs_matrix_type::local_inds_host_view_type;
353 using values_host_view_type =
typename block_crs_matrix_type::values_host_view_type;
354 using local_inds_device_view_type =
typename block_crs_matrix_type::local_inds_device_view_type;
355 using values_device_view_type =
typename block_crs_matrix_type::values_device_view_type;
357 void allocate_L_and_U_blocks();
358 void initAllValues ();
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:347
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:146
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_Experimental_RBILUK_def.hpp:136
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:193
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:136
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:1010
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:153
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:245
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:150
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:165
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:142
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:179
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:595
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:141
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:165
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:159
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:1205
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:128