KarhunenLoeveQuadratureAlgorithm

(Source code, png, hires.png, pdf)

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

Computation of Karhunen-Loeve decomposition using Quadrature approximation.

Available constructors:

KarhunenLoeveQuadratureAlgorithm(domain, covariance, experiment, basis, basisSize, mustScale, s)

KarhunenLoeveQuadratureAlgorithm(domain, covariance, marginalDegree, s)

Parameters:

domain : Domain

The domain on which the covariance model and the Karhunen-Loeve eigenfunctions (modes) are discretized.

covariance : CovarianceModel

The covariance function to decompose.

experiment : WeightedExperiment

The points and weights used in the quadrature approximation.

basis : Basis

The basis in wich the eigenfunctions are projected.

basisSize : int

The number of elements of the basis considered.

marginalDegree : int

The maximum degree to take into account in the tensorized Legendre basis.

mustScale : boolean

Flag to tell if the bounding box of the weighted experiment and the domain have to be maped or not.

s : float, \geq0

The threshold used to select the most significant eigenmodes, defined in KarhunenLoeveAlgorithm.

Notes

The Karhunen-Loeve quadrature algorithm solves the Fredholm problem associated to the covariance function C: see KarhunenLoeveAlgorithm to get the notations.

The Karhunen-Loeve quadrature approximation consists in replacing the integral by a quadrature approximation: if (\omega_\ell, \vect{\xi}_\ell)_{1 \leq \ell \leq L} is the weighted experiment (see WeightedExperiment) associated to the measure \mu, then for all functions measurable wrt \mu, we have:

\int_{\Rset^n} f(\vect{x}) \mu(\vect{x})\, d\vect{x} \simeq \sum_{\ell=1}^{\ell=L} \omega_\ell f(\vect{\xi}_\ell)

If we note \eta_{\ell}=\omega_{\ell}\dfrac{1_{\cD}(\vect{\xi}_{\ell})}{\mu(\vect{\xi}_{\ell})}, we build a more general quadrature approximation (\eta_\ell, \xi_\ell)_{1 \leq l \leq L} such that:

\int_{\cD} f(\vect{t})  \, d\vect{t} \simeq \sum_{\ell=1}^L \eta_\ell f(\vect{\xi}_\ell)

where only the points \vect{\xi}_\ell \in \cD are considered.

We introduce the matrices \mat{\Theta}=(\mat{\Theta}_{ij})_{1 \leq i \leq L, 1 \leq j \leq P} \in \cM_{Ld, Pd}(\Rset) such that \mat{\Theta}_{ij} = \theta_{j}(\vect{\xi}_i)\mat{I}_d, \mat{C}_{\ell, \ell'} =  \mat{C}(\vect{\xi}_{\ell},\vect{\xi}_{\ell'}) and \mat{W}= \mathrm{diag} (\mat{W}_{ii})_{1 \leq i \leq L} \in \cM_{Ld, Ld}(\Rset) such that \mat{W}_{ii} = \eta_i \mat{I}_d.

The normalisation constraint \|\vect{\varphi}_k\|^2_{L^2(\cD, \Rset^d)}=1 ang the orthogonality of the \vect{\varphi}_k in L^2(\cD, \Rset^d) leads to:

(1)\Tr{\vect{\Phi}} \,\Tr{\mat{\Theta}} \,\mat{W}\, \mat{\Theta} \,\vect{\Phi}=\mat{I}

The Galerkin approach leads to the following generalized eigenvalue problem:

(2)\Tr{\mat{\Theta}} \,\mat{W} \,\mat{C}  \,\mat{W} \,\mat{\Theta}\,\vect{\Phi}_k = \lambda_k  \Tr{\mat{\Theta}}\, \mat{W}\,\mat{\Theta}\,\vect{\Phi}_k

where \Tr{\mat{\Theta}}\, \mat{W} \,\mat{C} \, \mat{W}\, \mat{\Theta} and \Tr{\mat{\Theta}}\, \mat{W}\, \mat{\Theta}\in \cM_{Pd,Pd}(\Rset).

The collocation approach leads to the following generalized eigenvalue problem:

(3)\mat{C}\, \mat{W} \,\mat{\Theta}\, \vect{\Phi}_k = \lambda_k \mat{\Theta}\,\vect{\Phi}_k

Equations (2) and (3) are equivalent when \mat{\Theta} is invertible.

OpenTURNS solves the equation (2).

The second constructor is a short-hand to the first one, where basis is the tensorized Legendre basis (see OrthogonalProductPolynomialFactory and LegendreFactory), experiment is a tensorized Gauss-Legendre quadrature (see GaussProductExperiment), basisSize is equal to marginalDegree to the power the dimension of domain and mustScale is set to True.

Examples

Discretize the domain \cD and create a covariance model:

