1.3.2
CCC
 
krige.cpp File Reference
#include "geoslib_f.h"
#include "geoslib_old_f.h"
#include "geoslib_f_private.h"
#include "geoslib_define.h"
#include "Enum/EAnam.hpp"
#include "Enum/ECalcMember.hpp"
#include "Polynomials/Hermite.hpp"
#include "Db/Db.hpp"
#include "Db/DbGrid.hpp"
#include "Db/DbStringFormat.hpp"
#include "Model/Model.hpp"
#include "Model/CovInternal.hpp"
#include "Neigh/NeighMoving.hpp"
#include "Neigh/NeighImage.hpp"
#include "Neigh/NeighUnique.hpp"
#include "Neigh/ANeigh.hpp"
#include "Anamorphosis/AnamDiscreteDD.hpp"
#include "Anamorphosis/AnamDiscreteIR.hpp"
#include "Anamorphosis/AnamHermite.hpp"
#include "Basic/String.hpp"
#include "Basic/NamingConvention.hpp"
#include "Basic/Utilities.hpp"
#include "Basic/Law.hpp"
#include "Basic/File.hpp"
#include "Basic/OptDbg.hpp"
#include "Basic/OptCustom.hpp"
#include "Covariances/CovContext.hpp"
#include "Drifts/DriftList.hpp"
#include "Estimation/KrigingSystem.hpp"
#include "Space/SpaceRN.hpp"
#include "Matrix/MatrixFactory.hpp"
#include <math.h>
#include <string.h>
#include <stdio.h>
#include <limits.h>

Classes

struct  Disc_Structure
 

Functions

int is_flag_data_disc_defined (void)
 ‍****************************************************************************‍/ *! More...
 
void set_DBIN (Db *dbin)
 
void set_DBOUT (Db *dbout)
 
int krige_koption_manage (int mode, int flag_check, const EKrigOpt &calcul, int flag_rand, const VectorInt &ndiscs)
 
void krige_lhs_print (int nech, int neq, int nred, int *flagloc, const double *lhs)
 
void krige_rhs_print (int nvar, int nech, int neq, int nred, int *flag, double *rhs)
 
void krige_dual_print (int nech, int neq, int nred, int *flag, double *dual)
 
int _krigsim (Db *dbin, Db *dbout, const Model *model, ANeigh *neigh, bool flag_bayes, const VectorDouble &dmean, const MatrixSquareSymmetric &dcov, int icase, int nbsimu, bool flag_dgm)
 
int global_transitive (DbGrid *dbgrid, Model *model, int flag_verbose, int flag_regular, int ndisc, double *abundance, double *sse, double *cvtrans)
 
int anakexp_f (DbGrid *db, double *covdd, double *covd0, double top, double bot, int cov_radius, int neigh_radius, int flag_sym, int nfeq)
 
int anakexp_3D (DbGrid *db, double *cov_ref, int cov_radius, int neigh_ver, int neigh_hor, int flag_sym, Model *model, double nugget, int nfeq, int dbg_ix, int dbg_iy)
 
int bayes_simulate (Model *model, int nbsimu, const VectorDouble &rmean, const VectorDouble &rcov, VectorDouble &smean)
 
int krigsum (Db *dbin, Db *dbout, Model *model, ANeigh *neigh, bool flag_positive, const NamingConvention &namconv)
 
int st_krige_data (Db *db, Model *model, double beta, VectorInt &ranks1, VectorInt &ranks2, VectorInt &rother, int flag_abs, double *data_est, double *data_var)
 
int st_crit_global (Db *db, Model *model, VectorInt &ranks1, VectorInt &rother, double *crit)
 
int sampling_f (Db *db, Model *model, double beta, int method1, int nsize1_max, VectorInt &ranks1, int method2, int nsize2_max, VectorInt &ranks2, int verbose)
 
int krigsampling_f (Db *dbin, Db *dbout, Model *model, double beta, VectorInt &ranks1, VectorInt &ranks2, bool flag_std, int verbose)
 
int declustering (Db *dbin, Model *model, int method, ANeigh *neigh, DbGrid *dbgrid, const VectorDouble &radius, const VectorInt &ndiscs, int flag_sel, bool verbose)
 
int inhomogeneous_kriging (Db *dbdat, Db *dbsrc, Db *dbout, double power, int flag_source, Model *model_dat, Model *model_src)
 
void _image_smoother (DbGrid *dbgrid, const NeighImage *neigh, int type, double range, int iptr0)
 

Function Documentation

◆ _image_smoother()

void _image_smoother ( DbGrid dbgrid,
const NeighImage neigh,
int  type,
double  range,
int  iptr0 
)

Smooth a regular grid

Parameters
[in]dbgridinput and output Db grid structure
[in]neighNeigh structure
[in]type1 for Uniform; 2 for Gaussian
[in]rangeRange (used for Gaussian only)
[in]iptr0Storage address
Remarks
Limited to the monovariate case

◆ _krigsim()

int _krigsim ( Db dbin,
Db dbout,
const Model model,
ANeigh neigh,
bool  flag_bayes,
const VectorDouble dmean,
const MatrixSquareSymmetric dcov,
int  icase,
int  nbsimu,
bool  flag_dgm 
)

