|
| RBFInterpolant (int maxPoints, int dim) |
| Constructor. More...
|
|
| RBFInterpolant (int maxPoints, int dim, double eta) |
| Constructor. More...
|
|
| RBFInterpolant (int maxPoints, int dim, vec xLow, vec xUp, double eta) |
| Constructor. More...
|
|
| RBFInterpolant (int maxPoints, int dim, vec xLow, vec xUp) |
| Constructor with default eta. More...
|
|
int | dim () const |
| Method for getting the number of dimensions. More...
|
|
void | reset () |
| Method for resetting the surrogate model.
|
|
int | numPoints () const |
| Method for getting the current number of points. More...
|
|
mat | X () const |
| Method for getting the current points. More...
|
|
vec | X (int i) const |
| Method for getting current point number i (0 is the first) More...
|
|
vec | fX () const |
| Method for getting the values of the current points. More...
|
|
double | fX (int i) const |
| Method for getting the value of current point number i (0 is the first) More...
|
|
vec | coeffs () |
| Method for getting the RBF interpolation coefficients. More...
|
|
void | fit () |
| Method for fitting the surrogate model.
|
|
void | addPoint (const vec &ppoint, double funVal) |
| Method for adding one points with known value. More...
|
|
void | addPoints (const mat &ppoints, const vec &funVals) |
| Method for adding multiple points with known values. More...
|
|
double | eval (const vec &ppoint) const |
| Method for evaluating the surrogate model at a point. More...
|
|
double | eval (const vec &ppoint, const vec &dists) const |
| Method for evaluating the surrogate model at a point with known distances. More...
|
|
vec | evals (const mat &ppoints) const |
| Method for evaluating the surrogate model at multiple points. More...
|
|
vec | evals (const mat &ppoints, const mat &dists) const |
| Method for evaluating the surrogate model at multiple points with known distances. More...
|
|
vec | deriv (const vec &ppoint) const |
| Method for evaluating the derivative of the surrogate model at a point. More...
|
|
template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
class sot::RBFInterpolant< RBFKernel, PolyTail >
Radial basis function.
A radial basis function (RBF) interpolant is a weighted sum of radial basis functions. It is common to add a polynomial tail as well to assure that the interpolant can exactly reproduce polynomial of that degree. This leads to an interpolant of the form
\( s(y) = \displaystyle\sum_{i=1}^n \lambda_i \varphi(\|y-x_i\|) + \displaystyle\sum_{i=1}^m c_i \pi_i(y) \)
where \( y_i\) are the n centers and \( \{\pi_i\}_{i=1}^m\) is a basis of the polynomial space of the tail. Given a set of points \( X=\{x_1,\ldots,x_n\}\) with values \( f_X = \{f(x_1),\ldots,f(x_n)\}\) the interpolation conditions are
\( s(x_i) = f(x_i)\) for \( i=1,\ldots,n\)
and in order to get a unique interpolant one usually adds the conditions
\( \displaystyle\sum_{i=1}^n \lambda_i \pi_j(x_i) = 0\) for \( j=1,\ldots,m\)
which leads to a system of equations of size \((m+n) x (m+n)\). With the notation \( \Phi_{i,j}=\varphi(\|x_i-x_j\|)\) and \( P_{ij} = \pi_j(x_i)\) we can write in a more compact form
\( \left(\begin{array}{cc} 0 & P^T \\ P & \Phi \end{array}\right) \left(\begin{array}{c} c \\ \lambda \end{array}\right) = \left(\begin{array}{c} 0 \\ f_X \end{array}\right) \)
We can see that adding one more point corresponds to adding a column and a row to this matrix which is why this ordering of \(\Phi\) and \(P\) is convenient. We store the LU decomposition of this matrix and use the fact that the Schur complement is symmetric and positive definite so we can update the LU-decomposition by computing the Cholesky decomposition of the Schur complement. We can hence add \(k\) new points using roughly \(2n^2k\) flops under the assumption that \(n\gg k\). Solving for the new coefficients is then a matter of back and forward subsitution which will take roughly \(2(m+n)^2\) flops. This is better than solving the system from scratch which takes roughly \((2/3)(m+n)^3\) flops.
In order to compute the initial decomposition we need at least m initial points that serve as "outer" points that are needed to uniquely fit the polynomial tail and to construct the initial LU-decomposition with pivoting. No pivoting is needed from that point.
The domain is automatically scaled to the unit box to avoid scaling issues since the kernel and polynomial tail scale differently.
- Template Parameters
-
RBFKernel | Radial kernel (Cubic is default) |
PolyTail | Polynomial tail (Linear is default) |
- Author
- David Eriksson, dme65.nosp@m.@cor.nosp@m.nell..nosp@m.edu