English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Modeler: TGPBuilder.h Source File
Uranie / Modeler  v4.10.0
/* @license-end */
TGPBuilder.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 
56 #ifndef TGPBUILDER_H
57 #define TGPBUILDER_H
58 
60 #include "TObjArray.h"
61 
63 #include "TCorrelationFunction.h"
64 #include "TShareEval.h"
65 #include "TDataServer.h"
66 #include "TAttribute.h"
67 #include "TKriging.h"
68 
69 #include <vector>
70 
72 extern "C"
73 {
74 int gpMod0(double *Kin, double *yin, int n, double *gamma, double *iL);
75 int gpModh(double *Hin, double *Kin, double *yin, int n, int p, double *hbeta, double *vbeta, double *gamma, double *iL, double *iR, double *M);
76 int gpModhBayes(double *BetaPrior, double *QPrior, double *Hin, double *Kin, double *yin, int n, int p, double *BetaPost, double *QPost, double *gamma, double *iL, double *iR, double *M);
77 
78 int gpLoo0Error(double *Kin, double *yin, int n, double *errorLoo, double *varLoo, double *eqmLoo);
79 int gpLoohError(double *Hin, double *Kin, double *yin, int n, int p, double *errorLoo, double *varLoo, double *eqmLoo);
80 int gpLoohBayesError(double *BetaPrior, double *Qprior, double *Hin, double *Kin, double *yin, int n, int p, double *errorLoo, double *varLoo, double *eqmLoo);
81 }
82 ;
83 
84 namespace URANIE
85 {
86  namespace Modeler
87  {
88  class TGPBuilder
89  {
90 
91  // Attributes
92  private:
93  Bool_t _blog;
94 
95  void LSort(Double_t *costValues, Int_t *minIndexes, Int_t nbthreads, Int_t n);
96 
97  public:
98  URANIE::DataServer::TDataServer *_tds;
99 
102  Int_t _nX;
103  Int_t _nY;
104 
109 
110  Int_t _nS;
111  Int_t _nP;
112  Double_t reg;
113  Double_t tolerance;
114  Double_t sigma2;
115  Double_t alpha;
116  std::vector<double> Vmes;
117  std::vector<double> H;
118  std::vector<double> C;
119  std::vector<double> K;
120  std::vector<double> xObs;
121  std::vector<double> yObs;
122  std::vector<double> xNormParams;
123 
124  Double_t _modlengthMax;
125  Double_t _modlengthMin;
127 
128  TString trendString;
129  std::vector<TFormula> trendCoefList;
130  Bool_t bHas_trend;
132 
133  std::vector<double> BetaPrior;
134  std::vector<double> QPrior;
135  Bool_t bUse_prior;
136  Int_t _nSeed;
137 
138  Int_t _retOptim;
139 
141 
142  // Operations
143  public:
144  //---------------------------------------------
148  TGPBuilder();
150 
152 
161  TGPBuilder(URANIE::DataServer::TDataServer *tdsObs,
162  const char *varexpinput,
163  const char *varexpoutput,
164  TCorrelationFunction* corrFunc,
165  const char *trend = "",
166  Option_t *option = "");
167 
169 
180  TGPBuilder(URANIE::DataServer::TDataServer *tdsObs,
181  const char *varexpinput,
182  const char *varexpoutput,
183  const char *corrFunc,
184  const char *trend = "",
185  Option_t *option = "");
186 
188 
200  TGPBuilder(URANIE::DataServer::TDataServer *tdsModelData);
201 
203 
205  void initVariables();
206 
208  void initMatrices();
209 
211  void createTrend();
212 
214  virtual ~TGPBuilder();
216 
217  //---------------------------------------------
221  Bool_t hasCorrFunction()
223  {
224  return bHas_corrFunction;
225  }
226 
229  {
230  return bVariance_in_factor;
231  }
232 
235  {
236  return bHas_measurement_error;
237  }
238 
240  Bool_t hasTrend()
241  {
242  return bHas_trend;
243  }
244 
246  Bool_t usePrior()
247  {
248  return bUse_prior;
249  }
250 
253  {
254  return _retOptim;
255  }
256 
259  {
260  return correlationFunction;
261  }
262 
264  Int_t getNInputs()
265  {
266  return _nX;
267  }
268 
270  Int_t getNObs()
271  {
272  return _nS;
273  }
274 
277  {
278  return _nP;
279  }
280 
282 
286  std::vector<double> const& getObsInputs()
287  {
288  return xObs;
289  }
290 
292  std::vector<double> const& getObsOutput()
293  {
294  return yObs;
295  }
296 
298 
303  std::vector<double> const& getObsNormParams()
304  {
305  return xNormParams;
306  }
307 
309  std::vector<double> const& getCovarianceMatrix()
310  {
311  return K;
312  }
313 
314  void setKandCandCorrFunc(const std::vector<double>& inK, const std::vector<double>& inC,
315  TCorrelationFunction* corrFunc)
316  {
317  K = inK;
318  C = inC;
319  correlationFunction = corrFunc;
321  }
322 
324  std::vector<double> const& getCorrelationMatrix()
325  {
326  return C;
327  }
328 
330  std::vector<double> const& getMeasurementErrorCovMatrix()
331  {
332  return Vmes;
333  }
334 
336  std::vector<double> const& getTrendMatrix()
337  {
338  return H;
339  }
340 
342  std::vector<TFormula> const& getTrendFormula()
343  {
344  return trendCoefList;
345  }
346 
348  TString getTrendString()
349  {
350  return trendString;
351  }
352 
354  std::vector<double> const& getTrendPriorMeans()
355  {
356  return BetaPrior;
357  }
358 
360  std::vector<double> const& getTrendPriorCovMat()
361  {
362  return QPrior;
363  }
364 
366  Double_t getVariance()
367  {
368  return sigma2;
369  }
370 
372  Double_t getTolerance()
373  {
374  return tolerance;
375  }
376 
378  Double_t getAlpha()
379  {
380  return alpha;
381  }
382 
384  Double_t getRegularisation()
385  {
386  return reg;
387  }
388 
389 
391  Int_t getSeed()
392  {
393  return _nSeed;
394  }
395 
397 
398  //---------------------------------------------
402 
410  void setLengthRange(double themin, double themax);
411 
413 
420  void setLengthMaxScreeningSF(double sf);
421 
423 
429  void setCorrFunction(TCorrelationFunction* corrFunc);
430 
432 
439  void setCorrFunction(const char* corrFuncType);
440 
442 
445  void setVariance(Double_t var)
446  {
447  sigma2 = var;
448  }
449 
451 
454  void setTolerance(Double_t tol)
455  {
456  tolerance = tol;
457  }
458 
460 
468  void setRegularisation(Double_t regularisation)
469  {
470  reg = regularisation;
471  }
472 
474 
484  void setHasMeasurementError(Bool_t has_measurement_error = kTRUE)
485  {
486  bHas_measurement_error = has_measurement_error;
487  }
488 
490 
498  void setAlpha(Double_t a)
499  {
500  alpha = a;
501  bHas_measurement_error = kTRUE;
502  bVariance_in_factor = kTRUE;
503  }
504 
506 
513  void setMeasurementErrorCovMatrix(Double_t* covMes);
514 
516 
523  void setMeasurementErrorVariance(Double_t varMes);
524 
526 
544  void setTrendString(const char* trend)
545  {
546  trendString = TString(trend);
547  }
548 
550 
559  void setPriorData(Double_t* betaPrior, Double_t* qPrior);
560 
562 
567  void setUsePrior(Bool_t use_prior = kTRUE);
568 
570  void computeCovarianceMatrix(bool print = false);
571 
573 
578  void setSeed(Int_t nSeed = 0);
579 
581 
582  //---------------------------------------------
586 
588 
593  void updateObservations();
594 
596 
609  void exportGPData(const char* fileName);
610 
611 
613 
677  void findOptimalParameters(TString criterion = "ML",
678  Int_t screeningSize = 100,
679  TString optimAlgo = "BFGS",
680  Int_t nbMaxOptimSteps = 1000,
681  Bool_t reset = kTRUE,
682  Option_t * option = "");
683 
685 
695  URANIE::Modeler::TKriging* buildGP(Bool_t computeLooErrors = kTRUE);
696 
698  void cleanAttributeList(TList *inputlist);
700 
701  //---------------------------------------------
707  void setLog()
708  {
709  _blog = kTRUE;
710  }
711 
712  void unsetLog()
713  {
714  _blog = kFALSE;
715  }
716 
717  void changeLog()
718  {
719  _blog = !_blog;
720  }
721 
722  Bool_t getLog()
723  {
724  return _blog;
725  }
726 
727  virtual void printLog(Option_t *option = "");
729 
730  ClassDef(URANIE::Modeler::TGPBuilder, ID_MODELER)
731  //Classe de
732  };
733  } // Fin du namespace Modeler
734 } // Fin du namespace URANIE
735 #endif
ROOT.
Definition: TAnisp.h:163
Interface de la classe URANIE::Modeler::TKriging.
void setPriorData(Double_t *betaPrior, Double_t *qPrior)
Set the a priori mean and variance of the trend coefficients.
void initMatrices()
Initialisation of the matrices.
std::vector< double > const & getCorrelationMatrix()
Return the array of the correlation matrix of the gaussian process.
Definition: TGPBuilder.h:324
void setVariance(Double_t var)
Set the variance of the gaussian process.
Definition: TGPBuilder.h:445
TString _sInputAttributes
input attributes names separated by ":"
Definition: TGPBuilder.h:100
std::vector< double > const & getTrendPriorMeans()
Return the array of the a priori means of the trend parameters.
Definition: TGPBuilder.h:354
TString getTrendString()
Return the deterministic trend character string.
Definition: TGPBuilder.h:348
Bool_t bUse_prior
kFALSE if no prior is used, kTRUE otherwise
Definition: TGPBuilder.h:135
URANIE::DataServer::TDataServer * _tds
Pointer to the data server containing the observations.
Definition: TGPBuilder.h:98
void createTrend()
Create the matrices for the deterministic trend.
std::vector< double > const & getMeasurementErrorCovMatrix()
Return the array of the covariance matrix of the measurement errors.
Definition: TGPBuilder.h:330
Int_t _retOptim
optimisation results
Definition: TGPBuilder.h:138
std::vector< double > yObs
Array representation of the outputs (size: _nS)
Definition: TGPBuilder.h:121
void computeCovarianceMatrix(bool print=false)
Computes the covariance matrix of the observed data.
Double_t sigma2
Variance of the gaussian process.
Definition: TGPBuilder.h:114
Bool_t hasVarianceInFactor()
Return kTRUE if the variance of the GP is determined analytically.
Definition: TGPBuilder.h:228
Int_t getNObs()
Return the number of observations used to build the GP model.
Definition: TGPBuilder.h:270
std::vector< double > QPrior
a priori variance matrix for the deterministic trend parameters (size: _nP x _nP) ...
Definition: TGPBuilder.h:134
void exportGPData(const char *fileName)
Export the minimum information necessary to build a new TKriging object.
void findOptimalParameters(TString criterion="ML", Int_t screeningSize=100, TString optimAlgo="BFGS", Int_t nbMaxOptimSteps=1000, Bool_t reset=kTRUE, Option_t *option="")
Search for the optimal parameters of the Gaussian Process.
Double_t alpha
Quotient of the variance of the measurment error and the variance of the gaussian process...
Definition: TGPBuilder.h:115
void setTolerance(Double_t tol)
Set the tolerance used to stop the optimisation process.
Definition: TGPBuilder.h:454
Double_t getVariance()
Return the variance of the gaussian process.
Definition: TGPBuilder.h:366
std::vector< TFormula > trendCoefList
List of the formulas of the trend coefficients (_nP elements)
Definition: TGPBuilder.h:129
Int_t _nSeed
number of seed for pseudo-random generator
Definition: TGPBuilder.h:136
int gpMod0(double *Kin, double *yin, int n, double *gamma, double *iL)
Double_t tolerance
Tolerance value used to stop the optimisation process.
Definition: TGPBuilder.h:113
Int_t _nX
Number of input attributes.
Definition: TGPBuilder.h:102
int gpLoohError(double *Hin, double *Kin, double *yin, int n, int p, double *errorLoo, double *varLoo, double *eqmLoo)
void setLengthRange(double themin, double themax)
Set the Minimun and maximum boundaries for the correlation length.
Bool_t hasTrend()
Return kTRUE if a deterministic trend has been defined.
Definition: TGPBuilder.h:240
std::vector< TFormula > const & getTrendFormula()
Return the deterministic trend list of formulas.
Definition: TGPBuilder.h:342
void initVariables()
Initialisation of the GPBuilder variables.
std::vector< double > BetaPrior
a priori mean vector for the deterministic trend parameters (size: _nP)
Definition: TGPBuilder.h:133
Bool_t _blog
Boolean for edit the log.
Definition: TGPBuilder.h:93
Double_t getTolerance()
Return the tolerance used to stop the optimisation process.
Definition: TGPBuilder.h:372
int getReturnOptim()
Return the return of the optimisation procedure.
Definition: TGPBuilder.h:252
TCorrelationFunction * correlationFunction
The correlation function.
Definition: TGPBuilder.h:105
Bool_t bUse_normalisation
if kFALSE, the matrix H will be computed with unormalised input data. This happens when the user is p...
Definition: TGPBuilder.h:131
std::vector< double > const & getTrendMatrix()
Return the array of the deterministic trend matrix.
Definition: TGPBuilder.h:336
Bool_t bIsbuildWithDuplicateKandCandCorrFunc
kTRUE if build with duplicate K, C and CorrFunc, otherwise kFALSE
Definition: TGPBuilder.h:140
virtual void printLog(Option_t *option="")
int gpLoo0Error(double *Kin, double *yin, int n, double *errorLoo, double *varLoo, double *eqmLoo)
void setLengthMaxScreeningSF(double sf)
Set the scale factor to reduce the maximum value of the correlation length during screening...
Description of the class TCorrelationFunction.
Definition: TCorrelationFunction.h:52
std::vector< double > Vmes
The covariance matrix of the measurement errors (dimension: _nS * _nS)
Definition: TGPBuilder.h:116
Bool_t bVariance_in_factor
if kTRUE, GP variance is calculated analytically. Otherwise, it is a parameter of the cost function...
Definition: TGPBuilder.h:107
void setAlpha(Double_t a)
Set the value of the alpha parameter.
Definition: TGPBuilder.h:498
std::vector< double > const & getObsOutput()
Return the observations outputs.
Definition: TGPBuilder.h:292
TString _sOutputAttributes
output attribute name
Definition: TGPBuilder.h:101
void setLog()
Definition: TGPBuilder.h:707
TString trendString
Character string containing the trend parameters separated by ":" or a trend type ("Const" or "Linear...
Definition: TGPBuilder.h:128
std::vector< double > const & getObsInputs()
Return the observations inputs as a flat array.
Definition: TGPBuilder.h:286
Double_t _modlengthMin
Minimum value of the correlation length if default behaviour is requested to be changed by the user...
Definition: TGPBuilder.h:125
Bool_t bHas_measurement_error
if kTRUE, the construction of the GP will take into account the existence of a measurement error vari...
Definition: TGPBuilder.h:108
virtual ~TGPBuilder()
Default destructor.
Double_t getAlpha()
Return the alpha parameter, quotient of the variance of the measurement error and the variance of the...
Definition: TGPBuilder.h:378
void cleanAttributeList(TList *inputlist)
Clean list of attribute.
std::vector< double > C
The correlation matrix (dimension: _nS * _nS)
Definition: TGPBuilder.h:118
Int_t _nS
Number of observations.
Definition: TGPBuilder.h:110
Int_t getNTrendParams()
Return the number of parameters of the deterministic trend.
Definition: TGPBuilder.h:276
void setRegularisation(Double_t regularisation)
Set the regularisation coefficient.
Definition: TGPBuilder.h:468
Bool_t getLog()
Definition: TGPBuilder.h:722
void LSort(Double_t *costValues, Int_t *minIndexes, Int_t nbthreads, Int_t n)
void setMeasurementErrorVariance(Double_t varMes)
Set the variance of the measurement errors.
Double_t reg
regularisation parameter, used to improve numerical stability of matrix operations ...
Definition: TGPBuilder.h:112
Bool_t hasMeasurementError()
Return kTRUE if a measurement error should be considered to build the GP.
Definition: TGPBuilder.h:234
Double_t getRegularisation()
Return the regularisation coefficient.
Definition: TGPBuilder.h:384
std::vector< double > const & getCovarianceMatrix()
Return the array of the covariance matrix of the gaussian process.
Definition: TGPBuilder.h:309
int gpModhBayes(double *BetaPrior, double *QPrior, double *Hin, double *Kin, double *yin, int n, int p, double *BetaPost, double *QPost, double *gamma, double *iL, double *iR, double *M)
std::vector< double > H
Deterministic trend matrix (dimension: _nS * _nP)
Definition: TGPBuilder.h:117
void setCorrFunction(TCorrelationFunction *corrFunc)
Set the correlation function of the gaussian process.
void updateObservations()
Update the GP builder matrices.
Double_t _modlengthMax
Maximum value of the correlation length if default behaviour is requested to be changed by the user...
Definition: TGPBuilder.h:124
int gpLoohBayesError(double *BetaPrior, double *Qprior, double *Hin, double *Kin, double *yin, int n, int p, double *errorLoo, double *varLoo, double *eqmLoo)
Interface de la classe URANIE::Modeler::TCorrelationFunction.
std::vector< double > xObs
The matrix of normalised inputs (dimension: _nX * _nS)
Definition: TGPBuilder.h:120
std::vector< double > const & getTrendPriorCovMat()
Return the array of the a priori covariance matrix of the trend parameters.
Definition: TGPBuilder.h:360
void unsetLog()
Definition: TGPBuilder.h:712
TGPBuilder()
Default Constructor.
void setKandCandCorrFunc(const std::vector< double > &inK, const std::vector< double > &inC, TCorrelationFunction *corrFunc)
Definition: TGPBuilder.h:314
Int_t _nP
Number of coefficients of the deterministic trend.
Definition: TGPBuilder.h:111
Description of the class TGPBuilder.
Definition: TGPBuilder.h:88
Bool_t bHas_corrFunction
kFALSE if correlationFunction == NULL, kTRUE otherwise
Definition: TGPBuilder.h:106
std::vector< double > const & getObsNormParams()
Return the normalisation parameters of the inputs as a flat array.
Definition: TGPBuilder.h:303
std::vector< double > K
The covariance matrix (dimension: _nS * _nS)
Definition: TGPBuilder.h:119
void setSeed(Int_t nSeed=0)
set a seed number
void setHasMeasurementError(Bool_t has_measurement_error=kTRUE)
Indicate if a measurement error should be taken into account or not.
Definition: TGPBuilder.h:484
void setUsePrior(Bool_t use_prior=kTRUE)
Set wether or not to use a priori knowledge on the deterministic trend parameters.
Int_t _nY
Number of output attributes (should be equal to 1)
Definition: TGPBuilder.h:103
Description of the class TKriging.
Definition: TKriging.h:81
Bool_t usePrior()
Return kTRUE if prior knowledge on the deterministic trend is used (bayesian approach).
Definition: TGPBuilder.h:246
Int_t getSeed()
Return the seed number.
Definition: TGPBuilder.h:391
Int_t getNInputs()
Return the number of input variables.
Definition: TGPBuilder.h:264
void setMeasurementErrorCovMatrix(Double_t *covMes)
Set the covariance matrix of the measurement errors.
void setTrendString(const char *trend)
Set the deterministic trend character string.
Definition: TGPBuilder.h:544
TCorrelationFunction * getCorrFunction()
Return a pointer to the correlation function of the gaussian process.
Definition: TGPBuilder.h:258
Double_t _lengthMaxScreeningSF
Scale factor to change the Maximum value of the correlation length during screening procedure: new le...
Definition: TGPBuilder.h:126
URANIE::Modeler::TKriging * buildGP(Bool_t computeLooErrors=kTRUE)
Build the gaussian process model.
Bool_t bHas_trend
kFALSE if no deterministic trend is defined, kTRUE otherwise.
Definition: TGPBuilder.h:130
Bool_t hasCorrFunction()
Return kTRUE if a correlation function is defined, kFALSE otherwise.
Definition: TGPBuilder.h:222
std::vector< double > xNormParams
The matrix of normalisation parameters. Contains min and max value of each input variable: [x1Min...
Definition: TGPBuilder.h:122
void changeLog()
Definition: TGPBuilder.h:717
int gpModh(double *Hin, double *Kin, double *yin, int n, int p, double *hbeta, double *vbeta, double *gamma, double *iL, double *iR, double *M)
TGPBuilder * buildWithDuplicateKandCandCorrFunc()