Welcome to DABI_UQToolbox’s documentation!

Contents:

Created on Wed Jul 11 10:40:30 2012

@author: Daniele Bigoni (dabi@imm.dtu.dk)

UncertaintyQuantificationToolbox is the collection of tools for Uncertainty Quantification.

@copyright: 2012-2014 The Technical University of Denmark

class UQToolbox.RandomSampling.Experiments(f, params, dists, marshal_f=True, store_file='')[source]

This class is devoted to the sampling from a multi dimensional distribution and the evaluation of a function f on it.

Parameters:
  • f (function) – the function representing the experiment
  • params (object) – parameters to be passed to the function f (pickable)
  • dists (list) – list of distributions, instance of scipy.stats
  • marshal_f (bool) – whether to marshal the function f or not.
  • store_file (str) – file path where to store the computed values
get_new_samples()
Returns:the list of samples for which the experiment have not been evaluated yet.
get_results()
Returns:the list of results of the evaluation of f on the samples.
get_samples()
Returns:the list of samples for which the experiment have been evaluated already.
run(maxprocs=None, store_freq=None)

Evaluate the function f on the samples.

Parameters:
  • maxprocs (int) – the maximum number of processors available
  • store_freq (int) – number of evaluations to be executed before storing. If maxprocs > 1 then the object is store every store_freq*maxprocs evaluations.
sample(size, method='mc', append=True)

Sample from the multidimensional distribution defined by the dists.

Parameters:
  • size (int) – the size of the sample set
  • method (str) – the sampling method to be used. Options are: ‘mc’ for Monte Carlo, ‘lhc’ for Latin Hyper Cube, ‘sobol’ for Quasi Monte Carlo with Sobol sequence.
  • append (bool) – if True, the new samples are appended to the old samples. Using run() will apply the experiments on only the new samples. If False, the old samples and old results are discarded. This is possible only for methods ‘mc’ and ‘lhc’.
Note :

If the Experiments object has been reloaded using pickle, the dists parameters should be reset using set_dists.

UQToolbox.RandomSampling.LatinHyperCube(dists, N, experiment, params, paramUpdate, postProc)[source]

Run Latin Hyper Cube Simulations

..deprecated:: 0.1.5

UQToolbox.RandomSampling.MonteCarlo(dists, N, experiment, params, paramUpdate, postProc)[source]

Run Monte Carlo Simulations

..deprecated:: 0.1.5

class UQToolbox.RandomSampling.MultiDimDistribution(dists)[source]

Tensor construction of a multidimensional distribution

Parameters:dists (list) – list of distributions along each dimension.
lhc(size)

Samples size realizations from the multidimensional distribution using the Latin Hyper Cube method.

Parameters:size (int) – number of samples
Returns:list of samples
rvs(size)[source]

Samples size realizations from the multidimensional distribution. Equivalent to Monte Carlo sampling.

Parameters:size (int) – number of samples
Returns:list of samples
sobol(size, skip=None)

Samples size realizations from the multidimensional distribution using the Sobol sequence for Quasi Monte Carlo method.

Parameters:size (int) – number of samples
Returns:list of samples
UQToolbox.RandomSampling.QuasiMonteCarlo(dists, N, experiment, params, paramUpdate, postProc, skip=None)[source]

Run Quasi Monte Carlo Simulations

..deprecated:: 0.1.5

UQToolbox.RandomSampling.lhc(N, d, dists=None)[source]

Latin Hyper Cube

..deprecated:: 0.1.5

UQToolbox.RandomSampling.run_experiments(f, samples, params, paramUpdate, action)

Compute the Experiments f on the samples. The implementation uses MPI for parallel computations.

Parameters:
  • f (function) – experiment function handle. Signature: f( params )
  • samples – nd-array with the set of samples grouped by the first dimension
  • params – set of parameters to be passed to the experiment
  • action (function) – post processing action
Returns:

Array of computed values, ordered by the first dimension of the array.

..deprecated:: 0.1.5

UQToolbox.gPC.PCM_Tensor(dist, poly, distPoly, N)[source]

Generates the full tensor grid using the list of distributions of the parameters space approximated by the selected polynomials.

Parameters:
  • dist – list of n distributions representing the n-Dimensional random space
  • poly – list of n polynomials (class Spectral1D) used for approximating each of the directions in the random space
  • distPoly – list of n standard distributions associated to the polynomial selected
  • Ns – list of n integers for the order of polynomial expansion in each direction
  • allocBasis – boolean. True if you want the tensor basis.
Returns:

dictionary with the following attributes ‘x’: tensor product of collocation points in the standard distributions ‘w’: tensor product of the weights for cubature/projection rules (sum to 1 or prod(gamma0) respectively) ‘vals’: tensor product of collocation points ‘V’: tensor basis functions

