Basix 0.10.0.0

Home     Installation     Demos     C++ docs     Python docs

basix::element Namespace Reference

Interfaces for creating finite elements. More...

Typedefs

template<typename T , std::size_t d>
using mdspan_t = impl::mdspan_t< T, d >
 Typedef for mdspan.
 

Enumerations

enum class  lagrange_variant {
  unset = 0 , equispaced = 1 , gll_warped = 2 , gll_isaac = 3 ,
  gll_centroid = 4 , chebyshev_warped = 5 , chebyshev_isaac = 6 , chebyshev_centroid = 7 ,
  gl_warped = 8 , gl_isaac = 9 , gl_centroid = 10 , legendre = 11 ,
  bernstein = 12
}
 Variants of a Lagrange space that can be created.
 
enum class  dpc_variant {
  unset = 0 , simplex_equispaced = 1 , simplex_gll = 2 , horizontal_equispaced = 3 ,
  horizontal_gll = 4 , diagonal_equispaced = 5 , diagonal_gll = 6 , legendre = 7
}
 
enum class  family {
  custom = 0 , P = 1 , RT = 2 , N1E = 3 ,
  BDM = 4 , N2E = 5 , CR = 6 , Regge = 7 ,
  DPC = 8 , bubble = 9 , serendipity = 10 , HHJ = 11 ,
  Hermite = 12 , iso = 13
}
 Available element families.
 

Functions

