basix

Basix is a finite element definition and tabulation library.

The core of the library is written in C++, but the majority of Basix’s functionality can be used via this Python interface.

class basix.CellType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Cell type.

hexahedron = basix._basixcpp.CellType.hexahedron
interval = basix._basixcpp.CellType.interval
point = basix._basixcpp.CellType.point
prism = basix._basixcpp.CellType.prism
pyramid = basix._basixcpp.CellType.pyramid
quadrilateral = basix._basixcpp.CellType.quadrilateral
tetrahedron = basix._basixcpp.CellType.tetrahedron
triangle = basix._basixcpp.CellType.triangle
class basix.DPCVariant(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

DPC variant.

diagonal_equispaced = basix._basixcpp.DPCVariant.diagonal_equispaced
diagonal_gll = basix._basixcpp.DPCVariant.diagonal_gll
horizontal_equispaced = basix._basixcpp.DPCVariant.horizontal_equispaced
horizontal_gll = basix._basixcpp.DPCVariant.horizontal_gll
legendre = basix._basixcpp.DPCVariant.legendre
simplex_equispaced = basix._basixcpp.DPCVariant.simplex_equispaced
simplex_gll = basix._basixcpp.DPCVariant.simplex_gll
unset = basix._basixcpp.DPCVariant.unset
class basix.ElementFamily(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Element family.

BDM = basix._basixcpp.ElementFamily.BDM
CR = basix._basixcpp.ElementFamily.CR
DPC = basix._basixcpp.ElementFamily.DPC
HHJ = basix._basixcpp.ElementFamily.HHJ
Hermite = basix._basixcpp.ElementFamily.Hermite
N1E = basix._basixcpp.ElementFamily.N1E
N2E = basix._basixcpp.ElementFamily.N2E
P = basix._basixcpp.ElementFamily.P
RT = basix._basixcpp.ElementFamily.RT
Regge = basix._basixcpp.ElementFamily.Regge
bubble = basix._basixcpp.ElementFamily.bubble
custom = basix._basixcpp.ElementFamily.custom
iso = basix._basixcpp.ElementFamily.iso
serendipity = basix._basixcpp.ElementFamily.serendipity
class basix.LagrangeVariant(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Lagrange variant.

bernstein = basix._basixcpp.LagrangeVariant.bernstein
chebyshev_centroid = basix._basixcpp.LagrangeVariant.chebyshev_centroid
chebyshev_isaac = basix._basixcpp.LagrangeVariant.chebyshev_isaac
chebyshev_warped = basix._basixcpp.LagrangeVariant.chebyshev_warped
equispaced = basix._basixcpp.LagrangeVariant.equispaced
gl_centroid = basix._basixcpp.LagrangeVariant.gl_centroid
gl_isaac = basix._basixcpp.LagrangeVariant.gl_isaac
gl_warped = basix._basixcpp.LagrangeVariant.gl_warped
gll_centroid = basix._basixcpp.LagrangeVariant.gll_centroid
gll_isaac = basix._basixcpp.LagrangeVariant.gll_isaac
gll_warped = basix._basixcpp.LagrangeVariant.gll_warped
legendre = basix._basixcpp.LagrangeVariant.legendre
unset = basix._basixcpp.LagrangeVariant.unset
class basix.LatticeSimplexMethod(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Lattice simplex method.

centroid = basix._basixcpp.LatticeSimplexMethod.centroid
isaac = basix._basixcpp.LatticeSimplexMethod.isaac
none = basix._basixcpp.LatticeSimplexMethod.none
warp = basix._basixcpp.LatticeSimplexMethod.warp
class basix.LatticeType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Lattice type.

chebyshev = basix._basixcpp.LatticeType.chebyshev
equispaced = basix._basixcpp.LatticeType.equispaced
gl = basix._basixcpp.LatticeType.gl
gll = basix._basixcpp.LatticeType.gll
class basix.MapType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Map type.

L2Piola = basix._basixcpp.MapType.L2Piola
contravariantPiola = basix._basixcpp.MapType.contravariantPiola
covariantPiola = basix._basixcpp.MapType.covariantPiola
doubleContravariantPiola = basix._basixcpp.MapType.doubleContravariantPiola
doubleCovariantPiola = basix._basixcpp.MapType.doubleCovariantPiola
identity = basix._basixcpp.MapType.identity
class basix.PolynomialType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Polynomial type.

bernstein = basix._basixcpp.PolynomialType.bernstein
legendre = basix._basixcpp.PolynomialType.legendre
class basix.PolysetType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Polyset type.

macroedge = basix._basixcpp.PolysetType.macroedge
standard = basix._basixcpp.PolysetType.standard
class basix.QuadratureType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Quadrature type.

Default = basix._basixcpp.QuadratureType.Default
gauss_jacobi = basix._basixcpp.QuadratureType.gauss_jacobi
gll = basix._basixcpp.QuadratureType.gll
xiao_gimbutas = basix._basixcpp.QuadratureType.xiao_gimbutas
class basix.SobolevSpace(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Sobolev space.

H1 = basix._basixcpp.SobolevSpace.H1
H2 = basix._basixcpp.SobolevSpace.H2
H3 = basix._basixcpp.SobolevSpace.H3
HCurl = basix._basixcpp.SobolevSpace.HCurl
HDiv = basix._basixcpp.SobolevSpace.HDiv
HDivDiv = basix._basixcpp.SobolevSpace.HDivDiv
HEin = basix._basixcpp.SobolevSpace.HEin
HInf = basix._basixcpp.SobolevSpace.HInf
L2 = basix._basixcpp.SobolevSpace.L2
basix.compute_interpolation_operator(e0: FiniteElement, e1: FiniteElement) ndarray[Any, dtype[_ScalarType_co]]

Compute a matrix that represents the interpolation between two elements.

If the two elements have the same value size, this function returns the interpolation between them.

If element_from has value size 1 and element_to has value size > 1, then this function returns a matrix to interpolate from a blocked element_from (ie multiple copies of element_from) into element_to.

If element_to has value size 1 and element_from has value size > 1, then this function returns a matrix that interpolates the components of element_from into copies of element_to.

Note

If the elements have different value sizes and both are greater than 1, this function throws a runtime error

In order to interpolate functions between finite element spaces on arbitrary cells, the functions must be pulled back to the reference element (this pull back includes applying DOF transformations). The matrix that this function returns can then be applied, then the result pushed forward to the cell. If element_from and element_to have the same map type, then only the DOF transformations need to be applied, as the pull back and push forward cancel each other out.

Parameters:
  • e0 – The element to interpolate from

  • e1 – The element to interpolate to

Returns:

Matrix operator that maps the ‘from’ degrees-of-freedom to the ‘to’ degrees-of-freedom. Shape is (ndofs(element_to), ndofs(element_from))

basix.create_custom_element(cell_type: ~basix.cell.CellType, value_shape: tuple[int, ...], wcoeffs: ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.floating]], x: list[list[~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.floating]]]], M: list[list[~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.floating]]]], interpolation_nderivs: int, map_type: ~basix.maps.MapType, sobolev_space: ~basix.sobolev_spaces.SobolevSpace, discontinuous: bool, embedded_subdegree: int, embedded_superdegree: int, poly_type: ~basix.polynomials.PolysetType, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>) FiniteElement

Create a custom finite element.

Parameters:
  • cell_type – Element cell type.

  • value_shape – Value shape of the element.

  • wcoeffs – Matrices for the k-th value index containing the expansion coefficients defining a polynomial basis spanning the polynomial space for this element. Shape is (dim(finite element polyset), dim(Legendre polynomials)).

  • x – Interpolation points. Indices are (tdim, entity index, point index, dim).

  • M – Interpolation matrices. Indices are (tdim, entity index, dof, vs, point_index, derivative).

  • interpolation_nderivs – Number of derivatives that need to be used during interpolation.

  • map_type – Type of map to be used to map values from the reference to a physical cell.

  • sobolev_space – Underlying Sobolev space for the element.

  • discontinuous – If True create the discontinuous version of the element.

  • embedded_subdegree – Highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element.

  • embedded_superdegree – Degree of a polynomial in this element’s polyset.

  • poly_type – Type of polyset to use for this element.

  • dtype – Element scalar type.

Returns:

A custom finite element.

basix.create_element(family: ~basix.finite_element.ElementFamily, celltype: ~basix.cell.CellType, degree: int, lagrange_variant: ~basix.finite_element.LagrangeVariant = LagrangeVariant.unset, dpc_variant: ~basix.finite_element.DPCVariant = DPCVariant.unset, discontinuous: bool = False, dof_ordering: list[int] | None = None, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>) FiniteElement

Create a finite element.

