Home Installation C++ docs Python docs
Functions
basix::moments Namespace Reference

Functions

xt::xtensor< double, 3 > create_dot_moment_dof_transformations (const FiniteElement &moment_space)
 
xt::xtensor< double, 3 > create_moment_dof_transformations (const FiniteElement &moment_space)
 
xt::xtensor< double, 3 > create_normal_moment_dof_transformations (const FiniteElement &moment_space)
 
xt::xtensor< double, 3 > create_tangent_moment_dof_transformations (const FiniteElement &moment_space)
 
std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > make_integral_moments (const FiniteElement &moment_space, cell::type celltype, std::size_t value_size, int q_deg)
 
std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > make_dot_integral_moments (const FiniteElement &V, cell::type celltype, std::size_t value_size, int q_deg)
 
std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > make_tangent_integral_moments (const FiniteElement &V, cell::type celltype, std::size_t value_size, int q_deg)
 
std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > make_normal_integral_moments (const FiniteElement &V, cell::type celltype, std::size_t value_size, int q_deg)
 

Detailed Description

Integral moments

These functions generate dual set matrices for integral moments against spaces on a subentity of the cell

Function Documentation

◆ create_dot_moment_dof_transformations()

xt::xtensor< double, 3 > basix::moments::create_dot_moment_dof_transformations ( const FiniteElement moment_space)

Create the dof transformations for the DOFs defined using a dot integral moment.

A dot integral moment is defined by

\[l_i(v) = \int v\cdot\phi_i,\]

where \(\phi_i\) is a basis function in the moment space, and \(v\) and \(\phi_i\) are either both scalars or are vectors of the same size.

If the moment space is an interval, this returns one matrix representing the reversal of the interval. If the moment space is a face, this returns two matrices: one representing a rotation, the other a reflection.

These matrices are computed by calculation the interpolation coefficients of a rotated/reflected basis into the original basis.

Parameters
[in]moment_spaceThe finite element space that the integral moment is taken against
Returns
A list of dof transformations

◆ create_moment_dof_transformations()

xt::xtensor< double, 3 > basix::moments::create_moment_dof_transformations ( const FiniteElement moment_space)

Create the DOF transformations for the DOFs defined using an integral moment.

An integral moment is defined by

\[l_{i,j}(v) = \int v\cdot e_j\phi_i,\]

where \(\phi_i\) is a basis function in the moment space, \(e_j\) is a coordinate direction (of the cell sub-entity the moment is taken on), \(v\) is a vector, and \(\phi_i\) is a scalar.

This will combine multiple copies of the result of create_dot_moment_dof_transformations to give the transformations for integral moments of each vector component against the moment space.

Parameters
[in]moment_spaceThe finite element space that the integral moment is taken against
Returns
A list of dof transformations

◆ create_normal_moment_dof_transformations()

xt::xtensor< double, 3 > basix::moments::create_normal_moment_dof_transformations ( const FiniteElement moment_space)

Create the dof transformations for the DOFs defined using a normal integral moment.

A normal integral moment is defined by

\[l_{i,j}(v) = \int v\cdot n\phi_i,\]

where \(\phi_i\) is a basis function in the moment space, \(n\) is normal to the cell sub-entity, \(v\) is a vector, and \(\phi_i\) is a scalar.

This does the same as create_dot_moment_dof_transformations with some additional factors of -1 to account for the changing of the normal direction when the entity is reflected.

Parameters
[in]moment_spaceThe finite element space that the integral moment is taken against
Returns
A list of dof transformations

◆ create_tangent_moment_dof_transformations()

xt::xtensor< double, 3 > basix::moments::create_tangent_moment_dof_transformations ( const FiniteElement moment_space)

Create the dof transformations for the DOFs defined using a tangential integral moment.

A tangential integral moment is defined by

\[l_{i,j}(v) = \int v\cdot t\phi_i,\]

where \(\phi_i\) is a basis function in the moment space, \(t\) is tangential to the edge, \(v\) is a vector, and \(\phi_i\) is a scalar.

This does the same as create_dot_moment_dof_transformations with some additional factors of -1 to account for the changing of the tangent direction when the edge is reflected.

Parameters
[in]moment_spaceThe finite element space that the integral moment is taken against
Returns
A list of dof transformations

◆ make_dot_integral_moments()

std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > basix::moments::make_dot_integral_moments ( const FiniteElement V,
cell::type  celltype,
std::size_t  value_size,
int  q_deg 
)

Make interpolation points and weights for dot product integral moments

These will represent the integral of each function in the moment space over each sub entity of the moment space's cell type in a cell with the given type. For example, if the input cell type is a triangle and the moment space is a P1 space on an edge, this will perform two integrals for each of the 3 edges of the triangle.

Todo:
Clarify what happens value size of the moment space is less than value_size.
Parameters
VThe space to compute the integral moments against
celltypeThe cell type of the cell on which the space is being defined
value_sizeThe value size of the space being defined
q_degThe quadrature degree used for the integrals

◆ make_integral_moments()

std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > basix::moments::make_integral_moments ( const FiniteElement moment_space,
cell::type  celltype,
std::size_t  value_size,
int  q_deg 
)

Make interpolation points and weights for simple integral moments

These will represent the integral of each function in the moment space over each sub entity of the moment space's cell type in a cell with the given type. For example, if the input cell type is a triangle, and the moment space is a P1 space on an edge, this will perform two integrals for each of the 3 edges of the triangle.

Parameters
moment_spaceThe space to compute the integral moments against
celltypeThe cell type of the cell on which the space is being defined
value_sizeThe value size of the space being defined
q_degThe quadrature degree used for the integrals

◆ make_normal_integral_moments()

std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > basix::moments::make_normal_integral_moments ( const FiniteElement V,
cell::type  celltype,
std::size_t  value_size,
int  q_deg 
)

Compute interpolation points and weights for normal integral moments

These can only be used when the moment space is defined on facets of the cell

Parameters
[in]VThe space to compute the integral moments against
[in]celltypeThe cell type of the cell on which the space is being defined
[in]value_sizeThe value size of the space being defined
[in]q_degThe quadrature degree used for the integrals
Returns
(interpolation points, interpolation matrix)

◆ make_tangent_integral_moments()

std::pair< std::vector< xt::xtensor< double, 2 > >, std::vector< xt::xtensor< double, 3 > > > basix::moments::make_tangent_integral_moments ( const FiniteElement V,
cell::type  celltype,
std::size_t  value_size,
int  q_deg 
)

Make interpolation points and weights for tangent integral moments

These can only be used when the moment space is defined on edges of the cell

Parameters
VThe space to compute the integral moments against
celltypeThe cell type of the cell on which the space is being defined
value_sizeThe value size of the space being defined the space
q_degThe quadrature degree used for the integrals