Basix
Loading...
Searching...
No Matches
precompute.h
1// Copyright (c) 2020 Matthew Scroggs
2// FEniCS Project
3// SPDX-License-Identifier: MIT
4
5#pragma once
6
7#include "math.h"
8#include "mdspan.hpp"
9#include <concepts>
10#include <cstdint>
11#include <span>
12#include <tuple>
13#include <type_traits>
14#include <vector>
15
18{
19namespace impl
20{
23template <typename T, typename = void>
24struct scalar_value_type
25{
27 typedef T value_type;
28};
30template <typename T>
31struct scalar_value_type<T, std::void_t<typename T::value_type>>
32{
33 typedef typename T::value_type value_type;
34};
36template <typename T>
37using scalar_value_type_t = typename scalar_value_type<T>::value_type;
38} // namespace impl
39
84void prepare_permutation(std::span<std::size_t> perm);
85
136template <typename E>
137void apply_permutation(std::span<const std::size_t> perm, std::span<E> data,
138 std::size_t offset = 0, std::size_t n = 1)
139{
140 for (std::size_t i = 0; i < perm.size(); ++i)
141 for (std::size_t b = 0; b < n; ++b)
142 std::swap(data[n * (offset + i) + b], data[n * (offset + perm[i]) + b]);
143}
144
146template <typename E>
147void apply_permutation_mapped(std::span<const std::size_t> perm,
148 std::span<E> data, std::span<const int> emap,
149 std::size_t n = 1)
150{
151 for (std::size_t i = 0; i < perm.size(); ++i)
152 for (std::size_t b = 0; b < n; ++b)
153 std::swap(data[n * emap[i] + b], data[n * emap[perm[i]] + b]);
154}
155
164template <typename E>
165void apply_inv_permutation_right(std::span<const std::size_t> perm,
166 std::span<E> data, std::size_t offset = 0,
167 std::size_t n = 1)
168{
169 const std::size_t dim = perm.size();
170 const std::size_t data_size = (data.size() + (dim < n ? n - dim : 0)) / n;
171 for (std::size_t b = 0; b < n; ++b)
172 {
173 for (std::size_t i = 0; i < dim; ++i)
174 {
175 std::swap(data[data_size * b + offset + i],
176 data[data_size * b + offset + perm[i]]);
177 }
178 }
179}
180
197template <std::floating_point T>
198std::vector<std::size_t>
199prepare_matrix(std::pair<std::vector<T>, std::array<std::size_t, 2>>& A)
200{
201 return math::transpose_lu<T>(A);
202}
203
235template <typename T, typename E>
236void apply_matrix(std::span<const std::size_t> v_size_t,
237 md::mdspan<const T, md::dextents<std::size_t, 2>> M,
238 std::span<E> data, std::size_t offset = 0, std::size_t n = 1)
239{
240 using U = typename impl::scalar_value_type_t<E>;
241
242 const std::size_t dim = v_size_t.size();
243 apply_permutation(v_size_t, data, offset, n);
244 for (std::size_t b = 0; b < n; ++b)
245 {
246 for (std::size_t i = 0; i < dim; ++i)
247 {
248 for (std::size_t j = i + 1; j < dim; ++j)
249 {
250 data[n * (offset + i) + b]
251 += static_cast<U>(M(i, j)) * data[n * (offset + j) + b];
252 }
253 }
254
255 for (std::size_t i = 1; i <= dim; ++i)
256 {
257 data[n * (offset + dim - i) + b] *= static_cast<U>(M(dim - i, dim - i));
258 for (std::size_t j = 0; j < dim - i; ++j)
259 {
260 data[n * (offset + dim - i) + b]
261 += static_cast<U>(M(dim - i, j)) * data[n * (offset + j) + b];
262 }
263 }
264 }
265}
266
275template <typename T, typename E>
277 std::span<const std::size_t> v_size_t,
278 md::mdspan<const T, md::dextents<std::size_t, 2>> M, std::span<E> data,
279 std::size_t offset = 0, std::size_t n = 1)
280{
281 using U = typename impl::scalar_value_type_t<E>;
282
283 const std::size_t dim = v_size_t.size();
284 const std::size_t data_size = (data.size() + (dim < n ? n - dim : 0)) / n;
285 apply_inv_permutation_right(v_size_t, data, offset, n);
286 for (std::size_t b = 0; b < n; ++b)
287 {
288 for (std::size_t i = 0; i < dim; ++i)
289 {
290 for (std::size_t j = i + 1; j < dim; ++j)
291 {
292 data[data_size * b + offset + i]
293 += static_cast<U>(M(i, j)) * data[data_size * b + offset + j];
294 }
295 }
296 for (std::size_t i = 1; i <= dim; ++i)
297 {
298 data[data_size * b + offset + dim - i]
299 *= static_cast<U>(M(dim - i, dim - i));
300 for (std::size_t j = 0; j < dim - i; ++j)
301 {
302 data[data_size * b + offset + dim - i]
303 += static_cast<U>(M(dim - i, j)) * data[data_size * b + offset + j];
304 }
305 }
306 }
307}
308
309} // namespace basix::precompute
std::vector< std::size_t > transpose_lu(std::pair< std::vector< T >, std::array< std::size_t, 2 > > &A)
Compute the LU decomposition of the transpose of a square matrix A.
Definition math.h:271
Matrix and permutation pre-computation.
Definition precompute.h:18
std::vector< std::size_t > prepare_matrix(std::pair< std::vector< T >, std::array< std::size_t, 2 > > &A)
Prepare a square matrix.
Definition precompute.h:199
void apply_permutation(std::span< const std::size_t > perm, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) permutation .
Definition precompute.h:137
void apply_inv_permutation_right(std::span< const std::size_t > perm, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) permutation to some transposed data.
Definition precompute.h:165
void prepare_permutation(std::span< std::size_t > perm)
Prepare a permutation.
Definition precompute.cpp:10
void apply_permutation_mapped(std::span< const std::size_t > perm, std::span< E > data, std::span< const int > emap, std::size_t n=1)
Permutation of mapped data.
Definition precompute.h:147
void apply_tranpose_matrix_right(std::span< const std::size_t > v_size_t, md::mdspan< const T, md::dextents< std::size_t, 2 > > M, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) matrix to some transposed data.
Definition precompute.h:276
void apply_matrix(std::span< const std::size_t > v_size_t, md::mdspan< const T, md::dextents< std::size_t, 2 > > M, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) matrix.
Definition precompute.h:236