46#ifndef IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
47#define IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
49#include <Tpetra_BlockCrsMatrix.hpp>
51#include <Ifpack2_RILUK.hpp>
127template<
class MatrixType>
129 typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_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;
153 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
172 template <
class NewMatrixType>
friend class RBILUK;
178 typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
179 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
180 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
181 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
182 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
183 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
184 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
185 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
186 <
typename lno_row_view_t::const_value_type,
typename lno_nonzero_view_t::const_value_type,
typename scalar_nonzero_view_t::value_type,
187 HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;
191 Teuchos::RCP<kk_handle_type> KernelHandle_;
200 RBILUK (
const Teuchos::RCP<const row_matrix_type>& A_in);
205 RBILUK (
const Teuchos::RCP<const block_crs_matrix_type>& A_in);
241 using RILUK<Tpetra::RowMatrix<
typename MatrixType::scalar_type,
242 typename MatrixType::local_ordinal_type,
243 typename MatrixType::global_ordinal_type,
269 setMatrix (
const Teuchos::RCP<const block_crs_matrix_type>& A);
312 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
313 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
314 Teuchos::ETransp mode = Teuchos::NO_TRANS,
315 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one (),
316 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero ())
const;
322 Teuchos::RCP<const block_crs_matrix_type>
getBlockMatrix ()
const;
325 const block_crs_matrix_type&
getLBlock ()
const;
328 const block_crs_matrix_type&
getDBlock ()
const;
331 const block_crs_matrix_type&
getUBlock ()
const;
334 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
335 typedef Teuchos::ScalarTraits<impl_scalar_type> STS;
336 typedef Teuchos::ScalarTraits<magnitude_type> STM;
337 typedef typename block_crs_matrix_type::little_block_type little_block_type;
338 typedef typename block_crs_matrix_type::little_block_host_type little_block_host_type;
339 typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
340 typedef typename block_crs_matrix_type::little_host_vec_type little_host_vec_type;
341 typedef typename block_crs_matrix_type::const_host_little_vec_type const_host_little_vec_type;
343 using local_inds_host_view_type =
typename block_crs_matrix_type::local_inds_host_view_type;
344 using values_host_view_type =
typename block_crs_matrix_type::values_host_view_type;
345 using local_inds_device_view_type =
typename block_crs_matrix_type::local_inds_device_view_type;
346 using values_device_view_type =
typename block_crs_matrix_type::values_device_view_type;
348 void allocate_L_and_U_blocks();
349 void initAllValues (
const block_crs_matrix_type& A);
352 Teuchos::RCP<const row_matrix_type> A_;
355 Teuchos::RCP<const block_crs_matrix_type> A_block_;
361 Teuchos::RCP<block_crs_matrix_type> L_block_;
363 Teuchos::RCP<block_crs_matrix_type> U_block_;
365 Teuchos::RCP<block_crs_matrix_type> D_block_;
368 Teuchos::RCP<block_crs_matrix_type> D_block_inverse_;
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:130
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:150
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:139
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
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:182
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:142
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:86
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:800
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:458
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
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:136
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:81
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:125
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:187
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:952
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:111
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
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:254
Ifpack2 features that are experimental. Use at your own risk.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:74