>>> import openturns as ot
>>> domain = ot.IntervalMesher([10]*2).build(ot.Interval([-1.0]*2, [1.0]*2))
>>> s = 0.01
>>> model = ot.AbsoluteExponential([1.0]*2)

Give the basis used to decompose the eigenfunctions:

here, the 10 first Legendre polynomials family:

>>> basis = ot.OrthogonalProductPolynomialFactory([ot.LegendreFactory()]*2)
>>> basisSize = 10

Create the weighted experiment of the quadrature approximation: here, a Monte Carlo experiment from the measure orthogonal wrt the Legendre polynomials family:

>>> experiment = ot.MonteCarloExperiment(basis.getMeasure(), 1000)

Create the Karhunen-Loeve Quadrature algorithm:

>>> algorithm = ot.KarhunenLoeveQuadratureAlgorithm(domain, model, experiment, basis, basisSize, True, s)

Run it!

>>> algorithm.run()
>>> result = algorithm.getResult()

Methods

getBasis() Accessor to the functional basis.
getBasisSize() Accessor to the number of elements in the functional basis.
getClassName() Accessor to the object’s name.
getCovarianceModel() Accessor to the covariance model.
getDomain() Accessor to the domain.
getExperiment() Accessor to the points and weights of the quadrature approximation.
getId() Accessor to the object’s id.
getMustScale() Accessor to scale option.
getName() Accessor to the object’s name.
getResult() Get the result structure.
getShadowedId() Accessor to the object’s shadowed id.
getThreshold() Accessor to the threshold used to select the most significant eigenmodes.
getVisibility() Accessor to the object’s visibility state.
hasName() Test if the object is named.
hasVisibleName() Test if the object has a distinguishable name.
run() Computation of the eigenvalues and eigenfunctions values at the quadrature points.
setCovarianceModel(covariance) Accessor to the covariance model.
setName(name) Accessor to the object’s name.
setShadowedId(id) Accessor to the object’s shadowed id.
setThreshold(threshold) Accessor to the limit ratio on eigenvalues.
setVisibility(visible) Accessor to the object’s visibility state.
__init__(*args)

x.__init__(…) initializes x; see help(type(x)) for signature

getBasis()

Accessor to the functional basis.

Returns:

basis : Basis

The basis in wich the eigenfunctions are projected.

getBasisSize()

Accessor to the number of elements in the functional basis.

Returns:

basisSize : int

The number of elements of the functional basis considered. The basis in wich the eigenfunctions are projected.

getClassName()

Accessor to the object’s name.

Returns:

class_name : str

The object class name (object.__class__.__name__).

getCovarianceModel()

Accessor to the covariance model.

Returns:

covModel : CovarianceModel

The covariance model.

getDomain()

Accessor to the domain.

Returns:

domain : Domain

The domain \cD_N that discretizes the domin \cD.

getExperiment()

Accessor to the points and weights of the quadrature approximation.

Returns:

experiment : WeightedExperiment

The points and weights used in the quadrature approximation.

getId()

Accessor to the object’s id.

Returns:

id : int

Internal unique identifier.

getMustScale()

Accessor to scale option.

Returns:

mustScale : boolean

Flag to tell if the bounding box of the weighted experiment and the domain have to be maped or not.

getName()

Accessor to the object’s name.

Returns:

name : str

The name of the object.

getResult()

Get the result structure.

Returns:

resKL : KarhunenLoeveResult

The structure containing all the results of the Fredholm problem.

Notes

The structure contains all the results of the Fredholm problem.

getShadowedId()

Accessor to the object’s shadowed id.

Returns:

id : int

Internal unique identifier.

getThreshold()

Accessor to the threshold used to select the most significant eigenmodes.

Returns:

s : float, positive

The threshold s.

Notes

OpenTURNS truncates the sequence (\lambda_k,  \vect{\varphi}_k)_{k \geq 1} at the index K defined in (3).

getVisibility()

Accessor to the object’s visibility state.

Returns:

visible : bool

Visibility flag.

hasName()

Test if the object is named.

Returns:

hasName : bool

True if the name is not empty.

hasVisibleName()

Test if the object has a distinguishable name.

Returns:

hasVisibleName : bool

True if the name is not empty and not the default one.

run()

Computation of the eigenvalues and eigenfunctions values at the quadrature points.

Notes

Runs the algorithm and creates the result structure KarhunenLoeveResult.

setCovarianceModel(covariance)

Accessor to the covariance model.

Parameters:

covModel : CovarianceModel

The covariance model.

setName(name)

Accessor to the object’s name.

Parameters:

name : str

The name of the object.

setShadowedId(id)

Accessor to the object’s shadowed id.

Parameters:

id : int

Internal unique identifier.

setThreshold(threshold)

Accessor to the limit ratio on eigenvalues.

Parameters:

s : float, \geq 0

The threshold s defined in (3).

setVisibility(visible)

Accessor to the object’s visibility state.

Parameters:

visible : bool

Visibility flag.