queens.example_simulator_functions package#

Example simulator functions.

Modules containing example simulator functions for testing and benchmarking.

example_simulator_function_by_name(function_name)[source]#

Get example simulator function by name.

Parameters:

function_name (str) – Name of the example simulator function

Returns:

func – Function

Submodules#

queens.example_simulator_functions.agawal09 module#

Agawal09 function.

[1]: Agarwal, N., & Aluru, N. R. (2009). A domain adaptive stochastic collocation approach for analysis of MEMS under uncertainties. Journal of Computational Physics, 228(20), 7662?7688. http://doi.org/10.1016/j.jcp.2009.07.014

agawal09a(x1, x2, a1=0.5, a2=0.5)[source]#

Compute the Agawal09a function.

Two dimensional benchmark function for UQ approaches proposed in [1].

The function is defined as follows:

\[\begin{split}f({\bf x}) = \begin{cases} 0 & \textrm{if } x_1 > \alpha_1 \textrm{ or } x_2 > \alpha_2 \\ \sin(\pi x_1 )\sin(\pi x_2 ) & \textrm{otherwise} \end{cases}\end{split}\]

Distribution of the input random variables is probably uniform on [0,1].

Parameters:
  • x1 (float) – First input parameter

  • x2 (float) – Second input parameter

  • a1 (float) – Coefficient (optional), with default value 0.5

  • a2 (float) – Coefficient (optional), with default value 0.5

Returns:

float – Value of agawal09a function

queens.example_simulator_functions.borehole83 module#

Borehole function.

borehole83_hifi(rw, r, tu, hu, tl, hl, l, kw)[source]#

High-fidelity version of Borehole benchmark function.

Very simple and quick to evaluate eight dimensional function, that models water flow through a borehole. Frequently used function for testing a wide variety of methods in computer experiments, see, e.g., [1]-[10].

The high-fidelity version is defined as:

\(f_{hifi}({\\bf x}) = \frac{2 \pi T_u(H_u-H_l)}{\ln(\frac{r}{r_w}) \left(1 + \frac{2LT_u}{\ln(\frac{r}{r_w})r_w^2K_w} \right)+ \frac{T_u}{T_l}}\)

For the purpose of multi-fidelity simulation, Xiong et al. (2013) [8] use the following function for the lower fidelity code:

\(f_{lofi}({\bf x}) = \frac{5 T_u (H_u-H_l)}{\ln(\frac{r}{r_w})(1.5) + \frac{2 L T_u}{\ln(\frac{r}{r_w})r_w^2 K_w} + \frac{T_u}{T_l}}\)

For the purposes of uncertainty quantification, the distributions of the input random variables are often chosen as:

rw ~ N(0.10,0.0161812)
r ~ Lognormal(7.71,1.0056)
tu ~ Uniform[63070, 115600]
hu ~ Uniform[990, 1110]
tl ~ Uniform[63.1, 116]
hl ~ Uniform[700, 820]
l ~ Uniform[1120, 1680]
kw ~ Uniform[9855, 12045]
Parameters:
  • rw (float) – Radius of borehole \((m)\) [0.05, 0.15]

  • r (float) – Radius of influence \((m)\) [100, 50000]

  • tu (float) – Transmissivity of upper aquifer \((\frac{m^2}{yr})\) [63070, 115600]

  • hu (float) – Potentiometric head of upper aquifer \((m)\) [990, 1110]

  • tl (float) – Transmissivity of lower aquifer \((\frac{m^2}{yr})\) [63.1, 116]

  • hl (float) – Potentiometric head of lower aquifer \((m)\) [700, 820]

  • l (float) – Length of borehole \((m)\) [1120, 1680]

  • kw (float) – Hydraulic conductivity of borehole \((\frac{m}{yr})\) [9855, 12045]

Returns:

float – The response is water flow rate, in \(\frac{m^3}{yr}\)

References

[1] An, J., & Owen, A. (2001). Quasi-regression. Journal of Complexity,

17(4), 588-607.

[2] Gramacy, R. B., & Lian, H. (2012). Gaussian process single-index models

as emulators for computer experiments. Technometrics, 54(1), 30-41.

[3] Harper, W. V., & Gupta, S. K. (1983). Sensitivity/uncertainty analysis

of a borehole scenario comparing Latin Hypercube Sampling and deterministic sensitivity approaches (No. BMI/ONWI-516). Battelle Memorial Inst., Columbus, OH (USA). Office of Nuclear Waste Isolation.

[4] Joseph, V. R., Hung, Y., & Sudjianto, A. (2008). Blind kriging:

A new method for developing metamodels. Journal of mechanical design, 130, 031102.

[5] Moon, H. (2010). Design and Analysis of Computer Experiments for

Screening Input Variables (Doctoral dissertation, Ohio State University).

[6] Moon, H., Dean, A. M., & Santner, T. J. (2012). Two-stage

sensitivity-based group screening in computer experiments. Technometrics, 54(4), 376-387.

[7] Morris, M. D., Mitchell, T. J., & Ylvisaker, D. (1993).

Bayesian design and analysis of computer experiments: use of derivatives in surface prediction. Technometrics, 35(3), 243-255.

[8] Xiong, S., Qian, P. Z., & Wu, C. J. (2013). Sequential design and

analysis of high-accuracy and low-accuracy computer codes. Technometrics, 55(1), 37-46.

[9] Worley, B. A. (1987). Deterministic uncertainty analysis

(No. CONF-871101-30). Oak Ridge National Lab., TN (USA).

[10] Zhou, Q., Qian, P. Z., & Zhou, S. (2011). A simple approach to

emulation for computer models with qualitative and quantitative factors. Technometrics, 53(3).

For further information, see also http://www.sfu.ca/~ssurjano/borehole.html

borehole83_lofi(rw, r, tu, hu, tl, hl, l, kw)[source]#

Low-fidelity version of Borehole benchmark function.

Very simple and quick to evaluate eight dimensional function that models water flow through a borehole. Frequently used function for testing a wide variety of methods in computer experiments.

The low-fidelity version is defined as in [1] as:

\(f_{lofi}({\bf x}) = \frac{5 T_u (H_u-H_l)}{\ln(\frac{r}{r_w})(1.5) + \frac{2 L T_u}{\ln(\frac{r}{r_w})r_w^2 K_w} + \frac{T_u}{T_l}}\)

For the purposes of uncertainty quantification, the distributions of the input random variables are often chosen as:

rw ~ N(0.10, 0.0161812)
r ~ Lognormal(7.71, 1.0056)
tu ~ Uniform[63070, 115600]
hu ~ Uniform[990, 1110]
tl ~ Uniform[63.1, 116]
hl ~ Uniform[700, 820]
l ~ Uniform[1120, 1680]
kw ~ Uniform[9855, 12045]
Parameters:
  • rw (float) – Radius of borehole \((m)\) [0.05, 0.15]

  • r (float) – Radius of influence \((m)\) [100, 50000]

  • tu (float) – Transmissivity of upper aquifer \((\frac{m^2}{yr})\) [63070, 115600]

  • hu (float) – Potentiometric head of upper aquifer \((m)\) [990, 1110]

  • tl (float) – Transmissivity of lower aquifer \((\frac{m^2}{yr})\) [63.1, 116]

  • hl (float) – Potentiometric head of lower aquifer \((m)\) [700, 820]

  • l (float) – Length of borehole \((m)\) [1120, 1680]

  • kw (float) – Hydraulic conductivity of borehole \((\frac{m}{yr})\) [9855, 12045]

Returns:

float – The response is water flow rate, in \((\frac{m^3}{yr})\)

References

[1] Xiong, S., Qian, P. Z., & Wu, C. J. (2013). Sequential design and

analysis of high-accuracy and low-accuracy computer codes. Technometrics, 55(1), 37-46.

queens.example_simulator_functions.branin78 module#

Branin function.

Originates back to [1].

[1]: Dixon, L. C. W., & Szego, G. P. (1978). The global optimization problem: an introduction. Towards global optimization, 2, 1-15.

branin78_hifi(x1, x2)[source]#

High-fidelity Branin function.

Compute value of high fidelity version of Branin function as described in [1]. The corresponding medium- and low-fidelity versions are implemented in branin_medfi and branin_lofi, respectively.

Parameters:
  • x1 (float) – First input parameter [−5, 10]

  • x2 (float) – Second input parameter [0, 15]

Returns:

float – Value of high-fidelity Branin function

References

[1] Perdikaris, P. et al., 2017. Nonlinear information fusion algorithms

for data-efficient multi-fidelity modelling. Proceedings of the Royal Society of London A: Mathematical, Physical and Engineering Sciences, 473(2198), pp.20160751–16.

branin78_lofi(x1, x2)[source]#

Compute the value of the low-fidelity Branin function.

This function computes the value of a low-fidelity version of the Branin function, using a medium-fidelity variant as described in [1]. The corresponding high- and medium-fidelity versions are implemented in branin_hifi and branin_medfi, respectively.

Parameters:
  • x1 (float) – First input parameter.

  • x2 (float) – Second input parameter.

Returns:

float – Value of the low-fidelity Branin function.

References

[1] Perdikaris, P. et al., 2017. Nonlinear information fusion algorithms for data-efficient

multi-fidelity modelling. Proceedings of the Royal Society of London A: Mathematical, Physical and Engineering Sciences, 473(2198), pp.20160751–16.

branin78_medfi(x1, x2)[source]#

Medium fidelity Branin function.

Compute the value of a medium-fidelity version of Branin function as described in [1]. The corresponding high- and low-fidelity versions are implemented in branin_hifi and branin_lofi, respectively.

Parameters:
  • x1 (float) – First input parameter

  • x2 (float) – Second input parameter

Returns:

float – Value of medium fidelity Branin function

References

[1] Perdikaris, P. et al., 2017. Nonlinear information fusion algorithms

for data-efficient multi-fidelity modelling. Proceedings of the Royal Society of London A: Mathematical, Physical and Engineering Sciences, 473(2198), pp.20160751–16.

queens.example_simulator_functions.currin88 module#

Currin functions.

currin88_hifi(x1, x2)[source]#

High-fidelity version of the Currin88 benchmark function.

Simple two-dimensional example which appears several times in the literature, see, e.g., [1]-[3].

The high-fidelity version is defined as follows:

\(f_{hifi}({\bf x}) = \left[1 - \exp(\frac{1}{2 x_2}) \right] \frac{2300x_1^3 + 1900 x_1^2 + 2092x_1 +60}{100x_1^3 + 500x_1^2 +4 x_1+20}\)

[3] proposed the following low-fidelity version of the function:

\(f_{lofi}({\bf x}) = \frac{1}{4}[f_{hifi}(x_1+0.05,x_2+0.05) + f_{hifi}(x_1+0.05,max(0,x_2-0.05))] + \frac{1}{4}[f_{hifi}(x_1-0.05,x_2+0.05) + f_{hifi}(x_1-0.05,max(0,x_2-0.05))]\)

Parameters:
  • x1 (float) – Input parameter 1 in [0,1]

  • x2 (float) – Input parameter 2 in [0,1]

Returns:

float – Value of the currin88 function

References

[1] Currin, C., Mitchell, T., Morris, M., & Ylvisaker, D. (1988).

A Bayesian approach to the design and analysis of computer experiments. Technical Report 6498. Oak Ridge National Laboratory.

[2] Currin, C., Mitchell, T., Morris, M., & Ylvisaker, D. (1991).

Bayesian prediction of deterministic functions, with applications to the design and analysis of computer experiments. Journal of the American Statistical Association, 86(416), 953-963.

[3] Xiong, S., Qian, P. Z., & Wu, C. J. (2013). Sequential design and

analysis of high-accuracy and low-accuracy computer codes. Technometrics, 55

currin88_lofi(x1, x2)[source]#

Low-fidelity version of the Currin88 benchmark function.

Simple two-dimensional example which appears several times in the literature, see e.g. [1]-[3].

The low-fidelity version is defined as follows [3]:

\(f_{lofi}({\bf x}) = \frac{1}{4} \left[f_{hifi}(x_1+0.05,x_2+0.05) + f_{hifi}(x_1+0.05,max(0,x_2-0.05)) \right] + \frac{1}{4} \left[f_{hifi}(x_1-0.05,x_2+0.05) + f_{hifi}(x_1-0.05,max(0,x_2-0.05)) \right]\)

Parameters:
  • x1 (float) – Input parameter 1 in [0,1]

  • x2 (float) – Input parameter 2 in [0,1]

Returns:

float – Value of the currin88 function

References

[1] Currin, C., Mitchell, T., Morris, M., & Ylvisaker, D. (1988).

A Bayesian approach to the design and analysis of computer experiments. Technical Report 6498. Oak Ridge National Laboratory.

[2] Currin, C., Mitchell, T., Morris, M., & Ylvisaker, D. (1991).

Bayesian prediction of deterministic functions, with applications to the design and analysis of computer experiments. Journal of the American Statistical Association, 86(416), 953-963.

[3] Xiong, S., Qian, P. Z., & Wu, C. J. (2013). Sequential design and

analysis of high-accuracy and low-accuracy computer codes. Technometrics, 55(1), 37-46.

queens.example_simulator_functions.executable_park91a_hifi_on_grid_with_gradients module#

High-fidelity Park91a with x3 and x4 as fixed coordinates as executable.

main(run_type, params)[source]#

Interface to Park91a test function.

Parameters:
  • run_type (str) –

    Run type for the main function:

    • ’s’: standard run without gradients

    • ’p’: provided gradient run

    • ’a’: adjoint run that solves the adjoint equation

  • params (dict) – Dictionary with parameters

Returns:
  • function_output (np.array) – Value of function input parameters

  • evaluated_gradient_expression (np.array) – Value of gradient input parameters Note: This can be different gradients depending on the run_type!

park91a_hifi_coords(x1, x2, x3, x4)[source]#

High-fidelity Park91a function.

High-fidelity Park91a function with x3 and x4 as fixed coordinates. Coordinates are prescribed in the main function of this module.

Simple four-dimensional benchmark function as proposed in [1] to mimic a computer model. For the purpose of multi-fidelity simulation, [3] defined a corresponding lower fidelity function, which is implemented in park91a_lofi.

The high-fidelity version is defined as:

\(f({\bf x}) = \frac{x_1}{2} \left[\sqrt{1+(x_2+x_3^2)\frac{x_4}{x_1^2}}-1 \right]+(x_1+3x_4)\exp[1-\sin(x_3)]\)

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

  • x3 (float) – Input parameter 3 [0,1)

  • x4 (float) – Input parameter 4 [0,1)

Returns:

float – Value of function at parameters

References

[1] Park, J.-S.(1991). Tuning complex computer codes to data and optimal

designs, Ph.D. Thesis

[2] Cox, D. D., Park, J.-S., & Singer, C. E. (2001). A statistical method

for tuning a computer code to a database. Computational Statistics & Data Analysis, 37(1), 77?92. http://doi.org/10.1016/S0167-9473(00)00057-8

[3] Xiong, S., Qian, P., & Wu, C. (2013). Sequential design and analysis of

high-accuracy and low-accuracy computer codes. Technometrics. http://doi.org/10.1080/00401706.2012.723572

read_input_file(input_file_path)[source]#

Read-in input from csv file.

run_type_adjoint(x3_vec, x4_vec, params)[source]#

Run that only solves the adjoint problem.

