UQ[PY]LAB USER MANUAL
POLYNOMIAL CHAOS KRIGING
R. Sch
¨
obi, S. Marelli, B. Sudret
CHAIR OF RISK, SAFETY AND UNCERTAINTY QUANTIFICATION
STEFANO-FRANSCINI-PLATZ 5
CH-8093 Z
¨
URICH
Risk, Safety &
Uncertainty Quantification
How to cite UQ[PY]LAB
C. Lataniotis, S. Marelli, B. Sudret, Uncertainty Quantification in the cloud with UQCloud, Proceedings of the 4th
International Conference on Uncertainty Quantification in Computational Sciences and Engineering (UNCECOMP
2021), Athens, Greece, June 27–30, 2021.
How to cite this manual
R. Sch
¨
obi, S. Marelli, B. Sudret, UQ[py]Lab user manual Polynomial chaos Kriging, Report UQ[py]Lab -V1.0-109,
Chair of Risk, Safety and Uncertainty Quantification, ETH Zurich, Switzerland, 2024
BIBT
E
X entry
@TechReport{UQdoc_10_109,
author = {Sch\"obi, R. and Marelli, S. and Sudret, B.},
title = {{UQ[py]Lab user manual -- Polynomial chaos Kriging }},
institution = {Chair of Risk, Safety and Uncertainty Quantification, ETH Zurich,
Switzerland},
year = {2024},
note = {Report UQ[py]Lab - V1.0-109}
}
List of contributors:
Name Contribution
A. Hlobilov
´
a Translation from the UQLab manual
Document Data Sheet
Document Ref. UQ[PY]LAB-V1.0-109
Title: UQ[PY]LAB user manual – Polynomial chaos Kriging
Authors: R. Sch
¨
obi, S. Marelli, B. Sudret
Chair of Risk, Safety and Uncertainty Quantification, ETH Zurich,
Switzerland
Date: 27/05/2024
Doc. Version Date Comments
V1.0 27/05/2024 Initial release
Abstract
Polynomial chaos Kriging (PC-Kriging) is a novel metamodelling technique which combines
the advantages of Kriging (Gaussian process modelling) and polynomial chaos expansions
(PCE). More specifically, PC-Kriging consists of a universal Kriging model, the trend of which
is modelled by a sparse set of orthogonal polynomials.
UQ[PY]LAB metamodelling tools provide an efficient, flexible and easy-to-use PC-Kriging
module that allows one to apply state-of-the-art algorithms for different variations of PC-
Kriging on a variety of applications. This manual for the PC-Kriging metamodelling module
is divided into three parts:
A short introduction to the main concepts and techniques behind PC-Kriging, with a
selection of references to the relevant literature;
A detailed example-based guide, with the explanation of most of the available options
and methods;
A comprehensive reference list detailing all the available functionalities in UQ[PY]LAB.
Keywords: UQ[PY]LAB, metamodelling, Kriging, PC-Kriging, Polynomial Chaos Expansions,
PCE, Sparse PCE
Contents
1 Theory 1
1.1 Introduction . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 1
1.2 Basics . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 1
1.2.1 Universal Kriging . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 1
1.2.2 Polynomial chaos expansions . . . . . . . . . . . . . . . . . . . . . . . 2
1.3 Polynomial-Chaos-Kriging . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 2
1.3.1 Framework . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 2
1.3.2 Sequential PC-Kriging . . . . . . . . . . . . . . . . . . . . . . . . . . . 3
1.3.3 Optimal PC-Kriging . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 3
1.3.4 A posteriori error estimation . . . . . . . . . . . . . . . . . . . . . . . . 3
2 Usage 5
2.1 Reference problem: one-dimensional function . . . . . . . . . . . . . . . . . . 5
2.2 Problem setup . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 5
2.2.1 Full model and probabilistic input model . . . . . . . . . . . . . . . . . 5
2.3 Setup of a PC-Kriging model . . . . . . . . . . . . . . . . . . . . . . . . . . . . 6
2.3.1 Sequential PC-Kriging model . . . . . . . . . . . . . . . . . . . . . . . 6
2.3.2 Using PC-Kriging as a model (predictor) . . . . . . . . . . . . . . . . . 7
2.4 Assessing the results . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 8
2.4.1 Experimental design . . . . . . . . . . . . . . . . . . . . . . . . . . . . 9
2.4.2 PCE-trend . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 9
2.4.3 Gaussian process . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 9
2.5 Advanced options . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 10
2.5.1 Optimal PC-Kriging model . . . . . . . . . . . . . . . . . . . . . . . . . 10
2.5.2 PCE options . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 10
2.5.3 Custom set A . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 10
2.5.4 Kriging options . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 10
2.5.5 Use of a validation set . . . . . . . . . . . . . . . . . . . . . . . . . . . 11
2.6 Vector-valued models . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 11
2.7 Using PC-Kriging with constant parameters . . . . . . . . . . . . . . . . . . . 12
3 Reference List 13
3.1 Create a PC-Kriging metamodel . . . . . . . . . . . . . . . . . . . . . . . . . . 15
3.1.1 Validation Set . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 16
3.2 Accessing the results . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 16
Chapter 1
Theory
1.1 Introduction
In modern engineering and applied sciences in general, computational simulations become
more and more complex and hence more time-consuming. In fact, the evaluation of a state-
of-the-art finite element model may take hours to days. In this context, metamodelling (a.k.a
surrogate modelling) attempts to reduce the computational costs and hence to allow for more
sophisticated analyses, such as reliability analysis and design optimizations.
Polynomial-Chaos-Kriging (PC-Kriging) is a state-of-the-art metamodelling algorithm which
was developed by the authors and which is based on the well-established metamodelling
techniques Polynomial Chaos Expansions (PCE, see also UQ[PY]LAB User Manual Polyno-
mial Chaos Expansions) and Kriging (see also UQ[PY]LAB User Manual Kriging (Gaussian
process modelling)). PC-Kriging makes use of the regression-type PCE to capture the global
behaviour of the computational model as well as the interpolation-type Kriging to capture
local variations. This combination results in a metamodelling technique that is more efficient
than PCE and Kriging separately.
The recent developments of PC-Kriging are implemented in UQ[PY]LAB. Hence, this part is
intended as an overview of the relevant theory and associated literature.
1.2 Basics
1.2.1 Universal Kriging
Kriging is a stochastic interpolation algorithm which assumes that the model output y = M(x)
is a realization of a Gaussian process indexed by x D
X
R
M
. A Kriging model is described
by the following equation (Santner et al., 2003; Lataniotis et al., 2021):
y M
K
(x) = β
T
f(x) + σ
2
Z(x, ω), (1.1)
where β
T
f(x) is the mean value of the Gaussian process (a.k.a. trend), σ
2
is the variance
of Gaussian process, and Z(x, ω) is a zero mean, unit variance, stationary Gaussian process
which is characterized by a correlation function R and its hyperparameters θ.
1
UQ[PY]LAB user manual
For more details, the reader is referred to UQ[PY]LAB User Manual Kriging (Gaussian
process modelling).
1.2.2 Polynomial chaos expansions
Consider a random vector with independent components X R
M
described by the joint
probability density function (PDF) f
X
. Polynomial Chaos Expansions (PCE) approximates
the computational model output Y = M(X) by a sum of orthonormal polynomials (Xiu and
Karniadakis, 2002; Sudret, 2007):
Y M
P C
=
X
α∈A
y
α
Ψ
α
(X), (1.2)
where Ψ
α
(X) are multivariate polynomials orthonormal with respect to the input distribu-
tions f
X
, α A N
M
are multi-indices, and y
α
are the corresponding coefficients. For more
details, the reader is referred to UQ[PY]LAB User Manual – Polynomial Chaos Expansions.
1.3 Polynomial-Chaos-Kriging
1.3.1 Framework
Kriging interpolates the local variations of Y as a function of the neighbouring experimental
design points, whereas PCE approximates well the global behaviour of Y . By combining the
global and local approximation of these techniques, a more accurate metamodel is achieved.
Polynomial-Chaos-Kriging (PC-Kriging) is defined as a universal Kriging model the trend of
which consists of a set of orthonormal polynomials (Sch
¨
obi et al., 2015, 2016; Kersaudy
et al., 2015):
y M
(PCK)
(x) =
X
α∈A
y
α
Ψ
α
(X) + σ
2
Z(x, ω), (1.3)
where
P
α∈A
y
α
Ψ
α
(X) is a weighted sum of orthonormal polynomials describing the trend
of the PC-Kriging model, σ
2
and Z(x, ω) denote the variance and the zero mean, unit vari-
ance, stationary Gaussian process, respectively, as introduced in Section 1.2.1. Hence, PC-
Kriging can be interpreted as a universal Kriging model with a specific trend. In other words,
the Kriging equations described in UQ[PY]LAB User Manual Kriging (Gaussian process
modelling) are valid.
Constructing a PC-Kriging model consists of two parts: the determination of the optimal set of
polynomials contained in the trend and the calibration of the Kriging model (i.e.determining
the parameters {θ, σ
2
, y
α
}). The two parts can be combined in various ways. In UQ[PY]LAB,
Sequential PC-Kriging and Optimal PC-Kriging are implemented and presented in the follow-
ing sections (Sch
¨
obi et al., 2015, 2016).
UQ[PY]LAB-V1.0-109 - 2 -
Polynomial chaos Kriging
1.3.2 Sequential PC-Kriging
In Sequential PC-Kriging (SPCK), the set of polynomials and the Kriging metamodel are
determined sequentially. In a first step the optimal set of polynomials A is determined by
sparse PCE based on least angle regression selection (LARS) (Blatman and Sudret, 2011).
Then, then entire set A is embedded into the PC-Kriging equation Eq. (1.3) as a trend,
consisting of P = |A| regressors. Finally, the PC-Kriging metamodel is calibrated as a usual
Kriging model, including the computation of the coefficients y
α
.
1.3.3 Optimal PC-Kriging
In Optimal PC-Kriging (OPCK), the PC-Kriging model is obtained iteratively. As in SPCK, the
optimal set of polynomials is determined by LARS. The LARS algorithm results in a sparse
set of polynomials which are ranked according to their correlation to the current residual at
each LARS iterations (in decreasing order). Each polynomial is then added individually to
the trend of a PC-Kriging model. In each iteration, a new PC-Kriging model is calibrated.
In the end, P = |A| PC-Kriging models are available, which are compared by means of
their leave-one-out (LOO) error estimator (see UQ[PY]LAB User Manual – Kriging (Gaussian
process modelling)for more details). The Optimal PC-Kriging metamodel is then chosen as
the PC-Kriging model that minimizes the LOO error.
1.3.4 A posteriori error estimation
After the PC-Kriging metamodel is set up, its predictive accuracy on new data can be assessed
by using the so-called validation error. It is calculated as the relative generalization error on
an independent set of inputs and outputs [X
Val
, Y
Val
= M(X
Val
)]:
ϵ
Val
=
N 1
N
N
P
i=1
M(x
(i)
Val
) M
P C K
(x
(i)
Val
)
2
N
P
i=1
M(x
(i)
Val
) ˆµ
Y
Val
2
(1.4)
where ˆµ
Y
Val
=
1
N
N
P
i=1
M(x
(i)
Val
) is the sample mean of the validation set response. This error
measure is useful to compare the performance of different surrogate models when evaluated
on the same validation set.
UQ[PY]LAB-V1.0-109 - 3 -
Chapter 2
Usage
In this section a reference problem will be set up to showcase how each of the techniques
described in Part 1 can be deployed in UQ[PY]LAB.
2.1 Reference problem: one-dimensional function
In the context of this manual, the following one-dimensional function is used as the reference
computational model:
f(x) = x sin(x). (2.1)
The input random vector consists of a uniform random variable X U(0, 15). An example
UQ[PY]LAB script that showcases several of the currently available PC-Kriging techniques on
this function can be found in the example file:
01_XsinX.ipynb
2.2 Problem setup
2.2.1 Full model and probabilistic input model
The first step is to initialize UQ[PY]LAB and fixing a random seed for reproducibility:
# Package imports
from uqpylab import sessions
import numpy as np
# Start the session
mySession = sessions.cloud()
# (Optional) Get a convenient handle to the command line interface
uq = mySession.cli
# Reset the session
mySession.reset()
# Set the random seed for reproducibility
uq.rng(100, 'twister');
To surrogate it using UQ[PY]LAB, we need to first configure a basic MODEL object:
MOpts = {
'Type': 'Model',
5
UQ[PY]LAB user manual
'ModelFun': 'xsinx.model'
}
myModel = uq.createModel(MOpts)
The model is implemented as a Python file in xsinx.py which you can download at https:
//uqpylab.uq-cloud.io/examples.
For more details about the configuration options available for a model, please refer to the
UQ[PY]LAB User Manual – the MODEL module.
Compared to plain Kriging, PC-Kriging always requires an input model due to the existence
of PCE, even when using a given experimental design. The input variable can be defined as:
InputOpts = {
'Marginals': {
'Type': 'Uniform',
'Parameters': [0, 15]
}
}
myInput = uq.createInput(InputOpts)
For more details about the configuration options available for an INPUT object, please refer
to the UQ[PY]LAB User Manual – the INPUT module.
2.3 Setup of a PC-Kriging model
The PC-Kriging module creates a MODEL object that can be used as any other model. Its
configuration options, however, are generally more complex than for a basic Kriging or PCE
metamodel.
The basic options common to any PC-Kriging metamodelling MODEL read:
MetaOpts = {
'Type': 'Metamodel',
'MetaType': 'PCK',
}
The additional configuration options needed to properly create a PCK object in UQ[PY]LAB
are given in the following subsections.
2.3.1 Sequential PC-Kriging model
Computing a Sequential PC-Kriging (SPCK) requires then the specification of an experimen-
tal design. Here, we specify an experimental design consisting of 9 samples from a Sobol’
sequence:
X = uq.getSample(myInput, 9, 'Sobol')
Y = uq.evalModel(myModel, X)
MetaOpts['ExpDesign'] = {
'X': X.tolist(),
'Y': Y.tolist(),
}
UQ[PY]LAB-V1.0-109 - 6 -
Polynomial chaos Kriging
Further, we specify the type of PC-Kriging model and compute the metamodel:
MetaOpts['Mode'] = 'sequential'
mySPCK = uq.createModel(MetaOpts)
Once the model is created, a report with basic information about the PC-Kriging model can
be printed as follows:
uq.print(mySPCK)
which produces the following output in the Python workspace:
%--------------- PC-Kriging metamodel ---------------%
Object Name: Model 2
Input Dimension: 1
Experimental Design
Sampling: User
X size: [9x1]
Y size: [9x1]
Combination
Mode: sequential
Trend
Type: orthogonal polynomials
No. polys: 1
Gaussian Process
Corr. Handle: uq_eval_Kernel
Hyperparameters
theta: [ 8.50790 ]
Optim. method: Hybrid Genetic Algorithm
Leave-one-out error: 8.3328982e-01
%--------------- PC-Kriging metamodel ---------------%
Similarly, a visual representation of the PC-Kriging model can be obtained as follows:
uq.display(mySPCK)
which produces the image in Figure 1. Please note that such a representation is only available
for 1D models.
Further results are stored in the object mySPCK that was created by UQ[PY]LAB and the
default configurations for SPCK. We can directly access this object in order to obtain various
results and information on the metamodel. The content of PC-Kriging objects is summarized
in Section 3.2.
2.3.2 Using PC-Kriging as a model (predictor)
Regardless of the specifications of the metamodel, a PC-Kriging model can be used to predict
an output of new points in the input domain. Indeed, after a PC-Kriging MODEL object is
created in UQ[PY]LAB, it can be used like an ordinary model (for details, see the UQ[PY]LAB
User Manual – the MODEL module).
UQ[PY]LAB-V1.0-109 - 7 -
UQ[PY]LAB user manual
Figure 1: The figure created by uq.display of a SPCK-Kriging MODEL object which has
a one-dimensional input.
Consider the point x = 3, which is part of the input domain D
X
= [0, 15]. To evaluate the
SPCK metamodel on the point x, one can write:
x = 3
[ySPCK, ySPCKs2] = uq.evalModel(mySPCK, x, nargout=2)
where ySPCK and ySPCKs2 are the prediction mean value and variance of the fitted Gaussian
process at point x = 3, respectively. In this example, we obtain:
ySPCK = 0.15823773
ySPCKs2 = 0.15677765
2.4 Assessing the results
Apart from uq.display and uq.print, all results are accessible through the dictionary
mySPCK:
{
'Internal': {...},
'Name': 'Model 2',
'Type': 'uq_metamodel',
'core_component': 'model',
'Options': {...},
'ExpDesign': {
'X': [1x9 list of type float],
'Y': [1x9 list of type float]
},
'MetaType': 'PCK',
'PCK': {...},
'Error': {...}
}
UQ[PY]LAB-V1.0-109 - 8 -
Polynomial chaos Kriging
2.4.1 Experimental design
The experimental design used to calculate the PC-Kriging model can be accessed in
mySPCK['ExpDesign']:
{
'Sampling': 'Sobol',
'NSamples': 9,
'ED_Input': 'Input 2',
'X': [1x9 list of type float],
'Y': [1x9 list of type float],
'U': [1x9 list of type float],
}
Note that mySPCK['ExpDesign']['Sampling'] key has the value 'Sobol' due to the
Sobol sequence used to generate the experimental design.
2.4.2 PCE-trend
As the trend of PC-Kriging is based on orthogonal polynomials, the associated PCE model can
be accessed by the command:
myPCE = uq.extractFromModel(parentName=mySPCK['Name'],objPath='Internal.PCE')
The myPCE dictionary includes the following keys:
{
'Name': 'Model 2_xxx',
'Internal': {...},
'Type': 'uq_metamodel',
'Class': 'uq_model',
'PCE': {...},
'ExpDesign': {...},
'Error': {...},
'Options': {...},
'MetaType': 'PCE'
}
Note that the experimental design key myPCE['ExpDesign'] is identical to the one of the
PC-Kriging model in key mySPCK['ExpDesign'].
2.4.3 Gaussian process
The properties of the Gaussian process can be accessed via the key mySPCK['PCK']. As an
example, mySPCK['PCK'] includes the following keys:
{
'beta': 1.1263e+03,
'sigmaSQ': 4.1660e+07,
'theta': 8.5079
}
Note that mySPCK['Internal']['Kriging'] is of type string containing the name of the
model that can be extracted from UQCLOUD by the command:
UQ[PY]LAB-V1.0-109 - 9 -
UQ[PY]LAB user manual
myKriging =
uq.extractFromModel(parentName=mySPCK['Name'],objPath='Internal.Kriging')
and hence be used as such in prediction, too.
2.5 Advanced options
2.5.1 Optimal PC-Kriging model
An Optimal PC-Kriging model can be obtained by specifying the following option:
MetaOpts['Mode'] = 'optimal'
myOPCK = uq.createModel(MetaOpts)
The resulting metamodel can be summarized by uq.print(myOPCK) and visualized by
uq.display(myOPCK). Its usage is identical to that of the Sequential PCK.
2.5.2 PCE options
The trend of a PC-Kriging model is obtained by default by LARS and a polynomial de-
gree p = 3. However, it is possible to change the options for LARS by specifying them in
MetaOpts['PCE']. Every option related to LARS (described in UQ[PY]LAB User Manual
Polynomial Chaos Expansions) may be set to a non-default value in the corresponding key.
In the above example, we considered a degree-adaptive LARS of degrees one to ten as fol-
lows:
MetaOpts['PCE'] = {'Degree': np.arange(1,11).tolist()}
2.5.3 Custom set A
By default, the trend of the PC-Kriging is obtained by LARS. However, if a suitable set of
polynomials is available, it can be specified in MetaOpts by giving the polynomial types and
the corresponding index set. Then, no LARS will be computed.
As an example, the polynomial types are set to Legendre polynomials and index set is α =
(1, 4, 6)
T
(see also Section 2 in UQ[PY]LAB User Manual – Polynomial Chaos Expansions):
MetaOpts['PolyTypes'] = ['Legendre']
MetaOpts['PolyIndices'] = [1, 4, 6]
In case of multidimensional inputs MetaOpts['PolyIndices'] is a regular nested list of
lists of size P × M where P is the number of α-vectors and M is the size of the input vector.
2.5.4 Kriging options
By default, the Gaussian process of a PC-Kriging model is computed with the default Kriging
settings. However, it is possible to change the options by specifying them in MetaOpts['Kriging']
UQ[PY]LAB-V1.0-109 - 10 -
Polynomial chaos Kriging
with the same syntax as in a normal Kriging model (see also UQ[PY]LAB User Manual – Krig-
ing (Gaussian process modelling)). Every option related to the Kriging model (other than
those associated to the trend) may the set manually.
As an example, we consider a separable Gaussian correlation function as follows:
MetaOpts['Kriging'] = {
'Corr': {
'Type': 'separable',
'Family': 'Gaussian'
}
}
2.5.5 Use of a validation set
If a validation set is provided (see Table 3), UQ[PY]LAB automatically computes the valida-
tion error given in Eq. (1.4). To provide a validation set, the following command shall be
used:
MetaOpts['ValidationSet'] = {
'X': XVal,
'Y': YVal
}
The value of the validation error is stored in myPCK['Error']['Val'] next to the other
error measures (see Table 6) and will also be displayed when typing uq.print(myPCK).
2.6 Vector-valued models
The examples presented so far in this chapter dealt with scalar-valued models. In case the
model (or the experimental design, if manually specified) has multi-component outputs (de-
noted by N
out
), UQ[PY]LAB performs an independent PC-Kriging for each output component
on the shared experimental design. No additional configuration is needed to enable this be-
haviour.
A PC-Kriging with multi-component outputs can be found in the UQ[PY]LAB examples in:
03_MultipleOutputs.ipynb.
Running a PC-Kriging calculation on a multi-component output model will result in a multi-
component output dictionary. As an example, a model with four outputs will produce the
following list of dictionaries myPCK['PCK']:
[
{'beta': [...], 'siqmaSQ': ..., 'theta': [...]},
{'beta': [...], 'siqmaSQ': ..., 'theta': [...]},
{'beta': [...], 'siqmaSQ': ..., 'theta': [...]},
{'beta': [...], 'siqmaSQ': ..., 'theta': [...]},
]
Each element of the myPCK dictionary is functionally identical to its scalar counterpart in
Section 2.4. Similarly, the myPCE['Error'] dictionary becomes a nested list of dictionaries
myPCK['Error']:
UQ[PY]LAB-V1.0-109 - 11 -
UQ[PY]LAB user manual
[
{'LOO': ...},
{'LOO': ...},
{'LOO': ...},
{'LOO': ...}
]
2.7 Using PC-Kriging with constant parameters
In some analyses, one may need to assign a constant value to one or to a set of parameters.
When this is the case, the PC-Kriging metamodel is built by internally removing the constant
parameters from the inputs. This process is transparent to the users as they shall still eval-
uate the model using the full set of parameters (including those which were set constant).
UQ[PY]LAB will automatically and appropriately account for the set of input parameters
which were declared constant.
To set a parameter to constant, the following command can be used (See UQ[PY]LAB User
Manual – the INPUT module):
InputOpts['Marginals']['Type'] = 'Constant'
InputOpts['Marginals']['Parameters'] = value
Furthermore, when the standard deviation of a parameter is set to zero, UQ[PY]LAB auto-
matically sets this parameter’s marginal to the type Constant. For example, the following
uniformly distributed variable whose upper and lower bounds are identical is automatically
set to a constant with value 1:
InputOpts['Marginals']['Type'] = 'Uniform'
InputOpts['Marginals']['Parameters'] = [1, 1]
UQ[PY]LAB-V1.0-109 - 12 -
Chapter 3
Reference List
How to read the reference list
Python dictionaries play an important role throughout the UQ[PY]LAB syntax. They offer
a natural way to semantically group configuration options and output quantities. Due to
the complexity of the algorithms implemented, it is not uncommon to employ nested dictio-
naries to fine-tune the inputs and outputs. Throughout this reference guide, a table-based
description of the configuration dictionaries is adopted.
The simplest case is given when a value of a dictionary key is a simple value or a list:
Table X: Input
Name String A description of the field is put here
which corresponds to the following syntax:
Input = {
'Name' : 'My Input'
}
The columns, from left to right, correspond to the name, the data type and a brief description
of each key-value pair. At the beginning of each row a symbol is given to inform as to whether
the corresponding key is mandatory, optional, mutually exclusive, etc. The comprehensive
list of symbols is given in the following table:
Mandatory
Optional
Mandatory, mutually exclusive (only one of
the keys can be set)
Optional, mutually exclusive (one of them
can be set, if at least one of the group is set,
otherwise none is necessary)
13
UQ[PY]LAB user manual
When the value of one of the keys of a dictionary is a dictionary itself, a link to a table that
describes the structure of that nested dictionary is provided, as in the case of the Options
key in the following example:
Table X: Input
Name String Description
Options Table Y Description of the Options
dictionary
Table Y: Input['Options']
Key1 String Description of Key1
Key2 Double Description of Key2
In some cases, an option value gives the possibility to define further options related to that
value. The general syntax would be:
Input = {
'Option1' : 'VALUE1',
'VALUE1' : {
'Val1Opt1' : ... ,
'Val1Opt2' : ...
}
}
This is illustrated as follows:
Table X: Input
Option1 String Short description
'VALUE1' Description of 'VALUE1'
'VALUE2' Description of 'VALUE2'
VALUE1 Table Y Options for 'VALUE1'
VALUE2 Table Z Options for 'VALUE2'
Table Y: Input['VALUE1']
Val1Opt1 String Description
Val1Opt2 Float Description
Table Z: Input['VALUE2']
Val2Opt1 String Description
Val2Opt2 Float Description
UQ[PY]LAB-V1.0-109 - 14 -
Polynomial chaos Kriging
3.1 Create a PC-Kriging metamodel
Syntax
myPCK = uq.createModel(MetaOpts)
Input
The dictionary variable MetaOpts contains the following keys:
Table 1: MetaOpts
Type 'Metamodel' Select the metamodelling tool
MetaType 'PCK' Select PC-Kriging
Input INPUT object Probabilistic input model
Name String Unique identifier for the metamodel
Mode String
default: 'sequential'
Combination strategy
'sequential' Computation of Sequential
PC-Kriging model
'optimal' Computation of Optimal PC-Kriging
model
PCE Dictionary Options related to the polynomial
trend of PC-Kriging (see also
UQ[PY]LAB User Manual –
Polynomial Chaos Expansions). By
default
MetaOpts['PCE']['Degree']=
3.
PolyTypes 1 × M List of strings List of polynomial families to be used
to build the PCE basis and hence A.
If used, then also PolyIndices is
required.
PolyIndices P × M List of lists of type
Float
Set of multi-indices α to define the
polynomial set A. If used, then also
PolyTypes is required.
Kriging Dictionary Options related to the polynomial
trend of PC-Kriging (see also
UQ[PY]LAB User Manual – Kriging
(Gaussian process modelling))
FullModel MODEL object UQ[PY]LAB model used to create an
experimental design (Section 2.3)
ExpDesign Table 2 Experimental design-specific options
(Section 2.3)
UQ[PY]LAB-V1.0-109 - 15 -
UQ[PY]LAB user manual
ValidationSet Table 3 Validation set components
(Section 2.5.5)
If a model is specified, UQ[PY]LAB can automatically create an experimental design for
PC-Kriging. The available options are listed in Table 2.
Table 2: MetaOpts['ExpDesign']
Sampling String
default: 'MC'
Sampling type
'MC' Monte Carlo sampling
'LHS' Latin Hypercube sampling
'Sobol' Sobol sequence sampling
'Halton' Halton sequence sampling
Nsamples Integer The number of samples to draw. It is
required when Sampling is
specified.
X N × M List of lists of type
Float
User-defined experimental design X.
If specified, Sampling is ignored.
Y N × N
Out
List of lists of
type Float
User-defined model response Y. If
specified, Sampling is ignored.
3.1.1 Validation Set
If a validation set is provided, UQ[PY]LAB automatically calculates the validation error of
the created PCK. The required information is listed in Table 3.
Table 3: MetaOpts['ValidationSet']
X N × M List of lists of type
Float
User-specified validation set X
V al
Y N × N
Out
List of lists of
type Float
User-specified validation set response
Y
V al
3.2 Accessing the results
Syntax
myPCK = uq.createModel(MetaOpts)
Output
Regardless on the configuration options given at creation time in the MetaOpts dictionary,
all PC-Kriging metamodels share the same output dictionary, given in Table 4. Note that
UQ[PY]LAB-V1.0-109 - 16 -
Polynomial chaos Kriging
because of the similarities of PC-Kriging and Kriging, many entries are the same as in the ref-
erence list of the UQ[PY]LAB User Manual Kriging (Gaussian process modelling) and hence
some results are not explicitly elaborated here. For further details, the reader is referred to
the UQ[PY]LAB User Manual – Kriging (Gaussian process modelling), Chapter 3.
Table 4: myPCK
Name String Unique name of the PC-Kriging metamodel
Options Table 1 Copy of the MetaOpts dictionary used to
create the metamodel
PCK Table 5 Information about the final Kriging model
ExpDesign Table 7 Experimental design used for calculating the
coefficients
Error Table 6 Error estimates of the metamodels’s accuracy
Internal Table 8 Internal state of the MODEL object (useful for
debug/diagnostics)
Table 5: myPCK['PCK']
beta Float or List of
type Float and
length P
The value of β(θ) in Eq. (1.1)
sigmaSQ Float The value of σ
2
(θ) in Eq. (1.1)
theta Variable size
Float
The value of θ in Eq. (1.1)
Table 6: myPCK['Error']
LOO Float The Leave-One-Out error
Val Float Validation error (see Eq. (1.4) and
Section 2.5.5). Only available if a validation
set is provided (see Table 3).
Note: In general the fields myPCK['PCK'] and myPCK['Error'] are list of dictionaries
with length equal to the number of outputs N
out
of the metamodel.
Table 7: myPCK['ExpDesign']
NSamples Float The number of samples
Sampling String The sampling method
ED_Input INPUT object The input module that was used to generate
the experimental design (X)
UQ[PY]LAB-V1.0-109 - 17 -
UQ[PY]LAB user manual
X N × M List of
lists of type
Float
The experimental design values
U N × M List of
lists of type
Float
The experimental design values in the reduced
space
Y N × N
out
List of
lists of type
Float
The output Y that corresponds to the input X
Table 8: myPCK['Internal']
Runtime Dictionary Variables that are used during the calculation
of the PC-Kriging metamodel
Error Dictionary Internal keys related to the error of the Kriging
part of PC-Kriging (see also UQ[PY]LAB User
Manual – Kriging (Gaussian process
modelling))
Kriging MODEL object Kriging model generated with a PCE trend (see
also UQ[PY]LAB User Manual – Kriging
(Gaussian process modelling))
PCE MODEL object The PCE model generated to model the trend
of the PC-Kriging model. More details can be
found in UQ[PY]LAB User Manual –
Polynomial Chaos Expansions
NumberOfPoly Integer Number of polynomials in the trend of the
PC-Kriging model
AuxSpace INPUT object Auxiliary space where PC-Kriging is performed
TrendMethod String Method that defines the trend, either 'pce' or
'user'
Note: The internal keys of the PC-Kriging module are not intended to be accessed or
changed for most typical usage scenarios of the module.
UQ[PY]LAB-V1.0-109 - 18 -
References
Blatman, G. and B. Sudret (2011). Adaptive sparse polynomial chaos expansion based on
Least Angle Regression. Journal of Computational Physics 230, 2345–2367. 3
Kersaudy, P., B. Sudret, N. Varsier, O. Picon, and J. Wiart (2015). A new surrogate model-
ing technique combining Kriging and polynomial chaos expansions – Application to uncer-
tainty analysis in computational dosimetry. Journal of Computational Physics 286, 103–117.
2
Lataniotis, C., D. Wicaksono, S. Marelli, and B. Sudret (2021). UQLab user manual – Kriging.
Technical report, Chair of Risk, Safety & Uncertainty Quantification, ETH Zurich. Report
# UQLab-V1.4-105. 1
Santner, T., B. Williams, and W. Notz (2003). The design and analysis of computer experiments.
Springer series in Statistics. Springer. 1
Sch
¨
obi, R., B. Sudret, and S. Marelli (2016). Rare event estimation using Polynomial-Chaos-
Kriging. ASCE-ASME Journal of Risk and Uncertainty in Engineering Systems, Part A: Civil
Engineering, D4016002. 2
Sch
¨
obi, R., B. Sudret, and J. Wiart (2015). Polynomial-chaos-based Kriging. International
Journal for Uncertainty Quantification 5(2), 171–193. 2
Sudret, B. (2007). Uncertainty propagation and sensitivity analysis in mechanical models -
Contributions to structural reliability and stochastic spectral methods. Habilitation thesis,
Universit
´
e Blaise Pascal, Clermont-Ferrand, France. 2
Xiu, D. and G. E. Karniadakis (2002). The Wiener-Askey polynomial chaos for stochastic
differential equations. SIAM Journal of Scientific Computing 24(2), 619–644. 2
19