Placeholder. More...
Namespaces | |
cell | |
doftransforms | |
maps | |
Information about finite element maps. | |
moments | |
polyset | |
precompute | |
quadrature | |
basix | |
Classes | |
class | FiniteElement |
Functions | |
FiniteElement | create_bdm (cell::type celltype, int degree) |
FiniteElement | create_bubble (cell::type celltype, int degree) |
FiniteElement | create_cr (cell::type celltype, int degree) |
xt::xtensor< double, 3 > | compute_expansion_coefficients (cell::type cell_type, const xt::xtensor< double, 2 > &B, const std::vector< std::vector< xt::xtensor< double, 3 >>> &M, const std::vector< std::vector< xt::xtensor< double, 2 >>> &x, int degree, double kappa_tol=0.0) |
FiniteElement | create_element (std::string family, std::string cell, int degree) |
Create an element by name. | |
FiniteElement | create_element (element::family family, cell::type cell, int degree) |
Create an element by name. | |
std::string | version () |
constexpr int | idx (int p) |
constexpr int | idx (int p, int q) |
constexpr int | idx (int p, int q, int r) |
FiniteElement | create_lagrange (cell::type celltype, int degree) |
FiniteElement | create_dlagrange (cell::type celltype, int degree) |
FiniteElement | create_dpc (cell::type celltype, int degree) |
FiniteElement | create_rtc (cell::type celltype, int degree) |
FiniteElement | create_nce (cell::type celltype, int degree) |
FiniteElement | create_nedelec (cell::type celltype, int degree) |
FiniteElement | create_nedelec2 (cell::type celltype, int degree) |
FiniteElement | create_rt (cell::type celltype, int degree) |
FiniteElement | create_regge (cell::type celltype, int degree) |
FiniteElement | create_serendipity (cell::type celltype, int degree) |
FiniteElement | create_serendipity_div (cell::type celltype, int degree) |
FiniteElement | create_serendipity_curl (cell::type celltype, int degree) |
Placeholder.
xt::xtensor< double, 3 > basix::compute_expansion_coefficients | ( | cell::type | cell_type, |
const xt::xtensor< double, 2 > & | B, | ||
const std::vector< std::vector< xt::xtensor< double, 3 >>> & | M, | ||
const std::vector< std::vector< xt::xtensor< double, 2 >>> & | x, | ||
int | degree, | ||
double | kappa_tol = 0.0 |
||
) |
Calculates the basis functions of the finite element, in terms of the polynomial basis.
The below explanation uses Einstein notation.
The basis functions \({\phi_i}\) of a finite element are represented as a linear combination of polynomials \(\{p_j\}\) in an underlying polynomial basis that span the space of all d-dimensional polynomials up to order \(k \ (P_k^d)\):
\[ \phi_i = c_{ij} p_j \]
In some cases, the basis functions \(\{\phi_i\}\) do not span the full space \(P_k\), in which case we denote space spanned by the basis functions by \(\{q_k\}\), which can be represented by:
\[ q_i = b_{ij} p_j. \]
This leads to
\[ \phi_i = c^{\prime}_{ij} q_j = c^{\prime}_{ij} b_{jk} p_k, \]
and in matrix form:
\[ \phi = C^{\prime} B p \]
If the basis functions span the full space, then \( B \) is simply the identity.
The basis functions \(\phi_i\) are defined by a dual set of functionals \(\{f_i\}\). The basis functions are the functions in span{ \(q_k\)} such that
\[ f_i(\phi_j) = \delta_{ij} \]
and inserting the expression for \(\phi_{j}\):
\[ f_i(c^{\prime}_{jk}b_{kl}p_{l}) = c^{\prime}_{jk} b_{kl} f_i \left( p_{l} \right) \]
Defining a matrix D given by applying the functionals to each polynomial \(p_j\):
\[ [D] = d_{ij},\mbox{ where } d_{ij} = f_i(p_j), \]
we have:
\[ C^{\prime} B D^{T} = I \]
and
\[ C^{\prime} = (B D^{T})^{-1}. \]
Recalling that \(C = C^{\prime} B\), where \(C\) is the matrix form of \(c_{ij}\),
\[ C = (B D^{T})^{-1} B \]
This function takes the matrices B (span_coeffs) and D (dual) as inputs and returns the matrix C.
On a triangle, the scalar expansion basis is:
\[ p_0 = \sqrt{2}/2 \qquad p_1 = \sqrt{3}(2x + y - 1) \qquad p_2 = 3y - 1 \]
These span the space \(P_1\).
Lagrange order 1 elements span the space P_1, so in this example, B (span_coeffs) is the identity matrix:
\[ B = \begin{bmatrix} 1 & 0 & 0 \\ 0 & 1 & 0 \\ 0 & 0 & 1 \end{bmatrix} \]
The functionals defining the Lagrange order 1 space are point evaluations at the three vertices of the triangle. The matrix D (dual) given by applying these to p_0 to p_2 is:
\[ \mbox{dual} = \begin{bmatrix} \sqrt{2}/2 & -\sqrt{3} & -1 \\ \sqrt{2}/2 & \sqrt{3} & -1 \\ \sqrt{2}/2 & 0 & 2 \end{bmatrix} \]
For this example, this function outputs the matrix:
\[ C = \begin{bmatrix} \sqrt{2}/3 & -\sqrt{3}/6 & -1/6 \\ \sqrt{2}/3 & \sqrt{3}/6 & -1/6 \\ \sqrt{2}/3 & 0 & 1/3 \end{bmatrix} \]
The basis functions of the finite element can be obtained by applying the matrix C to the vector \([p_0, p_1, p_2]\), giving:
\[ \begin{bmatrix} 1 - x - y \\ x \\ y \end{bmatrix} \]
On a triangle, the 2D vector expansion basis is:
\[ \begin{matrix} p_0 & = & (\sqrt{2}/2, 0) \\ p_1 & = & (\sqrt{3}(2x + y - 1), 0) \\ p_2 & = & (3y - 1, 0) \\ p_3 & = & (0, \sqrt{2}/2) \\ p_4 & = & (0, \sqrt{3}(2x + y - 1)) \\ p_5 & = & (0, 3y - 1) \end{matrix} \]
These span the space \( P_1^2 \).
Raviart-Thomas order 1 elements span a space smaller than \( P_1^2 \), so B (span_coeffs) is not the identity. It is given by:
\[ B = \begin{bmatrix} 1 & 0 & 0 & 0 & 0 & 0 \\ 0 & 0 & 0 & 1 & 0 & 0 \\ 1/12 & \sqrt{6}/48 & -\sqrt{2}/48 & 1/12 & 0 & \sqrt{2}/24 \end{bmatrix} \]
Applying the matrix B to the vector \([p_0, p_1, ..., p_5]\) gives the basis of the polynomial space for Raviart-Thomas:
\[ \begin{bmatrix} \sqrt{2}/2 & 0 \\ 0 & \sqrt{2}/2 \\ \sqrt{2}x/8 & \sqrt{2}y/8 \end{bmatrix} \]
The functionals defining the Raviart-Thomas order 1 space are integral of the normal components along each edge. The matrix D (dual) given by applying these to \(p_0\) to \(p_5\) is:
\[ D = \begin{bmatrix} -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 & -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 \\ -\sqrt{2}/2 & \sqrt{3}/2 & -1/2 & 0 & 0 & 0 \\ 0 & 0 & 0 & \sqrt{2}/2 & 0 & -1 \end{bmatrix} \]
In this example, this function outputs the matrix:
\[ C = \begin{bmatrix} -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 & -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 \\ -\sqrt{2}/2 & \sqrt{3}/2 & -1/2 & 0 & 0 & 0 \\ 0 & 0 & 0 & \sqrt{2}/2 & 0 & -1 \end{bmatrix} \]
The basis functions of the finite element can be obtained by applying the matrix C to the vector \([p_0, p_1, ..., p_5]\), giving:
\[ \begin{bmatrix} -x & -y \\ x - 1 & y \\ -x & 1 - y \end{bmatrix} \]
[in] | cell_type | The cells shape |
[in] | B | Matrices for the kth value index containing the expansion coefficients defining a polynomial basis spanning the polynomial space for this element |
[in] | M | The interpolation tensor, such that the dual matrix \(D\) is computed by \(D = MP\) |
[in] | x | The interpolation points. The vector index is for points on entities of the same dimension, ordered with the lowest topological dimension being first. Each 3D tensor hold the points on cell entities of a common dimension. The shape of the 3d tensors is (num_entities, num_points_per_entity, tdim). |
[in] | degree | The degree of the polynomial basis P used to create the element (before applying B) |
[in] | kappa_tol | If positive, the condition number is computed and an error thrown if the condition number of \(B D^{T}\) is greater than kappa_tol . If kappa_tol is less than 1 the condition number is not checked. |
Flatten D and take transpose
FiniteElement basix::create_bdm | ( | cell::type | celltype, |
int | degree | ||
) |
Create BDM element
celltype | |
degree |
FiniteElement basix::create_bubble | ( | cell::type | celltype, |
int | degree | ||
) |
Create a bubble element on cell with given degree
[in] | celltype | interval, triangle, tetrahedral, quadrilateral or hexahedral celltype |
[in] | degree |
FiniteElement basix::create_cr | ( | cell::type | celltype, |
int | degree | ||
) |
Crouzeix-Raviart element
celltype | |
degree |
FiniteElement basix::create_dlagrange | ( | cell::type | celltype, |
int | degree | ||
) |
Create a Discontinuous Lagrange element on cell with given degree
celltype | interval, triangle, quadrilateral, tetrahedral, or hexahedral celltype | |
[in] | degree |
FiniteElement basix::create_dpc | ( | cell::type | celltype, |
int | degree | ||
) |
Create a DPC element on cell with given degree
celltype | interval, quadrilateral or hexahedral celltype | |
[in] | degree |
FiniteElement basix::create_lagrange | ( | cell::type | celltype, |
int | degree | ||
) |
Create a Lagrange element on cell with given degree
celltype | interval, triangle, quadrilateral, tetrahedral, or hexahedral celltype | |
[in] | degree |
FiniteElement basix::create_nce | ( | cell::type | celltype, |
int | degree | ||
) |
Create NC H(curl) element
celltype | |
degree |
FiniteElement basix::create_nedelec | ( | cell::type | celltype, |
int | degree | ||
) |
Create Nedelec element (first kind)
celltype | |
degree |
FiniteElement basix::create_nedelec2 | ( | cell::type | celltype, |
int | degree | ||
) |
Create Nedelec element (second kind)
celltype | |
degree |
FiniteElement basix::create_regge | ( | cell::type | celltype, |
int | degree | ||
) |
Create Regge element
celltype | |
degree |
FiniteElement basix::create_rt | ( | cell::type | celltype, |
int | degree | ||
) |
Create Raviart-Thomas element
celltype | |
degree |
FiniteElement basix::create_rtc | ( | cell::type | celltype, |
int | degree | ||
) |
Create RTC H(div) element
celltype | |
degree |
FiniteElement basix::create_serendipity | ( | cell::type | celltype, |
int | degree | ||
) |
Create a serendipity element on cell with given degree
[in] | celltype | quadrilateral or hexahedral celltype |
[in] | degree |
FiniteElement basix::create_serendipity_curl | ( | cell::type | celltype, |
int | degree | ||
) |
Create a serendipity H(curl) element on cell with given degree
[in] | celltype | quadrilateral or hexahedral celltype |
[in] | degree |
FiniteElement basix::create_serendipity_div | ( | cell::type | celltype, |
int | degree | ||
) |
Create a serendipity H(div) element on cell with given degree
[in] | celltype | quadrilateral or hexahedral celltype |
[in] | degree |
|
constexpr |
Compute trivial indexing in a 1D array (for completeness)
p | Index in x |
|
constexpr |
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 |
|
constexpr |
Compute indexing in a 3D tetrahedral array compressed into a 1D array
p | Index in x |
q | Index in y |
r | Index in z |
std::string basix::version | ( | ) |
Return the version number of basix across projects