This chapter describes routines for performing numerical integration
(quadrature) of a function in one dimension. There are routines for
adaptive and non-adaptive integration of general functions, with
specialised routines for specific cases. These include integration over
infinite and semi-infinite ranges, singular integrals, including
logarithmic singularities, computation of Cauchy principal values and
oscillatory integrals. The library reimplements the algorithms used in
QUADPACK, a numerical integration package written by Piessens,
Doncker-Kapenga, Uberhuber and Kahaner. Fortran code for QUADPACK is
available on Netlib.
The functions described in this chapter are declared in the header file
`gsl_integration.h'.
Each algorithm computes an approximation to a definite integral of the
form,
I = \int_a^b f(x) w(x) dx
where w(x) is a weight function (for general integrands w(x)=1).
The user provides absolute and relative error bounds
(epsabs, epsrel) which specify the following accuracy requirement,
|RESULT - I| <= max(epsabs, epsrel |I|)
where
RESULT is the numerical approximation obtained by the
algorithm. The algorithms attempt to estimate the absolute error
ABSERR = |RESULT - I| in such a way that the following inequality
holds,
|RESULT - I| <= ABSERR <= max(epsabs, epsrel |I|)
The routines will fail to converge if the error bounds are too
stringent, but always return the best approximation obtained up to that
stage.
The algorithms in QUADPACK use a naming convention based on the
following letters,
Q - quadrature routine
N - non-adaptive integrator
A - adaptive integrator
G - general integrand (user-defined)
W - weight function with integrand
S - singularities can be more readily integrated
P - points of special difficulty can be supplied
I - infinite range of integration
O - oscillatory weight function, cos or sin
F - Fourier integral
C - Cauchy principal value
The algorithms are built on pairs of quadrature rules, a higher order
rule and a lower order rule. The higher order rule is used to compute
the best approximation to an integral over a small range. The
difference between the results of the higher order rule and the lower
order rule gives an estimate of the error in the approximation.
The algorithms for general functions (without a weight function) are
based on Gauss-Kronrod rules. A Gauss-Kronrod rule begins with a
classical Gaussian quadrature rule of order m. This is extended
with additional points between each of the abscissae to give a higher
order Kronrod rule of order 2m+1. The Kronrod rule is efficient
because it reuses existing function evaluations from the Gaussian rule.
The higher order Kronrod rule is used as the best approximation to the
integral, and the difference between the two rules is used as an
estimate of the error in the approximation.
For integrands with weight functions the algorithms use Clenshaw-Curtis
quadrature rules. A Clenshaw-Curtis rule begins with an n-th
order Chebyschev polynomial approximation to the integrand. This
polynomial can be integrated exactly to give an approximation to the
integral of the original function. The Chebyschev expansion can be
extended to higher orders to improve the approximation. The presence of
singularities (or other behavior) in the integrand can cause slow
convergence in the Chebyschev approximation. The modified
Clenshaw-Curtis rules used in QUADPACK separate out several common
weight functions which cause slow convergence. These weight functions
are integrated analytically against the Chebyschev polynomials to
precompute modified Chebyschev moments. Combining the moments with
the Chebyschev approximation to the function gives the desired
integral. The use of analytic integration for the singular part of the
function allows exact cancellations and substantially improves the
overall convergence behavior of the integration.
The QNG algorithm is non-adaptive procedure which uses fixed
Gauss-Kronrod abscissae to sample the integrand at a maximum of 87
points. It is provided for fast integration of smooth functions.
This function applies the Gauss-Kronrod 10-point, 21-point, 43-point and
87-point integration rules in succession until an estimate of the
integral of f over (a,b) is achieved within the desired
absolute and relative error limits, epsabs and epsrel. The
function returns the final approximation, result, an estimate of
the absolute error, abserr and the number of function evaluations
used, neval. The Gauss-Kronrod rules are designed in such a way
that each rule uses all the results of its predecessors, in order to
minimize the total number of function evaluations.
The QAG algorithm is a simple adaptive integration procedure. The
integration region is divided into subintervals, and on each iteration
the subinterval with the largest estimated error is bisected. This
reduces the overall error rapidly, as the subintervals become
concentrated around local difficulties in the integrand. These
subintervals are managed by a gsl_integration_workspace struct,
which handles the memory for the subinterval ranges, results and error
estimates.
This function frees the memory associated with the workspace w.
Function: int gsl_integration_qag(const gsl_function *f, double a, double b, double epsabs, double epsrel, size_t limit, int key, gsl_integration_workspace * workspace, double * result, double * abserr)
This function applies an integration rule adaptively until an estimate
of the integral of f over (a,b) is achieved within the
desired absolute and relative error limits, epsabs and
epsrel. The function returns the final approximation,
result, and an estimate of the absolute error, abserr. The
integration rule is determined by the value of key, which should
be chosen from the following symbolic names,
corresponding to the 15, 21, 31, 41, 51 and 61 point Gauss-Kronrod
rules. The higher-order rules give better accuracy for smooth functions,
while lower-order rules save time when the function contains local
difficulties, such as discontinuities.
On each iteration the adaptive integration strategy bisects the interval
with the largest error estimate. The subintervals and their results are
stored in the memory provided by workspace. The maximum number of
subintervals is given by limit, which may not exceed the allocated
size of the workspace.
The presence of an integrable singularity in the integration region
causes an adaptive routine to concentrate new subintervals around the
singularity. As the subintervals decrease in size the successive
approximations to the integral converge in a limiting fashion. This
approach to the limit can be accelerated using an extrapolation
procedure. The QAGS algorithm combines adaptive bisection with the Wynn
epsilon-algorithm to speed up the integration of many types of
integrable singularities.
This function applies the Gauss-Kronrod 21-point integration rule
adaptively until an estimate of the integral of f over
(a,b) is achieved within the desired absolute and relative error
limits, epsabs and epsrel. The results are extrapolated
using the epsilon-algorithm, which accelerates the convergence of the
integral in the presence of discontinuities and integrable
singularities. The function returns the final approximation from the
extrapolation, result, and an estimate of the absolute error,
abserr. The subintervals and their results are stored in the
memory provided by workspace. The maximum number of subintervals
is given by limit, which may not exceed the allocated size of the
workspace.
This function applies the adaptive integration algorithm QAGS taking
account of the user-supplied locations of singular points. The array
pts of length npts should contain the endpoints of the
integration ranges defined by the integration region and locations of
the singularities. For example, to integrate over the region
(a,b) with break-points at x_1, x_2, x_3 (where
a < x_1 < x_2 < x_3 < b) the following pts array should be used
pts[0] = a
pts[1] = x_1
pts[2] = x_2
pts[3] = x_3
pts[4] = b
with npts = 5.
If you know the locations of the singular points in the integration
region then this routine will be faster than QAGS.
This function computes the integral of the function f over the
infinite interval (-\infty,+\infty). The integral is mapped onto the
interval (0,1] using the transformation x = (1-t)/t,
It is then integrated using the QAGS algorithm. The normal 21-point
Gauss-Kronrod rule of QAGS is replaced by a 15-point rule, because the
transformation can generate an integrable singularity at the origin. In
this case a lower-order rule is more efficient.
Function: int gsl_integration_qagiu(gsl_function * f, double a, double epsabs, double epsrel, size_t limit, gsl_integration_workspace * workspace, double *result, double *abserr)
This function computes the integral of the function f over the
semi-infinite interval (a,+\infty). The integral is mapped onto the
interval (0,1] using the transformation x = a + (1-t)/t,
This function computes the integral of the function f over the
semi-infinite interval (-\infty,b). The integral is mapped onto the
region (0,1] using the transformation x = b - (1-t)/t,
This function computes the Cauchy principal value of the integral of
f over (a,b), with a singularity at c,
I = \int_a^b dx f(x) / (x - c)
The adaptive bisection algorithm of QAG is used, with modifications to
ensure that subdivisions do not occur at the singular point x = c.
When a subinterval contains the point x = c or is close to
it then a special 25-point modified Clenshaw-Curtis rule is used to control
the singularity. Further away from the
singularity the algorithm uses an ordinary 15-point Gauss-Kronrod
integration rule.
The QAWS algorithm is designed for integrands with algebraic-logarithmic
singularities at the end-points of an integration region. In order to
work efficiently the algorithm requires a precomputed table of
Chebyschev moments.
Function: gsl_integration_qaws_table * gsl_integration_qaws_table_alloc(double alpha, double beta, int mu, int nu)
This function allocates space for a gsl_integration_qaws_table
struct and associated workspace describing a singular weight function
W(x) with the parameters (\alpha, \beta, \mu, \nu),
This function computes the integral of the function f(x) over the
interval (a,b) with the singular weight function
(x-a)^\alpha (b-x)^\beta \log^\mu (x-a) \log^\nu (b-x). The parameters
of the weight function (\alpha, \beta, \mu, \nu) are taken from the
table t. The integral is,
The adaptive bisection algorithm of QAG is used. When a subinterval
contains one of the endpoints then a special 25-point modified
Clenshaw-Curtis rule is used to control the singularities. For
subintervals which do not include the endpoints an ordinary 15-point
Gauss-Kronrod integration rule is used.
The QAWO algorithm is designed for integrands with an oscillatory
factor, \sin(\omega x) or \cos(\omega x). In order to
work efficiently the algorithm requires a table of Chebyschev moments
which must be pre-computed with calls to the functions below.
This function allocates space for a gsl_integration_qawo_table
struct and its associated workspace describing a sine or cosine weight
function W(x) with the parameters (\omega, L),
W(x) = sin(omega x)
W(x) = cos(omega x)
The parameter L must be the length of the interval over which the
function will be integrated L = b - a. The choice of sine or
cosine is made with the parameter sine which should be chosen from
one of the two following symbolic values:
GSL_INTEG_COSINE
GSL_INTEG_SINE
The gsl_integration_qawo_table is a table of the trigonometric
coefficients required in the integration process. The parameter n
determines the number of levels of coefficients that are computed. Each
level corresponds to one bisection of the interval L, so that
n levels are sufficient for subintervals down to the length
L/2^n. The integration routine gsl_integration_qawo
returns the error GSL_ETABLE if the number of levels is
insufficient for the requested accuracy.
This function uses an adaptive algorithm to compute the integral of
f over (a,b) with the weight function
\sin(\omega x) or \cos(\omega x) defined
by the table wf.
I = \int_a^b dx f(x) sin(omega x)
I = \int_a^b dx f(x) cos(omega x)
The results are extrapolated using the epsilon-algorithm to accelerate
the convergence of the integral. The function returns the final
approximation from the extrapolation, result, and an estimate of
the absolute error, abserr. The subintervals and their results are
stored in the memory provided by workspace. The maximum number of
subintervals is given by limit, which may not exceed the allocated
size of the workspace.
Those subintervals with "large" widths d, d\omega > 4 are
computed using a 25-point Clenshaw-Curtis integration rule, which handles the
oscillatory behavior. Subintervals with a "small" width
d\omega < 4 are computed using a 15-point Gauss-Kronrod integration.
This function attempts to compute a Fourier integral of the function
f over the semi-infinite interval [a,+\infty).
I = \int_a^{+\infty} dx f(x) sin(omega x)
I = \int_a^{+\infty} dx f(x) cos(omega x)
The parameter \omega is taken from the table wf (the length
L can take any value, since it is overridden by this function to a
value appropriate for the fourier integration). The integral is computed
using the QAWO algorithm over each of the subintervals,
C_1 = [a, a + c]
C_2 = [a + c, a + 2 c]
... = ...
C_k = [a + (k-1) c, a + k c]
where
c = (2 floor(|\omega|) + 1) \pi/|\omega|. The width c is
chosen to cover an odd number of periods so that the contributions from
the intervals alternate in sign and are monotonically decreasing when
f is positive and monotonically decreasing. The sum of this
sequence of contributions is accelerated using the epsilon-algorithm.
This function works to an overall absolute tolerance of
abserr. The following strategy is used: on each interval
C_k the algorithm tries to achieve the tolerance
TOL_k = u_k abserr
where
u_k = (1 - p)p^{k-1} and p = 9/10.
The sum of the geometric series of contributions from each interval
gives an overall tolerance of abserr.
If the integration of a subinterval leads to difficulties then the
accuracy requirement for subsequent intervals is relaxed,
TOL_k = u_k max(abserr, max_{i<k}{E_i})
where E_k is the estimated error on the interval C_k.
The subintervals and their results are stored in the memory provided by
workspace. The maximum number of subintervals is given by
limit, which may not exceed the allocated size of the workspace.
The integration over each subinterval uses the memory provided by
cycle_workspace as workspace for the QAWO algorithm.
The integrator QAGS will handle a large class of definite
integrals. For example, consider the following integral, which has a
algebraic-logarithmic singularity at the origin,
\int_0^1 x^{-1/2} log(x) dx = -4
The program below computes this integral to a relative accuracy bound of
1e-7.
The results below show that the desired accuracy is achieved after 8
subdivisions.
bash$ ./a.out
result = -3.999999999999973799
exact result = -4.000000000000000000
estimated error = 0.000000000000246025
actual error = 0.000000000000026201
intervals = 8
In fact, the extrapolation procedure used by QAGS produces an
accuracy of almost twice as many digits. The error estimate returned by
the extrapolation procedure is larger than the actual error, giving a
margin of safety of one order of magnitude.
The following book is the definitive reference for QUADPACK, and was
written by the original authors. It provides descriptions of the
algorithms, program listings, test programs and examples. It also
includes useful advice on numerical integration and many references to
the numerical integration literature used in developing QUADPACK.
R. Piessens, E. de Doncker-Kapenga, C.W. Uberhuber, D.K. Kahaner.
QUADPACK A subroutine package for automatic integration
Springer Verlag, 1983.