Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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"
21 #include "Ifpack2_ScalingType.hpp"
22 #include "Ifpack2_IlukGraph.hpp"
23 #include "Ifpack2_LocalSparseTriangularSolver_decl.hpp"
24 
25 #include <memory>
26 #include <type_traits>
27 
28 namespace Teuchos {
29 class ParameterList; // forward declaration
30 }
31 
32 namespace Ifpack2 {
33 
212 template <class MatrixType>
213 class 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 
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 
276 
277  typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
278  typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
279  typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
280  typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
281  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
282  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
283  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
284  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,
285  HandleExecSpace, TemporaryMemorySpace, PersistentMemorySpace>
286  kk_handle_type;
288 
293 
302 
303  private:
306  RILUK(const RILUK<MatrixType>& src);
307 
308  public:
310  virtual ~RILUK();
311 
320  void setParameters(const Teuchos::ParameterList& params);
321 
323  void initialize();
324 
333  void compute();
334 
336  bool isInitialized() const {
337  return isInitialized_;
338  }
340  bool isComputed() const {
341  return isComputed_;
342  }
343 
345  int getNumInitialize() const {
346  return numInitialize_;
347  }
349  int getNumCompute() const {
350  return numCompute_;
351  }
353  int getNumApply() const {
354  return numApply_;
355  }
356 
358  double getInitializeTime() const {
359  return initializeTime_;
360  }
362  double getComputeTime() const {
363  return computeTime_;
364  }
366  double getApplyTime() const {
367  return applyTime_;
368  }
369 
371  size_t getNodeSmootherComplexity() const;
372 
374 
375 
398  virtual void
400 
402 
404 
406  std::string description() const;
407 
409 
411 
414  getDomainMap() const;
415 
418  getRangeMap() const;
419 
449  void
450  apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
451  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
456 
457  private:
479  void
480  multiply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
481  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
482  const Teuchos::ETransp mode = Teuchos::NO_TRANS) const;
483 
484  public:
487 
488  // Attribute access functions
489 
491  magnitude_type getRelaxValue() const { return RelaxValue_; }
492 
494  magnitude_type getAbsoluteThreshold() const { return Athresh_; }
495 
497  magnitude_type getRelativeThreshold() const { return Rthresh_; }
498 
500  int getLevelOfFill() const { return LevelOfFill_; }
501 
503  Tpetra::CombineMode getOverlapMode() {
505  true, std::logic_error,
506  "Ifpack2::RILUK::SetOverlapMode: "
507  "RILUK no longer implements overlap on its own. "
508  "Use RILUK with AdditiveSchwarz if you want overlap.");
509  }
510 
512  Tpetra::global_size_t getGlobalNumEntries() const {
513  return getL().getGlobalNumEntries() + getU().getGlobalNumEntries();
514  }
515 
519  node_type>,
520  kk_handle_type> >
521  getGraph() const {
522  return Graph_;
523  }
524 
526  const crs_matrix_type& getL() const;
527 
529  const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>&
530  getD() const;
531 
533  const crs_matrix_type& getU() const;
534 
537 
545 
546  private:
547  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
550 
551  void allocateSolvers();
552  void allocate_L_and_U();
553  static void checkOrderingConsistency(const row_matrix_type& A);
554  void initAllValues(const row_matrix_type& A);
555 
556  void compute_serial();
557  void compute_kkspiluk();
558 // Workaround Cuda limitation of KOKKOS_LAMBDA in private/protected member functions
559 #ifdef KOKKOS_ENABLE_CUDA
560  public:
561 #endif
562  void compute_kkspiluk_stream();
563 #ifdef KOKKOS_ENABLE_CUDA
564  private:
565 #endif
566 
567  protected:
568  typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vec_type;
569 
572 
575  std::vector<Teuchos::RCP<iluk_graph_type> > Graph_v_;
581  Teuchos::RCP<crs_matrix_type> A_local_crs_nc_;
582  std::vector<local_matrix_device_type> A_local_diagblks;
583  std::vector<lno_row_view_t> A_local_diagblks_rowmap_v_;
584  std::vector<lno_nonzero_view_t> A_local_diagblks_entries_v_;
585  std::vector<scalar_nonzero_view_t> A_local_diagblks_values_v_;
586 
589  std::vector<Teuchos::RCP<crs_matrix_type> > L_v_;
594  std::vector<Teuchos::RCP<crs_matrix_type> > U_v_;
597 
600 
601  int LevelOfFill_;
602  double Overalloc_;
603 
604  bool isAllocated_;
605  bool isInitialized_;
606  bool isComputed_;
607 
608  int numInitialize_;
609  int numCompute_;
610  mutable int numApply_;
611 
612  double initializeTime_;
613  double computeTime_;
614  mutable double applyTime_;
615 
616  magnitude_type RelaxValue_;
617  magnitude_type Athresh_;
618  magnitude_type Rthresh_;
619 
622  Teuchos::RCP<kk_handle_type> KernelHandle_;
623  std::vector<Teuchos::RCP<kk_handle_type> > KernelHandle_v_;
624  bool isKokkosKernelsStream_;
625  int num_streams_;
626  std::vector<execution_space> exec_space_instances_;
627  bool hasStreamReordered_;
628  std::vector<typename lno_nonzero_view_t::non_const_type> perm_v_;
629  std::vector<typename lno_nonzero_view_t::non_const_type> reverse_perm_v_;
630  mutable std::unique_ptr<MV> Y_tmp_;
631  mutable std::unique_ptr<MV> reordered_x_;
632  mutable std::unique_ptr<MV> reordered_y_;
633 };
634 
635 // NOTE (mfh 11 Feb 2015) This used to exist in order to deal with
636 // different behavior of Tpetra::Crs{Graph,Matrix} for
637 // KokkosClassic::ThrustGPUNode. In particular, fillComplete on a
638 // CrsMatrix used to make the graph go away by default, so we had to
639 // pass in a parameter to keep a host copy of the graph. With the new
640 // (Kokkos refactor) version of Tpetra, this problem has gone away.
641 namespace detail {
642 template <class MatrixType, class NodeType>
643 struct setLocalSolveParams {
645  setParams(const Teuchos::RCP<Teuchos::ParameterList>& param) {
646  return param;
647  }
648 };
649 } // namespace detail
650 
651 } // namespace Ifpack2
652 
653 #endif /* IFPACK2_CRSRILUK_DECL_HPP */
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:491
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack2_RILUK_decl.hpp:494
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:521
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:228
node_type::execution_space execution_space
The Kokkos execution space of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:241
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:420
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:229
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:588
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:166
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:450
#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:248
Teuchos::RCP< iluk_graph_type > Graph_
The ILU(k) graph.
Definition: Ifpack2_RILUK_decl.hpp:574
double getApplyTime() const
Total time in seconds taken by all successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:366
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_RILUK_def.hpp:117
std::string description() const
A one-line description of this object.
Definition: Ifpack2_RILUK_def.hpp:1380
void setParameters(const Teuchos::ParameterList &params)
Definition: Ifpack2_RILUK_def.hpp:290
bool isKokkosKernelsSpiluk_
Optional KokkosKernels implementation.
Definition: Ifpack2_RILUK_decl.hpp:621
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:191
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition: Ifpack2_RILUK_decl.hpp:503
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:414
bool isComputed() const
Whether compute() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:340
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:250
node_type::device_type device_type
The Kokkos device type of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:238
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:179
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:232
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:208
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:226
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_RILUK_def.hpp:94
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition: Ifpack2_RILUK_decl.hpp:596
int getLevelOfFill() const
Get level of fill (the &quot;k&quot; in ILU(k)).
Definition: Ifpack2_RILUK_decl.hpp:500
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:1054
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition: Ifpack2_RILUK_decl.hpp:497
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:1131
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_RILUK_decl.hpp:235
crs_matrix_type::impl_scalar_type impl_scalar_type
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector)
Definition: Ifpack2_RILUK_decl.hpp:260
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:336
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
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:571
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition: Ifpack2_RILUK_def.hpp:408
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:362
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition: Ifpack2_RILUK_decl.hpp:599
int getNumCompute() const
Number of successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:349
int getNumApply() const
Number of successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:353
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack2_RILUK_decl.hpp:512
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object. ...
Definition: Ifpack2_RILUK_decl.hpp:358
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:150
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:579
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:593
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:345
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:223