Note
Go to the end to download the full example code
Adaptive Leja Sequences
This tutorial describes how to construct a polynomial chaos expansion (PCE) of a function with uncertain parameters using Leja sequences.
First lets import necessary modules and define a function useful for estimating the error in the PCE. We will also set the random seed for reproductibility
import numpy as np
import matplotlib.pyplot as plt
from pyapprox.surrogates.interp.sparse_grid import plot_sparse_grid_2d
from pyapprox.benchmarks import setup_benchmark
def compute_l2_error(validation_samples, validation_values, pce,
relative=True):
pce_values = pce(validation_samples)
error = np.linalg.norm(pce_values-validation_values, axis=0)
if not relative:
error /= np.sqrt(validation_samples.shape[1])
else:
error /= np.linalg.norm(validation_values, axis=0)
return error
np.random.seed(1)
Our goal is to demonstrate how to use a polynomial chaos expansion (PCE) to approximate a function
c = np.array([[10, 0.01]]).T
w = np.full((2, 1), 0.25)
benchmark = setup_benchmark('genz', nvars=2, test_name='oscillatory',
coeff=(c, w))
model = benchmark.fun
variable = benchmark.variable
# benchmark = setup_benchmark('ishigami', a=7, b=0.1)
# variable = benchmark.variable
# model = benchmark.fun
Here we have intentionally set the coefficients
PCE represent the model output
where
Before starting the adaptive algorithm we will generate some test data to estimate the error in the PCE as the adaptive algorithm evolves. We will compute the error at each step using a callback function.
from pyapprox import variables
var_trans = variables.AffineTransform(variable)
validation_samples = variable.rvs(int(1e3))
validation_values = model(validation_samples)
errors = []
num_samples = []
def callback(pce):
error = compute_l2_error(validation_samples, validation_values, pce)
errors.append(error)
num_samples.append(pce.samples.shape[1])
Now we setup the options for the adaptive algorithm.
from pyapprox import surrogates
max_num_samples = 200
opts = {"method": "leja",
"options": {"max_nsamples": 100, "tol": 1e-10, "callback": callback}}
The AdaptiveLejaPCE object is used to build an adaptive Leja sequence. Before building the sequence, let us first introduce the basic concepts of Leja sequences.
A Leja sequence (LS) is essentially a doubly-greedy computation of a determinant maximization procedure. Given an existing set of nodes
In one dimension, a weighted LS can be understood without linear algebra: Let
- We omit notation indicating the dependence of
on . By iterating the above equation, one can progressively build up the Leja sequence
by recomputing and maximizing the objective function for increasing .
Traditionally Leja sequences were developed with
which is the square-root of the Christoffel function.
Note univaraite weighted Leja sequence were intially developed setting
In multiple dimensions, formulating a generalization of the univariate procedure is challenging. The following linear algebra formulation greedily maximizes the weighted Vandermonde-like determinant
The above procedure is an optimization with no known explicit solution, so constructing a Leja sequence is challenging. In [NJ2014], gradient based optimization was used to construct weighted Leja sequences. However a simpler procedure based upon LU factorization can also be used [JFNMP2019]. The simpler approach comes at a cost of slight degradation in the achieved determinant of the LS. We adopt the LU-based approach here due to its ease of implementation.
The algorithm for generating weighted Leja sequences using LU factorization is outlined in :ref:` Algorithm 1 <algo1>`. The algorithm consists of 5 steps. First a polynomial basis must be specified. The number of polynomial basis elements must be greater than or equal to the number of desired samples in the Leja sequence, i.e.
for some polynomial degree
Algorithm 1:
Require number of desired samples
, preconditioning function , basis
Choose the index set
such that Specifying an ordering of the basis
Generate set of
candidate samples Build
, , , Compute preconditioning matrix
, Compute first M pivots of LU factorization,
Once a Leja sequence
where the matrices
These two steps are carried out at each iteration of the adaptive algorithm. The PCE coefficients are used to guide refinement of the polynomial index set
In the following we use an adaptive algorithm first developed for generalized sparse grid approximation (this is discussed in another tutorial). At each iteration the algorithm identifies a number of different sets
We end this section by noting that (approximate) Fekete points are an alternative determinant-maximizing choice for interpolation points. We opt to use Leja sequences here because they are indeed a sequence, whereas a Fekete point construction is not nested.
Now we are in a position to start the adaptive process
pce = surrogates.adaptive_approximate(
benchmark.fun, benchmark.variable, "polynomial_chaos", opts).approx
[1 2] 4.1460014785669786e-15
Accuracy misleadingly appears reached because admissibility criterion is preventing new subspaces from being added to the active set
[0 3] 4.1460014785669786e-15
Accuracy misleadingly appears reached because admissibility criterion is preventing new subspaces from being added to the active set
And finally we plot the final polynomial index set
plot_sparse_grid_2d(
pce.samples, np.ones(pce.samples.shape[1]),
pce.pce.indices, pce.subspace_indices)
plt.figure()
plt.loglog(num_samples, errors, 'o-')
plt.show()
References
Total running time of the script: ( 0 minutes 0.995 seconds)