Integration

Measure

Bcube.MeasureType

A Measure is geometrical domain of integration associated to a way to integrate on it (i.e a quadrature rule).

Q is the quadrature type used to integrate expressions using this measure.

source
Bcube.MeasureMethod
Measure(domain::AbstractDomain, degree::Integer)
Measure(domain::AbstractDomain, ::Val{degree}) where {degree}

Build a Measure on the designated AbstractDomain with a default quadrature of degree degree.

Arguments

  • domain::AbstractDomain : the domain to integrate over
  • degree : the degree of the quadrature rule (Legendre quadrature type by default)

Examples

julia> mesh = line_mesh(10)
julia> Ω = CellDomain(mesh)
julia> dΩ = Measure(Ω, 2)
source

Integration methods

Base.:*Method
*(a::Number, b::Integration)
*(a::Integration, b::Number)

Multiplication of an Integration is based on a rewriting rule following the linearity rules of integration : k*∫(f(x))dx => ∫(k*f(x))dx

source
Base.:-Method
-(a::Integrand)

Soustraction on an Integrand is treated as a multiplication by "(-1)" : -a ≡ ((-1)*a)

source
Base.:-Method
-(a::Integration)

Soustraction on an Integration is treated as a multiplication by "(-1)" : -a ≡ ((-1)*a)

source
Bcube.apply_quadratureMethod
apply_quadrature(
    g_ref,
    shape::AbstractShape,
    quadrature::AbstractQuadrature,
    ::MapComputeQuadratureStyle,
)

Apply quadrature rule to function g_ref expressed on reference shape shape. Computation is optimized according to the given concrete type T<:AbstractComputeQuadratureStyle.

source
Bcube.apply_quadrature_v2Method
apply_quadrature_v2(
    g_ref,
    shape::AbstractShape,
    quadrature::AbstractQuadrature,
    ::MapComputeQuadratureStyle,
)

Alternative version of apply_quadrature thats seems to be more efficient for face integration (this observation is not really understood)

source
Bcube.integrate_on_ref_elementMethod
integrate_on_ref_element(g, cellinfo::CellInfo, quadrature::AbstractQuadrature, [::T]) where {N,[T<:AbstractComputeQuadratureStyle]}
integrate_on_ref_element(g, faceinfo::FaceInfo, quadrature::AbstractQuadrature, [::T]) where {N,[T<:AbstractComputeQuadratureStyle]}

Integrate a function g over a cell/face described by cellinfo/faceinfo.

The function g can be expressed in the reference or the physical space corresponding to the cell, both cases are automatically handled by applying necessary mapping when needed.

If the last argument is given, computation is optimized according to the given concrete type T<:AbstractComputeQuadratureStyle.

source

Quadrature rules

Bcube.AbstractQuadratureNodeType
AbstractQuadratureNode{S,Q}

Abstract type representing a quadrature node for a shape S and a quadrature Q. This type is used to represent and identify easily a quadrature node in a quadrature rules.

Derived types must implement the following method:

- get_index(quadnode::AbstractQuadratureNode{S,Q})
- get_coords(quadnode::AbstractQuadratureNode)
source
Bcube.AbstractQuadratureRuleType
AbstractQuadratureRule{S,Q}

Abstract type representing a quadrature rule for a shape S and quadrature Q.

Derived types must implement the following method: - [get_weights(qr::AbstractQuadratureRule)] - [get_nodes(qr::AbstractQuadratureRule)] - [Base.length(qr::AbstractQuadratureRule)]

source
Bcube.QuadratureNodeType
QuadratureNode{S,Q}

Type representing a quadrature node for a shape S and a quadrature Q. This type can be used to represent and identify easily a quadrature node in the corresponding parent quadrature rule.

source
Bcube.QuadratureRuleMethod
QuadratureRule(shape::AbstractShape, q::AbstractQuadrature)

Return the quadrature rule corresponding to the given Shape according to the selected quadrature 'q'.

source
Base.lengthMethod
length(qr::AbstractQuadratureRule)

Returns the number of quadrature nodes of qr.

source
Bcube._gausslegendre1DMethod
_gausslegendre1D(::Val{N}) where N
_gausslobatto1D(::Val{N}) where N

Return N-point Gauss quadrature weights and nodes on the domain [-1:1].

