Quadrature

Ferrite.QuadratureRuleType
QuadratureRule{dim,shape}([quad_rule_type::Symbol], order::Int)

Create a QuadratureRule used for integration. dim is the space dimension, shape an AbstractRefShape and order 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 a CellValues or FaceValues object.

Common methods:

Example:

julia> QuadratureRule{2, RefTetrahedron}(1)
Ferrite.QuadratureRule{2,Ferrite.RefTetrahedron,Float64}([0.5], Tensors.Tensor{1,2,Float64,2}[[0.333333, 0.333333]])

julia> QuadratureRule{1, RefCube}(:lobatto, 2)
Ferrite.QuadratureRule{1,Ferrite.RefCube,Float64}([1.0, 1.0], Tensors.Tensor{1,1,Float64,1}[[-1.0], [1.0]])
source
Ferrite.AbstractRefShapeType

Represents a reference shape which quadrature rules and interpolations are defined on. Currently, the only concrete types that subtype this type are RefCube in 1, 2 and 3 dimensions, and RefTetrahedron in 2 and 3 dimensions.

source
Ferrite.getpointsFunction
getpoints(qr::QuadratureRule)

Return the points of the quadrature rule.

Examples

julia> qr = QuadratureRule{2, RefTetrahedron}(:legendre, 2);

julia> getpoints(qr)
3-element Array{Tensors.Tensor{1,2,Float64,2},1}:
 [0.166667, 0.166667]
 [0.166667, 0.666667]
 [0.666667, 0.166667]
source
Ferrite.getweightsFunction
getweights(qr::QuadratureRule)

Return the weights of the quadrature rule.

Examples

julia> qr = QuadratureRule{2, RefTetrahedron}(:legendre, 2);

julia> getweights(qr)
3-element Array{Float64,1}:
 0.166667
 0.166667
 0.166667
source