1.4.0
CCC
 
mlayers.cpp File Reference
#include "geoslib_old_f.h"
#include "Calculators/CalcMigrate.hpp"
#include "Variogram/Vario.hpp"
#include "Basic/Utilities.hpp"
#include "Basic/OptDbg.hpp"
#include "Model/Model.hpp"
#include "Neigh/ANeigh.hpp"
#include "Db/Db.hpp"
#include "Db/DbGrid.hpp"
#include "Basic/Memory.hpp"
#include "Core/Keypair.hpp"
#include <math.h>

Classes

struct  LMlayers
 

Functions

static LMlayerslmlayers_free (LMlayers *lmlayers)
 
static int st_get_number_drift (int irf_rank, int flag_ext)
 
static LMlayerslmlayers_alloc (int flag_same, int flag_vel, int flag_cumul, int flag_ext, int flag_z, int colrefd, int colreft, int colrefb, int irf_rank, int match_time, int nlayers)
 
static void lmlayers_print (LMlayers *lmlayers)
 
static int st_locate_sample_in_output (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, int iech, int *igrid)
 
static void st_check_layer (const char *string, LMlayers *lmlayers, int ilayer0)
 
static int st_get_props_result (LMlayers *lmlayers, Db *dbout, int iech, int ilayer0, VectorDouble &props)
 
static int st_get_props_data (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, int iech, int ilayer0, VectorDouble &props)
 
static double st_get_drift_result (LMlayers *lmlayers, Db *dbout, int iech, int ilayer0)
 
static double st_get_drift_data (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, int iech, int ilayer0)
 
static void st_covariance_c00 (LMlayers *lmlayers, Model *model, const VectorDouble &prop1, MatrixSquareGeneral &covtab, double *c00)
 
static double st_cij (LMlayers *lmlayers, Model *model, int ilayer, const VectorDouble &prop1, int jlayer, const VectorDouble &prop2, const double *dd, MatrixSquareGeneral &covtab)
 
static double st_ci0 (LMlayers *lmlayers, Model *model, int ilayer, const VectorDouble &prop1, int jlayer, const double *dd, MatrixSquareGeneral &covtab)
 
static int st_drift (LMlayers *lmlayers, const double *coor, double propval, double drext, int *ipos_loc, VectorDouble &b)
 
static int st_lhs_one (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, Model *model, VectorInt &seltab, int iech0, int ilayer0, double *coor, VectorDouble &prop0, VectorDouble &prop2, MatrixSquareGeneral &covtab, VectorDouble &b)
 
static int st_rhs (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, Model *model, double *coor, VectorInt &seltab, int iechout, int ilayer0, VectorDouble &prop0, VectorDouble &prop2, MatrixSquareGeneral &covtab, VectorDouble &b)
 
static int st_lhs (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, Model *model, VectorInt &seltab, VectorDouble &prop1, VectorDouble &prop2, MatrixSquareGeneral &covtab, double *a, double *acov)
 
static void st_data_vector (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, VectorInt &seltab, VectorDouble &zval)
 
static int st_subtract_optimal_drift (LMlayers *lmlayers, int verbose, Db *dbin, DbGrid *dbout, VectorInt &seltab, VectorDouble &zval)
 
static int st_get_close_sample (LMlayers *lmlayers, Db *dbin, int iech0, const double *coor)
 
static int st_collocated_prepare (LMlayers *lmlayers, int iechout, double *coor, Db *dbin, DbGrid *dbout, Model *model, VectorInt &seltab, double *a, VectorDouble &zval, VectorDouble &prop1, VectorDouble &prop2, MatrixSquareGeneral &covtab, double *b2, VectorDouble &baux, double *ratio)
 
static void st_estimate_regular (LMlayers *lmlayers, int flag_std, double c00, double *a, VectorDouble &b, double *dual, double *wgt, double *estim, double *stdev)
 
static void st_estimate_bayes (LMlayers *lmlayers, int flag_std, double c00, const double *acov, VectorDouble &zval, VectorDouble &b, double *wgt, double *post_mean, double *a0, double *cc, double *ss, const double *gs, double *estim, double *stdev)
 
