9 #include "ElementDofLayout.h"
12 #include <dolfinx/common/span.hpp>
13 #include <dolfinx/mesh/cell_types.h>
17 #include <unsupported/Eigen/CXX11/Tensor>
77 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>
79 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
80 Eigen::Dynamic, Eigen::RowMajor>>& X,
81 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
82 Eigen::Dynamic, Eigen::RowMajor>>&
88 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>& X,
89 Eigen::Tensor<double, 3, Eigen::RowMajor>& J, tcb::span<double> detJ,
90 Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
91 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
92 Eigen::Dynamic, Eigen::RowMajor>>& x,
93 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
94 Eigen::Dynamic, Eigen::RowMajor>>&
98 void permute_dofs(
int* dofs,
const uint32_t cell_perm)
const;
112 std::string _signature;
121 int _basix_element_handle;
124 bool _needs_permutation_data;
127 std::function<int(
int*,
const uint32_t)> _permute_dofs;
130 std::function<int(
int*,
const uint32_t)> _unpermute_dofs;
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:26
virtual ~CoordinateElement()=default
Destructor.
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:65
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 push_forward(Eigen::Ref< Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> x, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &cell_geometry) const
Compute physical coordinates x for points X in the reference configuration.
Definition: CoordinateElement.cpp:67
int non_affine_max_its
Maximum number of iterations for non-affine Newton solver.
Definition: CoordinateElement.h:68
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:242
void compute_reference_geometry(Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &X, Eigen::Tensor< double, 3, Eigen::RowMajor > &J, tcb::span< double > detJ, Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &x, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &cell_geometry) const
Compute reference coordinates X, and J, detJ and K for physical coordinates x.
Definition: CoordinateElement.cpp:85
void permute_dofs(int *dofs, const uint32_t cell_perm) const
Permutes a list of DOF numbers on a cell.
Definition: CoordinateElement.cpp:231
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:236
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:37
Finite element method functionality.
Definition: assemble_matrix_impl.h:24
CellType
Cell type identifier.
Definition: cell_types.h:21