1.1.0
CCC
 
Model Class Reference

Class containing the Model Information describing the formal Spatial (or Temporal) Characteristics of the (set of) random variable(s) under study. More...

#include <Model.hpp>

Inheritance diagram for Model:
AStringable ASerializable ICloneable

Public Member Functions

 Model (const CovContext &ctxt=CovContext())
 
 Model (int nvar, int ndim=2)
 
 Model (const Model &m)
 
Modeloperator= (const Model &m)
 
virtual ~Model ()
 
virtual String toString (const AStringFormat *strfmt=nullptr) const override
 ICloneable interface. More...
 
int resetFromDb (const Db *db)
 
void setCovList (const ACovAnisoList *covalist)
 
void addCov (const CovAniso *cov)
 
void addCovFromParam (const ECov &type, double range=0., double sill=1., double param=1., const VectorDouble &ranges=VectorDouble(), const VectorDouble &sills=VectorDouble(), const VectorDouble &angles=VectorDouble(), bool flagRange=true)
 
void delCova (int icov)
 
void delAllCovas ()
 
void setDriftList (const DriftList *driftlist)
 
void setDriftIRF (int order=0, int nfex=0)
 
void addDrift (const ADrift *drift)
 
void setDrifts (const VectorString &driftSymbols)
 
void delDrift (int rank)
 
void delAllDrifts ()
 
int setAnam (const AAnam *anam, const VectorInt &strcnt=VectorInt())
 
int unsetAnam ()
 
bool isFlagGradient () const
 
bool isFlagGradientNumerical () const
 
bool isFlagGradientFunctional () const
 
bool isFlagLinked () const
 
CovAniso extractCova (int icov) const
 
void switchToGradient ()
 
const ACovAnisoListgetCovAnisoList () const
 TODO : to be removed (encapsulation of ACovAnisoList) More...
 
const CovAnisogetCova (unsigned int icov) const
 
CovAnisogetCova (unsigned int icov)
 
int getCovaNumber (bool skipNugget=false) const
 
const ECov & getCovaType (int icov) const
 
const MatrixSquareSymmetric getSillValues (int icov) const
 
double getSill (int icov, int ivar, int jvar) const
 
double getParam (int icov) const
 
bool isCovaFiltered (int icov) const
 
bool isStationary () const
 
String getCovName (int icov) const
 
int getGradParamNumber (int icov) const
 
double getTotalSill (int ivar, int jvar) const
 
double getBallRadius () const
 
const AnamHermitegetAnamHermite () const
 
double getMaximumDistance () const
 
int getCovaMinIRFOrder () const
 
bool hasAnam () const
 
const AAnamgetAnam () const
 
bool isChangeSupportDefined () const
 
void normalize (double sill)
 
bool hasNugget () const
 
VectorInt getActiveCovList () const
 
VectorInt getAllActiveCovList () const
 
bool isAllActiveCovList () const
 
void setTapeRange (double range)
 
void setOptimEnabled (bool flagOptim)
 
bool isOptimEnabled () const
 
