GaussLegendre

(Source code, png)

../../_images/GaussLegendre.png
class GaussLegendre(*args)

Tensorized integration algorithm of Gauss-Legendre.

Available constructors:

GaussLegendre(dimension=1)

GaussLegendre(discretization)

Parameters:
dimensionint, dimension>0

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:

\int_{\vect{a}}^\vect{b} f(\vect{t})\di{\vect{t}}

with f: \Rset^d \mapsto \Rset^p, \vect{a}, \vect{b}\in\Rset^d using a fixed tensorized Gauss-Legendre approximation:

\int_{\vect{a}}^\vect{b} f(\vect{t})\di{\vect{t}} = \sum_{\vect{n}\in \cN}\left(\prod_{i=1}^d w^{N_i}_{n_i}\right)f(\xi^{N_1}_{n_1},\dots,\xi^{N_d}_{n_d})

where \xi^{N_i}_{n_i} is the n_i th node of the N_i -points Gauss-Legendre 1D integration rule and w^{N_i}_{n_i} the associated weight.

Examples

Create a Gauss-Legendre algorithm:

>>> import openturns as ot
>>> algo = ot.GaussLegendre(2)
>>> algo = ot.GaussLegendre([2, 4, 5])

Methods

getClassName()

Accessor to the object's name.

getDiscretization()

Accessor to the discretization of the tensorized rule.

getName()

Accessor to the object's name.

getNodes()

Accessor to the integration nodes.

getWeights()

Accessor to the integration weights.

hasName()

Test if the object is named.

integrate(*args)

Evaluation of the integral of f on an interval.

integrateWithNodes(function, interval)

Evaluation of the integral of f 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:
discretizationIndices

The number of integration point in each dimension.

getName()

Accessor to the object’s name.

Returns:
namestr

The name of the object.

getNodes()

Accessor to the integration nodes.

Returns:
nodesSample

The tensorized Gauss-Legendre integration nodes on [0,1]^d where d>0 is the dimension of the integration algorithm.

getWeights()

Accessor to the integration weights.

Returns:
weightsPoint

The tensorized Gauss-Legendre integration weights on [0,1]^d where d>0 is the dimension of the integration algorithm.

hasName()

Test if the object is named.

Returns:
hasNamebool

True if the name is not empty.

integrate(*args)

Evaluation of the integral of f on an interval.

Available usages:

integrate(f, interval)

Parameters:
fFunction, f: \Rset^d \mapsto \Rset^p

The integrand function.

intervalInterval, interval \in \Rset^d

The integration domain.

Returns:
valuePoint

Approximation of the integral.

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 f on an interval.

Parameters:
fFunction, f: \Rset^d \mapsto \Rset^p

The integrand function.

intervalInterval, interval \in \Rset^d

The integration domain.

Returns:
valuePoint

Approximation of the integral.

nodesSample.

The integration nodes.

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.