static void st_estimate (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, Model *model, VectorInt &seltab, int flag_bayes, int flag_std, double *a, VectorDouble &zval, double *dual, VectorDouble &prop1, VectorDouble &prop2, MatrixSquareGeneral &covtab, VectorDouble &b, double *b2, VectorDouble &baux, double *wgt, double *c00, VectorDouble &, double *a0, double *cc, double *ss, double *gs, double *, double *post_mean)
 
static int st_check_auxiliary_variables (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, VectorInt &seltab)
 
static void st_convert_results (LMlayers *lmlayers, Db *dbout, int flag_std)
 
static int st_drift_data (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, VectorInt &seltab, VectorDouble &prop1, VectorDouble &fftab)
 
static int st_drift_bayes (LMlayers *lmlayers, int verbose, double *prior_mean, const double *prior_vars, double *acov, VectorDouble &zval, VectorDouble &fftab, double *a0, double *cc, double *ss, double *gs, double *post_mean, double *post_S)
 
int multilayers_kriging (Db *dbin, DbGrid *dbout, Model *model, ANeigh *neigh, int flag_same, int flag_z, int flag_vel, int flag_cumul, int flag_ext, int flag_std, int flag_bayes, int irf_rank, int match_time, int dim_prior, double *prior_mean, double *prior_vars, int colrefd, int colreft, int colrefb, int verbose)
 
static int st_evaluate_lag (LMlayers *lmlayers, Db *dbin, DbGrid *dbout, Vario_Order *vorder, int nlayers, int ifirst, int ilast, VectorDouble &zval, int *nval, double *distsum, int *stat, VectorDouble &phia, VectorDouble &phib, double *atab, double *btab)
 
static int st_varioexp_chh (LMlayers *lmlayers, int verbose, Db *dbin, DbGrid *dbout, Vario_Order *vorder, VectorDouble &zval, int idir, Vario *vario)
 
int multilayers_vario (Db *dbin, DbGrid *dbout, Vario *vario, int nlayers, int flag_vel, int flag_ext, int irf_rank, int match_time, int colrefd, int colreft, int verbose)
 
static int st_get_prior (int nech, int npar, VectorDouble &zval, VectorDouble &fftab, double *mean, double *vars)
 
int multilayers_get_prior (Db *dbin, DbGrid *dbout, Model *model, int flag_same, int flag_vel, int flag_ext, int irf_rank, int match_time, int colrefd, int colreft, int colrefb, int verbose, int *npar_arg, double **mean, double **vars)
 

Function Documentation

◆ lmlayers_alloc()

static LMlayers* lmlayers_alloc ( int  flag_same,
int  flag_vel,
int  flag_cumul,
int  flag_ext,
int  flag_z,
int  colrefd,
int  colreft,
int  colrefb,
int  irf_rank,
int  match_time,
int  nlayers 
)
static

Allocate the Multi-Layers internal structure

Returns
Pointer to the allocated structure
Parameters
[in]flag_same1 if input and output files coincide
[in]flag_vel1 if work is performed in Velocity, 0 for Depth
[in]flag_cumul1 if work must be done on Depth; 0 on Thickness
[in]flag_ext1 if external drift must be used; 0 otherwise
[in]flag_z1 if the output must be provided in depth
[in]colrefdRank of the reference depth variable
[in]colreftRank of the reference Time variable in Dbout
[in]colrefbRank of the reference Bottom Depth variable
[in]irf_rankRank of the Intrinsic Random Function (0 or 1)
[in]match_timePointer to the Time pointer (1 if defined as ELoc::F or 0 for ELoc::TIME)
[in]nlayersNumber of layers

◆ lmlayers_free()

static LMlayers* lmlayers_free ( LMlayers lmlayers)
static

Free the Multi-Layers internal structure

Returns
Pointer to the freed structure
Parameters
[in]lmlayersPointer to the LMlayers structure to be freed

◆ lmlayers_print()

static void lmlayers_print ( LMlayers lmlayers)
static

Print the Multi-Layers internal structure

Parameters
[in]lmlayersPointer to the LMlayers structure

