|
Tabulation Project basix
|
basix More...
Namespaces | |
| cell | |
| cr | |
| Crouzeix-Raviart element. | |
| dofperms | |
| moments | |
| polyset | |
| quadrature | |
Classes | |
| class | FiniteElement |
Functions | |
| int | register_element (const char *family_name, const char *cell_type, int degree) |
| Create element in global registry and return handle. | |
| void | release_element (int handle) |
| Delete from global registry. | |
| std::vector< Eigen::ArrayXXd > | tabulate (int handle, int nd, const Eigen::ArrayXXd &x) |
| Tabulate. | |
| const char * | cell_type (int handle) |
| Cell type. | |
| int | degree (int handle) |
| Degree. | |
| int | value_size (int handle) |
| Value size. | |
| const std::vector< int > & | value_shape (int handle) |
| Value shape. | |
| int | dim (int handle) |
| Finite Element dimension. | |
| const char * | family_name (int handle) |
| Family name. | |
| const char * | mapping_name (int handle) |
| Mapping name (affine, piola etc.) | |
| const std::vector< std::vector< int > > & | entity_dofs (int handle) |
| Number of dofs per entity, ordered from vertex, edge, facet, cell. | |
| const std::vector< Eigen::MatrixXd > & | base_permutations (int handle) |
| Base permutations. | |
| const Eigen::ArrayXXd & | points (int handle) |
| Interpolation points. | |
| const Eigen::MatrixXd & | interpolation_matrix (int handle) |
| Interpolation matrix. | |
| Eigen::ArrayXXd | geometry (const char *cell_type) |
| Cell geometry. | |
| std::vector< std::vector< std::vector< int > > > | topology (const char *cell_type) |
| Cell topology. | |
| std::pair< Eigen::ArrayXXd, Eigen::ArrayXd > | make_quadrature (const char *rule, const char *cell_type, int order) |
| Create quadrature points and weights. | |
| FiniteElement | create_bdm (cell::type celltype, int degree, const std::string &name=std::string()) |
| Eigen::MatrixXd | compute_expansion_coefficients (const Eigen::MatrixXd &span_coeffs, const Eigen::MatrixXd &dual, bool condition_check=false) |
| FiniteElement | create_element (std::string family, std::string cell, int degree) |
| Create an element by name. | |
| const 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, const std::string &name=std::string()) |
| FiniteElement | create_dlagrange (cell::type celltype, int degree, const std::string &name=std::string()) |
| FiniteElement | create_nedelec (cell::type celltype, int degree, const std::string &name=std::string()) |
| FiniteElement | create_nedelec2 (cell::type celltype, int degree, const std::string &name=std::string()) |
| FiniteElement | create_rt (cell::type celltype, int degree, const std::string &=std::string()) |
| FiniteElement | create_regge (cell::type celltype, int degree, const std::string &name=std::string()) |
basix
| Eigen::MatrixXd basix::compute_expansion_coefficients | ( | const Eigen::MatrixXd & | span_coeffs, |
| const Eigen::MatrixXd & | dual, | ||
| bool | condition_check = false |
||
| ) |
Calculates the basis functions of the finite element, in terms of the polynomial basis.
The basis functions \((\phi_i)\) of a finite element can be 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 = \sum_j c_{ij} p_j \]
This function computed the matrix \(C = (c_{ij})\).
In some cases, the basis functions \((\phi_i)\) do not span the full space \(P_k\). In these cases, we represent the space spanned by the basis functions as the span of some polynomials \((q_k)\). These can be represented in terms of the underlying polynomial basis:
\[ q_k = \sum_j b_{kj} p_j \]
If the basis functions span the full space, then \(B = (b_{kj})\) is simply the identity.
The basis functions \(\phi_i\) are defined by a dual set of functionals \((f_l)\). The basis functions are the functions in span{ \(q_k\)} such that:
\[ f_l(\phi_i) = 1 \mbox{ if } i=l \mbox{ else } 0 \]
We can define a matrix D given by applying the functionals to each polynomial p_j:
\[ D = (d_{lj}),\mbox{ where } d_{lj} = f_l(p_j) \]
This function takes the matrices B (span_coeffs) and D (dual) as inputs and returns the matrix C. It computed C using:
\[ C = (B D^T)^{-1} B \]
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: dual =
\[ \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] | span_coeffs | The matrix B containing the expansion coefficients defining a polynomial basis spanning the polynomial space for this element. |
| [in] | dual | The matrix D of values obtained by applying each functional in the dual set to each expansion polynomial |
| [in] | condition_check | If set, checks the condition of the matrix B.D^T and throws an error if it is ill-conditioned. |
| FiniteElement basix::create_bdm | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create BDM element
| celltype | |
| degree | |
| name |
| FiniteElement basix::create_dlagrange | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create a Discontinuous Lagrange element on cell with given degree
| celltype | interval, triangle or tetrahedral celltype | |
| [in] | degree | |
| [in] | name | Identifier string (optional) |
| FiniteElement basix::create_lagrange | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create a Lagrange element on cell with given degree
| [in] | celltype | interval, triangle or tetrahedral celltype |
| [in] | degree | |
| [in] | name | Identifier string (optional) |
| FiniteElement basix::create_nedelec | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create Nedelec element (first kind)
| celltype | |
| degree | |
| name |
| FiniteElement basix::create_nedelec2 | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create Nedelec element (second kind)
| celltype | |
| degree | |
| name |
| FiniteElement basix::create_regge | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create Regge element
| celltype | |
| degree | |
| name |
| FiniteElement basix::create_rt | ( | cell::type | celltype, |
| int | degree, | ||
| const std::string & | name = std::string() |
||
| ) |
Create Raviart-Thomas element
| celltype | |
| degree | |
| name |
|
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 |
| const std::string & basix::version | ( | ) |
Return the version number of basix across projects