1.5.1
CCC
 
model.cpp File Reference

Functions

Modelmodel_duplicate_for_gradient (const Model *model, double ball_radius)
 
void model_covupdt (Model *model, const double *c0, int flag_verbose, int *flag_nugget, double *nugget)
 
void model_cova_characteristics (const ECov &type, char cov_name[STRING_LENGTH], int *flag_range, int *flag_param, int *min_order, int *max_ndim, int *flag_int_1d, int *flag_int_2d, int *flag_aniso, int *flag_rotation, double *scale, double *parmax)
 
Modelmodel_combine (const Model *model1, const Model *model2, double r)
 
int model_covmat_inchol (int verbose, Db *db, Model *model, double eta, int npivot_max, int nsize1, const int *ranks1, const double *center, int flag_sort, int *npivot_arg, int **Pret, double **Gret, const CovCalcMode *mode)
 

Variables

int NDIM_LOCAL = 0
 
VectorDouble X1_LOCAL = VectorDouble()
 
VectorDouble X2_LOCAL = VectorDouble()
 

Function Documentation

◆ model_combine()

Model* model_combine ( const Model model1,
const Model model2,
double  r 
)

Combine two monovariate models into a bivariate model (residuals model)

Returns
Pointer to the newly created Model structure
Parameters
[in]model1First input Model
[in]model2Second input Model
[in]rCorrelation coefficient
Remarks
: The drift is not copied into the new model
: It has been extended to the case where only one model is defined

◆ model_cova_characteristics()

void model_cova_characteristics ( const ECov &  type,
char  cov_name[STRING_LENGTH],
int *  flag_range,
int *  flag_param,
int *  min_order,
int *  max_ndim,
int *  flag_int_1d,
int *  flag_int_2d,
int *  flag_aniso,
int *  flag_rotation,
double *  scale,
double *  parmax 
)

Returns the characteristics of the covariance

Parameters
[in]typeType of the covariance
[out]cov_nameName of the covariance
[out]flag_rangerange definition
  • +1 if the range is defined
  • -1 if the range is redundant with the sill
[out]flag_param1 if the third parameter is defined
[out]min_orderMinimum IRF order for validity
[out]max_ndimMaximum dimension for validity
[out]flag_int_1dIntegral range in 1-D
[out]flag_int_2dIntegral range in 2-D
[out]flag_aniso1 if anisotropy is meaningful
[out]flag_rotation1 if an anisotropy rotation is meaningful
[out]scaleScaling parameter
[out]parmaxMaximum value for the third parameter

◆ model_covmat_inchol()

int model_covmat_inchol ( int  verbose,
Db db,
Model model,
double  eta,
int  npivot_max,
int  nsize1,
const int *  ranks1,
const double *  center,
int  flag_sort,
int *  npivot_arg,
int **  Pret,
double **  Gret,
const CovCalcMode mode 
)

Establish and invert a covariance matrix using Incomplete Cholesky method

Parameters
[in]verboseVerbose option
[in]dbDb structure
[in]modelModel structure
[in]npivot_maxMaximum number of pivots (or 0)
[in]etaPrecision (or TEST)
[in]nsize1Number of pivots already selected
[in]ranks1Ranks of pivots already selected
[in]centerOptional Centering point (for increments)
[in]flag_sortReordering flag (see remarks)
[in]modeCovCalcMode structure
[out]npivot_argNumber of pivots
[out]PretArray of indices of the retained samples (from 1) Dimension: nech
[out]GretRectangular matrix Dimension: nech * npivot_arg
Remarks
The output arrays Pret and Gret should be freed by calling function
The array G contains as many lines as there are samples
If flag_sort = FALSE, the first lines concentrate on pivots,
and the other points are located afterwards
If flag_sort = TRUE, the lines are sorted in the same order as the
initial set of samples
The incomplete Cholsky algorithm stops when either the next pivot
value is below 'eta' or when maximum number of pivots 'npivot_max'
has been reached
If the center point is provided in 'center', the calculations
of covariance of increments are calculated instead. Then 'center'
must provide the coordinates of the origin point.

◆ model_covupdt()

void model_covupdt ( Model model,
const double *  c0,
int  flag_verbose,
int *  flag_nugget,
double *  nugget 
)

Update the model for fitting Covariance or Covariogram

Parameters
[in]modelModel structure
[in]c0Array of variance values at the origin
[in]flag_verbose1 for verbose output
[out]flag_nugget1 if a nugget component must be added
[out]nuggetArray of sills for the nugget component

TODO : dead code ?

◆ model_duplicate_for_gradient()

Model* model_duplicate_for_gradient ( const Model model,
double  ball_radius 
)

Duplicates a Model from another Model for Gradients

Returns
The modified Model structure
Parameters
[in]modelInput Model
[in]ball_radiusRadius for Gradient calculation

Variable Documentation

◆ NDIM_LOCAL

int NDIM_LOCAL = 0

◆ X1_LOCAL

VectorDouble X1_LOCAL = VectorDouble()

◆ X2_LOCAL

VectorDouble X2_LOCAL = VectorDouble()