13 #ifndef IFPACK2_CRSRILUK_DECL_HPP
14 #define IFPACK2_CRSRILUK_DECL_HPP
16 #include "KokkosSparse_spiluk.hpp"
20 #include "Tpetra_CrsMatrix_decl.hpp"
23 #include "Ifpack2_LocalSparseTriangularSolver_decl.hpp"
26 #include <type_traits>
212 template<
class MatrixType>
215 typename MatrixType::local_ordinal_type,
216 typename MatrixType::global_ordinal_type,
217 typename MatrixType::node_type>,
219 typename MatrixType::local_ordinal_type,
220 typename MatrixType::global_ordinal_type,
221 typename MatrixType::node_type> >
252 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
"Ifpack2::RILUK: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.");
263 template <
class NewMatrixType>
friend class RILUK;
265 typedef typename crs_matrix_type::global_inds_host_view_type global_inds_host_view_type;
266 typedef typename crs_matrix_type::local_inds_host_view_type local_inds_host_view_type;
267 typedef typename crs_matrix_type::values_host_view_type values_host_view_type;
270 typedef typename crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
271 typedef typename crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
272 typedef typename crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
279 typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
280 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
281 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
282 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
283 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
284 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
285 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
286 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
287 <
typename lno_row_view_t::const_value_type,
typename lno_nonzero_view_t::const_value_type,
typename scalar_nonzero_view_t::value_type,
288 HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;
339 return isInitialized_;
348 return numInitialize_;
361 return initializeTime_;
454 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
455 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
484 multiply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
485 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
509 true, std::logic_error,
"Ifpack2::RILUK::SetOverlapMode: "
510 "RILUK no longer implements overlap on its own. "
511 "Use RILUK with AdditiveSchwarz if you want overlap.");
516 return getL ().getGlobalNumEntries () +
getU ().getGlobalNumEntries ();
530 const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>&
548 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
552 void allocateSolvers ();
553 void allocate_L_and_U ();
557 void compute_serial();
558 void compute_kkspiluk();
560 #ifdef KOKKOS_ENABLE_CUDA
563 void compute_kkspiluk_stream();
564 #ifdef KOKKOS_ENABLE_CUDA
569 typedef Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> vec_type;
576 std::vector< Teuchos::RCP<iluk_graph_type> > Graph_v_;
583 std::vector<local_matrix_device_type> A_local_diagblks;
584 std::vector< lno_row_view_t > A_local_diagblks_rowmap_v_;
585 std::vector< lno_nonzero_view_t > A_local_diagblks_entries_v_;
586 std::vector< scalar_nonzero_view_t > A_local_diagblks_values_v_;
590 std::vector< Teuchos::RCP<crs_matrix_type> > L_v_;
595 std::vector< Teuchos::RCP<crs_matrix_type> > U_v_;
611 mutable int numApply_;
613 double initializeTime_;
615 mutable double applyTime_;
624 std::vector< Teuchos::RCP<kk_handle_type> > KernelHandle_v_;
625 bool isKokkosKernelsStream_;
627 std::vector<execution_space> exec_space_instances_;
628 bool hasStreamReordered_;
629 std::vector<typename lno_nonzero_view_t::non_const_type> perm_v_;
630 std::vector<typename lno_nonzero_view_t::non_const_type> reverse_perm_v_;
631 mutable std::unique_ptr<MV> Y_tmp_;
632 mutable std::unique_ptr<MV> reordered_x_;
633 mutable std::unique_ptr<MV> reordered_y_;
643 template<
class MatrixType,
class NodeType>
644 struct setLocalSolveParams{
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
Declaration and definition of IlukGraph.
magnitude_type getRelaxValue() const
Get RILU(k) relaxation parameter.
Definition: Ifpack2_RILUK_decl.hpp:495
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack2_RILUK_decl.hpp:498
Ifpack2::ScalingType enumerable type.
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type >, kk_handle_type > > getGraph() const
Return the Ifpack2::IlukGraph associated with this factored matrix.
Definition: Ifpack2_RILUK_decl.hpp:522
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_RILUK_def.hpp:241
node_type::execution_space execution_space
The Kokkos execution space of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:243
Ifpack2::Preconditioner< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::scalar_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::local_ordinal_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::global_ordinal_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::node_type >::magnitude_type Teuchos::ScalarTraits< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Preconditioner.hpp:78
static Teuchos::RCP< const row_matrix_type > makeLocalFilter(const Teuchos::RCP< const row_matrix_type > &A)
Return A, wrapped in a LocalFilter, if necessary.
Definition: Ifpack2_RILUK_def.hpp:448
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:231
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:589
const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > & getD() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:178
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:481
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_RILUK_decl.hpp:249
Teuchos::RCP< iluk_graph_type > Graph_
The ILU(k) graph.
Definition: Ifpack2_RILUK_decl.hpp:575
double getApplyTime() const
Total time in seconds taken by all successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:368
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_RILUK_def.hpp:125
std::string description() const
A one-line description of this object.
Definition: Ifpack2_RILUK_def.hpp:1434
void setParameters(const Teuchos::ParameterList ¶ms)
Definition: Ifpack2_RILUK_def.hpp:306
bool isKokkosKernelsSpiluk_
Optional KokkosKernels implementation.
Definition: Ifpack2_RILUK_decl.hpp:622
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:213
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_RILUK_def.hpp:205
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition: Ifpack2_RILUK_decl.hpp:507
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
Definition: Ifpack2_RILUK_def.hpp:441
bool isComputed() const
Whether compute() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:342
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_RILUK_decl.hpp:252
node_type::device_type device_type
The Kokkos device type of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:240
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:192
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:234
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_RILUK_def.hpp:222
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:228
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_RILUK_def.hpp:98
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition: Ifpack2_RILUK_decl.hpp:597
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition: Ifpack2_RILUK_decl.hpp:504
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:1099
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition: Ifpack2_RILUK_decl.hpp:501
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_RILUK_def.hpp:1183
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_RILUK_decl.hpp:237
crs_matrix_type::impl_scalar_type impl_scalar_type
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector)
Definition: Ifpack2_RILUK_decl.hpp:261
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:67
Declaration of interface for preconditioners that can change their matrix after construction.
bool isInitialized() const
Whether initialize() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:338
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition: Ifpack2_RILUK_decl.hpp:592
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:572
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition: Ifpack2_RILUK_def.hpp:434
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:364
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition: Ifpack2_RILUK_decl.hpp:600
int getNumCompute() const
Number of successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:351
int getNumApply() const
Number of successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:355
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack2_RILUK_decl.hpp:515
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object. ...
Definition: Ifpack2_RILUK_decl.hpp:360
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:161
Teuchos::RCP< const row_matrix_type > A_local_
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_t...
Definition: Ifpack2_RILUK_decl.hpp:580
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:594
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:347
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:225