Sampling

Sampling

In UQit, different types of samples can be taken from the parameter space. From one point of view, the parameter samples are divided into training and test. To construct a surrogate or perform a UQ forward problem, we need to take training samples from a mapped or standardized space Γ and then map them onto the parameter admissible space Q. In contrast, the test samples which are, for instance, used to evaluate the constructed surrogates, are taken from Q and then are mapped onto Γ.

Available types of training samples:

  • GQ: Gauss-Quadrature nodes

    can be used with distributions Unif, Norm

  • GLL: Gauss-Lobatto-Legendre nodes

  • unifSpaced: Uniformly-spaced samples

  • unifRand: Uniformly distributed random samples

  • normRand: Gaussian distributed random samples

  • Clenshaw: Clenshaw points

  • Clenshaw-Curtis: Clenshaw-Curtis points

Available types of test samples:

  • GLL: Gauss-Lobatto-Legendre nodes

  • unifSpaced: Uniformly-spaced points

  • unifRand: Uniformly distributed random

  • normRand: Gaussian distributed random

Note that the argument qInfo appearing in sampling methods:
  • qInfo=[a,b], if the parameter is Unif over range [a,b], i.e. qU[a,b]

  • qInfo=[m,s] contains the mean m and standard-deviation s, if the parameter is Norm, i.e. qN(m,s2)

Example

tr_=trainSample(sampleType='GQ',GQdistType='Unif',qInfo=[2,3],nSamp=10)
tr_=trainSample(sampleType='NormRand',qInfo=[2,3],nSamp=10)
tr_=trainSample(sampleType='GLL',qInfo=[2,3],nSamp=10)
ts_=testSample(sampleType='unifRand',GQdistType='Unif',qBound=[-1,3],nSamp=10)
ts_=testSample(sampleType='unifRand',qBound=[-1,3],nSamp=10)
ts_=testSample(sampleType='normRand',GQdistType='Norm',qBound=[-1,3],qInfo=[0.5,2],nSamp=10)
ts_=testSample(sampleType='unifSpaced',GQdistType='Norm',qBound=[-1,3],qInfo=[0.5,2],nSamp=10)
ts_=testSample(sampleType='unifSpaced',GQdistType='Unif',qBound=[-1,3],nSamp=10)
ts_=testSample(sampleType='GLL',qBound=[-1,3],nSamp=10)

Implementation

sampling.LHS_sampling(n, xBound)

LHS (Latin Hypercube) sampler from a p-D random variable distributed uniformly. Credits: https://zmurchok.github.io/2019/03/15/Latin-Hypercube-Sampling.html

Args:
n: int

Number of samples to be taken

xBound: list of length p

=[[min(x1),max(x1)],…[min(xp),max(xp)]], where [min(xi),max(xi)] specifies the range of the i-th parameter

Returns:
x: 2D numpy array of size (n,p)

Samples taken from the p-D space with ranges xBound

class sampling.testSample(sampleType, qBound, nSamp, GQdistType='Unif', qInfo=[])

Generating test samples from a 1D paramter space using different methods. Samples of q in parameter space Q are drawn and then mapped to xi in the mapped space Gamma.

Args:
sampleType: string
Type of sample, chosen from the following list:
  • ‘GLL’: Gauss-Lobatto-Legendre nodes

  • ‘unifSpaced’: Uniformly-spaced

  • ‘unifRand’: Uniformly distributed random

  • ‘normRand’: Gaussian distributed random

GQdistType: string
Type of standard distribution in gPCE; default is ‘Unif’
  • ‘Unif’: Uniform distribution, Gamma=[-1,1]

  • ‘Norm’: Gaussian distribution, Gamma=[-infty,infty]

qInfo: List of length 2 (optional)

qInfo=[mean,sdev] only if GQdistType==’Norm’

qBound: List of length 2

Admissible range of q

nSamp: int

Number of samples to draw

Attributes:
xi: 1D numpy array of size nSamp

Samples on the mapped space Gamma

xiBound: List of length 2

Admissible range of xi

q: 1D numpy array of size nSamp

Samples q from the mapped space Q

qBound: List of length 2

Admissible range of q. It will be the same as the argument qBound if GQdistType==’Unif’

Examples:

