FEValues

Main types

CellValues and FaceValues are the most common subtypes of Ferrite.AbstractValues. For more details about how these work, please see the related topic guide.

Ferrite.CellValuesType
CellValues([::Type{T},] quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A CellValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. in the finite element cell.

Arguments:

  • T: an optional argument (default to Float64) to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of a Interpolation which is used to interpolate the geometry. By default linear Lagrange interpolation is used. For embedded elements the geometric interpolations should be vectorized to the spatial dimension.

Common methods:

source
Ferrite.FaceValuesType
FaceValues([::Type{T}], quad_rule::FaceQuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A FaceValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. on the faces of finite elements.

Arguments:

  • T: an optional argument (default to Float64) to determine the type the internal data is stored as.
  • quad_rule: an instance of a FaceQuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of an Interpolation which is used to interpolate the geometry. By default linear Lagrange interpolation is used.

Common methods:

source

Applicable functions

The following functions are applicable to both CellValues and FaceValues.

Ferrite.reinit!Function
reinit!(cv::CellValues, cell::AbstractCell, x::Vector)
reinit!(cv::CellValues, x::Vector)
reinit!(fv::FaceValues, cell::AbstractCell, x::Vector, face::Int)
reinit!(fv::FaceValues, x::Vector, face::Int)

Update the CellValues/FaceValues object for a cell or face with coordinates x. The derivatives of the shape functions, and the new integration weights are computed. For interpolations with non-identity mappings, the current cell is also required.

source
Ferrite.getnquadpointsFunction
getnquadpoints(fe_v::AbstractValues)

Return the number of quadrature points. For FaceValues, this is the number for the current face.

source
Ferrite.getdetJdVFunction
getdetJdV(fe_v::AbstractValues, q_point::Int)

Return the product between the determinant of the Jacobian and the quadrature point weight for the given quadrature point: $\det(J(\mathbf{x})) w_q$.

This value is typically used when one integrates a function on a finite element cell or face as

$\int\limits_\Omega f(\mathbf{x}) d \Omega \approx \sum\limits_{q = 1}^{n_q} f(\mathbf{x}_q) \det(J(\mathbf{x})) w_q$ $\int\limits_\Gamma f(\mathbf{x}) d \Gamma \approx \sum\limits_{q = 1}^{n_q} f(\mathbf{x}_q) \det(J(\mathbf{x})) w_q$

source
Ferrite.shape_valueMethod
shape_value(fe_v::AbstractValues, q_point::Int, base_function::Int)

Return the value of shape function base_function evaluated in quadrature point q_point.

source
Ferrite.shape_gradientMethod
shape_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)

Return the gradient of shape function base_function evaluated in quadrature point q_point.

source
Ferrite.shape_symmetric_gradientFunction
shape_symmetric_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)

Return the symmetric gradient of shape function base_function evaluated in quadrature point q_point.

source
Ferrite.shape_divergenceFunction
shape_divergence(fe_v::AbstractValues, q_point::Int, base_function::Int)

Return the divergence of shape function base_function evaluated in quadrature point q_point.

source
Ferrite.shape_curlFunction
shape_curl(fe_v::AbstractValues, q_point::Int, base_function::Int)

Return the curl of shape function base_function evaluated in quadrature point q_point.

source
Ferrite.geometric_valueFunction
geometric_value(fe_v::AbstractValues, q_point, base_function::Int)

Return the value of the geometric shape function base_function evaluated in quadrature point q_point.

source
Ferrite.function_valueFunction
function_value(iv::InterfaceValues, q_point::Int, u; here::Bool)
function_value(iv::InterfaceValues, q_point::Int, u, dof_range_here, dof_range_there; here::Bool)

Compute the value of the function in quadrature point q_point on the "here" (here=true) or "there" (here=false) side of the interface. u_here and u_there are the values of the degrees of freedom for the respective element.

u is a vector of scalar values for the degrees of freedom. This function can be used with a single u vector containing the dofs of both elements of the interface or two vectors (u_here and u_there) which contain the dofs of each cell of the interface respectively.

here determines which element to use for calculating function value. true uses the value on the first element's side of the interface, while false uses the value on the second element's side.

The value of a scalar valued function is computed as $u(\mathbf{x}) = \sum\limits_{i = 1}^n N_i (\mathbf{x}) u_i$ where $u_i$ are the value of $u$ in the nodes. For a vector valued function the value is calculated as $\mathbf{u}(\mathbf{x}) = \sum\limits_{i = 1}^n N_i (\mathbf{x}) \mathbf{u}_i$ where $\mathbf{u}_i$ are the nodal values of $\mathbf{u}$.

source
function_value(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])

Compute the value of the function in a quadrature point. u is a vector with values for the degrees of freedom. For a scalar valued function, u contains scalars. For a vector valued function, u can be a vector of scalars (for use of VectorValues) or u can be a vector of Vecs (for use with ScalarValues).

