PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_quad_1d_gauss_vec (d01ua)
Purpose
nag_quad_1d_gauss_vec (d01ua) computes an estimate of the definite integral of a function of known analytical form, using a Gaussian quadrature formula with a specified number of abscissae. Formulae are provided for a finite interval (Gauss–Legendre), a semi-infinite interval (Gauss–Laguerre, rational Gauss), and an infinite interval (Gauss–Hermite).
Syntax
Description
General
nag_quad_1d_gauss_vec (d01ua) evaluates an estimate of the definite integral of a function
, over a finite or infinite range, by
-point Gaussian quadrature (see
Davis and Rabinowitz (1975),
Fröberg (1970),
Ralston (1965) or
Stroud and Secrest (1966)). The integral is approximated by a summation of the product of a set of weights and a set of function evaluations at a corresponding set of abscissae
. For adjusted weights, the function values correspond to the values of the integrand
, and hence the sum will be
where the
are called the weights, and the
the abscissae. A selection of values of
is available. (See
Arguments.)
Where applicable, normal weights may instead be used, in which case the corresponding weight function
is factored out of the integrand as
and hence the sum will be
where the normal weights
are computed internally.
nag_quad_1d_gauss_vec (d01ua) uses a vectorized
f to evaluate the integrand or normalized integrand at a set of abscissae,
, for
. If adjusted weights are used, the integrand
must be evaluated otherwise the normalized integrand
must be evaluated.
Both Limits Finite
The Gauss–Legendre weights and abscissae are used, and the formula is exact for any function of the form:
The formula is appropriate for functions which can be well approximated by such a polynomial over
. It is inappropriate for functions with algebraic singularities at one or both ends of the interval, such as
on
.
One Limit Infinite
Two quadrature formulae are available for these integrals.
(a) |
The Gauss–Laguerre formula is exact for any function of the form:
This formula is appropriate for functions decaying exponentially at infinity; the argument should be chosen if possible to match the decay rate of the function.
If the adjusted weights are selected, the complete integrand should be provided through f.
If the normal form is selected, the contribution of is accounted for internally, and f should only return , where .
If is supplied, the interval of integration will be . Otherwise if is supplied, the interval of integration will be . |
(b) |
The rational Gauss formula is exact for any function of the form:
This formula is likely to be more accurate for functions having only an inverse power rate of decay for large . Here the choice of a suitable value of may be more difficult; unfortunately a poor choice of can make a large difference to the accuracy of the computed integral.
Only the adjusted form of the rational Gauss formula is available, and as such, the complete integrand must be supplied in f.
If , the interval of integration will be . Otherwise if , the interval of integration will be . |
Both Limits Infinite
The Gauss–Hermite weights and abscissae are used, and the formula is exact for any function of the form:
where
. Again, for general functions not of this exact form, the argument
should be chosen to match if possible the decay rate at
.
If the adjusted weights are selected, the complete integrand
should be provided through
f.
If the normal form is selected, the contribution of
is accounted for internally, and
f should only return
, where
.
References
Davis P J and Rabinowitz P (1975) Methods of Numerical Integration Academic Press
Fröberg C E (1970) Introduction to Numerical Analysis Addison–Wesley
Ralston A (1965) A First Course in Numerical Analysis pp. 87–90 McGraw–Hill
Stroud A H and Secrest D (1966) Gaussian Quadrature Formulas Prentice–Hall
Parameters
Compulsory Input Parameters
- 1:
– int64int32nag_int scalar
-
Indicates the quadrature formula.
- Gauss–Legendre quadrature on a finite interval, using normal weights.
- Gauss–Laguerre quadrature on a semi-infinite interval, using normal weights.
- Gauss–Laguerre quadrature on a semi-infinite interval, using adjusted weights.
- Gauss–Hermite quadrature on an infinite interval, using normal weights.
- Gauss–Hermite quadrature on an infinite interval, using adjusted weights.
- Rational Gauss quadrature on a semi-infinite interval, using adjusted weights.
Constraint:
, , , , or .
- 2:
– double scalar
- 3:
– double scalar
-
The quantities
and
as described in the appropriate subsection of
Description.
Constraints:
- Rational Gauss: ;
- Gauss–Laguerre: ;
- Gauss–Hermite: .
- 4:
– int64int32nag_int scalar
-
, the number of abscissae to be used.
Constraint:
,
,
,
,
,
,
,
,
,
,
,
,
,
,
or
.
If the soft fail option is used, the answer is evaluated for the largest valid value of
n less than the requested value.
- 5:
– function handle or string containing name of m-file
-
f must return the value of the integrand
, or the normalized integrand
, at a specified point.
[fv, iflag, user] = f(x, nx, iflag, user)
Input Parameters
- 1:
– double array
-
The abscissae,
, for at which function values are required.
- 2:
– int64int32nag_int scalar
-
, the number of abscissae.
- 3:
– int64int32nag_int scalar
-
.
- 4:
– Any MATLAB object
f is called from
nag_quad_1d_gauss_vec (d01ua) with the object supplied to
nag_quad_1d_gauss_vec (d01ua).
Output Parameters
- 1:
– double array
-
If adjusted weights are used, the values of the integrand
.
, for
.
Otherwise the values of the normalized integrand .
, for .
- 2:
– int64int32nag_int scalar
-
Set if you wish to force an immediate exit from nag_quad_1d_gauss_vec (d01ua) with .
- 3:
– Any MATLAB object
Some points to bear in mind when coding
f are mentioned in
Accuracy.
Optional Input Parameters
- 1:
– Any MATLAB object
user is not used by
nag_quad_1d_gauss_vec (d01ua), but is passed to
f. Note that for large objects it may be more efficient to use a global variable which is accessible from the m-files than to use
user.
Output Parameters
- 1:
– double scalar
-
The estimate of the definite integral.
- 2:
– Any MATLAB object
- 3:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
Note: nag_quad_1d_gauss_vec (d01ua) may return useful information for one or more of the following detected errors or warnings.
Errors or warnings detected by the function:
Cases prefixed with W are classified as warnings and
do not generate an error of type NAG:error_n. See nag_issue_warnings.
- W
-
The -point rule is not among those stored.
- W
-
Underflow occurred in calculation of normal weights.
- W
-
No nonzero weights were generated for the provided parameters.
-
-
Constraint: , , , , or .
-
The value of
a and/or
b is invalid for the chosen
key. Either:
- Constraint: .
- Constraint: .
- Constraint: .
-
-
Constraint: .
- W
-
Exit requested from
f with
.
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
The accuracy depends on the behaviour of the integrand, and on the number of abscissae used. No tests are carried out in nag_quad_1d_gauss_vec (d01ua) to estimate the accuracy of the result. If such an estimate is required, the function may be called more than once, with a different number of abscissae each time, and the answers compared. It is to be expected that for sufficiently smooth functions a larger number of abscissae will give improved accuracy.
Alternatively, the range of integration may be subdivided, the integral estimated separately for each sub-interval, and the sum of these estimates compared with the estimate over the whole range.
The coding of
f may also have a bearing on the accuracy. For example, if a high-order Gauss–Laguerre formula is used, and the integrand is of the form
it is possible that the exponential term may underflow for some large abscissae. Depending on the machine, this may produce an error, or simply be assumed to be zero. In any case, it would be better to evaluate the expression with
Another situation requiring care is exemplified by
The integrand here assumes very large values; for example, when , the peak value exceeds . Now, if the machine holds floating-point numbers to an accuracy of significant decimal digits, we could not expect such terms to cancel in the summation leaving an answer of much less than (the weights being of order unity); that is, instead of zero we obtain a rather large answer through rounding error. Such situations are characterised by great variability in the answers returned by formulae with different values of .
In general, you should be aware of the order of magnitude of the integrand, and should judge the answer in that light.
Further Comments
The time taken by nag_quad_1d_gauss_vec (d01ua) depends on the complexity of the expression for the integrand and on the number of abscissae required.
Example
This example evaluates the integrals
by Gauss–Legendre quadrature;
by rational Gauss quadrature with
;
by Gauss–Laguerre quadrature with
; and
by Gauss–Hermite quadrature with
and
.
The formulae with and are used in each case. Both adjusted and normal weights are used for Gauss–Laguerre and Gauss–Hermite quadrature.
Open in the MATLAB editor:
d01ua_example
function d01ua_example
fprintf('d01ua example results\n\n');
for funid=1:6
switch funid
case 1
fprintf('\nGauss-Legendre example\n');
a = 0.0;
b = 1.0;
key = int64(0);
case 2
fprintf('\nRational Gauss example\n');
a = 2.0;
b = 0.0;
key = int64(-5);
case 3
fprintf('\nGauss-Laguerre example (adjusted weights)\n');
a = 2.0;
b = 1.0;
key = int64(-3);
case 4
fprintf('\nGauss-Laguerre example (normal weights)\n');
a = 2.0;
b = 1.0;
key = int64(3);
case 5
fprintf('\nGauss-Hermite example (adjusted weights)\n');
a = -1.0;
b = 3.0;
key = int64(-4);
case 6
fprintf('\nGauss-Hermite example (normal weights)\n');
a = -1.0;
b = 3.0;
key = int64(4);
end
for i=1:6
n = int64(2^i);
[dinest, user, ifail] = d01ua(key, a, b, n, @f, 'user', funid);
if ifail == 0 || ifail == 1
fprintf('%2d Points Answer = %10.5f\n', n, dinest);
end
end
end
function [fv, iflag, user] = f(x, nx, iflag, user)
switch user
case 1
fv = 4./(1+x.*x);
case 2
fv = 1./(x.*x.*log(x));
case 3
fv = exp(-x)./x;
case 4
fv = 1./x;
case 5
fv = exp(-3.*x.*x-4.*x-1);
case 6
fv = exp(2*x+2);
otherwise
fv = zeros(nx, 1);
iflag = -1;
end
d01ua example results
Gauss-Legendre example
2 Points Answer = 3.14754
4 Points Answer = 3.14161
8 Points Answer = 3.14159
16 Points Answer = 3.14159
32 Points Answer = 3.14159
64 Points Answer = 3.14159
Rational Gauss example
2 Points Answer = 0.37989
4 Points Answer = 0.37910
8 Points Answer = 0.37876
16 Points Answer = 0.37869
32 Points Answer = 0.37867
64 Points Answer = 0.37867
Gauss-Laguerre example (adjusted weights)
2 Points Answer = 0.04833
4 Points Answer = 0.04887
8 Points Answer = 0.04890
16 Points Answer = 0.04890
32 Points Answer = 0.04890
64 Points Answer = 0.04890
Gauss-Laguerre example (normal weights)
2 Points Answer = 0.04833
4 Points Answer = 0.04887
8 Points Answer = 0.04890
16 Points Answer = 0.04890
32 Points Answer = 0.04890
64 Points Answer = 0.04890
Gauss-Hermite example (adjusted weights)
2 Points Answer = 1.38381
4 Points Answer = 1.42803
8 Points Answer = 1.42817
16 Points Answer = 1.42817
32 Points Answer = 1.42817
64 Points Answer = 1.42817
Gauss-Hermite example (normal weights)
2 Points Answer = 1.38381
4 Points Answer = 1.42803
8 Points Answer = 1.42817
16 Points Answer = 1.42817
32 Points Answer = 1.42817
64 Points Answer = 1.42817
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015