Home Installation C++ docs Python docs
Public Member Functions | Public Attributes | List of all members
basix::FiniteElement Class Reference

#include <finite-element.h>

Public Member Functions

 FiniteElement (element::family family, cell::type cell_type, int degree, const std::vector< std::size_t > &value_shape, const xt::xtensor< double, 3 > &coeffs, const std::map< cell::type, xt::xtensor< double, 3 >> &entity_transformations, const std::array< std::vector< xt::xtensor< double, 2 >>, 4 > &x, const std::array< std::vector< xt::xtensor< double, 3 >>, 4 > &M, maps::type map_type=maps::type::identity)
 
 FiniteElement (const FiniteElement &element)=default
 Copy constructor.
 
 FiniteElement (FiniteElement &&element)=default
 Move constructor.
 
 ~FiniteElement ()=default
 Destructor.
 
FiniteElementoperator= (const FiniteElement &element)=default
 Assignment operator.
 
FiniteElementoperator= (FiniteElement &&element)=default
 Move assignment operator.
 
xt::xtensor< double, 4 > tabulate (int nd, const xt::xarray< double > &x) const
 
void tabulate (int nd, const xt::xarray< double > &x, xt::xtensor< double, 4 > &basis_data) const
 
cell::type cell_type () const
 
int degree () const
 
int value_size () const
 
const std::vector< int > & value_shape () const
 
int dim () const
 
element::family family () const
 
maps::type mapping_type () const
 
bool dof_transformations_are_permutations () const
 
bool dof_transformations_are_identity () const
 
xt::xtensor< double, 3 > map_push_forward (const xt::xtensor< double, 3 > &U, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K) const
 
template<typename O , typename P , typename Q , typename S , typename T >
void map_push_forward_m (const O &U, const P &J, const Q &detJ, const S &K, T &&u) const
 
xt::xtensor< double, 3 > map_pull_back (const xt::xtensor< double, 3 > &u, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K) const
 
template<typename O , typename P , typename Q , typename S , typename T >
void map_pull_back_m (const O &u, const P &J, const Q &detJ, const S &K, T &&U) const
 
const std::vector< std::vector< int > > & num_entity_dofs () const
 
const std::vector< std::vector< int > > & num_entity_closure_dofs () const
 
const std::vector< std::vector< std::set< int > > > & entity_dofs () const
 
const std::vector< std::vector< std::set< int > > > & entity_closure_dofs () const
 
xt::xtensor< double, 3 > base_transformations () const
 
std::map< cell::type, xt::xtensor< double, 3 > > entity_transformations () const
 Return the entity dof transformation matricess.
 
void permute_dofs (const xtl::span< std::int32_t > &dofs, std::uint32_t cell_info) const
 
void unpermute_dofs (const xtl::span< std::int32_t > &dofs, std::uint32_t cell_info) const
 