ts1=testSample(sampleType=’unifRand’,GQdistType=’Unif’,qBound=[-1,3],nSamp=10) ts2=testSample(sampleType=’unifRand’,qBound=[-1,3],nSamp=10) ts3=testSample(sampleType=’unifSpaced’,GQdistType=’Norm’,qBound=[-1,3],qInfo=[0.5,2],nSamp=10) ts4=testSample(sampleType=’normRand’,GQdistType=’Norm’,qBound=[-1,3],qInfo=[0.5,2],nSamp=10) ts5=testSample(sampleType=’unifSpaced’,GQdistType=’Unif’,qBound=[-1,3],nSamp=10) ts6=testSample(sampleType=’GLL’,qBound=[-1,3],nSamp=10)

q2xi_map()

Linearly map q in Q to xi in Gamma

class sampling.trainSample(sampleType='', GQdistType='', qInfo=[], nSamp=0)

Generating training samples from a 1D paramter space using different methods. Samples of xi are drawn from the mapped space Gamma and are then mapped to the parameter space Q.

Args:
sampleType: string
Sample type, chosen from the following list:
  • ‘GQ’: Gauss-Quadrature nodes

  • ‘GLL’: Gauss-Lobatto-Legendre nodes

  • ‘unifSpaced’: Uniformly-spaced

  • ‘unifRand’: Uniformly distributed random

  • ‘normRand’: Gaussian distributed random

  • ‘Clenshaw’: Clenshaw points

  • ‘Clenshaw-Curtis’: Clenshaw-Curtis points

GQdistType: string (optional)
Specifies type of gPCE standard distribution if sampleType==’GQ’
  • ‘Unif’: Uniform distribution, Gamma=[-1,1]

  • ‘Norm’: Gaussian distribution, Gamma=[-infty,infty]

qInfo: List of length 2 (optional)
Information about the parameter range or distribution.
  • If q is Gaussian (‘Norm’ or ‘normRand’) => qInfo=[mean,sdev]

  • Otherwise, qInfo`=[min(q),max(q)]=admissible range of `q

nSamp: Integer

Number of samples to draw

Attributes:
xi: 1D numpy array of size nSamp

Samples drawn from the mapped space Gamma

xiBound: List of length 2

Admissible range of xi

q: 1D numpy array of size nSamp

Samples over the parameter space Q

qBound: List of length 2

Admissible range of q

w: 1D numpy array of size nSamp

Weights in Gauss-Quadrature rule only if sampleType==’GQ’

Examples:

ts1=trainSample(sampleType=’GQ’,GQdistType=’Unif’,qInfo=[2,3],nSamp=10) ts2=trainSample(sampleType=’NormRand’,qInfo=[2,3],nSamp=10) ts3=trainSample(sampleType=’GLL’,qInfo=[2,3],nSamp=10)

gqPtsWts()

Gauss-Quadrature nodes and weights according to the gPCE rule

xi2q_map()

Linearly map xi in Gamma to q in Q

Nodes

Some of the sampling methods rely on generating nodes from mathematical polynomials, for instance see [Canuto87]. The associated methods are implemented in nodes.py.

nodes.ClenshawCurtis_pts(l)

Generates Clenshaw-Curtis nodes at level l over range [0,1]

Args:
l: int

Level

Returns:
x: 1D numpy array of size n

Contains Clenshaw-Curtis points

nodes.Clenshaw_pts(n)

Generates Clenshaw points over range [-1,1]

Args:
n: int

Number of nodes

Returns
x_: 1D numpy array of size n

Clenshaw points

nodes.gllPts(n, eps=1e-08, maxIter=1000)

Generating Gauss-Lobatto-Legendre (GLL) nodes of order n using the Newton-Raphson iteration.

Args:
n: int

Number of GLL nodes

eps: float (optional)

Min error to keep the iteration running

maxIter: float (optional)

Max number of iterations

Outputs:
xi: 1D numpy array of size n

GLL nodes

w: 1D numpy array of size n

GLL weights

Reference:

Canuto C., Hussaini M. Y., Quarteroni A., Tang T. A., “Spectral Methods in Fluid Dynamics,” Section 2.3. Springer-Verlag 1987. https://link.springer.com/book/10.1007/978-3-642-84108-8