GaussLegendre¶
(Source code
, png
)
- class GaussLegendre(*args)¶
Tensorized integration algorithm of Gauss-Legendre.
- Available constructors:
GaussLegendre(dimension=1)
GaussLegendre(discretization)
- Parameters:
- dimensionint,
The input dimension of the functions to integrate.
- discretizationsequence of int
The number of nodes in each dimension. The sequence must be non-empty and must contain only positive values.
The default discretization is stored as the GaussLegendre- DefaultMarginalIntegrationPointsNumber
ResourceMap
key, which is used for each dimension.
Methods
Accessor to the object's name.
Accessor to the discretization of the tensorized rule.
getName
()Accessor to the object's name.
getNodes
()Accessor to the integration nodes in .
Accessor to the integration weights for nodes in .
hasName
()Test if the object is named.
integrate
(*args)Evaluation of the integral of on an interval.
integrateWithNodes
(function, interval)Evaluation of the integral of on an interval with nodes.
setName
(name)Accessor to the object's name.
Notes
We consider a function and a domain where . The Gauss-Legendre algorithm approximates the definite integral:
(1)¶
For , let be the number of marginal nodes of the univariate Gauss-Legendre quadrature. The algorithm uses a fixed tensorized Gauss-Legendre approximation based on the tensor-product Gauss quadrature using marginal nodes. In other words, this algorithm uses an anisotropic quadrature rule, i.e. the number of nodes in each marginal dimension are not necessarily equal.
In order to get a more straightforward access to the nodes and weights, please use the
GaussProductExperiment
class. In order to create a more general tensor product quadrature, please use theTensorProductExperiment
class.Quadrature rule
We detail here the quadrature rule. As integral (1) is computed for each marginal of the integrand function, we consider now that . The approximated integral is:
(2)¶
where is the -th node computed from the -points Gauss- Legendre 1D integration rule and the associated weight.
The total number of nodes is the product of marginal number of nodes:
To use that approximation, the integration variables are scaled to the domain with the mapping function defined by:
(3)¶
where each component is for .
Let be the volume of the domain . Then, we have:
(4)¶
where we introduced the scaled function on defined by:
(5)¶
Then, the Gauss-Legendre quadrature rule is applied to the function defined on . We denote by the -th node of the -points Gauss-Legendre 1D integration rule and the associated weight. Then we have:
(6)¶
(7)¶
Then, the nodes used to compute approximation (2) are defined as follows:
(8)¶
and the weights are:
(9)¶
Polynomial exactness
The Gauss-Legendre quadrature rule is exact for polynomials up to some degree. More precisely, for any , let be the set of mono-variable polynomials of degree lower or equal to . Consider the tensor product of 1D polynomials:
Therefore the Gauss-Legendre tensorized quadrature is exact for all polynomials of the vector space:
Examples
Create a Gauss-Legendre algorithm in dimension 2.
>>> import openturns as ot >>> algo = ot.GaussLegendre(2)
Create a Gauss-Legendre algorithm in dimension 3 with 2 nodes for , 4 nodes for and 5 nodes for .
>>> algo = ot.GaussLegendre([2, 4, 5])
We show how to use this method to evaluate the following integral:
First, we use the default number of nodes as defined in the
ResourceMap
.>>> dimension = 3 >>> bounds = ot.Interval([0.0] * dimension, [1.0] * dimension) >>> polynomial = ot.SymbolicFunction(['x1', 'x2', 'x3'], ['x1^5 * x2^3 * x3^7']) >>> algo = ot.GaussLegendre(dimension) >>> computedIntegral = algo.integrate(polynomial, bounds) >>> print(computedIntegral) [0.00520833]
In the previous example, we used the default number of nodes in each dimension, as defined by the
ResourceMap
. In the next example, we set the number of marginal nodes: we use 3 nodes for , 2 nodes for , 4 nodes for . In this case, the vector space for which the Gauss-Legendre quadrature rule is exact is . In other words, we test the quadrature rule against the polynomial with maximum possible marginal degrees.>>> dimension = 3 >>> bounds = ot.Interval([0.0] * dimension, [1.0] * dimension) >>> polynomial = ot.SymbolicFunction(['x1', 'x2', 'x3'], ['x1^5 * x2^3 * x3^7']) >>> algo = ot.GaussLegendre([3, 2, 4]) >>> computedIntegral = algo.integrate(polynomial, bounds) >>> print(computedIntegral) [0.00520833]
- __init__(*args)¶
- getClassName()¶
Accessor to the object’s name.
- Returns:
- class_namestr
The object class name (object.__class__.__name__).
- getDiscretization()¶
Accessor to the discretization of the tensorized rule.
- Returns:
- discretization
Indices
The number of integration point in each dimension.
- discretization
- getName()¶
Accessor to the object’s name.
- Returns:
- namestr
The name of the object.
- getNodes()¶
Accessor to the integration nodes in .
- Returns:
- nodes
Sample
The tensorized Gauss-Legendre integration nodes on where is the dimension of the integration algorithm.
- nodes
Notes
The nodes are those associated to the scaled function: .
- getWeights()¶
Accessor to the integration weights for nodes in .
- Returns:
- weights
Point
The tensorized Gauss-Legendre integration weights on where is the dimension of the integration algorithm.
- weights
Notes
The weights are those associated to the scaled function: .
- hasName()¶
Test if the object is named.
- Returns:
- hasNamebool
True if the name is not empty.
- integrate(*args)¶
Evaluation of the integral of on an interval.
- Available usages:
integrate(f, interval)
- Parameters:
- Returns:
- value
Point
Approximation of the integral.
- value
Examples
>>> import openturns as ot >>> f = ot.SymbolicFunction(['x'], ['sin(x)']) >>> a = -2.5 >>> b = 4.5 >>> algoGL = ot.GaussLegendre([10]) >>> value = algoGL.integrate(f, ot.Interval(a, b))[0] >>> print(value) -0.590...
- integrateWithNodes(function, interval)¶
Evaluation of the integral of on an interval with nodes.
- Parameters:
- Returns:
Notes
The nodes are those associated to the function: .
Examples
>>> import openturns as ot >>> f = ot.SymbolicFunction(['x'], ['sin(x)']) >>> a = -2.5 >>> b = 4.5 >>> algoGL = ot.GaussLegendre([10]) >>> value, nodes = algoGL.integrateWithNodes(f, ot.Interval(a, b)) >>> print(value[0]) -0.590... >>> print(nodes) 0 : [ -2.40867 ] 1 : [ -2.02772 ] 2 : [ -1.37793 ] 3 : [ -0.516884 ] 4 : [ 0.47894 ] 5 : [ 1.52106 ] 6 : [ 2.51688 ] 7 : [ 3.37793 ] 8 : [ 4.02772 ] 9 : [ 4.40867 ]
- setName(name)¶
Accessor to the object’s name.
- Parameters:
- namestr
The name of the object.