surrogate.kriging

surrogate.kriging

Kriging (Gaussian Process) surrogate model for SpotOptim.

Adapted from spotpython.surrogate.kriging_basic for compatibility with SpotOptim. This implementation follows Forrester et al. (2008) “Engineering Design via Surrogate Modelling”.

Specific references: - Section 2.4 “Kriging”: Core implementation of the Kriging predictor and likelihood. - Section 6 “Surrogate Modeling of Noisy Data”: Implementation of “regression” and “reinterpolation” methods. - Validated against the book’s Matlab implementation: - likelihood.m: Concentrated log-likelihood calculation. - pred.m: Prediction and error estimation.

Classes

Name Description
Kriging A scikit-learn compatible Kriging model class for regression tasks.

Kriging

surrogate.kriging.Kriging(
    noise=None,
    penalty=10000.0,
    method='regression',
    var_type=None,
    name='Kriging',
    seed=124,
    model_fun_evals=None,
    n_theta=None,
    min_theta=-3.0,
    max_theta=2.0,
    theta_init_zero=False,
    p_val=2.0,
    n_p=1,
    optim_p=False,
    min_Lambda=-9.0,
    max_Lambda=0.0,
    metric_factorial='canberra',
    isotropic=False,
    theta=None,
    **kwargs,
)

A scikit-learn compatible Kriging model class for regression tasks.

This class provides Ordinary Kriging with support for: - Mixed variable types (continuous, integer, factor) - Gaussian/RBF correlation function - Three fitting methods (Forrester (2008), Section 6): - Isotropic or anisotropic length scales

Compatible with SpotOptim’s variable type conventions: - ‘float’: continuous numeric variables - ‘int’: integer variables - ‘factor’: categorical/unordered variables

Parameters

Name Type Description Default
noise float Small regularization term for numerical stability (nugget effect). If None, defaults to sqrt(machine epsilon). Only used for “interpolation” method. For “regression” and “reinterpolation”, this is replaced by the Lambda parameter. Defaults to None. None
penalty float Large negative log-likelihood assigned if correlation matrix is not positive-definite. Defaults to 1e4. 10000.0
method str Fitting method (Forrester (2008), Section 6). Options: - “interpolation”: Pure Kriging interpolation (Eq 2.X). Fits exact data points. Uses a small noise (nugget) for numerical stability. - “regression”: Regression Kriging (Section 6.2). Optimizes a regularization parameter Lambda (nugget) along with theta. Suitable for noisy data. - “reinterpolation”: Re-interpolation (Section 6.3). Fits hyperparameters using regression (with Lambda), but predicts using the “noise-free” correlation matrix (removing Lambda). This creates a surrogate that glosses over noise but passes closer to the underlying trend (interpolating the regression model). Defaults to “regression”. 'regression'
var_type List[str] Variable types for each dimension. SpotOptim uses: ‘float’ (continuous), ‘int’ (integer), ‘factor’ (categorical). Defaults to [“float”]. None
name str Name of the Kriging instance. Defaults to “Kriging”. 'Kriging'
seed int Random seed for reproducibility. Defaults to 124. 124
model_fun_evals int Maximum function evaluations for hyperparameter optimization. Defaults to 100. None
n_theta int Number of theta parameters. If None, set during fit. Defaults to None. None
min_theta float Minimum log10(theta) bound. Defaults to -3.0. -3.0
max_theta float Maximum log10(theta) bound. Defaults to 2.0. 2.0
theta_init_zero bool Initialize theta to zero. Defaults to False. False
p_val float Power parameter for correlation (fixed at 2.0 for Gaussian). Defaults to 2.0. 2.0
n_p int Number of p parameters (currently not optimized). Defaults to 1. 1
optim_p bool Optimize p parameters (currently not supported). Defaults to False. False
min_Lambda float Minimum log10(Lambda) bound. Defaults to -9.0. -9.0
max_Lambda float Maximum log10(Lambda) bound. Defaults to 0.0. 0.0
metric_factorial str Distance metric for factor variables. Defaults to “canberra”. 'canberra'
isotropic bool Use single theta for all dimensions. Defaults to False. False
theta np.ndarray Initial theta values (log10 scale). Defaults to None. None

