Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_LocalCrsMatrixOperator_decl.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP
11 #define TPETRA_LOCALCRSMATRIXOPERATOR_DECL_HPP
12 
13 #include "Tpetra_LocalCrsMatrixOperator_fwd.hpp"
14 #include "Tpetra_LocalOperator.hpp"
15 #include "KokkosSparse_CrsMatrix.hpp"
16 #include <memory> // std::shared_ptr
17 
18 namespace Tpetra {
19 
29  template<class MultiVectorScalar, class MatrixScalar, class Device>
31  public LocalOperator<MultiVectorScalar, Device> {
32  private:
33  using mv_scalar_type =
35  using matrix_scalar_type =
36  typename LocalOperator<MatrixScalar, Device>::scalar_type;
37  using array_layout =
39  using device_type =
41  using local_ordinal_type =
43  using execution_space = typename Device::execution_space;
44  public:
45  using local_matrix_device_type =
46  KokkosSparse::CrsMatrix<matrix_scalar_type,
47  local_ordinal_type,
48  device_type,
49  void,
50  size_t>;
51  private:
52  //The type of a matrix with offset=ordinal, but otherwise the same as local_matrix_device_type
53  using local_cusparse_matrix_type =
54  KokkosSparse::CrsMatrix<matrix_scalar_type,
55  local_ordinal_type,
56  device_type,
57  void,
58  local_ordinal_type>;
59  using local_graph_device_type = typename local_matrix_device_type::StaticCrsGraphType;
60 
61  public:
62  using ordinal_view_type = typename local_graph_device_type::entries_type::non_const_type;
63 
64  LocalCrsMatrixOperator (const std::shared_ptr<local_matrix_device_type>& A);
65  LocalCrsMatrixOperator (const std::shared_ptr<local_matrix_device_type>& A, const ordinal_view_type& A_ordinal_rowptrs);
66  ~LocalCrsMatrixOperator () override = default;
67 
68  void
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;
76 
77  void
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;
86 
87  bool hasTransposeApply () const override;
88 
89  const local_matrix_device_type& getLocalMatrixDevice () const;
90 
91  private:
92  std::shared_ptr<local_matrix_device_type> A_;
93  local_cusparse_matrix_type A_cusparse;
94  const bool have_A_cusparse;
95  };
96 
97 } // namespace Tpetra
98 
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...