template<typename T >
void apply_dof_transformation (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_transpose_dof_transformation (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_transpose_dof_transformation (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_dof_transformation (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_dof_transformation_to_transpose (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_transpose_dof_transformation_to_transpose (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_transpose_dof_transformation_to_transpose (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_dof_transformation_to_transpose (const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
 
const xt::xtensor< double, 2 > & points () const
 
int num_points () const
 Return the number of interpolation points.
 
const xt::xtensor< double, 2 > & interpolation_matrix () const
 

Public Attributes

maps::type map_type
 Element map type.
 

Detailed Description

Finite Element The basis is stored as a set of coefficients, which are applied to the underlying expansion set for that cell type, when tabulating.

Constructor & Destructor Documentation

◆ FiniteElement()

FiniteElement::FiniteElement ( element::family  family,
cell::type  cell_type,
int  degree,
const std::vector< std::size_t > &  value_shape,
const xt::xtensor< double, 3 > &  coeffs,
const std::map< cell::type, xt::xtensor< double, 3 >> &  entity_transformations,
const std::array< std::vector< xt::xtensor< double, 2 >>, 4 > &  x,
const std::array< std::vector< xt::xtensor< double, 3 >>, 4 > &  M,
maps::type  map_type = maps::type::identity 
)
Todo:
Document A finite element
Parameters
[in]family
[in]cell_type
[in]degree
[in]value_shape
[in]coeffsExpansion coefficients. The shape is (num_dofs, value_size, basis_dim)
[in]entity_transformationsEntity transformations
[in]xInterpolation points. Shape is (tdim, entity index, point index, dim)
[in]MThe interpolation matrices. Indices are (tdim, entity index, dof, vs, point_index)
[in]map_type

Member Function Documentation

◆ apply_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_dof_transformation ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply DOF transformations to some data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_dof_transformation_to_transpose ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply DOF transformations to some transposed data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_inverse_dof_transformation ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse DOF transformations to some data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_inverse_dof_transformation_to_transpose ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse DOF transformations to some transposed data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_transpose_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_inverse_transpose_dof_transformation ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse transpose DOF transformations to some data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_transpose_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_inverse_transpose_dof_transformation_to_transpose ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse transpose DOF transformations to some transposed data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_transpose_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_transpose_dof_transformation ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply transpose DOF transformations to some data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_transpose_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_transpose_dof_transformation_to_transpose ( const xtl::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply transpose DOF transformations to some transposed data

Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ base_transformations()

xt::xtensor< double, 3 > FiniteElement::base_transformations ( ) const

Get the base transformations The base transformations represent the effect of rotating or reflecting a subentity of the cell on the numbering and orientation of the DOFs. This returns a list of matrices with one matrix for each subentity permutation in the following order: Reversing edge 0, reversing edge 1, ... Rotate face 0, reflect face 0, rotate face 1, reflect face 1, ...

Example: Order 3 Lagrange on a triangle

This space has 10 dofs arranged like:

2
|\
6 4
| \
5 9 3
| \
0-7-8-1

For this element, the base transformations are: [Matrix swapping 3 and 4, Matrix swapping 5 and 6, Matrix swapping 7 and 8] The first row shows the effect of reversing the diagonal edge. The second row shows the effect of reversing the vertical edge. The third row shows the effect of reversing the horizontal edge.

Example: Order 1 Raviart-Thomas on a triangle

This space has 3 dofs arranged like:

|\
| \
| \
<-1 0
| / \
| L ^ \
| | \
---2---

These DOFs are integrals of normal components over the edges: DOFs 0 and 2 are oriented inward, DOF 1 is oriented outwards. For this element, the base transformation matrices are:

0: [[-1, 0, 0],
[ 0, 1, 0],
[ 0, 0, 1]]
1: [[1, 0, 0],
[0, -1, 0],
[0, 0, 1]]
2: [[1, 0, 0],
[0, 1, 0],
[0, 0, -1]]

The first matrix reverses DOF 0 (as this is on the first edge). The second matrix reverses DOF 1 (as this is on the second edge). The third matrix reverses DOF 2 (as this is on the third edge).

Example: DOFs on the face of Order 2 Nedelec first kind on a tetrahedron

On a face of this tetrahedron, this space has two face tangent DOFs:

|\ |\
| \ | \
| \ | ^\
| \ | | \
| 0->\ | 1 \
| \ | \
------ ------

For these DOFs, the subblocks of the base transformation matrices are:

rotation: [[-1, 1],
[ 1, 0]]
reflection: [[0, 1],
[1, 0]]

◆ cell_type()

cell::type FiniteElement::cell_type ( ) const

Get the element cell type

Returns
The cell type

◆ degree()

int FiniteElement::degree ( ) const

Get the element polynomial degree

Returns
Polynomial degree

◆ dim()

int FiniteElement::dim ( ) const

Dimension of the finite element space (number of degrees of freedom for the element)

Returns
Number of degrees of freedom

◆ dof_transformations_are_identity()

bool FiniteElement::dof_transformations_are_identity ( ) const

Indicates whether the dof transformations are all the identity

Returns
True or False

◆ dof_transformations_are_permutations()

bool FiniteElement::dof_transformations_are_permutations ( ) const

Indicates whether the dof transformations are all permutations

Returns
True or False

◆ entity_closure_dofs()

const std::vector< std::vector< std::set< int > > > & FiniteElement::entity_closure_dofs ( ) const

Get the dofs on the closure of each topological entity: (vertices, edges, faces, cell) in that order. For example, Lagrange degree 2 on a triangle has vertices: [[0], [1], [2]], edges: [[1, 2, 3], [0, 2, 4], [0, 1, 5]], cell: [[0, 1, 2, 3, 4, 5]]

Returns
Dofs associated with the closre of an entity of a given topological dimension. The shape is (tdim + 1, num_entities, num_dofs).

◆ entity_dofs()

const std::vector< std::vector< std::set< int > > > & FiniteElement::entity_dofs ( ) const

Get the dofs on each topological entity: (vertices, edges, faces, cell) in that order. For example, Lagrange degree 2 on a triangle has vertices: [[0], [1], [2]], edges: [[3], [4], [5]], cell: [[]]

Returns
Dofs associated with an entity of a given topological dimension. The shape is (tdim + 1, num_entities, num_dofs).

◆ family()

element::family FiniteElement::family ( ) const

Get the finite element family

Returns
The family

◆ interpolation_matrix()

const xt::xtensor< double, 2 > & FiniteElement::interpolation_matrix ( ) const

Return a matrix of weights interpolation To interpolate a function in this finite element, the functions should be evaluated at each point given by FiniteElement::points(). These function values should then be multiplied by the weight matrix to give the coefficients of the interpolated function.

◆ map_pull_back()

xt::xtensor< double, 3 > FiniteElement::map_pull_back ( const xt::xtensor< double, 3 > &  u,
const xt::xtensor< double, 3 > &  J,
const xtl::span< const double > &  detJ,
const xt::xtensor< double, 3 > &  K 
) const

Map function values from a physical cell to the reference

Parameters
[in]uThe function values on the cell
[in]JThe Jacobian of the mapping
[in]detJThe determinant of the Jacobian of the mapping
[in]KThe inverse of the Jacobian of the mapping
Returns
The function values on the reference

◆ map_pull_back_m()

template<typename O , typename P , typename Q , typename S , typename T >
void basix::FiniteElement::map_pull_back_m ( const O &  u,
const P &  J,
const Q &  detJ,
const S &  K,
T &&  U 
) const
inline

Map function values from a physical cell back to to the reference

Parameters
[in]uData defined on the physical element. It must have dimension 3. The first index is for the geometric/map data, the second is the point index for points that share map data, and the third index is (vector) component, e.g. u[i,:,:] are points that are mapped by J[i,:,:].
[in]JThe Jacobians. It must have dimension 3. The first index is for the ith Jacobian, i.e. J[i,:,:] is the ith Jacobian.
[in]detJThe determinant of J. detJ[i] is equal to det(J[i,:,:]). It must have dimension 1.
[in]KThe inverse of J, K[i,:,:] = J[i,:,:]^-1. It must have dimension 3.
[out]UThe input u mapped to the reference element. It must have dimension 3.

◆ map_push_forward()

xt::xtensor< double, 3 > FiniteElement::map_push_forward ( const xt::xtensor< double, 3 > &  U,
const xt::xtensor< double, 3 > &  J,
const xtl::span< const double > &  detJ,
const xt::xtensor< double, 3 > &  K 
) const

Map function values from the reference to a physical cell. This function can perform the mapping for multiple points, grouped by points that share a common Jacobian.

Parameters
UThe function values on the reference. The indices are [Jacobian index, point index, components].
JThe Jacobian of the mapping. The indices are [Jacobian index, J_i, J_j].
detJThe determinant of the Jacobian of the mapping. It has length J.shape(0)
KThe inverse of the Jacobian of the mapping. The indices are [Jacobian index, K_i, K_j].
Returns
The function values on the cell. The indices are [Jacobian index, point index, components].

◆ map_push_forward_m()

template<typename O , typename P , typename Q , typename S , typename T >
void basix::FiniteElement::map_push_forward_m ( const O &  U,
const P &  J,
const Q &  detJ,
const S &  K,
T &&  u 
) const
inline

Direct to memory push forward

Parameters
[in]UData defined on the reference element. It must have dimension 3. The first index is for the geometric/map data, the second is the point index for points that share map data, and the third index is (vector) component, e.g. u[i,:,:] are points that are mapped by J[i,:,:].
[in]JThe Jacobians. It must have dimension 3. The first index is for the ith Jacobian, i.e. J[i,:,:] is the ith Jacobian.
[in]detJThe determinant of J. detJ[i] is equal to det(J[i,:,:]). It must have dimension 1.
[in]KThe inverse of J, K[i,:,:] = J[i,:,:]^-1. It must have dimension 3.
[out]uThe input U mapped to the physical. It must have dimension 3.

◆ mapping_type()

maps::type FiniteElement::mapping_type ( ) const

Get the mapping type used for this element

Returns
The mapping

◆ num_entity_closure_dofs()

const std::vector< std::vector< int > > & FiniteElement::num_entity_closure_dofs ( ) const

Get the number of dofs on the closure of each topological entity: (vertices, edges, faces, cell) in that order. For example, Lagrange degree 2 on a triangle has vertices: [1, 1, 1], edges: [3, 3, 3], cell: [6]

Returns
Number of dofs associated with the closure of an entity of a given topological dimension. The shape is (tdim + 1, num_entities).

◆ num_entity_dofs()

const std::vector< std::vector< int > > & FiniteElement::num_entity_dofs ( ) const

Get the number of dofs on each topological entity: (vertices, edges, faces, cell) in that order. For example, Lagrange degree 2 on a triangle has vertices: [1, 1, 1], edges: [1, 1, 1], cell: [0] The sum of the entity dofs must match the total number of dofs reported by FiniteElement::dim,

const std::vector<std::vector<int>>& dofs = e.entity_dofs();
int num_dofs0 = dofs[1][3]; // Number of dofs associated with edge 3
int num_dofs1 = dofs[2][0]; // Number of dofs associated with face 0
Returns
Number of dofs associated with an entity of a given topological dimension. The shape is (tdim + 1, num_entities).

◆ permute_dofs()

void FiniteElement::permute_dofs ( const xtl::span< std::int32_t > &  dofs,
std::uint32_t  cell_info 
) const

Permute the dof numbering on a cell

Parameters
[in,out]dofsThe dof numbering for the cell
cell_infoThe permutation info for the cell

◆ points()

const xt::xtensor< double, 2 > & FiniteElement::points ( ) const

Return the interpolation points, i.e. the coordinates on the reference element where a function need to be evaluated in order to interpolate it in the finite element space.

Returns
Array of coordinate with shape (num_points, tdim)

◆ tabulate() [1/2]

xt::xtensor< double, 4 > FiniteElement::tabulate ( int  nd,
const xt::xarray< double > &  x 
) const

Compute 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).
Returns
The basis functions (and derivatives). The shape is (derivative, point, basis fn index, value index).
  • The first index is the derivative, with higher derivatives are stored in triangular (2D) or tetrahedral (3D) ordering, i.e. for the (x,y) derivatives in 2D: (0,0), (1,0), (0,1), (2,0), (1,1), (0,2), (3,0)... The function basix::idx can be used to find the appropriate derivative.
  • The second index is the point index
  • The third index is the basis function index
  • The fourth index is the basis function component. Its has size one for scalar basis functions.

◆ tabulate() [2/2]

void FiniteElement::tabulate ( int  nd,
const xt::xarray< double > &  x,
xt::xtensor< double, 4 > &  basis_data 
) const

Direct to memory block tabulation

Parameters
ndNumber of derivatives
xPoints
basis_dataMemory location to fill

◆ unpermute_dofs()

void FiniteElement::unpermute_dofs ( const xtl::span< std::int32_t > &  dofs,
std::uint32_t  cell_info 
) const

Unpermute the dof numbering on a cell

Parameters
[in,out]dofsThe dof numbering for the cell
cell_infoThe permutation info for the cell

◆ value_shape()

const std::vector< int > & FiniteElement::value_shape ( ) const

Get the element value tensor shape, e.g. returning [1] for scalars.

Returns
Value shape

◆ value_size()

int FiniteElement::value_size ( ) const

Get the element value size This is just a convenience function returning product(value_shape)

Returns
Value size

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