Class interp_krige (o2scl)¶
-
template<class vec_t, class vec2_t = vec_t, class covar_func_t = std::function<double(double, double)>, class covar_integ_t = std::function<double(double, double, double)>, class mat_t = boost::numeric::ublas::matrix<double>, class mat_inv_t = o2scl_linalg::matrix_invert_det_cholesky<mat_t>>
class interp_krige : public o2scl::interp_base<vec_t, vec_t>¶ Interpolation by Kriging with a user-specified covariance function.
See also the Interpolation section of the O2scl User’s guide.
Note
The function set() stores a pointer to the covariance function and its derivatives and integrals so they cannot go out of scope before any of the interpolation functions are called.
Note
This class is experimental.
Subclassed by o2scl::interp_krige_optim< vec_t, vec2_t, func_t, mat_t, mat_inv_t, vec_vec_t >
Public Types
-
typedef boost::numeric::ublas::vector<double> ubvector¶
Public Functions
-
inline interp_krige()¶
-
inline virtual ~interp_krige()¶
-
inline virtual void set(size_t size, const vec_t &x, const vec2_t &y)¶
Initialize interpolation routine.
-
inline virtual int set_covar_di(size_t n_dim, const vec_t &x, const vec2_t &y, covar_func_t &fcovar, covar_func_t &fderiv, covar_func_t &fderiv2, covar_integ_t &finteg, bool rescale = false)¶
Initialize interpolation routine with covariance function, derivatives and integrals.
-
inline virtual int set_covar(size_t n_dim, const vec_t &x, const vec2_t &y, covar_func_t &fcovar, bool rescale = false)¶
Initialize interpolation routine with covariance function, but without derivatives and integrals.
-
inline virtual double eval(double x0) const¶
Give the value of the function \( y(x=x_0) \) .
-
template<class covar_func2_t>
inline double eval_covar(double x0, covar_func2_t &user_f) const¶ Evaluate the interpolation at \( x=x_0 \) using an alternate covariance function.
Note
This function only works if the data is not rescaled.
-
inline virtual double deriv(double x0) const¶
Give the value of the derivative \( y^{\prime}(x=x_0) \) .
-
inline virtual double deriv2(double x0) const¶
Give the value of the second derivative \( y^{\prime \prime}(x=x_0) \).
-
inline virtual double integ(double a, double b) const¶
Give the value of the integral \( \int_a^{b}y(x)~dx \) .
-
inline double sigma(double x0) const¶
Return the interpolation uncertainty from the Gaussian process.
-
inline prob_dens_gaussian gen_dist(double x0) const¶
Generate a probability distribution for the interpolation at a specified point.
This function implements Eq. 2.19 of R&W
-
inline double sample(double x0) const¶
Sample the probability distribution for the interpolation at a specified point.
This creates a new distribution at the point \x c0 and then samples that distribution. If one wants to perform several samples at the same point, it is much more efficient to use gen_dist() instead.
-
template<class vec3_t, class vec4_t>
inline void sample_vec(vec3_t &x, vec4_t &y) const¶ Generate a probability distribution for the interpolation at a vector of points.
-
inline virtual const char *type() const¶
Return the type,
"interp_krige"
.
Public Members
-
bool err_nonconv¶
If true, call the error handler for convergence errors (default true)
-
bool keep_matrix¶
If true, keep \( K^{-1} \) (default true)
This must be set to true in order to use sigma() or gen_dist() .
Protected Functions
-
inline virtual int set_covar_di_internal(size_t n_dim, const vec_t &x, const vec2_t &y, covar_func_t &fcovar, covar_func_t *fderiv, covar_func_t *fderiv2, covar_integ_t *finteg, bool rescale = false)¶
Initialize interpolation routine, internal version with pointers for derivative and integral functions.
Protected Attributes
-
covar_func_t *f¶
Pointer to user-specified covariance function.
-
covar_func_t *fd¶
Pointer to user-specified derivative.
-
covar_func_t *fd2¶
Pointer to user-specified second derivative.
-
covar_integ_t *fi¶
Pointer to user-specified integral.
-
bool rescaled¶
If true, then the data was rescaled to zero mean and unit standard deviation.
-
double mean_y¶
Mean before rescaling.
-
double std_y¶
Standard deviation before rescaling.
Private Functions
-
interp_krige(const interp_krige<vec_t, vec2_t, covar_func_t, covar_integ_t, mat_t, mat_inv_t>&)¶
-
interp_krige<vec_t, vec2_t, covar_func_t, covar_integ_t, mat_t, mat_inv_t> &operator=(const interp_krige<vec_t, vec2_t, covar_func_t, covar_integ_t, mat_t, mat_inv_t>&)¶
-
typedef boost::numeric::ublas::vector<double> ubvector¶