CeedBasis¶
A CeedBasis defines the discrete finite element basis and associated quadrature rule.
Discrete element bases and quadrature¶
-
typedef struct CeedBasis_private *CeedBasis¶
Handle for object describing discrete finite element evaluations.
-
const CeedBasis CEED_BASIS_NONE = &ceed_basis_none¶
Argument for CeedOperatorSetField indicating that the field does not require a CeedBasis.
-
const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_none¶
This feature will be removed. Use CEED_BASIS_NONE.
-
int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis)¶
Create a tensor-product basis for H^1 discretizations.
User Functions
- Parameters
ceed – [in] Ceed object where the CeedBasis will be created
dim – [in] Topological dimension
num_comp – [in] Number of field components (1 for scalar fields)
P_1d – [in] Number of nodes in one dimension
Q_1d – [in] Number of quadrature points in one dimension
interp_1d – [in] Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
grad_1d – [in] Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
q_ref_1d – [in] Array of length Q_1d holding the locations of quadrature points on the 1D reference element [-1, 1]
q_weight_1d – [in] Array of length Q_1d holding the quadrature weights on the reference element
basis – [out] Address of the variable where the newly created CeedBasis will be stored.
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis)¶
Create a tensor-product Lagrange basis.
User Functions
- Parameters
ceed – [in] Ceed object where the CeedBasis will be created
dim – [in] Topological dimension of element
num_comp – [in] Number of field components (1 for scalar fields)
P – [in] Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the resulting Q_k element is k=P-1.
Q – [in] Number of quadrature points in one dimension.
quad_mode – [in] Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
basis – [out] Address of the variable where the newly created CeedBasis will be stored.
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis)¶
Create a non tensor-product basis for H^1 discretizations.
User Functions
- Parameters
ceed – [in] Ceed object where the CeedBasis will be created
topo – [in] Topology of element, e.g. hypercube, simplex, ect
num_comp – [in] Number of field components (1 for scalar fields)
num_nodes – [in] Total number of nodes
num_qpts – [in] Total number of quadrature points
interp – [in] Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
grad – [in] Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
q_ref – [in] Array of length num_qpts * dim holding the locations of quadrature points on the reference element
q_weight – [in] Array of length num_qpts holding the quadrature weights on the reference element
basis – [out] Address of the variable where the newly created CeedBasis will be stored.
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis)¶
Create a non tensor-product basis for \(H(\mathrm{div})\) discretizations.
User Functions
- Parameters
ceed – [in] Ceed object where the CeedBasis will be created
topo – [in] Topology of element (
CEED_TOPOLOGY_QUAD
,CEED_TOPOLOGY_PRISM
, etc.), dimension of which is used in some array sizes belownum_comp – [in] Number of components (usually 1 for vectors in H(div) bases)
num_nodes – [in] Total number of nodes (dofs per element)
num_qpts – [in] Total number of quadrature points
interp – [in] Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
div – [in] Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points
q_ref – [in] Array of length num_qpts * dim holding the locations of quadrature points on the reference element
q_weight – [in] Array of length num_qpts holding the quadrature weights on the reference element
basis – [out] Address of the variable where the newly created CeedBasis will be stored.
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis)¶
Create a non tensor-product basis for \(H(\mathrm{curl})\) discretizations.
User Functions
- Parameters
ceed – [in] Ceed object where the CeedBasis will be created
topo – [in] Topology of element (
CEED_TOPOLOGY_QUAD
,CEED_TOPOLOGY_PRISM
, etc.), dimension of which is used in some array sizes belownum_comp – [in] Number of components (usually 1 for vectors in H(curl) bases)
num_nodes – [in] Total number of nodes (dofs per element)
num_qpts – [in] Total number of quadrature points
interp – [in] Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
curl – [in] Row-major (curl_comp * num_qpts * num_nodes, curl_comp = 1 if dim < 3 else dim) matrix expressing curl of basis functions at quadrature points
q_ref – [in] Array of length num_qpts * dim holding the locations of quadrature points on the reference element
q_weight – [in] Array of length num_qpts holding the quadrature weights on the reference element
basis – [out] Address of the variable where the newly created CeedBasis will be stored.
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project)¶
Create a CeedBasis for projection from the nodes of
basis_from
to the nodes ofbasis_to
.Only
CEED_EVAL_INTERP
will be valid for the new basis,basis_project
. For H^1 spaces,CEED_EVAL_GRAD
will also be valid. The interpolation is given byinterp_project = interp_to^+ * interp_from
, where the pseudoinverseinterp_to^+
is given by QR factorization. The gradient (for the H^1 case) is given bygrad_project = interp_to^+ * grad_from
.Note:
basis_from
andbasis_to
must have compatible quadrature spaces.Note:
basis_project
will have the same number of components asbasis_from
, regardless of the number of components thatbasis_to
has. Ifbasis_from
has 3 components andbasis_to
has 5 components, thenbasis_project
will have 3 components.User Functions
- Parameters
basis_from – [in] CeedBasis to prolong from
basis_to – [in] CeedBasis to prolong to
basis_project – [out] Address of the variable where the newly created CeedBasis will be stored.
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy)¶
Copy the pointer to a CeedBasis.
Note: If the value of
basis_copy
passed into this function is non-NULL, then it is assumed thatbasis_copy
is a pointer to a CeedBasis. This CeedBasis will be destroyed ifbasis_copy
is the only reference to this CeedBasis.User Functions
- Parameters
basis – [in] CeedBasis to copy reference to
basis_copy – [inout] Variable to store copied reference
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisView(CeedBasis basis, FILE *stream)¶
View a CeedBasis.
User Functions
- Parameters
basis – [in] CeedBasis to view
stream – [in] Stream to view to, e.g., stdout
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v)¶
Apply basis evaluation from nodes to quadrature points or vice versa.
User Functions
- Parameters
basis – [in] CeedBasis to evaluate
num_elem – [in] The number of elements to apply the basis evaluation to; the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
t_mode – [in] CEED_NOTRANSPOSE to evaluate from nodes to quadrature points; CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
eval_mode – [in] CEED_EVAL_NONE to use values directly, CEED_EVAL_INTERP to use interpolated values, CEED_EVAL_GRAD to use gradients, CEED_EVAL_DIV to use divergence, CEED_EVAL_CURL to use curl, CEED_EVAL_WEIGHT to use quadrature weights.
u – [in] Input CeedVector
v – [out] Output CeedVector
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v)¶
Apply basis evaluation from nodes to arbitrary points.
User Functions
- Parameters
basis – [in] CeedBasis to evaluate
num_points – [in] The number of points to apply the basis evaluation to
t_mode – [in] CEED_NOTRANSPOSE to evaluate from nodes to points; CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
eval_mode – [in] CEED_EVAL_INTERP to use interpolated values, CEED_EVAL_GRAD to use gradients
x_ref – [in] CeedVector holding reference coordinates of each point
u – [in] Input CeedVector, of length
num_nodes * num_comp
forCEED_NOTRANSPOSE
v – [out] Output CeedVector, of length
num_points * num_q_comp
forCEED_NOTRANSPOSE
withCEED_EVAL_INTERP
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed)¶
Get Ceed associated with a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
ceed – [out] Variable to store Ceed
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim)¶
Get dimension for given CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
dim – [out] Variable to store dimension of basis
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo)¶
Get topology for given CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
topo – [out] Variable to store topology of basis
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp)¶
Get number of components for given CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
num_comp – [out] Variable to store number of components of basis
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P)¶
Get total number of nodes (in dim dimensions) of a CeedBasis.
Utility Functions
- Parameters
basis – [in] CeedBasis
P – [out] Variable to store number of nodes
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d)¶
Get total number of nodes (in 1 dimension) of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
P_1d – [out] Variable to store number of nodes
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q)¶
Get total number of quadrature points (in dim dimensions) of a CeedBasis.
Utility Functions
- Parameters
basis – [in] CeedBasis
Q – [out] Variable to store number of quadrature points
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d)¶
Get total number of quadrature points (in 1 dimension) of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
Q_1d – [out] Variable to store number of quadrature points
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref)¶
Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
q_ref – [out] Variable to store reference coordinates of quadrature points
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight)¶
Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
q_weight – [out] Variable to store quadrature weights
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp)¶
Get interpolation matrix of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
interp – [out] Variable to store interpolation matrix
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d)¶
Get 1D interpolation matrix of a tensor product CeedBasis.
Backend Developer Functions
- Parameters
basis – [in] CeedBasis
interp_1d – [out] Variable to store interpolation matrix
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad)¶
Get gradient matrix of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
grad – [out] Variable to store gradient matrix
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d)¶
Get 1D gradient matrix of a tensor product CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
grad_1d – [out] Variable to store gradient matrix
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div)¶
Get divergence matrix of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
div – [out] Variable to store divergence matrix
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl)¶
Get curl matrix of a CeedBasis.
Advanced Functions
- Parameters
basis – [in] CeedBasis
curl – [out] Variable to store curl matrix
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedBasisDestroy(CeedBasis *basis)¶
Destroy a CeedBasis.
User Functions
- Parameters
basis – [inout] CeedBasis to destroy
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d)¶
Construct a Gauss-Legendre quadrature.
Utility Functions
- Parameters
Q – [in] Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
q_ref_1d – [out] Array of length Q to hold the abscissa on [-1, 1]
q_weight_1d – [out] Array of length Q to hold the weights
- Returns
An error code: 0 - success, otherwise - failure
-
int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d)¶
Construct a Gauss-Legendre-Lobatto quadrature.
Utility Functions
- Parameters
Q – [in] Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
q_ref_1d – [out] Array of length Q to hold the abscissa on [-1, 1]
q_weight_1d – [out] Array of length Q to hold the weights
- Returns
An error code: 0 - success, otherwise - failure
Typedefs and Enumerations¶
-
enum CeedTransposeMode¶
Denotes whether a linear transformation or its transpose should be applied.
Values:
-
enumerator CEED_NOTRANSPOSE¶
Apply the linear transformation.
-
enumerator CEED_TRANSPOSE¶
Apply the transpose.
-
enumerator CEED_NOTRANSPOSE¶
-
enum CeedEvalMode¶
Basis evaluation mode.
Values:
-
enumerator CEED_EVAL_NONE¶
Perform no evaluation (either because there is no data or it is already at quadrature points)
-
enumerator CEED_EVAL_INTERP¶
Interpolate from nodes to quadrature points.
-
enumerator CEED_EVAL_GRAD¶
Evaluate gradients at quadrature points from input in the basis.
-
enumerator CEED_EVAL_DIV¶
Evaluate divergence at quadrature points from input in the basis.
-
enumerator CEED_EVAL_CURL¶
Evaluate curl at quadrature points from input in the basis.
-
enumerator CEED_EVAL_WEIGHT¶
Using no input, evaluate quadrature weights on the reference element.
-
enumerator CEED_EVAL_NONE¶
-
enum CeedQuadMode¶
Type of quadrature; also used for location of nodes.
Values:
-
enumerator CEED_GAUSS¶
Gauss-Legendre quadrature.
-
enumerator CEED_GAUSS_LOBATTO¶
Gauss-Legendre-Lobatto quadrature.
-
enumerator CEED_GAUSS¶
-
enum CeedElemTopology¶
Type of basis shape to create non-tensor element basis.
Dimension can be extracted with bitwise AND (CeedElemTopology & 2**(dim + 2)) == TRUE
Values:
-
enumerator CEED_TOPOLOGY_LINE¶
Line.
-
enumerator CEED_TOPOLOGY_TRIANGLE¶
Triangle - 2D shape.
-
enumerator CEED_TOPOLOGY_QUAD¶
Quadralateral - 2D shape.
-
enumerator CEED_TOPOLOGY_TET¶
Tetrahedron - 3D shape.
-
enumerator CEED_TOPOLOGY_PYRAMID¶
Pyramid - 3D shape.
-
enumerator CEED_TOPOLOGY_PRISM¶
Prism - 3D shape.
-
enumerator CEED_TOPOLOGY_HEX¶
Hexehedron - 3D shape.
-
enumerator CEED_TOPOLOGY_LINE¶