double eval0 (int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
MatrixSquareGeneral eval0Nvar (const CovCalcMode *mode=nullptr) const
 
int isNoStat () const
 
const ANoStatgetNoStat () const
 
ANoStatgetNoStatModify () const
 
void eval0MatInPlace (MatrixSquareGeneral &mat, const CovCalcMode *mode=nullptr) const
 
double eval (const SpacePoint &p1, const SpacePoint &p2, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
MatrixSquareGeneral evalNvarIpas (double step, const VectorDouble &dir, const CovCalcMode *mode=nullptr) const
 
MatrixSquareGeneral evalMat (const SpacePoint &p1, const SpacePoint &p2, const CovCalcMode *mode=nullptr) const
 
void evalMatInPlace (const SpacePoint &p1, const SpacePoint &p2, MatrixSquareGeneral &mat, const CovCalcMode *mode=nullptr) const
 
MatrixSquareGeneral evalNvarIpasIncr (const VectorDouble &dincr, const CovCalcMode *mode=nullptr) const
 
VectorDouble evalIvarNpas (const VectorDouble &vec_step, const VectorDouble &dir=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
double evalIvarIpas (double step, const VectorDouble &dir=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
double evalCvv (const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
double evalCvvShift (const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &shift, const VectorDouble &angles=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
MatrixSquareGeneral evalCvvM (const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const CovCalcMode *mode=nullptr)
 
double evalCxv (const SpacePoint &p1, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
MatrixSquareGeneral evalCxvM (const SpacePoint &p1, const VectorDouble &ext, const VectorInt &ndisc, const VectorDouble &angles=VectorDouble(), const VectorDouble &x0=VectorDouble(), const CovCalcMode *mode=nullptr)
 
VectorDouble evalPointToDb (const SpacePoint &p1, const Db *db2, int ivar=0, int jvar=0, bool useSel=true, const VectorInt &nbgh2=VectorInt(), const CovCalcMode *mode=nullptr)
 
VectorDouble evalPointToDbAsSP (const std::vector< SpacePoint > &p1s, const SpacePoint &p2, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
double evalAverageDbToDb (const Db *db1, const Db *db2, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
double evalAverageIncrToIncr (const VectorVectorDouble &d1, const VectorVectorDouble &d2, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr) const
 
double evalAveragePointToDb (const SpacePoint &p1, const Db *db2, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
MatrixRectangular evalCovMatrix (Db *db1, Db *db2=nullptr, int ivar=0, int jvar=0, const VectorInt &nbgh1=VectorInt(), const VectorInt &nbgh2=VectorInt(), const CovCalcMode *mode=nullptr)
 
void evalMatOptimInPlace (int icas1, int iech1, int icas2, int iech2, MatrixSquareGeneral &mat, const CovCalcMode *mode=nullptr) const
 
VectorVectorDouble evalCovMatrixOptim (const Db *db1, const Db *db2=nullptr, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
double 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 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 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 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 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 evalZAndGradients (const SpacePoint &p1, const SpacePoint &p2, double &covVal, VectorDouble &covGp, VectorDouble &covGG, const CovCalcMode *mode=nullptr, bool flagGrad=false) const
 
void evalZAndGradients (const VectorDouble &vec, double &covVal, VectorDouble &covGp, VectorDouble &covGG, const CovCalcMode *mode=nullptr, bool flagGrad=false) const
 
double evalCov (const VectorDouble &incr, int icov=0, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
void setSill (int icov, int ivar, int jvar, double value)
 
void setCovaFiltered (int icov, bool filtered)
 
void updateCovByPoints (int icas1, int iech1, int icas2, int iech2)
 
void updateCovByMesh (int imesh)
 
void setActiveFactor (int iclass)
 
int getActiveFactor () const
 
int getAnamNClass () const
 
const DriftListgetDriftList () const
 TODO : to be removed (encapsulation of DriftList) More...
 
const ADriftgetDrift (int il) const
 
ADriftgetDrift (int il)
 
int getDriftNumber () const
 
int getExternalDriftNumber () const
 
int getRankFext (int il) const
 
const VectorDoublegetDriftCLs () const
 
double getDriftCL (int ivar, int il, int ib) const
 
int getDriftEquationNumber () const
 
bool isDriftFiltered (unsigned int il) const
 
bool isDriftDefined (const VectorInt &powers, int rank_fex=0) const
 
bool isDriftDifferentDefined (const VectorInt &powers, int rank_fex=-1) const
 
int getDriftMaxIRFOrder (void) const
 
void resetDriftCoef ()
 
void setDriftCoef (int ivar, int il, int ib, double coeff)
 
void setDriftFiltered (int il, bool filtered)
 
VectorDouble getDriftByColumn (const Db *db, int ib, bool useSel=true)
 
VectorVectorDouble getDrifts (const Db *db, bool useSel=true)
 
double evalDrift (const Db *db, int iech, int il, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
VectorDouble evalDriftVec (const Db *db, int iech, const ECalcMember &member=ECalcMember::fromKey("LHS")) const
 
VectorDouble evalDriftCoefVec (const Db *db, const VectorDouble &coeffs, int ivar=0, bool useSel=false) const
 
void evalDriftVecInPlace (const Db *db, int iech, const ECalcMember &member, VectorDouble &drftab) const
 
double evalDriftCoef (const Db *db, int iech, int ivar, const VectorDouble &coeffs) const
 
const CovContextgetContext () const
 TODO : to be removed (encapsulation of Context) More...
 
const VectorDoublegetMeans () const
 
double getMean (int ivar) const
 
const VectorDoublegetCovar0s () const
 
double getCovar0 (int ivar, int jvar) const
 
double getField () const
 
int getDimensionNumber () const
 
void setMeans (const VectorDouble &mean)
 
void setMean (double mean, int ivar=0)
 
void setCovar0s (const VectorDouble &covar0)
 
void setCovar0 (int ivar, int jvar, double covar0)
 
void setField (double field)
 
int addNoStat (const ANoStat *anostat)
 Shortcut for Non-stationary. More...
 
int getNoStatElemNumber () const
 
int addNoStatElem (int igrf, int icov, const EConsElem &type, int iv1, int iv2)
 
int addNoStatElems (const VectorString &codes)
 
CovParamId getCovParamId (int ipar) const
 
const EModelProperty & getCovMode () const
 
Modelduplicate () const
 
ModelcreateReduce (const VectorInt &validVars) const
 
int getVariableNumber () const
 
int hasExternalCov () const
 
MatrixSquareSymmetric covMatrixM (Db *db1, Db *db2=nullptr, int ivar=-1, int jvar=-1, const CovCalcMode *mode=nullptr)
 
VectorDouble covMatrixV (Db *db1, Db *db2=nullptr, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
void covMatrix (VectorDouble &covmat, Db *db1, Db *db2=nullptr, int ivar=0, int jvar=0, const CovCalcMode *mode=nullptr)
 
VectorDouble sample (const VectorDouble &hh, int ivar=0, int jvar=0, VectorDouble codir=VectorDouble(), const CovCalcMode *mode=nullptr)
 
VectorDouble sampleUnitary (const VectorDouble &hh, int ivar=0, int jvar=0, VectorDouble codir=VectorDouble(), const CovCalcMode *mode=nullptr)
 
VectorDouble envelop (const VectorDouble &hh, int ivar=0, int jvar=0, int isign=1, VectorDouble codir=VectorDouble(), const CovCalcMode *mode=nullptr)
 
int fitFromCovIndices (Vario *vario, const VectorECov &types=ECov::fromKeys({"EXPONENTIAL"}), const Constraints &constraints=Constraints(), Option_VarioFit optvar=Option_VarioFit(), Option_AutoFit mauto=Option_AutoFit(), bool verbose=false)
 
int fit (Vario *vario, const VectorECov &types=ECov::fromKeys({"SPHERICAL"}), const Constraints &constraints=Constraints(), Option_VarioFit optvar=Option_VarioFit(), Option_AutoFit mauto=Option_AutoFit(), bool verbose=false)
 
int fitFromVMap (DbGrid *dbmap, const VectorECov &types=ECov::fromKeys({"SPHERICAL"}), const Constraints &constraints=Constraints(), Option_VarioFit optvar=Option_VarioFit(), Option_AutoFit mauto=Option_AutoFit(), bool verbose=false)
 
int buildVmapOnDbGrid (DbGrid *dbgrid, const NamingConvention &namconv=NamingConvention("VMAP"))
 
int stabilize (double percent, bool verbose=false)
 
int standardize (bool verbose=false)
 
double gofToVario (const Vario *vario, bool verbose=true)
 
void gofDisplay (double gof, bool byValue=true, const VectorDouble &thresholds={2., 5., 10., 100})
 
VectorECov initCovList (const VectorInt &covranks)
 
bool isValid () const
 
- Public Member Functions inherited from AStringable
 AStringable ()
 
 AStringable (const AStringable &r)
 
AStringableoperator= (const AStringable &r)
 
virtual ~AStringable ()
 
virtual void display (const AStringFormat *strfmt=nullptr) const final
 
virtual void display (int level) const final
 
- Public Member Functions inherited from ASerializable
 ASerializable ()
 
 ASerializable (const ASerializable &r)
 
ASerializableoperator= (const ASerializable &r)
 
virtual ~ASerializable ()
 
bool deserialize (std::istream &is, bool verbose=true)
 
bool serialize (std::ostream &os, bool verbose=true) const
 
bool dumpToNF (const String &neutralFilename, bool verbose=false) const
 
- Public Member Functions inherited from ICloneable
 ICloneable ()
 
virtual ~ICloneable ()
 
virtual ICloneableclone () const =0
 

Static Public Member Functions

static Modelcreate (const CovContext &ctxt=CovContext())
 
static ModelcreateFromEnvironment (int nvar, int ndim=2)
 
static ModelcreateNugget (int nvar, int ndim=2, double sill=1.)
 
static ModelcreateFromParam (const ECov &type=ECov::fromKey("NUGGET"), double range=1., double sill=1., double param=1., const VectorDouble &ranges=VectorDouble(), const VectorDouble &sills=VectorDouble(), const VectorDouble &angles=VectorDouble(), const ASpace *space=nullptr, bool flagRange=true)
 
static ModelcreateFromDb (const Db *db)
 
static ModelcreateFromNF (const String &neutralFilename, bool verbose=true)
 
- Static Public Member Functions inherited from ASerializable
static String buildFileName (int status, const String &filename, bool ensureDirExist=false)
 
static String getHomeDirectory (const String &sub="")
 
static String getWorkingDirectory ()
 
static String getTestData (const String &subdir, const String &filename)
 
static String getFileIdentity (const String &filename, bool verbose=false)
 
static void setContainerName (bool useDefault, const String &containerName="", bool verbose=false)
 
static void unsetContainerName ()
 
static void setPrefixName (const String &prefixName)
 
static void unsetPrefixName ()
 
static const StringgetContainerName ()
 
static const StringgetPrefixName ()
 
static bool createDirectory (const String &dir)
 
static String getExecDirectory ()
 
static String getDirectory (const String &path)
 

Detailed Description

Class containing the Model Information describing the formal Spatial (or Temporal) Characteristics of the (set of) random variable(s) under study.

The Model is essentially a container with two main contents:

  • the covariance part: see ACov.hpp for more information
  • the drift part: see DriftList.hpp for more information

The additional member CovContext only serves in carrying the following information:

  • the number of variables: if more than 1, the Model becomes multivariate
  • the field extension: this information is needed to get a stationary version to any covariance
  • the experimental mean vector and the variance-covariance matrix (used to calibrate the Model)

Constructor & Destructor Documentation

Model::Model ( const CovContext ctxt = CovContext())
Model::Model ( int  nvar,
int  ndim = 2 
)
Model::Model ( const Model m)
Model::~Model ( )
virtual

Member Function Documentation

void Model::addCov ( const CovAniso cov)
void Model::addCovFromParam ( const ECov &  type,
double  range = 0.,
double  sill = 1.,
double  param = 1.,
const VectorDouble ranges = VectorDouble(),
const VectorDouble sills = VectorDouble(),
const VectorDouble angles = VectorDouble(),
bool  flagRange = true 
)
void Model::addDrift ( const ADrift drift)
int Model::addNoStat ( const ANoStat anostat)

Shortcut for Non-stationary.

Define Non-stationary parameters

Parameters
anostatANoStat pointer will be duplicated
Returns
Error return code
int Model::addNoStatElem ( int  igrf,
int  icov,
const EConsElem &  type,
int  iv1,
int  iv2 
)
int Model::addNoStatElems ( const VectorString codes)
int Model::buildVmapOnDbGrid ( DbGrid dbgrid,
const NamingConvention namconv = NamingConvention("VMAP") 
)

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

Returns
Error return code
Parameters
[in]dbgridGrid structure
[in]namconvNaming convention
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
void Model::covMatrix ( VectorDouble covmat,
Db db1,
Db db2 = nullptr,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)

Calculate the covariance matrix between active samples of Db1 and active samples of Db2.

Parameters
covmatReturned matrix (returned as a vector).
db1First Data Base
db2Second Data Base (if not provided, the first Db is provided instead)
ivarRank of the first variable (all variables if not defined)
jvarRank of the second variable (all variables if not defined)
modeCovCalcMode structure
Remarks
The returned argument must have been dimensioned beforehand to (nvar * nechA)^2 where:
-nvar stands for the number of (active) variables
-nechA stands for the number of active samples
MatrixSquareSymmetric Model::covMatrixM ( Db db1,
Db db2 = nullptr,
int  ivar = -1,
int  jvar = -1,
const CovCalcMode mode = nullptr 
)
VectorDouble Model::covMatrixV ( Db db1,
Db db2 = nullptr,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)
Model * Model::create ( const CovContext ctxt = CovContext())
static
Model * Model::createFromDb ( const Db db)
static
Model * Model::createFromEnvironment ( int  nvar,
int  ndim = 2 
)
static
Model * Model::createFromNF ( const String neutralFilename,
bool  verbose = true 
)
static
Model * Model::createFromParam ( const ECov &  type = ECov::fromKey("NUGGET"),
double  range = 1.,
double  sill = 1.,
double  param = 1.,
const VectorDouble ranges = VectorDouble(),
const VectorDouble sills = VectorDouble(),
const VectorDouble angles = VectorDouble(),
const ASpace space = nullptr,
bool  flagRange = true 
)
static
Model * Model::createNugget ( int  nvar,
int  ndim = 2,
double  sill = 1. 
)
static
Model * Model::createReduce ( const VectorInt validVars) const
void Model::delAllCovas ( )
void Model::delAllDrifts ( )
void Model::delCova ( int  icov)
void Model::delDrift ( int  rank)
Model * Model::duplicate ( ) const
VectorDouble Model::envelop ( const VectorDouble hh,
int  ivar = 0,
int  jvar = 0,
int  isign = 1,
VectorDouble  codir = VectorDouble(),
const CovCalcMode mode = nullptr 
)
double Model::eval ( const SpacePoint p1,
const SpacePoint p2,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
double Model::eval0 ( int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
void Model::eval0MatInPlace ( MatrixSquareGeneral mat,
const CovCalcMode mode = nullptr 
) const
inline

Calculate the Matrix of covariance for zero distance

Parameters
matCovariance matrix (Dimension: nvar * nvar)
modeCalculation Options
Remarks
: Matrix 'mat' should be dimensioned and initialized beforehand
MatrixSquareGeneral Model::eval0Nvar ( const CovCalcMode mode = nullptr) const
inline
double Model::evalAverageDbToDb ( const Db db1,
const Db db2,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
double Model::evalAverageIncrToIncr ( const VectorVectorDouble d1,
const VectorVectorDouble d2,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
double Model::evalAveragePointToDb ( const SpacePoint p1,
const Db db2,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)
inline
double Model::evalCov ( const VectorDouble incr,
int  icov = 0,
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const
MatrixRectangular Model::evalCovMatrix ( Db db1,
Db db2 = nullptr,
int  ivar = 0,
int  jvar = 0,
const VectorInt nbgh1 = VectorInt(),
const VectorInt nbgh2 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline
VectorVectorDouble Model::evalCovMatrixOptim ( const Db db1,
const Db db2 = nullptr,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)
double Model::evalCvv ( const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
MatrixSquareGeneral Model::evalCvvM ( const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const CovCalcMode mode = nullptr 
)
inline
double Model::evalCvvShift ( const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble shift,
const VectorDouble angles = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
double Model::evalCxv ( const SpacePoint p1,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
)
inline
MatrixSquareGeneral Model::evalCxvM ( const SpacePoint p1,
const VectorDouble ext,
const VectorInt ndisc,
const VectorDouble angles = VectorDouble(),
const VectorDouble x0 = VectorDouble(),
const CovCalcMode mode = nullptr 
)
inline
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
double Model::evalDriftCoef ( const Db db,
int  iech,
int  ivar,
const VectorDouble coeffs 
) const

Evaluate the drift with a given set of coefficients

Parameters
[in]dbDb structure
[in]iechRank of the sample
[in]ivarRank of the variable
[in]coeffsVector of coefficients
VectorDouble Model::evalDriftCoefVec ( const Db db,
const VectorDouble coeffs,
int  ivar = 0,
bool  useSel = false 
) const

A vector of the drift evaluation

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
VectorDouble Model::evalDriftVec ( const Db db,
int  iech,
const ECalcMember &  member = ECalcMember::fromKey("LHS") 
) const
void Model::evalDriftVecInPlace ( const Db db,
int  iech,
const ECalcMember &  member,
VectorDouble drftab 
) const
double Model::evalIvarIpas ( double  step,
const VectorDouble dir = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
VectorDouble Model::evalIvarNpas ( const VectorDouble vec_step,
const VectorDouble dir = VectorDouble(),
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
MatrixSquareGeneral Model::evalMat ( const SpacePoint p1,
const SpacePoint p2,
const CovCalcMode mode = nullptr 
) const
inline
void Model::evalMatInPlace ( const SpacePoint p1,
const SpacePoint p2,
MatrixSquareGeneral mat,
const CovCalcMode mode = nullptr 
) const
inline

Calculate the Matrix of covariance between two space points

Parameters
p1Reference of the first space point
p2Reference of the second space point
matCovariance matrix (Dimension: nvar * nvar)
modeCalculation Options
Remarks
: Matrix 'mat' should be dimensioned and initialized beforehand
void Model::evalMatOptimInPlace ( int  icas1,
int  iech1,
int  icas2,
int  iech2,
MatrixSquareGeneral mat,
const CovCalcMode mode = nullptr 
) const
inline

Calculate the Matrix of covariance between two elements of two Dbs (defined beforehand)

Parameters
icas1Origin of the Db containing the first point
iech1Rank of the first point
icas2Origin of the Db containing the second point
iech2Rank of the second point
matCovariance matrix (Dimension: nvar * nvar)
modeCalculation Options
Remarks
: Matrix 'mat' should be dimensioned and initialized beforehand
MatrixSquareGeneral Model::evalNvarIpas ( double  step,
const VectorDouble dir,
const CovCalcMode mode = nullptr 
) const
inline
MatrixSquareGeneral Model::evalNvarIpasIncr ( const VectorDouble dincr,
const CovCalcMode mode = nullptr 
) const
inline
VectorDouble Model::evalPointToDb ( const SpacePoint p1,
const Db db2,
int  ivar = 0,
int  jvar = 0,
bool  useSel = true,
const VectorInt nbgh2 = VectorInt(),
const CovCalcMode mode = nullptr 
)
inline
VectorDouble Model::evalPointToDbAsSP ( const std::vector< SpacePoint > &  p1s,
const SpacePoint p2,
int  ivar = 0,
int  jvar = 0,
const CovCalcMode mode = nullptr 
) const
inline
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::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
CovAniso Model::extractCova ( int  icov) const
int Model::fit ( Vario vario,
const VectorECov types = ECov::fromKeys({"SPHERICAL"}),
const Constraints constraints = Constraints(),
Option_VarioFit  optvar = Option_VarioFit(),
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 ?

int Model::fitFromCovIndices ( Vario vario,
const VectorECov types = ECov::fromKeys({"EXPONENTIAL"}),
const Constraints constraints = Constraints(),
Option_VarioFit  optvar = Option_VarioFit(),
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 ?

int Model::fitFromVMap ( DbGrid dbmap,
const VectorECov types = ECov::fromKeys({"SPHERICAL"}),
const Constraints constraints = Constraints(),
Option_VarioFit  optvar = Option_VarioFit(),
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
VectorInt Model::getActiveCovList ( ) const
int Model::getActiveFactor ( ) const
VectorInt Model::getAllActiveCovList ( ) const
const AAnam * Model::getAnam ( ) const
const AnamHermite * Model::getAnamHermite ( ) const
int Model::getAnamNClass ( ) const
double Model::getBallRadius ( ) const

Returns the Ball radius (from the first covariance of _covaList)

Returns
Value of the Ball Radius (if defined, i.e. for Numerical Gradient calculation)
const CovContext& Model::getContext ( ) const
inline

TODO : to be removed (encapsulation of Context)

const CovAniso * Model::getCova ( unsigned int  icov) const
CovAniso * Model::getCova ( unsigned int  icov)
int Model::getCovaMinIRFOrder ( ) const
const ACovAnisoList * Model::getCovAnisoList ( ) const

TODO : to be removed (encapsulation of ACovAnisoList)

int Model::getCovaNumber ( bool  skipNugget = false) const
double Model::getCovar0 ( int  ivar,
int  jvar 
) const
inline
const VectorDouble& Model::getCovar0s ( ) const
inline
const ECov & Model::getCovaType ( int  icov) const
const EModelProperty & Model::getCovMode ( ) const
String Model::getCovName ( int  icov) const
CovParamId Model::getCovParamId ( int  ipar) const
int Model::getDimensionNumber ( ) const
inline
const ADrift * Model::getDrift ( int  il) const
ADrift * Model::getDrift ( int  il)
VectorDouble Model::getDriftByColumn ( const Db db,
int  ib,
bool  useSel = true 
)
double Model::getDriftCL ( int  ivar,
int  il,
int  ib 
) const
const VectorDouble & Model::getDriftCLs ( ) const
int Model::getDriftEquationNumber ( ) const
const DriftList * Model::getDriftList ( ) const

TODO : to be removed (encapsulation of DriftList)

int Model::getDriftMaxIRFOrder ( void  ) const
int Model::getDriftNumber ( ) const
VectorVectorDouble Model::getDrifts ( const Db db,
bool  useSel = true 
)
int Model::getExternalDriftNumber ( ) const
double Model::getField ( ) const
inline
int Model::getGradParamNumber ( int  icov) const
double Model::getMaximumDistance ( ) const
double Model::getMean ( int  ivar) const
inline
const VectorDouble& Model::getMeans ( ) const
inline
const ANoStat* Model::getNoStat ( ) const
inline
int Model::getNoStatElemNumber ( ) const
ANoStat* Model::getNoStatModify ( ) const
inline
double Model::getParam ( int  icov) const
int Model::getRankFext ( int  il) const
double Model::getSill ( int  icov,
int  ivar,
int  jvar 
) const
const MatrixSquareSymmetric Model::getSillValues ( int  icov) const
double Model::getTotalSill ( int  ivar,
int  jvar 
) const
int Model::getVariableNumber ( ) const
inline
void Model::gofDisplay ( double  gof,
bool  byValue = true,
const VectorDouble thresholds = {2., 5., 10., 100} 
)

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
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)
bool Model::hasAnam ( ) const
int Model::hasExternalCov ( ) const
bool Model::hasNugget ( ) const
VectorECov Model::initCovList ( const VectorInt covranks)
bool Model::isAllActiveCovList ( ) const
bool Model::isChangeSupportDefined ( ) const
bool Model::isCovaFiltered ( int  icov) 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::isDriftFiltered ( unsigned int  il) const
bool Model::isFlagGradient ( ) const
bool Model::isFlagGradientFunctional ( ) const
bool Model::isFlagGradientNumerical ( ) const
bool Model::isFlagLinked ( ) const
int Model::isNoStat ( ) const
inline
bool Model::isOptimEnabled ( ) const
inline
bool Model::isStationary ( ) const
bool Model::isValid ( ) const
void Model::normalize ( double  sill)
Model & Model::operator= ( const Model m)
void Model::resetDriftCoef ( )
inline
int Model::resetFromDb ( const Db db)
VectorDouble Model::sample ( const VectorDouble hh,
int  ivar = 0,
int  jvar = 0,
VectorDouble  codir = VectorDouble(),
const CovCalcMode mode = nullptr 
)

Sample a Model for given variable(s) and given direction

Parameters
hhVector of distances
ivarRank of the first variable
jvarRank of the second variable
codirVector of direction coefficients
modeCovCalMode structure
Returns
The array of variogram evaluated at discretized positions
Note that its dimension is 'nh' (if 'addZero' is false and 'nh+1' otherwise)
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
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
void Model::setActiveFactor ( int  iclass)
int Model::setAnam ( const AAnam anam,
const VectorInt strcnt = VectorInt() 
)

Defining an Anamorphosis information for the Model (in fact, this is added to ACovAnisoList part and transforms it from CovLMC to CovLMCAnamorphosis

Parameters
anamPointer to the anamorphosis
strcntArray of covariance description used for IR case
Returns
void Model::setCovaFiltered ( int  icov,
bool  filtered 
)
void Model::setCovar0 ( int  ivar,
int  jvar,
double  covar0 
)
void Model::setCovar0s ( const VectorDouble covar0)
void Model::setCovList ( const ACovAnisoList covalist)

Add a list of Covariances. This operation cleans any previously stored covariance

Parameters
covalistList of Covariances to be added
void Model::setDriftCoef ( int  ivar,
int  il,
int  ib,
double  coeff 
)
void Model::setDriftFiltered ( int  il,
bool  filtered 
)
void Model::setDriftIRF ( int  order = 0,
int  nfex = 0 
)

Define the list of drift functions for:

  • a given degree of the IRF
  • a given number of external drifts
    Parameters
    orderOrder of the IRF
    nfexNumber of External Drifts
    Remarks
    This method deletes any pre-existing drift functions and replaces them by the new definition
    This replacement is performed accounting for information stored in 'model', such as:
  • the space dimension
  • the number of variables
void Model::setDriftList ( const DriftList driftlist)

Add a list of Drifts. This operation cleans any previously stored drift function

Parameters
driftlistList of Drifts to be added
Remarks
This method deletes any pre-existing drift functions
void Model::setDrifts ( const VectorString driftSymbols)
void Model::setField ( double  field)
void Model::setMean ( double  mean,
int  ivar = 0 
)
void Model::setMeans ( const VectorDouble mean)
void Model::setOptimEnabled ( bool  flagOptim)
inline
void Model::setSill ( int  icov,
int  ivar,
int  jvar,
double  value 
)
void Model::setTapeRange ( double  range)
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
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
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
int Model::standardize ( bool  verbose = false)

Normalize the model

Parameters
[in]verbosetrue for a verbose output
void Model::switchToGradient ( )

Switch to a Model dedicated to Gradients (transforms it from CovLMC to CovLMGradient)

String Model::toString ( const AStringFormat strfmt = nullptr) const
overridevirtual

ICloneable interface.

AStringable Interface

Reimplemented from AStringable.

int Model::unsetAnam ( )
void Model::updateCovByMesh ( int  imesh)
void Model::updateCovByPoints ( int  icas1,
int  iech1,
int  icas2,
int  iech2 
)

The documentation for this class was generated from the following files: