hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_quad_1d_gauss_wset (d01bb)

Purpose

nag_quad_1d_gauss_wset (d01bb) returns the weights and abscissae appropriate to a Gaussian quadrature formula with a specified number of abscissae. The formulae provided are Gauss–Legendre, rational Gauss, Gauss–Laguerre and Gauss–Hermite.
Note: this function is scheduled to be withdrawn, please see d01bb in Advice on Replacement Calls for Withdrawn/Superseded Routines..

Syntax

[weight, abscis, ifail] = d01bb(d01xxx, a, b, itype, n)
[weight, abscis, ifail] = nag_quad_1d_gauss_wset(d01xxx, a, b, itype, n)

Description

nag_quad_1d_gauss_wset (d01bb) returns the weights and abscissae for use in the Gaussian quadrature of a function f(x)f(x). The quadrature takes the form
n
S = wif(xi)
i = 1
S=i=1nwif(xi)
where wiwi are the weights and xixi are the abscissae (see Davis and Rabinowitz (1975), Fröberg (1970), Ralston (1965) or Stroud and Secrest (1966)).
Weights and abscissae are available for Gauss–Legendre, rational Gauss, Gauss–Laguerre and Gauss–Hermite quadrature, and for a selection of values of nn (see Section [Parameters]).
(a) Gauss–Legendre Quadrature:
b
Sf(x)dx
a
Sabf(x)dx
where aa and bb are finite and it will be exact for any function of the form
2n1
f(x) = cixi.
i = 0
f(x)=i=0 2n-1cixi.
(b) Rational Gauss quadrature, adjusted weights:
a
Sf(x)dx(a + b > 0)  or  Sf(x)dx(a + b < 0)
a
Saf(x) dx (a+b> 0)   or   S-a f(x) dx (a+b< 0)
and will be exact for any function of the form
2n + 1
f(x) = (ci)/((x + b)i) = (i = 02n1c2n + 1i(x + b)i)/((x + b)2n + 1).
i = 2
f(x)=i=2 2n+1ci(x+b)i=i=0 2n-1c2n+1-i(x+b)i(x+b)2n+1.
(c) Gauss–Laguerre quadrature, adjusted weights:
a
Sf(x)dx(b > 0)  or  Sf(x)dx(b < 0)
a
Saf(x) dx (b> 0)   or   S-a f(x) dx (b< 0)
and will be exact for any function of the form
2n1
f(x) = ebxcixi.
i = 0
f(x)=e-bxi=0 2n-1cixi.
(d) Gauss–Hermite quadrature, adjusted weights:
+
Sf(x)dx
S- + f(x) dx
and will be exact for any function of the form
2n1
f(x) = eb(xa)2cixi(b > 0).
i = 0
f(x)=e-b (x-a) 2i=0 2n-1cixi(b>0).
(e) Gauss–Laguerre quadrature, normal weights:
a
Sebxf(x)dx(b > 0)  or  Sebxf(x)dx(b < 0)
a
Sae-bxf(x) dx (b> 0)   or   S-a e-bxf(x) dx (b< 0)
and will be exact for any function of the form
2n1
f(x) = cixi.
i = 0
f(x)=i=0 2n-1cixi.
(f) Gauss–Hermite quadrature, normal weights:
+
Seb(xa)2f(x)dx
S- + e-b (x-a) 2f(x) dx
and will be exact for any function of the form
2n1
f(x) = cixi.
i = 0
f(x)=i=0 2n-1cixi.
Note:  the Gauss–Legendre abscissae, with a = 1a=-1, b = + 1b=+1, are the zeros of the Legendre polynomials; the Gauss–Laguerre abscissae, with a = 0a=0, b = 1b=1, are the zeros of the Laguerre polynomials; and the Gauss–Hermite abscissae, with a = 0a=0, b = 1b=1, are the zeros of the Hermite polynomials.

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:     d01xxx – string
String specifying the quadrature formula to be used:
  • 'd01baz', for Gauss–Legendre weights and abscissae;
  • 'd01bay', for rational Gauss weights and abscissae;
  • 'd01bax', for Gauss–Laguerre weights and abscissae;
  • 'd01baw', for Gauss–Hermite weights and abscissae.
2:     a – double scalar
3:     b – double scalar
The quantities aa and bb as described in the appropriate sub-section of Section [Description].
4:     itype – int64int32nag_int scalar
Indicates the type of weights for Gauss–Laguerre or Gauss–Hermite quadrature (see Section [Description]).
itype = 1itype=1
Adjusted weights will be returned.
itype = 0itype=0
Normal weights will be returned.
Constraint: itype = 0itype=0 or 11.
For Gauss–Legendre or rational Gauss quadrature, this parameter is not used.
5:     n – int64int32nag_int scalar
nn, the number of weights and abscissae to be returned.
Constraint: n = 1n=1, 22, 33, 44, 55, 66, 88, 1010, 1212, 1414, 1616, 2020, 2424, 3232, 4848 or 6464.

Optional Input Parameters

None.

Input Parameters Omitted from the MATLAB Interface

None.

Output Parameters

1:     weight(n) – double array
The n weights. For Gauss–Laguerre and Gauss–Hermite quadrature, these will be the adjusted weights if itype = 1itype=1, and the normal weights if itype = 0itype=0.
2:     abscis(n) – double array
The n abscissae.
3:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and 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 ifail = 1ifail=1
The N-point rule is not among those stored. If the soft fail option is used, the weights and abscissae returned will be those for the largest valid value of n less than the requested value, and the excess elements of weight and abscis (i.e., up to the requested n) will be filled with zeros.
  ifail = 2ifail=2
The value of a and/or b is invalid.
Rational Gauss: a + b = 0.0a+b=0.0 
Gauss–Laguerre: b = 0.0b=0.0 
Gauss–Hermite: b0.0b0.0
If the soft fail option is used the weights and abscissae are returned as zero.
W ifail = 3ifail=3
Laguerre and Hermite normal weights only: underflow is occurring in evaluating one or more of the normal weights. If the soft fail option is used, the underflowing weights are returned as zero. A smaller value of n must be used; or adjusted weights should be used (itype = 1itype=1). In the latter case, take care that underflow does not occur when evaluating the integrand appropriate for adjusted weights.

Accuracy

The weights and abscissae are stored for standard values of a and b to full machine accuracy.

Further Comments

Timing is negligible.

Example

function nag_quad_1d_gauss_wset_example
a = 0;
b = 1;
itype = int64(1);
n = int64(6);
[weight, abscis, ifail] = nag_quad_1d_gauss_wset('nag_quad_1d_gauss_laguerre', a, b, itype, n)
 

weight =

    0.5735
    1.3693
    2.2607
    3.3505
    4.8868
    7.8490


abscis =

    0.2228
    1.1889
    2.9927
    5.7751
    9.8375
   15.9829


ifail =

                    0


function d01bb_example
a = 0;
b = 1;
itype = int64(1);
n = int64(6);
[weight, abscis, ifail] = d01bb('d01bax', a, b, itype, n)
 

weight =

    0.5735
    1.3693
    2.2607
    3.3505
    4.8868
    7.8490


abscis =

    0.2228
    1.1889
    2.9927
    5.7751
    9.8375
   15.9829


ifail =

                    0



PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013