Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_RILUK_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
12
13#ifndef IFPACK2_CRSRILUK_DECL_HPP
14#define IFPACK2_CRSRILUK_DECL_HPP
15
16#include "KokkosSparse_spiluk.hpp"
17
20#include "Tpetra_CrsMatrix_decl.hpp"
22#include "Ifpack2_IlukGraph.hpp"
23#include "Ifpack2_LocalSparseTriangularSolver_decl.hpp"
24
25#include <memory>
26#include <type_traits>
27
28namespace Teuchos {
29class ParameterList; // forward declaration
30}
31
32namespace Ifpack2 {
33
212template <class MatrixType>
213class RILUK : virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
214 typename MatrixType::local_ordinal_type,
215 typename MatrixType::global_ordinal_type,
216 typename MatrixType::node_type>,
217 virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
218 typename MatrixType::local_ordinal_type,
219 typename MatrixType::global_ordinal_type,
220 typename MatrixType::node_type> > {
221 public:
223 typedef typename MatrixType::scalar_type scalar_type;
224
226 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
227
229 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
230
232 typedef typename MatrixType::node_type node_type;
233
235 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
236
238 typedef typename node_type::device_type device_type;
239
241 typedef typename node_type::execution_space execution_space;
242
244 typedef Tpetra::RowMatrix<scalar_type,
247 node_type>
249
250 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.");
251
253 typedef Tpetra::CrsMatrix<scalar_type,
256 node_type>
258
260 typedef typename crs_matrix_type::impl_scalar_type impl_scalar_type;
261
262 template <class NewMatrixType>
263 friend class RILUK;
264
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;
268
269 typedef typename crs_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
270 typedef typename crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
271 typedef typename crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
272
274 typedef Tpetra::MultiVector<magnitude_type, local_ordinal_type,
277
279
281
282 typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
283 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
284 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
285 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
286 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
287 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
288 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
289 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,
290 HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace>
291 kk_handle_type;
293 typedef Kokkos::View<magnitude_type**, Kokkos::LayoutLeft, device_type> coors_view_t;
294 typedef Kokkos::View<local_ordinal_type*, Kokkos::LayoutLeft, device_type> perm_view_t;
295
299 RILUK(const Teuchos::RCP<const row_matrix_type>& A_in, const Teuchos::RCP<const coord_type>& A_in_coordinates = Teuchos::null);
300
308 RILUK(const Teuchos::RCP<const crs_matrix_type>& A_in, const Teuchos::RCP<const coord_type>& A_in_coordinates = Teuchos::null);
309
310 private:
313 RILUK(const RILUK<MatrixType>& src);
314
315 public:
317 virtual ~RILUK();
318
327 void setParameters(const Teuchos::ParameterList& params);
328
330 void initialize();
331
340 void compute();
341
343 bool isInitialized() const {
344 return isInitialized_;
345 }
347 bool isComputed() const {
348 return isComputed_;
349 }
350
352 int getNumInitialize() const {
353 return numInitialize_;
354 }
356 int getNumCompute() const {
357 return numCompute_;
358 }
360 int getNumApply() const {
361 return numApply_;
362 }
363
365 double getInitializeTime() const {
366 return initializeTime_;
367 }
369 double getComputeTime() const {
370 return computeTime_;
371 }
373 double getApplyTime() const {
374 return applyTime_;
375 }
376
378 size_t getNodeSmootherComplexity() const;
379
381
382
405 virtual void
406 setMatrix(const Teuchos::RCP<const row_matrix_type>& A);
407
411 void
412 setCoord(const Teuchos::RCP<const coord_type>& A_coordinates);
413
415
417
419 std::string description() const;
420
422
424
426 Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
427 getDomainMap() const;
428
430 Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
431 getRangeMap() const;
432
462 void
463 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
464 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
465 Teuchos::ETransp mode = Teuchos::NO_TRANS,
466 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
467 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
469
470 private:
492 void
493 multiply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
494 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
495 const Teuchos::ETransp mode = Teuchos::NO_TRANS) const;
496
497 public:
499 Teuchos::RCP<const row_matrix_type> getMatrix() const;
500
502 Teuchos::RCP<const coord_type> getCoord() const;
503
504 // Attribute access functions
505
507 magnitude_type getRelaxValue() const { return RelaxValue_; }
508
510 magnitude_type getAbsoluteThreshold() const { return Athresh_; }
511
513 magnitude_type getRelativeThreshold() const { return Rthresh_; }
514
516 int getLevelOfFill() const { return LevelOfFill_; }
517
519 Tpetra::CombineMode getOverlapMode() {
520 TEUCHOS_TEST_FOR_EXCEPTION(
521 true, std::logic_error,
522 "Ifpack2::RILUK::SetOverlapMode: "
523 "RILUK no longer implements overlap on its own. "
524 "Use RILUK with AdditiveSchwarz if you want overlap.");
525 }
526
528 Tpetra::global_size_t getGlobalNumEntries() const {
529 return getL().getGlobalNumEntries() + getU().getGlobalNumEntries();
530 }
531
533 Teuchos::RCP<Ifpack2::IlukGraph<Tpetra::CrsGraph<local_ordinal_type,
535 node_type>,
536 kk_handle_type> >
537 getGraph() const {
538 return Graph_;
539 }
540
542 const crs_matrix_type& getL() const;
543
545 const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>&
546 getD() const;
547
549 const crs_matrix_type& getU() const;
550
552 Teuchos::RCP<const crs_matrix_type> getCrsMatrix() const;
553
559 static Teuchos::RCP<const row_matrix_type>
560 makeLocalFilter(const Teuchos::RCP<const row_matrix_type>& A);
561
562 private:
563 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
564 typedef Teuchos::ScalarTraits<scalar_type> STS;
565 typedef Teuchos::ScalarTraits<magnitude_type> STM;
566
567 void allocateSolvers();
568 void allocate_L_and_U();
569 static void checkOrderingConsistency(const row_matrix_type& A);
570 void initAllValues(const row_matrix_type& A);
571
572 void compute_serial();
573 void compute_kkspiluk();
574// Workaround Cuda limitation of KOKKOS_LAMBDA in private/protected member functions
575#ifdef KOKKOS_ENABLE_CUDA
576 public:
577#endif
578 void compute_kkspiluk_stream();
579#ifdef KOKKOS_ENABLE_CUDA
580 private:
581#endif
582
583 protected:
584 typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vec_type;
585
587 Teuchos::RCP<const row_matrix_type> A_;
588
591 Teuchos::RCP<const coord_type> A_coordinates_;
592
594 Teuchos::RCP<iluk_graph_type> Graph_;
595 std::vector<Teuchos::RCP<iluk_graph_type> > Graph_v_;
599 Teuchos::RCP<const row_matrix_type> A_local_;
600 Teuchos::RCP<const crs_matrix_type> A_local_crs_;
601 Teuchos::RCP<crs_matrix_type> A_local_crs_nc_;
602 std::vector<local_matrix_device_type> A_local_diagblks_v_;
603 std::vector<lno_row_view_t> A_local_diagblks_rowmap_v_;
604 std::vector<lno_nonzero_view_t> A_local_diagblks_entries_v_;
605 std::vector<scalar_nonzero_view_t> A_local_diagblks_values_v_;
606
608 Teuchos::RCP<crs_matrix_type> L_;
609 std::vector<Teuchos::RCP<crs_matrix_type> > L_v_;
611 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> > L_solver_;
613 Teuchos::RCP<crs_matrix_type> U_;
614 std::vector<Teuchos::RCP<crs_matrix_type> > U_v_;
616 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> > U_solver_;
617
619 Teuchos::RCP<vec_type> D_;
620
621 int LevelOfFill_;
622 double Overalloc_;
623
624 bool isAllocated_;
625 bool isInitialized_;
626 bool isComputed_;
627
628 int numInitialize_;
629 int numCompute_;
630 mutable int numApply_;
631
632 double initializeTime_;
633 double computeTime_;
634 mutable double applyTime_;
635
636 magnitude_type RelaxValue_;
637 magnitude_type Athresh_;
638 magnitude_type Rthresh_;
639
642 Teuchos::RCP<kk_handle_type> KernelHandle_;
643 std::vector<Teuchos::RCP<kk_handle_type> > KernelHandle_v_;
644 bool isKokkosKernelsStream_;
645 int num_streams_;
646 std::vector<execution_space> exec_space_instances_;
647 bool hasStreamReordered_;
648 bool hasStreamsWithRCB_;
649 std::vector<typename lno_nonzero_view_t::non_const_type> perm_v_;
650 std::vector<typename lno_nonzero_view_t::non_const_type> reverse_perm_v_;
651 mutable std::unique_ptr<MV> Y_tmp_;
652 mutable std::unique_ptr<MV> reordered_x_;
653 mutable std::unique_ptr<MV> reordered_y_;
654 perm_view_t perm_rcb_;
655 coors_view_t coors_rcb_;
656};
657
658// NOTE (mfh 11 Feb 2015) This used to exist in order to deal with
659// different behavior of Tpetra::Crs{Graph,Matrix} for
660// KokkosClassic::ThrustGPUNode. In particular, fillComplete on a
661// CrsMatrix used to make the graph go away by default, so we had to
662// pass in a parameter to keep a host copy of the graph. With the new
663// (Kokkos refactor) version of Tpetra, this problem has gone away.
664namespace detail {
665template <class MatrixType, class NodeType>
666struct setLocalSolveParams {
667 static Teuchos::RCP<Teuchos::ParameterList>
668 setParams(const Teuchos::RCP<Teuchos::ParameterList>& param) {
669 return param;
670 }
671};
672} // namespace detail
673
674} // namespace Ifpack2
675
676#endif /* IFPACK2_CRSRILUK_DECL_HPP */
Declaration of interface for preconditioners that can change their matrix after construction.
Declaration and definition of IlukGraph.
Ifpack2::ScalingType enumerable type.
Mix-in interface for preconditioners that can change their matrix after construction.
Definition Ifpack2_Details_CanChangeMatrix.hpp:60
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition Ifpack2_IlukGraph.hpp:67
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition Ifpack2_RILUK_decl.hpp:220
crs_matrix_type::impl_scalar_type impl_scalar_type
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector)
Definition Ifpack2_RILUK_decl.hpp:260
Tpetra::MultiVector< magnitude_type, local_ordinal_type, global_ordinal_type, node_type > coord_type
Tpetra::MultiVector specialization used for containing coordinates.
Definition Ifpack2_RILUK_decl.hpp:276
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:467
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:257
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition Ifpack2_RILUK_def.hpp:98
double getApplyTime() const
Total time in seconds taken by all successful apply() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:373
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:1171
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_RILUK_decl.hpp:235
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:369
int getNumApply() const
Number of successful apply() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:360
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:239
magnitude_type getRelaxValue() const
Get RILU(k) relaxation parameter.
Definition Ifpack2_RILUK_decl.hpp:507
int getNumCompute() const
Number of successful compute() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:356
Teuchos::RCP< const coord_type > A_coordinates_
Definition Ifpack2_RILUK_decl.hpp:591
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition Ifpack2_RILUK_decl.hpp:587
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition Ifpack2_RILUK_decl.hpp:519
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition Ifpack2_RILUK_decl.hpp:613
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:190
std::string description() const
A one-line description of this object.
Definition Ifpack2_RILUK_def.hpp:1438
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition Ifpack2_RILUK_decl.hpp:528
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_RILUK_def.hpp:121
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:223
bool isComputed() const
Whether compute() has been called on this object.
Definition Ifpack2_RILUK_decl.hpp:347
Teuchos::RCP< const coord_type > getCoord() const
Get the coordinates associated with the input matrix's rows.
Definition Ifpack2_RILUK_def.hpp:431
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:177
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:232
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:599
void setCoord(const Teuchos::RCP< const coord_type > &A_coordinates)
Set the matrix rows' coordinates.
Definition Ifpack2_RILUK_def.hpp:153
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition Ifpack2_RILUK_def.hpp:419
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition Ifpack2_RILUK_decl.hpp:510
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition Ifpack2_RILUK_decl.hpp:611
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition Ifpack2_RILUK_decl.hpp:516
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:425
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:219
void setParameters(const Teuchos::ParameterList &params)
Definition Ifpack2_RILUK_def.hpp:298
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition Ifpack2_RILUK_decl.hpp:513
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:352
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:226
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition Ifpack2_RILUK_decl.hpp:608
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:229
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:161
node_type::device_type device_type
The Kokkos device type of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:238
node_type::execution_space execution_space
The Kokkos execution space of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:241
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition Ifpack2_RILUK_decl.hpp:616
Teuchos::RCP< iluk_graph_type > Graph_
The ILU(k) graph.
Definition Ifpack2_RILUK_decl.hpp:594
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_RILUK_def.hpp:202
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:1084
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:437
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition Ifpack2_RILUK_decl.hpp:619
bool isKokkosKernelsSpiluk_
Optional KokkosKernels implementation.
Definition Ifpack2_RILUK_decl.hpp:641
bool isInitialized() const
Whether initialize() has been called on this object.
Definition Ifpack2_RILUK_decl.hpp:343
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:365
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:248
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:537
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40