UQ Forward Problem¶
In UQ forward (uncertainty propagation) problems, the aim is to estimate the propagation of uncertainties in the model response or QoIs, from known uncertain parameters/inputs.
In many situations, we are mostly interested in approximately estimating the statistical moments of the model outputs and QoIs.
In particular, among the moments, our main focus is on the expected value \(\mathbb{E}_\mathbf{q}[r]\) and variance \(\mathbb{V}_\mathbf{q}[r]\) of a QoI \(r\), with \(\mathbf{q}\) denoting the uncertain parameters.
Different approaches can be used for this purpose, see [Smith13], [Ghanem17].
The Monte Carlo approaches rely on taking independent samples from the parameters and running the model simulator as a blackbox.
However, the convergence rate of these methods can be low, see e.g. [Ghanem17].
As a more efficient technique, the spectral-based method non-intrusive generalized polynomial chaos expansion (gPCE) [Xiu02], [Xiu05], [Xiu07] is employed in UQit
for UQ forward problems.
For an overview of the technique with the same notations as in this document, refer to [Rezaeiravesh20].
Standard Polynomial Chaos Expansion¶
Theory¶
The generalized polynomial chaos expansion for \(\tilde{f}(\chi,\mathbf{q})\) is written as,
where the basis function for the multi-variate parameter \(\xi\in\Gamma\) is defined as \(\Psi_{k}(\xi)=\prod_{i=1}^p \psi_{k_i}(\xi_i)\). In the framework of generalized PCE, see [Xiu02], given the distribution of the single-variate parameter \(\xi_i\in \Gamma_i\), a set of orthogonal basis functions are chosen for \(\psi_{k_i}(\xi_i)\), for \(i=1,2,\ldots,p\). Note that, there is a one-to-one correspondence between any sample of \(\mathbf{q}\in \mathbb{Q}\) and \(\xi\in\Gamma\), where \(\mathbb{Q}=\bigotimes_{i=1}^p \mathbb{Q}_i\) and \(\Gamma=\bigotimes_{i=1}^p \Gamma_i\). The mapped space \(\Gamma_i\) is known based on the gPCE rule, see [Xiu02], [Eldred09].
Given a set of training data \(\mathcal{D}=\{(\mathbf{q}^{(i)},r^{(i)})\}_{i=1}^n\), there are two main steps to construct the above expansion.
First, a truncation scheme is needed to handle \(p\)-dimensional parameter space and determine \(K\).
Currently, tensor-product and total-order schemes are available in UQit
.
Second, the coefficients \(\{\hat{f}_k(\chi)\}_{k=0}^K\) have to be determined.
In UQit
, two different approaches can be used for this purpose: projection and regression method, see [Rezaeiravesh20] and the references therein.
In case the number of training data is less than \(K\), compressed sensing method can be adopted which is implemented in UQit
through the external Python library cvxpy
[Diamond16].
Once the coefficients \(\{\hat{f}_k(\chi)\}_{k=0}^K\) are obtained, the PCE can be used as a surrogate for the actual unobserved \(f(\chi,\mathbf{q})\). A main advantage of the PCE method is that the approximate estimation of the statistical moments of the \(f(\chi,\mathbf{q})\) or response \(r\) is a natural outcome of the surrogate construction. Using gPCE, the mean and variance of the simulator are estimated by,
where \(\gamma_k\) is the inner-product of the polynomial basis.
Example¶
In UQit
, to construct and estimate expected value and variance of \(f(\mathbf{q})\) for \(\mathbf{q}\in\mathbb{Q}\subset \mathbb{R}^p\), we have:
pce_=pce(fVal=fVal,xi=xiGrid,pceDict=pceDict,nQList=nQ)
fMean=pce_.fMean
fVar=pce_.fVar
pceCoefs=pce_.coefs
kSet=pce_.kSet
To evaluate a constructed PCE at a set of test parameter samples taken from \(\Gamma\), we write:
pcePred_=pce.pceEval(coefs=pceCoefs,xi=xiTest,distType=distType,kSet=kSet)
fPCE=pcePred_.pceVal
As for instance described in [Rezaeiravesh20], as an a-posteriori measure of the convergence of the PCE terms, we can evaluate the following indicator,
at different multi-indices \(\mathbf{k}\).
In UQit
this is done through running,
pce.convPlot(coefs=pceCoefs,distType=distType)
Implementation¶
In UQit
, the methods required for standard PCE are implemented in pce.py
.
-
class
pce.
convPlot
(coefs, distType, kSet=[], convPltOpts=[])¶ Computes and plots the convergence indicator of a PCE that is defined as,
\[\vartheta_k = ||f_k \Psi_k||/|f_0|\]versus \(|k|=\sum_{i=1}^p k_i\). The dimension of the parameter space, p, is arbitrary. But for p>1, kSet have to be provided.
- Args:
- coefs: 1D numpy array of length K
The PCE coefficients
- distType: List of length p of strings,
The i-th value specifies the distribution type of the i-th parameter (based on the gPCE rule)
- kSet: (optional, required only if p>1) List of length K
The index set produced when constructing the PCE with a specified truncation scheme. \(kSet=[[k_{1,1},k_{2,1},...k_{p,1}],...,[k_{1,K},k_{2,K},..,k_{p,K}]]\) with \(k_{i,j}\) being integer.
- convPltOpts: (optional) dict
- Containing the options to save the figure. It includes the following keys:
- ‘figDir’:
Path to the directory at which the figure is saved (if not exists, is created)
- ‘figName’:
Name of the figure
- Attributes:
- kMag: List of K integers
=|k|, sum of the PCE uni-directional indices
- pceConvIndic: 1D numpy array of size K
The PCE convergence indicator
- Methods:
- pceConv():
Computes the convergence indicator
- pceConvPlot():
Plots the PCE convergence indicator
-
pceConv
()¶ Computes the convergence indicator
-
pceConvPlot
()¶ Plots the PCE convergence indicator
-
class
pce.
pce
(fVal, xi, pceDict, nQList=[], verbose=True)¶ Constructs non-intrusive generalized Polynomial Chaos Expansion (PCE). The parameter space has dimension p. We have taken n samples from the p-D parameter space and evaluated the simulator at each sample. The samples can be taken using any approach, it is only enough to set the options correctly. The general aim is to estimate \(\hat{f}_k\) in
\[f(q)=\sum_{k=0}^K \hat{f}_k \psi_k(\xi)\]where
\[\psi_k(\xi)=\psi_{k_1}(\xi_1)\psi_{k_2}(\xi_2)\cdots\psi_{k_p}(\xi_p)\]- Args:
- fVal: 1D numpy array of size n
Simulator’s response values at n training samples.
- xi: 2D numpy array of shape (n,p)
Training parameter samples over the mapped space. NOTE: Always have to be provided unless ‘sampleType’:’GQ’ .
- nQList: (optional) List of length p,
- nQList=[nQ1,nQ2,…,nQp], where nQi: number of samples in i-th direction
nQList=[] (default) if p==1 or p>1 and ‘truncMethod’:’TO’
Required only if p>1 and ‘truncMethod’:’TP’
- pceDict: dict
- Contains settings required for constructing the PCE. The keys are:
- ‘p’: int
Dimension of the parameter
- ‘distType’: List of length p,
The i-th value specifies the distribution type of the i-th parameter (based on the gPCE rule)
- ‘sampleType’: string
- Type of parameter samples at which observations are made
‘GQ’ (Gauss quadrature nodes)
‘ ‘ (other nodal sets, see class trainSample in sampling.py)
- ‘pceSolveMethod’: string
- Method of solving for the PCE coefficients
‘Projection’: Projection method; samples have to be Gauss-quadrature nodes.
‘Regression’: Regression method for uniquely-, over-, and under-determined systems. If under-determined, compressed sensing with L1/L2 regularization is automatically used.
- ‘truncMethod’: string (mandatory only if p>1)
- Method of truncating the PCE
‘TP’: Tensor-Product method
‘TO’: Total-Order method
- ‘LMax’: int (optional)
- Maximum order of the PCE in each of the parameter dimensions.
- It is mandatory for p>1 and ‘TuncMethod’==’TO’
‘LMax’ can be used only with ‘pceSolveMethod’:’Regression’
If p==1 and ‘LMax’ is not provided, it will be assumed to be equal to n.
If p>1 and ‘LMax’ is not provided, it will be assumed to a default value.
- verbose: bool (optional)
If True (default), info is printed about the PCE being constructed
- Attributes:
- coefs: 1D numpy array of size K
Coefficients in the PCE
- fMean: scalar
PCE estimation for E[f(q)]
- fVar: scalar
PCE estimation for V[f(q)]
- kSet: List (size K) of p-D lists, p>1
Index set \([[k_{1,1},k_{2,1},...k_{p,1}],...,[k_{1,K},k_{2,K},..,k_{p,K}]]\) If p==1, then kSet=[]
-
classmethod
basis
(n_, xi_, distType_)¶ Evaluates gPCE polynomial basis of order n_ at xi_ points taken from the mapped space \(\Gamma\). The standard polynomials are chosen based on the gPCE rules.
- Args:
- n_: int
Order of the basis
- xi_: 1D numpy array of size m
Points taken from the mapped space
- distType_: string
Distribution type of the random parameter (based on the gPCE rule)
- Returns:
- psi: 1D numpy array of size m
Values of the gPCE basis at xi_
-
classmethod
basisNorm
(k_, distType_, nInteg=10000)¶ Evaluates the L2-norm of the gPCE polynomial basis of order k_ at nInteg points taken from the mapped space \(\Gamma\).
- Args:
- k_: int
Order of the gPCE basis
- distType_: string
Distribution type of the random parameter (according to the gPCE rule)
- nInteg: int (optional)
Number of points to evaluate the L2-norm integral
- Returns:
- xi_: scalar
L2-norm of the gPCE basis of order k_
-
cnstrct
()¶ Constructs a PCE over a p-D parameter space (p=1,2,3,…)
-
cnstrct_1d
()¶ Constructs a PCE over a 1D parameter space, i.e. p=1.
-
cnstrct_GQTP_pd
()¶ - Constructs a PCE over a pD parameter space (p>1) using the following settings:
‘sampType’:’GQ’ (Gauss-Quadrature nodes)
‘truncMethod’: ‘TP’ (Tensor-product)
‘pceSolveMethod’:’Projection’ or ‘Regression’
-
cnstrct_GQ_1d
()¶ Constructs a PCE over a 1D parameter space using Projection method with Gauss-quadrature nodes.
- Args:
- fVal: 1D numpy array of size n
Simulator’s response values at n training samples
- pceDict[‘distType’]: List of length 1
=[distType1], where distType1:string specifies the distribution type of the parameter based on the gPCE rule
-
cnstrct_nonGQTP_pd
()¶ - Constructs a PCE over a pD parameter space (p>1), for the following settings:
‘truncMethod’: ‘TO’ or ‘TP’
‘pceSolveMethod’: ‘Regression’ (only allowed method)
This method is used for any combination of ‘sampleType’ and ‘truncMethod’ but ‘GQ’`+’TP’`
-
cnstrct_nonGQ_1d
()¶ Constructs a PCE over a 1D parameter space using ‘Regression’ method with arbitrary truncation K=LMax to compute the PCE coefficients for arbitrarily chosen parameter samples.
- Args:
- fVal: 1D numpy array of size n
Simulator’s response values at n training samples
- pceDict[‘distType’]: List of length 1,
=[distType1], where distType1:string specifies the distribution type of the parameter based on the gPCE rule
- xi: 2D numpy array of size (n,1)
Training parameter samples over the mapped space
- pceDict[‘LMax’]: int
Maximum order of the PCE. LMax is required since ‘pceSolveMethod’==’Regression’.
-
cnstrct_pd
()¶ Construct a PCE over a p-D parameter space, where p>1.
-
classmethod
density
(xi_, distType_)¶ Evaluates the PDF of the standard gPCE random variables with distribution type distType_ at points xi_ taken from the mapped space \(\Gamma\).
- Args:
- xi_: 1D numpy array of size n
Samples taken from the mapped space
- distType_: string
Distribution type of the random parameter (according to the gPCE rule)
- Returns:
- pdf_: 1D numpy array of size n
PDF of the random parameter (according to the gPCE rule)
-
classmethod
gqPtsWts
(n_, type_)¶ Gauss quadrature nodes and weights associated to distribution type type_ based on the gPCE rule.
- Args:
- n_: int,
Order of the gPCE polynomial.
- type_: string,
Distribution of the random variable according to the gPCE rule
- Returns:
- quads: 1D numpy array of size n_
Gauss quadrature nodes
- weights: 1D numpy array of size n_
Gauss quadrature weights
-
classmethod
mapFromUnit
(xi_, xBound_)¶ Linearly maps xi_ in [-1,1] to x in xBound_
- Args:
xi_: Either a scalar or a 1D numpy array
xBound_: A list of length 2 specifying the range of x_
- Returns:
x_: Mapped value of xi_
-
classmethod
mapToUnit
(x_, xBound_)¶ Linearly maps x_ in xBound_ to xi_ in [-1,1]
- Args:
x_: Either a scalar or a 1D numpy array
xBound_: A list of length 2 specifying the range of x_
- Returns:
xi_: Mapped value of x_
-
classmethod
map_xi2q
(xi_, xAux_, distType_)¶ Maps xi_ in \(\Gamma\) to x_ in \(\mathbb{Q}\), where \(\Gamma\) is the mapped space corresponding to the standard gPCE.
- Args:
- xi_: 1D numpy array of size n
Samples taken from the mapped parameter space
- distType_: string
Distribution type of the parameter
- xAux: List of length 2
xAux==xBound (admissible range) if distType_==’Unif’, hence \(\Gamma=[-1,1]\)
xAux==[m,sdev] if disType_==’Norm’, where x_~N(m,sdev^2) and \(\Gamma=[-\infty,\infty]\)
- Returns:
- x_: 1D numpy array of size n
Mapped parameter value in \(\mathbb{Q}\)
-
class
pce.
pceEval
(coefs, xi, distType, kSet=[])¶ Evaluates a constructed PCE at test samples taken from the parameter space. The parameter space has dimension p. The number of test samples is m.
- Args:
- coefs: 1D numpy array of size K
PCE coefficients
- xi: A list of length p
xi=[xi_1,xi_2,..,xi_p], where xi_k is a 1D numpy array containing m_k test samples taken from the mapped space of the k-th parameter. Always a tensor-product grid of the test samples is constructed over the p-D space, therefore, m=m_1*m_2*…*m_p.
- distType: List of length p of strings,
The i-th value specifies the distribution type of the i-th parameter (based on the gPCE rule)
- kSet: (optional, required only if p>1) List of length K
The index set produced when constructing the PCE with a specified truncation scheme. \(kSet=[[k_{1,1},k_{2,1},...k_{p,1}],...,[k_{1,K},k_{2,K},..,k_{p,K}]]\) with \(k_{i,j}\) being integer
- Attribute:
- pceVal: 1D numpy array of size m,
Response values predicted (interpolated) by the PCE at the test samples
-
eval_1d
()¶ Evaluates a PCE over a 1D parameter space at a set of test samples xi taken from the mapped parameter space.
-
eval_pd
()¶ Evaluates a PCE over a p-D (p>1) parameter space at a set of test samples xi taken from the mapped parameter space.
Notebook¶
Try the PCE notebook to see how to use UQit
to perform standard polynomial chaos expansion (PCE). The provided examples can also be seen as a way to validate the implementation of the methods in UQit
.
Probabilistic Polynomial Chaos Expansion¶
Theory¶
The standard PCE (polynomial chaos expansion) and GPR (Gaussian process regression) are two powerful approaches for surrogate construction and metamodeling in UQ.
Combining these two approaches, probabilistic PCE is derived.
There are at least two different views to this derivation which can be found in Schobi et al. [Schobi15] and Owen [Owen17].
In UQit
, a generalization of the latter is implemented which is detailed in [Rezaeiravesh20].
Example¶
Given training parameter samples qTrain
and associated responses yTrain
with observation noise noiseSdev
, the ppce
is constructed as follows.
ppce_=ppce.ppce(qTrain,yTrain,noiseSdev,ppceDict)
optOut=ppce_.optOut
fMean_samples=ppce_.fMean_samps
fVar_samples=ppce_.fVar_samps
Implementation¶
-
class
ppce.
ppce
(qTrain, yTrain, noiseV, ppceDict)¶ Probabilistic Polynomial Chaos Expansion (PPCE) over a p-D parameter space, where p=1,2,… Model: y=f(q)+e
- Args:
- qTrain: 2D numpy array of shape (n,p)
Training samples for q
- yTrain: 1D numpy array of size n
Training observed values for y
- noiseV: 1D numpy array of size n
Standard-deviation of the noise in the training observations: e~N(0,noiseSdev)
- ppceDict: dict
- Dictionary containing controllers for PPCE, including:
- ‘nGQtest’: List of length p
Number of GQ test points in each of the p-directions
- ‘qInfo’: List of length p
- =[qInfo_1,…,qInfo_p], where qInfo_i is the information about distribution of q_i
if q_i~’Unif’, qInfo_i =[min(q_i),max(q_i)]
if q_i~’Norm’, qInfo_i =[m,v] for q~N(m,v^2)
- ‘nMC’: int
Number of independent samples drawn from the GPR to construct PCE
- ‘nIter_gpr’: int
Number of iterations for optimization of the GPR hyper-parameters
- ‘lr_gpr’: float
Learning rate for optimization of the GPR hyper-parameters
- ‘standardizeYTrain_gpr’: bool (optional, default: False)
If true, the training data are standardized by shifting by mean and scaling by sdev:
\(yStnd = (yTrain-mean(yTrain))/sdev(yTrain)\)
- ‘convPlot_gpr’: bool
If true, values of the hyper-parameters are plotted vs. iteration during the optimization.
- Attributes:
- fMean_list: 1D numpy array of size nMC
PCE estimates for the mean of f(q)
- fVar_list: 1D numpy array of size nMC
PCE estimates for the variance of f(q)
- optOut: dict
- Optional outputs for plotting using gprPlot, with the following keys:
post_f: Posterior density of f(q)
post_obs: Posterior density of y
qTest: A List of length p, =[qTest_1,qTest_2,…,qTest_p], where qTest_i is a 1D numpy array of size ppceDict[‘nGQtest’][i] containing the GQ test samples in the i-th direction of the parameter.
-
ppce_cnstrct_1d
()¶ Constructing a probabilistic PCE over a 1D parameter space
-
ppce_cnstrct_pd
()¶ Constructing a probabilistic PCE over a p-D parameter space, p>1
Notebook¶
Try this PPCE notebook to see how to use UQit
to perform probabilistic polynomial chaos expansion (PPCE).