Conditioning Kriging

Returns
Error return code
Parameters
[in]dbininput Db structure
[in]dboutoutput Db structure
[in]modelModel structure
[in]neighANeigh structure
[in]flag_bayes1 if Bayes option is switched ON
[in]dmeanArray giving the prior means for the drift terms
[in]dcovArray containing the prior covariance matrix for the drift terms
[in]icaseCase for PGS and GRF (or -1)
[in]nbsimuNumber of simulations
[in]flag_dgm1 if the DGM version of kriging should be used
Remarks
: The model contains an anamorphosis with a change of support
: coefficient as soon as flag_dgm is TRUE

◆ anakexp_3D()

int anakexp_3D ( DbGrid db,
double *  cov_ref,
int  cov_radius,
int  neigh_ver,
int  neigh_hor,
int  flag_sym,
Model model,
double  nugget,
int  nfeq,
int  dbg_ix,
int  dbg_iy 
)

Factorial Kriging analysis on a grid file using discretized covariances for the target variable. The discretized covariance of the total variable is calculated on the fly

Returns
Error return code
Parameters
[in]dbinput Db structure
[in]cov_refArray of discretized covariance for target variable
[in]cov_radiusRadius of the covariance array
[in]neigh_verRadius of the Neighborhood along Vertical
[in]neigh_horRadius of the Neighborhood along Horizontal
[in]flag_sym1 for symmetrized covariance
[in]modelModel structure (only used for horizontal)
[in]nuggetAdditional Nugget Effect component
[in]nfeq0 or 1 drift function(s)
[in]dbg_ixRank of the trace along X for variogram debug
[in]dbg_iyRank of the trace along Y for variogram debug
Remarks
The discretized covariance of the target variable is provided
in 1-D along the vertical. Its extension to the space dimension
is performed using the theoretical factorized model
If dbg_ix < -1 || dbg_iy < -1, no variogram debug file is created

◆ anakexp_f()

int anakexp_f ( DbGrid db,
double *  covdd,
double *  covd0,
double  top,
double  bot,
int  cov_radius,
int  neigh_radius,
int  flag_sym,
int  nfeq 
)

Factorial Kriging analysis on a 1-D grid file using discretized covariances for total and partial variables

Returns
Error return code
Parameters
[in]dbinput Db structure
[in]covddArray of discretized cov. for total variable
[in]covd0Array of discretized cov. for partial variable
[in]topElevation of the Top variable
[in]botElevation of the bottom variable
[in]cov_radiusRadius of the Covariance arrays
[in]neigh_radiusRadius of the Neighborhood
[in]flag_sym1 for symmetrized covariance
[in]nfeq0 or 1 drift function(s)

◆ bayes_simulate()

int bayes_simulate ( Model model,
int  nbsimu,
const VectorDouble rmean,
const VectorDouble rcov,
VectorDouble smean 
)

Simulate the drift coefficients from the posterior distributions

Returns
Error returned code
Parameters
[in]modelModel structure
[in]nbsimuNumber of simulation (0 for kriging)
[in]rmeanArray giving the posterior means for the drift terms
[in]rcovArray containing the posterior covariance matrix for the drift terms
[out]smeanArray for simulated posterior mean for the drift means

◆ declustering()

int declustering ( Db dbin,
Model model,
int  method,
ANeigh neigh,
DbGrid dbgrid,
const VectorDouble radius,
const VectorInt ndiscs,
int  flag_sel,
bool  verbose 
)

Perform the Declustering task

Returns
Error return code
Parameters
[in]dbininput Db structure
[in]modelModel structure
[in]methodMethod for declustering
[in]neighANeigh structure
[in]dbgridGrid auxiliary Db structure
[in]radiusArray of neighborhood radius
[in]ndiscsArray of discretization
[in]flag_sel1 to mask off samples with zero weight
[in]verboseVerbose option

◆ global_transitive()

int global_transitive ( DbGrid dbgrid,
Model model,
int  flag_verbose,
int  flag_regular,
int  ndisc,
double *  abundance,
double *  sse,
double *  cvtrans 
)

Estimation of the variance by transitive method

Returns
Error return code
Parameters
[in]dbgridDb structure containing the dicretization grid
[in]modelModel structure
[in]flag_verbose1 for a verbose output
[in]flag_regular1 for regular; 0 for stratified
[in]ndiscNumber of Discretization steps
[out]abundanceGlobal estimated abundance
[out]sseGlobal standard deviation
[out]cvtransCV transitive

◆ inhomogeneous_kriging()

int inhomogeneous_kriging ( Db dbdat,
Db dbsrc,
Db dbout,
double  power,
int  flag_source,
Model model_dat,
Model model_src 
)

Inhomogeneous Kriging with Sources

Returns
Error return code
Parameters
[in]dbdatDb structure containing Data
[in]dbsrcDb structure containing Sources
[in]dboutOutput Db structure
[in]powerPower of the Distance decay
[in]flag_sourceIf the result is the source, rather than diffusion
[in]model_datModel structure for the data
[in]model_srcModel structure for the sources

