Tabulation Project basix
Functions
basix::moments Namespace Reference

Functions

Eigen::MatrixXd make_integral_moments (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > make_integral_moments_interpolation (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
Eigen::MatrixXd make_dot_integral_moments (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > make_dot_integral_moments_interpolation (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
Eigen::MatrixXd make_tangent_integral_moments (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > make_tangent_integral_moments_interpolation (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
Eigen::MatrixXd make_normal_integral_moments (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const int q_deg)
 
std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > make_normal_integral_moments_interpolation (const FiniteElement &moment_space, const cell::type celltype, const int value_size, const int poly_deg, const 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

◆ make_dot_integral_moments()

Eigen::MatrixXd basix::moments::make_dot_integral_moments ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make dot product integral moments

This will compute 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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_dot_integral_moments_interpolation()

std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > basix::moments::make_dot_integral_moments_interpolation ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make interpolation points and weights for dot product integral moments

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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_integral_moments()

Eigen::MatrixXd basix::moments::make_integral_moments ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make simple integral moments

This will compute 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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_integral_moments_interpolation()

std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > basix::moments::make_integral_moments_interpolation ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make interpolation points and weights for simple integral moments

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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_normal_integral_moments()

Eigen::MatrixXd basix::moments::make_normal_integral_moments ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make normal integral moments

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

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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_normal_integral_moments_interpolation()

std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > basix::moments::make_normal_integral_moments_interpolation ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make 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
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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_tangent_integral_moments()

Eigen::MatrixXd basix::moments::make_tangent_integral_moments ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const int  q_deg 
)

Make tangential integral moments

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

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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals

◆ make_tangent_integral_moments_interpolation()

std::pair< Eigen::ArrayXXd, Eigen::MatrixXd > basix::moments::make_tangent_integral_moments_interpolation ( const FiniteElement moment_space,
const cell::type  celltype,
const int  value_size,
const int  poly_deg,
const 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
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
poly_degThe polynomial degree of the poly set that defines the space
q_degThe quadrature degree used for the integrals