9 #include "ElementDofLayout.h"
10 #include <Eigen/Dense>
12 #include <dolfinx/mesh/cell_types.h>
15 #include <unsupported/Eigen/CXX11/Tensor>
37 const std::function<
int(
double*,
int,
int,
const double*)>&
38 evaluate_basis_derivatives);
67 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>
69 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
70 Eigen::Dynamic, Eigen::RowMajor>>& X,
71 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
72 Eigen::Dynamic, Eigen::RowMajor>>&
78 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>& X,
79 Eigen::Tensor<double, 3, Eigen::RowMajor>& J,
80 Eigen::Ref<Eigen::Array<double, Eigen::Dynamic, 1>> detJ,
81 Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
82 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
83 Eigen::Dynamic, Eigen::RowMajor>>& x,
84 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
85 Eigen::Dynamic, Eigen::RowMajor>>&
87 double eps = 1.0e-16)
const;
97 std::string _signature;
110 std::function<int(
double*,
int,
int,
const double*)>
111 _evaluate_basis_derivatives;
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:24
virtual ~CoordinateElement()=default
Destructor.
int topological_dimension() const
Return the topological dimension of the cell shape.
Definition: CoordinateElement.cpp:31
void compute_reference_geometry(Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &X, Eigen::Tensor< double, 3, Eigen::RowMajor > &J, Eigen::Ref< Eigen::Array< double, Eigen::Dynamic, 1 >> 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, double eps=1.0e-16) const
Compute reference coordinates X, and J, detJ and K for physical coordinates x.
Definition: CoordinateElement.cpp:61
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:40
std::string signature() const
String identifying the finite element.
Definition: CoordinateElement.cpp:27
int geometric_dimension() const
Return the geometric dimension of the cell shape.
Definition: CoordinateElement.cpp:33
const ElementDofLayout & dof_layout() const
Return the dof layout.
Definition: CoordinateElement.cpp:35
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:29
CoordinateElement(mesh::CellType cell_type, int topological_dimension, int geometric_dimension, const std::string &signature, const ElementDofLayout &dof_layout, bool is_affine, const std::function< int(double *, int, int, const double *)> &evaluate_basis_derivatives)
Create a coordinate element.
Definition: CoordinateElement.cpp:14
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:23
CellType
Cell type identifier.
Definition: cell_types.h:21