◆ is_flag_data_disc_defined()

int is_flag_data_disc_defined ( void  )

‍****************************************************************************‍/ *!

‍****************************************************************************‍/ *!

◆ krige_dual_print()

void krige_dual_print ( int  nech,
int  neq,
int  nred,
int *  flag,
double *  dual 
)

Print the Dual matrix

Parameters
[in]nechNumber of active points (optional)
[in]neqNumber of equations
[in]nredReduced number of equations
[in]flagFlag array (optional)
[in]dualKriging Dual matrix

◆ krige_koption_manage()

int krige_koption_manage ( int  mode,
int  flag_check,
const EKrigOpt &  calcul,
int  flag_rand,
const VectorInt ndiscs 
)

Management of Kriging option

Returns
Error return code
Parameters
[in]mode1 for allocation; -1 for deallocation
[in]flag_check1 if the file should be checked
[in]calculType of calculation (EKrigOpt)
[in]flag_rand0 if the second discretization is regular 1 if the second point must be randomized
[in]ndiscsDiscretization parameters (or NULL)
Remarks
This function manages the global structure KOPTION

◆ krige_lhs_print()

void krige_lhs_print ( int  nech,
int  neq,
int  nred,
int *  flagloc,
const double *  lhs 
)

Print the L.H.S. matrix

Parameters
[in]nechNumber of active points (optional)
[in]neqNumber of equations
[in]nredReduced number of equations
[in]flaglocFlag array (optional)
[in]lhsKriging L.H.S

◆ krige_rhs_print()

void krige_rhs_print ( int  nvar,
int  nech,
int  neq,
int  nred,
int *  flag,
double *  rhs 
)

Print the R.H.S. matrix

Parameters
[in]nvarNumber of variables
[in]nechNumber of active points (optional)
[in]neqNumber of equations
[in]nredReduced number of equations
[in]flagFlag array (optional)
[in]rhsKriging R.H.S. matrix

◆ krigsampling_f()

int krigsampling_f ( Db dbin,
Db dbout,
Model model,
double  beta,
VectorInt ranks1,
VectorInt ranks2,
bool  flag_std,
int  verbose 
)

Perform the Kriging procedure using the parcimonious search within the whole input data set

Returns
Error retun code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel structure
[in]betaThresholding value
[in]ranks1Ranks of exact pivots
[in]ranks2Ranks of ACP pivots
[in]flag_stdOption for storing the standard deviation
[in]verboseVerbose flag

◆ krigsum()

int krigsum ( Db dbin,
Db dbout,
Model model,
ANeigh neigh,
bool  flag_positive,
const NamingConvention namconv 
)

Punctual Multivariate Kriging under a constraint

Returns
Error return code
Parameters
[in]dbininput Db structure
[in]dboutoutput Db structure
[in]modelModel structure (univariate)
[in]neighANeigh structure
[in]flag_positive1 for a positive constraints
[in]namconvNaming convention
Remarks
All the variables are estimated using the same model
In this procedure, we assume that:
- the problem is multivariate ("z" variables)
- the constraints is stored in "sum" (only used in dbout)

◆ sampling_f()

int sampling_f ( Db db,
Model model,
double  beta,
int  method1,
int  nsize1_max,
VectorInt ranks1,
int  method2,
int  nsize2_max,
VectorInt ranks2,
int  verbose 
)

Optimize the sampling design

Returns
Error retun code
Parameters
[in]dbDb structure
[in]modelModel structure
[in]betaThresholding value
[in]method1Criterion for choosing exact pivots 1 : Local evaluation 2 : Global evaluation
[in]nsize1_maxMaximum number of exact pivots
[in]ranks1Ranks of exact pivots
[in]method2Criterion for choosing ACP pivots 1 : Local evaluation 2 : Global evaluation
[in]nsize2_maxMaximum number of ACP pivots
[in]ranks2Ranks of ACP pivots
[in]verbose1 for a verbose output

◆ set_DBIN()

void set_DBIN ( Db dbin)

◆ set_DBOUT()

void set_DBOUT ( Db dbout)

◆ st_crit_global()

int st_crit_global ( Db db,
Model model,
VectorInt ranks1,
VectorInt rother,
double *  crit 
)

Evaluate the improvement in adding a new pivot on the global score

Returns
Error retun code
Parameters
[in]dbDb structure
[in]modelModel structure
[in]ranks1Ranks of exact pivots
[in]rotherRanks of the idle samples
[out]critArray of criterion

◆ st_krige_data()

int st_krige_data ( Db db,
Model model,
double  beta,
VectorInt ranks1,
VectorInt ranks2,
VectorInt rother,
int  flag_abs,
double *  data_est,
double *  data_var 
)

Perform the estimation at the data points

Returns
Error return code
Parameters
[in]dbDb structure
[in]modelModel structure
[in]betaThresholding value
[in]ranks1Ranks of exact pivots
[in]ranks2Ranks of ACP pivots
[in]rotherRanks of the idle samples
[in]flag_abs1 Modify 'daata_est' to store the estimation error
[out]data_estArray of estimation at samples
[out]data_varArray of estimation variance at samples