basix

Basix is a finite element definition and tabulation runtime 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(self: basix._basixcpp.CellType, value: int)

Bases: pybind11_object

Members:

point

interval

triangle

tetrahedron

quadrilateral

hexahedron

prism

pyramid

hexahedron = <CellType.hexahedron: 5>
interval = <CellType.interval: 1>
property name
point = <CellType.point: 0>
prism = <CellType.prism: 6>
pyramid = <CellType.pyramid: 7>
quadrilateral = <CellType.quadrilateral: 4>
tetrahedron = <CellType.tetrahedron: 3>
triangle = <CellType.triangle: 2>
property value
class basix.DPCVariant(self: basix._basixcpp.DPCVariant, value: int)

Bases: pybind11_object

Members:

unset

simplex_equispaced

simplex_gll

horizontal_equispaced

horizontal_gll

diagonal_equispaced

diagonal_gll

legendre

diagonal_equispaced = <DPCVariant.diagonal_equispaced: 4>
diagonal_gll = <DPCVariant.diagonal_gll: 5>
horizontal_equispaced = <DPCVariant.horizontal_equispaced: 2>
horizontal_gll = <DPCVariant.horizontal_gll: 3>
legendre = <DPCVariant.legendre: 6>
property name
simplex_equispaced = <DPCVariant.simplex_equispaced: 0>
simplex_gll = <DPCVariant.simplex_gll: 1>
unset = <DPCVariant.unset: -1>
property value
class basix.ElementFamily(self: basix._basixcpp.ElementFamily, value: int)

Bases: pybind11_object

Members:

custom

P

BDM

RT

N1E

N2E

Regge

HHJ

bubble

serendipity

DPC

CR

Hermite

iso

BDM = <ElementFamily.BDM: 4>
CR = <ElementFamily.CR: 6>
DPC = <ElementFamily.DPC: 8>
HHJ = <ElementFamily.HHJ: 11>
Hermite = <ElementFamily.Hermite: 12>
N1E = <ElementFamily.N1E: 3>
N2E = <ElementFamily.N2E: 5>
P = <ElementFamily.P: 1>
RT = <ElementFamily.RT: 2>
Regge = <ElementFamily.Regge: 7>
bubble = <ElementFamily.bubble: 9>
custom = <ElementFamily.custom: 0>
iso = <ElementFamily.iso: 13>
property name
serendipity = <ElementFamily.serendipity: 10>
property value
class basix.LagrangeVariant(self: basix._basixcpp.LagrangeVariant, value: int)

Bases: pybind11_object

Members:

unset

equispaced

gll_warped

gll_isaac

gll_centroid

chebyshev_warped

chebyshev_isaac

chebyshev_centroid

gl_warped

gl_isaac

gl_centroid

legendre

bernstein

vtk

bernstein = <LagrangeVariant.bernstein: 11>
chebyshev_centroid = <LagrangeVariant.chebyshev_centroid: 6>
chebyshev_isaac = <LagrangeVariant.chebyshev_isaac: 5>
chebyshev_warped = <LagrangeVariant.chebyshev_warped: 4>
equispaced = <LagrangeVariant.equispaced: 0>
gl_centroid = <LagrangeVariant.gl_centroid: 9>
gl_isaac = <LagrangeVariant.gl_isaac: 8>
gl_warped = <LagrangeVariant.gl_warped: 7>
gll_centroid = <LagrangeVariant.gll_centroid: 3>
gll_isaac = <LagrangeVariant.gll_isaac: 2>
gll_warped = <LagrangeVariant.gll_warped: 1>
legendre = <LagrangeVariant.legendre: 10>
property name
unset = <LagrangeVariant.unset: -1>
property value
vtk = <LagrangeVariant.vtk: 20>
class basix.LatticeSimplexMethod(self: basix._basixcpp.LatticeSimplexMethod, value: int)

Bases: pybind11_object

Members:

none

warp

isaac

centroid

centroid = <LatticeSimplexMethod.centroid: 3>
isaac = <LatticeSimplexMethod.isaac: 2>
property name
none = <LatticeSimplexMethod.none: 0>
property value
warp = <LatticeSimplexMethod.warp: 1>
class basix.LatticeType(self: basix._basixcpp.LatticeType, value: int)

Bases: pybind11_object

Members:

equispaced

gll

chebyshev

gl

chebyshev = <LatticeType.chebyshev: 2>
equispaced = <LatticeType.equispaced: 0>
gl = <LatticeType.gl: 4>
gll = <LatticeType.gll: 1>
property name
property value
class basix.MapType(self: basix._basixcpp.MapType, value: int)

Bases: pybind11_object

Members:

identity

L2Piola

covariantPiola

contravariantPiola

doubleCovariantPiola

doubleContravariantPiola

