ndelement/traits.rs
1//! Traits.
2use crate::types::DofTransformation;
3use rlst::{Array, MutableArrayImpl, RlstScalar, ValueArrayImpl};
4use std::fmt::Debug;
5use std::hash::Hash;
6
7/// A finite element.
8pub trait FiniteElement {
9 /// The scalar type
10 type T: RlstScalar;
11
12 /// Cell type
13 ///
14 /// The cell type is typically defined through [ReferenceCellType](crate::types::ReferenceCellType).
15 type CellType: Debug + PartialEq + Eq + Clone + Copy + Hash;
16
17 /// The reference cell type, eg one of `Point`, `Interval`, `Triangle`, etc.
18 fn cell_type(&self) -> Self::CellType;
19
20 /// The number of basis functions.
21 fn dim(&self) -> usize;
22
23 /// The shape of the values returned by functions in $\mathcal{V}$.
24 ///
25 /// If the values are scalar an empty slice is returned.
26 fn value_shape(&self) -> &[usize];
27
28 /// The number of values returned.
29 ///
30 /// If (for example) `value_shape` is `[3, 4]` then `value_size` is $3\times4 = 12$.
31 /// If `value_shape` returns an empty array (ie the shape functions are scalar) the
32 // convention is used that the product of the elements of an empty array is 1.
33 fn value_size(&self) -> usize;
34
35 /// Tabulate the values of the basis functions and their derivatives at a set of points
36 ///
37 /// - `points` is a two-dimensional array where each column contains the coordinates of one point.
38 /// - `nderivs` is the desired number of derivatives (0 for no derivatives, 1 for all first order derivatives, etc.).
39 /// - `data` is the 4-dimensional output array. The first dimension is the total number of partial derivatives,
40 /// the second dimension is the number of evaluation points, the third dimension is the number of
41 /// basis functions of the element, and the fourth dimension is the value size of the basis function output.
42 /// For example, `data[3, 2, 1, 0]` returns the 0th value of the third partial derivative on the point with index 2 for the
43 /// basis function with index 1.
44 ///
45 /// ## Remark
46 /// Let $d^{i + k} = dx^{i}dy^{j}$ be a derivative with respect to $x$, $y$ in two dimensions and
47 /// $d^{i + k + j} = dx^{i}dy^{j}dz^{k}$ be a derivative with respect to $x$, $y$, and $z$ in three dimensions.
48 /// Then the corresponding index $\ell$ in the first dimension of the `data` array is computed as follows.
49 ///
50 /// - Triangle: $l = (i + j + 1) * (i + j) / 2 + j$
51 /// - Quadrilateral: $l = i * (n + 1) + j$
52 /// - Tetrahedron: $l = (i + j + k) * (i + j + k + 1) * (i + j + k + 2) / 6 + (j + k) * (j + k + 1) / 2 + k$
53 /// - Hexahedron $l = i * (n + 1) * (n + 1) + j * (n + 1) + k$.
54 ///
55 /// For the quadrilaterals and hexahedra, the parameter $n$ denotes the Lagrange superdegree.
56 fn tabulate<
57 TGeo: RlstScalar,
58 Array2Impl: ValueArrayImpl<TGeo, 2>,
59 Array4MutImpl: MutableArrayImpl<Self::T, 4>,
60 >(
61 &self,
62 points: &Array<Array2Impl, 2>,
63 nderivs: usize,
64 data: &mut Array<Array4MutImpl, 4>,
65 );
66
67 /// Get the required shape for a tabulation array.
68 fn tabulate_array_shape(&self, nderivs: usize, npoints: usize) -> [usize; 4];
69
70 /// Return the dof indices that are associated with the subentity with index `entity_number` and dimension `entity_dim`.
71 ///
72 /// - For `entity_dim = 0` this returns the degrees of freedom (dofs) associated with the corresponding point.
73 /// - For `entity_dim = 1` this returns the dofs associated with the corresponding edge.
74 /// - For `entity_dim = 2` this returns the dofs associated with the corresponding face.
75 ///
76 /// Note that this does not return dofs on the boundary of an entity, that means that (for example)
77 /// for an edge the dofs associated with the two vertices at the boundary of the edge are not returned.
78 /// To return also the boundary dofs use [FiniteElement::entity_closure_dofs].
79 fn entity_dofs(&self, entity_dim: usize, entity_number: usize) -> Option<&[usize]>;
80
81 /// The DOFs that are associated with a closure of a subentity of the reference cell.
82 ///
83 /// This method is similar to [FiniteElement::entity_dofs], but it includes the dofs
84 /// associated with the boundary of an entity. For an edge (for example) it returns the dofs associated
85 /// with the vertices at the boundary of the edge (as well as the dofs associated with the edge itself).
86 fn entity_closure_dofs(&self, entity_dim: usize, entity_number: usize) -> Option<&[usize]>;
87}
88
89/// A finite element that is mapped from a reference cell.
90pub trait MappedFiniteElement: FiniteElement {
91 /// Transformation type
92 ///
93 /// The Transformation type specifies possible transformations of the dofs on the reference element.
94 /// In most cases these will be rotations and reflections as defined in [Transformation](crate::types::Transformation).
95 type TransformationType: Debug + PartialEq + Eq + Clone + Copy + Hash;
96
97 /// The smallest degree Lagrange space that contains all possible polynomials of the finite element's polynomial space.
98 ///
99 /// Details on the definition of the degree of Lagrange spaces of finite elements are
100 /// given [here](https://defelement.org/ciarlet.html#The+degree+of+a+finite+element).
101 fn lagrange_superdegree(&self) -> usize;
102
103 /// Push function values forward to a physical cell.
104 ///
105 /// For Lagrange elements, this is an identity map. For many other element types, the function values
106 /// on the reference cell differ from those on the physical cell: for example Nedlec elements use a covariant
107 /// Piola transform. This method implements the appropriate transformation for the element.
108 ///
109 /// - `reference_values`: The values on the reference cell. The shape of this input is the same as the `data` input to the function
110 /// [[FiniteElement::tabulate]].
111 /// - `nderivs`: The number of derivatives.
112 /// - `jacobians` should have shape [geometry_dimension, entity_topology_dimension, npts] and use column-major ordering;
113 /// - `jacobian_determinants` should have shape \[npts\];
114 /// - `inverse_jacobians` should have shape [entity_topology_dimension, geometry_dimension, npts] and use column-major ordering;
115 /// - `physical_values`: The output array of the push operation. This shape of this array is the same as the `reference_values`
116 /// but with reference value size replaced by the physical value size
117 fn push_forward<
118 TGeo: RlstScalar,
119 Array3GeoImpl: ValueArrayImpl<TGeo, 3>,
120 Array4Impl: ValueArrayImpl<Self::T, 4>,
121 Array4MutImpl: MutableArrayImpl<Self::T, 4>,
122 >(
123 &self,
124 reference_values: &Array<Array4Impl, 4>,
125 nderivs: usize,
126 jacobians: &Array<Array3GeoImpl, 3>,
127 jacobian_determinants: &[TGeo],
128 inverse_jacobians: &Array<Array3GeoImpl, 3>,
129 physical_values: &mut Array<Array4MutImpl, 4>,
130 );
131
132 /// Pull function values back to the reference cell.
133 ///
134 /// This is the inverse operation to [MappedFiniteElement::push_forward].
135 fn pull_back<
136 TGeo: RlstScalar,
137 Array3GeoImpl: ValueArrayImpl<TGeo, 3>,
138 Array4Impl: ValueArrayImpl<Self::T, 4>,
139 Array4MutImpl: MutableArrayImpl<Self::T, 4>,
140 >(
141 &self,
142 physical_values: &Array<Array4Impl, 4>,
143 nderivs: usize,
144 jacobians: &Array<Array3GeoImpl, 3>,
145 jacobian_determinants: &[TGeo],
146 inverse_jacobians: &Array<Array3GeoImpl, 3>,
147 reference_values: &mut Array<Array4MutImpl, 4>,
148 );
149
150 /// The value shape on a physical cell
151 fn physical_value_shape(&self, gdim: usize) -> Vec<usize>;
152
153 /// The value size on a physical cell
154 fn physical_value_size(&self, gdim: usize) -> usize {
155 let mut vs = 1;
156 for i in self.physical_value_shape(gdim) {
157 vs *= i;
158 }
159 vs
160 }
161
162 /// The DOF transformation for a sub-entity.
163 fn dof_transformation(
164 &self,
165 entity: Self::CellType,
166 transformation: Self::TransformationType,
167 ) -> Option<&DofTransformation<Self::T>>;
168
169 /// Apply permutation parts of DOF transformations.
170 fn apply_dof_permutations<T>(&self, data: &mut [T], cell_orientation: i32);
171
172 /// Apply non-permutation parts of DOF transformations.
173 fn apply_dof_transformations(&self, data: &mut [Self::T], cell_orientation: i32);
174
175 /// Apply DOF transformations.
176 fn apply_dof_permutations_and_transformations(
177 &self,
178 data: &mut [Self::T],
179 cell_orientation: i32,
180 ) {
181 self.apply_dof_permutations(data, cell_orientation);
182 self.apply_dof_transformations(data, cell_orientation);
183 }
184}
185
186/// A family of finite elements defined on various cell types, with appropriate continuity
187/// between different cell types.
188pub trait ElementFamily {
189 /// The scalar type
190 type T: RlstScalar;
191 /// Cell type
192 type CellType: Debug + PartialEq + Eq + Clone + Copy + Hash;
193 /// The finite element type
194 type FiniteElement: FiniteElement<T = Self::T, CellType = Self::CellType> + 'static;
195
196 /// Create an element for the given cell type.
197 fn element(&self, cell_type: Self::CellType) -> Self::FiniteElement;
198}
199
200/// A map from the reference cell to physical cells.
201pub trait Map {
202 /// Push values from the reference cell to physical cells.
203 fn push_forward<
204 T: RlstScalar,
205 TGeo: RlstScalar,
206 Array3GeoImpl: ValueArrayImpl<TGeo, 3>,
207 Array4Impl: ValueArrayImpl<T, 4>,
208 Array4MutImpl: MutableArrayImpl<T, 4>,
209 >(
210 &self,
211 reference_values: &Array<Array4Impl, 4>,
212 nderivs: usize,
213 jacobians: &Array<Array3GeoImpl, 3>,
214 jacobian_determinants: &[TGeo],
215 inverse_jacobians: &Array<Array3GeoImpl, 3>,
216 physical_values: &mut Array<Array4MutImpl, 4>,
217 );
218
219 /// Pull values back to the reference cell.
220 fn pull_back<
221 T: RlstScalar,
222 TGeo: RlstScalar,
223 Array3GeoImpl: ValueArrayImpl<TGeo, 3>,
224 Array4Impl: ValueArrayImpl<T, 4>,
225 Array4MutImpl: MutableArrayImpl<T, 4>,
226 >(
227 &self,
228 physical_values: &Array<Array4Impl, 4>,
229 nderivs: usize,
230 jacobians: &Array<Array3GeoImpl, 3>,
231 jacobian_determinants: &[TGeo],
232 inverse_jacobians: &Array<Array3GeoImpl, 3>,
233 reference_values: &mut Array<Array4MutImpl, 4>,
234 );
235
236 /// The value shape on a physical cell.
237 fn physical_value_shape(&self, gdim: usize) -> Vec<usize>;
238}