The value of a scalar valued function is computed as $u(\mathbf{x}) = \sum\limits_{i = 1}^n N_i (\mathbf{x}) u_i$ where $u_i$ are the value of $u$ in the nodes. For a vector valued function the value is calculated as $\mathbf{u}(\mathbf{x}) = \sum\limits_{i = 1}^n N_i (\mathbf{x}) \mathbf{u}_i$ where $\mathbf{u}_i$ are the nodal values of $\mathbf{u}$.

source
Ferrite.function_gradientFunction
function_gradient(iv::InterfaceValues, q_point::Int, u; here::Bool)
function_gradient(iv::InterfaceValues, q_point::Int, u, dof_range_here, dof_range_there; here::Bool)

Compute the gradient of the function in a quadrature point. u is a vector of scalar values for the degrees of freedom. This function can be used with a single u vector containing the dofs of both elements of the interface or two vectors (u_here and u_there) which contain the dofs of each cell of the interface respectively.

here determines which element to use for calculating function value. true uses the value on the first element's side of the interface, while false uses the value on the second element's side.

The gradient of a scalar function or a vector valued function with use of VectorValues is computed as $\mathbf{\nabla} u(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{\nabla} N_i (\mathbf{x}) u_i$ or $\mathbf{\nabla} u(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{\nabla} \mathbf{N}_i (\mathbf{x}) u_i$ respectively, where $u_i$ are the nodal values of the function. For a vector valued function with use of ScalarValues the gradient is computed as $\mathbf{\nabla} \mathbf{u}(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{u}_i \otimes \mathbf{\nabla} N_i (\mathbf{x})$ where $\mathbf{u}_i$ are the nodal values of $\mathbf{u}$.

source
function_gradient(fe_v::AbstractValues{dim}, q_point::Int, u::AbstractVector, [dof_range])

Compute the gradient of the function in a quadrature point. u is a vector with values for the degrees of freedom. For a scalar valued function, u contains scalars. For a vector valued function, u can be a vector of scalars (for use of VectorValues) or u can be a vector of Vecs (for use with ScalarValues).

The gradient of a scalar function or a vector valued function with use of VectorValues is computed as $\mathbf{\nabla} u(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{\nabla} N_i (\mathbf{x}) u_i$ or $\mathbf{\nabla} u(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{\nabla} \mathbf{N}_i (\mathbf{x}) u_i$ respectively, where $u_i$ are the nodal values of the function. For a vector valued function with use of ScalarValues the gradient is computed as $\mathbf{\nabla} \mathbf{u}(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{u}_i \otimes \mathbf{\nabla} N_i (\mathbf{x})$ where $\mathbf{u}_i$ are the nodal values of $\mathbf{u}$.

source
Ferrite.function_symmetric_gradientFunction
function_symmetric_gradient(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])

Compute the symmetric gradient of the function, see function_gradient. Return a SymmetricTensor.

The symmetric gradient of a scalar function is computed as $\left[ \mathbf{\nabla} \mathbf{u}(\mathbf{x_q}) \right]^\text{sym} = \sum\limits_{i = 1}^n \frac{1}{2} \left[ \mathbf{\nabla} N_i (\mathbf{x}_q) \otimes \mathbf{u}_i + \mathbf{u}_i \otimes \mathbf{\nabla} N_i (\mathbf{x}_q) \right]$ where $\mathbf{u}_i$ are the nodal values of the function.

source
Ferrite.function_divergenceFunction
function_divergence(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])

Compute the divergence of the vector valued function in a quadrature point.

The divergence of a vector valued functions in the quadrature point $\mathbf{x}_q)$ is computed as $\mathbf{\nabla} \cdot \mathbf{u}(\mathbf{x_q}) = \sum\limits_{i = 1}^n \mathbf{\nabla} N_i (\mathbf{x_q}) \cdot \mathbf{u}_i$ where $\mathbf{u}_i$ are the nodal values of the function.

source
Ferrite.function_curlFunction
function_curl(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])

Compute the curl of the vector valued function in a quadrature point.

The curl of a vector valued functions in the quadrature point $\mathbf{x}_q)$ is computed as $\mathbf{\nabla} \times \mathbf{u}(\mathbf{x_q}) = \sum\limits_{i = 1}^n \mathbf{\nabla} N_i \times (\mathbf{x_q}) \cdot \mathbf{u}_i$ where $\mathbf{u}_i$ are the nodal values of the function.

source
Ferrite.spatial_coordinateFunction
spatial_coordinate(fe_v::AbstractValues, q_point::Int, x::AbstractVector)

Compute the spatial coordinate in a quadrature point. x contains the nodal coordinates of the cell.

The coordinate is computed, using the geometric interpolation, as $\mathbf{x} = \sum\limits_{i = 1}^n M_i (\mathbf{x}) \mathbf{\hat{x}}_i$

source

In addition, there are some methods that are unique for FaceValues.

Ferrite.getcurrentfaceFunction
getcurrentface(fv::FaceValues)

Return the current active face of the FaceValues object (from last reinit!).

source
Ferrite.getnormalFunction
getnormal(fv::FaceValues, qp::Int)

Return the normal at the quadrature point qp for the active face of the FaceValues object(from last reinit!).

source
getnormal(iv::InterfaceValues, qp::Int; here::Bool=true)

Return the normal vector in the quadrature point qp on the interface. If here = true (default) the outward normal to the "here" element is returned, otherwise the outward normal to the "there" element.

source

InterfaceValues

All of the methods for FaceValues apply for InterfaceValues as well. In addition, there are some methods that are unique for InterfaceValues:

Ferrite.InterfaceValuesType
InterfaceValues

An InterfaceValues object facilitates the process of evaluating values, averages, jumps and gradients of shape functions and function on the interfaces between elements.

The first element of the interface is denoted "here" and the second element "there".

Constructors

  • InterfaceValues(qr::FaceQuadratureRule, ip::Interpolation): same quadrature rule and interpolation on both sides, default linear Lagrange geometric interpolation.
  • InterfaceValues(qr::FaceQuadratureRule, ip::Interpolation, ip_geo::Interpolation): same as above but with given geometric interpolation.
  • InterfaceValues(qr_here::FaceQuadratureRule, ip_here::Interpolation, qr_there::FaceQuadratureRule, ip_there::Interpolation): different quadrature rule and interpolation on the two sides, default linear Lagrange geometric interpolation.
  • InterfaceValues(qr_here::FaceQuadratureRule, ip_here::Interpolation, ip_geo_here::Interpolation, qr_there::FaceQuadratureRule, ip_there::Interpolation, ip_geo_there::Interpolation): same as above but with given geometric interpolation.
  • InterfaceValues(fv::FaceValues): quadrature rule and interpolations from face values (same on both sides).
  • InterfaceValues(fv_here::FaceValues, fv_there::FaceValues): quadrature rule and interpolations from the face values.

Associated methods:

Common methods:

source
Ferrite.shape_value_averageFunction
shape_value_average(iv::InterfaceValues, qp::Int, i::Int)

Compute the average of the value of shape function i at quadrature point qp across the interface.

source
Ferrite.shape_value_jumpFunction
shape_value_jump(iv::InterfaceValues, qp::Int, i::Int)

Compute the jump of the value of shape function i at quadrature point qp across the interface.

This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{here} -\vec{v}^\text{there}$. to obtain the form $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} ⋅ \vec{n}^\text{there} + \vec{v}^\text{here} ⋅ \vec{n}^\text{here}$ multiply by the outward facing normal to the first element's side of the interface (which is the default normal for getnormal with InterfaceValues).

