40 #ifndef TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP 41 #define TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP 43 #include "Tpetra_LocalCrsMatrixOperator_fwd.hpp" 44 #include "Tpetra_LocalOperator.hpp" 45 #include "KokkosSparse_CrsMatrix.hpp" 59 template<
class MultiVectorScalar,
class MatrixScalar,
class Device>
63 using mv_scalar_type =
65 using matrix_scalar_type =
66 typename LocalOperator<MatrixScalar, Device>::scalar_type;
71 using local_ordinal_type =
73 using execution_space =
typename Device::execution_space;
75 using local_matrix_type =
76 KokkosSparse::CrsMatrix<matrix_scalar_type,
80 using local_graph_type =
typename local_matrix_type::StaticCrsGraphType;
87 apply (Kokkos::View<
const mv_scalar_type**, array_layout,
88 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > X,
89 Kokkos::View<mv_scalar_type**, array_layout,
90 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Y,
91 const Teuchos::ETransp mode,
92 const mv_scalar_type alpha,
93 const mv_scalar_type beta)
const override;
97 Kokkos::View<
const mv_scalar_type**, array_layout,
98 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > X,
99 Kokkos::View<mv_scalar_type**, array_layout,
100 device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged> > Y,
101 const Teuchos::ETransp mode,
102 const mv_scalar_type alpha,
103 const mv_scalar_type beta)
const;
105 bool hasTransposeApply ()
const override;
107 const local_matrix_type& getLocalMatrix ()
const;
110 std::shared_ptr<local_matrix_type> A_;
115 #endif // TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP Namespace Tpetra contains the class and methods constituting the Tpetra library.
Abstract interface for local operators (e.g., matrices and preconditioners).
Abstract interface for local operators (e.g., matrices and preconditioners).
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...
int local_ordinal_type
Default value of Scalar template parameter.