DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
sort.h
1// Copyright (C) 2021 Igor Baratta
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 <algorithm>
10#include <cassert>
11#include <concepts>
12#include <cstdint>
13#include <functional>
14#include <iterator>
15#include <numeric>
16#include <span>
17#include <type_traits>
18#include <utility>
19#include <vector>
20
21namespace dolfinx
22{
23struct __radix_sort
24{
49 template <
50 std::ranges::random_access_range R, typename P = std::identity,
51 std::remove_cvref_t<std::invoke_result_t<P, std::iter_value_t<R>>> BITS
52 = 8>
53 requires std::integral<decltype(BITS)>
54 constexpr void operator()(R&& range, P proj = {}) const
55 {
56 // value type
57 using T = std::iter_value_t<R>;
58
59 // index type (if no projection is provided it holds I == T)
60 using I = std::remove_cvref_t<std::invoke_result_t<P, T>>;
61
62 if (range.size() <= 1)
63 return;
64
65 T max_value = proj(*std::ranges::max_element(range, std::less{}, proj));
66
67 // Sort N bits at a time
68 constexpr I bucket_size = 1 << BITS;
69 T mask = (T(1) << BITS) - 1;
70
71 // Compute number of iterations, most significant digit (N bits) of
72 // maxvalue
73 I its = 0;
74 while (max_value)
75 {
76 max_value >>= BITS;
77 its++;
78 }
79
80 // Adjacency list arrays for computing insertion position
81 std::array<I, bucket_size> counter;
82 std::array<I, bucket_size + 1> offset;
83
84 I mask_offset = 0;
85 std::vector<T> buffer(range.size());
86 std::span<T> current_perm = range;
87 std::span<T> next_perm = buffer;
88 for (I i = 0; i < its; i++)
89 {
90 // Zero counter array
91 std::ranges::fill(counter, 0);
92
93 // Count number of elements per bucket
94 for (const auto& c : current_perm)
95 counter[(proj(c) & mask) >> mask_offset]++;
96
97 // Prefix sum to get the inserting position
98 offset[0] = 0;
99 std::partial_sum(counter.begin(), counter.end(),
100 std::next(offset.begin()));
101 for (const auto& c : current_perm)
102 {
103 I bucket = (proj(c) & mask) >> mask_offset;
104 I new_pos = offset[bucket + 1] - counter[bucket];
105 next_perm[new_pos] = c;
106 counter[bucket]--;
107 }
108
109 mask = mask << BITS;
110 mask_offset += BITS;
111
112 std::swap(current_perm, next_perm);
113 }
114
115 // Copy data back to array
116 if (its % 2 != 0)
117 std::ranges::copy(buffer, range.begin());
118 }
119};
120
122inline constexpr __radix_sort radix_sort{};
123
134template <typename T, int BITS = 16>
135std::vector<std::int32_t> sort_by_perm(std::span<const T> x, std::size_t shape1)
136{
137 static_assert(std::is_integral_v<T>, "Integral required.");
138
139 if (x.empty())
140 return std::vector<std::int32_t>{};
141
142 assert(shape1 > 0);
143 assert(x.size() % shape1 == 0);
144 const std::size_t shape0 = x.size() / shape1;
145 std::vector<std::int32_t> perm(shape0);
146 std::iota(perm.begin(), perm.end(), 0);
147
148 // Sort by each column, right to left. Col 0 has the most significant
149 // "digit".
150 std::vector<T> column(shape0);
151 for (std::size_t i = 0; i < shape1; ++i)
152 {
153 std::size_t col = shape1 - 1 - i;
154 for (std::size_t j = 0; j < shape0; ++j)
155 column[j] = x[j * shape1 + col];
156
157 radix_sort(perm, [&column](auto index) { return column[index]; });
158 }
159
160 return perm;
161}
162
163} // namespace dolfinx
Top-level namespace.
Definition defines.h:12
std::vector< std::int32_t > sort_by_perm(std::span< const T > x, std::size_t shape1)
Compute the permutation array that sorts a 2D array by row.
Definition sort.h:135
constexpr __radix_sort radix_sort
Radix sort.
Definition sort.h:122