FEValues
Main types
CellValues
and FacetValues
are the most common subtypes of Ferrite.AbstractValues
. For more details about how these work, please see the related topic guide.
Ferrite.CellValues
— TypeCellValues([::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 toFloat64
) to determine the type the internal data is stored as.quad_rule
: an instance of aQuadratureRule
func_interpol
: an instance of anInterpolation
used to interpolate the approximated functiongeom_interpol
: an optional instance of aInterpolation
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.
Keyword arguments: The following keyword arguments are experimental and may change in future minor releases
update_gradients
: Specifies if the gradients of the shape functions should be updated (default true)update_hessians
: Specifies if the hessians of the shape functions should be updated (default false)update_detJdV
: Specifies if the volume associated with each quadrature point should be updated (default true)
Common methods:
Ferrite.FacetValues
— TypeFacetValues([::Type{T}], quad_rule::FacetQuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
A FacetValues
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 facets of finite elements.
Arguments:
T
: an optional argument (default toFloat64
) to determine the type the internal data is stored as.quad_rule
: an instance of aFacetQuadratureRule
func_interpol
: an instance of anInterpolation
used to interpolate the approximated functiongeom_interpol
: an optional instance of anInterpolation
which is used to interpolate the geometry. By default linear Lagrange interpolation is used.
Keyword arguments: The following keyword arguments are experimental and may change in future minor releases
update_gradients
: Specifies if the gradients of the shape functions should be updated (default true)update_hessians
: Specifies if the hessians of the shape functions should be updated (default false)
Common methods:
Currently, embedded FEValues
returns SArray
s, which behave differently from the Tensor
s for normal value. In the future, we expect to return an AbstractTensor
, this change may happen in a minor release, and the API for embedded FEValues
should therefore be considered experimental.
Applicable functions
The following functions are applicable to both CellValues
and FacetValues
.
Ferrite.reinit!
— Functionreinit!(cv::CellValues, cell::AbstractCell, x::AbstractVector)
reinit!(cv::CellValues, x::AbstractVector)
reinit!(fv::FacetValues, cell::AbstractCell, x::AbstractVector, facet::Int)
reinit!(fv::FacetValues, x::AbstractVector, function_gradient::Int)
Update the CellValues
/FacetValues
object for a cell or facet with cell 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.
Ferrite.getnquadpoints
— Functiongetnquadpoints(fe_v::AbstractValues)
Return the number of quadrature points. For FacetValues
, this is the number for the current facet.
Ferrite.getdetJdV
— FunctiongetdetJdV(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 facet 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$
Ferrite.shape_value
— Methodshape_value(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the value of shape function base_function
evaluated in quadrature point q_point
.
Ferrite.shape_gradient
— Methodshape_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the gradient of shape function base_function
evaluated in quadrature point q_point
.
Ferrite.shape_symmetric_gradient
— Functionshape_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
.
Ferrite.shape_divergence
— Functionshape_divergence(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the divergence of shape function base_function
evaluated in quadrature point q_point
.
Ferrite.shape_curl
— Functionshape_curl(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the curl of shape function base_function
evaluated in quadrature point q_point
.
Ferrite.geometric_value
— Functiongeometric_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
.
Ferrite.function_value
— Functionfunction_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}$.
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 Vec
s (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}$.
Ferrite.function_gradient
— Functionfunction_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}$.
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 Vec
s (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}$.
Ferrite.function_symmetric_gradient
— Functionfunction_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.
Ferrite.function_divergence
— Functionfunction_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.
Ferrite.function_curl
— Functionfunction_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.
Ferrite.spatial_coordinate
— Functionspatial_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{\xi}) \mathbf{\hat{x}}_i$.
where $\xi$is the coordinate of the given quadrature point q_point
of the associated quadrature rule.
spatial_coordinate(ip::ScalarInterpolation, ξ::Vec, x::AbstractVector{<:Vec{sdim, T}})
Compute the spatial coordinate in a given 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{\xi}) \mathbf{\hat{x}}_i$
In addition, there are some methods that are unique for FacetValues
.
Ferrite.getcurrentfacet
— Functiongetcurrentfacet(fv::FacetValues)
Return the current active facet of the FacetValues
object (from last reinit!
).
Ferrite.getnormal
— Functiongetnormal(fv::FacetValues, qp::Int)
Return the normal at the quadrature point qp
for the active facet of the FacetValues
object(from last reinit!
).
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.
InterfaceValues
All of the methods for FacetValues
apply for InterfaceValues
as well. In addition, there are some methods that are unique for InterfaceValues
:
Ferrite.InterfaceValues
— TypeInterfaceValues
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::FacetQuadratureRule, ip::Interpolation)
: same quadrature rule and interpolation on both sides, default linear Lagrange geometric interpolation.InterfaceValues(qr::FacetQuadratureRule, ip::Interpolation, ip_geo::Interpolation)
: same as above but with given geometric interpolation.InterfaceValues(qr_here::FacetQuadratureRule, ip_here::Interpolation, qr_there::FacetQuadratureRule, ip_there::Interpolation)
: different quadrature rule and interpolation on the two sides, default linear Lagrange geometric interpolation.InterfaceValues(qr_here::FacetQuadratureRule, ip_here::Interpolation, ip_geo_here::Interpolation, qr_there::FacetQuadratureRule, ip_there::Interpolation, ip_geo_there::Interpolation)
: same as above but with given geometric interpolation.InterfaceValues(fv::FacetValues)
: quadrature rule and interpolations from face values (same on both sides).InterfaceValues(fv_here::FacetValues, fv_there::FacetValues)
: quadrature rule and interpolations from the face values.
Associated methods:
Common methods:
Ferrite.shape_value_average
— Functionshape_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.
Ferrite.shape_value_jump
— Functionshape_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 in the default normal direction.
This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} -\vec{v}^\text{here}$. To obtain the form, $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} \cdot \vec{n}^\text{there} + \vec{v}^\text{here} \cdot \vec{n}^\text{here}$, multiply by minus the outward facing normal to the first element's side of the interface (which is the default normal for getnormal
with InterfaceValues
).
Ferrite.shape_gradient_average
— Functionshape_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.
Ferrite.shape_gradient_jump
— Functionshape_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 in the default normal direction.
This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} -\vec{v}^\text{here}$. 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 minus the outward facing normal to the first element's side of the interface (which is the default normal for getnormal
with InterfaceValues
).
Ferrite.function_value_average
— Functionfunction_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.
Ferrite.function_value_jump
— Functionfunction_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 along the default normal direction.
This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} -\vec{v}^\text{here}$. 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 minus the outward facing normal to the first element's side of the interface (which is the default normal for getnormal
with InterfaceValues
).
Ferrite.function_gradient_average
— Functionfunction_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.
Ferrite.function_gradient_jump
— Functionfunction_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 along the default normal direction.
This function uses the definition $\llbracket \vec{v} \rrbracket=\vec{v}^\text{there} -\vec{v}^\text{here}$. 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 minus the outward facing normal to the first element's side of the interface (which is the default normal for getnormal
with InterfaceValues
).