source
Bcube._get_num_nodesMethod

Gauss-Legendre formula with $n$ nodes has degree of exactness $2n-1$. Then, to obtain a given degree $D$, the number of nodes must satisfy: $2n-1 ≥ D$ or equivalently $n ≥ (D+1)/2$

source
Bcube._get_num_nodesMethod

Gauss-Lobatto formula with $n+1$ nodes has degree of exactness $2n-1$, which equivalent to a degree of $2n-3$ with $n$ nodes. Then, to obtain a given degree $D$, the number of nodes must satisfy: $2n-3 ≥ D$ or equivalently $n ≥ (D+3)/2$

source
Bcube._get_num_nodes_per_dimMethod
get_num_nodes_per_dim(quadrule::AbstractQuadratureRule{S}) where S<:Shape

Returns the number of nodes per dimension. This function is defined for shapes for which quadratures are based on a cartesian product : Line, Square, Cube

Remark : Here we assume that the same degree is used along each dimension (no anisotropy for now!)

source
Bcube.evalquadnodeMethod
evalquadnode(f, quadnode::AbstractQuadratureNode)

Evaluate the function f at the coordinates of quadnode.

Basically, it computes:

f(get_coords(quadnode))

Remark:

Optimization could be applied if f is a function based on a nodal basis such as one of the DoF and quadnode are collocated.

source
Bcube.get_coordsMethod
get_coords(quadnode::AbstractQuadratureNode)

Returns the coordinates of quadnode.

source
Bcube.get_indexMethod
get_index(quadnode::AbstractQuadratureNode{S,Q})

Returns the index of quadnode in the parent quadrature rule AbstractQuadRules{S,Q}

source
Bcube.get_nodesMethod
get_nodes(qr::AbstractQuadratureRule)

Returns an array containing the coordinates of all quadrature nodes of qr.

source
Bcube.get_quadnodesMethod
get_quadnodes(qr::QuadratureRule{S,Q,N}) where {S,Q,N}

Returns an vector containing each QuadratureNode of qr

source
Bcube.get_quadrature_points_gausslegendreMethod

Gauss-Legendre quadrature, 12 point rule on triangle.

Ref: Witherden, F. D.; Vincent, P. E. On the identification of symmetric quadrature rules for finite element methods. Comput. Math. Appl. 69 (2015), no. 10, 1232–1241

source
Bcube.get_quadrature_points_gausslegendreMethod

Gauss-Legendre quadrature, 16 point rule on triangle, degree 8.

Ref: Witherden, F. D.; Vincent, P. E. On the identification of symmetric quadrature rules for finite element methods. Comput. Math. Appl. 69 (2015), no. 10, 1232–1241

Note : quadrature is rescale to match our reference triangular shape which is defined in [0:1]² instead of [-1:1]²

source
Bcube.get_weightsMethod
get_weights(qr::AbstractQuadratureRule)

Returns an array containing the weights of all quadrature nodes of qr.

source
Bcube.quadrature_pointsMethod

ref : https://www.math.umd.edu/~tadmor/references/files/Chen%20&%20Shu%20entropy%20stable%20DG%20JCP2017.pdf

source
Bcube.quadrature_ruleMethod
quadrature_rule(::AbstractShape, ::Val{D}, ::AbstractQuadratureType) where {D}

Return the weights and points associated to the quadrature rule of degree D for the desired shape and quadrature type.

source
Bcube.quadrature_ruleMethod
quadrature_rule(iside::Int, shape::AbstractShape, degree::Val{N}) where N

Return the quadrature rule, computed with barycentric coefficients, corresponding to the given boundary of a shape and the given degree.

source
Bcube.quadrature_rule_baryMethod
quadrature_rule_bary(::Int, ::AbstractShape, degree::Val{N}) where N

Return the quadrature rule, computed with barycentric coefficients, corresponding to the given boundary of a shape and the given degree.

This function returns the quadrature weights and the barycentric weights to apply to each vertex of the reference shape. Hence, to apply the quadrature using this function, one needs to do : for (weight, l) in quadrature_rule_bary(iside, shape(etype), degree) xp = zeros(SVector{td}) for i=1:nvertices xp += l[i]*vertices[i] end # weight, xp is the quadrature couple (weight, node) end

source