◆ multilayers_get_prior()

int multilayers_get_prior ( Db dbin,
DbGrid dbout,
Model model,
int  flag_same,
int  flag_vel,
int  flag_ext,
int  irf_rank,
int  match_time,
int  colrefd,
int  colreft,
int  colrefb,
int  verbose,
int *  npar_arg,
double **  mean,
double **  vars 
)

Multi-layers get the mean and prior matrices for Bayesian prior

Returns
Error return code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel structure
[in]flag_same1 if input and output files coincide
[in]flag_vel1 if work is performed in Velocity, 0 for Depth
[in]flag_ext1 if external drift must be used; 0 otherwise
[in]irf_rankRank of the Intrinsic Random Function (0 or 1)
[in]match_time1 if external drift matches time; 0 otherwise
[in]colrefdRank of the reference Depth variable in Dbout
[in]colreftRank of the reference Time variable in Dbout
[in]colrefbRank of the Bottom Depth variable in Dbout (or -1)
[in]verboseVerbose option
[out]npar_argNumber of drift terms
[out]meanArray of means
[out]varsArray of variances

◆ multilayers_kriging()

int multilayers_kriging ( Db dbin,
DbGrid dbout,
Model model,
ANeigh neigh,
int  flag_same,
int  flag_z,
int  flag_vel,
int  flag_cumul,
int  flag_ext,
int  flag_std,
int  flag_bayes,
int  irf_rank,
int  match_time,
int  dim_prior,
double *  prior_mean,
double *  prior_vars,
int  colrefd,
int  colreft,
int  colrefb,
int  verbose 
)

Multi-layers architecture estimation

Returns
Error return code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel structure
[in]neighANeigh structure
[in]flag_same1 if input and output files coincide
[in]flag_z1 if the output must be converted back into depth
[in]flag_vel1 if work is performed in Velocity, 0 for Depth
[in]flag_cumul1 if work is performed in Depth; 0 in Thickness
[in]flag_ext1 if external drift must be used; 0 otherwise
[in]flag_std1 if the estimation error must be calculated
[in]flag_bayes1 if the Bayesian hypothesis is used on drift coeffs
[in]irf_rankRank of the Intrinsic Random Function (0 or 1)
[in]match_time1 if external drift matches time; 0 otherwise
[in]dim_priorDimension of the prior information (for verification)
[in]prior_meanVector of prior means for drift coefficients
[in]prior_varsVector of prior variances for drift coefficients
[in]colrefdRank of the reference Depth variable in Dbout
[in]colreftRank of the reference Time variable in Dbout
[in]colrefbRank of the Bottom Depth variable in Dbout (or -1)
[in]verboseVerbose option

◆ multilayers_vario()

int multilayers_vario ( Db dbin,
DbGrid dbout,
Vario vario,
int  nlayers,
int  flag_vel,
int  flag_ext,
int  irf_rank,
int  match_time,
int  colrefd,
int  colreft,
int  verbose 
)

Multi-layers architecture experimental variogram

Returns
Error return code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]varioVario structure
[in]nlayersNumber of layers
[in]flag_vel1 if work is performed in Velocity, 0 for Depth
[in]flag_ext1 if external drift must be used; 0 otherwise
[in]irf_rankRank of the Intrinsic Random Function (0 or 1)
[in]match_time1 if external drift matches time; 0 otherwise
[in]colrefdRank of the reference Depth variable in Dbout
[in]colreftRank of the reference Time variable in Dbout
[in]verbose1 for a verbose option

◆ st_check_auxiliary_variables()

static int st_check_auxiliary_variables ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
VectorInt seltab 
)
static

Check all the auxiliary variables

Returns
The number of active samples
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in,out]seltabNumber of sample definition (0, 1 or 2)

◆ st_check_layer()

static void st_check_layer ( const char *  string,
LMlayers lmlayers,
int  ilayer0 
)
static

Check if the target layer rank is consistent

Parameters
[in]stringName of the calling procedure
[in]lmlayersLMlayers structure
[in]ilayer0Rank of the target layer (starting from 1)
Remarks
If this target layer rank is not correct, mes_abort() is called
and the program is interrupted as this must never happen.