Parameters:
  • x3_vec (np.array) – Vector of x3 values from grid points

  • x4_vec (np.array) – Vector of x4 values from grid points

  • params (dict) – Dictionary with input parameters

Returns:
  • y_vec (np.array) – Empty dummy vector of function values at grid points

  • do_dx (np.array) – Vector of gradient values of the objective function at grid points

run_type_provided_gradient(x3_vec, x4_vec, params)[source]#

Run with provided gradients.

Parameters:
  • x3_vec (np.array) – Vector of x3 values from grid points

  • x4_vec (np.array) – Vector of x4 values from grid points

  • params (dict) – Dictionary with input parameters

Returns:
  • y_vec (np.array) – Vector of function values at grid points

  • y_grad (np.array) – Vector of gradient values at grid points

run_type_standard(x3_vec, x4_vec, params)[source]#

Run standard function without gradients.

Parameters:
  • x3_vec (np.array) – Vector of x3 values from grid points

  • x4_vec (np.array) – Vector of x4 values from grid points

  • params (dict) – Dictionary with input parameters

Returns:
  • y_vec (np.array) – Vector of function values at grid points

  • y_grad (np.array) – Empty dummy vector of gradient values at grid points

write_results(output, output_path)[source]#

Write solution to csv files.

queens.example_simulator_functions.gardner14a module#

Gradner2014a function.

This is a two-dimensional benchmark function for constraint Bayesian optimization.

gardner14a(x1, x2)[source]#

Gradner2014a function.

Two-dimensional benchmark function for constraint Bayesian optimization [1]:

\(f({\bf x}) = \cos(2x_1)\cos(x_2)+ \sin(x_1)\)

with the corresponding constraint function:

\(c({\bf x}) = \cos(x_1)\cos(x_2) - \sin(x_1)\sin(x_2)\)

with

\(c({\bf x}) \leq 0.5\)

Parameters:
  • x1 (float) – Input parameter 1 in [0, 6]

  • x2 (float) – Input parameter 2 in [0, 6]

Returns:
  • np.ndarray – Value of the gardner2014a function, value of corresponding

  • constraint function

References

[1] Gardner, Jacob R., Matt J. Kusner, Zhixiang Eddie Xu, Kilian Q.

Weinberger, and John P. Cunningham. “Bayesian Optimization with Inequality Constraints.” In ICML, pp. 937-945. 2014

queens.example_simulator_functions.gaussian_logpdf module#

Gaussian distributions.

gaussian_1d_logpdf(x)[source]#

1D Gaussian likelihood model.

Used as a basic test function for MCMC methods.

Returns:

float – The logpdf evaluated at x

gaussian_2d_logpdf(samples)[source]#

2D Gaussian logpdf.

Parameters:

samples (np.ndarray) – Samples to be evaluated

Returns:

np.ndarray – logpdf

gaussian_4d_logpdf(samples)[source]#

4D Gaussian logpdf.

Parameters:

samples (np.ndarray) – Samples to be evaluated

Returns:

np.ndarray – logpdf

queens.example_simulator_functions.gaussian_mixture_logpdf module#

Weighted mixture of 2 multivariate Gaussian distributions.

Example adapted from [1], section 3.4.

[1]: Minson, S. E., Simons, M. and Beck, J. L. (2013) ‘Bayesian inversion for finite fault earthquake source models I-theory and algorithm’, Geophysical Journal International, 194(3), pp. 1701–1726. doi: 10.1093/gji/ggt180.

gaussian_mixture_4d_logpdf(samples)[source]#

Multivariate Gaussian Mixture likelihood model.

Used as a basic test function for MCMC and SMC methods.

The log likelihood is defined as (see [2]):

\(f({x}) =\log \left( w_1 \frac{1}{\sqrt{(2 \pi)^k |\Sigma_1|}} \exp \left[-\frac{1}{2}(x-\mu_1)^T\Sigma_1^{-1}(x-\mu_1) \right]+ w_2 \frac{1}{\sqrt{(2 \pi)^k |\Sigma_2|}} \exp \left[-\frac{1}{2}(x-\mu_2)^T\Sigma_2^{-1}(x-\mu_2) \right] \right)\)

Returns:

numpy.array – The logpdf evaluated at x

References

[2] https://en.wikipedia.org/wiki/Multivariate_normal_distribution

queens.example_simulator_functions.ishigami90 module#

Ishigami function.

Ishigami function is a three-dimensional test function for sensitivity analysis and UQ.

It is nonlinear and non-monotonous.

first_effect_variance(p1=7, p2=0.1)[source]#

Total variance of the Ishigami test function.

According to (50)-(53) in [1].

[1] Homma, T., & Saltelli, A. (1996). Importance measures in global

sensitivity analysis of nonlinear models. Reliability Engineering & System Safety, 52(1), 1–17. https://doi.org/10.1016/0951-8320(96)00002-6

Parameters:
  • p1 (float) – Coefficient (optional), with default value 7

  • p2 (float) – Coefficient (optional), with default value 0.1

Returns:

float – Value of first effect (conditional) variance of ishigami function

first_order_indices(p1=7, p2=0.1)[source]#

First order Sobol’ indices of Ishigami test function.

According to (50)-(53) in [1].

[1] Homma, T., & Saltelli, A. (1996). Importance measures in global

sensitivity analysis of nonlinear models. Reliability Engineering & System Safety, 52(1), 1–17. https://doi.org/10.1016/0951-8320(96)00002-6

Parameters:
  • p1 (float) – Coefficient (optional), with default value 7

  • p2 (float) – Coefficient (optional), with default value 0.1

Returns:

float – Analytical values of first order Sobol indices of the ishigami function

ishigami90(x1, x2, x3, p1=7, p2=0.1)[source]#

Three-dimensional benchmark function.

Three-dimensional benchmark function from [2], used for UQ because it exhibits strong nonlinearity and non-monotonicity. It also has a peculiar dependence on \(x_3\), as described [5]). The values of \(a\) and \(b\) used in [1] and [3] are: \(p_1 = a = 7\) and \(p_2 = b = 0.1\). In [5] the values \(a = 7\) and \(b = 0.05\) are used.

The function is defined as:

\(f({\bf x}) = \sin(x_1) + p_1 \sin^2(x_2 +p_2x_3^4\sin(x_1))\)

Typically, distributions of the input random variables are: \(x_i \sim\) Uniform[\(-\pi, \pi\)], for all i = 1, 2, 3.

Parameters:
  • x1 (float) – Input parameter 1

  • x2 (float) – Input parameter 2

  • x3 (float) – Input parameter 3

  • p1 (float) – Coefficient (optional), with default value 7

  • p2 (float) – Coefficient (optional), with default value 0.1

Returns:

float – Value of the ishigami function at [x1, x2, x3]

References

[1] Crestaux, T., Martinez, J.-M., Le Maitre, O., & Lafitte, O. (2007).

Polynomial chaos expansion for uncertainties quantification and sensitivity analysis [PowerPoint slides]. Retrieved from SAMO 2007 website: http://samo2007.chem.elte.hu/lectures/Crestaux.pdf.

[2] Ishigami, T., & Homma, T. (1990, December). An importance

quantification technique in uncertainty analysis for computer models. In Uncertainty Modeling and Analysis, 1990. Proceedings.,

[3] Marrel, A., Iooss, B., Laurent, B., & Roustant, O. (2009).

Calculations of sobol indices for the gaussian process metamodel. Reliability Engineering & System Safety, 94(3), 742-751.

[4] Saltelli, A., Chan, K., & Scott, E. M. (Eds.). (2000).

Sensitivity analysis (Vol. 134). New York: Wiley.

[5] Sobol, I. M., & Levitan, Y. L. (1999). On the use of variance

reducing multipliers in Monte Carlo computations of a global sensitivity index. Computer Physics Communications, 117(1), 52-61.

total_order_indices(p1=7, p2=0.1)[source]#

Total order Sobol’ indices of Ishigami test function.

According to (50)-(57) in [1].

[1] Homma, T., & Saltelli, A. (1996). Importance measures in global

sensitivity analysis of nonlinear models. Reliability Engineering & System Safety, 52(1), 1–17. https://doi.org/10.1016/0951-8320(96)00002-6

Parameters:
  • p1 (float) – Coefficient (optional), with default value 7

  • p2 (float) – Coefficient (optional), with default value 0.1

Returns:

float – Analytical values of total order Sobol indices of the ishigami function

variance(p1=7, p2=0.1)[source]#

Variance of Ishigami test function.

According to (50) in [1].

[1] Homma, T., & Saltelli, A. (1996). Importance measures in global

sensitivity analysis of nonlinear models. Reliability Engineering & System Safety, 52(1), 1–17. https://doi.org/10.1016/0951-8320(96)00002-6

Parameters:
  • p1 (float) – Coefficient (optional), with default value 7

  • p2 (float) – Coefficient (optional), with default value 0.1

Returns:

float – Value of variance of the ishigami function

queens.example_simulator_functions.ma09 module#

Ma test function.

ma09(x1, x2)[source]#

Ma09 function: Two dimensional benchmark for UQ defined in [1].

\(f({\bf x}) = \frac{1}{|0.3-x_1^2 - x_2^2|+0.1}\)

Parameters:
  • x1 (float) – Input parameter 1 in [0, 1]

  • x2 (float) – Input parameter 2 in [0, 1]

Returns:

float – Value of the ma09 function

References

[1] Ma, X., & Zabaras, N. (2009). An adaptive hierarchical sparse grid

collocation algorithm for the solution of stochastic differential equations. Journal of Computational Physics, 228(8), 3084?3113.

queens.example_simulator_functions.oakley_ohagan04 module#

Oakley O’Hagan 15D function.

oakley_ohagan04(x1, x2, x3, x4, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15)[source]#

Oakley O’Hagan 2004 function, 15 dimensional benchmark function in [1].

\(f({\bf x})= {\bf a_1}^T{\bf x}+{\bf a_2}^T\sin({\bf x})+{\bf a_3}^T\cos({\bf x})+{\bf x^T M x}\)

The a-coefficients are chosen so that 5 of the input variables contribute significantly to the output variance, 5 have a much smaller effect, and the remaining 5 have almost no effect on the output variance. Values of the coefficient vectors a1, a2 and a3, and the matrix M, can be found at: http://www.jeremy-oakley.staff.shef.ac.uk/psa_example.txt.

In [1] \(x_i \sim N(\mu=0, \sigma=1), \textrm{for all i = 1, ..., 15}.\)

Parameters:
  • x1 (float) – Input parameter 1

  • x2 (float) – Input parameter 2

  • x3 (float) – Input parameter 3

  • x4 (float) – Input parameter 4

  • x5 (float) – Input parameter 5

  • x6 (float) – Input parameter 6

  • x7 (float) – Input parameter 7

  • x8 (float) – Input parameter 8

  • x9 (float) – Input parameter 9

  • x10 (float) – Input parameter 10

  • x11 (float) – Input parameter 11

  • x12 (float) – Input parameter 12

  • x13 (float) – Input parameter 13

  • x14 (float) – Input parameter 14

  • x15 (float) – Input parameter 15

Returns:

float – Value of the function at the parameters

References

[1] Oakley, J. E., & O’Hagan, A. (2004). Probabilistic sensitivity analysis

of complex models: a Bayesian approach. Journal of the Royal Statistical Society: Series B (Statistical Methodology), 66(3), 751-769.

queens.example_simulator_functions.parabola_residual module#

Residual of a parabola.

parabola_residual(x1)[source]#

Residual formulation of a parabola.

Parameters:

x1 (float) – Input parameter 1

Returns:

ndarray – Vector of residuals of the parabola

queens.example_simulator_functions.paraboloid module#

2D paraboloid.

paraboloid(x1, x2)[source]#

A paraboloid.

See https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html

Parameters:
  • x1 (float) – Input parameter 1

  • x2 (float) – Input parameter 2

Returns:

float – Value of the function

queens.example_simulator_functions.park91a module#

Park91a function.

park91a_hifi(x1, x2, x3, x4, gradient_bool=False)[source]#

High-fidelity Park91a function.

Simple four-dimensional benchmark function as proposed in [1] to mimic a computer model. For the purpose of multi-fidelity simulation, [3] defined a corresponding lower fidelity function, which is implemented in park91a_lofi.

The high-fidelity version is defined as:

\[f({\bf x}) = \frac{x_1}{2}\left[ \sqrt{1+(x_2+x_3^2)\frac{x_4}{x_1^2}}-1 \right] +(x_1+3x_4)\exp \left[1-\sin(x_3) \right]\]
Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

  • x3 (float) – Input parameter 3 [0,1)

  • x4 (float) – Input parameter 4 [0,1)

  • gradient_bool (bool, optional) – Switch to return gradient value w.r.t. x1 and x2

Returns:
  • y (float) – Value of function at parameters

  • dy_dx1 (float) – Gradient of the function w.r.t. x1

  • dy_dx2 (float) – Gradient of the function w.r.t. x2

References

[1] Park, J.-S.(1991). Tuning complex computer codes to data and optimal

designs, Ph.D. Thesis

[2] Cox, D. D., Park, J.-S., & Singer, C. E. (2001). A statistical method

for tuning a computer code to a database. Computational Statistics & Data Analysis, 37(1), 77?92. http://doi.org/10.1016/S0167-9473(00)00057-8

[3] Xiong, S., Qian, P., & Wu, C. (2013). Sequential design and analysis of

high-accuracy and low-accuracy computer codes. Technometrics. http://doi.org/10.1080/00401706.2012.723572

park91a_hifi_on_grid(x1, x2)[source]#

High-fidelity Park91a function on x3 and x4 grid.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

Returns:

np.ndarray – Value of the function at the parameters

park91a_hifi_on_grid_with_gradients(x1, x2)[source]#

High-fidelity Park91a function on x3 and x4 grid with gradients.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

Returns:
  • y (float) – Value of function at parameters

  • dy_dx1 (float) – Gradient of the function w.r.t. x1

  • dy_dx2 (float) – Gradient of the function w.r.t. x2

park91a_lofi(x1, x2, x3, x4, gradient_bool=False)[source]#

Low-fidelity Park91a function.

Simple four-dimensional benchmark function as proposed in [1] to mimic a computer model. For the purpose of multi-fidelity simulation, [3] defined a corresponding lower fidelity function, which is defined as:

\(f_{lofi}({\bf x})= \left[1 + \frac{\sin(x_1)}{10} \right] f_{hifi}({\bf x}) - 2x_1 + x_2^2 + x_3^2 + 0.5\)

The high-fidelity version is defined as is implemented in park91a_hifi.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

  • x3 (float) – Input parameter 3 [0,1)

  • x4 (float) – Input parameter 4 [0,1)

  • gradient_bool (bool, optional) – Switch to return gradient value w.r.t. x1 and x2

Returns:
  • y (float) – Value of function at parameters

  • dy_dx1 (float) – Gradient of the function w.r.t. x1

  • dy_dx2 (float) – Gradient of the function w.r.t. x2

References

[1] Park, J.-S.(1991). Tuning complex computer codes to data and optimal

designs, Ph.D Thesis

[2] Cox, D. D., Park, J.-S., & Singer, C. E. (2001). A statistical method

for tuning a computer code to a data base. Computational Statistics & Data Analysis, 37(1), 77-92. http://doi.org/10.1016/S0167-9473(00)00057-8

[3] Xiong, S., Qian, P., & Wu, C. (2013). Sequential design and analysis of

high-accuracy and low-accuracy computer codes. Technometrics. http://doi.org/10.1080/00401706.2012.723572

park91a_lofi_on_grid(x1, x2)[source]#

Low-fidelity Park91a function on fixed x3 and x4 grid.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

Returns:

np.ndarray – Value of the function at the parameters

park91a_lofi_on_grid_with_gradients(x1, x2)[source]#

Low-fidelity Park91a function on x3 and x4 grid with gradients.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

Returns:
  • y (float) – Value of function at parameters

  • dy_dx1 (float) – Gradient of the function w.r.t. x1

  • dy_dx2 (float) – Gradient of the function w.r.t. x2

unit_bounding(*args)[source]#

Bounding function.

This should be avoided…

x3_x4_grid_eval(park_function, x1, x2, gradient_bool=False)[source]#

Evaluate a park function on a x3 and x4 grid.

Parameters:
  • park_function (func) – Function to be evaluated

  • x1 (int) – Input parameter 1

  • x2 (int) – Input parameter 2

  • gradient_bool (bool, optional) – Switch to return gradient value w.r.t. x1 and x2

Returns:

np.ndarraypark_function evaluated on the grid

queens.example_simulator_functions.park91b module#

Park91b functions.

park91b_hifi(x1, x2, x3, x4)[source]#

High-fidelity version of Park91b benchmark function.

Simple four dimensional benchmark function as proposed in [1] to mimic a computer model. The high-fidelity version is defined as:

\(f_{hifi}({\bf x})= \frac{2}{3} \exp(x_1 + x_2) - x_4 \sin(x_3) + x_3\)

For the purpose of multi-fidelity simulation, [2] defined a corresponding lower fidelity function, which is implemented in park91b_lofi.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

  • x3 (float) – Input parameter 3 [0,1)

  • x4 (float) – Input parameter 4 [0,1)

Returns:

float – Value of function at parameters

References

[1] Park, J.-S.(1991). Tuning complex computer codes to data and optimal

designs, Ph.D Thesis

[2] Xiong, S., Qian, P., & Wu, C. (2013). Sequential design and analysis

of high-accuracy and low-accuracy computer codes. Technometrics. http://doi.org/10.1080/00401706.2012.723572

park91b_lofi(x1, x2, x3, x4)[source]#

Low-fidelity version of the Park91b benchmark function.

Simple four-dimensional benchmark function as proposed in [1], to mimic a computer model. The low-fidelity version is defined as:

\(f_{lofi}({\bf x})=1.2 f_{hifi}({\bf x})-1\)

The corresponding high-fidelity function is implemented in park91b_hifi.

Parameters:
  • x1 (float) – Input parameter 1 [0,1)

  • x2 (float) – Input parameter 2 [0,1)

  • x3 (float) – Input parameter 3 [0,1)

  • x4 (float) – Input parameter 4 [0,1)

Returns:

float – Value of the function at the parameters

References

[1] Park, J.-S.(1991). Tuning complex computer codes to data and optimal

designs, Ph.D Thesis

[2] Xiong, S., Qian, P., & Wu, C. (2013). Sequential design and analysis

of high-accuracy and low-accuracy computer codes. Technometrics. http://doi.org/10.1080/00401706.2012.723572

queens.example_simulator_functions.perdikaris17 module#

Perdikaris test functions.

perdikaris17_hifi(x)[source]#

High-fidelity version of simple 1D test function.

High-fidelity version of simple 1-dimensional benchmark function as proposed in [1] and defined as:

\(f_{hifi}(x)= (x-\sqrt{2})(f_{lofi}(x))^2\)

The low-fidelity version of the function was also proposed in [1] and is in implemented in perdikaris_1dsin_lofi.

Parameters:

x (float) – Input parameter [0,1]

Returns:

float – Value of function at x

References

[1] Perdikaris, P. et al., 2017. Nonlinear information fusion algorithms

for data-efficient multi-fidelity modelling. Proceedings of the Royal Society of London A: Mathematical, Physical and Engineering Sciences, 473(2198), pp.20160751?16.

perdikaris17_lofi(x)[source]#

Low-fidelity version of simple 1D test function.

Low-fidelity version of a simple 1-dimensional benchmark function as proposed in [1] and defined as:

\(f_{lofi}({\bf x}) = \sin(8.0\pi x)\)

The high-fidelity version of the function was also proposed in [1] and is in implemented in perdikaris_1dsin_hifi.

Parameters:

x (float) – Input parameter

Returns:

float – Value of function at x

References

[1] Perdikaris, P. et al., 2017. Nonlinear information fusion algorithms

for data-efficient multi-fidelity modelling. Proceedings of the Royal Society of London A: Mathematical, Physical and Engineering Sciences, 473(2198), pp.20160751?16.

queens.example_simulator_functions.rezende15 module#

Probabilistic models in [1].

[1]: Rezende, D. J., & Mohamed, S. (2016). Variational Inference with Normalizing Flows. ArXiv:1505. 05770 [Cs, Stat]. http://arxiv.org/abs/1505.05770

rezende15_potential1(x, theta=None, as_logpdf=False)[source]#

First potential in [1].

The unnormalized probabilistic model used is proportional to \(p(\theta)\propto \exp(-U(\theta))\) where \(U(\theta)\) is a potential. Hence the log_posterior_unnormalized is given by \(-U(\theta)\).