source
Ferrite.shape_gradient_averageFunction
shape_gradient_average(iv::InterfaceValues, qp::Int, i::Int)

Compute the average of the gradient of shape function i at quadrature point qp across the interface.

source
Ferrite.shape_gradient_jumpFunction
shape_gradient_jump(iv::InterfaceValues, qp::Int, i::Int)

Compute the jump of the gradient of shape function i at quadrature point qp across the interface.

This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{here} -\vec{v}^\text{there}$. to obtain the form $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} ⋅ \vec{n}^\text{there} + \vec{v}^\text{here} ⋅ \vec{n}^\text{here}$ multiply by the outward facing normal to the first element's side of the interface (which is the default normal for getnormal with InterfaceValues).

source
Ferrite.function_value_averageFunction
function_value_average(iv::InterfaceValues, q_point::Int, u)
function_value_average(iv::InterfaceValues, q_point::Int, u, dof_range_here, dof_range_there)

Compute the average of the function value at the quadrature point on the interface.

source
Ferrite.function_value_jumpFunction
function_value_jump(iv::InterfaceValues, q_point::Int, u)
function_value_jump(iv::InterfaceValues, q_point::Int, u, dof_range_here, dof_range_there)

Compute the jump of the function value at the quadrature point over the interface.

This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{here} -\vec{v}^\text{there}$. to obtain the form $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} ⋅ \vec{n}^\text{there} + \vec{v}^\text{here} ⋅ \vec{n}^\text{here}$ multiply by the outward facing normal to the first element's side of the interface (which is the default normal for getnormal with InterfaceValues).

source
Ferrite.function_gradient_averageFunction
function_gradient_average(iv::InterfaceValues, q_point::Int, u)
function_gradient_average(iv::InterfaceValues, q_point::Int, u, dof_range_here, dof_range_there)

Compute the average of the function gradient at the quadrature point on the interface.

source
Ferrite.function_gradient_jumpFunction
function_gradient_jump(iv::InterfaceValues, q_point::Int, u)
function_gradient_jump(iv::InterfaceValues, q_point::Int, u, dof_range_here, dof_range_there)

Compute the jump of the function gradient at the quadrature point over the interface.

This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{here} -\vec{v}^\text{there}$. to obtain the form $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} ⋅ \vec{n}^\text{there} + \vec{v}^\text{here} ⋅ \vec{n}^\text{here}$ multiply by the outward facing normal to the first element's side of the interface (which is the default normal for getnormal with InterfaceValues).

source