DOLFIN-X
DOLFIN-X C++ interface
|
9 #include <dolfinx/common/types.h>
10 #include <dolfinx/mesh/cell_types.h>
14 #include <unsupported/Eigen/CXX11/Tensor>
17 struct ufc_coordinate_mapping;
18 struct ufc_finite_element;
70 std::string
family()
const;
75 Eigen::Tensor<double, 3, Eigen::RowMajor>& reference_values,
76 const Eigen::Ref<
const Eigen::Array<
77 double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>& X)
const;
81 Eigen::Tensor<double, 3, Eigen::RowMajor>& values,
82 const Eigen::Tensor<double, 3, Eigen::RowMajor>& reference_values,
83 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
84 Eigen::Dynamic, Eigen::RowMajor>>& X,
85 const Eigen::Tensor<double, 3, Eigen::RowMajor>& J,
86 const Eigen::Ref<
const Eigen::Array<double, Eigen::Dynamic, 1>>& detJ,
87 const Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
88 const std::uint32_t permutation_info)
const;
92 Eigen::Tensor<double, 4, Eigen::RowMajor>& values, std::size_t order,
93 const Eigen::Tensor<double, 4, Eigen::RowMajor>& reference_values,
94 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
95 Eigen::Dynamic, Eigen::RowMajor>>& X,
96 const Eigen::Tensor<double, 3, Eigen::RowMajor>& J,
97 const Eigen::Ref<
const Eigen::Array<double, Eigen::Dynamic, 1>>& detJ,
98 const Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
99 const std::uint32_t permutation_info)
const;
107 const Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>&
113 PetscScalar* reference_values,
114 const Eigen::Ref<
const Eigen::Array<PetscScalar, Eigen::Dynamic,
115 Eigen::Dynamic, Eigen::RowMajor>>&
117 const Eigen::Ref<
const Eigen::Array<
double, Eigen::Dynamic,
118 Eigen::Dynamic, Eigen::RowMajor>>&
119 coordinate_dofs)
const;
125 std::size_t
hash()
const;
128 std::shared_ptr<const FiniteElement>
132 std::string _signature, _family;
136 int _tdim, _space_dim, _value_size, _reference_value_size, _degree;
140 Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> _refX;
143 std::vector<std::shared_ptr<const FiniteElement>> _sub_elements;
146 static std::shared_ptr<const FiniteElement>
148 const std::vector<int>& component);
154 std::vector<int> _value_dimension;
157 std::function<int(
double*,
int,
const double*)> _evaluate_reference_basis;
159 std::function<int(
double*,
int,
int,
const double*)>
160 _evaluate_reference_basis_derivatives;
162 std::function<int(
double*,
int,
int,
const double*,
const double*,
163 const double*,
const double*,
const double*,
164 const std::uint32_t)>
165 _transform_reference_basis_derivatives;
167 std::function<int(ufc_scalar_t*,
const ufc_scalar_t*,
const double*,
168 const ufc_coordinate_mapping*)>
FiniteElement(const ufc_finite_element &ufc_element)
Create finite element from UFC finite element.
Definition: FiniteElement.cpp:19
std::shared_ptr< const FiniteElement > extract_sub_element(const std::vector< int > &component) const
Extract sub finite element for component.
Definition: FiniteElement.cpp:207
int value_dimension(int i) const
Return the dimension of the value space for axis i.
Definition: FiniteElement.cpp:99
void evaluate_reference_basis(Eigen::Tensor< double, 3, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X) const
Evaluate all basis functions at given point in reference cell.
Definition: FiniteElement.cpp:110
int num_sub_elements() const
Return the number of sub elements (for a mixed element)
Definition: FiniteElement.cpp:202
const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > & dof_reference_coordinates() const
Tabulate the reference coordinates of all dofs on an element.
Definition: FiniteElement.cpp:176
CellType
Cell type identifier.
Definition: cell_types.h:22
mesh::CellType cell_shape() const
Cell shape.
Definition: FiniteElement.cpp:84
void transform_values(PetscScalar *reference_values, const Eigen::Ref< const Eigen::Array< PetscScalar, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &physical_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &coordinate_dofs) const
Map values of field from physical to reference space which has been evaluated at points given by dof_...
Definition: FiniteElement.cpp:187
Finite Element, containing the dof layout on a reference element, and various methods for evaluating ...
Definition: FiniteElement.h:25
void transform_reference_basis(Eigen::Tensor< double, 3, Eigen::RowMajor > &values, const Eigen::Tensor< double, 3, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Tensor< double, 3, Eigen::RowMajor > &J, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, 1 >> &detJ, const Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const std::uint32_t permutation_info) const
Push basis functions forward to physical element.
Definition: FiniteElement.cpp:126
int value_rank() const
Rank of the value space.
Definition: FiniteElement.cpp:97
std::size_t hash() const
Return simple hash of the signature string.
Definition: FiniteElement.cpp:204
virtual ~FiniteElement()=default
Destructor.
void transform_reference_basis_derivatives(Eigen::Tensor< double, 4, Eigen::RowMajor > &values, std::size_t order, const Eigen::Tensor< double, 4, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Tensor< double, 3, Eigen::RowMajor > &J, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, 1 >> &detJ, const Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const std::uint32_t permutation_info) const
Push basis function (derivatives) forward to physical element.
Definition: FiniteElement.cpp:148
int space_dimension() const
Dimension of the finite element function space.
Definition: FiniteElement.cpp:88
std::string signature() const
String identifying the finite element.
Definition: FiniteElement.cpp:82
std::string family() const
The finite element family.
Definition: FiniteElement.cpp:108
Finite element method functionality.
Definition: assemble_matrix_impl.h:34
int value_size() const
The value size, e.g. 1 for a scalar function, 2 for a 2D vector.
Definition: FiniteElement.cpp:90
int degree() const
Return the maximum polynomial degree.
Definition: FiniteElement.cpp:106
int reference_value_size() const
The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element.
Definition: FiniteElement.cpp:92
bool has_dof_reference_coordinates() const noexcept
Check if reference coordinates for dofs are defined.
Definition: FiniteElement.cpp:170