FiniteElementQuadratureRules
Documentation for FiniteElementQuadratureRules.
Base.DictFiniteElementQuadratureRules.AffineGeometryFiniteElementQuadratureRules.AffineGeometryFiniteElementQuadratureRules.BarycentricMonomialsFiniteElementQuadratureRules.BarycentricMonomialsFiniteElementQuadratureRules.CompactQuadratureRuleFiniteElementQuadratureRules.CompactQuadratureRuleFiniteElementQuadratureRules.CompactQuadratureRuleFiniteElementQuadratureRules.CompactQuadratureRuleWithWeightsFiniteElementQuadratureRules.CompactQuadratureRuleWithWeightsFiniteElementQuadratureRules.CompactQuadratureRuleWithWeightsFiniteElementQuadratureRules.JacobiPolySetFiniteElementQuadratureRules.JacobiPolySetFiniteElementQuadratureRules.JacobiPolySetFiniteElementQuadratureRules.LagrangePolySetFiniteElementQuadratureRules.LagrangePolySetFiniteElementQuadratureRules.LagrangePolySetFiniteElementQuadratureRules.MonomialPolySetFiniteElementQuadratureRules.MonomialPolySetFiniteElementQuadratureRules.MonomialPolySetFiniteElementQuadratureRules.MultiLinearGeometryFiniteElementQuadratureRules.MultiLinearGeometryFiniteElementQuadratureRules.MultiLinearGeometryFiniteElementQuadratureRules.QuadratureRuleFiniteElementQuadratureRules.QuadratureRuleFiniteElementQuadratureRules.QuadratureRuleFiniteElementQuadratureRules.QuadratureRuleFiniteElementQuadratureRules.QuadratureRuleFiniteElementQuadratureRules.ReferenceElementBase.positionFiniteElementQuadratureRules.barycentricCoordinatesFiniteElementQuadratureRules.checkInsideFiniteElementQuadratureRules.checkStrictlyInsideFiniteElementQuadratureRules.detectSymmetryOrbitsFiniteElementQuadratureRules.expandFiniteElementQuadratureRules.expandFiniteElementQuadratureRules.getCompactWeightsFiniteElementQuadratureRules.getCompactWeightsFiniteElementQuadratureRules.getCompactWeightsFiniteElementQuadratureRules.getOrbitWeightsFiniteElementQuadratureRules.getOrbitWeightsFiniteElementQuadratureRules.getPropertiesFiniteElementQuadratureRules.getPropertiesFiniteElementQuadratureRules.getQualityFiniteElementQuadratureRules.getWeightsFiniteElementQuadratureRules.getWeightsFiniteElementQuadratureRules.getWeightsFiniteElementQuadratureRules.getWeightsFiniteElementQuadratureRules.integrateFiniteElementQuadratureRules.integrateFiniteElementQuadratureRules.optimize
Functions exported from FiniteElementQuadratureRules:
Base.position — Method
position(ref::ReferenceElement, i::Integer, codim::Integer)Compute the center of the ith sub-entity of codimension codim of the reference element ref.
FiniteElementQuadratureRules.barycentricCoordinates — Method
barycentricCoordinates(domain::AbstractDomain, x::AbstractVector)
barycentricCoordinates(ref::ReferenceElement, x::AbstractVector)Transform reference element coordinates in a given domain into barycentric coordinates. This is in particular useful for simplex domains where these two coordinates differ. Also for prisms and pyramids where some coordinates parts are associated to a conical extension.
Parameters:
domain::AbstractDomain: The domain of the original coordinatesx::AbstractVector: The coordinates in the reference domain
Result: Barycentric coordinates associated to the domain.
FiniteElementQuadratureRules.checkInside — Function
checkInside(ref::ReferenceElement, x::AbstractVector, tol::Real)
checkInside(ref::ReferenceElement, x::AbstractVector{T})Check whether a point x lies inside the reference element ref or on its boundary. This check is performed with a given tolerance tol. If no tolerance is given, use the tolerance tol=eps(T).
FiniteElementQuadratureRules.checkStrictlyInside — Function
checkStrictlyInside(ref::ReferenceElement, x::AbstractVector, tol::Real)
checkStrictlyInside(ref::ReferenceElement, x::AbstractVector{T})Check whether a point x lies strictly inside the reference element ref. This check is performed with a given tolerance tol. If no tolerance is given, use the tolerance tol=eps(T).
FiniteElementQuadratureRules.detectSymmetryOrbits — Method
detectSymmetryOrbits(domain::AbstractDomain, points::AbstractVector;
weights=nothing,
candidate_orbits=collect(eachindex(symmetryOrbits(T, domain))),
point_tol=sqrt(eps(T)),
weight_tol=point_tol,
canonicalize=identity)Infer a compact orbit representation from points, which must already be expressed in the coordinate system expected by the symmetry orbits of domain (for example, barycentric coordinates on triangles and tetrahedra).
If weights are provided, only points with matching weights are grouped into the same orbit, and one representative weight per detected orbit is returned unchanged. The result is a named tuple with fields orbits, positions, and, when applicable, weights.
FiniteElementQuadratureRules.expand — Method
expand(cqr::CompactQuadratureRuleWithWeights{Ω,T})Expand a compact rule into a full quadrature rule. This first expands the symmetry orbits to generate a list of quadrature points and combines it with the given quadrature weights.
FiniteElementQuadratureRules.expand — Method
expand(cqr::CompactQuadratureRule{Ω,T})Expand a compact rule into a full quadrature rule. This first expands the symmetry orbits to generate a list of quadrature points and then computes the associated quadrature weights.
FiniteElementQuadratureRules.getCompactWeights — Method
getCompactWeights(qr::QuadratureRule, orbits::AbstractVector)Extract one compact-format weight per symmetry orbit from an expanded quadrature rule. The returned weights can be stored directly in a CompactQuadratureRuleWithWeights.
FiniteElementQuadratureRules.getCompactWeights — Method
getCompactWeights(cqr::CompactQuadratureRule)Compute one compact-format weight per symmetry orbit for the given compact rule.
FiniteElementQuadratureRules.getCompactWeights — Method
getCompactWeights(::Type{T}, domain::AbstractDomain, orbit_weights::AbstractVector)
getCompactWeights(domain::AbstractDomain, orbit_weights::AbstractVector)
getCompactWeights(::Type{T}, domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getCompactWeights(domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getCompactWeights(::Type{T}, ref::ReferenceElement, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getCompactWeights(ref::ReferenceElement, degree::Integer, points::AbstractVector, orbits::AbstractVector)Convert per-orbit reference-element weights into the compact-rule scaling, or compute compact-format per-orbit weights directly from expanded quadrature points grouped by symmetry orbits.
FiniteElementQuadratureRules.getOrbitWeights — Method
getOrbitWeights(domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getOrbitWeights(ref::ReferenceElement, degree::Integer, points::AbstractVector, orbits::AbstractVector)Compute one quadrature weight per symmetry orbit, inferring the weight type from the point type.
FiniteElementQuadratureRules.getOrbitWeights — Method
getOrbitWeights(::Type{T}, domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getOrbitWeights(::Type{T}, ref::ReferenceElement, degree::Integer, points::AbstractVector, orbits::AbstractVector)Compute one quadrature weight per symmetry orbit for the given expanded quadrature points. In contrast to getWeights, the returned vector contains exactly one weight per orbit rather than one weight per quadrature point.
FiniteElementQuadratureRules.getProperties — Method
getProperties(qr::QuadratureRule)Collect properties of a quadrature rule.
FiniteElementQuadratureRules.getProperties — Method
getProperties(domain::AbstractDomain, points::AbstractVector, weights::AbstractVector)Collect properties of a quadrature given as vector of points and weights. Currently, three properties are checked:
- positive weights: property
:positive - points strictly inside the domain: property
:inside - points are inside or on the boundary of the domain: property
:boundary
The function returns a vector of properties as corresponding Symbols.
FiniteElementQuadratureRules.getQuality — Method
getQuality(properties::AbstractVector{Symbol})
getQuality(ref::ReferenceElement, points::AbstractVector, weights::AbstractVector)
getQuality(qr::QuadratureRule)Return the short quality marker associated with a quadrature rule.
For real-valued rules, the following markers are used:
:PI: all weights are positive and all points are strictly inside the domain:PB: all weights are positive and all points are inside or on the boundary:PO: all weights are positive and at least one point is outside the domain:NI: at least one weight is non-positive and all points are strictly inside:NB: at least one weight is non-positive and all points are inside or on the boundary:NO: at least one weight is non-positive and at least one point is outside the domain
The literature also uses :PC, :NC, and :CC for rules with complex-valued coordinates or weights. These markers are not produced by QuadratureRule, which stores real-valued rules only.
FiniteElementQuadratureRules.getWeights — Method
getWeights(domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getWeights(ref::ReferenceElement, domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)Compute the quadrature weight for the given set of quadrature points, using the symmetry orbits. The data type of the weights is identical to the element type of the points.
FiniteElementQuadratureRules.getWeights — Method
getWeights(domain::AbstractDomain, degree::Integer, points::AbstractVector)
getWeights(ref::ReferenceElement, degree::Integer, points::AbstractVector)Compute the quadrature weight for the given set of quadrature points. The data type of the weights is identical to the element type of the points.
FiniteElementQuadratureRules.getWeights — Method
getWeights(::Type{T}, domain::AbstractDomain, degree::Integer, points::AbstractVector, orbits::AbstractVector)
getWeights(::Type{T}, ref::ReferenceElement, degree::Integer, points::AbstractVector, orbits::AbstractVector)Compute the quadrature weight associated to given quadrature points for quadrature in a given reference element ref, such that polynomial up to degree are integrated exactly. This overload takes into account that the points are associated to symmetry orbits and thus points in the same orbit share a quadrature weight. This makes the computation more stable and faster. The parameter T represents the number type to be used for the computation.
FiniteElementQuadratureRules.getWeights — Method
getWeights(::Type{T}, domain::AbstractDomain, degree::Integer, points::AbstractVector)
getWeights(::Type{T}, ref::ReferenceElement, degree::Integer, points::AbstractVector)Compute the quadrature weight associated to given quadrature points for quadrature in a given reference element ref, such that polynomial up to degree are integrated exactly. The parameter T represents the number type to be used for the computation.
FiniteElementQuadratureRules.integrate — Method
integrate(f::Function, qr::QuadratureRule)Compute the integral ∫f(x)dx over a domain Ω=domain(qr) using a quadrature rule.
FiniteElementQuadratureRules.integrate — Method
integrate(m::BarycentricMonomials, Ω::AbstractSimplex)Compute the integral of a barycentric monomial m over a generic simplex of volume 1. For a specific simplex, this integral needs to be scaled with its volume.
FiniteElementQuadratureRules.optimize — Method
optimize(qr::CompactQuadratureRule)Optimize the position of the quadrature points given as symmetric orbits in the compact rule, by minimizing the quadrature residual on a set of polynomial basis functions. We use basis functions and their integrals given by a JacobiPolySet.
Allowed solvers:
:lm: Levenberg-Marquardt:tr: TrustRegion:gn: Gauss-Newton:str: SimpleTrustRegion
Additional keyword arguments are forwarded to the selected solver constructor.
Types in module FiniteElementQuadratureRules:
Base.Dict — Method
Dict(qr::QuadratureRule, reference::String="unknown", precision::Int=32)Convert the given QuadratureRule into a Dict for exporting into a YAML file. The optional parameter reference refers to a bibtex key used to reference a publication where the quadrature rule is extracted from.
FiniteElementQuadratureRules.AffineGeometry — Type
AffineGeometry{Ω,mydim,cdom,T} <: AbstractGeometryAn AffineGeometry is a parametrization using an affine map of the form x → A*x + b. As an AbstractGeometry it is itself a mapping.
Example:
geo = AffineGeometry(ReferenceElement(Triangle()), [x0,x1,x2])
@assert domaintype(geo) == Triangle()
x = geo(λ) # given λ in coordinates of the reference elementFiniteElementQuadratureRules.AffineGeometry — Method
(geo::AffineGeometry)(λ::AbstractVector)Evaluate the affine geometry mapping geo in the reference element coordinate λ.
FiniteElementQuadratureRules.BarycentricMonomials — Type
BarycentricMonomials <: FunctionA representation of a multi-dimensional monomial basis function of the form λ[1]^α[1] * λ[2]^α[2] * ... with λ the barycentric coordinate vector and α the exponents associated to the components of the coordinates.
FiniteElementQuadratureRules.BarycentricMonomials — Method
(m::BarycentricMonomials)(λ::AbstractVector)Evaluate the monomial basis function m in a given barycentric coordinate λ, i.e., compute λ -> λ[1]^α[1] * λ[2]^α[2] * ...
FiniteElementQuadratureRules.CompactQuadratureRule — Type
CompactQuadratureRule{T,Ω}A compact (quadrature) rule is defined in terms of a vector of symmetry orbits and their corresponding arguments, which allows to generate a list of quadrature points in the reference element and also their corresponding quadrature weights.
The rule is defined on a domain, e.g., a triangle or quad, has an expected degree, which means the maximal polynomial degree up to which the quadrature rule is exact.
FiniteElementQuadratureRules.CompactQuadratureRule — Method
CompactQuadratureRule(data::Dict)Construct a CompactQuadratureRule from parsed YAML/Dict data using Float64.
FiniteElementQuadratureRules.CompactQuadratureRule — Method
CompactQuadratureRule(::Type{T}, data::Dict)Construct a CompactQuadratureRuleWithWeights from a YAML/Dict of strings and string arrays. This is a convenience constructor typically used when reading a quadrature rule from a YAML file. All information is encoded in the Dict, in particular the fields
dimandregion: characterizing the domaindegree: the quadrature degreeorbits: the symmetry orbitspositions: representing the arguments to the symmetry orbits
FiniteElementQuadratureRules.CompactQuadratureRuleWithWeights — Type
CompactQuadratureRuleWithWeights{T,Ω}A compact (quadrature) rule with weights is defined in terms of a vector of symmetry orbits, their corresponding arguments, and a list of associated weights. This allows to generate a list of quadrature points in the reference element.
The rule is defined on a domain, e.g., a triangle or quad, has an expected degree, which means the maximal polynomial degree up to which the quadrature rule is exact.
FiniteElementQuadratureRules.CompactQuadratureRuleWithWeights — Method
CompactQuadratureRuleWithWeights(data::Dict)Construct a CompactQuadratureRuleWithWeights from parsed YAML/Dict data using Float64.
FiniteElementQuadratureRules.CompactQuadratureRuleWithWeights — Method
CompactQuadratureRuleWithWeights(::Type{T}, data::Dict)
Construct a CompactQuadratureRuleWithWeights from a YAML/Dict of strings and string arrays. This is a convenience constructor typically used when reading a quadrature rule from a YAML file. All information is encoded in the Dict, in particular the fields
dimandregion: characterizing the domaindegree: the quadrature degreeorbits: the symmetry orbitspositions: representing the arguments to the symmetry orbitsweights: for the quadrature weights
FiniteElementQuadratureRules.JacobiPolySet — Type
JacobiPolySet{Ω<:AbstractDomain,R<:Real}A JacobiPolySet represents a set of polynomials on a domain Ω in terms of a set of basis polynomials, hereby given as Jacobi polynomials. As an AbstractPolySet it provides also the integral values over the domain of all basis functions.
FiniteElementQuadratureRules.JacobiPolySet — Method
JacobiPolySet(domain::AbstractDomain, degree::Integer)Construct a JacobiPolySet on the given domain of given polynomial degree, with Float64 as data type used for the integral values.
FiniteElementQuadratureRules.JacobiPolySet — Method
JacobiPolySet(::Type, domain::AbstractDomain, degree::Integer)Construct a JacobiPolySet on the given domain of given polynomial degree, with T the data type used for the integral values.
FiniteElementQuadratureRules.LagrangePolySet — Type
LagrangePolySet{Ω<:AbstractDomain,R<:Real}A LagrangePolySet represents a set of polynomials on a domain Ω in terms of a set of basis polynomials, hereby given as Lagrange polynomials. As an AbstractPolySet it provides also the integral values over the domain of all basis functions.
FiniteElementQuadratureRules.LagrangePolySet — Method
LagrangePolySet(domain::AbstractDomain, degree::Integer)Construct a LagrangePolySet on the given domain of given polynomial degree, with Float64 as data type used for the integral values.
FiniteElementQuadratureRules.LagrangePolySet — Method
LagrangePolySet(::Type{T}, domain::AbstractDomain, degree::Integer)Construct a LagrangePolySet on a given domain of given polynomial degree, with T the data type used for the integral values.
FiniteElementQuadratureRules.MonomialPolySet — Type
MonomialPolySet{Ω<:AbstractDomain,R<:Real}A MonomialPolySet represents a set of polynomials on a domain Ω in terms of a set of basis polynomials, hereby given as a monomial basis. As an AbstractPolySet it provides also the integral values over the domain of all basis functions.
FiniteElementQuadratureRules.MonomialPolySet — Method
MonomialPolySet(domain::AbstractDomain, degree::Integer)Construct a MonomialPolySet on the given domain of given polynomial degree, with Float64 as data type used for the integral values.
FiniteElementQuadratureRules.MonomialPolySet — Method
MonomialPolySet(::Type{T}, domain::AbstractSimplex, degree::Integer)Construct a MonomialPolySet on a simplex domain of given polynomial degree, with T the data type used for the integral values.
FiniteElementQuadratureRules.MultiLinearGeometry — Type
MultiLinearGeometry{T,cdim,Ω} <: AbstractGeometryA MultiLinearGeometry is a geometry mapping defined in terms of linear Lagrange basis function associated to the domain Ω. It can be constructed from a reference element and a vector of corner coordinates.
FiniteElementQuadratureRules.MultiLinearGeometry — Method
(geo::MultiLinearGeometry)(λ::AbstractVector)Evaluate the multilinear geometry mapping geo in the reference element coordinates λ, by linear combination of the associated Lagrange basis functions and corner vertices.
FiniteElementQuadratureRules.MultiLinearGeometry — Method
MultiLinearGeometry(ref::ReferenceElement, coordVector::AbstractVector)Construct a MultiLinearGeometry as an affine geometry if the domain is a simplex.
FiniteElementQuadratureRules.QuadratureRule — Type
QuadratureRule{Ω,T,Point}A QuadratureRule on a given domain Ω is a collection of points {xᵢ} in the associated reference element and weights {wᵢ}, such that for a given polynomial p the quadrature formula ∑ᵢ p(xᵢ) wᵢ = ∫p(x) dΩ is exact, up to polynomials of a certain degree.
FiniteElementQuadratureRules.QuadratureRule — Method
QuadratureRule(data::Dict)Construct a QuadratureRule from parsed YAML/Dict data using Float64.
FiniteElementQuadratureRules.QuadratureRule — Method
QuadratureRule(domain::AbstractDomain, degree::Integer, points::Vector)Construct a new QuadratureRule and compute the weights and properties of the rule automatically.
FiniteElementQuadratureRules.QuadratureRule — Method
QuadratureRule(domain::AbstractDomain, degree::Integer, points::Vector, weights::Vector)Construct a new QuadratureRule and compute the properties of the rule automatically.
FiniteElementQuadratureRules.QuadratureRule — Method
QuadratureRule(::Type{T}, data::Dict)Construct a QuadratureRule from a YAML/Dict of strings and string arrays. This is a convenience constructor typically used when reading a quadrature rule from a YAML file. All information is encoded in the Dict, in particular the fields
dimandregion: characterizing the domaindegree: the quadrature degreecoordinates: representing the quadrature pointsweights: for the quadrature weightsproperties: characterizing properties of the quadrature rule
FiniteElementQuadratureRules.ReferenceElement — Type
ReferenceElement{Ω,T,P}A reference representation of a given domain Ω in terms of corner coordinates and a list of facets connecting the corners.