blackboxopt.rbf module

Radial Basis Function model.

class blackboxopt.rbf.MedianLpfFilter

Bases: RbfFilter

Filter values by replacing large function values by the median of all.

This strategy was proposed by [1] based on results from [2]. Use this strategy to reduce oscillations of the interpolator, especially if the range target function is large. This filter may reduce the quality of the approximation by the surrogate.

References

__call__(x) ndarray

Call self as a function.

Return type:

ndarray

class blackboxopt.rbf.RbfFilter

Bases: object

Base filter class for the RBF target training set. Trivial identity filter.

__call__(x) ndarray

Call self as a function.

Return type:

ndarray

class blackboxopt.rbf.RbfKernel(value)

Bases: Enum

RBF kernel tags.

CUBIC = 3

Cubic kernel, i.e., \(\phi(r) = r^3\)

LINEAR = 1

Linear kernel, i.e., \(\phi(r) = r\)

THINPLATE = 2

Thinplate spline: \(\phi(r) = r^2 \log(r)\)

class blackboxopt.rbf.RbfModel(kernel: RbfKernel = RbfKernel.CUBIC, iindex: tuple[int, ...] = (), filter: RbfFilter | None = None)

Bases: object

Radial Basis Function model.

\[f(x) = \sum_{i=1}^{m} \beta_i \phi(\|x - x_i\|) + \sum_{i=1}^{n} \beta_{m+i} p_i(x),\]

where:

  • \(m\) is the number of sampled points.

  • \(x_i\) are the sampled points.

  • \(\beta_i\) are the coefficients of the RBF model.

  • \(\phi\) is the kernel function.

  • \(p_i\) are the basis functions of the polynomial tail.

  • \(n\) is the dimension of the polynomial tail.

This implementation focuses on quick successive updates of the model, which is essential for the good performance of active learning processes.

Parameters:
  • kernel (RbfKernel) – Kernel function \(\phi\) used in the RBF model. (default: <RbfKernel.CUBIC: 3>)

  • iindex (tuple[int, ...]) – Indices of integer variables in the feature space. (default: ())

  • filter (Optional[RbfFilter]) – Filter to be used in the target (image) space. (default: None)

kernel

Kernel function \(\phi\) used in the RBF model.

iindex

Indices of integer variables in the feature (domain) space.

filter

Filter to be used in the target (image) space.

__call__(x: ndarray) tuple[ndarray, ndarray]

Evaluates the model at one or multiple points.

Parameters:

x (ndarray) – m-by-d matrix with m point coordinates in a d-dimensional space.

Return type:

tuple[ndarray, ndarray]

Returns:

  • Value for the RBF model on each of the input points.

  • Matrix D where D[i, j] is the distance between the i-th input point and the j-th sampled point.

check_initial_design(sample: ndarray) bool

Check if the sample is able to generate a valid surrogate.

Parameters:

sample (ndarray) – m-by-d matrix with m training points in a d-dimensional space.

Return type:

bool

create_initial_design(dim: int, bounds, minm: int = 0, maxm: int = 0) None

Creates an initial sample for the RBF model.

The points are generated using a symmetric Latin hypercube design.

Parameters:
  • dim (int) – Dimension of the domain space.

  • bounds – List with the limits [x_min,x_max] of each direction x in the domain space.

  • minm (int) – Minimum number of points to generate. If not provided, the initial design will have min(2 * pdim(),maxm) points. (default: 0)

  • maxm (int) – Maximum number of points to generate. If not provided, the initial design will have max(2 * pdim(),minm) points. (default: 0)

Return type:

None

ddpbasis(x: ndarray, p: ndarray) ndarray

Computes the second derivative of the polynomial tail matrix for a given x and direction p.

Parameters:
  • x (ndarray) – Point in a d-dimensional space.

  • p (ndarray) – Direction in which the second derivative is evaluated.

Return type:

ndarray

ddphi(r)

Second derivative of the kernel function phi at the distance(s) r.

Parameters:

r – Vector with distance(s).

dim() int

Get the dimension of the domain space.

Return type:

int

dpbasis(x: ndarray) ndarray

Computes the derivative of the polynomial tail matrix for a given x.

Parameters:

x (ndarray) – Point in a d-dimensional space.

Return type:

ndarray

dphi(r)

Derivative of the kernel function phi at the distance(s) r.

Parameters:

r – Vector with distance(s).

dphiOverR(r)

Derivative of the kernel function phi divided by r at the distance(s) r.

This routine may avoid excessive numerical accumulation errors when phi(r)/r is needed.

Parameters:

r – Vector with distance(s).

get_RBFmatrix() ndarray

Get the complete matrix used to compute the RBF weights.

This is a blocked matrix \([[\Phi, P],[P^T, 0]]\), where \(\Phi\) is the kernel matrix, and \(P\) is the polynomial tail basis matrix.

Return type:

ndarray

Returns:

(m+pdim)-by-(m+pdim) matrix used to compute the RBF weights.

get_iindex() tuple[int, ...]

Return iindex, the sequence of integer variable indexes.

Return type:

tuple[int, ...]

get_matrixP() ndarray

Get the m-by-pdim matrix with the polynomial tail.

Return type:

ndarray

hessp(x: ndarray, p: ndarray) ndarray

Evaluates the Hessian of the model at x in the direction of p.

\[H(f)(x) v = \sum_{i=1}^{m} \beta_i \left( \phi''(\|x - x_i\|)\frac{(x^Tv)x}{\|x - x_i\|^2} + \frac{\phi'(\|x - x_i\|)}{\|x - x_i\|} \left(v - \frac{(x^Tv)x}{\|x - x_i\|^2}\right) \right) + \sum_{i=1}^{n} \beta_{m+i} H(p_i)(x) v.\]
Parameters:
  • x (ndarray) – Point in a d-dimensional space.

  • p (ndarray) – Direction in which the Hessian is evaluated.

Return type:

ndarray

jac(x: ndarray) ndarray

Evaluates the derivative of the model at one point.

\[\nabla f(x) = \sum_{i=1}^{m} \beta_i \frac{\phi'(\|x - x_i\|)}{\|x - x_i\|} x + \sum_{i=1}^{n} \beta_{m+i} \nabla p_i(x).\]
Parameters:

x (ndarray) – Point in a d-dimensional space.

Return type:

ndarray

min_design_space_size(dim: int) int

Return the minimum design space size for a given space dimension.

Return type:

int

ntrain() int

Get the number of sampled points.

Return type:

int

pbasis(x: ndarray) ndarray

Computes the polynomial tail matrix for a given set of points.

Parameters:

x (ndarray) – m-by-d matrix with m point coordinates in a d-dimensional space.

Return type:

ndarray

pdim() int

Get the dimension of the polynomial tail.

Return type:

int

phi(r)

Applies the kernel function phi to the distance(s) r.

Parameters:

r – Vector with distance(s).

reserve(maxeval: int, dim: int) None

Reserve space for the RBF model.

This routine avoids successive dynamic memory allocations with successive calls of update(). If the input maxeval is smaller than the current number of sample points, nothing is done.

Parameters:
  • maxeval (int) – Maximum number of function evaluations.

  • dim (int) – Dimension of the domain space.

Return type:

None

reset() None

Resets the RBF model.

Return type:

None

sample(i: int) ndarray

Get the i-th sampled point.

Parameters:

i (int) – Index of the sampled point.

Return type:

ndarray

update(Xnew, ynew) None

Updates the model with new pairs of data (x,y).

Parameters:
  • Xnew – m-by-d matrix with m point coordinates in a d-dimensional space.

  • ynew – Function values on the sampled points.

Return type:

None

update_coefficients(fx, filter: RbfFilter | None = None) None

Updates the coefficients of the RBF model.

Parameters:
  • fx – Function values on the sampled points.

  • filter (Optional[RbfFilter]) – Filter used with the function values. The default is None, which means the filter used in the initialization of the RBF model is used. (default: None)

Return type:

None

update_xtrain(xNew: ndarray, distNew=None) None

Updates the RBF model with new points, without retraining.

For training the model one still needs to call update_coefficients().

Parameters:
  • xNew (ndarray) – m-by-d matrix with m point coordinates in a d-dimensional space.

  • distNew – m-by-(self.ntrain() + m) matrix with distances between points in xNew and points in (self.xtrain(), xNew). If not provided, the distances are computed. (default: None)

Return type:

None

xtrain() ndarray

Get the training data points.

Return type:

ndarray

Returns:

m-by-d matrix with m training points in a d-dimensional space.

ytrain() ndarray

Get f(x) for the sampled points.

Return type:

ndarray