NAG Library Routine Document

d01rgf (dim1_fin_gonnet_vec)

1
Purpose

d01rgf is a general purpose integrator which calculates an approximation to the integral of a function fx over a finite interval a,b:
I= ab fx dx .  
The routine is suitable as a general purpose integrator, and can be used when the integrand has singularities and infinities. In particular, the routine can continue if the subroutine f explicitly returns a quiet or signalling NaN or a signed infinity.

2
Specification

Fortran Interface
Subroutine d01rgf ( a, b, f, epsabs, epsrel, dinest, errest, nevals, iuser, ruser, ifail)
Integer, Intent (Inout):: iuser(*), ifail
Integer, Intent (Out):: nevals
Real (Kind=nag_wp), Intent (In):: a, b, epsabs, epsrel
Real (Kind=nag_wp), Intent (Inout):: ruser(*)
Real (Kind=nag_wp), Intent (Out):: dinest, errest
External:: f
C Header Interface
#include <nagmk26.h>
void  d01rgf_ (const double *a, const double *b,
void (NAG_CALL *f)(const double x[], const Integer *nx, double fv[], Integer *iflag, Integer iuser[], double ruser[]),
const double *epsabs, const double *epsrel, double *dinest, double *errest, Integer *nevals, Integer iuser[], double ruser[], Integer *ifail)

3
Description

d01rgf uses the algorithm described in Gonnet (2010). It is an adaptive algorithm, similar to the QUADPACK routine QAGS (see Piessens et al. (1983), see also d01raf) but includes significant differences regarding how the integrand is represented, how the integration error is estimated and how singularities and divergent integrals are treated. The local error estimation is described in Gonnet (2010).
d01rgf requires a subroutine to evaluate the integrand at an array of different points and is therefore amenable to parallel execution.

4
References

Gonnet P (2010) Increasing the reliability of adaptive quadrature using explicit interpolants ACM Trans. Math. software 37 26
Piessens R, de Doncker–Kapenga E, Überhuber C and Kahaner D (1983) QUADPACK, A Subroutine Package for Automatic Integration Springer–Verlag

5
Arguments

1:     a – Real (Kind=nag_wp)Input
On entry: a, the lower limit of integration.
2:     b – Real (Kind=nag_wp)Input
On entry: b, the upper limit of integration. It is not necessary that a<b.
Note: if a=b, the routine will immediately return dinest=0.0, errest=0.0 and nevals=0.
3:     f – Subroutine, supplied by the user.External Procedure
f must return the value of the integrand f at a set of points.
The specification of f is:
Fortran Interface
Subroutine f ( x, nx, fv, iflag, iuser, ruser)
Integer, Intent (In):: nx
Integer, Intent (Inout):: iflag, iuser(*)
Real (Kind=nag_wp), Intent (In):: x(nx)
Real (Kind=nag_wp), Intent (Inout):: ruser(*)
Real (Kind=nag_wp), Intent (Out):: fv(nx)
C Header Interface
#include <nagmk26.h>
void  f (const double x[], const Integer *nx, double fv[], Integer *iflag, Integer iuser[], double ruser[])
1:     xnx – Real (Kind=nag_wp) arrayInput
On entry: the abscissae, xi, for i=1,2,,nx, at which function values are required.
2:     nx – IntegerInput
On entry: the number of abscissae at which a function value is required.
3:     fvnx – Real (Kind=nag_wp) arrayOutput
On exit: fv must contain the values of the integrand f. fvi=fxi for all i=1,2,,nx.
4:     iflag – IntegerInput/Output
On entry: iflag=0.
On exit: set iflag<0 to force an immediate exit with ifail=-1.
5:     iuser* – Integer arrayUser Workspace
6:     ruser* – Real (Kind=nag_wp) arrayUser Workspace
f is called with the arguments iuser and ruser as supplied to d01rgf. You should use the arrays iuser and ruser to supply information to f.
f must either be a module subprogram USEd by, or declared as EXTERNAL in, the (sub)program from which d01rgf is called. Arguments denoted as Input must not be changed by this procedure.
4:     epsabs – Real (Kind=nag_wp)Input
On entry: the absolute accuracy required.
If epsabs is negative, epsabs is used. See Section 7.
If epsabs=0.0, only the relative error will be used.
5:     epsrel – Real (Kind=nag_wp)Input
On entry: the relative accuracy required.
If epsrel is negative, epsrel is used. See Section 7.
If epsrel=0.0, only the absolute error will be used otherwise the actual value of epsrel used by d01rgf is maxmachine precision,epsrel.
Constraint: at least one of epsabs and epsrel must be nonzero.
6:     dinest – Real (Kind=nag_wp)Output
On exit: the estimate of the definite integral f.
7:     errest – Real (Kind=nag_wp)Output
On exit: the error estimate of the definite integral f.
8:     nevals – IntegerOutput
On exit: the total number of points at which the integrand, f, has been evaluated.
9:     iuser* – Integer arrayUser Workspace
10:   ruser* – Real (Kind=nag_wp) arrayUser Workspace
iuser and ruser are not used by d01rgf, but are passed directly to f and may be used to pass information to this routine.
11:   ifail – IntegerInput/Output
On entry: ifail must be set to 0, -1 or 1. If you are unfamiliar with this argument you should refer to Section 3.4 in How to Use the NAG Library and its Documentation for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1 or 1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, because for this routine the values of the output arguments may be useful even if ifail0 on exit, the recommended value is -1. When the value -1 or 1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6
Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Note: d01rgf may return useful information for one or more of the following detected errors or warnings.
Errors or warnings detected by the routine:
ifail=1
The requested accuracy was not achieved. Consider using larger values of epsabs and epsrel.
ifail=2
The integral is probably divergent or slowly convergent.
ifail=14
Both epsabs=0.0 and epsrel=0.0.
ifail=-1
Exit requested from f with iflag=value.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.9 in How to Use the NAG Library and its Documentation for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.8 in How to Use the NAG Library and its Documentation for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 3.7 in How to Use the NAG Library and its Documentation for further information.

7
Accuracy

d01rgf cannot guarantee, but in practice usually achieves, the following accuracy:
I-dinest tol ,  
where
tol = max epsabs , epsrel × I ,  
and epsabs and epsrel are user-specified absolute and relative error tolerances. Moreover, it returns the quantity errest which, in normal circumstances, satisfies
I-dinest errest tol .  

8
Parallelism and Performance

d01rgf is currently neither directly nor indirectly threaded. In particular, the user-supplied subroutine f is not called from within a parallel region initialized inside d01rgf.
The user-supplied subroutine 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 subroutine f.

9
Further Comments

The time taken by d01rgf depends on the integrand and the accuracy required.
d01rgf is suitable for evaluating integrals that have singularities within the requested interval.
In particular, d01rgf accepts non-finite values on return from the user-supplied subroutine f, and will adapt the integration rule accordingly to eliminate such points. Non-finite values include NaNs and infinities.

10
Example

This example computes
-1 1 sinx x ln101 - x .  

10.1
Program Text

Program Text (d01rgfe.f90)

10.2
Program Data

None.

10.3
Program Results

Program Results (d01rgfe.r)