DOLFINx 0.8.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{
20
26template <typename T>
28{
29public:
33 explicit AdjacencyList(const std::int32_t n) : _array(n), _offsets(n + 1)
34 {
35 std::iota(_array.begin(), _array.end(), 0);
36 std::iota(_offsets.begin(), _offsets.end(), 0);
37 }
38
43 template <typename U, typename V>
44 requires std::is_convertible_v<std::remove_cvref_t<U>, std::vector<T>>
45 and std::is_convertible_v<std::remove_cvref_t<V>,
46 std::vector<std::int32_t>>
47 AdjacencyList(U&& data, V&& offsets)
48 : _array(std::forward<U>(data)), _offsets(std::forward<V>(offsets))
49 {
50 _array.reserve(_offsets.back());
51 assert(_offsets.back() == (std::int32_t)_array.size());
52 }
53
59 template <typename X>
60 explicit AdjacencyList(const std::vector<X>& data)
61 {
62 // Initialize offsets and compute total size
63 _offsets.reserve(data.size() + 1);
64 _offsets.push_back(0);
65 for (auto& row : data)
66 _offsets.push_back(_offsets.back() + row.size());
67
68 _array.reserve(_offsets.back());
69 for (auto& e : data)
70 _array.insert(_array.end(), e.begin(), e.end());
71 }
72
74 AdjacencyList(const AdjacencyList& list) = default;
75
77 AdjacencyList(AdjacencyList&& list) = default;
78
80 ~AdjacencyList() = default;
81
83 AdjacencyList& operator=(const AdjacencyList& list) = default;
84
87
90 bool operator==(const AdjacencyList& list) const
91 {
92 return this->_array == list._array and this->_offsets == list._offsets;
93 }
94
97 std::int32_t num_nodes() const { return _offsets.size() - 1; }
98
102 int num_links(std::size_t node) const
103 {
104 assert((node + 1) < _offsets.size());
105 return _offsets[node + 1] - _offsets[node];
106 }
107
112 std::span<T> links(std::size_t node)
113 {
114 return std::span<T>(_array.data() + _offsets[node],
115 _offsets[node + 1] - _offsets[node]);
116 }
117
122 std::span<const T> links(std::size_t node) const
123 {
124 return std::span<const T>(_array.data() + _offsets[node],
125 _offsets[node + 1] - _offsets[node]);
126 }
127
129 const std::vector<T>& array() const { return _array; }
130
132 std::vector<T>& array() { return _array; }
133
135 const std::vector<std::int32_t>& offsets() const { return _offsets; }
136
138 std::vector<std::int32_t>& offsets() { return _offsets; }
139
142 std::string str() const
143 {
144 std::stringstream s;
145 s << "<AdjacencyList> with " + std::to_string(this->num_nodes()) + " nodes"
146 << std::endl;
147 for (std::size_t e = 0; e < _offsets.size() - 1; ++e)
148 {
149 s << " " << e << ": [";
150 for (auto link : this->links(e))
151 s << link << " ";
152 s << "]" << std::endl;
153 }
154 return s.str();
155 }
156
157private:
158 // Connections for all entities stored as a contiguous array
159 std::vector<T> _array;
160
161 // Position of first connection for each entity (using local index)
162 std::vector<std::int32_t> _offsets;
163};
164
166template <typename T, typename U>
167AdjacencyList(T, U) -> AdjacencyList<typename T::value_type>;
168
176template <typename U>
177 requires requires {
178 typename std::decay_t<U>::value_type;
179 requires std::convertible_to<
180 U, std::vector<typename std::decay_t<U>::value_type>>;
181 }
182AdjacencyList<typename std::decay_t<U>::value_type>
183regular_adjacency_list(U&& data, int degree)
184{
185 if (degree == 0 and !data.empty())
186 {
187 throw std::runtime_error("Degree is zero but data is not empty for "
188 "constant degree AdjacencyList");
189 }
190
191 if (degree > 0 and data.size() % degree != 0)
192 {
193 throw std::runtime_error(
194 "Incompatible data size and degree for constant degree AdjacencyList");
195 }
196
197 std::int32_t num_nodes = degree == 0 ? data.size() : data.size() / degree;
198 std::vector<std::int32_t> offsets(num_nodes + 1, 0);
199 for (std::size_t i = 1; i < offsets.size(); ++i)
200 offsets[i] = offsets[i - 1] + degree;
202 std::forward<U>(data), std::move(offsets));
203}
204
205} // namespace dolfinx::graph
Definition AdjacencyList.h:28
int num_links(std::size_t node) const
Definition AdjacencyList.h:102
std::span< T > links(std::size_t node)
Definition AdjacencyList.h:112
bool operator==(const AdjacencyList &list) const
Definition AdjacencyList.h:90
std::vector< std::int32_t > & offsets()
Offset for each node in array()
Definition AdjacencyList.h:138
const std::vector< T > & array() const
Return contiguous array of links for all nodes (const version)
Definition AdjacencyList.h:129
AdjacencyList & operator=(AdjacencyList &&list)=default
Move assignment operator.
std::span< const T > links(std::size_t node) const
Definition AdjacencyList.h:122
AdjacencyList(const std::int32_t n)
Definition AdjacencyList.h:33
AdjacencyList(AdjacencyList &&list)=default
Move constructor.
AdjacencyList & operator=(const AdjacencyList &list)=default
Assignment operator.
AdjacencyList(const std::vector< X > &data)
Definition AdjacencyList.h:60
std::int32_t num_nodes() const
Definition AdjacencyList.h:97
AdjacencyList(const AdjacencyList &list)=default
Copy constructor.
std::vector< T > & array()
Return contiguous array of links for all nodes.
Definition AdjacencyList.h:132
AdjacencyList(U &&data, V &&offsets)
Definition AdjacencyList.h:47
const std::vector< std::int32_t > & offsets() const
Offset for each node in array() (const version)
Definition AdjacencyList.h:135
std::string str() const
Definition AdjacencyList.h:142
~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:183