◆ st_ci0()

static double st_ci0 ( LMlayers lmlayers,
Model model,
int  ilayer,
const VectorDouble prop1,
int  jlayer,
const double *  dd,
MatrixSquareGeneral covtab 
)
static

Calculate the covariance between data and target

Returns
The covariance terms or TEST
Parameters
[in]lmlayersLMlayers structure
[in]modelModel
[in]ilayerLayer of interest (data point). Starting from 1
[in]prop1Working array at data point (Dimension: nlayers)
[in]jlayerLayer of interest (target point). Starting from 1
[in]ddDistance vector (or NULL for zero-distance)
[out]covtabWorking array (Dimension = nlayers * nlayers)
Remarks
: As this function may return TEST, TEST value should be tested

◆ st_cij()

static double st_cij ( LMlayers lmlayers,
Model model,
int  ilayer,
const VectorDouble prop1,
int  jlayer,
const VectorDouble prop2,
const double *  dd,
MatrixSquareGeneral covtab 
)
static

Calculate the covariance between data and data

Returns
The covariance terms or TEST
Parameters
[in]lmlayersPointer to the LMlayers structure to be freed
[in]modelModel
[in]ilayerLayer of interest (first point). Starting from 1
[in]prop1Working array at first point (Dimension: nlayers)
[in]jlayerLayer of interest (second point). Starting from 1
[in]prop2Working array at second point (Dimension: nlayers)
[in]ddDistance vector (or NULL for zero-distance)
[out]covtabWorking array (Dimension = nlayers * nlayers)
Remarks
: As this function may return TEST, TEST value should be tested

◆ st_collocated_prepare()

static int st_collocated_prepare ( LMlayers lmlayers,
int  iechout,
double *  coor,
Db dbin,
DbGrid dbout,
Model model,
VectorInt seltab,
double *  a,
VectorDouble zval,
VectorDouble prop1,
VectorDouble prop2,
MatrixSquareGeneral covtab,
double *  b2,
VectorDouble baux,
double *  ratio 
)
static

Perform the per-calculation for estimation with collocated option

Parameters
[in]lmlayersLMlayers structure
[in]iechoutRank of the target
[in]coorCoordinates of the target
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel
[in]seltabNumber of sample definition (0, 1 or 2)
[in]aL.H.S. (square) inverted matrix
[in]zvalData vector (extended)
[out]prop1Working array (Dimension: nlayers)
[out]prop2Working array (Dimension: nlayers)
[out]covtabWorking array (Dimension = nlayers * nlayers)
[out]b2Working vector (Dimension = neq)
[out]bauxWorking vector (Dimension = neq)
[out]ratioRatio value

◆ st_convert_results()

static void st_convert_results ( LMlayers lmlayers,
Db dbout,
int  flag_std 
)
static

Convert the results in the Depth scale

Parameters
[in]lmlayersLMlayers structure
[in]dboutOutput Db structure
[in]flag_std1 if the estimation error must be calculated
Remarks
The conversion is performed:
- if the calculations have been performed in Velocity or Thickness
- if the calculations have been performed in cumulative or not
The standard deviation is also transformed

◆ st_covariance_c00()

static void st_covariance_c00 ( LMlayers lmlayers,
Model model,
const VectorDouble prop1,
MatrixSquareGeneral covtab,
double *  c00 
)
static

Calculate the array of covariances for zero distance

Parameters
[in]lmlayersPointer to the LMlayers structure to be freed
[in]modelModel
[in]prop1Working array at first point (Dimension: nlayers)
[out]covtabWorking array (Dimension = nlayers * nlayers)
[out]c00Returned array (Dimension = nlayers)
Remarks
: This array depends on the target location through proportions

◆ st_data_vector()

static void st_data_vector ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
VectorInt seltab,
VectorDouble zval 
)
static

Establish the vector of data

Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]seltabNumber of sample definition (0, 1 or 2)
[out]zvalThe data vector (Dimension: neq)

◆ st_drift()

