Quadrature
Ferrite.QuadratureRule — TypeQuadratureRule{shape}([::Type{T},] [quad_rule_type::Symbol,] order::Int)
QuadratureRule{shape}(weights::AbstractVector{T}, points::AbstractVector{Vec{rdim, T}})Create a QuadratureRule used for integration on the refshape shape (of type AbstractRefShape). order is the order of the quadrature rule. quad_rule_type is an optional argument determining the type of quadrature rule, currently the :legendre and :lobatto rules are implemented for hypercubes. For triangles up to order 8 the default rule is the one by :dunavant (see [7]) and for tetrahedra the default rule is keast_minimal (see [8]). Wedges and pyramids default to :polyquad (see [9]). Furthermore we have implemented
:gaussjacobifor triangles (order 9-15):keast_minimal(see [8]) for tetrahedra (order 1-5), containing negative weights:keast_positive(see [8]) for tetrahedra (order 1-5), containing only positive weights
A QuadratureRule is used to approximate an integral on a domain by a weighted sum of function values at specific points:
$\int\limits_\Omega f(\mathbf{x}) \text{d} \Omega \approx \sum\limits_{q = 1}^{n_q} f(\mathbf{x}_q) w_q$
The quadrature rule consists of $n_q$ points in space $\mathbf{x}_q$ with corresponding weights $w_q$.
In Ferrite, the QuadratureRule type is mostly used as one of the components to create CellValues.
Common methods:
getpoints: the points of the quadrature rulegetweights: the weights of the quadrature rule
Example:
julia> qr = QuadratureRule{RefTriangle}(1)
QuadratureRule{RefTriangle, Float64, 2}([0.5], Vec{2, Float64}[[0.33333333333333, 0.33333333333333]])
julia> getpoints(qr)
1-element Vector{Vec{2, Float64}}:
[0.33333333333333, 0.33333333333333]Ferrite.FacetQuadratureRule — TypeFacetQuadratureRule{shape}([::Type{T},] [quad_rule_type::Symbol,] order::Int)
FacetQuadratureRule{shape}(facet_rules::NTuple{<:Any, <:QuadratureRule{shape}})
FacetQuadratureRule{shape}(facet_rules::AbstractVector{<:QuadratureRule{shape}})Create a FacetQuadratureRule used for integration of the facets of the refshape shape (of type AbstractRefShape). order is the order of the quadrature rule. If no symbol is provided, the default quad_rule_type for each facet's reference shape is used (see QuadratureRule). For non-default quad_rule_types on cells with mixed facet types (e.g. RefPrism and RefPyramid), the facet_rules must be provided explicitly.
FacetQuadratureRule is used as one of the components to create FacetValues.
Ferrite.getnquadpoints — Methodgetnquadpoints(qr::QuadratureRule)Return the number of quadrature points in qr.
Ferrite.getnquadpoints — Methodgetnquadpoints(qr::FacetQuadratureRule, facet::Int)Return the number of quadrature points in qr for local facet index facet.
Ferrite.getpoints — Functiongetpoints(qr::QuadratureRule)
getpoints(qr::FacetQuadratureRule, facet::Int)Return the points of the quadrature rule.
Examples
julia> qr = QuadratureRule{RefTriangle}(:legendre, 2);
julia> getpoints(qr)
3-element Vector{Vec{2, Float64}}:
[0.16666666666667, 0.16666666666667]
[0.16666666666667, 0.66666666666667]
[0.66666666666667, 0.16666666666667]Ferrite.getweights — Functiongetweights(qr::QuadratureRule)
getweights(qr::FacetQuadratureRule, facet::Int)Return the weights of the quadrature rule.
Examples
julia> qr = QuadratureRule{RefTriangle}(:legendre, 2);
julia> getweights(qr)
3-element Array{Float64,1}:
0.166667
0.166667
0.166667