DOLFIN-X
DOLFIN-X C++ interface
|
9 #include "ElementDofLayout.h"
10 #include <Eigen/Dense>
12 #include <dolfinx/mesh/cell_types.h>
15 #include <unsupported/Eigen/CXX11/Tensor>
40 std::function<
void(
double*,
int,
const double*,
const double*)>
41 compute_physical_coordinates,
42 std::function<
void(
double*,
double*,
double*,
double*,
int,
const double*,
73 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>
75 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
76 Eigen::Dynamic, Eigen::RowMajor>>& X,
77 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
78 Eigen::Dynamic, Eigen::RowMajor>>&
84 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>& X,
85 Eigen::Tensor<double, 3, Eigen::RowMajor>& J,
86 Eigen::Ref<Eigen::Array<double, Eigen::Dynamic, 1>> detJ,
87 Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
88 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
89 Eigen::Dynamic, Eigen::RowMajor>>& x,
90 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
91 Eigen::Dynamic, Eigen::RowMajor>>&
99 std::string _signature;
103 std::function<void(
double*,
int,
const double*,
const double*)>
104 _compute_physical_coordinates;
106 std::function<void(
double*,
double*,
double*,
double*,
int,
const double*,
108 _compute_reference_geometry;
CellType
Cell type identifier.
Definition: cell_types.h:22
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:43
virtual ~CoordinateElement()=default
Destructor.
const ElementDofLayout & dof_layout() const
Return the dof layout.
Definition: CoordinateElement.cpp:38
int geometric_dimension() const
Return the geometric dimension of the cell shape.
Definition: CoordinateElement.cpp:36
CoordinateElement(mesh::CellType cell_type, int topological_dimension, int geometric_dimension, const std::string &signature, const ElementDofLayout &dof_layout, std::function< void(double *, int, const double *, const double *)> compute_physical_coordinates, std::function< void(double *, double *, double *, double *, int, const double *, const double *)> compute_reference_geometry)
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:36
Finite element method functionality.
Definition: assemble_matrix_impl.h:34
std::string signature() const
String identifying the finite element.
Definition: CoordinateElement.cpp:30
int topological_dimension() const
Return the topological dimension of the cell shape.
Definition: CoordinateElement.cpp:34
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) const
Compute reference coordinates X, and J, detJ and K for physical coordinates x.
Definition: CoordinateElement.cpp:60
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:32
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:23