Main Page | Namespace List | Class List | File List | Namespace Members | Class Members | File Members

CGeostat Class Reference

variogram & universal kriging class More...

#include <geostat.h>

Collaboration diagram for CGeostat:

Collaboration graph
[legend]
List of all members.

Public Member Functions

 CGeostat (void)
 default constructor
virtual ~CGeostat (void)
 destructor
bool initialize (size_t smpl, size_t dim)
 initializer
bool isActive (void) const
size_t samples (void) const
 get number of samples
size_t dimension (void) const
 get dimension of control space
double maxDist (void) const
 get maximum distance in control space
bool isDomain (const gsl_vector *c) const
 domain validation
bool getCoordinate (gsl_vector *c, size_t s) const
 get coordinate of sample s
double getData (size_t s) const
 get sample s
bool getModelParameters (double &sill, double &range, double &nugget, double &power) const
 get parameters of estimated variogram
bool getWeightVector (gsl_vector *w, const gsl_vector *c) const
 estimation of optimum weights
bool getPredictData (double &pred, double &var, const gsl_vector *c) const
 prediction
bool getPredictData (double &pred, const gsl_vector *coord, const gsl_vector *w) const
 prediction using weights
bool getDomain (gsl_vector *lower, gsl_vector *upper) const
 get control space domain
double getTrend (const gsl_vector *c) const
 get trend component
double getResidual (size_t s) const
 get residual component
bool setDomain (const gsl_vector *lower, const gsl_vector *upper)
 set control space domain
bool setCoordinate (size_t s, const gsl_vector *c)
 set coordinate of sample s
bool setData (size_t s, double d)
 set sample data
bool estimate (int model, double power, double step)
 trend surface approximation, theoretical variogram estimate, and kriging system precomputation
bool removeRedundantSample (size_t s)
 removal of redundant sample
double evalModelValidity (int model, double power, double step)
 cross validation of theoretical variogram
bool getPredictionErrorExcluded (double &pred, double &var, size_t smpl) const
 leave-one-out cross validation
double selectEliminatee (size_t &candidate) const
 greedy search of eliminatee

Protected Member Functions

bool approximateTrendSurface (void)
 trend surface approximation by pseudo-inverse matrix
std::vector< CVariogram::VarioItmcomputeVariogramCloud (void) const
 compute variogram cloud
std::vector< CVariogram::VarioItmcomputeExperimentalVariogram (const std::vector< CVariogram::VarioItm > &vcloud, double step) const
 compute experimental variogram
bool precomputeKrigingSystem (void)
 precompute kriging system

Detailed Description

variogram & universal kriging class

Author:
Tomohiko Mukai

Definition at line 22 of file geostat.h.


Constructor & Destructor Documentation

CGeostat::CGeostat void   ) 
 

default constructor

Definition at line 40 of file geostat.cpp.

CGeostat::~CGeostat void   )  [virtual]
 

destructor

Definition at line 55 of file geostat.cpp.


Member Function Documentation

bool CGeostat::approximateTrendSurface void   )  [protected]
 

trend surface approximation by pseudo-inverse matrix

Return values:
true success
false fail

Definition at line 607 of file geostat.cpp.

std::vector< CVariogram::VarioItm > CGeostat::computeExperimentalVariogram const std::vector< CVariogram::VarioItm > &  vcloud,
double  step
const [protected]
 

compute experimental variogram

Parameters:
vcloud [in] variogram cloud
step [in] step width of equally devided region
Returns:
experimental variogram

Definition at line 335 of file geostat.cpp.

std::vector< CVariogram::VarioItm > CGeostat::computeVariogramCloud void   )  const [protected]
 

compute variogram cloud

Returns:
variogram cloud

Definition at line 292 of file geostat.cpp.

size_t CGeostat::dimension void   )  const [inline]
 

get dimension of control space

Returns:
dimension of control space

Definition at line 80 of file geostat.h.

bool CGeostat::estimate int  model,
double  power,
double  step
 

trend surface approximation, theoretical variogram estimate, and kriging system precomputation

Parameters:
model [in] type of theoretical variogram
power [in] power coefficient of stable variograms
step [in] step width
Return values:
true success
false fail

Definition at line 374 of file geostat.cpp.

double CGeostat::evalModelValidity int  model,
double  power,
double  step
 

cross validation of theoretical variogram

Parameters:
model [in] type of theoretical variogram
power [in] power coefficient of stable variograms
step [in] step width
Returns:
invalidity

Definition at line 898 of file geostat.cpp.

bool CGeostat::getCoordinate gsl_vector *  c,
size_t  s
const
 

get coordinate of sample s