L2Piola = <MapType.L2Piola: 1>
contravariantPiola = <MapType.contravariantPiola: 3>
covariantPiola = <MapType.covariantPiola: 2>
doubleContravariantPiola = <MapType.doubleContravariantPiola: 5>
doubleCovariantPiola = <MapType.doubleCovariantPiola: 4>
identity = <MapType.identity: 0>
property name
property value
class basix.PolynomialType(self: basix._basixcpp.PolynomialType, value: int)

Bases: pybind11_object

Members:

legendre

bernstein

bernstein = <PolynomialType.bernstein: 1>
legendre = <PolynomialType.legendre: 0>
property name
property value
class basix.PolysetType(self: basix._basixcpp.PolysetType, value: int)

Bases: pybind11_object

Members:

standard

macroedge

macroedge = <PolysetType.macroedge: 1>
property name
standard = <PolysetType.standard: 0>
property value
class basix.QuadratureType(self: basix._basixcpp.QuadratureType, value: int)

Bases: pybind11_object

Members:

Default

gauss_jacobi

gll

xiao_gimbutas

Default = <QuadratureType.Default: 0>
gauss_jacobi = <QuadratureType.gauss_jacobi: 1>
gll = <QuadratureType.gll: 2>
property name
property value
xiao_gimbutas = <QuadratureType.xiao_gimbutas: 3>
class basix.SobolevSpace(self: basix._basixcpp.SobolevSpace, value: int)

Bases: pybind11_object

Members:

L2

H1

H2

H3

HInf

HDiv

HCurl

HEin

HDivDiv

H1 = <SobolevSpace.H1: 1>
H2 = <SobolevSpace.H2: 2>
H3 = <SobolevSpace.H3: 3>
HCurl = <SobolevSpace.HCurl: 11>
HDiv = <SobolevSpace.HDiv: 10>
HDivDiv = <SobolevSpace.HDivDiv: 13>
HEin = <SobolevSpace.HEin: 12>
HInf = <SobolevSpace.HInf: 8>
L2 = <SobolevSpace.L2: 0>
property name
property value
basix.compute_interpolation_operator(arg0: basix._basixcpp.FiniteElement, arg1: basix._basixcpp.FiniteElement) numpy.ndarray[numpy.float64]

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:
  • element_from – The element to interpolate from

  • element_to – 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(arg0: basix._basixcpp.CellType, arg1: List[int], arg2: numpy.ndarray[numpy.float64], arg3: List[List[numpy.ndarray[numpy.float64]]], arg4: List[List[numpy.ndarray[numpy.float64]]], arg5: int, arg6: basix._basixcpp.MapType, arg7: basix._basixcpp.SobolevSpace, arg8: bool, arg9: int, arg10: int, arg11: basix::polyset::type) basix._basixcpp.FiniteElement

Create a custom finite element

Parameters:
  • cell_type – The cell type

  • value_shape – The value shape of the element

  • wcoeffs – Matrices for the kth 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 – The interpolation matrices. Indices are (tdim, entity index, dof, vs, point_index, derivative)

  • interpolation_nderivs – The number of derivatives that need to be used during interpolation

  • map_type – The type of map to be used to map values from the reference to a cell

  • sobolev_space – The underlying Sobolev space for the element

  • discontinuous – Indicates whether or not this is the discontinuous version of the element

  • highest_complete_degree – The highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element

  • highest_degree – The degree of a polynomial in this element’s polyset

  • poly_type – The type of polyset to use for this element

Returns:

A custom finite element

basix.create_element(family_name: basix._basixcpp.ElementFamily, cell_name: basix._basixcpp.CellType, degree: int, lagrange_variant: basix._basixcpp.LagrangeVariant = <LagrangeVariant.unset: -1>, dpc_variant: basix._basixcpp.DPCVariant = <DPCVariant.unset: -1>, discontinuous: bool = False, dof_ordering: List[int] = []) basix._basixcpp.FiniteElement

Create an element using a given Lagrange variant and a given DPC variant

Parameters:
  • family – The element family

  • cell – The reference cell type that the element is defined on

  • degree – The degree of the element

  • lvariant – The variant of Lagrange to use

  • dvariant – The variant of DPC to use

  • discontinuous – Indicates whether the element is discontinuous between cells points of the element. The discontinuous element will have the same DOFs, but they will all be associated with the interior of the cell.

  • dof_ordering – Ordering of dofs for ElementDofLayout

Returns:

A finite element

basix.create_lattice(*args, **kwargs)

Overloaded function.

  1. create_lattice(arg0: basix::cell::type, arg1: int, arg2: basix._basixcpp.LatticeType, arg3: bool) -> numpy.ndarray[numpy.float64]

Create a lattice of points on a reference cell optionally including the outer surface points.

For a given celltype, this creates a set of points on a regular grid which covers the cell, eg for a quadrilateral, with n=2, the points are: [0,0], [0.5,0], [1,0], [0,0.5], [0.5,0.5], [1,0.5], [0,1], [0.5,1], [1,1]. If the parameter exterior is set to false, the points lying on the external boundary are omitted, in this case for a quadrilateral with n == 2, the points are: [0.5, 0.5]. The lattice type can be chosen as type::equispaced or type::gll. The type::gll lattice has points spaced along each edge at the Gauss-Lobatto-Legendre quadrature points. These are the same as type::equispaced when n < 3.

Parameters:
  • celltype – The cell type

  • n – Size in each direction. There are n + 1 points along each edge of the cell

  • type – A lattice type

  • exterior – If set, includes outer boundaries

Returns:

Set of points. Shape is (npoints, tdim) and storage is row-major

  1. create_lattice(arg0: basix::cell::type, arg1: int, arg2: basix._basixcpp.LatticeType, arg3: bool, arg4: basix._basixcpp.LatticeSimplexMethod) -> numpy.ndarray[numpy.float64]

Create a lattice of points on a reference cell optionally including the outer surface points.

For a given celltype, this creates a set of points on a regular grid which covers the cell, eg for a quadrilateral, with n=2, the points are: [0,0], [0.5,0], [1,0], [0,0.5], [0.5,0.5], [1,0.5], [0,1], [0.5,1], [1,1]. If the parameter exterior is set to false, the points lying on the external boundary are omitted, in this case for a quadrilateral with n == 2, the points are: [0.5, 0.5]. The lattice type can be chosen as type::equispaced or type::gll. The type::gll lattice has points spaced along each edge at the Gauss-Lobatto-Legendre quadrature points. These are the same as type::equispaced when n < 3.

Parameters:
  • celltype – The cell type

  • n – Size in each direction. There are n + 1 points along each edge of the cell

  • type – A lattice type

  • exterior – If set, includes outer boundaries

  • simplex_method – The method used to generate points on simplices

Returns:

Set of points. Shape is (npoints, tdim) and storage is row-major

basix.geometry(arg0: basix::cell::type) numpy.ndarray[numpy.float64]

Cell geometry

Parameters:

celltype – Cell Type

Returns::

(0) Vertex point data of the cell and (1) the shape of the data array. The points are stored in row-major format and the shape is is (npoints, gdim)

basix.index(*args, **kwargs)

Overloaded function.

  1. index(arg0: int) -> int

Compute trivial indexing in a 1D array (for completeness)

Parameters:

p – Index in x

Returns:

1D Index

  1. index(arg0: int, arg1: int) -> int

Compute indexing in a 2D triangular array compressed into a 1D array. This can be used to find the index of a derivative returned by FiniteElement::tabulate. For instance to find d2N/dx2, use FiniteElement::tabulate(2, points)[idx(2, 0)];

Parameters:
  • p – Index in x

  • q – Index in y

Returns:

1D Index

  1. index(arg0: int, arg1: int, arg2: int) -> int

Compute indexing in a 3D tetrahedral array compressed into a 1D array

Parameters:
  • p – Index in x

  • q – Index in y

  • r – Index in z

Returns:

1D Index

basix.make_quadrature(cell: ~basix._basixcpp.CellType, degree: int, rule: ~basix._basixcpp.QuadratureType = <QuadratureType.Default: 0>, polyset_type: ~basix._basixcpp.PolysetType = <PolysetType.standard: 0>) 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:

The quadrature points and weights.

basix.polyset_restriction()

restriction(arg0: basix._basixcpp.PolysetType, arg1: basix._basixcpp.CellType, arg2: basix._basixcpp.CellType) -> basix._basixcpp.PolysetType

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

Parameters:
  • ptype – The polyset type

  • cell – The cell type

  • restriction_cell – The cell type of the subentity

Returns::

The restricted polyset type

basix.polyset_superset()

superset(arg0: basix._basixcpp.CellType, arg1: basix._basixcpp.PolysetType, arg2: basix._basixcpp.PolysetType) -> basix._basixcpp.PolysetType

Get the polyset types 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(arg0: basix._basixcpp.PolynomialType, arg1: basix::cell::type, arg2: int, arg3: numpy.ndarray[numpy.float64]) numpy.ndarray[numpy.float64]

Tabulate a set of polynomials.

Parameters:
  • polytype – Polynomial type

  • celltype – Cell type

  • d – Polynomial degree

  • x – Points at which to evaluate the basis. The shape is (number of points, geometric dimension).

Returns:

Polynomial sets, for each derivative, tabulated at points. The shape is (basis index, number of points).

basix.topology(arg0: basix::cell::type) List[List[List[int]]]

Cell topology

Parameters:

celltype – Cell Type

Returns::

List of topology (vertex indices) for each dimension (0..tdim)

Modules

basix.cell

Functions to get cell geometry information and manipulate cell types.

basix.finite_element

Functions for creating finite elements.

basix.lattice

Functions to manipulate lattice types.

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.variants

Functions to manipulate variant types.