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 {
29 class ParameterList; // forward declaration
30}
31
32namespace Ifpack2 {
33
212template<class MatrixType>
213class RILUK:
214 virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
215 typename MatrixType::local_ordinal_type,
216 typename MatrixType::global_ordinal_type,
217 typename MatrixType::node_type>,
218 virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
219 typename MatrixType::local_ordinal_type,
220 typename MatrixType::global_ordinal_type,
221 typename MatrixType::node_type> >
222{
223 public:
225 typedef typename MatrixType::scalar_type scalar_type;
226
228 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
229
231 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
232
234 typedef typename MatrixType::node_type node_type;
235
237 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
238
240 typedef typename node_type::device_type device_type;
241
243 typedef typename node_type::execution_space execution_space;
244
246 typedef Tpetra::RowMatrix<scalar_type,
250
251
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.");
253
255 typedef Tpetra::CrsMatrix<scalar_type,
259
261 typedef typename crs_matrix_type::impl_scalar_type impl_scalar_type;
262
263 template <class NewMatrixType> 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
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;
273
274
276
278
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;
290
294 RILUK (const Teuchos::RCP<const row_matrix_type>& A_in);
295
303 RILUK (const Teuchos::RCP<const crs_matrix_type>& A_in);
304
305 private:
308 RILUK (const RILUK<MatrixType> & src);
309
310 public:
312 virtual ~RILUK ();
313
322 void setParameters (const Teuchos::ParameterList& params);
323
325 void initialize ();
326
335 void compute ();
336
338 bool isInitialized () const {
339 return isInitialized_;
340 }
341
342 bool isComputed () const {
343 return isComputed_;
344 }
345
347 int getNumInitialize () const {
348 return numInitialize_;
349 }
350
351 int getNumCompute () const {
352 return numCompute_;
353 }
354
355 int getNumApply () const {
356 return numApply_;
357 }
358
360 double getInitializeTime () const {
361 return initializeTime_;
362 }
363
364 double getComputeTime () const {
365 return computeTime_;
366 }
367
368 double getApplyTime () const {
369 return applyTime_;
370 }
371
373 size_t getNodeSmootherComplexity() const;
374
375
376
378
379
402 virtual void
403 setMatrix (const Teuchos::RCP<const row_matrix_type>& A);
404
406
408
410 std::string description () const;
411
413
415
417 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
418 getDomainMap () const;
419
421 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
422 getRangeMap () const;
423
453 void
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,
456 Teuchos::ETransp mode = Teuchos::NO_TRANS,
457 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one (),
458 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero ()) const;
460
461private:
483 void
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,
486 const Teuchos::ETransp mode = Teuchos::NO_TRANS) const;
487public:
488
490 Teuchos::RCP<const row_matrix_type> getMatrix () const;
491
492 // Attribute access functions
493
495 magnitude_type getRelaxValue () const { return RelaxValue_; }
496
498 magnitude_type getAbsoluteThreshold () const { return Athresh_; }
499
501 magnitude_type getRelativeThreshold () const {return Rthresh_;}
502
504 int getLevelOfFill () const { return LevelOfFill_; }
505
507 Tpetra::CombineMode getOverlapMode () {
508 TEUCHOS_TEST_FOR_EXCEPTION(
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.");
512 }
513
515 Tpetra::global_size_t getGlobalNumEntries () const {
516 return getL ().getGlobalNumEntries () + getU ().getGlobalNumEntries ();
517 }
518
520 Teuchos::RCP<Ifpack2::IlukGraph<Tpetra::CrsGraph<local_ordinal_type,
522 node_type>, kk_handle_type> > getGraph () const {
523 return Graph_;
524 }
525
527 const crs_matrix_type& getL () const;
528
530 const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>&
531 getD () const;
532
534 const crs_matrix_type& getU () const;
535
537 Teuchos::RCP<const crs_matrix_type> getCrsMatrix () const;
538
544 static Teuchos::RCP<const row_matrix_type>
545 makeLocalFilter (const Teuchos::RCP<const row_matrix_type>& A);
546
547private:
548 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
549 typedef Teuchos::ScalarTraits<scalar_type> STS;
550 typedef Teuchos::ScalarTraits<magnitude_type> STM;
551
552 void allocateSolvers ();
553 void allocate_L_and_U ();
554 static void checkOrderingConsistency (const row_matrix_type& A);
555 void initAllValues (const row_matrix_type& A);
556
557 void compute_serial();
558 void compute_kkspiluk();
559// Workaround Cuda limitation of KOKKOS_LAMBDA in private/protected member functions
560#ifdef KOKKOS_ENABLE_CUDA
561public:
562#endif
563 void compute_kkspiluk_stream();
564#ifdef KOKKOS_ENABLE_CUDA
565private:
566#endif
567
568protected:
569 typedef Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> vec_type;
570
572 Teuchos::RCP<const row_matrix_type> A_;
573
575 Teuchos::RCP<iluk_graph_type> Graph_;
576 std::vector< Teuchos::RCP<iluk_graph_type> > Graph_v_;
580 Teuchos::RCP<const row_matrix_type> A_local_;
581 Teuchos::RCP<const crs_matrix_type> A_local_crs_;
582 Teuchos::RCP<crs_matrix_type> A_local_crs_nc_;
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_;
587
589 Teuchos::RCP<crs_matrix_type> L_;
590 std::vector< Teuchos::RCP<crs_matrix_type> > L_v_;
592 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> > L_solver_;
594 Teuchos::RCP<crs_matrix_type> U_;
595 std::vector< Teuchos::RCP<crs_matrix_type> > U_v_;
597 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> > U_solver_;
598
600 Teuchos::RCP<vec_type> D_;
601
602 int LevelOfFill_;
603 double Overalloc_;
604
605 bool isAllocated_;
606 bool isInitialized_;
607 bool isComputed_;
608
609 int numInitialize_;
610 int numCompute_;
611 mutable int numApply_;
612
613 double initializeTime_;
614 double computeTime_;
615 mutable double applyTime_;
616
617 magnitude_type RelaxValue_;
618 magnitude_type Athresh_;
619 magnitude_type Rthresh_;
620
623 Teuchos::RCP<kk_handle_type> KernelHandle_;
624 std::vector< Teuchos::RCP<kk_handle_type> > KernelHandle_v_;
625 bool isKokkosKernelsStream_;
626 int num_streams_;
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_;
634};
635
636// NOTE (mfh 11 Feb 2015) This used to exist in order to deal with
637// different behavior of Tpetra::Crs{Graph,Matrix} for
638// KokkosClassic::ThrustGPUNode. In particular, fillComplete on a
639// CrsMatrix used to make the graph go away by default, so we had to
640// pass in a parameter to keep a host copy of the graph. With the new
641// (Kokkos refactor) version of Tpetra, this problem has gone away.
642namespace detail {
643 template<class MatrixType, class NodeType>
644 struct setLocalSolveParams{
645 static Teuchos::RCP<Teuchos::ParameterList>
646 setParams (const Teuchos::RCP<Teuchos::ParameterList>& param) {
647 return param;
648 }
649 };
650} // namespace detail
651
652} // namespace Ifpack2
653
654#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:75
crs_matrix_type::impl_scalar_type impl_scalar_type
Scalar type stored in Kokkos::Views (CrsMatrix and MultiVector).
Definition Ifpack2_RILUK_decl.hpp:261
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:481
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:258
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:368
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
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:364
int getNumApply() const
Number of successful apply() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:355
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
magnitude_type getRelaxValue() const
Get RILU(k) relaxation parameter.
Definition Ifpack2_RILUK_decl.hpp:495
int getNumCompute() const
Number of successful compute() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:351
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition Ifpack2_RILUK_decl.hpp:572
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition Ifpack2_RILUK_decl.hpp:507
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition Ifpack2_RILUK_decl.hpp:594
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition Ifpack2_RILUK_def.hpp:192
std::string description() const
A one-line description of this object.
Definition Ifpack2_RILUK_def.hpp:1434
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition Ifpack2_RILUK_decl.hpp:515
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_RILUK_def.hpp:125
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:225
bool isComputed() const
Whether compute() has been called on this object.
Definition Ifpack2_RILUK_decl.hpp:342
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
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:234
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< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition Ifpack2_RILUK_def.hpp:434
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition Ifpack2_RILUK_decl.hpp:498
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition Ifpack2_RILUK_decl.hpp:592
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition Ifpack2_RILUK_decl.hpp:504
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
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
void setParameters(const Teuchos::ParameterList &params)
Definition Ifpack2_RILUK_def.hpp:306
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition Ifpack2_RILUK_decl.hpp:501
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:347
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:228
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition Ifpack2_RILUK_decl.hpp:589
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:231
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:240
node_type::execution_space execution_space
The Kokkos execution space of the input MatrixType.
Definition Ifpack2_RILUK_decl.hpp:243
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition Ifpack2_RILUK_decl.hpp:597
Teuchos::RCP< iluk_graph_type > Graph_
The ILU(k) graph.
Definition Ifpack2_RILUK_decl.hpp:575
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_RILUK_def.hpp:205
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_RILUK_def.hpp:1099
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
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition Ifpack2_RILUK_decl.hpp:600
bool isKokkosKernelsSpiluk_
Optional KokkosKernels implementation.
Definition Ifpack2_RILUK_decl.hpp:622
bool isInitialized() const
Whether initialize() has been called on this object.
Definition Ifpack2_RILUK_decl.hpp:338
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object.
Definition Ifpack2_RILUK_decl.hpp:360
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< 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
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41