basix.interpolation

Interpolation.

Functions

compute_interpolation_operator(e0, e1)

Compute a matrix that represents the interpolation between two elements.

class basix.interpolation.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[Any, dtype[floating]]]]

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() ndarray[Any, dtype[floating]]

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: ndarray[Any, dtype[floating]]

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: dtype[Any] | None | type[Any] | _SupportsDType[dtype[Any]] | str | tuple[Any, int] | tuple[Any, SupportsIndex | Sequence[SupportsIndex]] | list[Any] | _DTypeDict | tuple[Any, Any]

Element float type.

property dual_matrix: ndarray[Any, dtype[floating]]

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.

property interpolation_is_identity: bool

True if interpolation matrix for this element is the identity.

property interpolation_matrix: ndarray[Any, dtype[floating]]

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.

property points: ndarray[Any, dtype[floating]]

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[Any, dtype[_ScalarType_co]], J: ndarray[Any, dtype[_ScalarType_co]], detJ: ndarray[Any, dtype[_ScalarType_co]], K: ndarray[Any, dtype[_ScalarType_co]]) ndarray[Any, dtype[floating]]

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) ndarray[Any, dtype[floating]]

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[Any, dtype[_ScalarType_co]]) ndarray[Any, dtype[floating]]

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: ndarray[Any, dtype[floating]]

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

See C++ documentation for details.

property x: list[list[ndarray[Any, dtype[floating]]]]

Interpolation points for each sub-entity.

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

basix.interpolation.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))