Parameters:
  • x (np.ndarray) – Samples at which to evaluate the potential (2 \(\times\) n_samples)

  • theta (float) – Angle in radiants in which to rotate the potential

  • as_logpdf (bool,optional) – True if \(-U\) is to be returned

Returns:

np.ndarray – Potential or unnormalized logpdf

queens.example_simulator_functions.rosenbrock60 module#

Rosenbrock function.

[1]: Rosenbrock, H. H. (1960). An Automatic Method for Finding the Greatest or Least Value of a Function. The Computer Journal, 3(3), 175–184. doi:10.1093/comjnl/3.3.175

rosenbrock60(x1, x2)[source]#

Rosenbrocks banana function.

Parameters:
  • x1 (float) – Input parameter 1

  • x2 (float) – Input parameter 2

Returns:

float – Value of the Rosenbrock function

rosenbrock60_residual(x1, x2)[source]#

Residuals of the Rosenbrock banana function.

Parameters:
  • x1 (float) – Input parameter 1

  • x2 (float) – Input parameter 2

Returns:

ndarray – Vector of residuals of the Rosenbrock function

rosenbrock60_residual_1d(x1)[source]#

Residuals of the Rosenbrock banana function.

Parameters:

x1 (float) – Input parameter 1

Returns:

ndarray – Vector of residuals of the Rosenbrock function

rosenbrock60_residual_3d(x1, x2, x3)[source]#

Residuals of the Rosenbrock banana function.

Parameters:
  • x1 (float) – Input parameter 1

  • x2 (float) – Input parameter 2

  • x3 (float) – Input parameter 3

Returns:

ndarray – Vector of residuals of the Rosenbrock function

queens.example_simulator_functions.sinus module#

Sine test function.

gradient_sinus_test_fun(x1)[source]#

Evaluate sine and its gradient.

Parameters:

x1 (float) – Input of the sinus in RAD

Returns:

result (float) – Value of the sinus function

sinus_test_fun(x1)[source]#

Evaluate standard sine as a test function.

Parameters:

x1 (float) – Input of the sinus in RAD

Returns:

result (float) – Value of the sinus function

queens.example_simulator_functions.sobol_g_function module#

Sobol’s G function.

A test function for sensitivity analysis.

References

[1]: Saltelli, A., Annoni, P., Azzini, I., Campolongo, F., Ratto, M., & Tarantola, S. (2010). Variance based sensitivity analysis of model output. Design and estimator for the total sensitivity index. Computer Physics Communications, 181(2), 259–270. https://doi.org/10.1016/j.cpc.2009.09.018

first_effect_variance(a=array([0., 0.1, 0.2, 0.3, 0.4, 0.8, 1., 2., 3., 4.]), alpha=array([2, 2, 2, 2, 2, 2, 2, 2, 2, 2]))[source]#

Calculate first effect variance \(V_{x_{i}}[E_{x \sim i}[Y|x_i]]\).

Parameters:
  • a (ndarray) – Vector of a parameter values

  • alpha (ndarray) – Vector of alpha parameter values

Returns:

ndarray – Vector of first order variances

first_order_indices(a=array([0., 0.1, 0.2, 0.3, 0.4, 0.8, 1., 2., 3., 4.]), alpha=array([2, 2, 2, 2, 2, 2, 2, 2, 2, 2]))[source]#

Compute first order indices of the Sobol test function.

See (32) in [1].

Parameters:
  • a (ndarray) – Vector of a parameter values

  • alpha (ndarray) – Vector of alpha parameter values

Returns:

ndarray – Vector of first order indices

sobol_g_function(a=array([0., 0.1, 0.2, 0.3, 0.4, 0.8, 1., 2., 3., 4.]), alpha=array([2, 2, 2, 2, 2, 2, 2, 2, 2, 2]), delta=array([0, 0, 0, 0, 0, 0, 0, 0, 0, 0]), **kwargs)[source]#

Compute generalized Sobol’s G function.

With variable dimension. See (33) in [1]. Default is a 10-dimensional version as defined in [1].

Parameters:
  • a (ndarray) – Vector of a parameter values

  • alpha (ndarray) – Vector of alpha parameter values

  • delta (ndarray) – Vector of delta parameter values

  • kwargs (dict) – Contains the input x_i ~ U(0,1) (uniformly from [0,1])

Returns:

double – Value of the Sobol function

total_order_indices(a=array([0., 0.1, 0.2, 0.3, 0.4, 0.8, 1., 2., 3., 4.]), alpha=array([2, 2, 2, 2, 2, 2, 2, 2, 2, 2]))[source]#

Compute total indices of Sobol test function.

See (31)-(32) in [1].

Parameters:
  • a (ndarray) – Vector of a parameter values

  • alpha (ndarray) – Vector of alpha parameter values

Returns:

ndarray – Vector of total order indices

variance(variance_i=array([0.8, 0.66115702, 0.55555556, 0.47337278, 0.40816327, 0.24691358, 0.2, 0.08888889, 0.05, 0.032]))[source]#

Calculate variance of Sobol function \(V_x[Y]\).

Parameters:

variance_i (ndarray) – Vector of first effect variances

Returns:

double – Variance of Sobol’s G function