12#ifndef __IFPACK2_FIC_DECL_HPP__
13#define __IFPACK2_FIC_DECL_HPP__
15#include <Ifpack2_Details_FastILU_Base.hpp>
18template<
typename LocalOrdinal,
typename Scalar,
typename execution_space>
28template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
36 typedef FastICPrec<LocalOrdinal, ImplScalar, typename Base::execution_space> LocalFIC;
39 Fic(Teuchos::RCP<const TRowMatrix> mat_);
54 mutable Teuchos::RCP<LocalFIC> localPrec_;
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type ImplScalar
Kokkos scalar type.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:48
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TRowMatrix
Tpetra row matrix.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:50
FastILU_Base(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_FastILU_Base_def.hpp:31
Kokkos::View< ImplScalar *, execution_space > ImplScalarArray
Array of Scalar on device.
Definition Ifpack2_Details_FastILU_Base_decl.hpp:62
int getNTrisol() const
Get the number of triangular solves ("nTrisol") from localPrec_.
Definition Ifpack2_Details_Fic_def.hpp:46
void checkLocalIC() const
Verify and print debug info about the internal IC preconditioner.
Definition Ifpack2_Details_Fic_def.hpp:53
int getSweeps() const
Get the sweeps ("nFact") from localPrec_.
Definition Ifpack2_Details_Fic_def.hpp:32
void computeLocalPrec()
Get values array from the matrix and then call compute() on the underlying preconditioner.
Definition Ifpack2_Details_Fic_def.hpp:72
std::string getSpTrsvType() const
Get the name of triangular solve algorithm.
Definition Ifpack2_Details_Fic_def.hpp:39
void applyLocalPrec(ImplScalarArray x, ImplScalarArray y) const
Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only).
Definition Ifpack2_Details_Fic_def.hpp:82
Fic(Teuchos::RCP< const TRowMatrix > mat_)
Constructor.
Definition Ifpack2_Details_Fic_def.hpp:27
void initLocalPrec()
Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->ini...
Definition Ifpack2_Details_Fic_def.hpp:60
std::string getName() const
Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic").
Definition Ifpack2_Details_Fic_def.hpp:91
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41