Loading [MathJax]/extensions/tex2jax.js
DOLFINx 0.10.0.0
DOLFINx C++ interface
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Modules Pages Concepts
CoordinateElement.h
1// Copyright (C) 2018-2020 Garth N. Wells and Chris N. Richardson
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 "ElementDofLayout.h"
10#include <algorithm>
11#include <array>
12#include <basix/element-families.h>
13#include <basix/mdspan.hpp>
14#include <cmath>
15#include <concepts>
16#include <cstdint>
17#include <dolfinx/common/math.h>
18#include <dolfinx/mesh/cell_types.h>
19#include <limits>
20#include <memory>
21#include <span>
22
23namespace basix
24{
25template <std::floating_point T>
27}
28
29namespace dolfinx::fem
30{
36template <std::floating_point T>
38{
39public:
42 explicit CoordinateElement(
43 std::shared_ptr<const basix::FiniteElement<T>> element);
44
51 basix::element::lagrange_variant type
52 = basix::element::lagrange_variant::unset);
53
55 virtual ~CoordinateElement() = default;
56
60
63 int degree() const;
64
71 int dim() const;
72
74 basix::element::lagrange_variant variant() const;
75
79 std::uint64_t hash() const;
80
89 std::array<std::size_t, 4> tabulate_shape(std::size_t nd,
90 std::size_t num_points) const;
91
100 void tabulate(int nd, std::span<const T> X, std::array<std::size_t, 2> shape,
101 std::span<T> basis) const;
102
119 void permute_subentity_closure(std::span<std::int32_t> d,
120 std::uint32_t cell_info,
121 mesh::CellType entity_type,
122 int entity_index) const;
123
132 template <typename U, typename V, typename W>
133 static void compute_jacobian(const U& dphi, const V& cell_geometry, W&& J)
134 {
135 math::dot(cell_geometry, dphi, J, true);
136 }
137
141 template <typename U, typename V>
142 static void compute_jacobian_inverse(const U& J, V&& K)
143 {
144 int gdim = J.extent(0);
145 int tdim = K.extent(0);
146 if (gdim == tdim)
147 math::inv(J, K);
148 else
149 math::pinv(J, K);
150 }
151
157 template <typename U>
158 static double
159 compute_jacobian_determinant(const U& J, std::span<typename U::value_type> w)
160 {
161 static_assert(U::rank() == 2, "Must be rank 2");
162 if (J.extent(0) == J.extent(1))
163 return math::det(J);
164 else
165 {
166 assert(w.size() >= 2 * J.extent(0) * J.extent(1));
167 using X = typename U::element_type;
168 using mdspan2_t = md::mdspan<X, md::dextents<std::size_t, 2>>;
169 mdspan2_t B(w.data(), J.extent(1), J.extent(0));
170 mdspan2_t BA(w.data() + J.extent(0) * J.extent(1), B.extent(0),
171 J.extent(1));
172 for (std::size_t i = 0; i < B.extent(0); ++i)
173 for (std::size_t j = 0; j < B.extent(1); ++j)
174 B(i, j) = J(j, i);
175
176 // Zero working memory of BA
177 std::fill_n(BA.data_handle(), BA.size(), 0);
178 math::dot(B, J, BA);
179 return std::sqrt(math::det(BA));
180 }
181 }
182
185
193 template <typename U, typename V, typename W>
194 static void push_forward(U&& x, const V& cell_geometry, const W& phi)
195 {
196 for (std::size_t i = 0; i < x.extent(0); ++i)
197 for (std::size_t j = 0; j < x.extent(1); ++j)
198 x(i, j) = 0;
199
200 // Compute x = phi * cell_geometry;
201 math::dot(phi, cell_geometry, x);
202 }
203
214 template <typename U, typename V, typename W>
215 static void pull_back_affine(U&& X, const V& K, const std::array<T, 3>& x0,
216 const W& x)
217 {
218 assert(X.extent(0) == x.extent(0));
219 assert(X.extent(1) == K.extent(0));
220 assert(x.extent(1) == K.extent(1));
221 for (std::size_t i = 0; i < X.extent(0); ++i)
222 for (std::size_t j = 0; j < X.extent(1); ++j)
223 X(i, j) = 0;
224
225 // Calculate X for each point
226 for (std::size_t p = 0; p < x.extent(0); ++p)
227 for (std::size_t i = 0; i < K.extent(0); ++i)
228 for (std::size_t j = 0; j < K.extent(1); ++j)
229 X(p, i) += K(i, j) * (x(p, j) - x0[j]);
230 }
231
233 template <typename X>
234 using mdspan2_t = md::mdspan<X, md::dextents<std::size_t, 2>>;
235
248 mdspan2_t<const T> cell_geometry,
249 double tol = 1.0e-6, int maxit = 15) const;
250
252 void permute(std::span<std::int32_t> dofs, std::uint32_t cell_perm) const;
253
255 void permute_inv(std::span<std::int32_t> dofs, std::uint32_t cell_perm) const;
256
262 bool needs_dof_permutations() const;
263
266 bool is_affine() const noexcept { return _is_affine; }
267
268private:
269 // Flag denoting affine map
270 bool _is_affine;
271
272 // Basix Element
273 std::shared_ptr<const basix::FiniteElement<T>> _element;
274};
275} // namespace dolfinx::fem
Definition CoordinateElement.h:26
ElementDofLayout create_dof_layout() const
Compute and return the dof layout.
Definition CoordinateElement.cpp:75
void permute_subentity_closure(std::span< std::int32_t > d, std::uint32_t cell_info, mesh::CellType entity_type, int entity_index) const
Given the closure DOFs of a cell sub-entity in reference ordering, this function computes the permut...
Definition CoordinateElement.cpp:64
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:133
static void pull_back_affine(U &&X, const V &K, const std::array< T, 3 > &x0, const W &x)
Compute reference coordinates X for physical coordinates x for an affine map. For the affine case,...
Definition CoordinateElement.h:215
void tabulate(int nd, std::span< const T > X, std::array< std::size_t, 2 > shape, std::span< T > basis) const
Evaluate basis values and derivatives at set of points.
Definition CoordinateElement.cpp:55
CoordinateElement(std::shared_ptr< const basix::FiniteElement< T > > element)
Create a coordinate element from a Basix element.
Definition CoordinateElement.cpp:19
void permute_inv(std::span< std::int32_t > dofs, std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition CoordinateElement.cpp:182
basix::element::lagrange_variant variant() const
Variant of the element.
Definition CoordinateElement.cpp:212
mesh::CellType cell_shape() const
Cell shape.
Definition CoordinateElement.cpp:41
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:142
std::uint64_t hash() const
Element hash.
Definition CoordinateElement.cpp:219
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Shape of array to fill when calling tabulate.
Definition CoordinateElement.cpp:48
int degree() const
The polynomial degree of the element.
Definition CoordinateElement.cpp:198
void permute(std::span< std::int32_t > dofs, std::uint32_t cell_perm) const
Permute a list of DOF numbers on a cell.
Definition CoordinateElement.cpp:174
md::mdspan< X, md::dextents< std::size_t, 2 > > mdspan2_t
mdspan typedef
Definition CoordinateElement.h:234
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:205
void pull_back_nonaffine(mdspan2_t< T > X, mdspan2_t< const T > x, mdspan2_t< const T > cell_geometry, double tol=1.0e-6, int maxit=15) const
Compute reference coordinates X for physical coordinates x for a non-affine map.
Definition CoordinateElement.cpp:83
virtual ~CoordinateElement()=default
Destructor.
bool is_affine() const noexcept
Definition CoordinateElement.h:266
static void push_forward(U &&x, const V &cell_geometry, const W &phi)
Compute physical coordinates x for points X in the reference configuration.
Definition CoordinateElement.h:194
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:159
bool needs_dof_permutations() const
Indicates whether the geometry DOF numbers on each cell need permuting.
Definition CoordinateElement.cpp:190
Definition ElementDofLayout.h:30
Finite element method functionality.
Definition assemble_expression_impl.h:23
CellType
Cell type identifier.
Definition cell_types.h:22