QuadGK.jl
This package provides support for one-dimensional numerical integration in Julia using adaptive Gauss-Kronrod quadrature. The code was originally part of Base Julia.
Functions
QuadGK.quadgk
— Function.quadgk(f, a,b,c...; rtol=sqrt(eps), atol=0, maxevals=10^7, order=7, norm=norm)
Numerically integrate the function f(x)
from a
to b
, and optionally over additional intervals b
to c
and so on. Keyword options include a relative error tolerance rtol
(defaults to sqrt(eps)
in the precision of the endpoints), an absolute error tolerance atol
(defaults to 0), a maximum number of function evaluations maxevals
(defaults to 10^7
), and the order
of the integration rule (defaults to 7).
Returns a pair (I,E)
of the estimated integral I
and an estimated upper bound on the absolute error E
. If maxevals
is not exceeded then E <= max(atol, rtol*norm(I))
will hold. (Note that it is useful to specify a positive atol
in cases where norm(I)
may be zero.)
The endpoints a
et cetera can also be complex (in which case the integral is performed over straight-line segments in the complex plane). If the endpoints are BigFloat
, then the integration will be performed in BigFloat
precision as well.
It is advisable to increase the integration order
in rough proportion to the precision, for smooth integrands.
More generally, the precision is set by the precision of the integration endpoints (promoted to floating-point types).
The integrand f(x)
can return any numeric scalar, vector, or matrix type, or in fact any type supporting +
, -
, multiplication by real values, and a norm
(i.e., any normed vector space). Alternatively, a different norm can be specified by passing a norm
-like function as the norm
keyword argument (which defaults to norm
).
Only one-dimensional integrals are provided by this function. For multi-dimensional integration (cubature), there are many different algorithms (often much better than simple nested 1d integrals) and the optimal choice tends to be very problem-dependent. See the Julia external-package listing for available algorithms for multidimensional integration or other specialized tasks (such as integrals of highly oscillatory or singular functions).
The algorithm is an adaptive Gauss-Kronrod integration technique: the integral in each interval is estimated using a Kronrod rule (2*order+1
points) and the error is estimated using an embedded Gauss rule (order
points). The interval with the largest error is then subdivided into two intervals and the process is repeated until the desired error tolerance is achieved.
These quadrature rules work best for smooth functions within each interval, so if your function has a known discontinuity or other singularity, it is best to subdivide your interval to put the singularity at an endpoint. For example, if f
has a discontinuity at x=0.7
and you want to integrate from 0 to 1, you should use quadgk(f, 0,0.7,1)
to subdivide the interval at the point of discontinuity. The integrand is never evaluated exactly at the endpoints of the intervals, so it is possible to integrate functions that diverge at the endpoints as long as the singularity is integrable (for example, a log(x)
or 1/sqrt(x)
singularity).
For real-valued endpoints, the starting and/or ending points may be infinite. (A coordinate transformation is performed internally to map the infinite interval to a finite one.)
QuadGK.gauss
— Function.gauss([T,] N)
Return a pair (x, w)
of N
quadrature points x[i]
and weights w[i]
to integrate functions on the interval (-1, 1)
, i.e. sum(w .* f.(x))
approximates the integral. Uses the method described in Trefethen & Bau, Numerical Linear Algebra, to find the N
-point Gaussian quadrature in O(N
²) operations.
T
is an optional parameter specifying the floating-point type, defaulting to Float64
. Arbitrary precision (BigFloat
) is also supported.
QuadGK.kronrod
— Function.kronrod([T,] n)
Compute 2n+1
Kronrod points x
and weights w
based on the description in Laurie (1997), appendix A, simplified for a=0
, for integrating on [-1,1]
. Since the rule is symmetric, this only returns the n+1
points with x <= 0
. The function Also computes the embedded n
-point Gauss quadrature weights gw
(again for x <= 0
), corresponding to the points x[2:2:end]
. Returns (x,w,wg)
in O(n
²) operations.
T
is an optional parameter specifying the floating-point type, defaulting to Float64
. Arbitrary precision (BigFloat
) is also supported.
Given these points and weights, the estimated integral I
and error E
can be computed for an integrand f(x)
as follows:
x, w, wg = kronrod(n)
fx⁰ = f(x[end]) # f(0)
x⁻ = x[1:end-1] # the x < 0 Kronrod points
fx = f.(x⁻) .+ f.((-).(x⁻)) # f(x < 0) + f(x > 0)
I = sum(fx .* w[1:end-1]) + fx⁰ * w[end]
if isodd(n)
E = abs(sum(fx[2:2:end] .* wg[1:end-1]) + fx⁰*wg[end] - I)
else
E = abs(sum(fx[2:2:end] .* wg[1:end])- I)
end