Quadrature
Ferrite.QuadratureRule
— TypeQuadratureRule{shape}([quad_rule_type::Symbol], order::Int)
QuadratureRule{shape, T}([quad_rule_type::Symbol], order::Int)
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.
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.FaceQuadratureRule
— TypeFaceQuadratureRule{shape}([quad_rule_type::Symbol], order::Int)
FaceQuadratureRule{shape, T}([quad_rule_type::Symbol], order::Int)
Create a FaceQuadratureRule
used for integration of the faces of 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.
FaceQuadratureRule
is used as one of the components to create FaceValues
.
Ferrite.getnquadpoints
— Methodgetnquadpoints(qr::QuadratureRule)
Return the number of quadrature points in qr
.
Ferrite.getnquadpoints
— Methodgetnquadpoints(qr::FaceQuadratureRule, face::Int)
Return the number of quadrature points in qr
for local face index face
.
Ferrite.getpoints
— Functiongetpoints(qr::QuadratureRule)
getpoints(qr::FaceQuadratureRule, face::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::FaceQuadratureRule, face::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