UQToolbox.gPC.UQ_MEPCM_Tensor(domain, solve, params, solvePostProcessing, paramUpdate, refOrder, minRefOrd=2, maxRefOrd=2, iMax=10, gamma=0.5, theta1=0.01, theta2=0.75)[source]
DESCRIPTION:
UQ_MEPCM_Tensor(solve,params,solvePostProcessing,dists,paramUpdate,polys,distsPoly,N,method) Performs Uncertainty Quantification analysis of the system described by rhs using the Probabilistic Collocation Method by Tensor Product of the random space.
INPUT:
domain: (MEPCM_Forest object) initial domain on which to evaluate the functions solve: function handle to the experiment to be quantified. Signature: sol = solve(sample,params) params: parameters to be passed to the ODE system (these contain the uncertain parameter as well) solvePostProcessing: post processing function called after all the simulation have been run. Signature: solsEl = solvePostProcessing(solsEl). It returns properly aligned solutions per element, meaning that for the solutions in each element solsEl[i] contains an [Ncoll x NDOF] array. paramUpdate: function handle for updating the set of uncertain parameters. Signature: params = paramUpdate(params,sample) refinement: 0 = No refinement, 1 = refinement + partial limiting, 2 = refinement + final limiting iMax: maximum number of refinement iterations gamma: parameter for adaptivity refinement criteria theta1: parameter for adaptivity refinement criteria theta2: parameter for the dimensional refinement criteria
OUTPUT: dictionary with the following attributes
domain: MEPCM_Forest object solsEl: (list #<multi-dim El>) solutions at collocation points per each element meanSol: estimated mean varSol: estimated variance uHatEl: (list #<multi-dim El>) Galerkin coefficients (available only if method == ‘projection’) fEval: total number of solve evalutations
UQToolbox.gPC.UQ_PCM_Tensor(solve, params, solvePostProcessing, dists, paramUpdate, polys, distsPoly, Ns, method)[source]
DESCRIPTION:
UQ_PCM_Tensor(solve,params,solvePostProcessing,dists,paramUpdate,polys,distsPoly,N,method) Performs Uncertainty Quantification analysis of the system described by rhs using the Probabilistic Collocation Method by Tensor Product of the random space.
INPUT:
  • solve: function handle to the experiment to be quantified. Signature: sol = solve(sample,params)
  • params: parameters to be passed to the ODE system (these contain the uncertain parameter as well)
  • solvePostProcessing: post processing function called after all the simulation have been run. Signature: sols = solvePostProcessing(sols)
  • dists: list of n scipy.stats distribution functions for each uncertain parameter
  • paramUpdate: function handle for updating the set of uncertain parameters. Signature: params = paramUpdate(params,sample)
  • polys: list of n Spectral1D polynomials using for approximating the distribution of the random parameters
  • distsPoly: list of n standard distributions associated to the selected approximating polynomials
  • Ns: list of n integers for the order of polynomial expansion in each random direction
  • method: ‘projection’ use projection rule in order to compute statistics, ‘cubature’ use cubature rule in order to compute statistics
OUTPUT: dictionary with the following attributes
  • ‘Samples’: list of collocation points (samples) used for UQ
  • ‘Y’: solutions at collocation points
  • ‘YHat’: Galerkin coefficients (available only if method == ‘projection’)
  • ‘meanSol’: estimated mean
  • ‘varSol’: estimated variance
UQToolbox.gPC.UQ_SparseGrid(solve, params, solvePostProcessing, dists, paramUpdate, method, k, sym=1)[source]
DESCRIPTION:
UQ_PCM_Tensor(solve,params,solvePostProcessing,dists,paramUpdate,polys,distsPoly,N,method) Performs Uncertainty Quantification analysis of the system described by rhs using the Probabilistic Collocation Method by Tensor Product of the random space.
INPUT:
  • solve: function handle to the experiment to be quantified. Signature: sol = solve(sample,params)
  • params: parameters to be passed to the ODE system (these contain the uncertain parameter as well)
  • solvePostProcessing: post processing function called after all the simulation have been run. Signature: sols = solvePostProcessing(sols)
  • dists: list of n scipy.stats distribution functions for each uncertain parameter
  • paramUpdate: function handle for updating the set of uncertain parameters. Signature: params = paramUpdate(params,sample)
  • method: Sparse Grid method selected. See description below.
  • k: accuracy of the Sparse Grid rule
  • sym: [optional, default = 1] symmetry parameter for Sparse Grid package
OUTPUT: dictionary with the following attributes
  • ‘Samples’: list of collocation points (samples) used for UQ
  • ‘Y’: solutions at collocation points
  • ‘meanSol’: estimated mean
  • ‘varSol’: estimated variance
Description:
The methods available for sparse grid are:
  • KPU = Nested rule for unweighted integral over [0,1]
  • KPN = Nested rule for integral with Gaussian weight
  • GQU = Gaussian quadrature for unweighted integral over [0,1] (Gauss-Legendre)
  • GQN = Gaussian quadrature for integral with Gaussian weight (Gauss-Hermite)
  • func = any function provided by the user that accept level l and returns nodes n and weights w for univariate quadrature rule with polynomial exactness 2l-1 as [n w] = feval(func,level)
UQToolbox.gPC.gPC_MultiIndex(dist, poly, distPoly, N)[source]

Generates the Multi Index basis Vandermonde matrix using the list of distributions of the parameters space approximated by the selected polynomials.

Parameters:
  • dist – list of n distributions representing the n-Dimensional random space
  • poly – list of n polynomials (class Spectral1D) used for approximating each of the directions in the random space
  • distPoly – list of n standard distributions associated to the polynomial selected
  • N – maximum polynomial order in the basis
  • NI – (int) number of points in which to evaluate the basis functions for precise cubature rules (if NI=[], the interpolated vandermonde is not evaluated)
  • quadTypes – (list, Spectral1D.AVAIL_POINTS) list of n types of quadrature points among Gauss, Gauss-Lobatto and Gauss-Radau. If None, Gauss quadrature points are used.
  • left – (list,float) list of left values used by ORTHPOL for Gauss-Lobatto rules (the dimensions where the value is not used can be set to anything)
  • right – (list,float) list of left values used by ORTHPOL for Gauss-Lobatto rules (the dimensions where the value is not used can be set to anything)
  • end – (list,float) list of left values used by ORTHPOL for Gauss-Radau rules (the dimensions where the value is not used can be set to anything)
  • warnings – True if you want to read warning on memory allocation
Returns:

dictionary with the following attributes x: tensor product of collocation points in the standard distributions w: tensor product of the weights for cubature/projection rules (sum to 1 or prod(gamma0) respectively) vals: tensor product of collocation points V: Pascal’s simplex of the basis functions


Created on Wed Mar 13 09:03:41 2013

@author: Daniele Bigoni (dabi@dtu.dk)

Description

This module is used to construct High Dimensionar Model Representation using cut-HDMR based on spectral methods from the module SpectralToolbox.Spectral1D. Additionally the cut-HDMR can be used to compute the associated ANOVA-HDMR.

class UQToolbox.ModelReduction.GaussianField(C)[source]

Provided a Covariance function and discretization points, it generates the Gaussian fields.

Parameters:C (function) – Covariance function C(x1,x2)
fit(x)[source]

Given some discretization points, it finds the Cholesky factorization to be used for the simulation of the random field

Parameters:x (np.ndarray) – array containing the discretization points
transform(xi)[source]

Generates random fields from the input variables (Normally distributed)

Parameters:xi (np.ndarray) – 1 or 2 dimensional array containing the generating variables. If 2-dimensional, each row should contain a different sample.
Returns:the realizations of the random field
Return type:np.ndarray
class UQToolbox.ModelReduction.KLExpansion(C, d, spans)[source]

KLExpansion: Computes the 1D KL expansion of the covariance matrix C, using the Legendre-Gauss-Loabatto rule for the solution to the generalized eigenvalue problem A * u = lmb * M * u. The KLExpansion of d-dimensional random fields is obtained by tensor product.

Parameters:
  • C (function) – Covariance function C(x1,x2)
  • d (int) – dimension of the tensorized field
  • spans (list) – list of tuples containing the spans spans of the coordinates of the domain of C(x1,x2)
fit(N, targetVar=0.95)[source]

Sets up the KL-expansion eigenvalues and eigenvectors for the required tolerance.

Parameters:
  • N (list) – list of number of discretization points per dimension
  • targetVar (float) – total variance to be represented by the output KL-expansion
invscale(i, x)[source]

Scales values x in self.spans[i] to [-1,1].

scale(i, x)[source]

Scales values x in [-1,1] to self.spans[i].

transform(xs, xi)[source]
Params list xs:list of d coordinates for each direction
Params list xi:list of 1 dimensional arrays of self.nvars values xi
UQToolbox.ModelReduction.KLExpansion1D(C, N, span, targetVar=0.95)[source]
Note :Deprecated. Use the class KLExpansion instead
KLExpansion: Computes the 1D KL expansion of the covariance matrix C, using the
Legendre-Gauss-Loabatto rule for the solution to the generalized eigenvalue problem A * u = lmb * M * u
Parameters:
  • C (function) – Covariance function C(x1,x2)
  • N (int) – number of discretization points
  • span (tuple) – span of the domain of C(x1,x2)
  • targetVar (float) – total variance to be represented by the output KL-expansion
UQToolbox.sobol_lib.scramble(X)[source]

Scramble function as in Owen (1997)

Reference:

[1]Saltelli, A., Chan, K., Scott, E.M., “Sensitivity Analysis”
UQToolbox.sobol_lib.scrambled_sobol_generate(k, N, skip, leap)[source]

Scramble function as in Owen (1997)

Reference:

[1]Saltelli, A., Chan, K., Scott, E.M., “Sensitivity Analysis”
UQToolbox.sobol_lib.sobol_generate(k, N, skip, leap)[source]

Skip and leap sobol sequence

Reference:

[1]Saltelli, A., Chan, K., Scott, E.M., “Sensitivity Analysis”

Indices and tables

Table Of Contents

This Page