English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Modeler: TKriging.h Source File
Uranie / Modeler  v4.10.0
/* @license-end */
TKriging.h
Go to the documentation of this file.
1 // Copyright (C) 2013-2024 CEA/DES
3 //
4 // This program is free software: you can redistribute it and/or modify
5 // it under the terms of the GNU Lesser General Public License as published
6 // by the Free Software Foundation, either version 3 of the License, or any
7 // later version.
8 //
9 // This program is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 // GNU Lesser General Public License for more details.
13 //
14 // You should have received a copy of the GNU Lesser General Public License
15 // along with this program. If not, see <http://www.gnu.org/licenses/>.
18 // $Id$
19 // $Author$
20 // $Date$
21 // $Revision$
22 // $State$
24 
44 #ifndef TKRIGING_H
45 #define TKRIGING_H
46 
47 #include <cmath>
48 
50 #include "TObjArray.h"
51 #include "TFormula.h"
52 #include "TMath.h"
53 
55 #include "TCorrelationFunction.h"
56 #include "TStandardEval.h"
57 #include "TDataServer.h"
58 #include "TAttribute.h"
59 #include "TDataServer.h"
60 
61 #include <memory>
62 
64 extern "C"
65 {
66  void gpPred0(double *Kp, double variance, double *gamma, double *iL,
67  int nnew, int n, double *yp, double *vp);
68  void gpPred0Cov(double *Kp, double *Rp, double variance, double *gamma, double *iL,
69  int nnew, int n, double *yp, double *vp);
70  void gpPredh(double *Hp, double *Kp, double variance, double *hbeta, double *gamma, double *iL, double *iR, double *M,
71  int nnew, int n, int p, double *yp, double *vp);
72  void gpPredhCov(double *Hp, double *Kp, double *Rp, double variance, double *hbeta, double *gamma, double *iL, double *iR, double *M,
73  int nnew, int n, int p, double *yp, double *vp);
74 }
75 ;
76 
77 namespace URANIE
78 {
79 namespace Modeler
80 {
81 class TKriging: public URANIE::Relauncher::TDoubleEval
82 {
83 
84  // Attributes
85 private:
86  Bool_t _blog;
87 public:
88  TString _outputVarName;
89  TString _inputVarNames;
90  TObjArray* _inputNamesList;
91 
92  Int_t _nX;
93  Int_t _nS;
94  Int_t _nYDeclared;
95  Int_t _nXDeclared;
96  Double_t* _xObs;
97  Double_t* _xNormParams;
98  Double_t* _yObs;
99  Double_t* _iL;
100  Double_t* _gamma;
101 
102  Double_t _sigma2;
103  Double_t _vErrMes;
105  Double_t* _corrLengths;
106  Double_t* _normCorrLengths;
107 
109  Bool_t bHas_trend;
110  TString _trendString;
111  vector<TFormula*> _trendCoefList;
112  vector<vector<int>> _trendIndex;
113  std::unique_ptr<TObjArray> _inpname;
114  Int_t _nP;
115  Double_t* _hbeta;
116  Double_t* _vbeta;
117  Double_t* _M;
118  Double_t* _iR;
119 
120  Double_t* _errorLoo;
121  Double_t* _varLoo;
122  Double_t _mseLoo;
123 
124  // Operations
125 public:
126  //---------------------------------------------
130 
154  TKriging(Int_t nX, Int_t nS, Double_t* xObs, Double_t* xNormParams, Double_t* yObs, Double_t* gamma, Double_t* iL,
155  Double_t gpVar, Double_t vErrMes, URANIE::Modeler::TCorrelationFunction* corrFunc, TString varNames,
156  Int_t nP = 0, Double_t* hbeta = NULL, Double_t* vbeta= NULL, Double_t* M = NULL, Double_t* iR = NULL, TString trend = "");
157 
159 #ifdef WITH_LIBXML2
160 
166  TKriging(TString filename, TString modelName);
167 #endif
168  void initMatrices();
170 
172  virtual ~TKriging();
174 
175  //---------------------------------------------
179  const char* getInputNames()
181  {
182  return _inputVarNames.Data();
183  }
184 
186  const char* getOutputName()
187  {
188  return _outputVarName.Data();
189  }
190 
192  Bool_t hasTrend()
193  {
194  return bHas_trend;
195  }
196 
199  {
200  return _correlationFunction;
201  }
202 
204  Int_t getNInputs()
205  {
206  return _nX;
207  }
208 
210  Int_t getNObs()
211  {
212  return _nS;
213  }
214 
217  {
219  }
220 
223  {
225  }
226 
229  {
230  return _nP;
231  }
232 
234  Double_t* getObsInputs()
235  {
236  return _xObs;
237  }
238 
240  Double_t* getObsOutputs()
241  {
242  return _yObs;
243  }
244 
246  Double_t* getObsNormParams()
247  {
248  return _xNormParams;
249  }
250 
252  Double_t* getGPParams()
253  {
254  return _gamma;
255  }
256 
258  Double_t* getTrendParams()
259  {
260  return _hbeta;
261  }
262 
265  {
266  return _vbeta;
267  }
268 
270  vector<TFormula*> getTrendFormula()
271  {
272  return _trendCoefList;
273  }
274 
276  TString getTrendString()
277  {
278  return _trendString;
279  }
280 
282  Double_t getVariance()
283  {
284  return _sigma2;
285  }
286 
288 
294  {
295  return _vErrMes;
296  }
297 
299  Double_t* getiL()
300  {
301  return _iL;
302  }
303 
305  Double_t* getM()
306  {
307  return _M;
308  }
309 
311  Double_t* getiR()
312  {
313  return _iR;
314  }
315 
316 
318  Double_t* getLooErrors()
319  {
320  return _errorLoo;
321  }
322 
324  Double_t* getLooVariances()
325  {
326  return _varLoo;
327  }
328 
330  void getLooErrors(double *arr, int size){
331  try{
332  if(size!=_nS)
333  throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
334  Form("TKriging::getLooErrors method\n The size of provided array [%d] is not the size of the output to be dumped [%d]",
335  size, _nS));
336  } catch (URANIE::Exceptions::UErrorExceptions& ue)
337  {
338  ue.printMessage();
339  throw ue;
340  }
341 
342  for(int ival=0; ival<size; ival++)
343  arr[ival]=_errorLoo[ival];
344 
345  }
346 
348  void getLooVariances(double *arr, int size){
349  try{
350  if(size!=_nS)
351  throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
352  Form("TKriging::getLooErrors method\n The size of provided array [%d] is not the size of the output to be dumped [%d]",
353  size, _nS));
354  } catch (URANIE::Exceptions::UErrorExceptions& ue)
355  {
356  ue.printMessage();
357  throw ue;
358  }
359 
360  for(int ival=0; ival<size; ival++)
361  arr[ival]=_varLoo[ival];
362  }
364  Double_t getLooRMSE()
365  {
366  return sqrt(_mseLoo);
367  }
368 
370  Double_t getLooQ2()
371  {
372  Double_t mean = TMath::Mean(_nS, _yObs);
373  Double_t sumVar = 0.0;
374  for(Int_t i = 0; i < _nS; i++) sumVar += (_yObs[i]-mean)*(_yObs[i]-mean);
375  if(_nS > 1)
376  {
377  Double_t varYObs = sumVar/(_nS-1);
378  return 1.0 - _mseLoo/varYObs;
379  }
380  else
381  return 1.0 - _mseLoo/sumVar;
382  }
383 
385 
388  URANIE::DataServer::TDataServer* getLooData();
389 
391 
395  {
396  return _normCorrLengths;
397  }
398 
400 
404  Double_t* getCorrLengths()
405  {
406  return _corrLengths;
407  }
408 
410 
416  {
418  }
419 
420  //---------------------------------------------
424 
429  {
430  _correlationFunction = corrFunc;
431  }
432 
434 
437  void setVariance(Double_t var)
438  {
439  _sigma2 = var;
440  }
441 
443 
449  void setLooErrors(Double_t *err, Double_t *var, Double_t mse);
450 
452  void createTrend();
453 
454  static void StatMat(std::string title, int n1, int n2, Double_t *v);
455  static void StatMat(std::string title, int n1, int n2, Double_t *v,
456  std::vector<double> &tvmin, std::vector<double> &tvmax,
457  std::vector<double> &tmean, std::vector<double> &tvar,
458  std::vector<double> &tnorm);
459  static void StatVect(std::string title, int n, Double_t *v);
460  static void StatVect(std::string title, int n, Double_t *v, double &vmin, double &vmax, double &mean, double &var, double &norm);
461  void computeCovarianceMatrix(int ns, double *xObs,
462  bool bVariance_in_factor,
463  bool bHas_measurement_error, double alpha, double reg,
464  double sigma2, double *Vmes,
465  double *C, double *K);
466 
468 
469  //---------------------------------------------
473 
475 
485  Int_t eval(Double_t* x0, Double_t* y0, int=0);
486 
487 #ifdef WITH_GPUANN
488  void estimate(URANIE::DataServer::TDataServer *tdsQuery,
489  TString listOfInputs = "",
490  TString listOfOutputs = "",
491  Option_t * option = "",
492  Bool_t useGPU = true);
493  void estimate_GPU(URANIE::DataServer::TDataServer *tdsQuery,
494  TString listOfInputs = "",
495  TString listOfOutputs = "",
496  Option_t * option = "");
497 #else
498  void estimate(URANIE::DataServer::TDataServer *tdsQuery,
499  TString listOfInputs = "",
500  TString listOfOutputs = "",
501  Option_t * option = "",
502  Bool_t useGPU = false);
503 #endif
504  void estimate_CPU(URANIE::DataServer::TDataServer *tdsQuery,
505  TString listOfInputs = "",
506  TString listOfOutputs = "",
507  Option_t * option = "");
508 
509  void estimateWithCov(URANIE::DataServer::TDataServer *tdsQuery,
510  TString listOfInputs = "",
511  TString listOfOutputs = "",
512  TString realValue = "",
513  Option_t * option = "");
515 
516  //---------------------------------------------
524 
525  void exportFunction(const char* lang, const char* file, const char* name, Option_t *option="");
526 
528  virtual void exportModelCplusplus(std::ofstream *sourcefile, const char* name = "");
529 
531  virtual void exportModelFortran(std::ofstream *sourcefile)
532  {
533  UNUSED(sourcefile);
534  std::cout << "TKriging::exportModelFortran not yet implemented" << std::endl;
535  }
536  virtual void exportModelPMML(const char* file, const char* name, Option_t *option) const;
537 
538  void freeze();
540 
541  //---------------------------------------------
545  virtual void printLog(Option_t *option = "");
547 
548  ClassDef(URANIE::Modeler::TKriging, ID_MODELER)
549  //Classe de
550 };
551 
552 } // Fin du namespace Modeler
553 } // Fin du namespace URANIE
554 #endif
ROOT.
Definition: TAnisp.h:163
const char * getOutputName()
Return the output variable name.
Definition: TKriging.h:186
Double_t getLooQ2()
Return the Q2 of the model on the observation.
Definition: TKriging.h:370
Double_t getVariance()
Return the variance of the gaussian process.
Definition: TKriging.h:282
virtual void exportModelFortran(std::ofstream *sourcefile)
Export the model as a Fortran function (not implemented)
Definition: TKriging.h:531
void getLooVariances(double *arr, int size)
Return the vector of the leave one out prediction variances.
Definition: TKriging.h:348
Double_t * _gamma
Array of parameters of the gaussian process (size: _nS)
Definition: TKriging.h:100
virtual void printLog(Option_t *option="")
Bool_t bUse_normalisation
if kFALSE, the matrix H will be computed with unormalised input data. This happens when the user is p...
Definition: TKriging.h:108
void estimate(URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", Option_t *option="", Bool_t useGPU=false)
Double_t * getObsInputs()
Return the observations inputs as a flat array: [a1, b1, c1, a2, b2, c2,...].
Definition: TKriging.h:234
static void StatMat(std::string title, int n1, int n2, Double_t *v)
Int_t _nP
Number of coefficients of the deterministic trend.
Definition: TKriging.h:114
virtual ~TKriging()
Default destructor.
TString _inputVarNames
Name of the inputs variables separated by ":". This information can be used to remind the user of the...
Definition: TKriging.h:89
Double_t getMeasurementErrorVariance()
Return the variance of the measurement error.
Definition: TKriging.h:293
Int_t getNParameters()
Return the number of parameters of the correlation function.
Definition: TCorrelationFunction.h:85
Double_t * getM()
Return the M matrix as an array (cf. internal variable description)
Definition: TKriging.h:305
virtual void exportModelCplusplus(std::ofstream *sourcefile, const char *name="")
Export the model as a C function (not implemented)
Double_t * getiL()
Return the iL matrix as an array (cf. internal variable description)
Definition: TKriging.h:299
Int_t _nXDeclared
Number of attribute added to the tkiging object before running it.
Definition: TKriging.h:95
const char * getInputNames()
Return the input variable names.
Definition: TKriging.h:180
Int_t _nS
Number of observations.
Definition: TKriging.h:93
Double_t * getCorrLengths()
Return the correlation lengths.
Definition: TKriging.h:404
void gpPredh(double *Hp, double *Kp, double variance, double *hbeta, double *gamma, double *iL, double *iR, double *M, int nnew, int n, int p, double *yp, double *vp)
Double_t * _varLoo
If available, vector of the variances of the Leave One Out errors for the model.
Definition: TKriging.h:121
void gpPred0Cov(double *Kp, double *Rp, double variance, double *gamma, double *iL, int nnew, int n, double *yp, double *vp)
Bool_t hasTrend()
Return kTRUE if a deterministic trend has been defined.
Definition: TKriging.h:192
TString getTrendString()
Return the deterministic trend string.
Definition: TKriging.h:276
TCorrelationFunction * _correlationFunction
The correlation function.
Definition: TKriging.h:104
void estimate_CPU(URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", Option_t *option="")
Double_t * getParameters()
Returns an array containing the parameters of the function.
Definition: TCorrelationFunction.h:103
Double_t * getLooErrors()
Return the vector of the leave one out errors.
Definition: TKriging.h:318
Double_t _vErrMes
variance of the measurement error
Definition: TKriging.h:103
static void StatVect(std::string title, int n, Double_t *v)
Double_t * _hbeta
Estimates of the deterministic trend parameters (size: _nP)
Definition: TKriging.h:115
void setLooErrors(Double_t *err, Double_t *var, Double_t mse)
Set the information about the leave one out error of the model.
Int_t eval(Double_t *x0, Double_t *y0, int=0)
Evaluate the gaussian process at the new location x0.
void estimateWithCov(URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", TString realValue="", Option_t *option="")
Double_t _sigma2
Variance of the gaussian process.
Definition: TKriging.h:102
void getLooErrors(double *arr, int size)
Return the vector of the leave one out errors.
Definition: TKriging.h:330
Double_t * _iL
Inverse of the matrix L, Cholesky decomposition of K (size: _nS * _nS)
Definition: TKriging.h:99
Int_t getNInputs()
Return the number of input parameters.
Definition: TKriging.h:204
Description of the class TCorrelationFunction.
Definition: TCorrelationFunction.h:52
Double_t * getTrendParamsCovMatrix()
Return the array of the trend parameters covariance matrix.
Definition: TKriging.h:264
vector< TFormula * > getTrendFormula()
Return the deterministic trend list of formulas.
Definition: TKriging.h:270
TObjArray * _inputNamesList
list of the inputs&#39; names.
Definition: TKriging.h:90
Double_t * getCorrFunctionParams()
Return all the parameters of the correlation function.
Definition: TKriging.h:415
Double_t * _corrLengths
The correlation lengths in the real input space.
Definition: TKriging.h:105
Int_t getNCorrFunctionParams()
Return the number of parameters of the correlation function.
Definition: TKriging.h:216
void setVariance(Double_t var)
Set the variance of the gaussian process.
Definition: TKriging.h:437
TString _outputVarName
Name of the output variable. This information can be used to remind the user of the original names of...
Definition: TKriging.h:88
Double_t * _iR
Inverse of the matrix R, from the QR decomposition of M (size: _nP * _nP)
Definition: TKriging.h:118
void gpPredhCov(double *Hp, double *Kp, double *Rp, double variance, double *hbeta, double *gamma, double *iL, double *iR, double *M, int nnew, int n, int p, double *yp, double *vp)
Int_t _nYDeclared
Number of attribute added to the tkiging object before running it.
Definition: TKriging.h:94
void setCorrFunction(TCorrelationFunction *corrFunc)
Set the correlation function of the gaussian process.
Definition: TKriging.h:428
void gpPred0(double *Kp, double variance, double *gamma, double *iL, int nnew, int n, double *yp, double *vp)
void createTrend()
Create the list of trend coefficients.
TString _trendString
Character string containing the trend parameters separated by ":" or a trend type ("const" or "linear...
Definition: TKriging.h:110
Bool_t bHas_trend
kFALSE if no deterministic trend is defined, kTRUE otherwise
Definition: TKriging.h:109
URANIE::DataServer::TDataServer * getLooData()
Return a pointer to a TDataServer containing the observations and the Leave One Out results...
Double_t * _yObs
Array representation of the outputs (size: _nS)
Definition: TKriging.h:98
Int_t getNObs()
Return the number of observations used to build the GP model.
Definition: TKriging.h:210
Double_t * _errorLoo
If available, vector of the Leave One Out errors for the model.
Definition: TKriging.h:120
Int_t getNTrendParams()
Return the number of parameters of the deterministic trend.
Definition: TKriging.h:228
Interface de la classe URANIE::Modeler::TCorrelationFunction.
Double_t * getTrendParams()
Return the array of the trend parameters.
Definition: TKriging.h:258
Int_t getNCorrLengths()
Return the number of correlation lengths.
Definition: TCorrelationFunction.h:97
Double_t * _xObs
The matrix of inputs (dimension: _nX * _nS)
Definition: TKriging.h:96
Double_t * _vbeta
Covariance matrix of the trend parameters (size: _nP*_nP)
Definition: TKriging.h:116
std::unique_ptr< TObjArray > _inpname
Definition: TKriging.h:113
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.
Double_t * getiR()
Return the iR matrix as an array (cf. internal variable description)
Definition: TKriging.h:311
vector< TFormula * > _trendCoefList
List of the formulas of the trend coefficients (_nP elements)
Definition: TKriging.h:111
virtual void exportModelPMML(const char *file, const char *name, Option_t *option) const
TCorrelationFunction * getCorrFunction()
Return a pointer to the correlation function of the gaussian process.
Definition: TKriging.h:198
Double_t * _normCorrLengths
The correlation lengths in the normalised input space. This corresponds to the correlation lengths st...
Definition: TKriging.h:106
Description of the class TKriging.
Definition: TKriging.h:81
Double_t * getLooVariances()
Return the vector of the leave one out prediction variances.
Definition: TKriging.h:324
Double_t * getGPParams()
Return the array of the Gaussian Process parameters.
Definition: TKriging.h:252
Double_t * getObsNormParams()
Return the normalisation parameters of the inputs as a flat array: [aMin, aMax, bMin, bMax,...].
Definition: TKriging.h:246
void exportFunction(const char *lang, const char *file, const char *name, Option_t *option="")
void initMatrices()
Initialisation of the matrices.
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)
Double_t _mseLoo
if available, Leave one Out mean squared error of the model.
Definition: TKriging.h:122
Int_t getNCorrLengths()
Return the number of correlation lengths.
Definition: TKriging.h:222
Double_t * getObsOutputs()
Return the observations outputs.
Definition: TKriging.h:240
vector< vector< int > > _trendIndex
Definition: TKriging.h:112
Double_t * _xNormParams
The matrix of normalisation parameters. Contains min and max value of each input variable: [aMin...
Definition: TKriging.h:97
Int_t _nX
Number of input attributes.
Definition: TKriging.h:92
Bool_t _blog
Boolean for edit the log.
Definition: TKriging.h:86
Double_t * getCorrLengthsNormalised()
Return the normalised correlation lengths.
Definition: TKriging.h:394
Double_t getLooRMSE()
Return the leave one out mean squared error.
Definition: TKriging.h:364
Double_t * _M
Matrix M = L^{-T} H (size: _nS * _nP)
Definition: TKriging.h:117