basix.finite_element

Functions for creating finite elements.

Functions

create_custom_element(cell_type, ...[, dtype])

Create a custom finite element.

create_element(family, celltype, degree[, ...])

Create a finite element.

create_tp_element(family, celltype, degree)

Create a finite element with tensor product ordering.

lex_dof_ordering(family, celltype, degree[, ...])

Lexicographic DOF ordering for an element.

string_to_family(family, cell)

Basix ElementFamily enum representing the family type on the given cell.

tp_dof_ordering(family, celltype, degree[, ...])

Tensor product DOF ordering for an element.

tp_factors(family, celltype, degree[, ...])

Elements in the tensor product factorisation of an element.

Classes

FiniteElement(e)

Finite element class.

class basix.finite_element.FiniteElement(e: FiniteElement_float32 | FiniteElement_float64)

Bases: object

Finite element class.

Initialise a finite element wrapper.

Note

This initialiser is intended for internal library use.

property M: list[list[ndarray[tuple[Any, ...], dtype[_ScalarT]]]]

Interpolation matrices for each sub-entity.

See C++ documentation for details.

T_apply(data, block_size, cell_info) None

Apply DOF transformations to some data in-place.

Note

This function is designed to be called at runtime, so its performance is critical.

Parameters:
  • data – The data

  • block_size – The number of data points per DOF

  • cell_info – The permutation info for the cell

Tt_apply_right(data, block_size, cell_info) None

Post-apply DOF transformations to some transposed data in-place.

Note

This function is designed to be called at runtime, so its performance is critical.

Parameters:
  • data – The data.

  • block_size – The number of data points per DOF.

  • cell_info – The permutation info for the cell.

Tt_inv_apply(data, block_size, cell_info) None

Pre-apply inverse transpose DOF transformations to some data.

Note

This function is designed to be called at runtime, so its performance is critical.

Parameters:
  • data – The data.

  • block_size – The number of data points per DOF.

  • cell_info – The permutation info for the cell.

base_transformations() _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Get the base transformations.

The base transformations represent the effect of rotating or reflecting a subentity of the cell on the numbering and orientation of the DOFs. This returns a list of matrices with one matrix for each subentity permutation in the following order: Reversing edge 0, reversing edge 1, … Rotate face 0, reflect face 0, rotate face 1, reflect face 1, …

Example: Order 3 Lagrange on a triangle

This space has 10 dofs arranged like:

2
|\
6 4
|  \
5 9 3
|    \
0-7-8-1

For this element, the base transformations are: [Matrix swapping 3 and 4, Matrix swapping 5 and 6, Matrix swapping 7 and 8] The first row shows the effect of reversing the diagonal edge. The second row shows the effect of reversing the vertical edge. The third row shows the effect of reversing the horizontal edge.

Example: Order 1 Raviart-Thomas on a triangle

This space has 3 dofs arranged like:

  |\
  | \
  |  \
<-1   0
  |  / \
  | L ^ \
  |   |  \
   ---2---

These DOFs are integrals of normal components over the edges: DOFs 0 and 2 are oriented inward, DOF 1 is oriented outwards. For this element, the base transformation matrices are:

0: [[-1, 0, 0],
    [ 0, 1, 0],
    [ 0, 0, 1]]
1: [[1,  0, 0],
    [0, -1, 0],
    [0,  0, 1]]
2: [[1, 0,  0],
    [0, 1,  0],
    [0, 0, -1]]

The first matrix reverses DOF 0 (as this is on the first edge). The second matrix reverses DOF 1 (as this is on the second edge). The third matrix reverses DOF 2 (as this is on the third edge).

Example: DOFs on the face of Order 2 Nedelec first kind on a tetrahedron

On a face of this tetrahedron, this space has two face tangent DOFs:

|\        |\
| \       | \
|  \      | ^\
|   \     | | \
| 0->\    | 1  \
|     \   |     \
 ------    ------

For these DOFs, the subblocks of the base transformation matrices are:

rotation: [[-1, 1],
           [ 1, 0]]
reflection: [[0, 1],
             [1, 0]]
Returns:

The base transformations for this element. The shape is (ntranformations, ndofs, ndofs).

property cell_type: CellType

Element cell type.

property coefficient_matrix: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Matrix of coefficients.

property degree: int

Element polynomial degree.

property dim: int

Dimension of the finite element space.

This is the number of degrees-of-freedom for the element.

property discontinuous: bool

True is element is the discontinuous variant.

property dof_ordering: list[int]

DOF layout.

property dof_transformations_are_identity: bool

True if DOF transformations are all the identity.

property dof_transformations_are_permutations: bool

True if the dof transformations are all permutations.

property dpc_variant: DPCVariant

DPC variant of the element.

property dtype: type[Any] | dtype[Any] | _SupportsDType[dtype[Any]] | tuple[Any, Any] | list[Any] | _DTypeDict | str | None

Element float type.

property dual_matrix: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Matrix $BD^{T}$.

See C++ documentation.

property embedded_subdegree: int

Embedded polynomial sub-degree.

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

property embedded_superdegree: int

Embedded polynomial degree.

Lowest degree n such that the highest degree polynomial in this element is contained in a Lagrange (or vector Lagrange) element of degree n.

property entity_closure_dofs: list[list[list[int]]]

Get the dofs on the closure of each topological entity.

Data is in the order (vertices, edges, faces, cell). For example, Lagrange degree 2 on a triangle has vertices: [[0], [1], [2]], edges: [[1, 2, 3], [0, 2, 4], [0, 1, 5]], cell: [[0, 1, 2, 3, 4, 5]].

property entity_dofs: list[list[list[int]]]

Dofs on each topological entity.

Data is order (vertices, edges, faces, cell). For example, Lagrange degree 2 on a triangle has vertices: [[0], [1], [2]], edges: [[3], [4], [5]], cell: [[]].

entity_transformations() dict

Entity dof transformation matrices.

Returns:

The base transformations for this element. The shape is (ntranformations, ndofs, ndofs).

property family: ElementFamily

Finite element family.

get_tensor_product_representation() list[list[FiniteElement]]

Get the tensor product representation of this element.

Raises an exception if no such factorisation exists.

The tensor product representation will be a vector of tuples. Each tuple contains a vector of finite elements, and a vector of integers. The vector of finite elements gives the elements on an interval that appear in the tensor product representation. The vector of integers gives the permutation between the numbering of the tensor product DOFs and the number of the DOFs of this Basix element.

Returns:

The tensor product representation

property has_tensor_product_factorisation: bool

True if element has tensor-product structure.

hash() int

Hash.

property interpolation_is_identity: bool

True if interpolation matrix for this element is the identity.

property interpolation_matrix: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Interpolation points.

Coordinates on the reference element where a function need to be evaluated in order to interpolate it in the finite element space.

property interpolation_nderivs: int

Number of derivatives needed when interpolating.

property lagrange_variant: LagrangeVariant

Lagrange variant of the element.

property map_type: MapType

Map type for this element.

property num_entity_closure_dofs: list[list[int]]

Number of entity closure dofs.

Warning

This property may be removed.

property num_entity_dofs: list[list[int]]

Number of entity dofs.

Warning

This property may be removed.

permute_subentity_closure(indices: ndarray[tuple[Any, ...], dtype[_ScalarT]], cell_or_entity_info: int, entity_type: CellType, entity_index: int | None = None) ndarray[tuple[Any, ...], dtype[_ScalarT]]

Permute DOF indices on the closure of a sub-entity.

Parameters:
  • indices – The indices to permute

  • cell_or_entity_info – Bit packed entity info (if entity_index is None) or cell info (if entity_index is not None)

  • entity_type – The cell type of the entity

  • entity_index – The index of the entity

permute_subentity_closure_inv(indices: ndarray[tuple[Any, ...], dtype[_ScalarT]], cell_or_entity_info: int, entity_type: CellType, entity_index: int | None = None) ndarray[tuple[Any, ...], dtype[_ScalarT]]

Apply inverse permutation to DOF indices on the closure of a sub-entity.

Parameters:
  • indices – The indices to permute

  • cell_or_entity_info – Bit packed entity info (if entity_index is None) or cell info (if entity_index is not None)

  • entity_type – The cell type of the entity

  • entity_index – The index of the entity

property points: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Interpolation points.

Coordinates on the reference element where a function need to be evaluated in order to interpolate it in the finite element space. Shape is (num_points, tdim).

property polyset_type: PolysetType

Element polyset type.

pull_back(u: ndarray[tuple[Any, ...], dtype[_ScalarT]], J: ndarray[tuple[Any, ...], dtype[_ScalarT]], detJ: ndarray[tuple[Any, ...], dtype[_ScalarT]], K: ndarray[tuple[Any, ...], dtype[_ScalarT]]) _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Map function values from a physical cell to the reference.

Parameters:
  • u – The function values on the cell.

  • J – The Jacobian of the mapping.

  • detJ – The determinant of the Jacobian of the mapping.

  • K – The inverse of the Jacobian of the mapping.

Returns:

The function values on the reference. The indices are (Jacobian index, point index, components).

push_forward(U, J, detJ, K) _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Map function values from the reference to a physical cell.

This function can perform the mapping for multiple points, grouped by points that share a common Jacobian.

Parameters:
  • U – The function values on the reference cell. The indices are (Jacobian index, point index, components).

  • J – The Jacobian of the mapping. The indices are (Jacobian index, J_i, J_j).

  • detJ – The determinant of the Jacobian of the mapping. It has length J.shape(0).

  • K – The inverse of the Jacobian of the mapping. The indices are (Jacobian index, K_i, K_j).

Returns:

The function values on the cell. The indices are (Jacobian index, point index, components).

property sobolev_space: SobolevSpace

Underlying Sobolev space for this element.

tabulate(n: int, x: ndarray[tuple[Any, ...], dtype[_ScalarT]]) _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Compute basis values and derivatives at set of points.

Note

The version of FiniteElement::tabulate with the basis data as an out argument should be preferred for repeated call where performance is critical.

Parameters:
  • n – The order of derivatives, up to and including, to compute. Use 0 for the basis functions only.

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

Returns:

The basis functions (and derivatives). The shape is (derivative, point, basis fn index, value index).

  • The first index is the derivative, with higher derivatives

    are stored in triangular (2D) or tetrahedral (3D) ordering, i.e. for the (x,y) derivatives in 2D: (0,0), (1,0), (0,1), (2,0), (1,1), (0,2), (3,0).... The function basix::indexing::idx can be used to find the appropriate derivative.

  • The second index is the point index

  • The fourth index is the basis function component. Its has size

  • The third index is the basis function index one for scalar

    basis functions.

property value_shape: list[int]

Element value tensor shape.

E.g., returning (,) for scalars, (3,) for vectors in 3D, (2, 2) for a rank-2 tensor in 2D.

property value_size: int

Value size.

property wcoeffs: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]

Coefficients that define the polynomial set in terms of the orthonormal polynomials.

See C++ documentation for details.

property x: list[list[ndarray[tuple[Any, ...], dtype[_ScalarT]]]]

Interpolation points for each sub-entity.

The indices of this data are (tdim, entity index, point index, dim).

basix.finite_element.create_custom_element(cell_type: ~basix._basixcpp.CellType, value_shape: tuple[int, ...], wcoeffs: ~numpy.ndarray[tuple[~typing.Any, ...], ~numpy.dtype[~numpy.floating]], x: list[list[~numpy.ndarray[tuple[~typing.Any, ...], ~numpy.dtype[~numpy.floating]]]], M: list[list[~numpy.ndarray[tuple[~typing.Any, ...], ~numpy.dtype[~numpy.floating]]]], interpolation_nderivs: int, map_type: ~basix._basixcpp.MapType, sobolev_space: ~basix._basixcpp.SobolevSpace, discontinuous: bool, embedded_subdegree: int, embedded_superdegree: int, poly_type: ~basix._basixcpp.PolysetType, dtype: type[~typing.Any] | ~numpy.dtype[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | tuple[~typing.Any, ~typing.Any] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | str | None = <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.finite_element.create_element(family: ~basix._basixcpp.ElementFamily, celltype: ~basix._basixcpp.CellType, degree: int, lagrange_variant: ~basix._basixcpp.LagrangeVariant = LagrangeVariant.unset, dpc_variant: ~basix._basixcpp.DPCVariant = DPCVariant.unset, discontinuous: bool = False, dof_ordering: list[int] | None = None, dtype: type[~typing.Any] | ~numpy.dtype[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | tuple[~typing.Any, ~typing.Any] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | str | None = <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.finite_element.create_tp_element(family: ~basix._basixcpp.ElementFamily, celltype: ~basix._basixcpp.CellType, degree: int, lagrange_variant: ~basix._basixcpp.LagrangeVariant = LagrangeVariant.unset, dpc_variant: ~basix._basixcpp.DPCVariant = DPCVariant.unset, discontinuous: bool = False, dtype: type[~typing.Any] | ~numpy.dtype[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | tuple[~typing.Any, ~typing.Any] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | str | None = <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.finite_element.lex_dof_ordering(family: ElementFamily, celltype: CellType, degree: int, lagrange_variant: LagrangeVariant = LagrangeVariant.unset, dpc_variant: DPCVariant = DPCVariant.unset, discontinuous: bool = False) list[int]

Lexicographic DOF ordering for an element.

This DOF ordering can be passed into create_element to create the element with DOFs ordered in a lexicographic order.

If the element contains DOFs that are not point evaluations, raises a RuntimeError.

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.

Returns:

The DOF ordering.

basix.finite_element.string_to_family(family: str, cell: str) ElementFamily

Basix ElementFamily enum representing the family type on the given cell.

Parameters:
  • family – Element family as a string.

  • cell – Cell type as a string.

Returns:

Element family.

basix.finite_element.tp_dof_ordering(family: ElementFamily, celltype: CellType, degree: int, lagrange_variant: LagrangeVariant = LagrangeVariant.unset, dpc_variant: DPCVariant = DPCVariant.unset, discontinuous: bool = False) list[int] | None

Tensor product DOF ordering for an element.

This DOF ordering can be passed into create_element to create the element with DOFs ordered in a tensor product order.

If the element has no factorisation, raises a RuntimeError.

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.

Returns:

The DOF ordering.

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

Elements in the tensor product factorisation of an element.

If the element has no factorisation, raises a RuntimeError.

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 list of finite elements.