Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/ddd/CoordinateElement_8h_source.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
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{
31
37template <std::floating_point T>
39{
40public:
43 explicit CoordinateElement(
44 std::shared_ptr<const basix::FiniteElement<T>> element);
45
52 basix::element::lagrange_variant type
53 = basix::element::lagrange_variant::unset);
54
56 virtual ~CoordinateElement() = default;
57
61
64 int degree() const;
65
72 int dim() const;
73
75 basix::element::lagrange_variant variant() const;
76
83 std::array<std::size_t, 4> tabulate_shape(std::size_t nd,
84 std::size_t num_points) const;
85
94 void tabulate(int nd, std::span<const T> X, std::array<std::size_t, 2> shape,
95 std::span<T> basis) const;
96
105 template <typename U, typename V, typename W>
106 static void compute_jacobian(const U& dphi, const V& cell_geometry, W&& J)
107 {
108 math::dot(cell_geometry, dphi, J, true);
109 }
110
114 template <typename U, typename V>
115 static void compute_jacobian_inverse(const U& J, V&& K)
116 {
117 int gdim = J.extent(0);
118 int tdim = K.extent(0);
119 if (gdim == tdim)
120 math::inv(J, K);
121 else
122 math::pinv(J, K);
123 }
124
130 template <typename U>
131 static double
132 compute_jacobian_determinant(const U& J, std::span<typename U::value_type> w)
133 {
134 static_assert(U::rank() == 2, "Must be rank 2");
135 if (J.extent(0) == J.extent(1))
136 return math::det(J);
137 else
138 {
139 assert(w.size() >= 2 * J.extent(0) * J.extent(1));
140 using X = typename U::element_type;
141 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
142 X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
143 mdspan2_t B(w.data(), J.extent(1), J.extent(0));
144 mdspan2_t BA(w.data() + J.extent(0) * J.extent(1), B.extent(0),
145 J.extent(1));
146 for (std::size_t i = 0; i < B.extent(0); ++i)
147 for (std::size_t j = 0; j < B.extent(1); ++j)
148 B(i, j) = J(j, i);
149
150 // Zero working memory of BA
151 std::fill_n(BA.data_handle(), BA.size(), 0);
152 math::dot(B, J, BA);
153 return std::sqrt(math::det(BA));
154 }
155 }
156
159
167 template <typename U, typename V, typename W>
168 static void push_forward(U&& x, const V& cell_geometry, const W& phi)
169 {
170 for (std::size_t i = 0; i < x.extent(0); ++i)
171 for (std::size_t j = 0; j < x.extent(1); ++j)
172 x(i, j) = 0;
173
174 // Compute x = phi * cell_geometry;
175 math::dot(phi, cell_geometry, x);
176 }
177
188 template <typename U, typename V, typename W>
189 static void pull_back_affine(U&& X, const V& K, const std::array<T, 3>& x0,
190 const W& x)
191 {
192 assert(X.extent(0) == x.extent(0));
193 assert(X.extent(1) == K.extent(0));
194 assert(x.extent(1) == K.extent(1));
195 for (std::size_t i = 0; i < X.extent(0); ++i)
196 for (std::size_t j = 0; j < X.extent(1); ++j)
197 X(i, j) = 0;
198
199 // Calculate X for each point
200 for (std::size_t p = 0; p < x.extent(0); ++p)
201 for (std::size_t i = 0; i < K.extent(0); ++i)
202 for (std::size_t j = 0; j < K.extent(1); ++j)
203 X(p, i) += K(i, j) * (x(p, j) - x0[j]);
204 }
205
207 template <typename X>
208 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
209 X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
210
223 mdspan2_t<const T> cell_geometry,
224 double tol = 1.0e-6, int maxit = 15) const;
225
227 void permute_dofs(const std::span<std::int32_t>& dofs,
228 std::uint32_t cell_perm) const;
229
231 void unpermute_dofs(const std::span<std::int32_t>& dofs,
232 std::uint32_t cell_perm) const;
233
239 bool needs_dof_permutations() const;
240
243 bool is_affine() const noexcept { return _is_affine; }
244
245private:
246 // Flag denoting affine map
247 bool _is_affine;
248
249 // Basix Element
250 std::shared_ptr<const basix::FiniteElement<T>> _element;
251};
252} // namespace dolfinx::fem
Definition CoordinateElement.h:26
A CoordinateElement manages coordinate mappings for isoparametric cells.
Definition CoordinateElement.h:39
ElementDofLayout create_dof_layout() const
Compute and return the dof layout.
Definition CoordinateElement.cpp:62
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives...
Definition CoordinateElement.h:106
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:189
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:53
basix::element::lagrange_variant variant() const
Variant of the element.
Definition CoordinateElement.cpp:200
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > mdspan2_t
mdspan typedef
Definition CoordinateElement.h:209
mesh::CellType cell_shape() const
Cell shape.
Definition CoordinateElement.cpp:39
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:115
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:46
int degree() const
The polynomial degree of the element.
Definition CoordinateElement.cpp:186
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:193
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:70
virtual ~CoordinateElement()=default
Destructor.
bool is_affine() const noexcept
Check is geometry map is affine.
Definition CoordinateElement.h:243
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:168
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:132
void permute_dofs(const std::span< std::int32_t > &dofs, std::uint32_t cell_perm) const
Permute a list of DOF numbers on a cell.
Definition CoordinateElement.cpp:162
void unpermute_dofs(const std::span< std::int32_t > &dofs, std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition CoordinateElement.cpp:170
bool needs_dof_permutations() const
Indicates whether the geometry DOF numbers on each cell need permuting.
Definition CoordinateElement.cpp:178
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh ent...
Definition ElementDofLayout.h:31
Finite element method functionality.
Definition assemble_matrix_impl.h:25
CellType
Cell type identifier.
Definition cell_types.h:22