Basix 0.9.0

Home     Installation     Demos     C++ docs     Python docs

basix::moments Namespace Reference

Functions to create integral moment DOFs. More...

Functions

template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > make_integral_moments (const FiniteElement< T > &moment_space, cell::type celltype, polyset::type ptype, std::size_t value_size, int q_deg)
 Make interpolation points and weights for simple integral moments. More...
 
template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > make_dot_integral_moments (const FiniteElement< T > &V, cell::type celltype, polyset::type ptype, std::size_t value_size, int q_deg)
 Make interpolation points and weights for dot product integral moments. More...
 
template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > make_tangent_integral_moments (const FiniteElement< T > &V, cell::type celltype, polyset::type ptype, std::size_t value_size, int q_deg)
 Make interpolation points and weights for tangent integral moments. More...
 
template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > make_normal_integral_moments (const FiniteElement< T > &V, cell::type celltype, polyset::type ptype, std::size_t value_size, int q_deg)
 Compute interpolation points and weights for normal integral moments. More...
 

Detailed Description

Functions to create integral moment DOFs.

Function Documentation

◆ make_integral_moments()

template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > basix::moments::make_integral_moments ( const FiniteElement< T > &  moment_space,
cell::type  celltype,
polyset::type  ptype,
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
ptypeThe polyset type of the element this moment is being used to define
value_sizeThe value size of the space being defined
q_degThe quadrature degree used for the integrals
Returns
(interpolation points, interpolation matrix). The indices of the interpolation points are (number of entities, npoints, gdim). The indices on the interpolation matrix are (number of entities, ndofs, value_size, npoints, derivative)

◆ make_dot_integral_moments()

template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > basix::moments::make_dot_integral_moments ( const FiniteElement< T > &  V,
cell::type  celltype,
polyset::type  ptype,
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
ptypeThe polyset type of the element this moment is being used to define
value_sizeThe value size of the space being defined
q_degThe quadrature degree used for the integrals
Returns
(interpolation points, interpolation shape, interpolation matrix, interpolation shape). The indices of the interpolation points are (number of entities, npoints, gdim). The indices on the interpolation matrix are (number of entities, ndofs, value_size, npoints, derivative)

◆ make_tangent_integral_moments()

template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > basix::moments::make_tangent_integral_moments ( const FiniteElement< T > &  V,
cell::type  celltype,
polyset::type  ptype,
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
ptypeThe polyset type of the element this moment is being used to define
value_sizeThe value size of the space being defined the space
q_degThe quadrature degree used for the integrals
Returns
(interpolation points, interpolation shape, interpolation matrix, interpolation shape). The indices of the interpolation points are (number of entities, npoints, gdim). The indices on the interpolation matrix are (number of entities, ndofs, value_size, npoints, derivative)

◆ make_normal_integral_moments()

template<std::floating_point T>
std::tuple< std::vector< std::vector< T > >, std::array< std::size_t, 2 >, std::vector< std::vector< T > >, std::array< std::size_t, 4 > > basix::moments::make_normal_integral_moments ( const FiniteElement< T > &  V,
cell::type  celltype,
polyset::type  ptype,
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
ptypeThe polyset type of the element this moment is being used to define
[in]value_sizeThe value size of the space being defined
[in]q_degThe quadrature degree used for the integrals
Returns
(interpolation points, interpolation shape, interpolation matrix, interpolation shape). The indices of the interpolation points are (number of entities, npoints, gdim). The indices on the interpolation matrix are (number of entities, ndofs, value_size, npoints, derivative)