10 #ifndef TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP
11 #define TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP
13 #include "Tpetra_LocalCrsMatrixOperator_fwd.hpp"
14 #include "Tpetra_LocalOperator.hpp"
15 #include "KokkosSparse_CrsMatrix.hpp"
29 template<
class MultiVectorScalar,
class MatrixScalar,
class Device>
33 using mv_scalar_type =
35 using matrix_scalar_type =
36 typename LocalOperator<MatrixScalar, Device>::scalar_type;
41 using local_ordinal_type =
43 using execution_space =
typename Device::execution_space;
45 using local_matrix_device_type =
46 KokkosSparse::CrsMatrix<matrix_scalar_type,
53 using local_cusparse_matrix_type =
54 KokkosSparse::CrsMatrix<matrix_scalar_type,
59 using local_graph_device_type =
typename local_matrix_device_type::StaticCrsGraphType;
62 using ordinal_view_type =
typename local_graph_device_type::entries_type::non_const_type;
65 LocalCrsMatrixOperator (
const std::shared_ptr<local_matrix_device_type>& A,
const ordinal_view_type& A_ordinal_rowptrs);
69 apply (Kokkos::View<
const mv_scalar_type**, array_layout,
70 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > X,
71 Kokkos::View<mv_scalar_type**, array_layout,
72 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Y,
73 const Teuchos::ETransp mode,
74 const mv_scalar_type alpha,
75 const mv_scalar_type beta)
const override;
79 Kokkos::View<
const mv_scalar_type**, array_layout,
80 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > X,
81 Kokkos::View<mv_scalar_type**, array_layout,
82 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Y,
83 const Teuchos::ETransp mode,
84 const mv_scalar_type alpha,
85 const mv_scalar_type beta)
const;
87 bool hasTransposeApply ()
const override;
89 const local_matrix_device_type& getLocalMatrixDevice ()
const;
92 std::shared_ptr<local_matrix_device_type> A_;
93 local_cusparse_matrix_type A_cusparse;
94 const bool have_A_cusparse;
99 #endif // TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP
Abstract interface for local operators (e.g., matrices and preconditioners).
Abstract interface for local operators (e.g., matrices and preconditioners).
int local_ordinal_type
Default value of Scalar template parameter.
void applyImbalancedRows(Kokkos::View< const mv_scalar_type **, array_layout, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > X, Kokkos::View< mv_scalar_type **, array_layout, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > Y, const Teuchos::ETransp mode, const mv_scalar_type alpha, const mv_scalar_type beta) const
Same behavior as apply() above, except give KokkosKernels a hint to use an SPMV algorithm that can ef...