NAG Library Function Document
nag_quad_1d_gauss_vec (d01uac)
1
Purpose
nag_quad_1d_gauss_vec (d01uac) 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).
2
Specification
#include <nag.h> |
#include <nagd01.h> |
void |
nag_quad_1d_gauss_vec (Nag_QuadType quad_type,
double a,
double b,
Integer n,
void |
(*f)(const double x[],
Integer nx,
double fv[],
Integer *iflag,
Nag_Comm *comm),
|
|
double *dinest,
Nag_Comm *comm,
NagError *fail) |
|
3
Description
3.1
General
nag_quad_1d_gauss_vec (d01uac) 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
Section 5.)
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 (d01uac) 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.
3.2
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
.
3.3
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 . |
3.4
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
.
4
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
5
Arguments
- 1:
– Nag_QuadTypeInput
-
On entry: 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:
– doubleInput
- 3:
– doubleInput
-
On entry: the quantities
and
as described in the appropriate subsection of
Section 3.
Constraints:
- Rational Gauss: ;
- Gauss–Laguerre: ;
- Gauss–Hermite: .
- 4:
– IntegerInput
-
On entry: , 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, supplied by the userExternal Function
-
f must return the value of the integrand
, or the normalized integrand
, at a specified point.
The specification of
f is:
void |
f (const double x[],
Integer nx,
double fv[],
Integer *iflag,
Nag_Comm *comm)
|
|
- 1:
– const doubleInput
-
On entry: the abscissae,
, for at which function values are required.
- 2:
– IntegerInput
-
On entry: , the number of abscissae.
- 3:
– doubleOutput
-
On exit: if adjusted weights are used, the values of the integrand
.
, for
.
Otherwise the values of the normalized integrand .
, for .
- 4:
– Integer *Input/Output
-
On entry: .
On exit: set
if you wish to force an immediate exit from
nag_quad_1d_gauss_vec (d01uac) with
NE_USER_STOP.
- 5:
– Nag_Comm *
Pointer to structure of type Nag_Comm; the following members are relevant to
f.
- user – double *
- iuser – Integer *
- p – Pointer
The type Pointer will be
void *. Before calling
nag_quad_1d_gauss_vec (d01uac) you may allocate memory and initialize these pointers with various quantities for use by
f when called from
nag_quad_1d_gauss_vec (d01uac) (see
Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
Note: f should not return floating-point NaN (Not a Number) or infinity values, since these are not handled by
nag_quad_1d_gauss_vec (d01uac). If your code inadvertently
does return any NaNs or infinities,
nag_quad_1d_gauss_vec (d01uac) is likely to produce unexpected results.
Some points to bear in mind when coding
f are mentioned in
Section 7.
- 6:
– double *Output
-
On exit: the estimate of the definite integral.
- 7:
– Nag_Comm *
-
The NAG communication argument (see
Section 3.3.1.1 in How to Use the NAG Library and its Documentation).
- 8:
– NagError *Input/Output
-
The NAG error argument (see
Section 3.7 in How to Use the NAG Library and its Documentation).
6
Error Indicators and Warnings
- NE_ALLOC_FAIL
-
Dynamic memory allocation failed.
See
Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
- NE_BAD_PARAM
-
The value of
a and/or
b is invalid for the chosen
quad_type. Either:
- On entry, argument had an illegal value.
-
The value of a and/or b is invalid.
On entry, .
On entry, and .
Constraint: .
-
The value of a and/or b is invalid.
On entry, .
On entry, and .
Constraint: .
-
The value of a and/or b is invalid.
On entry, .
On entry, and .
Constraint: .
- NE_INT
-
On entry, .
Constraint: .
- NE_INTERNAL_ERROR
-
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact
NAG for assistance.
See
Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
- NE_NO_LICENCE
-
Your licence key may have expired or may not have been installed correctly.
See
Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.
- NE_QUAD_GAUSS_NPTS_RULE
-
The -point rule is not among those stored.
On entry: .
-point rule used: .
- NE_TOO_SMALL
-
Underflow occurred in calculation of normal weights.
Reduce
n or use adjusted weights:
.
- NE_USER_STOP
-
Exit requested from
f with
.
- NE_WEIGHT_ZERO
-
No nonzero weights were generated for the provided parameters.
7
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 (d01uac) 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.
8
Parallelism and Performance
nag_quad_1d_gauss_vec (d01uac) is currently neither directly nor indirectly threaded. In particular, the user-supplied argument
f is not called from within a parallel region initialized inside
nag_quad_1d_gauss_vec (d01uac).
The user-supplied argument
f uses a vectorized interface, allowing for the required vector of function values to be evaluated in parallel; for example by placing appropriate OpenMP directives in the code for the user-supplied argument
f.
The time taken by nag_quad_1d_gauss_vec (d01uac) depends on the complexity of the expression for the integrand and on the number of abscissae required.
10
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.
10.1
Program Text
Program Text (d01uace.c)
10.2
Program Data
None.
10.3
Program Results
Program Results (d01uace.r)