nag_quad_1d_gauss_wgen (d01tcc) (PDF version)
d01 Chapter Contents
d01 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_quad_1d_gauss_wgen (d01tcc)

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_quad_1d_gauss_wgen (d01tcc) returns the weights (normal or adjusted) and abscissae for a Gaussian integration rule with a specified number of abscissae. Six different types of Gauss rule are allowed.

2  Specification

#include <nag.h>
#include <nagd01.h>
void  nag_quad_1d_gauss_wgen (Nag_QuadType quad_type, double a, double b, double c, double d, Integer n, double weight[], double abscis[], NagError *fail)

3  Description

nag_quad_1d_gauss_wgen (d01tcc) returns the weights wi and abscissae xi for use in the summation
S=i=1nwifxi
which approximates a definite integral (see Davis and Rabinowitz (1975) or Stroud and Secrest (1966)). The following types are provided:
(a) Gauss–Legendre
Sab fxdx,   exact for ​fx=P2n- 1x.
Constraint: b>a.
(b) Gauss–Jacobi
normal weights:
Sabb-xcx-adfxdx,   exact for ​fx=P2n-1x,
adjusted weights:
Sab fxdx,   exact for ​fx=b-xcx-ad P2n- 1x.
Constraint: c>-1, d>-1, b>a.
(c) Exponential Gauss
normal weights:
S ab x- a+b 2 c fxdx,   exact for ​fx=P2n-1x,
adjusted weights:
S ab fxdx,   exact for ​fx = x-a+b2 c P2n- 1 x.
Constraint: c>-1, b>a.
(d) Gauss–Laguerre
normal weights:
S ax-ace-bxfxdxb>0, -ax-ace-bxfxdxb<0,   exact for ​fx=P2n-1x,
adjusted weights:
S a fx dx b> 0, -a fx dx b< 0,   exact for ​fx=x-ace-bxP2n- 1x.
Constraint: c>-1, b0.
(e) Gauss–Hermite
normal weights:
S- +x-ace-b x-a 2fxdx,   exact for ​fx=P2n-1x,
adjusted weights:
S- + fxdx,   exact for ​fx=x-ac e-b x-a 2 P2n- 1x.
Constraint: c>-1, b>0.
(f) Rational Gauss
normal weights:
S ax-acx+bdfxdxa+b>0, -ax-acx+bdfxdxa+b<0,   exact for ​fx=P2n-1 1x+b ,
adjusted weights:
S a fx dx a+b> 0, -a fx dx a+b< 0,   exact for ​fx=x-acx+bd P2n- 1 1x+b .
Constraint: c>-1, d>c+1, a+b0.
In the above formulae, P2n-1x stands for any polynomial of degree 2n-1 or less in x.
The method used to calculate the abscissae involves finding the eigenvalues of the appropriate tridiagonal matrix (see Golub and Welsch (1969)). The weights are then determined by the formula
wi= j=0 n-1Pj* xi 2 -1
where Pj*x is the jth orthogonal polynomial with respect to the weight function over the appropriate interval.
The weights and abscissae produced by nag_quad_1d_gauss_wgen (d01tcc) may be passed to nag_quad_md_gauss (d01fbc), which will evaluate the summations in one or more dimensions.

4  References

Davis P J and Rabinowitz P (1975) Methods of Numerical Integration Academic Press
Golub G H and Welsch J H (1969) Calculation of Gauss quadrature rules Math. Comput. 23 221–230
Stroud A H and Secrest D (1966) Gaussian Quadrature Formulas Prentice–Hall

5  Arguments

1:     quad_typeNag_QuadTypeInput
On entry: indicates the type of quadrature rule.
quad_type=Nag_Quad_Gauss_Legendre
Gauss–Legendre, with normal weights.
quad_type=Nag_Quad_Gauss_Jacobi
Gauss–Jacobi, with normal weights.
quad_type=Nag_Quad_Gauss_Jacobi_Adjusted
Gauss–Jacobi, with adjusted weights.
quad_type=Nag_Quad_Gauss_Exponential
Exponential Gauss, with normal weights.
quad_type=Nag_Quad_Gauss_Exponential_Adjusted
Exponential Gauss, with adjusted weights.
quad_type=Nag_Quad_Gauss_Laguerre
Gauss–Laguerre, with normal weights.
quad_type=Nag_Quad_Gauss_Laguerre_Adjusted
Gauss–Laguerre, with adjusted weights.
quad_type=Nag_Quad_Gauss_Hermite
Gauss–Hermite, with normal weights.
quad_type=Nag_Quad_Gauss_Hermite_Adjusted
Gauss–Hermite, with adjusted weights.
quad_type=Nag_Quad_Gauss_Rational
Rational Gauss, with normal weights.
quad_type=Nag_Quad_Gauss_Rational_Adjusted
Rational Gauss, with adjusted weights.
Constraint: quad_type=Nag_Quad_Gauss_Legendre, Nag_Quad_Gauss_Jacobi, Nag_Quad_Gauss_Jacobi_Adjusted, Nag_Quad_Gauss_Exponential, Nag_Quad_Gauss_Exponential_Adjusted, Nag_Quad_Gauss_Laguerre, Nag_Quad_Gauss_Laguerre_Adjusted, Nag_Quad_Gauss_Hermite, Nag_Quad_Gauss_Hermite_Adjusted, Nag_Quad_Gauss_Rational or Nag_Quad_Gauss_Rational_Adjusted.
2:     adoubleInput
3:     bdoubleInput
4:     cdoubleInput
5:     ddoubleInput
On entry: the parameters a, b, c and d which occur in the quadrature formulae. c is not used if quad_type=Nag_Quad_Gauss_Legendre; d is not used unless quad_type=Nag_Quad_Gauss_Jacobi, Nag_Quad_Gauss_Jacobi_Adjusted, Nag_Quad_Gauss_Rational or Nag_Quad_Gauss_Rational_Adjusted. For some rules c and d must not be too large (see Section 6).
Constraints:
  • if quad_type=Nag_Quad_Gauss_Legendre, a<b;
  • if quad_type=Nag_Quad_Gauss_Jacobi or Nag_Quad_Gauss_Jacobi_Adjusted, a<b and c>-1.0 and d>-1.0;
  • if quad_type=Nag_Quad_Gauss_Exponential or Nag_Quad_Gauss_Exponential_Adjusted, a<b and c>-1.0;
  • if quad_type=Nag_Quad_Gauss_Laguerre or Nag_Quad_Gauss_Laguerre_Adjusted, b0.0 and c>-1.0;
  • if quad_type=Nag_Quad_Gauss_Hermite or Nag_Quad_Gauss_Hermite_Adjusted, b>0.0 and c>-1.0;
  • if quad_type=Nag_Quad_Gauss_Rational or Nag_Quad_Gauss_Rational_Adjusted, a+b0.0 and c>-1.0 and d>c+1.0.
