NAG Library Routine Document

d01tdf  (dim1_gauss_wrec)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

d01tdf computes the weights and abscissae of a Gaussian quadrature rule using the method of Golub and Welsch.

2
Specification

Fortran Interface
Subroutine d01tdf ( n, a, b, c, muzero, weight, abscis, ifail)
Integer, Intent (In):: n
Integer, Intent (Inout):: ifail
Real (Kind=nag_wp), Intent (In):: a(n), muzero
Real (Kind=nag_wp), Intent (Inout):: b(n), c(n)
Real (Kind=nag_wp), Intent (Out):: weight(n), abscis(n)
C Header Interface
#include nagmk26.h
void  d01tdf_ ( const Integer *n, const double a[], double b[], double c[], const double *muzero, double weight[], double abscis[], Integer *ifail)

3
Description

A tri-diagonal system of equations is formed from the coefficients of an underlying three-term recurrence formula:
pjx=ajx+bjpj-1x-cjpj-2x  
for a set of othogonal polynomials pj induced by the quadrature. This is described in greater detail in the D01 Chapter Introduction. The user is required to specify the three-term recurrence and the value of the integral of the chosen weight function over the chosen interval.
As described in Golub and Welsch (1969) the abscissae are computed from the eigenvalues of this matrix and the weights from the first component of the eigenvectors.
LAPACK routines are used for the linear algebra to speed up computation.

4
References

Golub G H and Welsch J H (1969) Calculation of Gauss quadrature rules Math. Comput. 23 221–230

5
Arguments

1:     n – IntegerInput
On entry: n, the number of Gauss points required. The resulting quadrature rule will be exact for all polynomials of degree 2n-1.
Constraint: n>0.
2:     an – Real (Kind=nag_wp) arrayInput
On entry: a contains the coefficients aj.
3:     bn – Real (Kind=nag_wp) arrayInput/Output
On entry: b contains the coefficients bj.
On exit: elements of b are altered to make the underlying eigenvalue problem symmetric.
4:     cn – Real (Kind=nag_wp) arrayInput/Output
On entry: c contains the coefficients cj.
On exit: elements of c are altered to make the underlying eigenvalue problem symmetric.
5:     muzero – Real (Kind=nag_wp)Input
On entry: muzero contains the definite integral of the weight function for the interval of interest.
6:     weightn – Real (Kind=nag_wp) arrayOutput
On exit: weightj contains the weight corresponding to the jth abscissa.
7:     abscisn – Real (Kind=nag_wp) arrayOutput
On exit: abscisj the jth abscissa.
8:     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, if you are not familiar with this argument, the recommended value is 0. 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).
Errors or warnings detected by the routine:
ifail=1
The number of weights and abscissae requested (n) is less than 1: n=value.
ifail=4
Unexpected failure in eigenvalue computation. Please contact NAG.
ifail=5
The algorithm failed to converge. The ith diagonal was not zero: i=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

In general the computed weights and abscissae are accurate to a reasonable multiple of machine precision.

8
Parallelism and Performance

d01tdf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
d01tdf makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9
Further Comments

The weight function must be non-negative to obtain sensible results. This and the validity of muzero are not something that the routine can check, so please be particularly careful. If possible check the computed weights and abscissae by integrating a function with a function for which you already know the integral.

10
Example

This example program generates the weights and abscissae for the 4-point Gauss rules: Legendre, Chebyshev1, Chebyshev2, Jacobi, Laguerre and Hermite.

10.1
Program Text

Program Text (d01tdfe.f90)

10.2
Program Data

None.

10.3
Program Results

Program Results (d01tdfe.r)

© The Numerical Algorithms Group Ltd, Oxford, UK. 2017