46 #ifndef IFPACK2_CRSRILUK_DECL_HPP 47 #define IFPACK2_CRSRILUK_DECL_HPP 51 #include "Tpetra_CrsMatrix_decl.hpp" 54 #include "Ifpack2_LocalSparseTriangularSolver_decl.hpp" 56 #include <type_traits> 59 #ifdef IFPACK2_ILUK_EXPERIMENTAL 60 #include <Kokkos_Core.hpp> 61 #include <shylubasker_decl.hpp> 62 # ifdef IFPACK2_HTS_EXPERIMENTAL 63 # include <ShyLUHTS_config.h> 64 # include <shylu_hts_decl.hpp> 253 template<
class MatrixType>
256 typename MatrixType::local_ordinal_type,
257 typename MatrixType::global_ordinal_type,
258 typename MatrixType::node_type>,
260 typename MatrixType::local_ordinal_type,
261 typename MatrixType::global_ordinal_type,
262 typename MatrixType::node_type> >
278 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
287 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.");
295 template <
class NewMatrixType>
friend class RILUK;
300 RILUK (
const Teuchos::RCP<const row_matrix_type>& A_in);
309 RILUK (
const Teuchos::RCP<const crs_matrix_type>& A_in);
322 template <
typename NewMatrixType>
323 Teuchos::RCP< RILUK< NewMatrixType > >
324 clone (
const Teuchos::RCP<const NewMatrixType>& A_newnode)
const;
353 return isInitialized_;
362 return numInitialize_;
375 return initializeTime_;
412 setMatrix (
const Teuchos::RCP<const row_matrix_type>& A);
426 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
430 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
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;
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;
499 Teuchos::RCP<const row_matrix_type>
getMatrix ()
const;
517 TEUCHOS_TEST_FOR_EXCEPTION(
518 true, std::logic_error,
"Ifpack2::RILUK::SetOverlapMode: " 519 "RILUK no longer implements overlap on its own. " 520 "Use RILUK with AdditiveSchwarz if you want overlap.");
525 return getL ().getGlobalNumEntries () +
getU ().getGlobalNumEntries ();
529 Teuchos::RCP<Ifpack2::IlukGraph<Tpetra::CrsGraph<local_ordinal_type,global_ordinal_type,node_type> > >
getGraph ()
const {
537 const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>&
544 Teuchos::RCP<const crs_matrix_type>
getCrsMatrix ()
const;
547 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
548 typedef Teuchos::ScalarTraits<scalar_type> STS;
549 typedef Teuchos::ScalarTraits<magnitude_type> STM;
551 void allocate_L_and_U ();
560 static Teuchos::RCP<const row_matrix_type>
561 makeLocalFilter (
const Teuchos::RCP<const row_matrix_type>& A);
564 typedef Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> vec_type;
567 Teuchos::RCP<const row_matrix_type>
A_;
579 Teuchos::RCP<crs_matrix_type>
L_;
581 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> >
L_solver_;
583 Teuchos::RCP<crs_matrix_type>
U_;
585 Teuchos::RCP<LocalSparseTriangularSolver<row_matrix_type> >
U_solver_;
587 Teuchos::RCP<vec_type>
D_;
590 #ifdef IFPACK2_ILUK_EXPERIMENTAL 591 typedef typename node_type::device_type kokkos_device;
592 typedef typename kokkos_device::execution_space kokkos_exe;
594 static_assert( std::is_same< kokkos_exe,
595 Kokkos::OpenMP>::value,
596 "Kokkos node type not supported by exepertimentalthread basker RILUK decl");
599 Teuchos::RCP<const crs_matrix_type> A_local_crs_;
602 void initLocalCrs ();
604 Teuchos::RCP< BaskerNS::Basker<local_ordinal_type, scalar_type, Kokkos::OpenMP> >
610 # ifdef IFPACK2_HTS_EXPERIMENTAL 613 typedef ::Experimental::HTS<local_ordinal_type, local_ordinal_type, scalar_type> HTST;
614 struct HTSData :
public HTST::Deallocator {
620 virtual void free_CrsMatrix_data();
622 bool operator< (
const Entry& o)
const {
return j < o.j; }
628 typename HTST::Impl* Limpl, * Uimpl;
632 Teuchos::RCP<HTSManager> hts_mgr_;
641 bool isExperimental_;
645 mutable int numApply_;
647 double initializeTime_;
649 mutable double applyTime_;
664 template<
class MatrixType,
class NodeType>
665 struct setLocalSolveParams{
666 static Teuchos::RCP<Teuchos::ParameterList>
667 setParams (
const Teuchos::RCP<Teuchos::ParameterList>& param) {
673 template <
class MatrixType>
674 template <
typename NewMatrixType>
675 Teuchos::RCP<RILUK<NewMatrixType> >
677 clone (
const Teuchos::RCP<const NewMatrixType>& A_newnode)
const 679 using Teuchos::ParameterList;
683 typedef typename NewMatrixType::scalar_type new_scalar_type;
684 typedef typename NewMatrixType::local_ordinal_type new_local_ordinal_type;
685 typedef typename NewMatrixType::global_ordinal_type new_global_ordinal_type;
686 typedef typename NewMatrixType::node_type new_node_type;
687 typedef Tpetra::RowMatrix<new_scalar_type, new_local_ordinal_type,
688 new_global_ordinal_type, new_node_type> new_row_matrix_type;
691 RCP<new_riluk_type> new_riluk = rcp (
new new_riluk_type (A_newnode));
693 RCP<ParameterList> plClone = Teuchos::parameterList ();
694 plClone = detail::setLocalSolveParams<NewMatrixType, new_node_type>::setParams (plClone);
696 RCP<new_node_type> new_node = A_newnode->getNode ();
697 new_riluk->L_ = L_->clone (new_node, plClone);
698 new_riluk->U_ = U_->clone (new_node, plClone);
699 new_riluk->D_ = D_->clone (new_node);
701 new_riluk->LevelOfFill_ = LevelOfFill_;
703 new_riluk->isAllocated_ = isAllocated_;
704 new_riluk->isInitialized_ = isInitialized_;
705 new_riluk->isComputed_ = isComputed_;
707 new_riluk->numInitialize_ = numInitialize_;
708 new_riluk->numCompute_ = numCompute_;
709 new_riluk->numApply_ = numApply_;
711 new_riluk->RelaxValue_ = RelaxValue_;
712 new_riluk->Athresh_ = Athresh_;
713 new_riluk->Rthresh_ = Rthresh_;
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
Declaration and definition of IlukGraph.
Ifpack2::ScalingType enumerable type.
double getComputeTime() const
Total time in seconds taken by all successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:378
int getNumApply() const
Number of successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:369
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:111
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:272
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:579
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:502
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition: Ifpack2_RILUK_def.hpp:455
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:284
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:462
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_RILUK_def.hpp:106
void setParameters(const Teuchos::ParameterList ¶ms)
Definition: Ifpack2_RILUK_def.hpp:255
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:254
Tpetra::CombineMode getOverlapMode()
Get overlap mode type.
Definition: Ifpack2_RILUK_decl.hpp:516
bool isComputed() const
Whether compute() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:356
double getApplyTime() const
Total time in seconds taken by all successful apply() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:382
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:287
magnitude_type getAbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack2_RILUK_decl.hpp:507
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > getGraph() const
Return the Ifpack2::IlukGraph associated with this factored matrix.
Definition: Ifpack2_RILUK_decl.hpp:529
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:275
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:159
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:269
virtual ~RILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_RILUK_def.hpp:101
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition: Ifpack2_RILUK_decl.hpp:585
int getNumCompute() const
Number of successful compute() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:365
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:189
magnitude_type getRelaxValue() const
Get RILU(k) relaxation parameter.
Definition: Ifpack2_RILUK_decl.hpp:504
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > Graph_
The ILU(k) graph.
Definition: Ifpack2_RILUK_decl.hpp:572
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:107
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_RILUK_def.hpp:849
const crs_matrix_type & getU() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:173
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_RILUK_decl.hpp:278
Teuchos::RCP< RILUK< NewMatrixType > > clone(const Teuchos::RCP< const NewMatrixType > &A_newnode) const
Clone preconditioner to a new node type.
Definition: Ifpack2_RILUK_decl.hpp:677
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:97
Declaration of interface for preconditioners that can change their matrix after construction.
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition: Ifpack2_RILUK_decl.hpp:581
bool isInitialized() const
Whether initialize() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:352
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:567
Tpetra::global_size_t getGlobalNumEntries() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack2_RILUK_decl.hpp:524
Definition: Ifpack2_Container.hpp:761
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition: Ifpack2_RILUK_decl.hpp:587
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:1138
std::string description() const
A one-line description of this object.
Definition: Ifpack2_RILUK_def.hpp:1327
magnitude_type getRelativeThreshold() const
Get relative threshold value.
Definition: Ifpack2_RILUK_decl.hpp:510
int getNumInitialize() const
Number of successful initialize() calls for this object.
Definition: Ifpack2_RILUK_decl.hpp:361
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition: Ifpack2_RILUK_decl.hpp:513
double getInitializeTime() const
Total time in seconds taken by all successful initialize() calls for this object. ...
Definition: Ifpack2_RILUK_decl.hpp:374
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:576
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:583
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:208
const crs_matrix_type & getL() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_RILUK_def.hpp:142
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_RILUK_decl.hpp:266