9 #include "ElementDofLayout.h"
11 #include <dolfinx/common/array2d.h>
12 #include <dolfinx/common/span.hpp>
13 #include <dolfinx/mesh/cell_types.h>
80 tcb::span<double> detJ,
81 std::vector<double>& K,
86 void permute_dofs(
int* dofs,
const uint32_t cell_perm)
const;
100 std::string _signature;
109 int _basix_element_handle;
112 bool _needs_permutation_data;
115 std::function<int(
int*,
const uint32_t)> _permute_dofs;
118 std::function<int(
int*,
const uint32_t)> _unpermute_dofs;
This class provides a dynamic 2-dimensional row-wise array data structure.
Definition: array2d.h:21
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:25
virtual ~CoordinateElement()=default
Destructor.
void push_forward(array2d< double > &x, const array2d< double > &X, const array2d< double > &cell_geometry) const
Compute physical coordinates x for points X in the reference configuration.
Definition: CoordinateElement.cpp:67
int topological_dimension() const
Return the topological dimension of the cell shape.
Definition: CoordinateElement.cpp:54
double non_affine_atol
Absolute increment stopping criterium for non-affine Newton solver.
Definition: CoordinateElement.h:64
CoordinateElement(int basix_element_handle, int geometric_dimension, const std::string &signature, const ElementDofLayout &dof_layout, bool needs_permutation_data, std::function< int(int *, const uint32_t)> permute_dofs, std::function< int(int *, const uint32_t)> unpermute_dofs)
Create a coordinate element.
Definition: CoordinateElement.cpp:15
void compute_reference_geometry(array2d< double > &X, std::vector< double > &J, tcb::span< double > detJ, std::vector< double > &K, const array2d< double > &x, const array2d< double > &cell_geometry) const
Compute reference coordinates X, and J, detJ and K for physical coordinates x.
Definition: CoordinateElement.cpp:89
int non_affine_max_its
Maximum number of iterations for non-affine Newton solver.
Definition: CoordinateElement.h:67
std::string signature() const
String identifying the finite element.
Definition: CoordinateElement.cpp:34
int geometric_dimension() const
Return the geometric dimension of the cell shape.
Definition: CoordinateElement.cpp:60
bool needs_permutation_data() const
Indicates whether the coordinate map needs permutation data passing in (for higher order geometries)
Definition: CoordinateElement.cpp:247
void permute_dofs(int *dofs, const uint32_t cell_perm) const
Permutes a list of DOF numbers on a cell.
Definition: CoordinateElement.cpp:236
const ElementDofLayout & dof_layout() const
Return the dof layout.
Definition: CoordinateElement.cpp:62
void unpermute_dofs(int *dofs, const uint32_t cell_perm) const
Reverses a DOF permutation.
Definition: CoordinateElement.cpp:241
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:36
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh ent...
Definition: ElementDofLayout.h:36
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
CellType
Cell type identifier.
Definition: cell_types.h:22