Integration
Measure
Bcube.Measure — Type
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.
Bcube.Measure — Method
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 overdegree: the degree of the quadrature rule (Legendrequadrature type by default)
Examples
julia> mesh = line_mesh(10)
julia> Ω = CellDomain(mesh)
julia> dΩ = Measure(Ω, 2)Bcube.get_bcube_backend — Method
get_bcube_backend(m::AbstractMeasure)Retrieve the Bcube backend associated with the domain of the given measure.
Integration methods
Bcube.MultiIntegration — Type
MultiIntegration{N, I <: Tuple{Vararg{Integration, N}}}Container for multiple Integration objects. Stores a tuple integrations of length N, where each element is an Integration. Guarantees that all integrations share the same backend.
Bcube.apply_quadrature — Method
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.
Bcube.apply_quadrature_v2 — Method
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)
Bcube.get_bcube_backend — Method
get_bcube_backend(integration::Integration)Return the Bcube backend associated with the given Integration by delegating to its measure.
Bcube.get_bcube_backend — Method
get_bcube_backend(a::MultiIntegration)Return the Bcube backend associated with the given MultiIntegration by delegating to the first contained Integration. All integrations in a MultiIntegration are guaranteed to share the same backend.
Bcube.integrate_on_ref_element — Method
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.
Quadrature rules
Bcube.AbstractQuadrature — Type
AbstractQuadrature{T,D}Abstract type representing quadrature of type T and degree D
Bcube.AbstractQuadratureNode — Type
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)Bcube.AbstractQuadratureRule — Type
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)]
Bcube.Quadrature — Type
Quadrature{T,D}Quadrature of type T and degree D
Bcube.QuadratureNode — Type
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.
Bcube.QuadratureRule — Type
QuadratureRule{S,Q}Abstract type representing a quadrature rule for a shape S and quadrature Q
Bcube.QuadratureRule — Method
QuadratureRule(shape::AbstractShape, q::AbstractQuadrature)Return the quadrature rule corresponding to the given Shape according to the selected quadrature 'q'.
Base.length — Method
length(qr::AbstractQuadratureRule)Returns the number of quadrature nodes of qr.
Bcube._gausslegendre1D — Method
_gausslegendre1D(::Val{N}) where N
_gausslobatto1D(::Val{N}) where NReturn N-point Gauss quadrature weights and nodes on the domain [-1:1].
Bcube._get_num_nodes — Method
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$
Bcube._get_num_nodes — Method
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$
Bcube._get_num_nodes_per_dim — Method
get_num_nodes_per_dim(quadrule::AbstractQuadratureRule{S}) where S<:ShapeReturns 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!)
Bcube.evalquadnode — Method
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.
Bcube.get_coords — Method
get_coords(quadnode::AbstractQuadratureNode)Returns the coordinates of quadnode.
Bcube.get_index — Method
get_index(quadnode::AbstractQuadratureNode{S,Q})Returns the index of quadnode in the parent quadrature rule AbstractQuadRules{S,Q}
Bcube.get_nodes — Method
get_nodes(qr::AbstractQuadratureRule)Returns an array containing the coordinates of all quadrature nodes of qr.
Bcube.get_quadnodes — Method
get_quadnodes(qr::QuadratureRule{S,Q,N}) where {S,Q,N}Returns an vector containing each QuadratureNode of qr
Bcube.get_quadrature_points_gausslegendre — Method
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
Bcube.get_quadrature_points_gausslegendre — Method
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]²
Bcube.get_quadrature_points_gausslegendre — Method
Gauss-Legendre quadrature, 7 point rule on triangle.
Bcube.get_weights — Method
get_weights(qr::AbstractQuadratureRule)Returns an array containing the weights of all quadrature nodes of qr.
Bcube.quadrature_points — Method
ref : https://www.math.umd.edu/~tadmor/references/files/Chen%20&%20Shu%20entropy%20stable%20DG%20JCP2017.pdf
Bcube.quadrature_rule — Method
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.
Bcube.quadrature_rule — Method
quadrature_rule(iside::Int, shape::AbstractShape, degree::Val{N}) where NReturn the quadrature rule, computed with barycentric coefficients, corresponding to the given boundary of a shape and the given degree.
Bcube.quadrature_rule_bary — Method
quadrature_rule_bary(::Int, ::AbstractShape, degree::Val{N}) where NReturn 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