DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | List of all members
CoordinateElement< T > Class Template Reference

#include <CoordinateElement.h>

Public Types

template<typename X >
using mdspan2_t
 mdspan typedef
 

Public Member Functions

 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.
 
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
 

Static Public Member Functions

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}.
 

Detailed Description

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
TFloating point (real) type for the geometry and for the element basis.

Member Typedef Documentation

◆ mdspan2_t

template<std::floating_point T>
template<typename X >
using mdspan2_t
Initial value:
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>

mdspan typedef

Constructor & Destructor Documentation

◆ CoordinateElement() [1/2]

template<std::floating_point T>
CoordinateElement ( std::shared_ptr< const basix::FiniteElement< T > > element)
explicit

Create a coordinate element from a Basix element.

Parameters
[in]elementBasix finite element

◆ CoordinateElement() [2/2]

template<std::floating_point T>
CoordinateElement ( mesh::CellType celltype,
int degree,
basix::element::lagrange_variant type = basix::element::lagrange_variant::unset )

Create a Lagrange coordinate element.

Parameters
[in]celltypeCell shape.
[in]degreePolynomial degree of the map.
[in]typeType of Lagrange element (see Basix documentation for possible types).

Member Function Documentation

◆ cell_shape()

template<std::floating_point T>
mesh::CellType cell_shape ( ) const

Cell shape.

Returns
The cell shape

◆ compute_jacobian()

template<std::floating_point T>
template<typename U , typename V , typename W >
static void compute_jacobian ( const U & dphi,
const V & cell_geometry,
W && J )
inlinestatic

Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives.

Parameters
[in]dphiDerivatives of the basis functions (shape=(tdim, num geometry nodes))
[in]cell_geometryThe cell nodes coordinates (shape=(num geometry nodes, gdim))
[out]JThe Jacobian. It must have shape=(gdim, tdim) and must initialized to zero

◆ compute_jacobian_determinant()

template<std::floating_point T>
template<typename U >
static double compute_jacobian_determinant ( const U & J,
std::span< typename U::value_type > w )
inlinestatic

Compute the determinant of the Jacobian.

Parameters
[in]JJacobian (shape=(gdim, tdim)).
[in]wWorking memory, required when gdim != tdim. Size must be at least 2 * gdim * tdim.
Returns
Determinant of J.

◆ compute_jacobian_inverse()

template<std::floating_point T>
template<typename U , typename V >
static void compute_jacobian_inverse ( const U & J,
V && K )
inlinestatic

Compute the inverse of the Jacobian.

Parameters
[in]JJacobian (shape=(gdim, tdim)).
[out]KInverse Jacobian (shape=(tdim, gdim)).

◆ degree()

template<std::floating_point T>
int degree ( ) const

The polynomial degree of the element.

Returns
The degree

◆ dim()

template<std::floating_point T>
int dim ( ) const

The dimension of the coordinate element space.

The number of basis function is returned. E.g., for a linear triangle cell the dimension will be 3.

Returns
Dimension of the coordinate element space.

◆ is_affine()

template<std::floating_point T>
bool is_affine ( ) const
inlinenoexcept

Check is geometry map is affine

Returns
True is geometry map is affine

◆ needs_dof_permutations()

template<std::floating_point T>
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.

◆ pull_back_affine()

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]XReference coordinates to compute (shape=(num_points, tdim)),
[in]KInverse of the geometry Jacobian (shape=(tdim, gdim)).
[in]x0Physical coordinate of reference coordinate X0=(0, 0, 0).
[in]xPhysical coordinates (shape=(num_points, gdim)).

◆ pull_back_nonaffine()

template<std::floating_point T>
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.

Parameters
[in,out]XThe reference coordinates to compute (shape=(num_points, tdim)).
[in]xPhysical coordinates (shape=(num_points, gdim)).
[in]cell_geometryCell nodes coordinates (shape=(num geometry nodes, gdim)).
[in]tolTolerance for termination of Newton method.
[in]maxitMaximum number of Newton iterations
Note
If convergence is not achieved within maxit, the function throws a runtime error.

◆ push_forward()

template<std::floating_point T>
template<typename U , typename V , typename W >
static void push_forward ( U && x,
const V & cell_geometry,
const W & phi )
inlinestatic

Compute physical coordinates x for points X in the reference configuration.

Parameters
[in,out]xThe physical coordinates of the reference points X (rank 2).
[in]cell_geometryCell node physical coordinates (rank 2).
[in]phiTabulated basis functions at reference points X (rank 2).

◆ tabulate()

template<std::floating_point T>
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.

Parameters
[in]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only.
[in]XThe points at which to compute the basis functions. The shape of X is (number of points, geometric dimension).
[in]shapeThe shape of X.
[out]basisThe array to fill with the basis function values. The shape can be computed using tabulate_shape.

◆ tabulate_shape()

template<std::floating_point T>
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.

Parameters
[in]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only
[in]num_pointsNumber of points at which to evaluate the basis functions.
Returns
Shape of the array to be filled by tabulate.

The documentation for this class was generated from the following files: