9 #include "ElementDofLayout.h"
11 #include <dolfinx/mesh/cell_types.h>
14 #include <xtensor/xtensor.hpp>
15 #include <xtl/xspan.hpp>
57 xt::xtensor<double, 4>
tabulate(
int n,
const xt::xtensor<double, 2>& X)
const;
69 const xt::xtensor<double, 2>& cell_geometry,
70 xt::xtensor<double, 3>& J)
const;
80 xt::xtensor<double, 3>& K)
const;
89 xt::xtensor<double, 1>& detJ)
const;
100 const xt::xtensor<double, 2>& cell_geometry,
101 const xt::xtensor<double, 2>& phi);
105 void pull_back(xt::xtensor<double, 2>& X, xt::xtensor<double, 3>& J,
106 xt::xtensor<double, 1>& detJ, xt::xtensor<double, 3>& K,
107 const xt::xtensor<double, 2>& x,
108 const xt::xtensor<double, 2>& cell_geometry)
const;
112 const std::uint32_t cell_perm)
const;
116 const std::uint32_t cell_perm)
const;
136 std::shared_ptr<basix::FiniteElement> _element;
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:29
void unpermute_dofs(xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition: CoordinateElement.cpp:302
virtual ~CoordinateElement()=default
Destructor.
int topological_dimension() const
Return the topological dimension of the cell shape.
Definition: CoordinateElement.cpp:73
double non_affine_atol
Absolute increment stopping criterium for non-affine Newton solver.
Definition: CoordinateElement.h:126
bool needs_dof_permutations() const
Indicates whether the geometry DOF numbers on each cell need permuting.
Definition: CoordinateElement.cpp:309
void pull_back(xt::xtensor< double, 2 > &X, xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ, xt::xtensor< double, 3 > &K, const xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry) const
Compute reference coordinates X, and J, detJ and K for physical coordinates x.
Definition: CoordinateElement.cpp:208
void compute_jacobian(const xt::xtensor< double, 4 > &dphi, const xt::xtensor< double, 2 > &cell_geometry, xt::xtensor< double, 3 > &J) const
Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives...
Definition: CoordinateElement.cpp:84
void permute_dofs(xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const
Permutes a list of DOF numbers on a cell.
Definition: CoordinateElement.cpp:295
int non_affine_max_its
Maximum number of iterations for non-affine Newton solver.
Definition: CoordinateElement.h:129
ElementDofLayout dof_layout() const
Return the dof layout.
Definition: CoordinateElement.cpp:183
static void push_forward(xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry, const xt::xtensor< double, 2 > &phi)
Compute physical coordinates x for points X in the reference configuration.
Definition: CoordinateElement.cpp:193
void compute_jacobian_inverse(const xt::xtensor< double, 3 > &J, xt::xtensor< double, 3 > &K) const
Compute the inverse of the Jacobian. If the coordinate element is affine, it computes the inverse at ...
Definition: CoordinateElement.cpp:125
void compute_jacobian_determinant(const xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ) const
Compute the determinant of the Jacobian. If the coordinate element is affine, it computes the determi...
Definition: CoordinateElement.cpp:162
xt::xtensor< double, 4 > tabulate(int n, const xt::xtensor< double, 2 > &X) const
Compute basis values and derivatives at set of points.
Definition: CoordinateElement.cpp:79
CoordinateElement(std::shared_ptr< basix::FiniteElement > element)
Create a coordinate element from a Basix element.
Definition: CoordinateElement.cpp:36
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:55
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh ent...
Definition: ElementDofLayout.h:34
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
CellType
Cell type identifier.
Definition: cell_types.h:22