static int st_drift ( LMlayers lmlayers,
const double *  coor,
double  propval,
double  drext,
int *  ipos_loc,
VectorDouble b 
)
static

Calculate the drift terms

Returns
Error return code
Parameters
[in]lmlayersPointer to the LMlayers structure to be freed
[in]coorArray of coordinates
[in]propvalValue for the proportion (used if flag_cumul=TRUE)
[in]drextValue of the external drift
[in,out]ipos_locAddress for the first drift term. On output, address for the next term after the drift
[out]bArray for storing the drift

◆ st_drift_bayes()

static int st_drift_bayes ( LMlayers lmlayers,
int  verbose,
double *  prior_mean,
const double *  prior_vars,
double *  acov,
VectorDouble zval,
VectorDouble fftab,
double *  a0,
double *  cc,
double *  ss,
double *  gs,
double *  post_mean,
double *  post_S 
)
static

Calculate the posterior mean and variances in Bayesian case

Returns
Error return code
Remarks
At the end of this function,
invS contains the inverse of prior_S
post_S contains the inverse of post_S

◆ st_drift_data()

static int st_drift_data ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
VectorInt seltab,
VectorDouble prop1,
VectorDouble fftab 
)
static

Fill the array of drift values at data points

Returns
1 Error return code
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]seltabNumber of sample definition (0, 1 or 2)
[in]prop1Working array (Dimension: nlayers)
[out]fftabDrift array (Dimension: npar[nrow] * nech[ncol])

◆ st_estimate()

static void st_estimate ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
Model model,
VectorInt seltab,
int  flag_bayes,
int  flag_std,
double *  a,
VectorDouble zval,
double *  dual,
VectorDouble prop1,
VectorDouble prop2,
MatrixSquareGeneral covtab,
VectorDouble b,
double *  b2,
VectorDouble baux,
double *  wgt,
double *  c00,
VectorDouble ,
double *  a0,
double *  cc,
double *  ss,
double *  gs,
double *  ,
double *  post_mean 
)
static

Perform the estimation at the grid nodes

Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel
[in]seltabNumber of sample definition (0, 1 or 2)
[in]flag_bayes1 if the Bayesian hypothesis is used on drift coeffs
[in]flag_std1 if the estimation error must be calculated
[in]aL.H.S. (square) inverted matrix
[in]zvalData vector (extended)
[in]dualDual vector
[out]prop1Working array (Dimension: nlayers)
[out]prop2Working array (Dimension: nlayers)
[out]covtabWorking array (Dimension = nlayers * nlayers)
[out]bWorking vector (Dimension = neq)
[out]b2Working vector (Dimension = neq)
[out]bauxWorking vector (Dimension = neq)
[out]wgtWorking array (Dimension = neq)
[out]c00Working array (Dimension = nlayers)
[out]a0Constant term
[out]ccOutput value
[out]ssOutput value
[out]gsOutput value
[out]post_meanArray of posterior mean

◆ st_estimate_bayes()

static void st_estimate_bayes ( LMlayers lmlayers,
int  flag_std,
double  c00,
const double *  acov,
VectorDouble zval,
VectorDouble b,
double *  wgt,
double *  post_mean,
double *  a0,
double *  cc,
double *  ss,
const double *  gs,
double *  estim,
double *  stdev 
)
static

Perform the estimation at the grid nodes in the bayesian case

Parameters
[in]lmlayersLMlayers structure
[in]flag_std1 if the estimation error must be calculated
[in]c00Variance at target
[in]acovL.H.S. (square) inverted matrix
[in]zvalData vector
[in]bWorking vector (Dimension = neq)
[in]wgtWorking array (Dimension = neq)
[out]post_meanArray of posterior mean
[out]a0Constant term
[out]ccOutput value
[out]ssOutput value
[out]gsOutput value
[out]estimEstimated value
[out]stdevStandard deviation of estimation error

◆ st_estimate_regular()

static void st_estimate_regular ( LMlayers lmlayers,
int  flag_std,
double  c00,
double *  a,
VectorDouble b,
double *  dual,
double *  wgt,
double *  estim,
double *  stdev 
)
static

Perform the estimation at the grid nodes in regular case