Parameters:
c [out] coordinate of sample s
s [in] sample index
Return values:
true success
false fail

Definition at line 149 of file geostat.cpp.

double CGeostat::getData size_t  s  )  const
 

get sample s

Parameters:
s [in] sample index
Returns:
sample value

Definition at line 164 of file geostat.cpp.

bool CGeostat::getDomain gsl_vector *  lower,
gsl_vector *  upper
const
 

get control space domain

Parameters:
lower [out] lower bounds
upper [out] upper bounds
Return values:
true success
false fail

Definition at line 216 of file geostat.cpp.

bool CGeostat::getModelParameters double &  sill,
double &  range,
double &  nugget,
double &  power
const
 

get parameters of estimated variogram

Parameters:
sill [out] sill
range [out] range
nugget [out] nugget
power [out] power coefficient
Return values:
true success
false fail

Definition at line 408 of file geostat.cpp.

bool CGeostat::getPredictData double &  pred,
const gsl_vector *  c,
const gsl_vector *  weight
const
 

prediction using weights

Parameters:
pred [out] predictive value
c [in] target parameter
weight [in] blending weights
Return values:
true success
false fail

Definition at line 589 of file geostat.cpp.

bool CGeostat::getPredictData double &  pred,
double &  var,
const gsl_vector *  c
const
 

prediction

Parameters:
pred [out] predictive value
var [out] estimation variance
c [in] target parameter
Return values:
true success
false fail

Definition at line 544 of file geostat.cpp.

bool CGeostat::getPredictionErrorExcluded double &  pred,
double &  var,
size_t  s
const
 

leave-one-out cross validation

Parameters:
pred [out] predictive value
var [out] estimation variance
s [in] sample index
Return values:
true success
false fail

Definition at line 694 of file geostat.cpp.

double CGeostat::getResidual size_t  s  )  const
 

get residual component

Parameters:
s [in] sample index
Returns:
residual component

Definition at line 200 of file geostat.cpp.

double CGeostat::getTrend const gsl_vector *  c  )  const
 

get trend component

Parameters:
c [in] target parameter
Returns:
trend component

Definition at line 182 of file geostat.cpp.

bool CGeostat::getWeightVector gsl_vector *  weight,
const gsl_vector *  c
const
 

estimation of optimum weights

Parameters:
weight [out] blending weight
c [in] target parameter
Return values:
success 
fail 

Definition at line 507 of file geostat.cpp.

bool CGeostat::initialize size_t  smpl,
size_t  dim
 

initializer

Parameters:
smpl [in] number of samples
dim [in] dimension of control space
Return values:
true success
false fail

Definition at line 68 of file geostat.cpp.

bool CGeostat::isActive void   )  const [inline]
 

Return values:
true active
false inactive

Definition at line 61 of file geostat.h.

bool CGeostat::isDomain const gsl_vector *  c  )  const
 

domain validation

Parameters:
c [in] coordinate
Return values:
true success
false domaion error

Definition at line 488 of file geostat.cpp.

double CGeostat::maxDist void   )  const
 

get maximum distance in control space

Returns:
maximum distance in control space

Definition at line 472 of file geostat.cpp.

bool CGeostat::precomputeKrigingSystem void   )  [protected]
 

precompute kriging system

Return values:
success 
fail 

Definition at line 427 of file geostat.cpp.

bool CGeostat::removeRedundantSample size_t  smpl  ) 
 

removal of redundant sample

Parameters:
smpl [in] sample index
Return values:
true success
false fail

Definition at line 646 of file geostat.cpp.

size_t CGeostat::samples void   )  const [inline]
 

get number of samples

Returns:
number of samples

Definition at line 70 of file geostat.h.

double CGeostat::selectEliminatee size_t &  candidate  )  const
 

greedy search of eliminatee

Parameters:
candidate [out] candidate index
Returns:
prediction reliability after removal

Definition at line 791 of file geostat.cpp.

bool CGeostat::setCoordinate size_t  s,
const gsl_vector *  c
 

set coordinate of sample s

Parameters:
s [in] sample index
c [in] coordinate
Return values:
true success
false fail

Definition at line 260 of file geostat.cpp.

bool CGeostat::setData size_t  s,
double  d
 

set sample data

Parameters:
s [in] sample index
d [in] data
Return values:
true success
false fail

Definition at line 278 of file geostat.cpp.

bool CGeostat::setDomain const gsl_vector *  lower,
const gsl_vector *  upper
 

set control space domain

Parameters:
lower [in] lower bounds
upper [in] upper bounds
Return values:
true success
false fail

Definition at line 238 of file geostat.cpp.


Generated on Thu Aug 11 13:01:38 2005 for Kriging sample by  doxygen 1.4.3