template<std::floating_point T>
FiniteElement< T > create_bdm (cell::type celltype, int degree, lagrange_variant lvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_bubble (cell::type celltype, int degree, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_cr (cell::type celltype, int degree, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_hermite (cell::type celltype, int degree, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_hhj (cell::type celltype, int degree, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_lagrange (cell::type celltype, int degree, lagrange_variant variant, bool discontinuous, std::vector< int > dof_ordering={})
 Create a Lagrange(-like) element on cell with given degree. More...
 
template<std::floating_point T>
FiniteElement< T > create_iso (cell::type celltype, int degree, lagrange_variant variant, bool discontinuous)
 Create an iso macro element on cell with given degree. More...
 
template<std::floating_point T>
FiniteElement< T > create_rtc (cell::type celltype, int degree, element::lagrange_variant lvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_nce (cell::type celltype, int degree, element::lagrange_variant lvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_nedelec (cell::type celltype, int degree, lagrange_variant lvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_nedelec2 (cell::type celltype, int degree, lagrange_variant lvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_rt (cell::type celltype, int degree, element::lagrange_variant lvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_regge (cell::type celltype, int degree, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_serendipity (cell::type celltype, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_dpc (cell::type celltype, int degree, element::dpc_variant variant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_serendipity_div (cell::type celltype, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
 
template<std::floating_point T>
FiniteElement< T > create_serendipity_curl (cell::type celltype, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
 
template<std::floating_point T>
std::tuple< std::array< std::vector< std::vector< T > >, 4 >, std::array< std::vector< std::array< std::size_t, 2 > >, 4 >, std::array< std::vector< std::vector< T > >, 4 >, std::array< std::vector< std::array< std::size_t, 4 > >, 4 > > make_discontinuous (const std::array< std::vector< mdspan_t< const T, 2 >>, 4 > &x, const std::array< std::vector< mdspan_t< const T, 4 >>, 4 > &M, std::size_t tdim, std::size_t value_size)
 

Detailed Description

Interfaces for creating finite elements.

Enumeration Type Documentation

◆ dpc_variant

Variants of a DPC (discontinuous polynomial cubical) space that can be created. DPC spaces span the same set of polynomials as Lagrange spaces on simplices but are defined on tensor product cells.

Function Documentation

◆ create_bdm()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_bdm ( cell::type  celltype,
int  degree,
lagrange_variant  lvariant,
bool  discontinuous 
)

Create BDM element

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe Lagrange variant to use when defining the element to take integral moments against
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_bubble()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_bubble ( cell::type  celltype,
int  degree,
bool  discontinuous 
)

Create a bubble element on cell with given degree

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_cr()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_cr ( cell::type  celltype,
int  degree,
bool  discontinuous 
)

Crouzeix-Raviart element

Note
degree must be 1 for Crouzeix-Raviart
Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_hermite()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_hermite ( cell::type  celltype,
int  degree,
bool  discontinuous 
)

Create a Hermite element on cell with given degree

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_hhj()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_hhj ( cell::type  celltype,
int  degree,
bool  discontinuous 
)

Create Hellan-Herrmann-Johnson element

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_lagrange()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_lagrange ( cell::type  celltype,
int  degree,
lagrange_variant  variant,
bool  discontinuous,
std::vector< int >  dof_ordering = {} 
)

Create a Lagrange(-like) element on cell with given degree.

Parameters
[in]celltypeThe element cell type
[in]degreeThe degree of the element
[in]variantThe variant of the element to be created
[in]discontinuousTrue if the is discontinuous
[in]dof_orderingDOF reordering
Returns
A finite element

◆ create_iso()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_iso ( cell::type  celltype,
int  degree,
lagrange_variant  variant,
bool  discontinuous 
)

Create an iso macro element on cell with given degree.

Parameters
[in]celltypeThe element cell type
[in]degreeThe degree of the element
[in]variantThe variant of the element to be created
[in]discontinuousTrue if the is discontinuous
Returns
A finite element

◆ create_rtc()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_rtc ( cell::type  celltype,
int  degree,
element::lagrange_variant  lvariant,
bool  discontinuous 
)

Create RTC H(div) element

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe Lagrange variant to use when defining the element to take integral moments against
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_nce()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_nce ( cell::type  celltype,
int  degree,
element::lagrange_variant  lvariant,
bool  discontinuous 
)

Create NC H(curl) element

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe Lagrange variant to use when defining the element to take integral moments against
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_nedelec()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_nedelec ( cell::type  celltype,
int  degree,
lagrange_variant  lvariant,
bool  discontinuous 
)

Create Nedelec element (first kind)

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe Lagrange variant to use when defining the element to take integral moments against
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_nedelec2()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_nedelec2 ( cell::type  celltype,
int  degree,
lagrange_variant  lvariant,
bool  discontinuous 
)

Create Nedelec element (second kind)

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe Lagrange variant to use when defining the element to take integral moments against
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_rt()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_rt ( cell::type  celltype,
int  degree,
element::lagrange_variant  lvariant,
bool  discontinuous 
)

Create Raviart-Thomas element

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe Lagrange variant to use when defining the element to take integral moments against
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_regge()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_regge ( cell::type  celltype,
int  degree,
bool  discontinuous 
)

Create Regge element

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_serendipity()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_serendipity ( cell::type  celltype,
int  degree,
element::lagrange_variant  lvariant,
element::dpc_variant  dvariant,
bool  discontinuous 
)

Create a serendipity element on cell with given degree

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe variant of the Lagrange element to be used for integral moments on the edges of the cell
[in]dvariantThe variant of the DPC element to be used for integral moments on the interior of the cell (for quads and hexes). For elements on an interval element::dpc_variant::unset can be passed in
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_dpc()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_dpc ( cell::type  celltype,
int  degree,
element::dpc_variant  variant,
bool  discontinuous 
)

Create a DPC (discontinuous polynomial cubical) element on cell with given degree.

Note
DPC elements must be discontinuous
Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]variantThe variant of the element to be created
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_serendipity_div()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_serendipity_div ( cell::type  celltype,
int  degree,
element::lagrange_variant  lvariant,
element::dpc_variant  dvariant,
bool  discontinuous 
)

Create a serendipity H(div) element on cell with given degree

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe variant of the Lagrange element to be used for integral moments
[in]dvariantThe variant of the DPC element to be used for integral moments
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ create_serendipity_curl()

template<std::floating_point T>
template FiniteElement< double > basix::element::create_serendipity_curl ( cell::type  celltype,
int  degree,
element::lagrange_variant  lvariant,
element::dpc_variant  dvariant,
bool  discontinuous 
)

Create a serendipity H(curl) element on cell with given degree

Parameters
[in]celltypeThe cell type
[in]degreeThe degree of the element
[in]lvariantThe variant of the Lagrange element to be used for integral moments
[in]dvariantThe variant of the DPC element to be used for integral moments
[in]discontinuousControls whether the element is continuous or discontinuous
Returns
A finite element

◆ make_discontinuous()

template<std::floating_point T>
std::tuple< std::array< std::vector< std::vector< T > >, 4 >, std::array< std::vector< std::array< std::size_t, 2 > >, 4 >, std::array< std::vector< std::vector< T > >, 4 >, std::array< std::vector< std::array< std::size_t, 4 > >, 4 > > basix::element::make_discontinuous ( const std::array< std::vector< mdspan_t< const T, 2 >>, 4 > &  x,
const std::array< std::vector< mdspan_t< const T, 4 >>, 4 > &  M,
std::size_t  tdim,
std::size_t  value_size 
)

Create a version of the interpolation points, interpolation matrices and entity transformation that represent a discontinuous version of the element. This discontinuous version will have the same DOFs but they will all be associated with the interior of the reference cell.

Parameters
[in]xInterpolation points. Indices are (tdim, entity index, point index, dim)
[in]MThe interpolation matrices. Indices are (tdim, entity index, dof, vs, point_index, derivative)
[in]tdimThe topological dimension of the cell the element is defined on
[in]value_sizeThe value size of the element
Returns
(xdata, xshape, Mdata, Mshape), where the x and M data are for a discontinuous version of the element (with the same shapes as x and M)