9 #include "ElementDofLayout.h"
12 #include <basix/element-families.h>
13 #include <basix/mdspan.hpp>
15 #include <dolfinx/common/math.h>
16 #include <dolfinx/mesh/cell_types.h>
37 std::shared_ptr<const basix::FiniteElement> element);
45 basix::element::lagrange_variant type
46 = basix::element::lagrange_variant::equispaced);
67 basix::element::lagrange_variant
variant()
const;
76 std::size_t num_points)
const;
87 void tabulate(
int nd, std::span<const double> X,
88 std::array<std::size_t, 2> shape,
89 std::span<double> basis)
const;
99 template <
typename U,
typename V,
typename W>
102 math::dot(cell_geometry, dphi, J,
true);
108 template <
typename U,
typename V>
111 const int gdim = J.extent(1);
112 const int tdim = K.extent(1);
124 template <
typename U>
128 static_assert(U::rank() == 2,
"Must be rank 2");
129 if (J.extent(0) == J.extent(1))
133 assert(w.size() >= 2 * J.extent(0) * J.extent(1));
135 using T =
typename U::element_type;
136 namespace stdex = std::experimental;
137 using mdspan2_t = stdex::mdspan<T, stdex::dextents<std::size_t, 2>>;
138 mdspan2_t B(w.data(), J.extent(1), J.extent(0));
139 mdspan2_t BA(w.data() + J.extent(0) * J.extent(1), B.extent(0),
142 for (std::size_t i = 0; i < B.extent(0); ++i)
143 for (std::size_t j = 0; j < B.extent(1); ++j)
147 std::fill_n(BA.data_handle(), BA.size(), 0);
149 return std::sqrt(math::det(BA));
162 template <
typename U,
typename V,
typename W>
165 for (std::size_t i = 0; i < x.extent(0); ++i)
166 for (std::size_t j = 0; j < x.extent(1); ++j)
170 math::dot(phi, cell_geometry, x);
183 template <
typename U,
typename V,
typename W>
185 const std::array<double, 3>& x0,
const W& x)
187 assert(X.extent(0) == x.extent(0));
188 assert(X.extent(1) == K.extent(0));
189 assert(x.extent(1) == K.extent(1));
190 for (std::size_t i = 0; i < X.extent(0); ++i)
191 for (std::size_t j = 0; j < X.extent(1); ++j)
195 for (std::size_t p = 0; p < x.extent(0); ++p)
196 for (std::size_t i = 0; i < K.extent(0); ++i)
197 for (std::size_t j = 0; j < K.extent(1); ++j)
198 X(p, i) += K(i, j) * (x(p, j) - x0[j]);
203 = std::experimental::mdspan<double,
204 std::experimental::dextents<std::size_t, 2>>;
207 = std::experimental::mdspan<
const double,
208 std::experimental::dextents<std::size_t, 2>>;
222 double tol = 1.0e-8,
int maxit = 10)
const;
226 std::uint32_t cell_perm)
const;
230 std::uint32_t cell_perm)
const;
248 std::shared_ptr<const basix::FiniteElement> _element;
A CoordinateElement manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:32
ElementDofLayout create_dof_layout() const
Compute and return the dof layout.
Definition: CoordinateElement.cpp:55
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:100
static void pull_back_affine(U &&X, const V &K, const std::array< double, 3 > &x0, const W &x)
Compute reference coordinates X for physical coordinates x for an affine map. For the affine case,...
Definition: CoordinateElement.h:184
basix::element::lagrange_variant variant() const
The variant of the element.
Definition: CoordinateElement.cpp:189
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:36
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition: CoordinateElement.h:109
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 FiniteElement::tabulate
Definition: CoordinateElement.cpp:42
int degree() const
The polynomial degree of the element.
Definition: CoordinateElement.cpp:177
void pull_back_nonaffine(mdspan2_t X, cmdspan2_t x, cmdspan2_t cell_geometry, double tol=1.0e-8, int maxit=10) const
Compute reference coordinates X for physical coordinates x for a non-affine map.
Definition: CoordinateElement.cpp:62
std::experimental::mdspan< double, std::experimental::dextents< std::size_t, 2 > > mdspan2_t
mdspan typedef
Definition: CoordinateElement.h:204
void tabulate(int nd, std::span< const double > X, std::array< std::size_t, 2 > shape, std::span< double > basis) const
Evaluate basis values and derivatives at set of points.
Definition: CoordinateElement.cpp:47
int dim() const
The dimension of the geometry element space.
Definition: CoordinateElement.cpp:183
virtual ~CoordinateElement()=default
Destructor.
bool is_affine() const noexcept
Check is geometry map is affine.
Definition: CoordinateElement.h:241
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:163
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition: CoordinateElement.h:126
CoordinateElement(std::shared_ptr< const basix::FiniteElement > element)
Create a coordinate element from a Basix element.
Definition: CoordinateElement.cpp:17
void permute_dofs(const std::span< std::int32_t > &dofs, std::uint32_t cell_perm) const
Permutes a list of DOF numbers on a cell.
Definition: CoordinateElement.cpp:156
void unpermute_dofs(const std::span< std::int32_t > &dofs, std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition: CoordinateElement.cpp:163
bool needs_dof_permutations() const
Indicates whether the geometry DOF numbers on each cell need permuting.
Definition: CoordinateElement.cpp:170
std::experimental::mdspan< const double, std::experimental::dextents< std::size_t, 2 > > cmdspan2_t
mdspan typedef
Definition: CoordinateElement.h:208
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