naginterfaces.library.quad.dim1_gauss_vec¶
- naginterfaces.library.quad.dim1_gauss_vec(key, a, b, n, f, data=None)[source]¶
dim1_gauss_vec
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).For full information please refer to the NAG Library document for d01ua
https://support.nag.com/numeric/nl/nagdoc_30.2/flhtml/d01/d01uaf.html
- Parameters
- keyint
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.
- afloat
The quantities and as described in the appropriate subsection of Notes.
- bfloat
The quantities and as described in the appropriate subsection of Notes.
- nint
, the number of abscissae to be used.
- fcallable (fv, iflag) = f(x, iflag, data=None)
must return the value of the integrand , or the normalized integrand , at a set of points.
- Parameters
- xfloat, ndarray, shape
The abscissae, , for at which function values are required.
- iflagint
.
- dataarbitrary, optional, modifiable in place
User-communication data for callback functions.
- Returns
- fvfloat, array-like, shape
If adjusted weights are used, the values of the integrand . , for .
Otherwise the values of the normalized integrand . , for .
- iflagint
Set if you wish to force an immediate exit from
dim1_gauss_vec
with = -1.
- dataarbitrary, optional
User-communication data for callback functions.
- Returns
- dinestfloat
The estimate of the definite integral.
- Raises
- NagValueError
- (errno )
On entry, .
Constraint: , , , , or .
- (errno )
The value of and/or is invalid.
On entry, .
On entry, and .
Constraint: .
- (errno )
The value of and/or is invalid.
On entry, .
On entry, and .
Constraint: .
- (errno )
The value of and/or is invalid.
On entry, .
On entry, and .
Constraint: .
- (errno )
On entry, .
Constraint: .
- Warns
- NagAlgorithmicWarning
- (errno )
Exit requested from with .
- (errno )
The -point rule is not among those stored.
On entry: .
-point rule used: .
- (errno )
Underflow occurred in calculation of normal weights.
Reduce or use adjusted weights: .
- (errno )
No nonzero weights were generated for the provided parameters.
Notes
General
dim1_gauss_vec
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) and 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 bewhere the are called the weights, and the the abscissae. A selection of values of is available. (See Parameters.)
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.
dim1_gauss_vec
uses a vectorized 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.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 .
One Limit Infinite
Two quadrature formulae are available for these integrals.
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 .
If the normal form is selected, the contribution of is accounted for internally, and should only return , where .
If is supplied, the interval of integration will be . Otherwise if is supplied, the interval of integration will be .
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 .
If , the interval of integration will be . Otherwise if , the interval of integration will be .
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 .
If the normal form is selected, the contribution of is accounted for internally, and should only return , where .
- 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