Class interpm_krige (o2scl)¶
-
template<class
vec_t
= boost::numeric::ublas::vector<double>, classmat_t
= boost::numeric::ublas::vector<double>, classmat_row_t
= boost::numeric::ublas::matrix_row<boost::numeric::ublas::vector<double>>>
classo2scl
::
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.
Select matrix inversion method
-
size_t
np
¶ The number of points.
-
size_t
nd_in
¶ The number of dimensions of the inputs.
-
size_t
nd_out
¶ The number of dimensions of the outputs.
-
bool
data_set
¶ True if the data has been specified.
-
bool
rescaled
¶ True if the data needs to be rescaled.
-
size_t
matrix_mode
¶ Method for matrix inversion.
-
int
verbose
¶ Verbosity parameter (default 0)
-
const size_t
matrix_cholesky
= 0¶ Use Cholesky decomposition.
-
const size_t
matrix_LU
= 1¶ Use LU decomposition.
-
template<class
mat2_row_t
, classmat2_t
, classfunc_vec_t
>
intset_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.
-
template<class
mat2_row_t
, classmat2_t
, classfunc_vec_t
>
intset_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
, classvec3_t
, classvec_func_t
>
voideval
(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 Types
-
typedef boost::numeric::ublas::vector<double>
ubvector
¶
-
typedef boost::numeric::ublas::matrix<double>
ubmatrix
¶
-
typedef boost::numeric::ublas::vector<size_t>
ubvector_size_t
¶
Public Functions
-
interpm_krige
()¶