SOT
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
sot::RBFInterpolant< RBFKernel, PolyTail > Class Template Reference

Radial basis function. More...

#include <rbf.h>

Inheritance diagram for sot::RBFInterpolant< RBFKernel, PolyTail >:
sot::Surrogate sot::RBFInterpolantCap< RBFKernel, PolyTail >

Public Member Functions

 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...
 

Protected Member Functions

void setPoints (const mat &ppoints, const vec &funVals)
 Computes the initial LU decomposition. More...
 

Protected Attributes

RBFKernel mKernel
 
PolyTail mTail
 
mat mL
 
mat mU
 
uvec mp
 
vec mCoeffs
 
vec mfX
 
mat mCenters
 
int mMaxPoints
 
int mDimTail
 
int mNumPoints
 
bool mDirty
 
vec mxLow
 
vec mxUp
 
double mEta = 1e-8
 
int mDim
 

Detailed Description

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
RBFKernelRadial kernel (Cubic is default)
PolyTailPolynomial tail (Linear is default)
Author
David Eriksson, dme65.nosp@m.@cor.nosp@m.nell..nosp@m.edu

Constructor & Destructor Documentation

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
sot::RBFInterpolant< RBFKernel, PolyTail >::RBFInterpolant ( int  maxPoints,
int  dim 
)
inline

Constructor.

Parameters
maxPointsCapacity
dimNumber of dimensions
template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
sot::RBFInterpolant< RBFKernel, PolyTail >::RBFInterpolant ( int  maxPoints,
int  dim,
double  eta 
)
inline

Constructor.

Parameters
maxPointsCapacity
dimNumber of dimensions
etaDamping coefficient (non-negative)
template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
sot::RBFInterpolant< RBFKernel, PolyTail >::RBFInterpolant ( int  maxPoints,
int  dim,
vec  xLow,
vec  xUp,
double  eta 
)
inline

Constructor.

Parameters
maxPointsCapacity
dimNumber of dimensions
xLowLower variable bounds
xUpUpper variable bounds
etaDamping coefficient (non-negative)
template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
sot::RBFInterpolant< RBFKernel, PolyTail >::RBFInterpolant ( int  maxPoints,
int  dim,
vec  xLow,
vec  xUp 
)
inline

Constructor with default eta.

Parameters
maxPointsCapacity
dimNumber of dimensions
xLowLower variable bounds
xUpUpper variable bounds
Exceptions
std::logic_errorIf the polynomial tail degree is too low

Member Function Documentation

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
void sot::RBFInterpolant< RBFKernel, PolyTail >::addPoint ( const vec ppoint,
double  funVal 
)
inlinevirtual

Method for adding one points with known value.

Parameters
ppointPoint to be added
funValFunction value at the point
Exceptions
std::logic_errorif capacity is exceeded

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
void sot::RBFInterpolant< RBFKernel, PolyTail >::addPoints ( const mat ppoints,
const vec funVals 
)
inlinevirtual

Method for adding multiple points with known values.

Parameters
ppointsPoints to be added
funValsFunction values at the points
Exceptions
std::logic_errorif one point is supplied or if capacity is exceeded

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::coeffs ( )
inline

Method for getting the RBF interpolation coefficients.

Returns
RBF interpolation coefficients.
Exceptions
std::logic_errorIf the polynomial tail degree is too low
template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::deriv ( const vec ppoint) const
inlinevirtual

Method for evaluating the derivative of the surrogate model at a point.

Parameters
ppointPoint for which to evaluate the derivative of the surrogate model
Returns
Value of the derivative of the surrogate model at the point
Exceptions
std::logic_errorif coefficients aren't updated

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
int sot::RBFInterpolant< RBFKernel, PolyTail >::dim ( ) const
inlinevirtual

Method for getting the number of dimensions.

Returns
Number of dimensions

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
double sot::RBFInterpolant< RBFKernel, PolyTail >::eval ( const vec ppoint) const
inlinevirtual

Method for evaluating the surrogate model at a point.

