FEValues
CellValues
Ferrite.CellValues
— TypeCellScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
CellVectorValues([::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. There are two different types of CellValues
: CellScalarValues
and CellVectorValues
. As the names suggest, CellScalarValues
utilizes scalar shape functions and CellVectorValues
utilizes vectorial shape functions. For a scalar field, the CellScalarValues
type should be used. For vector field, both subtypes can be used.
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
Common methods:
Ferrite.reinit!
— Functionreinit!(cv::CellValues, x::Vector)
reinit!(bv::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.
Ferrite.getnquadpoints
— Functiongetnquadpoints(fe_v::Values)
Return the number of quadrature points for the Values
object.
Ferrite.getdetJdV
— FunctiongetdetJdV(fe_v::Values, 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$
Ferrite.shape_value
— Functionshape_value(fe_v::Values, q_point::Int, base_function::Int)
Return the value of shape function base_function
evaluated in quadrature point q_point
.
Ferrite.shape_gradient
— Functionshape_gradient(fe_v::Values, 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::Values, 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::Values, q_point::Int, base_function::Int)
Return the divergence of shape function base_function
evaluated in quadrature point q_point
.
Ferrite.function_value
— Functionfunction_value(fe_v::Values, q_point::Int, u::AbstractVector)
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(fe_v::Values{dim}, q_point::Int, u::AbstractVector)
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::Values, q_point::Int, u::AbstractVector)
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::Values, q_point::Int, u::AbstractVector)
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.spatial_coordinate
— Functionspatial_coordinate(fe_v::Values{dim}, 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$
FaceValues
All of the methods for CellValues
apply for FaceValues
as well. In addition, there are some methods that are unique for FaecValues
:
Ferrite.FaceValues
— TypeFaceScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
FaceVectorValues([::Type{T}], quad_rule::QuadratureRule, 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. There are two different types of FaceValues
: FaceScalarValues
and FaceVectorValues
. As the names suggest, FaceScalarValues
utilizes scalar shape functions and FaceVectorValues
utilizes vectorial shape functions. For a scalar field, the FaceScalarValues
type should be used. For vector field, both subtypes can be used.
The quadrature rule for the face should be given with one dimension lower. I.e. for a 3D case, the quadrature rule should be in 2D.
Arguments:
T
: an optional argument 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 anInterpolation
which is used to interpolate the geometry
Common methods:
Ferrite.getcurrentface
— Functiongetcurrentface(fv::FaceValues)
Return the current active face of the FaceValues
object (from last reinit!
).