Parameters
[in]lmlayersLMlayers structure
[in]flag_std1 if the estimation error must be calculated
[in]c00Variance for target
[in]aL.H.S. (square) inverted matrix
[in]bWorking vector (Dimension = neq)
[in]dualDual vector
[in]wgtWorking array (Dimension = neq)
[out]estimEstimated value
[out]stdevStandard deviation of estimation error

◆ st_evaluate_lag()

static int st_evaluate_lag ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
Vario_Order vorder,
int  nlayers,
int  ifirst,
int  ilast,
VectorDouble zval,
int *  nval,
double *  distsum,
int *  stat,
VectorDouble phia,
VectorDouble phib,
double *  atab,
double *  btab 
)
static

Evaluate the sill matrix for a given lag

Returns
Error return code (proportions not calculatable)
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]vorderVario_Order structure
[in]nlayersNumber of layers
[in]ifirstIndex of the first pair (included)
[in]ilastIndex of the last pair (excluded)
[in]zvalArray containing the sample values
[out]nvalNumber of relevant pairs
[out]distsumAverage distance
[out]statWorking array (Dimension: nlayers * nlayers)
[out]phiaWorking array for proportions (Dimension: nlayers)
[out]phibWorking array for proportions (Dimension: nlayers)
[out]atabWorking array (Dimension: nhalf * nhalf)
[out]btabWorking array (Dimension: nhalf)

◆ st_get_close_sample()

static int st_get_close_sample ( LMlayers lmlayers,
Db dbin,
int  iech0,
const double *  coor 
)
static

Check if an intercept with the bottom layer is located close enough to the current sample

Returns
1 if a duplicate must be generated; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]iech0Rank of sample to be discarded (or -1)
[in]coorCoordinates of the target

◆ st_get_drift_data()

static double st_get_drift_data ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
int  iech,
int  ilayer0 
)
static

Return the external drift value at an input location for a target layer

Returns
The external drift value of TEST
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]iechRank of the target sample (in output Db)
[in]ilayer0Rank of the target layer (Starting from 1)

◆ st_get_drift_result()

static double st_get_drift_result ( LMlayers lmlayers,
Db dbout,
int  iech,
int  ilayer0 
)
static

Return the external drift value at an output location for a target layer

Returns
The external drift value of TEST
Parameters
[in]lmlayersLMlayers structure
[in]dboutOutput Db structure
[in]iechRank of the target sample (in output Db)
[in]ilayer0Rank of the target layer (Starting from 1)

◆ st_get_number_drift()

static int st_get_number_drift ( int  irf_rank,
int  flag_ext 
)
static

Returns the number of drift functions

Returns
Number of drift conditions
Parameters
[in]irf_rankRank of the Intrinsic Random Function (0 or 1)
[in]flag_ext1 if external drift must be used; 0 otherwise

◆ st_get_prior()

static int st_get_prior ( int  nech,
int  npar,
VectorDouble zval,
VectorDouble fftab,
double *  mean,
double *  vars 
)
static

Determine the mean and variance of drift coefficients

Returns
Error return code
Parameters
[in]nechNumber of samples
[in]nparNumber of drift coefficients
[in]zvalThe data vector (Dimension: neq)
[in]fftabDrift array (Dimension: npar[nrow] * nech[ncol])
[out]meanArray of means
[out]varsArray of variances

◆ st_get_props_data()

static int st_get_props_data ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
int  iech,
int  ilayer0,
VectorDouble props 
)
static

Fill the proportion vector at a data location, up to the target layer

Returns
1 if the proportion vector cannot be defined; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]iechRank of the target sample
[in]ilayer0Rank of the target layer (starting from 1)
[out]propsWorking array (Dimension: nlayers)

◆ st_get_props_result()

static int st_get_props_result ( LMlayers lmlayers,
Db dbout,
int  iech,
int  ilayer0,
VectorDouble props 
)
static

Fill the proportion vector at a output location, up to the target layer

Returns
1 if the proportion vector cannot be defined; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dboutOutput Db structure
[in]iechRank of the target sample (in output Db)
[in]ilayer0Rank of the target layer (starting from 1)
[out]propsWorking array (Dimension: nlayers)

