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.
Bases: pybind11_object
Members:
point
interval
triangle
tetrahedron
quadrilateral
hexahedron
prism
pyramid
Bases: pybind11_object
Members:
unset
simplex_equispaced
simplex_gll
horizontal_equispaced
horizontal_gll
diagonal_equispaced
diagonal_gll
legendre
Bases: pybind11_object
Members:
custom
P
BDM
RT
N1E
N2E
Regge
HHJ
bubble
serendipity
DPC
CR
Hermite
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
Bases: pybind11_object
Members:
none
warp
isaac
centroid
Bases: pybind11_object
Members:
equispaced
gll
chebyshev
gl
Bases: pybind11_object
Members:
identity
L2Piola
covariantPiola
contravariantPiola
doubleCovariantPiola
doubleContravariantPiola
Bases: pybind11_object
Members:
legendre
bernstein
Bases: pybind11_object
Members:
Default
gauss_jacobi
gll
xiao_gimbutas
Bases: pybind11_object
Members:
L2
H1
H2
H3
HInf
HDiv
HCurl
HEin
HDivDiv
Computes 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.
element_from – The element to interpolate from
element_to – The element to interpolate to
Matrix operator that maps the ‘from’ degrees-of-freedom to the ‘to’ degrees-of-freedom. Shape is (ndofs(element_to), ndofs(element_from))
Create a custom finite element
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
A custom finite element
Overloaded function.
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: bool) -> basix._basixcpp.FiniteElement
Create an element
family – The element family
cell – The reference cell type that the element is defined on
degree – The degree of the element
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.
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: basix._basixcpp.LagrangeVariant, arg4: bool) -> basix._basixcpp.FiniteElement
Create an element using a given Lagrange variant
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
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.
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: basix._basixcpp.DPCVariant, arg4: bool) -> basix._basixcpp.FiniteElement
Create an element using a given DPC variant
family – The element family
cell – The reference cell type that the element is defined on
degree – The degree of the element
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.
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: basix._basixcpp.LagrangeVariant, arg4: basix._basixcpp.DPCVariant, arg5: bool) -> basix._basixcpp.FiniteElement
Create an element using a given Lagrange variant and a given DPC variant
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.
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int) -> basix._basixcpp.FiniteElement
Create a continuous element
family – The element family
cell – The reference cell type that the element is defined on
degree – The degree of the element
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: basix._basixcpp.LagrangeVariant) -> basix._basixcpp.FiniteElement
Create a continuous element using a given Lagrange variant
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
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: basix._basixcpp.DPCVariant) -> basix._basixcpp.FiniteElement
Create a continuous element using a given DPC variant
family – The element family
cell – The reference cell type that the element is defined on
degree – The degree of the element
dvariant – The variant of DPC to use
A finite element
create_element(arg0: basix._basixcpp.ElementFamily, arg1: basix._basixcpp.CellType, arg2: int, arg3: basix._basixcpp.LagrangeVariant, arg4: basix._basixcpp.DPCVariant) -> basix._basixcpp.FiniteElement
Create a continuous element using a given Lagrange variant and a given DPC variant
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
A finite element
Overloaded function.
create_lattice(arg0: basix::cell::type, arg1: int, arg2: basix._basixcpp.LatticeType, arg3: bool) -> numpy.ndarray[numpy.float64]
@brief 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.
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
Set of points. Shape is (npoints, tdim) and storage is row-major
create_lattice(arg0: basix::cell::type, arg1: int, arg2: basix._basixcpp.LatticeType, arg3: bool, arg4: basix._basixcpp.LatticeSimplexMethod) -> numpy.ndarray[numpy.float64]
@brief 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.
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
Set of points. Shape is (npoints, tdim) and storage is row-major
Cell geometry
celltype – Cell Type
(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)
Overloaded function.
index(arg0: int) -> int
Compute trivial indexing in a 1D array (for completeness)
p – Index in x
1D Index
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)];
p – Index in x
q – Index in y
1D Index
index(arg0: int, arg1: int, arg2: int) -> int
Compute indexing in a 3D tetrahedral array compressed into a 1D array
p – Index in x
q – Index in y
r – Index in z
1D Index
Overloaded function.
make_quadrature(arg0: basix._basixcpp.QuadratureType, arg1: basix._basixcpp.CellType, arg2: int) -> Tuple[numpy.ndarray[numpy.float64], numpy.ndarray[numpy.float64]]
Make a quadrature rule on a reference cell
rule – Type of quadrature rule (or use quadrature::Default)
celltype – The cell type
m – Maximum degree of polynomial that this quadrature rule will integrate exactly
List of points and list of weights. The number of points arrays has shape (num points, gdim)
make_quadrature(arg0: basix._basixcpp.CellType, arg1: int) -> Tuple[numpy.ndarray[numpy.float64], numpy.ndarray[numpy.float64]]
Make a default quadrature rule on reference cell
celltype – The cell type
m – Maximum degree of polynomial that this quadrature rule will integrate exactly
List of points and list of weights. The number of points arrays has shape (num points, gdim)
@brief Tabulate a set of polynomials.
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).
Polynomial sets, for each derivative, tabulated at points. The shape is (basis index, number of points).
Cell topology
celltype – Cell Type
List of topology (vertex indices) for each dimension (0..tdim)
Modules
Functions to get cell geometry information and manipulate cell types. |
|
Functions for creating finite elements. |
|
Functions to manipulate lattice types. |
|
Helper functions for writing DOLFINx custom kernels using Numba. |
|
Functions for working with polynomials. |
|
Functions to manipulate quadrature types. |
|
Functions for handling Sobolev spaces. |
|
Functions to directly wrap Basix elements in UFL. |
|
Functions to manipulate variant types. |