Ginkgo  Generated from pipelines/1330831941 branch based on master. Ginkgo version 1.8.0
A numerical linear algebra library targeting many-core architectures
sparsity_csr.hpp
1 // SPDX-FileCopyrightText: 2017 - 2024 The Ginkgo authors
2 //
3 // SPDX-License-Identifier: BSD-3-Clause
4 
5 #ifndef GKO_PUBLIC_CORE_MATRIX_SPARSITY_CSR_HPP_
6 #define GKO_PUBLIC_CORE_MATRIX_SPARSITY_CSR_HPP_
7 
8 
9 #include <vector>
10 
11 
12 #include <ginkgo/core/base/array.hpp>
13 #include <ginkgo/core/base/lin_op.hpp>
14 #include <ginkgo/core/base/polymorphic_object.hpp>
15 
16 
17 namespace gko {
18 namespace matrix {
19 
20 
21 template <typename ValueType, typename IndexType>
22 class Csr;
23 
24 
25 template <typename ValueType>
26 class Dense;
27 
28 
29 template <typename ValueType, typename IndexType>
30 class Fbcsr;
31 
32 
51 template <typename ValueType = default_precision, typename IndexType = int32>
52 class SparsityCsr : public EnableLinOp<SparsityCsr<ValueType, IndexType>>,
53  public ConvertibleTo<Csr<ValueType, IndexType>>,
54  public ConvertibleTo<Dense<ValueType>>,
55  public ReadableFromMatrixData<ValueType, IndexType>,
56  public WritableToMatrixData<ValueType, IndexType>,
57  public Transposable {
58  friend class EnablePolymorphicObject<SparsityCsr, LinOp>;
59  friend class Csr<ValueType, IndexType>;
60  friend class Dense<ValueType>;
61  friend class Fbcsr<ValueType, IndexType>;
62 
63 public:
66  using ConvertibleTo<Csr<ValueType, IndexType>>::convert_to;
67  using ConvertibleTo<Csr<ValueType, IndexType>>::move_to;
68  using ConvertibleTo<Dense<ValueType>>::convert_to;
69  using ConvertibleTo<Dense<ValueType>>::move_to;
71 
72  using value_type = ValueType;
73  using index_type = IndexType;
74  using transposed_type = SparsityCsr<IndexType, ValueType>;
75  using mat_data = matrix_data<ValueType, IndexType>;
76  using device_mat_data = device_matrix_data<ValueType, IndexType>;
77 
78  void convert_to(Csr<ValueType, IndexType>* result) const override;
79 
80  void move_to(Csr<ValueType, IndexType>* result) override;
81 
82  void convert_to(Dense<ValueType>* result) const override;
83 
84  void move_to(Dense<ValueType>* result) override;
85 
86  void read(const mat_data& data) override;
87 
88  void read(const device_mat_data& data) override;
89 
90  void read(device_mat_data&& data) override;
91 
92  void write(mat_data& data) const override;
93 
94  std::unique_ptr<LinOp> transpose() const override;
95 
96  std::unique_ptr<LinOp> conj_transpose() const override;
97 
107  std::unique_ptr<SparsityCsr> to_adjacency_matrix() const;
108 
112  void sort_by_column_index();
113 
114  /*
115  * Tests if all col_idxs are sorted by column index
116  *
117  * @returns True if all col_idxs are sorted.
118  */
119  bool is_sorted_by_column_index() const;
120 
126  index_type* get_col_idxs() noexcept { return col_idxs_.get_data(); }
127 
135  const index_type* get_const_col_idxs() const noexcept
136  {
137  return col_idxs_.get_const_data();
138  }
139 
145  index_type* get_row_ptrs() noexcept { return row_ptrs_.get_data(); }
146 
154  const index_type* get_const_row_ptrs() const noexcept
155  {
156  return row_ptrs_.get_const_data();
157  }
158 
164  value_type* get_value() noexcept { return value_.get_data(); }
165 
173  const value_type* get_const_value() const noexcept
174  {
175  return value_.get_const_data();
176  }
177 
183  size_type get_num_nonzeros() const noexcept { return col_idxs_.get_size(); }
184 
192  static std::unique_ptr<SparsityCsr> create(
193  std::shared_ptr<const Executor> exec, const dim<2>& size = dim<2>{},
194  size_type num_nonzeros = {});
195 
215  static std::unique_ptr<SparsityCsr> create(
216  std::shared_ptr<const Executor> exec, const dim<2>& size,
217  array<index_type> col_idxs, array<index_type> row_ptrs,
218  value_type value = one<ValueType>());
219 
225  template <typename ColIndexType, typename RowPtrType>
226  GKO_DEPRECATED(
227  "explicitly construct the gko::array argument instead of passing "
228  "initializer lists")
229  static std::unique_ptr<SparsityCsr> create(
230  std::shared_ptr<const Executor> exec, const dim<2>& size,
231  std::initializer_list<ColIndexType> col_idxs,
232  std::initializer_list<RowPtrType> row_ptrs,
233  value_type value = one<ValueType>())
234  {
235  return create(exec, size, array<index_type>{exec, std::move(col_idxs)},
236  array<index_type>{exec, std::move(row_ptrs)}, value);
237  }
238 
246  static std::unique_ptr<SparsityCsr> create(
247  std::shared_ptr<const Executor> exec,
248  std::shared_ptr<const LinOp> matrix);
249 
263  static std::unique_ptr<const SparsityCsr> create_const(
264  std::shared_ptr<const Executor> exec, const dim<2>& size,
265  gko::detail::const_array_view<IndexType>&& col_idxs,
266  gko::detail::const_array_view<IndexType>&& row_ptrs,
267  ValueType value = one<ValueType>())
268  {
269  // cast const-ness away, but return a const object afterwards,
270  // so we can ensure that no modifications take place.
271  return std::unique_ptr<const SparsityCsr>(new SparsityCsr{
272  exec, size, gko::detail::array_const_cast(std::move(col_idxs)),
273  gko::detail::array_const_cast(std::move(row_ptrs)), value});
274  }
275 
281 
288 
293  SparsityCsr(const SparsityCsr&);
294 
301 
302 protected:
303  SparsityCsr(std::shared_ptr<const Executor> exec,
304  const dim<2>& size = dim<2>{}, size_type num_nonzeros = {});
305 
306  SparsityCsr(std::shared_ptr<const Executor> exec, const dim<2>& size,
307  array<index_type> col_idxs, array<index_type> row_ptrs,
308  value_type value = one<ValueType>());
309 
310  SparsityCsr(std::shared_ptr<const Executor> exec,
311  std::shared_ptr<const LinOp> matrix);
312 
313  void apply_impl(const LinOp* b, LinOp* x) const override;
314 
315  void apply_impl(const LinOp* alpha, const LinOp* b, const LinOp* beta,
316  LinOp* x) const override;
317 
318 private:
319  array<index_type> col_idxs_;
320  array<index_type> row_ptrs_;
321  array<value_type> value_;
322 };
323 
324 
325 } // namespace matrix
326 } // namespace gko
327 
328 
329 #endif // GKO_PUBLIC_CORE_MATRIX_SPARSITY_CSR_HPP_
gko::EnablePolymorphicAssignment< ConcreteLinOp >::move_to
void move_to(result_type *result) override
Definition: polymorphic_object.hpp:732
gko::EnablePolymorphicAssignment< ConcreteLinOp >::convert_to
void convert_to(result_type *result) const override
Definition: polymorphic_object.hpp:730
gko::matrix::SparsityCsr::get_const_col_idxs
const index_type * get_const_col_idxs() const noexcept
Returns the column indices of the matrix.
Definition: sparsity_csr.hpp:135
gko::ReadableFromMatrixData::read
virtual void read(const matrix_data< ValueType, IndexType > &data)=0
Reads a matrix from a matrix_data structure.
gko::matrix::SparsityCsr
SparsityCsr is a matrix format which stores only the sparsity pattern of a sparse matrix by compressi...
Definition: csr.hpp:40
gko::size_type
std::size_t size_type
Integral type used for allocation quantities.
Definition: types.hpp:108
gko::matrix::SparsityCsr::transpose
std::unique_ptr< LinOp > transpose() const override
Returns a LinOp representing the transpose of the Transposable object.
gko::matrix::SparsityCsr::get_const_row_ptrs
const index_type * get_const_row_ptrs() const noexcept
Returns the row pointers of the matrix.
Definition: sparsity_csr.hpp:154
gko::matrix::SparsityCsr::get_col_idxs
index_type * get_col_idxs() noexcept
Returns the column indices of the matrix.
Definition: sparsity_csr.hpp:126
gko::matrix::SparsityCsr::create
static std::unique_ptr< SparsityCsr > create(std::shared_ptr< const Executor > exec, const dim< 2 > &size=dim< 2 >{}, size_type num_nonzeros={})
Creates an uninitialized SparsityCsr matrix of the specified size.
gko
The Ginkgo namespace.
Definition: abstract_factory.hpp:20
gko::matrix::SparsityCsr::get_value
value_type * get_value() noexcept
Returns the value stored in the matrix.
Definition: sparsity_csr.hpp:164
gko::array< index_type >
gko::matrix::SparsityCsr::read
void read(const mat_data &data) override
Reads a matrix from a matrix_data structure.
gko::dim< 2 >
gko::matrix::SparsityCsr::get_num_nonzeros
size_type get_num_nonzeros() const noexcept
Returns the number of elements explicitly stored in the matrix.
Definition: sparsity_csr.hpp:183
gko::matrix::SparsityCsr::sort_by_column_index
void sort_by_column_index()
Sorts each row by column index.
gko::array::get_data
value_type * get_data() noexcept
Returns a pointer to the block of memory used to store the elements of the array.
Definition: array.hpp:674
gko::matrix::SparsityCsr::create_const
static std::unique_ptr< const SparsityCsr > create_const(std::shared_ptr< const Executor > exec, const dim< 2 > &size, gko::detail::const_array_view< IndexType > &&col_idxs, gko::detail::const_array_view< IndexType > &&row_ptrs, ValueType value=one< ValueType >())
Creates a constant (immutable) SparsityCsr matrix from constant arrays.
Definition: sparsity_csr.hpp:263
gko::matrix::SparsityCsr::SparsityCsr
SparsityCsr(const SparsityCsr &)
Copy-constructs a SparsityCsr matrix.
gko::Executor
The first step in using the Ginkgo library consists of creating an executor.
Definition: executor.hpp:616
gko::array::get_const_data
const value_type * get_const_data() const noexcept
Returns a constant pointer to the block of memory used to store the elements of the array.
Definition: array.hpp:683
gko::matrix::SparsityCsr::get_row_ptrs
index_type * get_row_ptrs() noexcept
Returns the row pointers of the matrix.
Definition: sparsity_csr.hpp:145
gko::matrix::SparsityCsr::operator=
SparsityCsr & operator=(const SparsityCsr &)
Copy-assigns a SparsityCsr matrix.
gko::array::get_size
size_type get_size() const noexcept
Returns the number of elements in the array.
Definition: array.hpp:657
gko::matrix::SparsityCsr::conj_transpose
std::unique_ptr< LinOp > conj_transpose() const override
Returns a LinOp representing the conjugate transpose of the Transposable object.
gko::matrix::SparsityCsr::get_const_value
const value_type * get_const_value() const noexcept
Returns the value stored in the matrix.
Definition: sparsity_csr.hpp:173
gko::matrix::SparsityCsr::to_adjacency_matrix
std::unique_ptr< SparsityCsr > to_adjacency_matrix() const
Transforms the sparsity matrix to an adjacency matrix.
gko::matrix::SparsityCsr::write
void write(mat_data &data) const override
Writes a matrix to a matrix_data structure.
gko::LinOp::LinOp
LinOp(const LinOp &)=default
Copy-constructs a LinOp.
gko::one
constexpr T one()
Returns the multiplicative identity for T.
Definition: math.hpp:775