Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
DOLFINx C++ interface
Public Member Functions | Static Public Member Functions | Public Attributes | List of all members
dolfinx::fem::CoordinateElement Class Reference

This class manages coordinate mappings for isoparametric cells. More...

#include <CoordinateElement.h>

Public Member Functions

 CoordinateElement (std::shared_ptr< basix::FiniteElement > element)
 Create a coordinate element from a Basix element. More...
 
 CoordinateElement (mesh::CellType celltype, int degree)
 Create a Lagrage coordinate element. More...
 
virtual ~CoordinateElement ()=default
 Destructor.
 
mesh::CellType cell_shape () const
 Cell shape. More...
 
int topological_dimension () const
 Return the topological dimension of the cell shape.
 
xt::xtensor< double, 4 > tabulate (int n, const xt::xtensor< double, 2 > &X) const
 Compute basis values and derivatives at set of points. More...
 
void compute_jacobian (const xt::xtensor< double, 4 > &dphi, const xt::xtensor< double, 2 > &cell_geometry, xt::xtensor< double, 3 > &J) const
 Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives. More...
 
void compute_jacobian_inverse (const xt::xtensor< double, 3 > &J, xt::xtensor< double, 3 > &K) const
 Compute the inverse of the Jacobian. If the coordinate element is affine, it computes the inverse at only one point. More...
 
void compute_jacobian_determinant (const xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ) const
 Compute the determinant of the Jacobian. If the coordinate element is affine, it computes the determinant at only one point. More...
 
ElementDofLayout dof_layout () const
 Return the dof layout.
 
void pull_back (xt::xtensor< double, 2 > &X, xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ, xt::xtensor< double, 3 > &K, const xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry) const
 Compute reference coordinates X, and J, detJ and K for physical coordinates x.
 
void permute_dofs (xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const
 Permutes a list of DOF numbers on a cell.
 
void unpermute_dofs (xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const
 Reverses a DOF permutation.
 
bool needs_permutation_data () const
 Indicates whether the coordinate map needs permutation data passing in (for higher order geometries)
 

Static Public Member Functions

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

Public Attributes

double non_affine_atol = 1.0e-8
 Absolute increment stopping criterium for non-affine Newton solver.
 
int non_affine_max_its = 10
 Maximum number of iterations for non-affine Newton solver.
 

Detailed Description

This class manages coordinate mappings for isoparametric cells.

Constructor & Destructor Documentation

◆ CoordinateElement() [1/2]

CoordinateElement::CoordinateElement ( std::shared_ptr< basix::FiniteElement >  element)
explicit

Create a coordinate element from a Basix element.

Parameters
[in]elementElement from Basix

◆ CoordinateElement() [2/2]

CoordinateElement::CoordinateElement ( mesh::CellType  celltype,
int  degree 
)

Create a Lagrage coordinate element.

Parameters
[in]celltypeThe cell shape
[in]degreePolynomial degree of the map

Member Function Documentation

◆ cell_shape()

mesh::CellType CoordinateElement::cell_shape ( ) const

Cell shape.

Returns
The cell shape

◆ compute_jacobian()

void CoordinateElement::compute_jacobian ( const xt::xtensor< double, 4 > &  dphi,
const xt::xtensor< double, 2 > &  cell_geometry,
xt::xtensor< double, 3 > &  J 
) const

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

Parameters
[in]dphiPre-computed first order derivatives of basis functions at set of points. The shape of dphi is (tdim, number of points, number of basis fn, 1).
[in]cell_geometryCoordinates/geometry
[in,out]JThe Jacobian The shape of J is (number of points, geometric dimension, topological dimenson).

◆ compute_jacobian_determinant()

void CoordinateElement::compute_jacobian_determinant ( const xt::xtensor< double, 3 > &  J,
xt::xtensor< double, 1 > &  detJ 
) const

Compute the determinant of the Jacobian. If the coordinate element is affine, it computes the determinant at only one point.

Parameters
[in]JPolynomial degree of the map The shape of J is (number of points, geometric dimension, topological dimenson).
[in,out]detJJacobian The shape of detJ is (number of points)

◆ compute_jacobian_inverse()

void CoordinateElement::compute_jacobian_inverse ( const xt::xtensor< double, 3 > &  J,
xt::xtensor< double, 3 > &  K 
) const

Compute the inverse of the Jacobian. If the coordinate element is affine, it computes the inverse at only one point.

Parameters
[in]JThe Jacobian The shape of J is (number of points, geometric dimension, topological dimenson).
[in,out]KThe inverse of the Jacobian The shape of J is (number of points, tpological dimension, geometrical dimenson).

◆ push_forward()

void CoordinateElement::push_forward ( xt::xtensor< double, 2 > &  x,
const xt::xtensor< double, 2 > &  cell_geometry,
const xt::xtensor< double, 2 > &  phi 
)
static

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

Parameters
[in,out]xThe physical coordinates of the reference points X
[in]cell_geometryThe cell node coordinates (physical)
[in]phiTabulated basis functions at reference points X

◆ tabulate()

xt::xtensor< double, 4 > CoordinateElement::tabulate ( int  n,
const xt::xtensor< double, 2 > &  X 
) const

Compute basis values and derivatives at set of points.

Parameters
[in]nThe 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).
Returns
The basis functions (and derivatives). The shape is (derivative, number point, number of basis fn, value size).

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