GaussLegendre¶
(Source code
, png
)
- class GaussLegendre(*args)¶
Tensorized integration algorithm of Gauss-Legendre.
- Available constructors:
GaussLegendre(dimension=1)
GaussLegendre(discretization)
- Parameters:
- dimensionint,
The dimension of the functions to integrate. The default discretization is GaussLegendre-DefaultMarginalIntegrationPointsNumber in each dimension, see
ResourceMap
.- discretizationsequence of int
The number of nodes in each dimension. The sequence must be non-empty and must contain only positive values.
Notes
The Gauss-Legendre algorithm enables to approximate the definite integral:
with , using a fixed tensorized Gauss-Legendre approximation:
where is the th node of the -points Gauss-Legendre 1D integration rule and the associated weight.
Examples
Create a Gauss-Legendre algorithm:
>>> import openturns as ot >>> algo = ot.GaussLegendre(2) >>> algo = ot.GaussLegendre([2, 4, 5])
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.
Accessor to the integration weights.
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.
setName
(name)Accessor to the object's name.
- __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.
- Returns:
- nodes
Sample
The tensorized Gauss-Legendre integration nodes on where is the dimension of the integration algorithm.
- nodes
- getWeights()¶
Accessor to the integration weights.
- Returns:
- weights
Point
The tensorized Gauss-Legendre integration weights on where is the dimension of the integration algorithm.
- weights
- 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.
- Parameters:
- Returns:
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.