Skip to content

Quadratures

VortexPasta.Quadratures Module
julia
Quadratures

Module defining quadrature rules for numerical integration along filaments.

Quadrature rules are determined using the FastGaussQuadrature.jl package.

See the Wikipedia page on Gaussian quadratures for more details.

For convenience, quadrature rules in this module are defined on the [0,1] interval, as opposed to the standard [1,1] interval.

source

Quadrature rules

VortexPasta.Quadratures.AbstractQuadrature Type
julia
AbstractQuadrature

Abstract type defining a quadrature rule.

source

Fixed-size quadratures

These quadrature rules are meant to have a small size (typically less than 10). They have virtually zero creation cost, i.e. doing quad = GaussLegendre(4) is basically free.

VortexPasta.Quadratures.StaticSizeQuadrature Type
julia
StaticSizeQuadrature{N} <: AbstractQuadrature

Abstract type defining an N-point quadrature rule.

These quadrature rules can be created basically for free.

source
VortexPasta.Quadratures.GaussLegendre Type
julia
GaussLegendre(N) <: StaticSizeQuadrature{N}

N-point Gauss–Legendre quadrature.

source
VortexPasta.Quadratures.NoQuadrature Type
julia
NoQuadrature() <: StaticSizeQuadrature{1}

Inexpensive 1-point quadrature rule.

When integrating, evaluates approximated values at the segment midpoint. However, unlike a 1-point Gauss–Legendre quadrature, midpoint values are interpolated using a basic linear interpolation. In other words, segments are assumed to be straight.

source

Variable-size quadratures

These quadrature rules should be constructed just once, as they allocate vectors. These are usually adaptive quadratures.

VortexPasta.Quadratures.PreallocatedQuadrature Type
julia
PreallocatedQuadrature{T <: AbstractFloat} <: AbstractQuadrature

Abstract type defining a preallocated quadrature with element type T.

This is generally used for adaptive quadratures.

source
VortexPasta.Quadratures.AdaptiveTanhSinh Type
julia
AdaptiveTanhSinh([T = Float64]; nlevels = 10, rtol = sqrt(eps(T))) <: PreallocatedQuadrature{T}

Adaptive tanh-sinh quadrature.

Behaves well when there are singularities at (or near) the endpoints.

It can be easily made adaptive because it can be written as a simple trapezoidal rule after a change of variables.

Optional arguments

  • T = Float64: quadrature precision;

  • nlevels = 10: maximum number of adaptivity levels. Must be ≥ 2;

  • rtol = sqrt(eps(T)): relative tolerance.

Computations are stopped either when the maximum adaptivity level is reached, or when the difference between two levels falls below the relative tolerance rtol.

Note that the maximum number of function evaluations is 2^nlevels, so it can make sense to use a small number of levels (or a large tolerance) when function evaluations are expensive.

source

Estimating integrals

VortexPasta.Quadratures.integrate Function
julia
Quadratures.integrate(f::Function, quad::AbstractQuadrature, lims::NTuple{2,T})
Quadratures.integrate(f::Function, quad::AbstractQuadrature, ::Type{T})

Integrate f(x) using the chosen quadrature.

There are two variants:

  • The first one requires passing the limits (a, b), which should be floats of the desired accuracy;

  • the second one assumes the default limits (0, 1), but requires setting the float type T. It may avoid some operations when the limits are (0, 1).

In both cases T must be a subtype of AbstractFloat (e.g. Float64).

Note that T is simply ignored for adaptive quadratures such as AdaptiveTanhSinh, which preallocate the quadrature coefficients with a possibly different element type T′ (which can be chosen when creating the quadrature).

source

Computing quadrature rules

VortexPasta.Quadratures.quadrature Function
julia
quadrature([T = Float64], q::StaticSizeQuadrature{N}) -> (xs, ws)

Return N-point quadrature rule valid in [0,1] interval.

Quadrature nodes xs and weights ws are returned as tuples.

This only works for fixed-size (non-adaptive) quadrature rules such as GaussLegendre.

source

Other functions

Base.length Function
julia
Base.length(q::StaticSizeQuadrature{N}) -> N

Return the number of points associated to the quadrature.

source
julia
Base.length(s::SegmentIterator{<:AbstractFilament}) -> Int

Return the number of segments in a filament.

source