DOLFINx 0.9.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
AdjacencyList.h
1// Copyright (C) 2019-2022 Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include <cassert>
10#include <concepts>
11#include <cstdint>
12#include <numeric>
13#include <span>
14#include <sstream>
15#include <utility>
16#include <vector>
17
18namespace dolfinx::graph
19{
25template <typename T>
26class AdjacencyList
27{
28public:
32 explicit AdjacencyList(const std::int32_t n) : _array(n), _offsets(n + 1)
33 {
34 std::iota(_array.begin(), _array.end(), 0);
35 std::iota(_offsets.begin(), _offsets.end(), 0);
36 }
37
42 template <typename U, typename V>
43 requires std::is_convertible_v<std::remove_cvref_t<U>, std::vector<T>>
44 and std::is_convertible_v<std::remove_cvref_t<V>,
45 std::vector<std::int32_t>>
46 AdjacencyList(U&& data, V&& offsets)
47 : _array(std::forward<U>(data)), _offsets(std::forward<V>(offsets))
48 {
49 _array.reserve(_offsets.back());
50 assert(_offsets.back() == (std::int32_t)_array.size());
51 }
52
58 template <typename X>
59 explicit AdjacencyList(const std::vector<X>& data)
60 {
61 // Initialize offsets and compute total size
62 _offsets.reserve(data.size() + 1);
63 _offsets.push_back(0);
64 for (auto& row : data)
65 _offsets.push_back(_offsets.back() + row.size());
66
67 _array.reserve(_offsets.back());
68 for (auto& e : data)
69 _array.insert(_array.end(), e.begin(), e.end());
70 }
71
73 AdjacencyList(const AdjacencyList& list) = default;
74
76 AdjacencyList(AdjacencyList&& list) = default;
77
79 ~AdjacencyList() = default;
80
82 AdjacencyList& operator=(const AdjacencyList& list) = default;
83
86
89 bool operator==(const AdjacencyList& list) const
90 {
91 return this->_array == list._array and this->_offsets == list._offsets;
92 }
93
96 std::int32_t num_nodes() const { return _offsets.size() - 1; }
97
101 int num_links(std::size_t node) const
102 {
103 assert((node + 1) < _offsets.size());
104 return _offsets[node + 1] - _offsets[node];
105 }
106
111 std::span<T> links(std::size_t node)
112 {
113 return std::span<T>(_array.data() + _offsets[node],
114 _offsets[node + 1] - _offsets[node]);
115 }
116
121 std::span<const T> links(std::size_t node) const
122 {
123 return std::span<const T>(_array.data() + _offsets[node],
124 _offsets[node + 1] - _offsets[node]);
125 }
126
128 const std::vector<T>& array() const { return _array; }
129
131 std::vector<T>& array() { return _array; }
132
134 const std::vector<std::int32_t>& offsets() const { return _offsets; }
135
137 std::vector<std::int32_t>& offsets() { return _offsets; }
138
141 std::string str() const
142 {
143 std::stringstream s;
144 s << "<AdjacencyList> with " + std::to_string(this->num_nodes()) + " nodes"
145 << std::endl;
146 for (std::size_t e = 0; e < _offsets.size() - 1; ++e)
147 {
148 s << " " << e << ": [";
149 for (auto link : this->links(e))
150 s << link << " ";
151 s << "]" << std::endl;
152 }
153 return s.str();
154 }
155
156private:
157 // Connections for all entities stored as a contiguous array
158 std::vector<T> _array;
159
160 // Position of first connection for each entity (using local index)
161 std::vector<std::int32_t> _offsets;
162};
163
165template <typename T, typename U>
166AdjacencyList(T, U) -> AdjacencyList<typename T::value_type>;
167
175template <typename U>
176 requires requires {
177 typename std::decay_t<U>::value_type;
178 requires std::convertible_to<
179 U, std::vector<typename std::decay_t<U>::value_type>>;
180 }
181AdjacencyList<typename std::decay_t<U>::value_type>
182regular_adjacency_list(U&& data, int degree)
183{
184 if (degree == 0 and !data.empty())
185 {
186 throw std::runtime_error("Degree is zero but data is not empty for "
187 "constant degree AdjacencyList");
188 }
189
190 if (degree > 0 and data.size() % degree != 0)
191 {
192 throw std::runtime_error(
193 "Incompatible data size and degree for constant degree AdjacencyList");
194 }
195
196 std::int32_t num_nodes = degree == 0 ? data.size() : data.size() / degree;
197 std::vector<std::int32_t> offsets(num_nodes + 1, 0);
198 for (std::size_t i = 1; i < offsets.size(); ++i)
199 offsets[i] = offsets[i - 1] + degree;
201 std::forward<U>(data), std::move(offsets));
202}
203
204} // namespace dolfinx::graph
Definition topologycomputation.h:24
int num_links(std::size_t node) const
Definition AdjacencyList.h:101
std::span< T > links(std::size_t node)
Definition AdjacencyList.h:111
bool operator==(const AdjacencyList &list) const
Definition AdjacencyList.h:89
std::vector< std::int32_t > & offsets()
Offset for each node in array()
Definition AdjacencyList.h:137
const std::vector< T > & array() const
Return contiguous array of links for all nodes (const version)
Definition AdjacencyList.h:128
AdjacencyList & operator=(AdjacencyList &&list)=default
Move assignment operator.
std::span< const T > links(std::size_t node) const
Definition AdjacencyList.h:121
AdjacencyList(const std::int32_t n)
Definition AdjacencyList.h:32
AdjacencyList(AdjacencyList &&list)=default
Move constructor.
AdjacencyList & operator=(const AdjacencyList &list)=default
Assignment operator.
AdjacencyList(const std::vector< X > &data)
Definition AdjacencyList.h:59
std::int32_t num_nodes() const
Definition AdjacencyList.h:96
AdjacencyList(const AdjacencyList &list)=default
Copy constructor.
std::vector< T > & array()
Return contiguous array of links for all nodes.
Definition AdjacencyList.h:131
AdjacencyList(U &&data, V &&offsets)
Definition AdjacencyList.h:46
const std::vector< std::int32_t > & offsets() const
Offset for each node in array() (const version)
Definition AdjacencyList.h:134
std::string str() const
Definition AdjacencyList.h:141
~AdjacencyList()=default
Destructor.
Graph data structures and algorithms.
Definition dofmapbuilder.h:26
AdjacencyList< typename std::decay_t< U >::value_type > regular_adjacency_list(U &&data, int degree)
Construct a constant degree (valency) adjacency list.
Definition AdjacencyList.h:182