multifi_cokriging.py#

Implements the Multi-Fidelity Co-Kriging method described in [LeGratiet2013].

(Author: Remi Vauclin vauclin.remi@gmail.com)

This code was implemented using the package scikit-learn as basis. (Author: Vincent Dubourg, vincent.dubourg@gmail.com)

OpenMDAO adaptation. Regression and correlation functions were directly copied from scikit-learn package here to avoid scikit-learn dependency. (Author: Remi Lafage, remi.lafage@onera.fr)

ISAE/DMSM - ONERA/DCPS

class openmdao.surrogate_models.multifi_cokriging.MultiFiCoKriging(regr='constant', rho_regr='constant', normalize=True, theta=None, theta0=None, thetaL=None, thetaU=None)[source]

Bases: object

Implement the Multi-Fidelity Co-Kriging method described in [LeGratiet2013].

Parameters:
regrstr or callable, optional

A regression function returning an array of outputs of the linear regression functional basis for Universal Kriging purpose. regr is assumed to be the same for all levels of code. Default assumes a simple constant regression trend. Available built-in regression models are: ‘constant’, ‘linear’.

rho_regrstr or callable, optional

A regression function returning an array of outputs of the linear regression functional basis. Defines the regression function for the autoregressive parameter rho. rho_regr is assumed to be the same for all levels of code. Default assumes a simple constant regression trend. Available built-in regression models are: ‘constant’, ‘linear’.

normalizebool, optional

When true, normalize X and Y so that the mean is at zero.

thetadouble, array_like or list, optional

Value of correlation parameters if they are known; no optimization is run. Default is None, so that optimization is run. if double: value is replicated for all features and all levels. if array_like: an array with shape (n_features, ) for isotropic calculation. It is replicated for all levels. if list: a list of nlevel arrays specifying value for each level.

theta0double, array_like or list, optional

Starting point for the maximum likelihood estimation of the best set of parameters. Default is None and meaning use of the default 0.5*np.ones(n_features) if double: value is replicated for all features and all levels. if array_like: an array with shape (n_features, ) for isotropic calculation. It is replicated for all levels. if list: a list of nlevel arrays specifying value for each level.

thetaLdouble, array_like or list, optional

Lower bound on the autocorrelation parameters for maximum likelihood estimation. Default is None meaning use of the default 1e-5*np.ones(n_features). if double: value is replicated for all features and all levels. if array_like: An array with shape matching theta0’s. It is replicated for all levels of code. if list: a list of nlevel arrays specifying value for each level.

thetaUdouble, array_like or list, optional

Upper bound on the autocorrelation parameters for maximum likelihood estimation. Default is None meaning use of default value 50*np.ones(n_features). if double: value is replicated for all features and all levels. if array_like: An array with shape matching theta0’s. It is replicated for all levels of code. if list: a list of nlevel arrays specifying value for each level.

Notes

Implementation is based on the Package Scikit-Learn (Author: Vincent Dubourg, vincent.dubourg@gmail.com) which translates the DACE Matlab toolbox, see [Rafec0a633dc4-NLNS2002].

References

[NLNS2002]

H. B. Nielsen, S. N. Lophaven, and J. Sondergaard. DACE - A MATLAB Kriging Toolbox. (2002) http://www2.imm.dtu.dk/~hbn/dace/dace.pdf

[WBSWM1992]

W. J. Welch, R. J. Buck, J. Sacks, H. P. Wynn, T. J. Mitchell, and M. D. Morris (1992). “Screening, predicting, and computer experiments.” Technometrics, 34(1) 15–25. http://www.jstor.org/pss/1269548

[LeGratiet2013]

L. Le Gratiet (2013). “Multi-fidelity Gaussian process regression for computer experiments.” PhD thesis, Universite Paris-Diderot-Paris VII.

[TBKH2011]

Toal, D. J., Bressloff, N. W., Keane, A. J., & Holden, C. M. E. (2011). “The development of a hybridized particle swarm for kriging hyperparameter tuning.” Engineering optimization, 43(6), 675-699.

Attributes:
corrobject

Correlation function to use, default is squared_exponential_correlation.