Parameters:
  • family – Finite element family.

  • celltype – Reference cell type that the element is defined on.

  • degree – Polynomial degree of the element.

  • lagrange_variant – Lagrange variant type.

  • dpc_variant – DPC variant type.

  • discontinuous – If True element is discontinuous. The discontinuous element will have the same DOFs as a continuous element, but the DOFs will all be associated with the interior of the cell.

  • dof_ordering – Ordering of dofs for ElementDofLayout.

  • dtype – Element scalar type.

Returns:

A finite element.

basix.create_lattice(celltype: CellType, n: int, ltype: LatticeType, exterior: bool, method: LatticeSimplexMethod = LatticeSimplexMethod.none) ndarray[Any, dtype[_ScalarType_co]]

Create a lattice of points on a reference cell.

Parameters:
  • celltype – Cell type.

  • n – The size in each direction. There will be n+1 points along each edge of the cell.

  • ltype – Lattice type.

  • exterior – If True, the points on the edges will be included.

  • method – The simplex method used to generate points on simplices.

Returns:

Lattice points

basix.create_tp_element(family: ~basix.finite_element.ElementFamily, celltype: ~basix.cell.CellType, degree: int, lagrange_variant: ~basix.finite_element.LagrangeVariant = LagrangeVariant.unset, dpc_variant: ~basix.finite_element.DPCVariant = DPCVariant.unset, discontinuous: bool = False, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>) FiniteElement

Create a finite element with tensor product ordering.

Parameters:
  • family – Finite element family.

  • celltype – Reference cell type that the element is defined on

  • degree – Polynomial degree of the element.

  • lagrange_variant – Lagrange variant type.

  • dpc_variant – DPC variant type.

  • discontinuous – If True element is discontinuous. The discontinuous element will have the same DOFs as a continuous element, but the DOFs will all be associated with the interior of the cell.

  • dtype – Element scalar type.

Returns:

A finite element.

basix.geometry(celltype: CellType) ndarray[Any, dtype[_ScalarType_co]]

Cell geometry.

Parameters:

celltype – Cell type.

Returns:

Vertices of the cell.

basix.index(p: int, q: int | None = None, r: int | None = None) int

Compute the indexing in a 1D, 2D or 3D simplex.

Parameters:
  • p – Index in x.

  • q – Index in y.

  • r – Index in z.

Returns:

Index in a flattened 1D array.

basix.make_quadrature(cell: CellType, degree: int, rule: QuadratureType = QuadratureType.Default, polyset_type: PolysetType = PolysetType.standard) tuple[ndarray[Any, dtype[float64]], ndarray[Any, dtype[float64]]]

Create a quadrature rule.

Parameters:
  • cell – Cell type.

  • degree – Maximum polynomial degree that will be integrated exactly.

  • rule – Quadrature rule.

  • polyset_type – Type of polynomial that will be integrated exactly.

Returns:

Quadrature points and weights.

basix.polyset_restriction(ptype: PolysetType, cell: CellType, restriction_cell: CellType) PolysetType

Get the polyset type that represents the restrictions of a type on a subentity.

Parameters:
  • ptype – The polynomial type

  • cell – The cell type

  • restriction_cell – The cell type if the subentity

Returns:

The restricted polyset type

basix.polyset_superset(cell: CellType, type1: PolysetType, type2: PolysetType) PolysetType

Get the polyset type that is a superset of two types on the given cell.

Parameters:
  • cell – The cell type

  • type1 – The first polyset type

  • type2 – The second polyset type

Returns:

The superset type

basix.tabulate_polynomials(ptype: PolynomialType, celltype: CellType, degree: int, pts: ndarray[Any, dtype[_ScalarType_co]]) ndarray[Any, dtype[_ScalarType_co]]

Tabulate a set of polynomials on a reference cell.

Parameters:
  • ptype – The polynomial type

  • celltype – The cell type

  • degree – The polynomial degree

  • pts – The points

Returns:

Tabulated polynomials

basix.topology(celltype: CellType) list[list[list[int]]]

Cell topology.

Parameters:

celltype – Cell type.

Returns:

Vertex indices for each sub-entity of the cell.

Modules

basix.cell

Functions to get cell geometry information and manipulate cell types.

basix.finite_element

Functions for creating finite elements.

basix.interpolation

Interpolation.

basix.lattice

Functions to manipulate lattice types.

basix.maps

Maps.

basix.numba_helpers

Helper functions for writing DOLFINx custom kernels using Numba.

basix.polynomials

Functions for working with polynomials.

basix.quadrature

Functions to manipulate quadrature types.

basix.sobolev_spaces

Functions for handling Sobolev spaces.

basix.ufl

Functions to directly wrap Basix elements in UFL.

basix.utils

Utility funcitons.