|
DOLFINx
0.4.1
DOLFINx C++ interface
|
#include <CoordinateElement.h>
Public Member Functions | |
| CoordinateElement (std::shared_ptr< const basix::FiniteElement > element) | |
| Create a coordinate element from a Basix element. More... | |
| CoordinateElement (mesh::CellType celltype, int degree, basix::element::lagrange_variant type=basix::element::lagrange_variant::equispaced) | |
| Create a Lagrange coordinate element. More... | |
| virtual | ~CoordinateElement ()=default |
| Destructor. | |
| mesh::CellType | cell_shape () const |
| Cell shape. More... | |
| int | degree () const |
| The polynomial degree of the element. | |
| int | dim () const |
| The dimension of the geometry element space. More... | |
| basix::element::lagrange_variant | variant () const |
| The variant of the element. | |
| std::array< std::size_t, 4 > | tabulate_shape (std::size_t nd, std::size_t num_points) const |
Shape of array to fill when calling FiniteElement::tabulate More... | |
| xt::xtensor< double, 4 > | tabulate (int nd, const xt::xtensor< double, 2 > &X) const |
| Evaluate basis values and derivatives at set of points. More... | |
| void | tabulate (int nd, const xt::xtensor< double, 2 > &X, xt::xtensor< double, 4 > &basis) const |
| Evaluate basis values and derivatives at set of points. More... | |
| ElementDofLayout | create_dof_layout () const |
| Compute and return the dof layout. | |
| void | pull_back_nonaffine (xt::xtensor< double, 2 > &X, const xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry, double tol=1.0e-8, int maxit=10) const |
| Compute reference coordinates X for physical coordinates x for a non-affine map. More... | |
| void | permute_dofs (const xtl::span< std::int32_t > &dofs, std::uint32_t cell_perm) const |
| Permutes a list of DOF numbers on a cell. | |
| void | unpermute_dofs (const xtl::span< std::int32_t > &dofs, std::uint32_t cell_perm) const |
| Reverses a DOF permutation. | |
| bool | needs_dof_permutations () const |
| Indicates whether the geometry DOF numbers on each cell need permuting. More... | |
| bool | is_affine () const noexcept |
| Check is geometry map is affine. More... | |
Static Public Member Functions | |
| template<typename U , typename V , typename W > | |
| static void | compute_jacobian (const U &dphi, const V &cell_geometry, W &&J) |
| Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives. More... | |
| template<typename U , typename V > | |
| static void | compute_jacobian_inverse (const U &J, V &&K) |
| Compute the inverse of the Jacobian. More... | |
| template<typename U > | |
| static double | compute_jacobian_determinant (const U &J) |
| Compute the determinant of the Jacobian. More... | |
| static void | push_forward (xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry, const xt::xtensor< double, 2 > &phi) |
| Compute physical coordinates x for points X in the reference configuration. More... | |
| static std::array< double, 3 > | x0 (const xt::xtensor< double, 2 > &cell_geometry) |
| Compute the physical coordinate of the reference point X=(0 , 0, 0) More... | |
| static void | pull_back_affine (xt::xtensor< double, 2 > &X, const xt::xtensor< double, 2 > &K, const std::array< double, 3 > &x0, const xt::xtensor< double, 2 > &x) |
Compute reference coordinates X for physical coordinates x for an affine map. For the affine case, x = J X + x0, and this function computes X = K(x -x0) where K = J^{-1}. More... | |
A CoordinateElement manages coordinate mappings for isoparametric cells.
|
explicit |
Create a coordinate element from a Basix element.
| [in] | element | Element from Basix |
| CoordinateElement | ( | mesh::CellType | celltype, |
| int | degree, | ||
| basix::element::lagrange_variant | type = basix::element::lagrange_variant::equispaced |
||
| ) |
Create a Lagrange coordinate element.
| [in] | celltype | The cell shape |
| [in] | degree | Polynomial degree of the map |
| [in] | type | The type of Lagrange element (see Basix documentation for possible types) |
| mesh::CellType cell_shape | ( | ) | const |
Cell shape.
|
inlinestatic |
Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives.
| [in] | dphi | Derivatives of the basis functions (shape=(tdim, num geometry nodes)) |
| [in] | cell_geometry | The cell nodes coordinates (shape=(num geometry nodes, gdim)) |
| [out] | J | The Jacobian. It must have shape=(gdim, tdim) and must initialized to zero |
|
inlinestatic |
Compute the determinant of the Jacobian.
| [in] | J | Jacobian (shape=(gdim, tdim)) |
J
|
inlinestatic |
Compute the inverse of the Jacobian.
| [in] | J | The Jacobian (shape=(gdim, tdim)) |
| [out] | K | The Jacobian (shape=(tdim, gdim)) |
| int dim | ( | ) | const |
The dimension of the geometry element space.
The number of basis function is returned. E.g., for a linear triangle cell the dimension will be 3.
|
inlinenoexcept |
Check is geometry map is affine.
| bool needs_dof_permutations | ( | ) | const |
Indicates whether the geometry DOF numbers on each cell need permuting.
For higher order geometries (where there is more than one DOF on a subentity of the cell), this will be true.
|
static |
Compute reference coordinates X for physical coordinates x for an affine map. For the affine case, x = J X + x0, and this function computes X = K(x -x0) where K = J^{-1}.
| [out] | X | The reference coordinates to compute (shape=(num_points, tdim)) |
| [in] | K | The inverse of the geometry Jacobian (shape=(tdim, gdim)) |
| [in] | x0 | The physical coordinate of reference coordinate X0=(0, 0, 0). |
| [in] | x | The physical coordinates (shape=(num_points, gdim)) |
| void pull_back_nonaffine | ( | xt::xtensor< double, 2 > & | X, |
| const xt::xtensor< double, 2 > & | x, | ||
| const xt::xtensor< double, 2 > & | cell_geometry, | ||
| double | tol = 1.0e-8, |
||
| int | maxit = 10 |
||
| ) | const |
Compute reference coordinates X for physical coordinates x for a non-affine map.
| [in,out] | X | The reference coordinates to compute (shape=(num_points, tdim)) |
| [in] | x | The physical coordinates (shape=(num_points, gdim)) |
| [in] | cell_geometry | The cell nodes coordinates (shape=(num geometry nodes, gdim)) |
| [in] | tol | Tolerance for termination of Newton method. |
| [in] | maxit | Maximum number of Newton iterations |
|
static |
Compute physical coordinates x for points X in the reference configuration.
| [in,out] | x | The physical coordinates of the reference points X |
| [in] | cell_geometry | The cell node coordinates (physical) |
| [in] | phi | Tabulated basis functions at reference points X |
| xt::xtensor< double, 4 > tabulate | ( | int | nd, |
| const xt::xtensor< double, 2 > & | X | ||
| ) | const |
Evaluate basis values and derivatives at set of points.
| [in] | nd | The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. |
| [in] | X | The points at which to compute the basis functions. The shape of x is (number of points, geometric dimension). |
| void tabulate | ( | int | nd, |
| const xt::xtensor< double, 2 > & | X, | ||
| xt::xtensor< double, 4 > & | basis | ||
| ) | const |
Evaluate basis values and derivatives at set of points.
| [in] | nd | The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. |
| [in] | X | The points at which to compute the basis functions. The shape of X is (number of points, geometric dimension). |
| [out] | basis | The array to fill with the basis function values. The shape can be computed using FiniteElement::tabulate_shape |
| std::array< std::size_t, 4 > tabulate_shape | ( | std::size_t | nd, |
| std::size_t | num_points | ||
| ) | const |
Shape of array to fill when calling FiniteElement::tabulate
| [in] | nd | The order of derivatives, up to and including, to compute. Use 0 for the basis functions only |
| [in] | num_points | Number of points at which to evaluate the basis functions |
FiniteElement::tabulate
|
static |
Compute the physical coordinate of the reference point X=(0 , 0, 0)
| [in] | cell_geometry | The cell nodes coordinates (shape=(num geometry nodes, gdim)) |