Parameters
ppointPoint for which to evaluate the surrogate model
Returns
Value of the surrogate model at the point
Exceptions
std::logic_errorif coefficients aren't updated

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
double sot::RBFInterpolant< RBFKernel, PolyTail >::eval ( const vec ppoint,
const vec dists 
) const
inlinevirtual

Method for evaluating the surrogate model at a point with known distances.

Parameters
ppointPoint for which to evaluate the surrogate model
distsDistances between the interpolation nodes and the point
Returns
Value of the surrogate model at the point
Exceptions
std::logic_errorif coefficients aren't updated

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::evals ( const mat ppoints) const
inlinevirtual

Method for evaluating the surrogate model at multiple points.

Parameters
ppointsPoints for which to evaluate the surrogate model
Returns
Values of the surrogate model at the points
Exceptions
std::logic_errorif coefficients aren't updated

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::evals ( const mat ppoints,
const mat dists 
) const
inlinevirtual

Method for evaluating the surrogate model at multiple points with known distances.

Parameters
ppointsPoints for which to evaluate the surrogate model
distsDistances between the interpolation nodes and the points
Returns
Value of the surrogate model at the point
Exceptions
std::logic_errorif coefficients aren't updated

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::fX ( ) const
inlinevirtual

Method for getting the values of the current points.

Returns
Values of current points

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
double sot::RBFInterpolant< RBFKernel, PolyTail >::fX ( int  i) const
inlinevirtual

Method for getting the value of current point number i (0 is the first)

Returns
Value of point number i

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
int sot::RBFInterpolant< RBFKernel, PolyTail >::numPoints ( ) const
inlinevirtual

Method for getting the current number of points.

Returns
Current number of points

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
void sot::RBFInterpolant< RBFKernel, PolyTail >::setPoints ( const mat ppoints,
const vec funVals 
)
inlineprotected

Computes the initial LU decomposition.

Parameters
ppointsInitial points
funValsValues at the initial points
Exceptions
std::logic_errorif the number of points are less than the dimension of the polynomial space
template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
mat sot::RBFInterpolant< RBFKernel, PolyTail >::X ( ) const
inlinevirtual

Method for getting the current points.

Returns
Current points

Implements sot::Surrogate.

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::X ( int  i) const
inlinevirtual

Method for getting current point number i (0 is the first)

Returns
Point number i

Implements sot::Surrogate.

Member Data Documentation

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
mat sot::RBFInterpolant< RBFKernel, PolyTail >::mCenters
protected

Interpolation nodes (centers)

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::mCoeffs
protected

Coefficient vector

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
int sot::RBFInterpolant< RBFKernel, PolyTail >::mDim
protected

Number of dimensions

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
int sot::RBFInterpolant< RBFKernel, PolyTail >::mDimTail
protected

Dimensionality of the polynomial space

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
bool sot::RBFInterpolant< RBFKernel, PolyTail >::mDirty
protected

True if the coefficients need to be recomputed

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
double sot::RBFInterpolant< RBFKernel, PolyTail >::mEta = 1e-8
protected

Damping added to the kernel to avoid ill-conditioning

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::mfX
protected

Function values

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
RBFKernel sot::RBFInterpolant< RBFKernel, PolyTail >::mKernel
protected

The radial kernel

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
mat sot::RBFInterpolant< RBFKernel, PolyTail >::mL
protected

Lower triangular part of the LU decomposition

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
int sot::RBFInterpolant< RBFKernel, PolyTail >::mMaxPoints
protected

Capacity

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
int sot::RBFInterpolant< RBFKernel, PolyTail >::mNumPoints
protected

Current number of points

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
uvec sot::RBFInterpolant< RBFKernel, PolyTail >::mp
protected

Permutation vector

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
PolyTail sot::RBFInterpolant< RBFKernel, PolyTail >::mTail
protected

The polynomial tail

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
mat sot::RBFInterpolant< RBFKernel, PolyTail >::mU
protected

Upper triangular part of the LU decomposition

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::mxLow
protected

Lower variable bounds

template<class RBFKernel = CubicKernel, class PolyTail = LinearTail>
vec sot::RBFInterpolant< RBFKernel, PolyTail >::mxUp
protected

Upper variable bounds


The documentation for this class was generated from the following file: