Cell Values

FerriteInterfaceElements.InterfaceCellValuesType
InterfaceCellValues([::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 interface
  • here::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 base CellValues for each base function of the InterfaceCellValues
source
Ferrite.shape_valueMethod
shape_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".

source
Ferrite.shape_gradientMethod
shape_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".

source
Ferrite.function_valueFunction
function_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.

source
Ferrite.function_gradientFunction
function_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.

source
Ferrite.shape_value_averageFunction
shape_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.

source
Ferrite.shape_gradient_averageFunction
shape_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.

source
Ferrite.shape_value_jumpFunction
shape_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.

source
Ferrite.shape_gradient_jumpFunction
shape_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.

source
Ferrite.function_value_averageFunction
function_value_average(cv::InterfaceCellValues, qp::Int, u::AbstractVector)

Compute the average value of the function in a quadrature point.

source
Ferrite.function_gradient_averageFunction
function_gradient_average(cv::InterfaceCellValues, qp::Int, u::AbstractVector)

Compute the average gradient of the function in a quadrature point.

source
Ferrite.function_value_jumpFunction
function_value_jump(cv::InterfaceCellValues, qp::Int, u::AbstractVector)

Compute the jump of the function value in a quadrature point.

source
Ferrite.function_gradient_jumpFunction
function_gradient_jump(cv::InterfaceCellValues, qp::Int, u::AbstractVector)

Compute the jump of the function gradient in a quadrature point.

source
FerriteInterfaceElements.getdetJdV_averageFunction
getdetJdV_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.

source
FerriteInterfaceElements.get_side_and_baseindexFunction
get_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.

source