n_featuresndarray

Number of features for each fidelity level.

n_samplesndarray

Number of samples for each fidelity level.

nlevelint

Number of fidelity levels.

normalizebool, optional

When true, normalize X and Y so that the mean is at zero.

regrstr or callable

A regression function returning an array of outputs of the linear regression functional basis for Universal Kriging purpose. regr is assumed to be the same for all levels of code. Default assumes a simple constant regression trend. Available built-in regression models are: ‘constant’, ‘linear’

rho_regrstr or callable or None

A regression function returning an array of outputs of the linear regression functional basis. Defines the regression function for the autoregressive parameter rho. rho_regr is assumed to be the same for all levels of code. Default assumes a simple constant regression trend. Available built-in regression models are: ‘constant’, ‘linear’

thetadouble, array_like or list or None

Value of correlation parameters if they are known; no optimization is run. Default is None, so that optimization is run. if double: value is replicated for all features and all levels. if array_like: an array with shape (n_features, ) for isotropic calculation. It is replicated for all levels. if list: a list of nlevel arrays specifying value for each level

theta0double, array_like or list or None

Starting point for the maximum likelihood estimation of the best set of parameters. Default is None and meaning use of the default 0.5*np.ones(n_features) if double: value is replicated for all features and all levels. if array_like: an array with shape (n_features, ) for isotropic calculation. It is replicated for all levels. if list: a list of nlevel arrays specifying value for each level

thetaLdouble, array_like or list or None

Lower bound on the autocorrelation parameters for maximum likelihood estimation. Default is None meaning use of the default 1e-5*np.ones(n_features). if double: value is replicated for all features and all levels. if array_like: An array with shape matching theta0’s. It is replicated for all levels of code. if list: a list of nlevel arrays specifying value for each level

thetaUdouble, array_like or list or None

Upper bound on the autocorrelation parameters for maximum likelihood estimation. Default is None meaning use of default value 50*np.ones(n_features). if double: value is replicated for all features and all levels. if array_like: An array with shape matching theta0’s. It is replicated for all levels of code. if list: a list of nlevel arrays specifying value for each level

X_meanfloat

Mean of the low fidelity training data for X.

X_stdfloat

Standard deviation of the low fidelity training data for X.

y_meanfloat

Mean of the low fidelity training data for y.

y_stdfloat

Standard deviation of the low fidelity training data for y.

_nfevint

Number of function evaluations.

__init__(regr='constant', rho_regr='constant', normalize=True, theta=None, theta0=None, thetaL=None, thetaU=None)[source]

Initialize all attributes.

fit(X, y, initial_range=0.3, tol=1e-06)[source]

Implement the Multi-Fidelity co-kriging model fitting method.

Parameters:
Xlist of double array_like elements