Attributes

Name Type Description
X_ ndarray Training data, shape (n_samples, n_features).
y_ ndarray Training targets, shape (n_samples,).
theta_ ndarray Optimized log10(theta) parameters.
Lambda_ float or None Optimized log10(Lambda) for regression methods.
mu_ float Mean of Kriging predictor.
sigma2_ float Variance of Kriging predictor.
U_ ndarray Cholesky factor of correlation matrix.
Psi_ ndarray Correlation matrix.
negLnLike float Negative log-likelihood value.

Examples

Basic usage with SpotOptim:

>>> import numpy as np
>>> from spotoptim import SpotOptim
>>> from spotoptim.surrogate import Kriging
>>>
>>> # Define objective
>>> def objective(X):
...     return np.sum(X**2, axis=1)
>>>
>>> # Create Kriging surrogate
>>> kriging = Kriging(
...     noise=1e-10,
...     method='regression',
...     min_theta=-3.0,
...     max_theta=2.0,
...     seed=42
... )
>>>
>>> # Use with SpotOptim
>>> opt = SpotOptim(
...     fun=objective,
...     bounds=[(-5, 5), (-5, 5)],
...     surrogate=kriging,
...     max_iter=30,
...     n_initial=10,
...     seed=42
... )
>>> result = opt.optimize()

Direct usage (scikit-learn compatible):

>>> from spotoptim.surrogate import Kriging
>>> import numpy as np
>>>
>>> # Training data
>>> X_train = np.array([[0.0, 0.0], [0.5, 0.5], [1.0, 1.0]])
>>> y_train = np.array([0.1, 0.2, 0.3])
>>>
>>> # Fit model
>>> model = Kriging(method='regression', seed=42)
>>> model.fit(X_train, y_train)
>>>
>>> # Predict
>>> X_test = np.array([[0.25, 0.25], [0.75, 0.75]])
>>> y_pred = model.predict(X_test)
>>>
>>> # Predict with uncertainty
>>> y_pred, std = model.predict(X_test, return_std=True)

Methods

Name Description
build_correlation_matrix Build correlation matrix from training data.
build_psi_vector Build correlation vector between x and training points.
fit Fit the Kriging model to training data.
get_params Get parameters for this estimator (scikit-learn compatibility).
likelihood Compute negative concentrated log-likelihood.
objective Objective function for hyperparameter optimization.
predict Predict at new points.
predict_single Predict at a single point.
set_params Set parameters (scikit-learn compatibility).
build_correlation_matrix
surrogate.kriging.Kriging.build_correlation_matrix()

Build correlation matrix from training data.

Returns
Name Type Description
np.ndarray np.ndarray: Upper triangle of correlation matrix.
Examples
>>> import numpy as np
>>> from spotoptim.surrogate import Kriging
>>> X = np.array([[0.], [1.]])
>>> y = np.array([0., 1.])
>>> k = Kriging(seed=42).fit(X, y)
>>> Psi_upper = k.build_correlation_matrix()
>>> print(Psi_upper.shape)
(2, 2)
build_psi_vector
surrogate.kriging.Kriging.build_psi_vector(x)

Build correlation vector between x and training points.

Parameters
Name Type Description Default
x np.ndarray Single test point, shape (n_features,). required
Returns
Name Type Description
np.ndarray np.ndarray: Correlation vector, shape (n_samples,).
Reference

Calculates the vector small psi from Eq 2.15 in Forrester (2008).

