Intrepid2
Function space tools

Table of contents

Overview

Intrepid2::FunctionSpaceTools is a stateless class of expert methods for operations on finite element subspaces of $H(grad,\Omega)$, $H(curl,\Omega)$, $H(div,\Omega)$ and $L^2(\Omega)$. In Intrepid these spaces are referred to as HGRAD, HCURL, HDIV and HVOL. There are four basic groups of methods:

Pullbacks

Notation in this section follows the standard definition of a finite element space by Ciarlet; see The Finite Element Method for Elliptic Problems, Classics in Applied Mathematics, SIAM, 2002. Given a reference cell $\{\widehat{{\mathcal C}},\widehat{P},\widehat{\Lambda}\}$ with a basis $\{\widehat{u}_i\}_{i=0}^n$, the basis $\{{u}_i\}_{i=0}^n$ of $\{{\mathcal C},P,\Lambda\}$ is defined as follows:

\[ u_i = \sigma_i \Phi^*(\widehat{u}_i), \qquad i=1,\ldots,n \,. \]

In this formula $\{\sigma_i\}_{i=0}^n$, where $\sigma_i = \pm 1$, are the field signs, and $\Phi^*$ is the pullback ("change of variables") transformation. For scalar spaces such as HGRAD and HVOL the field signs are always equal to 1 and can be disregarded. For vector field spaces such as HCURL or HDIV, the field sign of a basis function can be +1 or -1, depending on the orientation of the physical edge or face, associated with the basis function.

The actual form of the pullback depends on which one of the four function spaces HGRAD, HCURL, HDIV and HVOL is being approximated and is computed as follows. Let $F_{\mathcal C}$ denote the reference-to-physical map (see Section Reference-to-physical cell mapping); $DF_{\mathcal C}$ is its Jacobian (see Section Jacobian of the reference-to-physical cell mapping) and $J_{\mathcal C} = \det(DF_{\mathcal C})$. Then,

\[ \begin{array}{ll} \Phi^*_G : HGRAD(\widehat{{\mathcal C}}) \mapsto HGRAD({\mathcal C})& \qquad \Phi^*_G(\widehat{u}) = \widehat{u}\circ F^{-1}_{\mathcal C} \\[2ex] \Phi^*_C : HCURL(\widehat{{\mathcal C}}) \mapsto HCURL({\mathcal C})& \qquad \Phi^*_C(\widehat{\bf u}) = \left((DF_{\mathcal C})^{-{\sf T}}\cdot\widehat{\bf u}\right)\circ F^{-1}_{\mathcal C} \\[2ex] \Phi^*_D : HDIV(\widehat{{\mathcal C}}) \mapsto HDIV({\mathcal C})& \qquad \Phi^*_D(\widehat{\bf u}) = \left(J^{-1}_{\mathcal C} DF_{\mathcal C}\cdot\widehat{\bf u}\right)\circ F^{-1}_{\mathcal C} \\[2ex] \Phi^*_S : HVOL(\widehat{{\mathcal C}}) \mapsto HVOL({\mathcal C})& \qquad \Phi^*_S(\widehat{u}) = \left(J^{-1}_{\mathcal C} \widehat{u}\right) \circ F^{-1}_{\mathcal C} \,. \end{array} \]

Intrepid supports pullbacks only for cell topologies that have reference cells; see Reference cells.

Measure

