Mapping
Mappings
Bcube.cell_normal
— Methodcell_normal(ctype::AbstractEntityType, cnodes, ξ) where {T, N}
Compute the cell normal vector of an entity of topology dimension equals to (n-1) in a n-D space, for instance a curve in a 2D space. This vector is expressed in the cell-reference coordinate system.
Do not confuse the cell normal vector with the cell-side (i.e face) normal vector.
Topology dimension 1
the curve direction vector, u, is J/||J||. Then n = [-u.y, u.x].
Bcube.center
— Methodcenter(ctype::AbstractEntityType, cnodes)
Return the center of the AbstractEntityType
by mapping the center of the corresponding Shape
.
Warning
Do not use this function on a face of a cell : since the face is of dimension "n-1", the mapping won't be appropriate.
Bcube.get_cell_centers
— Methodget_cell_centers(mesh::Mesh)
Get mesh cell centers coordinates (assuming perfectly flat cells)
Bcube.mapping
— Methodmapping(ctype::AbstractEntityType, cnodes, ξ)
mapping(cshape::AbstractShape, cnodes, ξ)
Map the reference shape on the local shape.
Implementation
This function must be implemented for all shape.
::Bar2_t
Map the reference 2-nodes bar [-1,1] on the local bar:
\[F(\xi) = \dfrac{x_r - x_l}{2} \xi + \dfrac{x_r + x_l}{2}\]
::Tri3_t
Map the reference 3-nodes Triangle [0,1] x [0,1] on the local triangle.
\[F(\xi \\ \eta) = (1 - \xi - \eta) M_1 + x M_2 + y M_3\]
::Quad4_t
Map the reference 4-nodes square [-1,1] x [-1,1] on the 4-quadrilateral.
::Tri6_t
Map the reference 6-nodes triangle [0,1] x [0,1] on the P2 curved-triangle. $F(\xi) = \sum \lambda_i(\xi) x_i$ where $\lambda_i$ are the Lagrange P2 shape functions and $x_i$ are the local curved-triangle vertices' coordinates.
::Quad9_t
Map the reference 4-nodes square [-1,1] x [-1,1] on the P2 curved-quadrilateral. $F(\xi) = \sum \lambda_i(\xi) x_i$ where $\lambda_i$ are the Lagrange P2 shape functions and $x_i$ are the local curved-quadrilateral vertices' coordinates.
::Hexa8_t
Map the reference 8-nodes cube [-1,1] x [-1,1] x [-1,1] on the 8-hexa.
::Hexa27_t
Map the reference 8-nodes cube [-1,1] x [-1,1] x [-1,1] on the 27-hexa.
::Penta6_t
Map the reference 6-nodes prism [0,1] x [0,1] x [-1,1] on the 6-penta (prism).
::Pyra5_t
Map the reference 5-nodes pyramid [-1,1] x [-1,1] x [0,1] on the 5-pyra. See https://www.math.u-bordeaux.fr/~durufle/montjoie/pyramid.php
Bcube.mapping
— Methodmapping(nodes, ::Tetra4_t, ξ)
Map the reference 4-nodes Tetraahedron [0,1] x [0,1] x [0,1] on the local triangle.
```
Bcube.mapping_det_jacobian
— Methodmapping_det_jacobian(ctype::AbstractEntityType, cnodes, ξ)
mapping_det_jacobian(::TopologyStyle, ctype::AbstractEntityType, cnodes, ξ)
Absolute value of the determinant of the mapping Jacobian matrix, expressed in the reference element.
Implementation
For a volumic cell (line in 1D, quad in 2D, cube in 3D), the default version uses mapping_jacobian
, but the function can be specified for each shape.
For a curvilinear cell (line in 2D, or in 3D), the formula is always J = ||F'|| where F is the mapping. Hence we always fallback to the "standard" version, like for the volumic case.
Finally, the surfacic cell (quad in 3D) needs a special treatment, see mapping_jacobian_hypersurface
.
::Bar2_t
Absolute value of the determinant of the mapping Jacobian matrix for the reference 2-nodes bar [-1,1] to the local bar mapping. $|det(J(\xi))| = \dfrac{|x_r - x_l|}{2}$
::Tri3_t
Absolute value of the determinant of the mapping Jacobian matrix for the the reference 3-nodes Triangle [0,1] x [0,1] to the local triangle mapping.
$|J| = |(x_2 - x_1) (y_3 - y_1) - (x_3 - x_1) (y_2 - y_1)|$
Bcube.mapping_face
— Methodmapping_face(cshape::AbstractShape, side)
mapping_face(cshape::AbstractShape, side, permutation)
Build a mapping from the face reference element (corresponding to the side
-th face of cshape
) to the cell reference element (i.e the cshape
).
Build a mapping from the face reference element (corresponding to the side
-th face of cshape
) to the cell reference element (i.e the cshape
). If permutation
is present, the mapping is built using this permutation.
Bcube.mapping_inv
— Methodmapping_inv(::AbstractEntityType, cnodes, x)
Map the local shape on the reference shape.
Implementation
This function does not have to be implemented for all shape.
::Bar2_t
Map the local bar on the reference 2-nodes bar [-1,1]: $F^{-1}(x) = \dfrac{2x - x_r - x_l}{x_r - x_l}$
::Tri3_t
Map the local triangle on the reference 3-nodes Triangle [0,1] x [0,1].
TODO: check this formulae with SYMPY
\[F^{-1} \begin{pmatrix} x \\ y \end{pmatrix} = \frac{1}{x_1(y_2-y_3) + x_2(y_3-y_1) + x_3(y_1-y_2)} \begin{pmatrix} (y_3-y_1)x + (x_1 - x_3)y + (x_3 y_1 - x_1 y_3) \\ (y_1-y_2)x + (x_2 - x_1)x + (x_1 y_2 - x_2 y_1) \end{pmatrix}\]
::Quad4_t
Map the PARALLELOGRAM quadrilateral on the reference 4-nodes square [-1,1] x [-1,1]. Warning : this mapping is only corrects for parallelogram quadrilateral, not for any quadrilateral.
TODO: check this formulae with SYMPY
\[F^{-1} \begin{pmatrix} x \\ y \end{pmatrix} = \begin{pmatrix} a_1 x + b_1 y + c_1 \\ a_2 x + b_2 y + c_2 \end{pmatrix}\]
with
\[\begin{aligned} a_1 & = \dfrac{-2 (y_3-y_2)}{\Delta} \\ b_1 & = \dfrac{2 (x_3-x_2)}{\Delta} \\ c_1 & = -1 - a_1 x_1 - b_1 y_1 \\ a_2 & = \dfrac{-2 (y_1-y_2)}{\Delta} \\ b_2 & = \dfrac{2 (x_1 - x_2)}{\Delta} \\ c_2 & = -1 - a_2 x_1 - b_2 y_1 \end{aligned}\]
where $\Delta = (x_1 - x_2)(y_3 - y_2) - (x_3 - x_2)(y_1 - y_2)$
Bcube.mapping_inv_jacobian
— Methodmapping_inv_jacobian(ctype::AbstractEntityType, cnodes, x)
Jacobian matrix of the inverse mapping : $\dfrac{\partial F_i^{-1}}{\partial x_j}$
Contrary to mapping_jacobian_inv
, this function is not always defined because the inverse mapping, F^-1, is not always defined.
Implementation
Default version using LinearAlgebra to inverse the matrix, but can be specified for each shape (if it exists).
::Bar2_t
Mapping's jacobian matrix for the local bar to the reference 2-nodes bar [-1, 1].
\[\dfrac{\partial F^{-1}}{\partial x} = \dfrac{2}{x_r - x_l}\]
::Tri3_t
Mapping's jacobian matrix for the local triangle to the reference 3-nodes Triangle [0,1] x [0,1] mapping.
TODO: check this formulae with SYMPY
\[\frac{\partial F_i^{-1}}{\partial x_j} = \frac{1}{x_1 (y_2 - y_3) + x_2 (y_3 - y_1) + x_3 (y_1 - y_2)} \begin{pmatrix} y_3 - y_1 & x_1 - x_3 \\ y_1 - y_2 & x_2 - x_1 \end{pmatrix}\]
Bcube.mapping_jacobian
— Methodmapping_jacobian(ctype::AbstractEntityType, cnodes, ξ)
Jacobian matrix of the mapping : $\dfrac{\partial F_i}{\partial \xi_j}$.
Implementation
Default version using ForwardDiff, but can be specified for each shape.
::Bar2_t
Mapping's jacobian matrix for the reference 2-nodes bar [-1, 1] to the local bar. $\dfrac{\partial F}{\partial \xi} = \dfrac{x_r - x_l}{2}$
::Bar3_t
Mapping's jacobian matrix for the reference 2-nodes bar [-1, 1] to the local bar. ``\dfrac{\partial F}{\partial \xi} = \frac{1}{2} \left( (2\xi - 1) M1 + (2\xi + 1)M2 - 4 \xi M_3\right)
::Tri3_t
Mapping's jacobian matrix for the reference 3-nodes Triangle [0,1] x [0,1] to the local triangle mapping.
\[\dfrac{\partial F_i}{\partial \xi_j} = \begin{pmatrix} M_2 - M_1 & M_3 - M_1 \end{pmatrix}\]
::Quad4_t
Mapping's jacobian matrix for the reference square [-1,1] x [-1,1] to the 4-quadrilateral
\[\frac{\partial F}{\partial \xi} = -M_1 + M_2 + M_3 - M_4 + \eta (M_1 - M_2 + M_3 - M_4)\]
\[\frac{\partial F}{\partial \eta} = -M_1 - M_2 + M_3 + M_4 + \xi (M_1 - M_2 + M_3 - M_4)\]
Bcube.mapping_jacobian_hypersurface
— Methodmapping_jacobian_hypersurface(ctype, cnodes, ξ)
"Augmented" jacobian matrix of the mapping.
Let's consider a $\mathbb{R}^2$ surface in $\mathbb{R}^3$. The mapping $F_\Gamma(\xi, \eta)$ maps the reference coordinate system to the physical coordinate system. It's jacobian $J_\Gamma$ is not squared. We can 'extend' this mapping to reach any point in $\mathbb{R}^3$ (and not only the surface) using
\[F(\xi, \eta, \zeta) = F_\Gamma(\xi, \eta) + \zeta \nu\]
where $\nu$ is the conormal. Then the restriction of the squared jacobian of $F$ to the surface is simply
\[J|_\Gamma = (J_\Gamma~~\nu)\]
Bcube.mapping_jacobian_inv
— Methodmapping_jacobian_inv(ctype::AbstractEntityType, cnodes, ξ)
Inverse of the mapping jacobian matrix. This is not exactly equivalent to the mapping_inv_jacobian
since this function is evaluated in the reference element.
Implementation
Default version using ForwardDiff, but can be specified for each shape.
::Bar2_t
Inverse of mapping's jacobian matrix for the reference 2-nodes bar [-1, 1] to the local bar.
\[\dfrac{\partial F}{\partial \xi}^{-1} = \dfrac{2}{x_r - x_l}\]
::Bar3_t
Inverse of mapping's jacobian matrix for the reference 2-nodes bar [-1, 1] to the local bar.
\[\dfrac{\partial F}{\partial \xi}^{-1} = \frac{2}{(2\xi - 1) M_1 + (2\xi + 1)M_2 - 4 \xi M_3}\]
::Tri3_t
Inverse of mapping's jacobian matrix for the reference 3-nodes Triangle [0,1] x [0,1] to the local triangle mapping.
\[\dfrac{\partial F_i}{\partial \xi_j}^{-1} = \frac{1}{(x_1 - x_2)(y_1 - y_3) - (x_1 - x_3)(y_1 - y_2)} \begin{pmatrix} -y_1 + y_3 & x_1 - x_3 \\ y_1 - y_2 & -x_1 + x_2 \end{pmatrix}\]
::Quad4_t
Inverse of mapping's jacobian matrix for the reference square [-1,1] x [-1,1] to the 4-quadrilateral
Bcube.normal
— Methodnormal(ctype::AbstractEntityType, cnodes, iside, ξ)
normal(::TopologyStyle, ctype::AbstractEntityType, cnodes, iside, ξ)
Normal vector of the iside
th face of a cell, evaluated at position ξ
in the face reference element. So for the normal vector of the face of triangle living in a 3D space, ξ
will be 1D (because the face is a line, which 1D).
Beware this function needs the nodes cnodes
and the type ctype
of the cell (and not of the face).
TODO: If iside
is positive, then the outward normal (with respect to the cell) is returned, otherwise the inward normal is returned.
::isCurvilinear
Note that the "face" normal vector of a curve is the "direction" vector at the given extremity.
::isVolumic
$n^{loc} = J^{-\intercal} n^{ref}$
Reference to physical
Bcube.interpolate
— Methodinterpolate(λ, q)
interpolate(λ, q, ncomps)
Create the interpolation function from a set of value on dofs and the shape functions, given by:
\[ f(x) = \sum_{i=1}^N q_i \lambda_i(x)\]
So q
is a vector whose size equals the number of dofs in the cell.
If ncomps
is present, create the interpolation function for a vector field given by a set of value on dofs and the shape functions.
The interpolation formulae is the same than interpolate(λ, q)
but the result is a vector function. Here q
is a vector whose size equals the total number of dofs in the cell (all components mixed).
Note that the result function is expressed in the same coordinate system as the input shape functions (i.e reference or local).
Bcube.∂fξ_∂x
— Method∂fξ_∂x(f, n::Val{N}, ctype::AbstractEntityType, cnodes, ξ) where N
Compute the gradient, with respect to the physical coordinates, of a function f
on a point in the reference domain. N
is the size of the codomain of f
.
Bcube.∂λξ_∂x
— Function∂λξ_∂x(::AbstractFunctionSpace, ::Val{N}, ctype::AbstractEntityType, cnodes, ξ) where N
Gradient, with respect to the physical coordinate system, of the shape functions associated to the FunctionSpace
.
Depending on the value of N
, the shape functions are interpreted as associated to a scalar FESpace (N = 1) or a vector FESpace. For a vector FESpace, the result gradient is an array of size (n*Nc, n, d)
where Nc
is the number of dofs of one component (i.e scalar case), n
is the size of the FESpace, and d
the number of spatial dimensions.
Default version : the gradient shape functions are "replicated".
Specialize with a given FESpace for a custom behaviour.
Implementation
We cannot use the topology_style
to dispatch because this style is too specific to integration methods. For instance for the integration it is important the consider any line as isCurvilinear
. However for the gradient computation we must distinguish a line in 1D, a line in 2D and a line in 3D...