10 #ifndef IFPACK2_LOCALSPARSETRIANGULARSOLVER_DECL_HPP
11 #define IFPACK2_LOCALSPARSETRIANGULARSOLVER_DECL_HPP
15 #include "Tpetra_CrsMatrix_fwd.hpp"
16 #include "Teuchos_FancyOStream.hpp"
17 #include <type_traits>
19 #include "KokkosSparse_sptrsv.hpp"
45 template<
class MatrixType>
48 typename MatrixType::local_ordinal_type,
49 typename MatrixType::global_ordinal_type,
50 typename MatrixType::node_type>,
52 typename MatrixType::local_ordinal_type,
53 typename MatrixType::global_ordinal_type,
54 typename MatrixType::node_type> >
69 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type>
map_type;
77 static_assert (std::is_same<MatrixType, row_matrix_type>::value,
78 "Ifpack2::LocalSparseTriangularSolver: The template parameter "
79 "MatrixType must be a Tpetra::RowMatrix specialization. "
80 "Please don't use Tpetra::CrsMatrix (a subclass of "
81 "Tpetra::RowMatrix) here anymore. The constructor can take "
82 "either a RowMatrix or a CrsMatrix just fine.");
85 using local_matrix_device_type =
typename crs_matrix_type::local_matrix_device_type;
86 using local_matrix_graph_device_type =
typename local_matrix_device_type::StaticCrsGraphType;
87 using lno_row_view_t =
typename local_matrix_graph_device_type::row_map_type;
88 using lno_nonzero_view_t =
typename local_matrix_graph_device_type::entries_type;
89 using scalar_nonzero_view_t =
typename local_matrix_device_type::values_type;
90 using TemporaryMemorySpace =
typename local_matrix_graph_device_type::device_type::memory_space;
91 using PersistentMemorySpace =
typename local_matrix_graph_device_type::device_type::memory_space;
92 using HandleExecSpace =
typename local_matrix_graph_device_type::device_type::execution_space;
93 using k_handle =
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, HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace >;
179 return isInitialized_;
212 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
213 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
312 void setStreamInfo (
const bool& isKokkosKernelsStream,
const int& num_streams,
const std::vector<HandleExecSpace>& exec_space_instances);
329 std::vector< Teuchos::RCP<crs_matrix_type> > A_crs_v_;
331 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
346 bool isInternallyChanged_;
347 bool reverseStorage_;
349 mutable int numInitialize_;
350 mutable int numCompute_;
351 mutable int numApply_;
353 double initializeTime_;
362 bool isKokkosKernelsSptrsv_;
364 std::vector< Teuchos::RCP<k_handle> > kh_v_;
366 bool isKokkosKernelsStream_;
368 std::vector<HandleExecSpace> exec_space_instances_;
397 localApply (
const MV& X,
405 localTriangularSolve (
const MV& Y,
409 void initializeState();
411 KokkosSparse::Experimental::SPTRSVAlgorithm kokkosKernelsAlgorithm()
const;
416 #endif // IFPACK2_LOCALSPARSETRIANGULARSOLVER_DECL_HPP
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
double getApplyFlops() const
Return the number of flops for the application of the preconditioner.
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 preconditioner to X, and put the result in Y.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:706
int getNumCompute() const
Return the number of calls to compute().
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1050
std::string description() const
A one-line description of this object.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1085
void compute()
"Numeric" phase of setup
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:602
Teuchos::RCP< const row_matrix_type > getMatrix() const
The original input matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:243
MatrixType::global_ordinal_type global_ordinal_type
Type of the global indices of the input matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:62
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:189
double getApplyTime() const
Return the time spent in apply().
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1078
double getComputeTime() const
Return the time spent in compute().
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1071
double getInitializeTime() const
Return the time spent in initialize().
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1064
double getComputeFlops() const
Return the number of flops in the computation phase.
MatrixType::node_type node_type
Node type of the input matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:64
void setMatrices(const std::vector< Teuchos::RCP< crs_matrix_type > > &A_crs_v)
Set this preconditioner's matrices (used by stream interface of triangular solve).
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1240
void applyMat(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) const
Apply the original input matrix.
Teuchos::RCP< const map_type > getDomainMap() const
The domain of this operator.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1158
"Preconditioner" that solves local sparse triangular systems.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:46
void initialize()
"Symbolic" phase of setup
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:355
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
int getNumInitialize() const
Return the number of calls to initialize().
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1043
void setStreamInfo(const bool &isKokkosKernelsStream, const int &num_streams, const std::vector< HandleExecSpace > &exec_space_instances)
Set this triangular solver's stream information.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1229
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Specialization of Tpetra::Map used by this class.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:69
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Set this preconditioner's matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1183
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
This operator's communicator.
Declaration of interface for preconditioners that can change their matrix after construction.
static const EVerbosityLevel verbLevel_default
int getNumApply() const
Return the number of calls to apply().
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1057
MatrixType::local_ordinal_type local_ordinal_type
Type of the local indices of the input matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:60
Teuchos::RCP< const map_type > getRangeMap() const
The range of this operator.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1171
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Specialization of Tpetra::CrsMatrix used by this class.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:75
void setParameters(const Teuchos::ParameterList ¶ms)
Set this object's parameters.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:308
MatrixType::scalar_type scalar_type
Type of the entries of the input matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:58
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print this object with given verbosity to the given output stream.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:1120
MatrixType::mag_type magnitude_type
Type of the absolute value (magnitude) of a scalar_type value.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:67
virtual ~LocalSparseTriangularSolver()
Destructor (virtual for memory safety).
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:289
LocalSparseTriangularSolver()
Constructor that takes no input matrix.
Definition: Ifpack2_LocalSparseTriangularSolver_def.hpp:250
bool isInitialized() const
Return true if the preconditioner has been successfully initialized.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:178
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Specialization of Tpetra::RowMatrix used by this class.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:72