In Intrepid integrals of finite element functions over cells, 2-subcells (faces) and 1-subcells (edges) are computed by change of variables to reference frame and require three different kinds of measures.

  1. The integral of a scalar function over a cell ${\mathcal C}$

    \[ \int_{{\mathcal C}} f(x) dx = \int_{\widehat{{\mathcal C}}} f(F(\widehat{x})) |J | d\widehat{x} \]

    requires the volume measure defined by the determinant of the Jacobian. This measure is computed by Intrepid2::FunctionSpaceTools::computeCellMeasure
  2. The integral of a scalar function over 2-subcell $\mathcal{F}$

    \[ \int_{\mathcal{F}} f(x) dx = \int_{R} f(\Phi(u,v)) \left\|\frac{\partial\Phi}{\partial u}\times \frac{\partial\Phi}{\partial v}\right\| du\,dv \]

    requires the surface measure defined by the norm of the vector product of the surface tangents. This measure is computed by Intrepid2::FunctionSpaceTools::computeFaceMeasure. In this formula R is the parametrization domain for the 2-subcell; see Section Parametrization of physical 1- and 2-subcells for details.
  3. The integral of a scalar function over a 1-subcell $\mathcal{E}$

    \[ \int_{\mathcal{E}} f(x) dx = \int_{R} f(\Phi(s)) \|\Phi'\| ds \]

    requires the arc measure defined by the norm of the arc tangent vector. This measure is computed by Intrepid2::FunctionSpaceTools::computeEdgeMeasure. In this formula R is the parametrization domain for the 1-subcell; see Section Parametrization of physical 1- and 2-subcells for details.

Evaluation of finite element fields

To make this example more specific, assume curl-conforming finite element spaces. Suppose that we have a physical cell $\{{\mathcal C},P,\Lambda\}$ with a basis $\{{\bf u}_i\}_{i=0}^n$. A finite element function on this cell is defined by a set of n coefficients $\{c_i\}_{i=0}^n$:

\[ {\bf u}^h(x) = \sum_{i=0}^n c_i {\bf u}_i(x) \,. \]

From Section Pullbacks it follows that

\[ {\bf u}^h(x) = \sum_{i=0}^n c_i \sigma_i \left((DF_{\mathcal C})^{-{\sf T}}\cdot\widehat{\bf u}_i\right)\circ F^{-1}_{\mathcal C}(x) = \sum_{i=0}^n c_i \sigma_i (DF_{\mathcal C}(\widehat{x}))^{-{\sf T}}\cdot\widehat{\bf u}_i(\widehat{x})\,, \]

where $ \widehat{x} = F^{-1}_{\mathcal C}(x) \in \widehat{\mathcal C} $ is the pre-image of x in the reference cell.

Consequently, evaluation of finite element functions at a given set of points $\{x_p\}_{p=0}^P \subset {\mathcal C}$ comprises of the following four steps:

  1. Application of the inverse map $F^{-1}_{\mathcal C}$ to obtain the pre-images $\{\widehat{x}_p\}_{p=0}^P$ of the evaluation points in the reference cell $\widehat{\mathcal{C}}$; see Intrepid2::CellTools::mapToReferenceFrame
  2. Evaluation of the appropriate reference basis set $\{\widehat{\bf u}_i\}_{i=1}^n$ at the pre-image set $\{\widehat{x}_p\}_{p=0}^P$; see Intrepid2::Basis::getValues
  3. Application of the appropriate transformation and field signs. In our example the finite element space is curl-conforming and the appropriate transformation is implemented in Intrepid2::FunctionSpaceTools::HCURLtransformVALUE. Application of the signs to the transformed functions is done by Intrepid2::FunctionSpaceTools::applyFieldSigns.
  4. The final step is to compute the sum of the transformed and signed basis function values multiplied by the coefficients of the finite element function using Intrepid2::FunctionSpaceTools::evaluate.

Evaluation of adimssible derivatives of finite element functions is completely analogous and follows the same four steps. Evaluation of scalar finite element functions is simpler because application of the signes can be skipped for these functions.

Evaluation of finite element operators and functionals

Assume the same setting as in Section Evaluation of finite element fields. A finite element operator defined by the finite element basis on the physical cell $\mathcal{C}$ is a matrix

\[ \mathbf{K}^{\mathcal{C}}_{i,j} = \int_{\mathcal C} {\mathcal L}_L {\bf u}_i(x)\, {\mathcal L}_R {\bf u}_j(x) \, dx \,. \]

where ${\mathcal L}_L$ and ${\mathcal L}_R $ are left and right operators acting on the basis functions. Typically, when the left and the right basis functions are from the same finite element basis (as in this example), the left and right operators are the same. If they are set to VALUE we get a mass matrix; if they are set to an admissible differential operator we get a stiffnesss matrix. Assume again that the basis is curl-conforming and the operators are set to VALUE. Using the basis definition from Section Pullbacks we have that

\[ \mathbf{K}^{\mathcal{C}}_{i,j} = \int_{\widehat{\mathcal C}} \sigma_i \sigma_j (DF_{\mathcal C}(\widehat{x}))^{-{\sf T}}\cdot\widehat{\bf u}_i(\widehat{x})\cdot (DF_{\mathcal C}(\widehat{x}))^{-{\sf T}}\cdot\widehat{\bf u}_i(\widehat{x})\,d\widehat{x} \]

It follows that

\[ \mathbf{K}^{\mathcal{C}}_{i,j} = \mbox{diag}(\sigma_0,\ldots,\sigma_n)\widehat{\mathbf{K}}^{\mathcal{C}}\mbox{diag}(\sigma_0,\ldots,\sigma_n) \]

where

\[ \widehat{\mathbf{K}}^{\mathcal{C}}_{i,j} = \int_{\widehat{\mathcal C}} (DF_{\mathcal C}(\widehat{x}))^{-{\sf T}}\cdot\widehat{\bf u}_i(\widehat{x})\cdot (DF_{\mathcal C}(\widehat{x}))^{-{\sf T}}\cdot\widehat{\bf u}_i(\widehat{x})\,d\widehat{x} \]

is the raw cell operator matrix. The methods Intrepid2::FunctionSpaceTools::applyLeftFieldSigns and Intrepid2::FunctionSpaceTools::applyRightFieldSigns apply the left and right diagonal sign matrices to the raw cell operator.

A finite element operator defined by the finite element basis on the physical cell is a vector

\[ \mathbf{f}^{\mathcal{C}}_{i} = \int_{\mathcal C} f(x) {\mathcal L}_R u_i(x) \, dx \,. \]

Assuming again operator VALUE and using the same arguments as above, we see that

\[ \mathbf{f}^{\mathcal{C}} = \mbox{diag}(\sigma_0,\ldots,\sigma_n)\widehat{\mathbf{f}}^{\mathcal{C}}\,, \]

where

\[ \widehat{\mathbf{f}}^{\mathcal{C}} = \int_{\widehat{\mathcal C}} \mathbf{f}\circ F_{\mathcal C}(\widehat{x}) (DF_{\mathcal C}(\widehat{x}))^{-{\sf T}}\cdot\widehat{\bf u}_i(\widehat{x})\,d\widehat{x} \]

is the raw cell functional.