Kernel functions¶
Here we provide a description of the different classes used to build
and evaluate the RBF and the GRBF surrogate models. Three different
kernel functions \(\phi(r)\) are defined in the code, the Exponential
kernel (RBF_Exponential,
GRBF_Exponential),
the Matérn kernel (RBF_Matern,
GRBF_Matern) and the Cubic kernel
(RBF_Cubic,
GRBF_Cubic).
The fit method of these classes solves the system \(A s = F\).
In the case of RBF interpolants, the vector \(F\) contains the values
of the function where the points have been evaluated and the matrix \(A\)
is defined as:
In the case of GRBF interpolants, the vector \(F\) contains both the values of the function and its gradient and the matrix \(A\) is defined as:
where \(\Phi\) is the RBF kernel matrix, \(P\) is a matrix with linear polynomial terms, \(\Phi_\mathrm{d}\) is the first derivative of the RBF kernel matrix \(\Phi\) and \(\Phi_\mathrm{dd}\) is the second derivative 1, 2.
The matrix \(\Phi\) is defined as:
where \(r_{i,j} = \left \| x^i-x^j \right\|\) is the Euclidean distance between the points \(x^i\) and \(x^j\). Tha matrices \(\Phi_\mathrm{d}\) and \(\Phi_\mathrm{dd}\) are defined as:
where \(x_k^i\) is the kth component of the ith point in the vector of evaluated points x.
Once the fit method has been used to build the surrogate model, it is possible
to evaluate points using the method evaluate. If needed, the internal parameters
can be updated using the update method.
- 1
Regis, R G and C A Shoemaker. 2013. Combining radial basis function surrogates and dynamic coordinate search in high-dimensional expensive black-box optimization. Engineering Optimization 45 (5): 529-555.
- 2
Laurent, L, R Le Riche, B Soulier and P A Boucard. 2019. An Overview of Gradient-Enhanced Metamodels with Applications. Archives of Computational Methods in Engineering 26 (1). 61-106.
Exponential Kernel¶
The Exponential Kernel is defined as follows:
where \(l>0\) and \(r=\|y-x\|\), where \(\|\cdot\|\) is the Euclidean norm.
The first derivative of the Exponential Kernel with respect to the distance is defined as:
The second derivative is defined as:
- class DyCors.kernels.exponential.GRBF_Exponential(l=1.0)[source]¶
GRBF Exponential kernel class.
- Parameters
- lfloat or ndarray, shape (d,), optional
Internal parameter. Width of the kernel.
- Attributes
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- sndarray, shape(m*(d+1),)
GRBF coefficients.
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.- dfndarray, shape (m*d,)
Array of function gradient values at
x.
- evaluate(y)[source]¶
Evaluate surrogate model at given points.
- Parameters
- yndarray, shape (n,d,)
Array of points where we want to evaluate the surrogate model.
- Returns
- fndarray, shape(n,)
Array of interpolated values.
- fit(x, f, df)[source]¶
Build surrogate model.
- Parameters
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.- dfndarray, shape (m*d,)
Array of function gradient values at
x.
- Returns
- Phindarray, shape(m,m,)
RBF matrix.
- Andarray, shape(m*(d+1),m*(d+1),)
GRBF matrix with gradient terms.
- class DyCors.kernels.exponential.RBF_Exponential(l=1.0)[source]¶
RBF Exponential kernel class.
- Parameters
- lfloat or ndarray, shape (d,), optional
Internal parameter. Width of the kernel.
- Attributes
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- sndarray, shape(m+d+1,)
RBF coefficients.
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.
- evaluate(y)[source]¶
Evaluate surrogate model at given points.
- Parameters
- yndarray, shape (n,d,)
Array of points where we want to evaluate the surrogate model.
- Returns
- fndarray, shape(n,)
Array of interpolated values.
- fit(x, f)[source]¶
Build surrogate model.
- Parameters
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.
- Returns
- Phindarray, shape(m,m,)
RBF matrix.
- Andarray, shape(m*(d+1),m*(d+1),)
RBF matrix with linear polynomial terms.
- optimize()[source]¶
Optimize internal parameters based on the Leave-One-Out error using a differential evolution algorithm [1], [2]. A constraint is applied to the condition number of the kernel matrix to ensure the smoothness of the function.
References
- 1
Rippa, S. 1999. An algorithm for selecting a good value for the parameter c in radial basis function interpolation. Advances in Computational Mathematics 11 (2). 193-210.
- 2
Bompard, M, J Peter and J A Desideri. 2010. Surrogate models based on function and derivative values for aerodynamic global optimization. V European Conference on Computational Fluid Dynamics ECCOMAS CFD 2010, ECCOMAS. Lisbon, Portugal.
Matérn Kernel¶
The half integer simplification of the Matérn Kernel is defined as follows:
where \(p=\left \lfloor{\nu - 1/2}\right \rfloor , l>0, \nu>=1/2\) and \(r=\|y-x\|\), where \(\|\cdot\|\) is the Euclidean norm.
The first derivative of the Matérn Kernel with respect to the distance is defined as:
The second derivative is defined as:
- class DyCors.kernels.matern.GRBF_Matern(l=1.0, nu=2.5)[source]¶
GRBF Matérn kernel class.
- Parameters
- lfloat or ndarray, shape (d,), optional
Internal parameter of the kernel.
- nufloat, optional
Internal parameter. Order of Bessel function of the kernel.
- Attributes
- lfloat or ndarray, shape (d,)
Internal parameter of the kernel.
- nufloat
Internal parameter. Order of Bessel function of the kernel.
- sndarray, shape(m*(d+1),)
GRBF coefficients.
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.- dfndarray, shape (m*d,)
Array of function gradient values at
x.
- evaluate(y)[source]¶
Evaluate surrogate model at given points.
- Parameters
- yndarray, shape (n,d,)
Array of points where we want to evaluate the surrogate model.
- Returns
- fndarray, shape(n,)
Array of interpolated values.
- fit(x, f, df)[source]¶
Build surrogate model.
- Parameters
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.- dfndarray, shape (m*d,)
Array of function gradient values at
x.
- Returns
- Phindarray, shape(m,m,)
RBF matrix.
- Andarray, shape(m*(d+1),m*(d+1),)
GRBF matrix with gradient terms.
- class DyCors.kernels.matern.RBF_Matern(l=1.0, nu=2.5)[source]¶
RBF Matérn kernel class.
- Parameters
- lfloat or ndarray, shape (d,), optional
Internal parameter. Width of the kernel.
- nufloat, optional
Internal parameter. Order of Bessel function of the kernel.
- Attributes
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- nufloat
Internal parameter. Order of Bessel function of the kernel.
- sndarray, shape(m+d+1,)
RBF coefficients.
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.
- evaluate(y)[source]¶
Evaluate surrogate model at given points.
- Parameters
- yndarray, shape (n,d,)
Array of points where we want to evaluate the surrogate model.
- Returns
- fndarray, shape(n,)
Array of interpolated values.
- fit(x, f)[source]¶
Build surrogate model.
- Parameters
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.
- Returns
- Phindarray, shape(m,m,)
RBF matrix.
- Andarray, shape(m*(d+1),m*(d+1),)
RBF matrix with linear polynomial terms.
- optimize()[source]¶
Optimize internal parameters based on the Leave-One-Out error using a differential evolution algorithm [1], [2]. A constraint is applied to the condition number of the kernel matrix to ensure the smoothness of the function.
References
- 1
Rippa, S. 1999. An algorithm for selecting a good value for the parameter c in radial basis function interpolation. Advances in Computational Mathematics 11 (2). 193-210.
- 2
Bompard, M, J Peter and J A Desideri. 2010. Surrogate models based on function and derivative values for aerodynamic global optimization. V European Conference on Computational Fluid Dynamics ECCOMAS CFD 2010, ECCOMAS. Lisbon, Portugal.
Cubic Kernel¶
The Cubic Kernel is defined as follows:
where \(l>0\), \(r=\|y-x\|\), where \(\|\cdot\|\) is the Euclidean norm.
The first derivative of the Cubic Kernel with respect to the distance is defined as:
The second derivative is defined as:
- class DyCors.kernels.cubic.GRBF_Cubic(l=1.0)[source]¶
GRBF Cubic kernel class.
- Parameters
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- Attributes
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- sndarray, shape(m*(d+1),)
GRBF coefficients.
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.- dfndarray, shape (m*d,)
Array of function gradient values at
x.
- evaluate(y)[source]¶
Evaluate surrogate model at given points.
- Parameters
- yndarray, shape (n,d,)
Array of points where we want to evaluate the surrogate model.
- Returns
- fndarray, shape(n,)
Array of interpolated values.
- fit(x, f, df)[source]¶
Build surrogate model.
- Parameters
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.- dfndarray, shape (m*d,)
Array of function gradient values at
x.
- Returns
- Phindarray, shape(m,m,)
RBF matrix.
- Andarray, shape(m*(d+1),m*(d+1),)
GRBF matrix with gradient terms.
- class DyCors.kernels.cubic.RBF_Cubic(l=1.0)[source]¶
RBF Cubic kernel class.
- Parameters
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- Attributes
- lfloat or ndarray, shape (d,)
Internal parameter. Width of the kernel.
- sndarray, shape(m+d+1,)
RBF coefficients.
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.
- evaluate(y)[source]¶
Evaluate surrogate model at given points.
- Parameters
- yndarray, shape (n,d,)
Array of points where we want to evaluate the surrogate model.
- Returns
- fndarray, shape(n,)
Array of interpolated values.
- fit(x, f)[source]¶
Build surrogate model.
- Parameters
- xndarray, shape (m,d,)
Array of points where function values are known. m is the number of sampling points and d is the number of dimensions.
- fndarray, shape (m,)
Array of function values at
x.
- Returns
- Phindarray, shape(m,m,)
RBF matrix.
- Andarray, shape(m*(d+1),m*(d+1),)
RBF matrix with linear polynomial terms.
- optimize()[source]¶
Optimize internal parameters based on the Leave-One-Out error using a differential evolution algorithm [1], [2]. A constraint is applied to the condition number of the kernel matrix to ensure the smoothness of the function.
References
- 1
Rippa, S. 1999. An algorithm for selecting a good value for the parameter c in radial basis function interpolation. Advances in Computational Mathematics 11 (2). 193-210.
- 2
Bompard, M, J Peter and J A Desideri. 2010. Surrogate models based on function and derivative values for aerodynamic global optimization. V European Conference on Computational Fluid Dynamics ECCOMAS CFD 2010, ECCOMAS. Lisbon, Portugal.