|
| CoordinateElement (std::shared_ptr< const basix::FiniteElement< T > > element) |
| Create a coordinate element from a Basix element.
|
|
| CoordinateElement (mesh::CellType celltype, int degree, basix::element::lagrange_variant type=basix::element::lagrange_variant::unset) |
| Create a Lagrange coordinate element.
|
|
virtual | ~CoordinateElement ()=default |
| Destructor.
|
|
mesh::CellType | cell_shape () const |
| Cell shape.
|
|
int | degree () const |
| The polynomial degree of the element.
|
|
int | dim () const |
| The dimension of the coordinate element space.
|
|
basix::element::lagrange_variant | variant () const |
| 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 tabulate .
|
|
void | tabulate (int nd, std::span< const T > X, std::array< std::size_t, 2 > shape, std::span< T > basis) const |
| Evaluate basis values and derivatives at set of points.
|
|
void | permute_subentity_closure (std::span< std::int32_t > d, std::uint32_t cell_info, mesh::CellType entity_type, int entity_index) const |
| Given the closure DOFs \(\tilde{d}\) of a cell sub-entity in reference ordering, this function computes the permuted degrees-of-freedom.
|
|
ElementDofLayout | create_dof_layout () const |
| Compute and return the dof layout.
|
|
void | pull_back_nonaffine (mdspan2_t< T > X, mdspan2_t< const T > x, mdspan2_t< const T > cell_geometry, double tol=1.0e-6, int maxit=15) const |
| Compute reference coordinates X for physical coordinates x for a non-affine map.
|
|
void | permute (std::span< std::int32_t > dofs, std::uint32_t cell_perm) const |
| Permute a list of DOF numbers on a cell.
|
|
void | permute_inv (std::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.
|
|
bool | is_affine () const noexcept |
|
|
template<typename U , typename V , typename W > |
static void | compute_jacobian (const U &dphi, const V &cell_geometry, W &&J) |
|
template<typename U , typename V > |
static void | compute_jacobian_inverse (const U &J, V &&K) |
| Compute the inverse of the Jacobian.
|
|
template<typename U > |
static double | compute_jacobian_determinant (const U &J, std::span< typename U::value_type > w) |
| Compute the determinant of the Jacobian.
|
|
template<typename U , typename V , typename W > |
static void | push_forward (U &&x, const V &cell_geometry, const W &phi) |
| Compute physical coordinates x for points X in the reference configuration.
|
|
template<typename U , typename V , typename W > |
static void | pull_back_affine (U &&X, const V &K, const std::array< T, 3 > &x0, const W &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} .
|
|
template<std::floating_point T>
class dolfinx::fem::CoordinateElement< T >
A CoordinateElement manages coordinate mappings for isoparametric cells.
- Todo
- A dof layout on a reference cell needs to be defined.
- Template Parameters
-
T | Floating point (real) type for the geometry and for the element basis. |
template<std::floating_point T>
void permute_subentity_closure |
( |
std::span< std::int32_t > | d, |
|
|
std::uint32_t | cell_info, |
|
|
mesh::CellType | entity_type, |
|
|
int | entity_index ) const |
Given the closure DOFs \(\tilde{d}\) of a cell sub-entity in reference ordering, this function computes the permuted degrees-of-freedom.
\[ d = P \tilde{d},\]
ordered to be consistent with the entity's mesh orientation, where \(P\) is a permutation matrix. This accounts for orientation discrepancies between the entity's cell and mesh orientation. All DOFs are rotated and reflected together, unlike permute
, which considered sub-entities independently.
- Parameters
-
[in,out] | d | Indices associated with the reference element degree-of-freedom (in). Indices associated with each physical element degree-of-freedom (out). |
[in] | cell_info | Permutation info for the cell |
[in] | entity_type | The cell type of the sub-entity |
[in] | entity_index | The local (with respect to the cell) index of the entity |
template<std::floating_point T>
template<typename U , typename V , typename W >
static void pull_back_affine |
( |
U && | X, |
|
|
const V & | K, |
|
|
const std::array< T, 3 > & | x0, |
|
|
const W & | x ) |
|
inlinestatic |
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}
.
- Parameters
-
[out] | X | Reference coordinates to compute (shape=(num_points, tdim) ), |
[in] | K | Inverse of the geometry Jacobian (shape=(tdim, / gdim) ). |
[in] | x0 | Physical coordinate of reference coordinate X0=(0, / 0, 0) . |
[in] | x | Physical coordinates (shape=(num_points, gdim)). |