6:     nIntegerInput
On entry: n, the number of weights and abscissae to be returned. If quad_type=Nag_Quad_Gauss_Exponential_Adjusted or Nag_Quad_Gauss_Hermite_Adjusted and c0.0, an odd value of n may raise problems (see fail.code= NE_INDETERMINATE).
Constraint: n>0.
7:     weight[n]doubleOutput
On exit: the n weights.
8:     abscis[n]doubleOutput
On exit: the n abscissae.
9:     failNagError *Input/Output
The NAG error argument (see Section 3.6 in the Essential Introduction).

6  Error Indicators and Warnings

NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_CONSTRAINT
On entry, a, b, c, or d is not in the allowed range: a=value, b=value c=value, d=value and quad_type=value.
NE_CONVERGENCE
The algorithm for computing eigenvalues of a tridiagonal matrix has failed to converge.
NE_INDETERMINATE
Exponential Gauss or Gauss–Hermite adjusted weights with n odd and c0.0. Theoretically, in these cases:
  • for c>0.0, the central adjusted weight is infinite, and the exact function fx is zero at the central abscissa;
  • for c<0.0, the central adjusted weight is zero, and the exact function fx is infinite at the central abscissa.
In either case, the contribution of the central abscissa to the summation is indeterminate.In practice, the central weight may not have overflowed or underflowed, if there is sufficient rounding error in the value of the central abscissa.The weights and abscissa returned may be usable; you must be particularly careful not to ‘round’ the central abscissa to its true value without simultaneously ‘rounding’ the central weight to zero or  as appropriate, or the summation will suffer. It would be preferable to use normal weights, if possible.Note: remember that, when switching from normal weights to adjusted weights or vice versa, redefinition of fx is involved.
NE_INT
On entry, n=value.
Constraint: n>0.
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.
NE_TOO_BIG
One or more of the weights are larger than rmax, the largest floating point number on this computer (see nag_real_largest_number (X02ALC)): rmax=value.
Possible solutions are to use a smaller value of n; or, if using adjusted weights to change to normal weights.
NE_TOO_SMALL
One or more of the weights are too small to be distinguished from zero on this machine.
The underflowing weights are returned as zero, which may be a usable approximation.
Possible solutions are to use a smaller value of n; or, if using normal weights, to change to adjusted weights.

7  Accuracy

The accuracy depends mainly on n, with increasing loss of accuracy for larger values of n. Typically, one or two decimal digits may be lost from machine accuracy with n20, and three or four decimal digits may be lost for n100.

8  Parallelism and Performance

nag_quad_1d_gauss_wgen (d01tcc) is not threaded by NAG in any implementation.
nag_quad_1d_gauss_wgen (d01tcc) 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 Users' Note for your implementation for any additional implementation-specific information.

9  Further Comments

The major portion of the time is taken up during the calculation of the eigenvalues of the appropriate tridiagonal matrix, where the time is roughly proportional to n3.

10  Example

This example returns the abscissae and (adjusted) weights for the seven-point Gauss–Laguerre formula.

10.1  Program Text

Program Text (d01tcce.c)

10.2  Program Data

Program Data (d01tcce.d)

10.3  Program Results

Program Results (d01tcce.r)

Produced by GNUPLOT 4.4 patchlevel 0 0 1 2 3 4 5 6 7 8 9 -5 0 5 10 15 20 25 Weights at Abscissae x Example Program Abscissae and Weights for the 7-point Gauss-Laguerre Formula (a=0, b=1)

nag_quad_1d_gauss_wgen (d01tcc) (PDF version)
d01 Chapter Contents
d01 Chapter Introduction
NAG Library Manual

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