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.9.0
/* @license-end */
TKriging.h
Go to the documentation of this file.
1
2// 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 "TMath.h"
52
55#include "TStandardEval.h"
56#include "TDataServer.h"
57#include "TAttribute.h"
58#include "TDataServer.h"
59
61extern "C"
62{
63 void gpPred0(double *Kp, double variance, double *gamma, double *iL,
64 int nnew, int n, double *yp, double *vp);
65 void gpPred0Cov(double *Kp, double *Rp, double variance, double *gamma, double *iL,
66 int nnew, int n, double *yp, double *vp);
67 void gpPredh(double *Hp, double *Kp, double variance, double *hbeta, double *gamma, double *iL, double *iR, double *M,
68 int nnew, int n, int p, double *yp, double *vp);
69 void gpPredhCov(double *Hp, double *Kp, double *Rp, double variance, double *hbeta, double *gamma, double *iL, double *iR, double *M,
70 int nnew, int n, int p, double *yp, double *vp);
71}
72;
73
74namespace URANIE
75{
76namespace Modeler
77{
78class TKriging: public URANIE::Relauncher::TDoubleEval
79{
80
81 // Attributes
82private:
83 Bool_t _blog;
84public:
87 TObjArray* _inputNamesList;
88
89 Int_t _nX;
90 Int_t _nS;
93 Double_t* _xObs;
94 Double_t* _xNormParams;
95 Double_t* _yObs;
96 Double_t* _iL;
97 Double_t* _gamma;
98
99 Double_t _sigma2;
100 Double_t _vErrMes;
102 Double_t* _corrLengths;
104
106 Bool_t bHas_trend;
107 TString _trendString;
108 TObjArray* _trendCoefList;
109 Int_t _nP;
110 Double_t* _hbeta;
111 Double_t* _vbeta;
112 Double_t* _M;
113 Double_t* _iR;
114
115 Double_t* _errorLoo;
116 Double_t* _varLoo;
117 Double_t _mseLoo;
118
119 // Operations
120public:
121 //---------------------------------------------
125
149 TKriging(Int_t nX, Int_t nS, Double_t* xObs, Double_t* xNormParams, Double_t* yObs, Double_t* gamma, Double_t* iL,
150 Double_t gpVar, Double_t vErrMes, URANIE::Modeler::TCorrelationFunction* corrFunc, TString varNames,
151 Int_t nP = 0, Double_t* hbeta = NULL, Double_t* vbeta= NULL, Double_t* M = NULL, Double_t* iR = NULL, TString trend = "");
152
154#ifdef WITH_LIBXML2
156
161 TKriging(TString filename, TString modelName);
162#endif
165
167 virtual ~TKriging();
169
170 //---------------------------------------------
174
175 const char* getInputNames()
176 {
177 return _inputVarNames.Data();
178 }
179
181 const char* getOutputName()
182 {
183 return _outputVarName.Data();
184 }
185
187 Bool_t hasTrend()
188 {
189 return bHas_trend;
190 }
191
197
200 {
201 return _nX;
202 }
203
205 Int_t getNObs()
206 {
207 return _nS;
208 }
209
215
218 {
220 }
221
224 {
225 return _nP;
226 }
227
229 Double_t* getObsInputs()
230 {
231 return _xObs;
232 }
233
235 Double_t* getObsOutputs()
236 {
237 return _yObs;
238 }
239
242 {
243 return _xNormParams;
244 }
245
247 Double_t* getGPParams()
248 {
249 return _gamma;
250 }
251
253 Double_t* getTrendParams()
254 {
255 return _hbeta;
256 }
257
260 {
261 return _vbeta;
262 }
263
265 TObjArray* getTrendFormula()
266 {
267 return _trendCoefList;
268 }
269
272 {
273 return _trendString;
274 }
275
277 Double_t getVariance()
278 {
279 return _sigma2;
280 }
281
283
289 {
290 return _vErrMes;
291 }
292
294 Double_t* getiL()
295 {
296 return _iL;
297 }
298
300 Double_t* getM()
301 {
302 return _M;
303 }
304
306 Double_t* getiR()
307 {
308 return _iR;
309 }
310
311
313 Double_t* getLooErrors()
314 {
315 return _errorLoo;
316 }
317
319 Double_t* getLooVariances()
320 {
321 return _varLoo;
322 }
323
325 void getLooErrors(double *arr, int size){
326 try{
327 if(size!=_nS)
328 throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
329 Form("TKriging::getLooErrors method\n The size of provided array [%d] is not the size of the output to be dumped [%d]",
330 size, _nS));
331 } catch (URANIE::Exceptions::UErrorExceptions& ue)
332 {
333 ue.printMessage();
334 throw ue;
335 }
336
337 for(int ival=0; ival<size; ival++)
338 arr[ival]=_errorLoo[ival];
339
340 }
341
343 void getLooVariances(double *arr, int size){
344 try{
345 if(size!=_nS)
346 throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
347 Form("TKriging::getLooErrors method\n The size of provided array [%d] is not the size of the output to be dumped [%d]",
348 size, _nS));
349 } catch (URANIE::Exceptions::UErrorExceptions& ue)
350 {
351 ue.printMessage();
352 throw ue;
353 }
354
355 for(int ival=0; ival<size; ival++)
356 arr[ival]=_varLoo[ival];
357 }
359 Double_t getLooRMSE()
360 {
361 return sqrt(_mseLoo);
362 }
363
365 Double_t getLooQ2()
366 {
367 Double_t mean = TMath::Mean(_nS, _yObs);
368 Double_t sumVar = 0.0;
369 for(Int_t i = 0; i < _nS; i++) sumVar += (_yObs[i]-mean)*(_yObs[i]-mean);
370 if(_nS > 1)
371 {
372 Double_t varYObs = sumVar/(_nS-1);
373 return 1.0 - _mseLoo/varYObs;
374 }
375 else
376 return 1.0 - _mseLoo/sumVar;
377 }
378
380
383 URANIE::DataServer::TDataServer* getLooData();
384
386
390 {
391 return _normCorrLengths;
392 }
393
395
399 Double_t* getCorrLengths()
400 {
401 return _corrLengths;
402 }
403
405
411 {
413 }
414
415 //---------------------------------------------
419
424 {
425 _correlationFunction = corrFunc;
426 }
427
429
432 void setVariance(Double_t var)
433 {
434 _sigma2 = var;
435 }
436
438
444 void setLooErrors(Double_t *err, Double_t *var, Double_t mse);
445
448
449 static void StatMat(std::string title, int n1, int n2, Double_t *v);
450 static void StatMat(std::string title, int n1, int n2, Double_t *v,
451 std::vector<double> &tvmin, std::vector<double> &tvmax,
452 std::vector<double> &tmean, std::vector<double> &tvar,
453 std::vector<double> &tnorm);
454 static void StatVect(std::string title, int n, Double_t *v);
455 static void StatVect(std::string title, int n, Double_t *v, double &vmin, double &vmax, double &mean, double &var, double &norm);
456 void computeCovarianceMatrix(int ns, double *xObs,
457 bool bVariance_in_factor,
458 bool bHas_measurement_error, double alpha, double reg,
459 double sigma2, double *Vmes,
460 double *C, double *K);
461
463
464 //---------------------------------------------
468
470
480 Int_t eval(Double_t* x0, Double_t* y0, int=0);
481
482#ifdef WITH_GPUANN
483 void estimate(URANIE::DataServer::TDataServer *tdsQuery,
484 TString listOfInputs = "",
485 TString listOfOutputs = "",
486 Option_t * option = "",
487 Bool_t useGPU = true);
488 void estimate_GPU(URANIE::DataServer::TDataServer *tdsQuery,
489 TString listOfInputs = "",
490 TString listOfOutputs = "",
491 Option_t * option = "");
492#else
493 void estimate(URANIE::DataServer::TDataServer *tdsQuery,
494 TString listOfInputs = "",
495 TString listOfOutputs = "",
496 Option_t * option = "",
497 Bool_t useGPU = false);
498#endif
499 void estimate_CPU(URANIE::DataServer::TDataServer *tdsQuery,
500 TString listOfInputs = "",
501 TString listOfOutputs = "",
502 Option_t * option = "");
503
504 void estimateWithCov(URANIE::DataServer::TDataServer *tdsQuery,
505 TString listOfInputs = "",
506 TString listOfOutputs = "",
507 TString realValue = "",
508 Option_t * option = "");
510
511 //---------------------------------------------
519
520 void exportFunction(const char* lang, const char* file, const char* name, Option_t *option="");
521
523 virtual void exportModelCplusplus(std::ofstream *sourcefile, const char* name = "");
524
526 virtual void exportModelFortran(std::ofstream *sourcefile)
527 {
528 UNUSED(sourcefile);
529 std::cout << "TKriging::exportModelFortran not yet implemented" << std::endl;
530 }
531 virtual void exportModelPMML(const char* file, const char* name, Option_t *option) const;
532
533 void freeze();
535
536 //---------------------------------------------
540 virtual void printLog(Option_t *option = "");
542
543 ClassDef(URANIE::Modeler::TKriging, ID_MODELER)
544 //Classe de
545};
546
547} // Fin du namespace Modeler
548} // Fin du namespace URANIE
549#endif
Interface de la classe URANIE::Modeler::TCorrelationFunction.
void gpPred0(double *Kp, double variance, double *gamma, double *iL, int nnew, int n, double *yp, double *vp)
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)
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)
void gpPred0Cov(double *Kp, double *Rp, double variance, double *gamma, double *iL, int nnew, int n, double *yp, double *vp)
Description of the class TCorrelationFunction.
Definition TCorrelationFunction.h:53
Int_t getNParameters()
Return the number of parameters of the correlation function.
Definition TCorrelationFunction.h:85
Double_t * getParameters()
Returns an array containing the parameters of the function.
Definition TCorrelationFunction.h:103
Int_t getNCorrLengths()
Return the number of correlation lengths.
Definition TCorrelationFunction.h:97
Description of the class TKriging.
Definition TKriging.h:79
Int_t _nYDeclared
Number of attribute added to the tkiging object before running it.
Definition TKriging.h:91
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:105
Double_t * _xNormParams
The matrix of normalisation parameters. Contains min and max value of each input variable: [aMin,...
Definition TKriging.h:94
void getLooErrors(double *arr, int size)
Return the vector of the leave one out errors.
Definition TKriging.h:325
void estimateWithCov(URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", TString realValue="", Option_t *option="")
Int_t _nS
Number of observations.
Definition TKriging.h:90
void initMatrices()
Initialisation of the matrices.
TObjArray * _trendCoefList
List of the formulas of the trend coefficients (_nP elements)
Definition TKriging.h:108
Int_t _nX
Number of input attributes
Definition TKriging.h:89
Double_t _mseLoo
if available, Leave one Out mean squared error of the model.
Definition TKriging.h:117
void estimate(URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", Option_t *option="", Bool_t useGPU=false)
Int_t getNObs()
Return the number of observations used to build the GP model.
Definition TKriging.h:205
Double_t getVariance()
Return the variance of the gaussian process.
Definition TKriging.h:277
virtual void exportModelFortran(std::ofstream *sourcefile)
Export the model as a Fortran function (not implemented)
Definition TKriging.h:526
TCorrelationFunction * getCorrFunction()
Return a pointer to the correlation function of the gaussian process.
Definition TKriging.h:193
URANIE::DataServer::TDataServer * getLooData()
Return a pointer to a TDataServer containing the observations and the Leave One Out results.
Double_t * _iL
Inverse of the matrix L, Cholesky decomposition of K (size: _nS * _nS)
Definition TKriging.h:96
void createTrend()
Create the list of trend coefficients.
Double_t * _gamma
Array of parameters of the gaussian process (size: _nS)
Definition TKriging.h:97
Bool_t bHas_trend
kFALSE if no deterministic trend is defined, kTRUE otherwise
Definition TKriging.h:106
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)
const char * getInputNames()
Return the input variable names.
Definition TKriging.h:175
Double_t * _errorLoo
If available, vector of the Leave One Out errors for the model.
Definition TKriging.h:115
virtual ~TKriging()
Default destructor.
Double_t * getTrendParamsCovMatrix()
Return the array of the trend parameters covariance matrix.
Definition TKriging.h:259
Double_t * getCorrLengthsNormalised()
Return the normalised correlation lengths.
Definition TKriging.h:389
Double_t getLooRMSE()
Return the leave one out mean squared error.
Definition TKriging.h:359
Double_t * getObsOutputs()
Return the observations outputs.
Definition TKriging.h:235
Double_t _vErrMes
variance of the measurement error
Definition TKriging.h:100
Double_t * getLooErrors()
Return the vector of the leave one out errors.
Definition TKriging.h:313
Double_t * getiL()
Return the iL matrix as an array (cf. internal variable description)
Definition TKriging.h:294
void setVariance(Double_t var)
Set the variance of the gaussian process.
Definition TKriging.h:432
Double_t * getiR()
Return the iR matrix as an array (cf. internal variable description)
Definition TKriging.h:306
Double_t * getLooVariances()
Return the vector of the leave one out prediction variances.
Definition TKriging.h:319
void estimate_CPU(URANIE::DataServer::TDataServer *tdsQuery, TString listOfInputs="", TString listOfOutputs="", Option_t *option="")
Double_t * _iR
Inverse of the matrix R, from the QR decomposition of M (size: _nP * _nP)
Definition TKriging.h:113
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)
TObjArray * _inputNamesList
list of the inputs' names.
Definition TKriging.h:87
Double_t * _yObs
Array representation of the outputs (size: _nS)
Definition TKriging.h:95
void getLooVariances(double *arr, int size)
Return the vector of the leave one out prediction variances.
Definition TKriging.h:343
TString getTrendString()
Return the deterministic trend string.
Definition TKriging.h:271
Double_t * getObsInputs()
Return the observations inputs as a flat array: [a1, b1, c1, a2, b2, c2,...].
Definition TKriging.h:229
Double_t * getObsNormParams()
Return the normalisation parameters of the inputs as a flat array: [aMin, aMax, bMin,...
Definition TKriging.h:241
static void StatMat(std::string title, int n1, int n2, Double_t *v)
Double_t * _corrLengths
The correlation lengths in the real input space.
Definition TKriging.h:102
Double_t * _M
Matrix M = L^{-T} H (size: _nS * _nP)
Definition TKriging.h:112
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.
Int_t getNCorrFunctionParams()
Return the number of parameters of the correlation function.
Definition TKriging.h:211
Double_t * _varLoo
If available, vector of the variances of the Leave One Out errors for the model.
Definition TKriging.h:116
Double_t * getCorrLengths()
Return the correlation lengths.
Definition TKriging.h:399
Int_t _nXDeclared
Number of attribute added to the tkiging object before running it.
Definition TKriging.h:92
TObjArray * getTrendFormula()
Return the deterministic trend list of formulas.
Definition TKriging.h:265
Double_t * getGPParams()
Return the array of the Gaussian Process parameters.
Definition TKriging.h:247
Double_t _sigma2
Variance of the gaussian process.
Definition TKriging.h:99
Double_t getLooQ2()
Return the Q2 of the model on the observation.
Definition TKriging.h:365
const char * getOutputName()
Return the output variable name.
Definition TKriging.h:181
Double_t * getCorrFunctionParams()
Return all the parameters of the correlation function.
Definition TKriging.h:410
TString _trendString
Character string containing the trend parameters separated by ":" or a trend type ("const" or "linear...
Definition TKriging.h:107
void setCorrFunction(TCorrelationFunction *corrFunc)
Set the correlation function of the gaussian process.
Definition TKriging.h:423
virtual void printLog(Option_t *option="")
Double_t * _hbeta
Estimates of the deterministic trend parameters (size: _nP)
Definition TKriging.h:110
TString _outputVarName
Name of the output variable. This information can be used to remind the user of the original names of...
Definition TKriging.h:85
static void StatVect(std::string title, int n, Double_t *v)
Int_t getNTrendParams()
Return the number of parameters of the deterministic trend.
Definition TKriging.h:223
Bool_t _blog
Boolean for edit the log.
Definition TKriging.h:83
Int_t getNInputs()
Return the number of input parameters.
Definition TKriging.h:199
Double_t * getTrendParams()
Return the array of the trend parameters.
Definition TKriging.h:253
Int_t getNCorrLengths()
Return the number of correlation lengths.
Definition TKriging.h:217
Int_t eval(Double_t *x0, Double_t *y0, int=0)
Evaluate the gaussian process at the new location x0.
static void StatVect(std::string title, int n, Double_t *v, double &vmin, double &vmax, double &mean, double &var, double &norm)
Double_t getMeasurementErrorVariance()
Return the variance of the measurement error.
Definition TKriging.h:288
Double_t * _normCorrLengths
The correlation lengths in the normalised input space. This corresponds to the correlation lengths st...
Definition TKriging.h:103
void setLooErrors(Double_t *err, Double_t *var, Double_t mse)
Set the information about the leave one out error of the model.
Double_t * getM()
Return the M matrix as an array (cf. internal variable description)
Definition TKriging.h:300
TString _inputVarNames
Name of the inputs variables separated by ":". This information can be used to remind the user of the...
Definition TKriging.h:86
Double_t * _vbeta
Covariance matrix of the trend parameters (size: _nP*_nP)
Definition TKriging.h:111
void exportFunction(const char *lang, const char *file, const char *name, Option_t *option="")
virtual void exportModelPMML(const char *file, const char *name, Option_t *option) const
virtual void exportModelCplusplus(std::ofstream *sourcefile, const char *name="")
Export the model as a C function (not implemented)
Double_t * _xObs
The matrix of inputs (dimension: _nX * _nS)
Definition TKriging.h:93
Bool_t hasTrend()
Return kTRUE if a deterministic trend has been defined.
Definition TKriging.h:187
TCorrelationFunction * _correlationFunction
The correlation function.
Definition TKriging.h:101
Int_t _nP
Number of coefficients of the deterministic trend.
Definition TKriging.h:109
ROOT.
Definition TAnisp.h:164