FejerAlgorithm

(Source code, png)

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

Fejer Integration algorithm

Available constructors:

FejerAlgorithm(dimension=1, method=FejerAlgorithm.FEJERTYPE1)

FejerAlgorithm(discretization, method=FejerAlgorithmFEJERTYPE1)

Parameters:
dimensionint, dimension>0

The dimension of the functions to integrate. The default discretization is FejerAlgorithm-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.

methodint, optional

Integer used to select the method of integration. (Amongst ot.FejerAlgorithm.FEJERTYPE1, ot.FejerAlgorithm.FEJERTYPE2 and ot.FejerAlgorithm.CLENSHAWCURTIS).

Default is ot.FejerAlgorithm.FEJERTYPE1

Notes

The FejerAlgorithm 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 the 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 and w^{N_i}_{n_i} are the associated weight.

For any k=0,1,...,n-1, let \theta_k = \dfrac{k\pi}{n}. The Clenshaw-Curtis nodes are:

x_k = \cos(\theta_k)

for any k=0,1,...,n-1 and its associated weights are:

w_k = \dfrac{c_k}{n}\left(1-\sum_{j=1}^{\lfloor n/2\rfloor}\dfrac{b_j}{4j^2-1}\cos\left(2j\theta_k\right)\right)

where:

b_j = 
\begin{cases}
2 & \textrm{ if } j < n/2, \\
1 & \textrm{ otherwise},
\end{cases}

and:

c_k = 
\begin{cases}
1 & \textrm{ if } k = 0 \textrm{ or } n - 1, \\
2 & \textrm{ otherwise}.
\end{cases}

The type-1 Fejer quadrature rule uses the nodes:

x_k = \cos(\theta_{k + 1/2})

for any k=0,1,...,n-1 and the associated weights are:

w_k = \dfrac{2}{n}\left(1-2\sum_{j=1}^{\lfloor n/2\rfloor}\dfrac{1}{4j^2-1}\cos\left(j\theta_{2k+1}\right)\right)

Finally, the type-2 Fejer quadrature rule is very close to the Clenshaw-Curtis rule. The two methods share the same nodes (except the endpoints that are set to 0 within the Fejer method). The weights of the type-2 Fejer quadrature rule are:

w_k=\dfrac{4}{n+1}\sin\theta_k\sum_{j=1}^{\lfloor n/2\rfloor}\dfrac{\sin\left((2j-1)\theta_k\right)}{2j-1}

for any k=0,1,...,n-1.

Examples

Create a FejerAlgorithm algorithm:

>>> import openturns as ot
>>> algo = ot.FejerAlgorithm(2)
>>> algo = ot.FejerAlgorithm([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.

setName(name)

Accessor to the object's name.

integrateWithNodes

__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 FejerAlgorithm 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 FejerAlgorithm 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)

integrate(f, interval, xi)

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

The integrand function.

intervalInterval, interval \in \Rset^d

The integration domain.

xiSample

The integration nodes.

Returns:
valuePoint

Approximation of the integral.

Examples

>>> import openturns as ot
>>> f = ot.SymbolicFunction(['x'], ['sin(x)'])
>>> a = -2.5
>>> b = 4.5
>>> algoF1 = ot.FejerAlgorithm([10])
>>> value = algoF1.integrate(f, ot.Interval(a, b))[0]
>>> print(value)
-0.590...
setName(name)

Accessor to the object’s name.

Parameters:
namestr

The name of the object.