Cell Values
FerriteInterfaceElements.InterfaceCellValues — TypeInterfaceCellValues([::Type{T},] qr::QuadratureRule, func_ip::InterfaceCellInterpolation, [geom_ip::InterfaceCellInterpolation]; use_same_cv=true)An InterfaceCellValues is based on two CellValues: one for each facet of an InterfaceCell. Since one can use the same CellValues for both sides, be default the same object is used for better performance. The keyword argument use_same_cv can be set to false to disable this behavior, if needed.
Fields
ip::InterfaceCellInterpolation: interpolation on the interfacehere::CellValues: values for facet "here"there::CellValues: values for facet "there"base_indices_here::Vector{Int}: base function indices on facet "here"base_indices_there::Vector{Int}: base function indices on facet "there"sides_and_baseindices::Tuple: side and base function for the baseCellValuesfor each base function of theInterfaceCellValues
Ferrite.shape_value — Methodshape_value(cv::InterfaceCellValues, qp::Int, i::Int, here::Bool)Return the value of shape function i evaluated in quadrature point qp on side here, where true means "here" and false means "there".
Ferrite.shape_gradient — Methodshape_gradient(cv::InterfaceCellValues, qp::Int, i::Int, here::Bool)Return the gradient of shape function i evaluated in quadrature point qp on side here, where true means "here" and false means "there".
Ferrite.function_value — Functionfunction_value(cv::InterfaceCellValues, qp::Int, u::AbstractVector, here::Bool)Compute the value of the function in a quadrature point on side here, where true means "here" and false means "there". u is a vector with values for the degrees of freedom.
Ferrite.function_gradient — Functionfunction_gradient(cv::InterfaceCellValues, qp::Int, u::AbstractVector, here::Bool)Compute the gradient of the function in a quadrature point on side here, where true means "here" and false means "there". u is a vector with values for the degrees of freedom.
Ferrite.shape_value_average — Functionshape_value_average(cv::InterfaceCellValues, qp::Int, i::Int)Return the value of shape function i evaluated in quadrature point qp for computing the average value on an interface.
Ferrite.shape_gradient_average — Functionshape_gradient_average(cv::InterfaceCellValues, qp::Int, i::Int)Return the gradient of shape function i evaluated in quadrature point qp for computing the average gradient on an interface.
Ferrite.shape_value_jump — Functionshape_value_jump(cv::InterfaceCellValues, qp::Int, i::Int)Return the value of shape function i evaluated in quadrature point qp for computing the value jump on an interface.
Ferrite.shape_gradient_jump — Functionshape_gradient_jump(cv::InterfaceCellValues, qp::Int, i::Int)Return the gradient of shape function i evaluated in quadrature point qp for computing the gradient jump on an interface.
Ferrite.function_value_average — Functionfunction_value_average(cv::InterfaceCellValues, qp::Int, u::AbstractVector)Compute the average value of the function in a quadrature point.
Ferrite.function_gradient_average — Functionfunction_gradient_average(cv::InterfaceCellValues, qp::Int, u::AbstractVector)Compute the average gradient of the function in a quadrature point.
Ferrite.function_value_jump — Functionfunction_value_jump(cv::InterfaceCellValues, qp::Int, u::AbstractVector)Compute the jump of the function value in a quadrature point.
Ferrite.function_gradient_jump — Functionfunction_gradient_jump(cv::InterfaceCellValues, qp::Int, u::AbstractVector)Compute the jump of the function gradient in a quadrature point.
FerriteInterfaceElements.getdetJdV_average — FunctiongetdetJdV_average(cv::InterfaceCellValues, qp::Int)Return the average of the product between the determinant of the Jacobian on each side of the interface and the quadrature point weight for the given quadrature point: $\det(J(\mathbf{x})) w_q$.
This value is typically used when integrating a function on the mid-plane of an interface element.
FerriteInterfaceElements.get_side_and_baseindex — Functionget_side_and_baseindex(cv::InterfaceCellInterpolation, i::Integer)
get_side_and_baseindex(cv::InterfaceCellValues, i::Integer)For an InterfaceCellInterpolation: given the base function index i return the side (:here or :there) and the base function index locally on that side.