Class interpm_krige (o2scl)¶
-
template<class vec_t, class mat_t, class mat_row_t, class mat_col_t, class mat2_t, class mat2_row_t, class mat3_t, class mat_inv_t = o2scl_linalg::matrix_invert_det_cholesky<mat3_t>>
class o2scl::interpm_krige¶ Multi-dimensional interpolation by kriging.
Note
The set data functions for this class uses a particular format, one different format than that in o2scl::interpm_idw . This design choice makes it easier to pass vector arguments to the covariance function and the linear algebra routines. The x and y objects should be of the form
x[n_points][n_in]
andy[n_out][n_points]
. A separate covariance function is required for each output.Note
This class assumes that the function specified in the call to set_data() is the same as that passed to the eval() functions. If this is not the case, the behavior of this class is undefined.
Note
Experimental.
Public Functions
-
inline interpm_krige()¶
-
template<class func_vec_t>
inline int set_data_noise(size_t n_in, size_t n_out, size_t n_points, mat_t &user_x, mat2_t &user_y, func_vec_t &fcovar, const vec_t &noise_var, bool rescale = false, bool err_on_fail = true)¶ Initialize the data for the interpolation.
Note
This function works differently than o2scl::interpm_idw::set_data() . See this class description for more details.
-
inline void unscale()¶
Remove the rescaling of the data.
-
template<class func_vec_t>
inline int set_data(size_t n_in, size_t n_out, size_t n_points, mat_t &user_x, mat2_t &user_y, func_vec_t &fcovar, bool rescale = false, bool err_on_fail = true)¶ Initialize the data for the interpolation.
Note
This function works differently than o2scl::interpm_idw::set_data() . See this class description for more details.
-
template<class vec2_t, class vec3_t, class vec_func_t>
inline void eval(const vec2_t &x0, vec3_t &y0, vec_func_t &fcovar)¶ Given covariance function
fcovar
and input vectorx
store the result of the interpolation iny
.
Public Members
-
bool keep_matrix¶
If true, keep \( K^{-1} \) (default true)
-
int verbose¶
Verbosity parameter (default 0)
-
inline interpm_krige()¶