setup_multi_index_advection_diffusion_benchmark

pyapprox.benchmarks.setup_multi_index_advection_diffusion_benchmark(kle_nvars=2, kle_length_scale=0.5, kle_stdev=1, max_eval_concurrency=1, time_scenario=None, functional=None, config_values=None, source_loc=[0.25, 0.75], source_scale=0.1, source_amp=100.0, vel_vec=[1.0, 0.0], kle_mean_field=0)[source]

This benchmark is used to test methods for forward propagation of uncertainty. The forward simulation model is the transient advection-diffusion model

\[\begin{split}\frac{\partial u}{\partial t}(x,t,\rv) = \nabla\cdot\left[k(x,\rv) \nabla u(x,t,\rv)\right] -\nabla \cdot (v u(x,t,\rv))+g(x,t) &(x,t,\rv)\in D\times [0,1]\times\rvdom\\ \mathcal{B}(x,t,\rv)=0 &(x,t,\rv)\in \partial D\times[0,1]\times\rvdom\\ u(x,t,\rv)=u_0(x,\rv) & (x,t,\rv)\in D\times\{t=0\}\times\rvdom\end{split}\]

where

\[g(x,t)=\frac{100}{2\pi 0.1^2}\exp\left(-\frac{\lvert x-[0.25,0.75]^\top\rvert^2}{2\cdot 0.1^2}\right)-\frac{s_\mathrm{sink}}{2\pi h_\mathrm{sink}^2}\exp\left(-\frac{\lvert x-x_\mathrm{sink}\rvert^2}{2h_\mathrm{sink}^2}\right)\]

and \(B(x,t,z)\) enforces Robin boundary conditions, i.e.

\[K(x,\rv)\nabla u(x,t,\rv)\cdot n -0.1 u(x,t,\rv)= 0 \quad\mathrm{on} \quad\partial D\]

As with the pyapprox.benchmarks.setup_advection_diffusion_kle_inversion_benchmark() we parameterize the uncertain diffusivity with a Karhunen Loeve Expansion (KLE)

\[k(x, \rv)=\exp\left(k_0+\sum_{d=1}^D \sqrt{\lambda_d}\psi_d(x)\rv_d\right).\]

If no initial condition is provided by the user then the governing equations in pyapprox.benchmarks.setup_advection_diffusion_kle_inversion_benchmark() is used to create an initial condition, where the forcing is set to be the first term of \(g\) here. I.e. the steady state solution before the second term of \(g\) is used to remove the concentration \(u\) from the domain.

The quantity of interest \(f(z)\) is the integral of the final solution in the subdomain \(S=[0.75, 1]\times[0, 0.25]\), i.e.

\[f(z)=\int_S u(x,T,z) dx\]

This model can be evaluated using different numerical discreizations that control the two spatial mesh resolutions and the timestep. The model is evaluated by specifying the random variables and the three numerical (configuration) variables.

If not time_scenario is provided. The QoI from the steady state solution is returned.

This benchmark can be modified by changing the default keyword arguments if necessary.

Parameters:
nvarsinteger

The number of variables of the KLE

kle_length_scalefloat

The correlation length \(L_c\) of the covariance kernel

kle_sigmafloat

The standard deviation of the KLE kernel

max_eval_concurrencyinteger

The maximum number of simulations that can be run in parallel. Should be no more than the maximum number of cores on the computer being used

time_scenariodict

Options defining the transient simulation. If None a steady state problem will be solved If True the default time scenario will be used which corresponds to specifying the dictionary

time_scenario = {
    "final_time": 0.2,
    "butcher_tableau": "im_crank2",
    "deltat": 0.1,  # default will be overwritten
    "init_sol_fun": None,
    "sink": None
    }

Respectively, the entries of sink are \(s_\mathrm{sink}, h_\mathrm{sink}, x_\mathrm{sink}\), e.g. [50, 0.1, [0.75, 0.75]]. If None then the sink will be turned off. init_sol is a callable function with signature init_sol_fun(x) -> np.ndarray (nx, 1) where x is np.ndarray (nphys_vars, nx) are physical coordinates in the mesh. butcher_tableau specifies the time-stepping scheme which can be either im_beuler1 or im_crank2. final_time specifies \(T\).

functionalcallable

Function used to compute the Quantities of interest with signature

functional(sol, z) -> float

Here sol: torch.tensor (ndof) is the solution at the mesh points and z -> np.ndarray(nkle_vars, 1) is the value of the KLE coefficients that produced sol. If None the subdomain intergral of sol at the final time will be used as defined above.

config_valueslist (np.ndarray)

List with three entries (two if time_scenario=None) The first two are the values of the degrees that can be used to construct the collocation mesh in each physical direction. The third is an array of the timestep sizes that can be used to integrate the PDE in time.

source_locnp.ndarray (2)

The center of the source

source_ampfloat

The source strength \(s\)

source_scalefloat

The source width \(h\)

vel_vec: iterable (2) default [1., 0.]

The spatially independent velocity field \(v\)

kle_mean_fieldfloat (default 0)

The spatially independent mean \(k_0\) of the KLE field in log space

Returns:
benchmarkBenchmark

Object containing the benchmark attributes documented below

funcallable

The quantity of interest \(f(w)\) with signature

fun(w) -> np.ndarray

where w is a 2D np.ndarray with shape (nvars+3,nsamples) and the output is a 2D np.ndarray with shape (nsamples,1). The first nvars rows of w are realizations of the random variables. The last 3 rows are configuration variables specifying the numerical discretization of the PDE model. See config_values documentation above. This is useful for testing multi-index multi-fidelity methods.

variableIndependentMarginalsVariable

Object containing information of the joint density of the inputs z which is the tensor product of independent and identically distributed Gaussian variables \(\mathcal{N}(0,1)\).

get_num_degrees_of_freedomcallable

Function that returns the number of mesh points multiplied by the number of timesteps, with signature

get_num_degrees_of_freedom(v) -> int

where v->np.ndarray(3) are the thre configuration values specifiying the numerical discretization

config_var_transConfigureVariableTransformation

A transform that maps the configuration values to and from a canonical space.

model_ensembleModelEnsemble

Function return the quantities of interest with the signature

fun(w) -> np.ndarray

where w is a 2D np.ndarray with shape (nvars+1, nsamples) and the output is a 2D np.ndarray with shape (nsamples, 1). The first nvars rows of w are realizations of the random variables. The last row is a model ID specifying a different numerical discretization. This is useful for testing multi-fidelity approximate control variate Monte Carlo estimators.

Examples

>>> from pyapprox_dev.benchmarks.benchmarks import setup_benchmark
>>> benchmark = setup_benchmark('multi_index_advection_diffusion', nvars=2)
>>> print(benchmark.keys())
dict_keys(['fun', 'variable'])