1.5.1
CCC
 

Functions

MatrixRectangular Model::evalCovMatrix (Db *db1, Db *db2=nullptr, int ivar0=-1, int jvar0=-1, const VectorInt &nbgh1=VectorInt(), const VectorInt &nbgh2=VectorInt(), const CovCalcMode *mode=nullptr)
 
MatrixSquareSymmetric Model::evalCovMatrixSymmetric (Db *db1, int ivar0=-1, const VectorInt &nbgh1=VectorInt(), const CovCalcMode *mode=nullptr)
 
MatrixSparseModel::evalCovMatrixSparse (Db *db1, Db *db2=nullptr, int ivar0=-1, int jvar0=-1, const VectorInt &nbgh1=VectorInt(), const VectorInt &nbgh2=VectorInt(), const CovCalcMode *mode=nullptr, double eps=EPSILON3)
 
VectorDouble Model::evalCovMatrixV (Db *db1, Db *db2=nullptr, int ivar0=-1, int jvar0=-1, const VectorInt &nbgh1=VectorInt(), const VectorInt &nbgh2=VectorInt(), const CovCalcMode *mode=nullptr)
 
MatrixRectangular Model::evalCovMatrixOptim (const Db *db1, const Db *db2=nullptr, int ivar0=-1, int jvar0=-1, const VectorInt &nbgh1=VectorInt(), const VectorInt &nbgh2=VectorInt(), const CovCalcMode *mode=nullptr)
 
MatrixSquareSymmetric Model::evalCovMatrixSymmetricOptim (const Db *db1, int ivar0=-1, const VectorInt &nbgh1=VectorInt(), const CovCalcMode *mode=nullptr)
 
double Model::extensionVariance (const Db *db, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), int ivar=0, int jvar=0)
 
double Model::samplingDensityVariance (const Db *db, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), int ivar=0, int jvar=0) const
 
double Model::specificVolume (const Db *db, double mean, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), int ivar=0, int jvar=0) const
 
double Model::coefficientOfVariation (const Db *db, double volume, double mean, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), int ivar=0, int jvar=0) const
 
double Model::specificVolumeFromCoV (Db *db, double cov, double mean, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), int ivar=0, int jvar=0) const
 
void Model::evalZAndGradients (const SpacePoint &p1, const SpacePoint &p2, double &covVal, VectorDouble &covGp, VectorDouble &covGG, const CovCalcMode *mode=nullptr, bool flagGrad=false) const
 
void Model::evalZAndGradients (const VectorDouble &vec, double &covVal, VectorDouble &covGp, VectorDouble &covGG, const CovCalcMode *mode=nullptr, bool flagGrad=false) const
 