Examples
>>> import numpy as np
>>> from spotoptim.surrogate import Kriging
>>> X = np.array([[0.], [1.]])
>>> y = np.array([0., 1.])
>>> k = Kriging(seed=42).fit(X, y)
>>> x_new = np.array([0.5])
>>> psi = k.build_psi_vector(x_new)
>>> print(psi.shape)
(2,)
fit
surrogate.kriging.Kriging.fit(X, y)

Fit the Kriging model to training data.

Optimizes hyperparameters (theta, Lambda) by maximizing the concentrated log-likelihood using differential evolution.

Parameters
Name Type Description Default
X np.ndarray Training inputs, shape (n_samples, n_features). required
y np.ndarray Training targets, shape (n_samples,). required
Returns
Name Type Description
Kriging Kriging Fitted model instance (self).
get_params
surrogate.kriging.Kriging.get_params(deep=True)

Get parameters for this estimator (scikit-learn compatibility).

Parameters
Name Type Description Default
deep bool Ignored, for compatibility. True
Returns
Name Type Description
dict Dict Parameter names mapped to their values.
likelihood
surrogate.kriging.Kriging.likelihood(params)

Compute negative concentrated log-likelihood.

Parameters
Name Type Description Default
params np.ndarray Hyperparameters [log10(theta), log10(Lambda)]. required
Returns
Name Type Description
tuple Tuple[float, np.ndarray, np.ndarray] (negLnLike, Psi, U) where U is Cholesky factor.
Reference

Forrester et al. (2008), Section 2.4. Matches implementation in likelihood.m from the book’s code. Concentrated Log-Likelihood approx: ln(L) ~ -(n/2)ln(sigma^2) - (1/2)ln|R|

Examples
>>> import numpy as np
>>> from spotoptim.surrogate import Kriging
>>> X = np.array([[0.], [1.]])
>>> y = np.array([0., 1.])
>>> k = Kriging(seed=42).fit(X, y)
>>> params = np.concatenate([k.theta_, [k.Lambda_]])
>>> nll, _, _ = k.likelihood(params)
objective
surrogate.kriging.Kriging.objective(params)

Objective function for hyperparameter optimization.

Parameters
Name Type Description Default
params np.ndarray Hyperparameters to evaluate. required
Returns
Name Type Description
float float Negative concentrated log-likelihood.
Examples
>>> import numpy as np
>>> from spotoptim.surrogate import Kriging
>>> X = np.array([[0.], [1.]])
>>> y = np.array([0., 1.])
>>> k = Kriging(seed=42).fit(X, y)
>>> # Evaluate objective at optimal parameters
>>> val = k.objective(np.concatenate([k.theta_, [k.Lambda_]]))
predict
surrogate.kriging.Kriging.predict(X, return_std=False)

Predict at new points.

Parameters
Name Type Description Default
X np.ndarray Test points, shape (n_samples, n_features). required
return_std bool Return standard deviations. Defaults to False. False
Returns
Name Type Description
np.ndarray np.ndarray: Predictions, shape (n_samples,).
tuple np.ndarray (predictions, std_devs) if return_std=True.
predict_single
surrogate.kriging.Kriging.predict_single(x)

Predict at a single point.

Parameters
Name Type Description Default
x np.ndarray Test point, shape (n_features,). required
Returns
Name Type Description
tuple Tuple[float, float] (prediction, std_dev).
Reference

Forrester et al. (2008), Eq 2.15 (Predictor) and Eq 2.19 (Error). Matches implementation in pred.m from the book’s code.

Examples
>>> import numpy as np
>>> from spotoptim.surrogate import Kriging
>>> X = np.array([[0.], [1.]])
>>> y = np.array([0., 1.])
>>> k = Kriging(seed=42).fit(X, y)
>>> x_new = np.array([0.5])
>>> y_pred, y_std = k.predict_single(x_new)
set_params
surrogate.kriging.Kriging.set_params(**params)

Set parameters (scikit-learn compatibility).

Parameters
Name Type Description Default
**params typing.Any Parameter names and values. {}
Returns
Name Type Description
Kriging Kriging Self.