◆ st_lhs()

static int st_lhs ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
Model model,
VectorInt seltab,
VectorDouble prop1,
VectorDouble prop2,
MatrixSquareGeneral covtab,
double *  a,
double *  acov 
)
static

Calculates the Kriging L.H.S.

Returns
: 1 if the L.H.S. vector cannot be calculated; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel
[in]seltabNumber of sample definition (0, 1 or 2)
[out]prop1Working array (Dimension: nlayers)
[out]prop2Working array (Dimension: nlayers)
[out]covtabWorking array (Dimension = nlayers * nlayers)
[out]aL.H.S. (square) matrix
[out]acovL.H.S. (square) covariance matrix

◆ st_lhs_one()

static int st_lhs_one ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
Model model,
VectorInt seltab,
int  iech0,
int  ilayer0,
double *  coor,
VectorDouble prop0,
VectorDouble prop2,
MatrixSquareGeneral covtab,
VectorDouble b 
)
static

Calculates the L.H.S. for one data

Returns
: 1 if the L.H.S. vector cannot be calculated; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel
[in]seltabNumber of sample definition (0, 1 or 2)
[in]iech0Rank of the target sample (used for ext. drift)
[in]ilayer0Rank of the layer of interest (Starting from 1)
[in]coorCoordinates of the data
[in]prop0Working array at first point (Dimension: nlayers)
[out]prop2Working array (Dimension: nlayers)
[out]covtabWorking array (Dimension = nlayers * nlayers)
[out]bR.H.S. vector (Dimension = neq)

◆ st_locate_sample_in_output()

static int st_locate_sample_in_output ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
int  iech,
int *  igrid 
)
static

Returns the absolute grid node absolute index which is the closest to a given sample of a Db In the case of same input and output file, simply return 'iech'

Returns
1 if the sample does not belong to the grid; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]iechRank in the input Db
[out]igridRank of the node in the output Db

◆ st_rhs()

static int st_rhs ( LMlayers lmlayers,
Db dbin,
DbGrid dbout,
Model model,
double *  coor,
VectorInt seltab,
int  iechout,
int  ilayer0,
VectorDouble prop0,
VectorDouble prop2,
MatrixSquareGeneral covtab,
VectorDouble b 
)
static

Calculates the Kriging R.H.S.

Returns
: 1 if the R.H.S. has not been calculated; 0 otherwise
Parameters
[in]lmlayersLMlayers structure
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel
[in]coorCoordinates of the target sample
[in]seltabNumber of sample definition (0, 1 or 2)
[in]iechoutRank of the target sample
[in]ilayer0Rank of the layer of interest (Starting from 1)
[out]prop0Working array (Dimension: nlayers)
[out]prop2Working array (Dimension: nlayers)
[out]covtabWorking array (Dimension = nlayers * nlayers)
[out]bR.H.S. vector (Dimension = neq)

◆ st_subtract_optimal_drift()

static int st_subtract_optimal_drift ( LMlayers lmlayers,
int  verbose,
Db dbin,
DbGrid dbout,
VectorInt seltab,
VectorDouble zval 
)
static

Calculate the Drift and subtract it from the Data

Parameters
[in]lmlayersLMlayers structure
[in]verbose1 for a verbose option
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]seltabNumber of sample definition (0, 1 or 2)
[out]zvalThe data vector (Dimension: neq)
Remarks
In the current version, the optimal coefficients of the Drift
are output using the set_keypair mechanism using the keyword:
"Optim_Drift_Coeffs" which returns 'ipos' values

◆ st_varioexp_chh()

static int st_varioexp_chh ( LMlayers lmlayers,
int  verbose,
Db dbin,
DbGrid dbout,
Vario_Order vorder,
VectorDouble zval,
int  idir,
Vario vario 
)
static

Calculate the variogram for each lag per direction

Returns
Error return code
Parameters
[in]lmlayersLMlayers structure
[in]verbose1 for a verbose option
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]vorderVario_Order structure
[in]zvalData vector
[in]idirRank of the Direction
[out]varioVario structure