double Model::evalCov (const VectorDouble &incr, int icov=0, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
void Model::setSill (int icov, int ivar, int jvar, double value)
 
void Model::setRangeIsotropic (int icov, double range)
 
void Model::setMarkovCoeffs (int icov, const VectorDouble &coeffs)
 
void Model::setCovaFiltered (int icov, bool filtered)
 
void Model::setActiveFactor (int iclass)
 
int Model::getActiveFactor () const
 
int Model::getAnamNClass () const
 
const DriftListModel::getDriftList () const
 TODO : to be removed (encapsulation of DriftList) More...
 
const ADriftModel::getDrift (int il) const
 
int Model::getDriftNumber () const
 
int Model::getExternalDriftNumber () const
 
int Model::getRankFext (int il) const
 
int Model::getDriftEquationNumber () const
 
bool Model::isDriftFiltered (unsigned int il) const
 
int Model::getDriftMaxIRFOrder (void) const
 
bool Model::isDriftDefined (const VectorInt &powers, int rank_fex=0) const
 
bool Model::isDriftDifferentDefined (const VectorInt &powers, int rank_fex=-1) const
 
bool Model::isDriftSampleDefined (const Db *db, int ib, int nech, const VectorInt &nbgh, const ELoc &loctype) const
 
void Model::setDriftFiltered (int il, bool filtered)
 
VectorVectorDouble Model::getDrifts (const Db *db, bool useSel=true)
 
void Model::setBetaHat (const VectorDouble &betaHat)
 
double Model::evalDrift (const Db *db, int iech, int il, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
double Model::evalDriftValue (const Db *db, int iech, int ivar, int ib, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
VectorDouble Model::evalDriftBySample (const Db *db, int iech, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
void Model::evalDriftBySampleInPlace (const Db *db, int iech, const ECalcMember &member, VectorDouble &drftab) const
 
MatrixRectangular Model::evalDriftMatrix (const Db *db, int ivar0=-1, const VectorInt &nbgh=VectorInt(), const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
double Model::evalDriftVarCoef (const Db *db, int iech, int ivar, const VectorDouble &coeffs) const
 
VectorDouble Model::evalDriftVarCoefs (const Db *db, const VectorDouble &coeffs, int ivar=0, bool useSel=false) const
 
const CovContextModel::getContext () const
 TODO : to be removed (encapsulation of Context) More...
 
const ASpaceModel::getASpace () const
 
const VectorDoubleModel::getMeans () const
 
double Model::getMean (int ivar) const
 
const VectorDoubleModel::getCovar0s () const
 
double Model::getCovar0 (int ivar, int jvar) const
 
double Model::getField () const
 
int Model::getDimensionNumber () const
 
void Model::setMeans (const VectorDouble &mean)
 
void Model::setMean (double mean, int ivar=0)
 
void Model::setCovar0s (const VectorDouble &covar0)
 
void Model::setCovar0 (int ivar, int jvar, double covar0)
 
void Model::setField (double field)
 
const EModelProperty & Model::getCovMode () const
 
ModelModel::duplicate () const
 
ModelModel::createReduce (const VectorInt &validVars) const
 
int Model::getVariableNumber () const
 
int Model::hasExternalCov () const
 
VectorDouble Model::sampleUnitary (const VectorDouble &hh, int ivar=0, int jvar=0, VectorDouble codir=VectorDouble(), const CovCalcMode *mode=nullptr)
 
VectorDouble Model::envelop (const VectorDouble &hh, int ivar=0, int jvar=0, int isign=1, VectorDouble codir=VectorDouble(), const CovCalcMode *mode=nullptr)
 
int Model::fitFromCovIndices (Vario *vario, const VectorECov &types=ECov::fromKeys({"EXPONENTIAL"}), const Constraints &constraints=Constraints(), const Option_VarioFit &optvar=Option_VarioFit(), const Option_AutoFit &mauto=Option_AutoFit(), bool verbose=false)
 
int Model::fit (Vario *vario, const VectorECov &types=ECov::fromKeys({"SPHERICAL"}), const Constraints &constraints=Constraints(), const Option_VarioFit &optvar=Option_VarioFit(), const Option_AutoFit &mauto=Option_AutoFit(), bool verbose=false)
 
int Model::fitFromVMap (DbGrid *dbmap, const VectorECov &types=ECov::fromKeys({"SPHERICAL"}), const Constraints &constraints=Constraints(), const Option_VarioFit &optvar=Option_VarioFit(), const Option_AutoFit &mauto=Option_AutoFit(), bool verbose=false)
 
int Model::buildVmapOnDbGrid (DbGrid *dbgrid, const NamingConvention &namconv=NamingConvention("VMAP")) const
 
int Model::stabilize (double percent, bool verbose=false)
 
int Model::standardize (bool verbose=false)
 
double Model::gofToVario (const Vario *vario, bool verbose=true)
 
static void Model::gofDisplay (double gof, bool byValue=true, const VectorDouble &thresholds={2., 5., 10., 100})
 
static VectorECov Model::initCovList (const VectorInt &covranks)
 
bool Model::isValid () const
 
VectorDouble Model::sample (const VectorDouble &h, const VectorDouble &codir=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr, const CovInternal *covint=nullptr)
 
double Model::evaluateOneIncr (double hh, const VectorDouble &codir=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
void Model::evaluateMatInPlace (const CovInternal *covint, const VectorDouble &d1, MatrixSquareGeneral &covtab, bool flag_init=false, double weight=1., const CovCalcMode *mode=nullptr)
 
double Model::evaluateOneGeneric (const CovInternal *covint, const VectorDouble &d1=VectorDouble(), double weight=1., const CovCalcMode *mode=nullptr)
 
VectorDouble Model::evaluateFromDb (Db *db, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
double Model::calculateStdev (Db *db1, int iech1, Db *db2, int iech2, bool verbose=false, double factor=1., const CovCalcMode *mode=nullptr)
 
double Model::computeLogLikelihood (Db *db, bool verbose=false)
 

Detailed Description

These functions are meant to calculate the covariance Matrix between two Dbs or between a Db and itself. They take into account the presence of a possible selection They also account for heterotopy (if Z-variables are defined in the Db(s)

Parameters
db1First Db
db2(Optional second Db)
ivar0Rank of the selected variable in the first Db (-1 for all variables)
jvar0Rank of the selected variable in the second Db (-1 for all variables)
nbgh1Vector of indices of active samples in first Db (optional)
nbgh2Vector of indices of active samples in second Db (optional)
modeCovCalcMode structure
Remarks
The returned matrix if dimension to nrows * ncols where
each term is the product of the number of active samples
by the number of samples where the variable is defined
Note
'dbin' and 'dbout' cannot be made 'const' as they can be updated
due to the presence of 'nostat'
Returns
A Matrix either in Dense or Sparse format

Function Documentation

◆ buildVmapOnDbGrid()

int Model::buildVmapOnDbGrid ( DbGrid dbgrid,
const NamingConvention namconv = NamingConvention("VMAP") 
) const

Calculate the variogram map from a Model (presented as Variogram, not Covariance)

Returns
Error return code
Parameters
[in]dbgridGrid structure
[in]namconvNaming convention

◆ calculateStdev()

double Model::calculateStdev ( Db db1,
int  iech1,
Db db2,
int  iech2,
bool  verbose = false,
double  factor = 1.,
const CovCalcMode mode = nullptr 
)

Returns the standard deviation at a given increment for a given model between two samples of two Dbs

Parameters
[in]db1First Db
[in]iech1Rank in the first Db
[in]db2Second Db
[in]iech2Rank in the second Db
[in]verboseVerbose flag
[in]factorMultiplicative factor for standard deviation
[in]modeCovCalcMode structure

◆ coefficientOfVariation()

double Model::coefficientOfVariation ( const Db db,
double  volume,
double  mean,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
int  ivar = 0,
int  jvar = 0 
) const
inline

◆ computeLogLikelihood()

double Model::computeLogLikelihood ( Db db,
bool  verbose = false 
)

Compute the log-likelihood (based on covariance)

Parameters
dbDb structure where variable are loaded from
verboseVerbose flag
Remarks
The calculation considers all the active samples.
It can work in multivariate case with or without drift conditions (linked or not)
The algorithm is stopped (with a message) in the heterotopic case // TODO; improve for heterotopic case

◆ createReduce()

Model * Model::createReduce ( const VectorInt validVars) const

◆ duplicate()

Model * Model::duplicate ( ) const

◆ envelop()

VectorDouble Model::envelop ( const VectorDouble hh,
int  ivar = 0,
int  jvar = 0,
int  isign = 1,
VectorDouble  codir = VectorDouble(),
const CovCalcMode mode = nullptr 
)

◆ evalCov()

double Model::evalCov ( const VectorDouble incr,
int  icov = 0,
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const

◆ evalCovMatrix()

MatrixRectangular Model::evalCovMatrix ( Db db1,
Db db2 = nullptr,
int  ivar0 = -1,
int  jvar0 = -1,
const VectorInt nbgh1 = VectorInt(),
const VectorInt nbgh2 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline

◆ evalCovMatrixOptim()

MatrixRectangular Model::evalCovMatrixOptim ( const Db db1,
const Db db2 = nullptr,
int  ivar0 = -1,
int  jvar0 = -1,
const VectorInt nbgh1 = VectorInt(),
const VectorInt nbgh2 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline

◆ evalCovMatrixSparse()

MatrixSparse* Model::evalCovMatrixSparse ( Db db1,
Db db2 = nullptr,
int  ivar0 = -1,
int  jvar0 = -1,
const VectorInt nbgh1 = VectorInt(),
const VectorInt nbgh2 = VectorInt(),
const CovCalcMode mode = nullptr,
double  eps = EPSILON3 
)
inline

◆ evalCovMatrixSymmetric()

MatrixSquareSymmetric Model::evalCovMatrixSymmetric ( Db db1,
int  ivar0 = -1,
const VectorInt nbgh1 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline

◆ evalCovMatrixSymmetricOptim()

MatrixSquareSymmetric Model::evalCovMatrixSymmetricOptim ( const Db db1,
int  ivar0 = -1,
const VectorInt nbgh1 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline

◆ evalCovMatrixV()

VectorDouble Model::evalCovMatrixV ( Db db1,
Db db2 = nullptr,
int  ivar0 = -1,
int  jvar0 = -1,
const VectorInt nbgh1 = VectorInt(),
const VectorInt nbgh2 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline

◆ evalDrift()

double Model::evalDrift ( const Db db,
int  iech,
int  il,
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const

Evaluate a given drift function for a given sample

Parameters
dbDb structure
iechRank of the target sample
ilRank of the drift function
memberMember type (used to check filtering)
Returns

◆ evalDriftBySample()

VectorDouble Model::evalDriftBySample ( const Db db,
int  iech,
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const

◆ evalDriftBySampleInPlace()

void Model::evalDriftBySampleInPlace ( const Db db,
int  iech,
const ECalcMember &  member,
VectorDouble drftab 
) const

◆ evalDriftMatrix()

MatrixRectangular Model::evalDriftMatrix ( const Db db,
int  ivar0 = -1,
const VectorInt nbgh = VectorInt(),
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const
inline

◆ evalDriftValue()

double Model::evalDriftValue ( const Db db,
int  iech,
int  ivar,
int  ib,
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const

◆ evalDriftVarCoef()

double Model::evalDriftVarCoef ( const Db db,
int  iech,
int  ivar,
const VectorDouble coeffs 
) const

Evaluate the drift with a given sample and a given variable The value is scaled by 'coeffs'

Parameters
[in]dbDb structure
[in]iechRank of the sample
[in]ivarRank of the variable
[in]coeffsVector of coefficients

◆ evalDriftVarCoefs()

VectorDouble Model::evalDriftVarCoefs ( const Db db,
const VectorDouble coeffs,
int  ivar = 0,
bool  useSel = false 
) const

A vector of the drift evaluation (for all samples)

Parameters
dbDb structure
coeffsVector of drift coefficients
ivarVariable rank (used for constant drift value)
useSelWhen TRUE, only non masked samples are returned
Returns
The vector of values
Remarks
When no drift is defined, a vector is returned filled with the variable mean

◆ evaluateFromDb()

VectorDouble Model::evaluateFromDb ( Db db,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)

Evaluate the model on a Db

Parameters
[in]dbDb structure
[in]ivarRank of the first variable
[in]jvarRank of the second variable
[in]modeCovCalcMode structure

◆ evaluateMatInPlace()

void Model::evaluateMatInPlace ( const CovInternal covint,
const VectorDouble d1,
MatrixSquareGeneral covtab,
bool  flag_init = false,
double  weight = 1.,
const CovCalcMode mode = nullptr 
)

Returns the covariances for an increment This is the generic internal function It can be called for stationary or non-stationary case

Parameters
[in]covintInternal structure for non-stationarityAddress for the next term after the drift or NULL (for stationary case)
[in]modeCovCalcMode structure
[in]flag_initInitialize the array beforehand
[in]weightMultiplicative weight
[in]d1Distance vector
[out]covtabCovariance array

◆ evaluateOneGeneric()

double Model::evaluateOneGeneric ( const CovInternal covint,
const VectorDouble d1 = VectorDouble(),
double  weight = 1.,
const CovCalcMode mode = nullptr 
)

Returns the covariance for an increment This is the generic internal function It can be called for stationary or non-stationary case

Parameters
[in]covintInternal structure for non-stationarityAddress for the next term after the drift or NULL (for stationary case)
[in]modeCovCalcMode structure
[in]weightMultiplicative weight
[in]d1Distance vector

◆ evaluateOneIncr()

double Model::evaluateOneIncr ( double  hh,
const VectorDouble codir = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)

Calculate the value of the model for a set of distances

Returns
The model value
Parameters
[in]ivarRank of the first variable
[in]jvarRank of the second variable
[in]modeCovCalcMode structure
[in]codirArray giving the direction coefficients (optional)
[in]hhVector of increments

◆ evalZAndGradients() [1/2]

void Model::evalZAndGradients ( const SpacePoint p1,
const SpacePoint p2,
double &  covVal,
VectorDouble covGp,
VectorDouble covGG,
const CovCalcMode mode = nullptr,
bool  flagGrad = false 
) const

◆ evalZAndGradients() [2/2]

void Model::evalZAndGradients ( const VectorDouble vec,
double &  covVal,
VectorDouble covGp,
VectorDouble covGG,
const CovCalcMode mode = nullptr,
bool  flagGrad = false 
) const

◆ extensionVariance()

double Model::extensionVariance ( const Db db,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
int  ivar = 0,
int  jvar = 0 
)
inline

◆ fit()

int Model::fit ( Vario vario,
const VectorECov types = ECov::fromKeys({"SPHERICAL"}),
const Constraints constraints = Constraints(),
const Option_VarioFit optvar = Option_VarioFit(),
const Option_AutoFit mauto = Option_AutoFit(),
bool  verbose = false 
)

Automatic Fitting procedure from an experimental Variogram

Parameters
varioExperimental variogram to be fitted
typesVector of ECov (see remarks)
constraintsSet of Constraints
optvarSet of options
mautoSpecial parameters for Automatic fitting procedure (instance of Option_AutoFit), for exemple wmode (type of weighting function)
verboseVerbose option
Remarks
If no list of specific basic structure is specified, the automatic fitting is performed using a single spherical structure by default.
Returns
0 if no error, 1 otherwise

TODO : What to do with that ?

◆ fitFromCovIndices()

int Model::fitFromCovIndices ( Vario vario,
const VectorECov types = ECov::fromKeys({"EXPONENTIAL"}),
const Constraints constraints = Constraints(),
const Option_VarioFit optvar = Option_VarioFit(),
const Option_AutoFit mauto = Option_AutoFit(),
bool  verbose = false 
)

Automatic Fitting procedure

Parameters
varioExperimental variogram to be fitted
typesVector of ECov integer values
constraintsSet of Constraints
optvarSet of options
mautoSpecial parameters for Automatic fitting procedure
verboseVerbose option
Returns
0 if no error, 1 otherwise

TODO : What to do with that ?

◆ fitFromVMap()

int Model::fitFromVMap ( DbGrid dbmap,
const VectorECov types = ECov::fromKeys({"SPHERICAL"}),
const Constraints constraints = Constraints(),
const Option_VarioFit optvar = Option_VarioFit(),
const Option_AutoFit mauto = Option_AutoFit(),
bool  verbose = false 
)

Automatic Fitting procedure from A Variogram Map stored on a DbGrid

Parameters
dbmapDbGrid containing the Variogram Map
typesVector of ECov
constraintsSet of Constraints
optvarSet of options
mautoSpecial parameters for Automatic fitting procedure (instance of Option_AutoFit), for exemple wmode (type of weighting function)
verboseVerbose option
Returns
0 if no error, 1 otherwise

◆ getActiveFactor()

int Model::getActiveFactor ( ) const

◆ getAnamNClass()

int Model::getAnamNClass ( ) const

◆ getASpace()

const ASpace* Model::getASpace ( ) const
inline

◆ getContext()

const CovContext& Model::getContext ( ) const
inline

TODO : to be removed (encapsulation of Context)

◆ getCovar0()

double Model::getCovar0 ( int  ivar,
int  jvar 
) const
inline

◆ getCovar0s()

const VectorDouble& Model::getCovar0s ( ) const
inline

◆ getCovMode()

const EModelProperty & Model::getCovMode ( ) const

◆ getDimensionNumber()

int Model::getDimensionNumber ( ) const
inline

◆ getDrift()

const ADrift * Model::getDrift ( int  il) const

◆ getDriftEquationNumber()

int Model::getDriftEquationNumber ( ) const

◆ getDriftList()

const DriftList * Model::getDriftList ( ) const

TODO : to be removed (encapsulation of DriftList)

◆ getDriftMaxIRFOrder()

int Model::getDriftMaxIRFOrder ( void  ) const

◆ getDriftNumber()

int Model::getDriftNumber ( ) const

◆ getDrifts()

VectorVectorDouble Model::getDrifts ( const Db db,
bool  useSel = true 
)

◆ getExternalDriftNumber()

int Model::getExternalDriftNumber ( ) const

◆ getField()

double Model::getField ( ) const
inline

◆ getMean()

double Model::getMean ( int  ivar) const
inline

◆ getMeans()

const VectorDouble& Model::getMeans ( ) const
inline

◆ getRankFext()

int Model::getRankFext ( int  il) const

◆ getVariableNumber()

int Model::getVariableNumber ( ) const
inline

◆ gofDisplay()

void Model::gofDisplay ( double  gof,
bool  byValue = true,
const VectorDouble thresholds = {2., 5., 10., 100} 
)
static

Printout of statement concerning the Quality of the GOF

Parameters
gofValue of the Gof
byValuetrue: display GOF value; false: print its quality level
thresholdsVector giving the Quality thresholds

◆ gofToVario()

double Model::gofToVario ( const Vario vario,
bool  verbose = true 
)

Evaluate the Goodness-of_fit of the Model on the Experimental Variogram It is expressed as the average departure between Model and Variogram scaled to the variance. As this variance may be poorly calculated (< gmax / 5), it may be replaced by the largest value (gmax) divided by 2 (highly non_stationary cases).

Parameters
varioExperimental variogram
verboseVerbose flag
Returns
Value for the Goodness-of_fit (as percentage of the total sill)

◆ hasExternalCov()

int Model::hasExternalCov ( ) const

◆ initCovList()

VectorECov Model::initCovList ( const VectorInt covranks)
static

◆ isDriftDefined()

bool Model::isDriftDefined ( const VectorInt powers,
int  rank_fex = 0 
) const

◆ isDriftDifferentDefined()

bool Model::isDriftDifferentDefined ( const VectorInt powers,
int  rank_fex = -1 
) const

◆ isDriftFiltered()

bool Model::isDriftFiltered ( unsigned int  il) const

◆ isDriftSampleDefined()

bool Model::isDriftSampleDefined ( const Db db,
int  ib,
int  nech,
const VectorInt nbgh,
const ELoc &  loctype 
) const

◆ isValid()

bool Model::isValid ( ) const

◆ sample()

VectorDouble Model::sample ( const VectorDouble h,
const VectorDouble codir = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr,
const CovInternal covint = nullptr 
)

Calculate the value of the model for a set of distances

Returns
Array containing the model values
Parameters
[in]ivarRank of the first variable
[in]jvarRank of the second variable
[in]codirArray giving the direction coefficients (optional)
[in]hVector of increments
[in]modeCovCalcMode structure
[in]covintNon-stationary parameters

◆ sampleUnitary()

VectorDouble Model::sampleUnitary ( const VectorDouble hh,
int  ivar = 0,
int  jvar = 0,
VectorDouble  codir = VectorDouble(),
const CovCalcMode mode = nullptr 
)

Returns the value of the normalized covariance (by the variance/covariance value) for a given pair of variables

Parameters
hhVector of distances
ivarRank of the first variable
jvarRank of the second variable
codirDirection coefficients
modeCovCalcMode structure
Returns

◆ samplingDensityVariance()

double Model::samplingDensityVariance ( const Db db,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
int  ivar = 0,
int  jvar = 0 
) const
inline

◆ setActiveFactor()

void Model::setActiveFactor ( int  iclass)

◆ setBetaHat()

void Model::setBetaHat ( const VectorDouble betaHat)

◆ setCovaFiltered()

void Model::setCovaFiltered ( int  icov,
bool  filtered 
)

◆ setCovar0()

void Model::setCovar0 ( int  ivar,
int  jvar,
double  covar0 
)

◆ setCovar0s()

void Model::setCovar0s ( const VectorDouble covar0)

◆ setDriftFiltered()

void Model::setDriftFiltered ( int  il,
bool  filtered 
)

◆ setField()

void Model::setField ( double  field)

◆ setMarkovCoeffs()

void Model::setMarkovCoeffs ( int  icov,
const VectorDouble coeffs 
)

◆ setMean()

void Model::setMean ( double  mean,
int  ivar = 0 
)

◆ setMeans()

void Model::setMeans ( const VectorDouble mean)

◆ setRangeIsotropic()

void Model::setRangeIsotropic ( int  icov,
double  range 
)

◆ setSill()

void Model::setSill ( int  icov,
int  ivar,
int  jvar,
double  value 
)

◆ specificVolume()

double Model::specificVolume ( const Db db,
double  mean,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
int  ivar = 0,
int  jvar = 0 
) const
inline

◆ specificVolumeFromCoV()

double Model::specificVolumeFromCoV ( Db db,
double  cov,
double  mean,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
int  ivar = 0,
int  jvar = 0 
) const
inline

◆ stabilize()

int Model::stabilize ( double  percent,
bool  verbose = false 
)

Stabilize the model (in the monovariate case)

Returns
Error returned code
Parameters
[in]percentPercentage of nugget effect added
[in]verbosetrue for a verbose output
Remarks
If the model only contains GAUSSIAN structures, add
a NUGGET EFFECT structure with a sill equal to a percentage
of the total sill of the GAUSSIAN component(s)
This function does not do anything in the multivariate case

◆ standardize()

int Model::standardize ( bool  verbose = false)

Normalize the model

Parameters
[in]verbosetrue for a verbose output