Documentation / Manuel développeur
Modules disponibles
Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,  
Uranie / Modeler v4.9.0
|
Description of the class TKriging. More...
#include <TKriging.h>
Public Member Functions | |
Constructor and Destructor | |
TKriging (Int_t nX, Int_t nS, Double_t *xObs, Double_t *xNormParams, Double_t *yObs, Double_t *gamma, Double_t *iL, Double_t gpVar, Double_t vErrMes, URANIE::Modeler::TCorrelationFunction *corrFunc, TString varNames, Int_t nP=0, Double_t *hbeta=NULL, Double_t *vbeta=NULL, Double_t *M=NULL, Double_t *iR=NULL, TString trend="") | |
Default constructor. | |
void | initMatrices () |
Initialisation of the matrices. | |
virtual | ~TKriging () |
Default destructor. | |
Get methods | |
const char * | getInputNames () |
Return the input variable names. | |
const char * | getOutputName () |
Return the output variable name. | |
Bool_t | hasTrend () |
Return kTRUE if a deterministic trend has been defined. | |
TCorrelationFunction * | getCorrFunction () |
Return a pointer to the correlation function of the gaussian process. | |
Int_t | getNInputs () |
Return the number of input parameters. | |
Int_t | getNObs () |
Return the number of observations used to build the GP model. | |
Int_t | getNCorrFunctionParams () |
Return the number of parameters of the correlation function. | |
Int_t | getNCorrLengths () |
Return the number of correlation lengths. | |
Int_t | getNTrendParams () |
Return the number of parameters of the deterministic trend. | |
Double_t * | getObsInputs () |
Return the observations inputs as a flat array: [a1, b1, c1, a2, b2, c2,...]. | |
Double_t * | getObsOutputs () |
Return the observations outputs. | |
Double_t * | getObsNormParams () |
Return the normalisation parameters of the inputs as a flat array: [aMin, aMax, bMin, bMax,...]. | |
Double_t * | getGPParams () |
Return the array of the Gaussian Process parameters. | |
Double_t * | getTrendParams () |
Return the array of the trend parameters. | |
Double_t * | getTrendParamsCovMatrix () |
Return the array of the trend parameters covariance matrix. | |
TObjArray * | getTrendFormula () |
Return the deterministic trend list of formulas. | |
TString | getTrendString () |
Return the deterministic trend string. | |
Double_t | getVariance () |
Return the variance of the gaussian process. | |
Double_t | getMeasurementErrorVariance () |
Return the variance of the measurement error. | |
Double_t * | getiL () |
Return the iL matrix as an array (cf. internal variable description) | |
Double_t * | getM () |
Return the M matrix as an array (cf. internal variable description) | |
Double_t * | getiR () |
Return the iR matrix as an array (cf. internal variable description) | |
Double_t * | getLooErrors () |
Return the vector of the leave one out errors. | |
Double_t * | getLooVariances () |
Return the vector of the leave one out prediction variances. | |
void | getLooErrors (double *arr, int size) |
Return the vector of the leave one out errors. | |
void | getLooVariances (double *arr, int size) |
Return the vector of the leave one out prediction variances. | |
Double_t | getLooRMSE () |
Return the leave one out mean squared error. | |
Double_t | getLooQ2 () |
Return the Q2 of the model on the observation. | |
URANIE::DataServer::TDataServer * | getLooData () |
Return a pointer to a TDataServer containing the observations and the Leave One Out results. | |
Double_t * | getCorrLengthsNormalised () |
Return the normalised correlation lengths. | |
Double_t * | getCorrLengths () |
Return the correlation lengths. | |
Double_t * | getCorrFunctionParams () |
Return all the parameters of the correlation function. | |
Evaluation | |
Int_t | eval (Double_t *x0, Double_t *y0, int=0) |
Evaluate the gaussian process at the new location x0. | |
void | estimate (URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", Option_t *option="", Bool_t useGPU=false) |
void | estimate_CPU (URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", Option_t *option="") |
void | estimateWithCov (URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", TString realValue="", Option_t *option="") |
Tools | |
Exporting the kriging model as a stand alone C function is difficult for various reasons. The main problem is that we would have to store potentially huge matrices as text, which is highly unpractical. It was thus decided to postpone this export as long as a workable solution isn't found. | |
void | exportFunction (const char *lang, const char *file, const char *name, Option_t *option="") |
virtual void | exportModelCplusplus (std::ofstream *sourcefile, const char *name="") |
Export the model as a C function (not implemented) | |
virtual void | exportModelFortran (std::ofstream *sourcefile) |
Export the model as a Fortran function (not implemented) | |
virtual void | exportModelPMML (const char *file, const char *name, Option_t *option) const |
void | freeze () |
Printing Log | |
virtual void | printLog (Option_t *option="") |
Public Attributes | |
TString | _outputVarName |
Name of the output variable. This information can be used to remind the user of the original names of the variables. | |
TString | _inputVarNames |
Name of the inputs variables separated by ":". This information can be used to remind the user of the original names of the variables. | |
TObjArray * | _inputNamesList |
list of the inputs' names. | |
Int_t | _nX |
Number of input attributes | |
Int_t | _nS |
Number of observations. | |
Int_t | _nYDeclared |
Number of attribute added to the tkiging object before running it. | |
Int_t | _nXDeclared |
Number of attribute added to the tkiging object before running it. | |
Double_t * | _xObs |
The matrix of inputs (dimension: _nX * _nS) | |
Double_t * | _xNormParams |
The matrix of normalisation parameters. Contains min and max value of each input variable: [aMin, aMax, bMin, bMax,...] (size: _nX * 2) | |
Double_t * | _yObs |
Array representation of the outputs (size: _nS) | |
Double_t * | _iL |
Inverse of the matrix L, Cholesky decomposition of K (size: _nS * _nS) | |
Double_t * | _gamma |
Array of parameters of the gaussian process (size: _nS) | |
Double_t | _sigma2 |
Variance of the gaussian process. | |
Double_t | _vErrMes |
variance of the measurement error | |
TCorrelationFunction * | _correlationFunction |
The correlation function. | |
Double_t * | _corrLengths |
The correlation lengths in the real input space. | |
Double_t * | _normCorrLengths |
The correlation lengths in the normalised input space. This corresponds to the correlation lengths stored in the correlation function. | |
Bool_t | bUse_normalisation |
if kFALSE, the matrix H will be computed with unormalised input data. This happens when the user is providing his own trend function. When the trend is constructed automatically from a trend type, this attribute is set to kTRUE. | |
Bool_t | bHas_trend |
kFALSE if no deterministic trend is defined, kTRUE otherwise | |
TString | _trendString |
Character string containing the trend parameters separated by ":" or a trend type ("const" or "linear") | |
TObjArray * | _trendCoefList |
List of the formulas of the trend coefficients (_nP elements) | |
Int_t | _nP |
Number of coefficients of the deterministic trend. | |
Double_t * | _hbeta |
Estimates of the deterministic trend parameters (size: _nP) | |
Double_t * | _vbeta |
Covariance matrix of the trend parameters (size: _nP*_nP) | |
Double_t * | _M |
Matrix M = L^{-T} H (size: _nS * _nP) | |
Double_t * | _iR |
Inverse of the matrix R, from the QR decomposition of M (size: _nP * _nP) | |
Double_t * | _errorLoo |
If available, vector of the Leave One Out errors for the model. | |
Double_t * | _varLoo |
If available, vector of the variances of the Leave One Out errors for the model. | |
Double_t | _mseLoo |
if available, Leave one Out mean squared error of the model. | |
Private Attributes | |
Bool_t | _blog |
Boolean for edit the log. | |
Set methods | |
void | setCorrFunction (TCorrelationFunction *corrFunc) |
Set the correlation function of the gaussian process. | |
void | setVariance (Double_t var) |
Set the variance of the gaussian process. | |
void | setLooErrors (Double_t *err, Double_t *var, Double_t mse) |
Set the information about the leave one out error of the model. | |
void | createTrend () |
Create the list of trend coefficients. | |
void | computeCovarianceMatrix (int ns, double *xObs, bool bVariance_in_factor, bool bHas_measurement_error, double alpha, double reg, double sigma2, double *Vmes, double *C, double *K) |
static void | StatMat (std::string title, int n1, int n2, Double_t *v) |
static void | StatMat (std::string title, int n1, int n2, Double_t *v, std::vector< double > &tvmin, std::vector< double > &tvmax, std::vector< double > &tmean, std::vector< double > &tvar, std::vector< double > &tnorm) |
static void | StatVect (std::string title, int n, Double_t *v) |
static void | StatVect (std::string title, int n, Double_t *v, double &vmin, double &vmax, double &mean, double &var, double &norm) |
Detailed Description
Description of the class TKriging.
This class defines a Kriging model. It can predict the output of new observations using the TKriging::eval function or the URANIE::Relauncher module formalism. It is created by TGPBuilder::buildGP. It is a wrapper to the gpPredX
functions from the "gaussian process Library" (gpL), by Jean-Marc Martinez.
Constructor & Destructor Documentation
◆ TKriging()
URANIE::Modeler::TKriging::TKriging | ( | Int_t | nX, |
Int_t | nS, | ||
Double_t * | xObs, | ||
Double_t * | xNormParams, | ||
Double_t * | yObs, | ||
Double_t * | gamma, | ||
Double_t * | iL, | ||
Double_t | gpVar, | ||
Double_t | vErrMes, | ||
URANIE::Modeler::TCorrelationFunction * | corrFunc, | ||
TString | varNames, | ||
Int_t | nP = 0 , |
||
Double_t * | hbeta = NULL , |
||
Double_t * | vbeta = NULL , |
||
Double_t * | M = NULL , |
||
Double_t * | iR = NULL , |
||
TString | trend = "" |
||
) |
Default constructor.
Constructs a TKriging from the raw matrices and other stuff produced by TGPBuilder. This constructor is not supposed to be used by a human being.
- Parameters
-
nX (Int_t): number of variables nS (Int_t): number of observations xObs (Double_t*): the normalised observations (size = nX * nS) xNormParams (Double_t*): the normalisation parameters of xObs (size = 2 * nX) yObs (Double_t*); the outputs (size = _nS) gamma (Double_t*): parameters of the Kriging model (size = nS) iL (Double_t*): inverse of matrix L (size = nS * nS) gpVar (Double_t): variance of the Gaussian Process vErrMes (Double_t): variance of the measurement error corrFunc (TCorrelationFunction*): correlation function of the Gaussian Process varNames (TString ): name of the variables separated by ":". The last name is the output variable. nP (Int_t): number of parameters of the trend (default = 0) hbeta (Double_t *): parameters of the trend (size = nP) (default = NULL) M (Double_t*): matrice M (size = nS * nP) (default = NULL) iR (Double_t*): inverse de la matrice R (size = nP * nP) (default = NULL) trend (TString): string listing the trend coefficients (default = "")
- Note
- the TKriging object is copying every information provided. This means that the source of this information can be safely destroyed after creating the TKriging object.
Referenced by ClassImp().
◆ ~TKriging()
|
virtual |
Default destructor.
Referenced by ClassImp().
Member Function Documentation
◆ computeCovarianceMatrix()
void URANIE::Modeler::TKriging::computeCovarianceMatrix | ( | int | ns, |
double * | xObs, | ||
bool | bVariance_in_factor, | ||
bool | bHas_measurement_error, | ||
double | alpha, | ||
double | reg, | ||
double | sigma2, | ||
double * | Vmes, | ||
double * | C, | ||
double * | K | ||
) |
Referenced by ClassImp().
◆ createTrend()
void URANIE::Modeler::TKriging::createTrend | ( | ) |
Create the list of trend coefficients.
Referenced by ClassImp().
◆ estimate()
void URANIE::Modeler::TKriging::estimate | ( | URANIE::DataServer::TDataServer * | tdsQuery, |
TString | listOfInputs = "" , |
||
TString | listOfOutputs = "" , |
||
Option_t * | option = "" , |
||
Bool_t | useGPU = false |
||
) |
Referenced by ClassImp().
◆ estimate_CPU()
void URANIE::Modeler::TKriging::estimate_CPU | ( | URANIE::DataServer::TDataServer * | tdsQuery, |
TString | listOfInputs = "" , |
||
TString | listOfOutputs = "" , |
||
Option_t * | option = "" |
||
) |
Referenced by ClassImp().
◆ estimateWithCov()
void URANIE::Modeler::TKriging::estimateWithCov | ( | URANIE::DataServer::TDataServer * | tdsQuery, |
TString | listOfInputs = "" , |
||
TString | listOfOutputs = "" , |
||
TString | realValue = "" , |
||
Option_t * | option = "" |
||
) |
Referenced by ClassImp().
◆ eval()
Int_t URANIE::Modeler::TKriging::eval | ( | Double_t * | x0, |
Double_t * | y0, | ||
int | = 0 |
||
) |
Evaluate the gaussian process at the new location x0.
This function aims to be compatible with the TEval prototype, which will allow to call it from a launcher or an optimizer.
- Parameters
-
x0 (Double_t*): a single set of input variables values. y0 (Double_t*): output of the gaussian process for x0.
- Warning
- the size of x0 is not checked. If it is different from the expected number of input variables, this will lead to unexpected results or an error.
- Note
- Calls the gpPred0 function from the gpL
Referenced by ClassImp().
◆ exportFunction()
void URANIE::Modeler::TKriging::exportFunction | ( | const char * | lang, |
const char * | file, | ||
const char * | name, | ||
Option_t * | option = "" |
||
) |
Referenced by ClassImp().
◆ exportModelCplusplus()
|
virtual |
Export the model as a C function (not implemented)
Referenced by ClassImp().
◆ exportModelFortran()
|
inlinevirtual |
Export the model as a Fortran function (not implemented)
Referenced by ClassImp().
◆ exportModelPMML()
|
virtual |
Referenced by ClassImp().
◆ freeze()
void URANIE::Modeler::TKriging::freeze | ( | ) |
Referenced by ClassImp().
◆ getCorrFunction()
|
inline |
Return a pointer to the correlation function of the gaussian process.
References _correlationFunction.
◆ getCorrFunctionParams()
|
inline |
Return all the parameters of the correlation function.
- Warning
- this is an alias to TCorrelationFunction::getParameters(). This means that the correlation lengths are stored at the end of the array, and are normalised.
References _correlationFunction, and URANIE::Modeler::TCorrelationFunction::getParameters().
Referenced by ClassImp().
◆ getCorrLengths()
|
inline |
Return the correlation lengths.
- Note
- these are NOT the correlation lengths stored in the correlation function. These values are denormalised in order to be interpretable in the "real" input space.
References _corrLengths.
Referenced by ClassImp().
◆ getCorrLengthsNormalised()
|
inline |
Return the normalised correlation lengths.
- Note
- this corresponds to the correlation lengths stored in the correlation function.
References _normCorrLengths.
Referenced by ClassImp().
◆ getGPParams()
|
inline |
Return the array of the Gaussian Process parameters.
References _gamma.
◆ getiL()
|
inline |
Return the iL matrix as an array (cf. internal variable description)
References _iL.
◆ getInputNames()
|
inline |
Return the input variable names.
References _inputVarNames.
◆ getiR()
|
inline |
Return the iR matrix as an array (cf. internal variable description)
References _iR.
◆ getLooData()
URANIE::DataServer::TDataServer * URANIE::Modeler::TKriging::getLooData | ( | ) |
Return a pointer to a TDataServer containing the observations and the Leave One Out results.
- Note
- the user is responsible of the memory allocated for the data server.
Referenced by ClassImp().
◆ getLooErrors() [1/2]
|
inline |
Return the vector of the leave one out errors.
References _errorLoo.
◆ getLooErrors() [2/2]
|
inline |
◆ getLooQ2()
|
inline |
Return the Q2 of the model on the observation.
References _mseLoo, _nS, and _yObs.
Referenced by ClassImp().
◆ getLooRMSE()
|
inline |
◆ getLooVariances() [1/2]
|
inline |
Return the vector of the leave one out prediction variances.
References _varLoo.
◆ getLooVariances() [2/2]
|
inline |
◆ getM()
|
inline |
Return the M matrix as an array (cf. internal variable description)
References _M.
◆ getMeasurementErrorVariance()
|
inline |
Return the variance of the measurement error.
- Note
- if a value of -1.0 is returned, it means that the measurement error used to create the model was user defined. It was not considered useful to saved it with the model.
References _vErrMes.
◆ getNCorrFunctionParams()
|
inline |
Return the number of parameters of the correlation function.
References _correlationFunction, and URANIE::Modeler::TCorrelationFunction::getNParameters().
Referenced by ClassImp().
◆ getNCorrLengths()
|
inline |
Return the number of correlation lengths.
References _correlationFunction, and URANIE::Modeler::TCorrelationFunction::getNCorrLengths().
Referenced by ClassImp().
◆ getNInputs()
|
inline |
Return the number of input parameters.
References _nX.
◆ getNObs()
|
inline |
Return the number of observations used to build the GP model.
References _nS.
◆ getNTrendParams()
|
inline |
Return the number of parameters of the deterministic trend.
References _nP.
◆ getObsInputs()
|
inline |
Return the observations inputs as a flat array: [a1, b1, c1, a2, b2, c2,...].
References _xObs.
◆ getObsNormParams()
|
inline |
Return the normalisation parameters of the inputs as a flat array: [aMin, aMax, bMin, bMax,...].
References _xNormParams.
◆ getObsOutputs()
|
inline |
Return the observations outputs.
References _yObs.
◆ getOutputName()
|
inline |
Return the output variable name.
References _outputVarName.
◆ getTrendFormula()
|
inline |
Return the deterministic trend list of formulas.
References _trendCoefList.
◆ getTrendParams()
|
inline |
Return the array of the trend parameters.
References _hbeta.
◆ getTrendParamsCovMatrix()
|
inline |
Return the array of the trend parameters covariance matrix.
References _vbeta.
◆ getTrendString()
|
inline |
Return the deterministic trend string.
References _trendString.
◆ getVariance()
|
inline |
Return the variance of the gaussian process.
References _sigma2.
◆ hasTrend()
|
inline |
Return kTRUE if a deterministic trend has been defined.
References bHas_trend.
◆ initMatrices()
void URANIE::Modeler::TKriging::initMatrices | ( | ) |
Initialisation of the matrices.
Referenced by ClassImp().
◆ printLog()
|
virtual |
Referenced by ClassImp().
◆ setCorrFunction()
|
inline |
Set the correlation function of the gaussian process.
- Parameters
-
corrFunc (URANIE::Modeler::TCorrelationFunction*): a pointer to the correlation function of the gaussian process.
References _correlationFunction.
◆ setLooErrors()
void URANIE::Modeler::TKriging::setLooErrors | ( | Double_t * | err, |
Double_t * | var, | ||
Double_t | mse | ||
) |
Set the information about the leave one out error of the model.
This function is called by the GPBuilder if the user wants the leave one out error of the model. The errors are computed thanks to the gpLooXError functions of the gpLib.
- Parameters
-
err (Double_t*): vector of size _nS of the errors from Leave One Out. var (Double_t*): vector of size _nS of the analytical variances of the leave one out errors. mse (Double_t): leave one out mean squared error.
Referenced by ClassImp(), and ClassImp().
◆ setVariance()
|
inline |
Set the variance of the gaussian process.
- Parameters
-
var (Double_t): value of the variance.
References _sigma2.
◆ StatMat() [1/2]
|
static |
Referenced by ClassImp().
◆ StatMat() [2/2]
|
static |
◆ StatVect() [1/2]
|
static |
Referenced by ClassImp().
◆ StatVect() [2/2]
|
static |
Member Data Documentation
◆ _blog
|
private |
Boolean for edit the log.
Referenced by ClassImp().
◆ _correlationFunction
TCorrelationFunction* URANIE::Modeler::TKriging::_correlationFunction |
The correlation function.
Referenced by ClassImp(), getCorrFunction(), getCorrFunctionParams(), getNCorrFunctionParams(), getNCorrLengths(), and setCorrFunction().
◆ _corrLengths
Double_t* URANIE::Modeler::TKriging::_corrLengths |
The correlation lengths in the real input space.
Referenced by ClassImp(), and getCorrLengths().
◆ _errorLoo
Double_t* URANIE::Modeler::TKriging::_errorLoo |
If available, vector of the Leave One Out errors for the model.
Referenced by ClassImp(), getLooErrors(), and getLooErrors().
◆ _gamma
Double_t* URANIE::Modeler::TKriging::_gamma |
Array of parameters of the gaussian process (size: _nS)
Referenced by ClassImp(), and getGPParams().
◆ _hbeta
Double_t* URANIE::Modeler::TKriging::_hbeta |
Estimates of the deterministic trend parameters (size: _nP)
Referenced by ClassImp(), and getTrendParams().
◆ _iL
Double_t* URANIE::Modeler::TKriging::_iL |
Inverse of the matrix L, Cholesky decomposition of K (size: _nS * _nS)
Referenced by ClassImp(), and getiL().
◆ _inputNamesList
TObjArray* URANIE::Modeler::TKriging::_inputNamesList |
list of the inputs' names.
Referenced by ClassImp().
◆ _inputVarNames
TString URANIE::Modeler::TKriging::_inputVarNames |
Name of the inputs variables separated by ":". This information can be used to remind the user of the original names of the variables.
Referenced by ClassImp(), and getInputNames().
◆ _iR
Double_t* URANIE::Modeler::TKriging::_iR |
Inverse of the matrix R, from the QR decomposition of M (size: _nP * _nP)
Referenced by ClassImp(), and getiR().
◆ _M
Double_t* URANIE::Modeler::TKriging::_M |
Matrix M = L^{-T} H (size: _nS * _nP)
Referenced by ClassImp(), and getM().
◆ _mseLoo
Double_t URANIE::Modeler::TKriging::_mseLoo |
if available, Leave one Out mean squared error of the model.
Referenced by ClassImp(), getLooQ2(), and getLooRMSE().
◆ _normCorrLengths
Double_t* URANIE::Modeler::TKriging::_normCorrLengths |
The correlation lengths in the normalised input space. This corresponds to the correlation lengths stored in the correlation function.
Referenced by ClassImp(), and getCorrLengthsNormalised().
◆ _nP
Int_t URANIE::Modeler::TKriging::_nP |
Number of coefficients of the deterministic trend.
Referenced by ClassImp(), and getNTrendParams().
◆ _nS
Int_t URANIE::Modeler::TKriging::_nS |
Number of observations.
Referenced by ClassImp(), getLooErrors(), getLooQ2(), getLooVariances(), and getNObs().
◆ _nX
Int_t URANIE::Modeler::TKriging::_nX |
Number of input attributes
Referenced by ClassImp(), and getNInputs().
◆ _nXDeclared
Int_t URANIE::Modeler::TKriging::_nXDeclared |
Number of attribute added to the tkiging object before running it.
Referenced by ClassImp().
◆ _nYDeclared
Int_t URANIE::Modeler::TKriging::_nYDeclared |
Number of attribute added to the tkiging object before running it.
Referenced by ClassImp().
◆ _outputVarName
TString URANIE::Modeler::TKriging::_outputVarName |
Name of the output variable. This information can be used to remind the user of the original names of the variables.
Referenced by ClassImp(), and getOutputName().
◆ _sigma2
Double_t URANIE::Modeler::TKriging::_sigma2 |
Variance of the gaussian process.
Referenced by ClassImp(), getVariance(), and setVariance().
◆ _trendCoefList
TObjArray* URANIE::Modeler::TKriging::_trendCoefList |
List of the formulas of the trend coefficients (_nP elements)
Referenced by ClassImp(), and getTrendFormula().
◆ _trendString
TString URANIE::Modeler::TKriging::_trendString |
Character string containing the trend parameters separated by ":" or a trend type ("const" or "linear")
Referenced by ClassImp(), and getTrendString().
◆ _varLoo
Double_t* URANIE::Modeler::TKriging::_varLoo |
If available, vector of the variances of the Leave One Out errors for the model.
Referenced by ClassImp(), getLooVariances(), and getLooVariances().
◆ _vbeta
Double_t* URANIE::Modeler::TKriging::_vbeta |
Covariance matrix of the trend parameters (size: _nP*_nP)
Referenced by ClassImp(), and getTrendParamsCovMatrix().
◆ _vErrMes
Double_t URANIE::Modeler::TKriging::_vErrMes |
variance of the measurement error
Referenced by ClassImp(), and getMeasurementErrorVariance().
◆ _xNormParams
Double_t* URANIE::Modeler::TKriging::_xNormParams |
The matrix of normalisation parameters. Contains min and max value of each input variable: [aMin, aMax, bMin, bMax,...] (size: _nX * 2)
Referenced by ClassImp(), and getObsNormParams().
◆ _xObs
Double_t* URANIE::Modeler::TKriging::_xObs |
The matrix of inputs (dimension: _nX * _nS)
Referenced by ClassImp(), and getObsInputs().
◆ _yObs
Double_t* URANIE::Modeler::TKriging::_yObs |
Array representation of the outputs (size: _nS)
Referenced by ClassImp(), getLooQ2(), and getObsOutputs().
◆ bHas_trend
Bool_t URANIE::Modeler::TKriging::bHas_trend |
kFALSE if no deterministic trend is defined, kTRUE otherwise
Referenced by ClassImp(), and hasTrend().
◆ bUse_normalisation
Bool_t URANIE::Modeler::TKriging::bUse_normalisation |
if kFALSE, the matrix H will be computed with unormalised input data. This happens when the user is providing his own trend function. When the trend is constructed automatically from a trend type, this attribute is set to kTRUE.
Referenced by ClassImp().