A list of arrays with the input at which observations were made, from lowest fidelity to highest fidelity. Designs must be nested with X[i] = np.vstack([…, X[i+1]).

ylist of double array_like elements

A list of arrays with the observations of the scalar output to be predicted, from lowest fidelity to highest fidelity.

initial_rangefloat

Initial range for the optimizer.

tolfloat

Optimizer terminates when the tolerance tol is reached.

predict(X, eval_MSE=True)[source]

Perform the predictions of the kriging model on X.

Parameters:
Xarray_like

An array with shape (n_eval, n_features) giving the point(s) at which the prediction(s) should be made.

eval_MSEbool, optional

A boolean specifying whether the Mean Squared Error should be evaluated or not. Default assumes evalMSE is True.

Returns:
array_like

An array with shape (n_eval, ) with the Best Linear Unbiased Prediction at X. If all_levels is set to True, an array with shape (n_eval, nlevel) giving the BLUP for all levels.

array_like, optional (if eval_MSE is True)

An array with shape (n_eval, ) with the Mean Squared Error at X. If all_levels is set to True, an array with shape (n_eval, nlevel) giving the MSE for all levels.

rlf(lvl, theta=None)[source]

Determine BLUP parameters and evaluate negative reduced likelihood function for theta.

Maximizing this function wrt the autocorrelation parameters theta is equivalent to maximizing the likelihood of the assumed joint Gaussian distribution of the observations y evaluated onto the design of experiments X.

Parameters:
lvlint

Level of fidelity.

thetaarray_like, optional

An array containing the autocorrelation parameters at which the Gaussian Process model parameters should be determined. Default uses the built-in autocorrelation parameters (ie theta = self.theta).

Returns:
double

The value of the negative concentrated reduced likelihood function associated to the given autocorrelation parameters theta.

class openmdao.surrogate_models.multifi_cokriging.MultiFiCoKrigingSurrogate(**kwargs)[source]

Bases: MultiFiSurrogateModel

OpenMDAO adapter of multi-fidelity recursive cokriging method described in [LeGratiet2013].

See MultiFiCoKriging class.

Parameters:
**kwargskeyword args

Some implementations of record_derivatives need additional args.

Attributes:
modelMultiFiCoKriging

Contains MultiFiCoKriging surrogate.

__init__(**kwargs)[source]

Initialize all attributes.

linearize(x)

Calculate the jacobian of the interpolant at the requested point.

Parameters:
xarray-like

Point at which the surrogate Jacobian is evaluated.

predict(new_x)[source]

Calculate a predicted value of the response based on the current trained model.

Parameters:
new_xarray_like

An array with shape (n_eval, n_features) giving the point(s) at which the prediction(s) should be made.

Returns:
array_like

An array with shape (n_eval, ) with the Best Linear Unbiased Prediction at X. If all_levels is set to True, an array with shape (n_eval, nlevel) giving the BLUP for all levels.

array_like

An array with shape (n_eval, ) with the square root of the Mean Squared Error at X.

train(x, y)

Train the surrogate model with the given set of inputs and outputs.

Parameters:
xarray-like

Point(s) at which the surrogate is evaluated.

yarray-like

Model responses at given inputs.

train_multifi(X, Y)[source]

Train the surrogate model with the given set of inputs and outputs.

Parameters:
Xlist of double array_like elements

A list of arrays with the input at which observations were made, from highest fidelity to lowest fidelity. Designs must be nested with X[i] = np.vstack([…, X[i+1]).

Ylist of double array_like elements

A list of arrays with the observations of the scalar output to be predicted, from highest fidelity to lowest fidelity.

vectorized_predict(x)

Calculate predicted values of the response based on the current trained model.

Parameters:
xarray-like

Vectorized point(s) at which the surrogate is evaluated.

openmdao.surrogate_models.multifi_cokriging.constant_regression(x)[source]

Zero order polynomial (constant, p = 1) regression model.

x –> f(x) = 1

Parameters:
xarray_like

Input data.

Returns:
array_like

Constant regression output.

openmdao.surrogate_models.multifi_cokriging.l1_cross_distances(X, Y=None)[source]

Compute the nonzero componentwise L1 cross-distances between the vectors in X and Y.

Parameters:
Xarray_like

An array with shape (n_samples_X, n_features).

Yarray_like

An array with shape (n_samples_Y, n_features).

Returns:
array with shape (n_samples * (n_samples - 1) / 2, n_features)

The array of componentwise L1 cross-distances.

openmdao.surrogate_models.multifi_cokriging.linear_regression(x)[source]

First order polynomial (linear, p = n+1) regression model.

x –> f(x) = [ 1, x_1, …, x_n ].T

Parameters:
xarray_like

Input data.

Returns:
array_like

Linear regression output.

openmdao.surrogate_models.multifi_cokriging.squared_exponential_correlation(theta, d)[source]

Squared exponential correlation model (Radial Basis Function).

(Infinitely differentiable stochastic process, very smooth):

                                    n
theta, dx --> r(theta, dx) = exp(  sum  - theta_i * (dx_i)^2 )
                                  i = 1
Parameters:
thetaarray_like

An array with shape 1 (isotropic) or n (anisotropic) giving the autocorrelation parameter(s).

darray_like

An array with shape (n_eval, n_features) giving the componentwise distances between locations x and x’ at which the correlation model should be evaluated.

Returns:
array_like

An array with shape (n_eval, ) containing the values of the autocorrelation model.