Basix 0.9.0.0

Home     Installation     Demos     C++ docs     Python docs

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 <span>
11 #include <tuple>
12 #include <type_traits>
13 #include <vector>
14 
17 {
18 namespace impl
19 {
22 template <typename T, typename = void>
23 struct scalar_value_type
24 {
26  typedef T value_type;
27 };
29 template <typename T>
30 struct scalar_value_type<T, std::void_t<typename T::value_type>>
31 {
32  typedef typename T::value_type value_type;
33 };
35 template <typename T>
36 using scalar_value_type_t = typename scalar_value_type<T>::value_type;
37 } // namespace impl
38 
88 void prepare_permutation(std::span<std::size_t> perm);
89 
143 template <typename E>
144 void apply_permutation(std::span<const std::size_t> perm, std::span<E> data,
145  std::size_t offset = 0, std::size_t n = 1)
146 {
147  for (std::size_t i = 0; i < perm.size(); ++i)
148  for (std::size_t b = 0; b < n; ++b)
149  std::swap(data[n * (offset + i) + b], data[n * (offset + perm[i]) + b]);
150 }
151 
153 template <typename E>
154 void apply_permutation_mapped(std::span<const std::size_t> perm,
155  std::span<E> data, std::span<const int> emap,
156  std::size_t n = 1)
157 {
158  for (std::size_t i = 0; i < perm.size(); ++i)
159  for (std::size_t b = 0; b < n; ++b)
160  std::swap(data[n * emap[i] + b], data[n * emap[perm[i]] + b]);
161 }
162 
171 template <typename E>
172 void apply_inv_permutation_right(std::span<const std::size_t> perm,
173  std::span<E> data, std::size_t offset = 0,
174  std::size_t n = 1)
175 {
176  const std::size_t dim = perm.size();
177  const std::size_t data_size = (data.size() + (dim < n ? n - dim : 0)) / n;
178  for (std::size_t b = 0; b < n; ++b)
179  {
180  for (std::size_t i = 0; i < dim; ++i)
181  {
182  std::swap(data[data_size * b + offset + i],
183  data[data_size * b + offset + perm[i]]);
184  }
185  }
186 }
187 
204 template <std::floating_point T>
205 std::vector<std::size_t>
206 prepare_matrix(std::pair<std::vector<T>, std::array<std::size_t, 2>>& A)
207 {
208  return math::transpose_lu<T>(A);
209 }
210 
242 template <typename T, typename E>
244  std::span<const std::size_t> v_size_t,
245  MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
246  const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
247  M,
248  std::span<E> data, std::size_t offset = 0, std::size_t n = 1)
249 {
250  using U = typename impl::scalar_value_type_t<E>;
251 
252  const std::size_t dim = v_size_t.size();
253  apply_permutation(v_size_t, data, offset, n);
254  for (std::size_t b = 0; b < n; ++b)
255  {
256  for (std::size_t i = 0; i < dim; ++i)
257  {
258  for (std::size_t j = i + 1; j < dim; ++j)
259  {
260  data[n * (offset + i) + b]
261  += static_cast<U>(M(i, j)) * data[n * (offset + j) + b];
262  }
263  }
264  for (std::size_t i = 1; i <= dim; ++i)
265  {
266  data[n * (offset + dim - i) + b] *= static_cast<U>(M(dim - i, dim - i));
267  for (std::size_t j = 0; j < dim - i; ++j)
268  {
269  data[n * (offset + dim - i) + b]
270  += static_cast<U>(M(dim - i, j)) * data[n * (offset + j) + b];
271  }
272  }
273  }
274 }
275 
284 template <typename T, typename E>
286  std::span<const std::size_t> v_size_t,
287  MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
288  const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
289  M,
290  std::span<E> data, std::size_t offset = 0, std::size_t n = 1)
291 {
292  using U = typename impl::scalar_value_type_t<E>;
293 
294  const std::size_t dim = v_size_t.size();
295  const std::size_t data_size = (data.size() + (dim < n ? n - dim : 0)) / n;
296  apply_inv_permutation_right(v_size_t, data, offset, n);
297  for (std::size_t b = 0; b < n; ++b)
298  {
299  for (std::size_t i = 0; i < dim; ++i)
300  {
301  for (std::size_t j = i + 1; j < dim; ++j)
302  {
303  data[data_size * b + offset + i]
304  += static_cast<U>(M(i, j)) * data[data_size * b + offset + j];
305  }
306  }
307  for (std::size_t i = 1; i <= dim; ++i)
308  {
309  data[data_size * b + offset + dim - i]
310  *= static_cast<U>(M(dim - i, dim - i));
311  for (std::size_t j = 0; j < dim - i; ++j)
312  {
313  data[data_size * b + offset + dim - i]
314  += static_cast<U>(M(dim - i, j)) * data[data_size * b + offset + j];
315  }
316  }
317  }
318 }
319 
320 } // namespace basix::precompute
Matrix and permutation pre-computation.
Definition: precompute.h:17
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:206
void apply_tranpose_matrix_right(std::span< const std::size_t > v_size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const T, MDSPAN_IMPL_STANDARD_NAMESPACE::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:285
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:144
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:172
void prepare_permutation(std::span< std::size_t > perm)
Prepare a permutation.
Definition: precompute.cpp:10
void apply_matrix(std::span< const std::size_t > v_size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const T, MDSPAN_IMPL_STANDARD_NAMESPACE::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:243
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:154