English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Modeler: TGPBuilder.h Source File
Uranie / Modeler  v4.11.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 #include "TOptimShare.h"
69 
70 #include <vector>
71 
73 extern "C"
74 {
75 int gpMod0(double *Kin, double *yin, int n, double *gamma, double *iL);
76 int gpModh(double *Hin, double *Kin, double *yin, int n, int p, double *hbeta, double *vbeta, double *gamma, double *iL, double *iR, double *M);
77 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);
78 
79 int gpLoo0Error(double *Kin, double *yin, int n, double *errorLoo, double *varLoo, double *eqmLoo);
80 int gpLoohError(double *Hin, double *Kin, double *yin, int n, int p, double *errorLoo, double *varLoo, double *eqmLoo);
81 int gpLoohBayesError(double *BetaPrior, double *Qprior, double *Hin, double *Kin, double *yin, int n, int p, double *errorLoo, double *varLoo, double *eqmLoo);
82 }
83 ;
84 
85 namespace URANIE
86 {
87  namespace Modeler
88  {
89  class TGPBuilder
90  {
91 
92  // Attributes
93  private:
94  Bool_t _blog;
95 
96  void LSort(Double_t *costValues, Int_t *minIndexes, Int_t nbthreads, Int_t n);
97 
98  public:
99  URANIE::DataServer::TDataServer *_tds;
100 
103  Int_t _nX;
104  Int_t _nY;
105 
110 
111  Int_t _nS;
112  Int_t _nP;
113  Double_t reg;
114  Double_t tolerance;
115  Double_t sigma2;
116  Double_t alpha;
117  std::vector<double> Vmes;
118  std::vector<double> H;
119  std::vector<double> C;
120  std::vector<double> K;
121  std::vector<double> xObs;
122  std::vector<double> yObs;
123  std::vector<double> xNormParams;
124 
125  Double_t _modlengthMax;
126  Double_t _modlengthMin;
128 
129  TString trendString;
130 // std::vector<TFormula> trendCoefList; ///< List of the formulas of the trend coefficients (_nP elements)
131  Bool_t bHas_trend;
133 
134  std::vector<double> BetaPrior;
135  std::vector<double> QPrior;
136  Bool_t bUse_prior;
137  Int_t _nSeed;
138 
139  Int_t _retOptim;
140 
142 
143  // Operations
144  public:
145  //---------------------------------------------
149  TGPBuilder();
151 
153 
162  TGPBuilder(URANIE::DataServer::TDataServer *tdsObs,
163  const char *varexpinput,
164  const char *varexpoutput,
165  TCorrelationFunction* corrFunc,
166  const char *trend = "",
167  Option_t *option = "");
168 
170 
181  TGPBuilder(URANIE::DataServer::TDataServer *tdsObs,
182  const char *varexpinput,
183  const char *varexpoutput,
184  const char *corrFunc,
185  const char *trend = "",
186  Option_t *option = "");
187 
189 
201  TGPBuilder(URANIE::DataServer::TDataServer *tdsModelData);
202 
204 
206  void initVariables();
207 
209  void initMatrices();
210 
212  void createTrend();
213 
215  virtual ~TGPBuilder();
217 
218  //---------------------------------------------
222  Bool_t hasCorrFunction()
224  {
225  return bHas_corrFunction;
226  }
227 
230  {
231  return bVariance_in_factor;
232  }
233 
236  {
237  return bHas_measurement_error;
238  }
239 
241  Bool_t hasTrend()
242  {
243  return bHas_trend;
244  }
245 
247  Bool_t usePrior()
248  {
249  return bUse_prior;
250  }
251 
254  {
255  return _retOptim;
256  }
257 
260  {
261  return correlationFunction;
262  }
263 
265  Int_t getNInputs()
266  {
267  return _nX;
268  }
269 
271  Int_t getNObs()
272  {
273  return _nS;
274  }
275 
278  {
279  return _nP;
280  }
281 
283 
287  std::vector<double> const& getObsInputs()
288  {
289  return xObs;
290  }
291 
293  std::vector<double> const& getObsOutput()
294  {
295  return yObs;
296  }
297 
299 
304  std::vector<double> const& getObsNormParams()
305  {
306  return xNormParams;
307  }
308 
310  std::vector<double> const& getCovarianceMatrix()
311  {
312  return K;
313  }
314 
316  {
317  correlationFunction = corrFunc;
319  }
320 
322  std::vector<double> const& getCorrelationMatrix()
323  {
324  return C;
325  }
326 
328  std::vector<double> const& getMeasurementErrorCovMatrix()
329  {
330  return Vmes;
331  }
332 
334  std::vector<double> const& getTrendMatrix()
335  {
336  return H;
337  }
338 
340  /* std::vector<TFormula> const& getTrendFormula()
341  {
342  return trendCoefList;
343  } */
344 
346  TString getTrendString()
347  {
348  return trendString;
349  }
350 
352  std::vector<double> const& getTrendPriorMeans()
353  {
354  return BetaPrior;
355  }
356 
358  std::vector<double> const& getTrendPriorCovMat()
359  {
360  return QPrior;
361  }
362 
364  Double_t getVariance()
365  {
366  return sigma2;
367  }
368 
370  Double_t getTolerance()
371  {
372  return tolerance;
373  }
374 
376  Double_t getAlpha()
377  {
378  return alpha;
379  }
380 
382  Double_t getRegularisation()
383  {
384  return reg;
385  }
386 
387 
389  Int_t getSeed()
390  {
391  return _nSeed;
392  }
393 
395 
396  //---------------------------------------------
400 
408  void setLengthRange(double themin, double themax);
409 
411 
418  void setLengthMaxScreeningSF(double sf);
419 
421 
427  void setCorrFunction(TCorrelationFunction* corrFunc);
428 
430 
437  void setCorrFunction(const char* corrFuncType);
438 
440 
443  void setVariance(Double_t var)
444  {
445  sigma2 = var;
446  }
447 
449 
452  void setTolerance(Double_t tol)
453  {
454  tolerance = tol;
455  }
456 
458 
466  void setRegularisation(Double_t regularisation)
467  {
468  reg = regularisation;
469  }
470 
472 
482  void setHasMeasurementError(Bool_t has_measurement_error = kTRUE)
483  {
484  bHas_measurement_error = has_measurement_error;
485  }
486 
488 
496  void setAlpha(Double_t a)
497  {
498  alpha = a;
499  bHas_measurement_error = kTRUE;
500  bVariance_in_factor = kTRUE;
501  }
502 
504 
511  void setMeasurementErrorCovMatrix(Double_t* covMes);
512 
514 
521  void setMeasurementErrorVariance(Double_t varMes);
522 
524 
542  void setTrendString(const char* trend)
543  {
544  trendString = TString(trend);
545  }
546 
548 
557  void setPriorData(Double_t* betaPrior, Double_t* qPrior);
558 
560 
565  void setUsePrior(Bool_t use_prior = kTRUE);
566 
568  void computeCovarianceMatrix(bool print = false);
569 
571 
576  void setSeed(Int_t nSeed = 0);
577 
579 
580  //---------------------------------------------
584 
586 
591  void updateObservations();
592 
594 
607  void exportGPData(const char* fileName);
608 
609 
611 
675  void findOptimalParameters(TString criterion = "ML",
676  Int_t screeningSize = 100,
677  TString optimAlgo = "BFGS",
678  Int_t nbMaxOptimSteps = 1000,
679  Bool_t reset = kTRUE,
680  Option_t * option = "");
681  void findOptimalParameters(TString criterion,
682  Int_t screeningSize,
683  URANIE::Reoptimizer::TOptimSolver *solv,
684  Int_t nbMaxOptimSteps = 1000,
685  Bool_t reset = kTRUE,
686  Option_t * option = "");
687 
689 
699  URANIE::Modeler::TKriging* buildGP(Bool_t computeLooErrors = kTRUE);
700 
702  void cleanAttributeList(TList *inputlist);
704 
705  //---------------------------------------------
711  void setLog()
712  {
713  _blog = kTRUE;
714  }
715 
716  void unsetLog()
717  {
718  _blog = kFALSE;
719  }
720 
721  void changeLog()
722  {
723  _blog = !_blog;
724  }
725 
726  Bool_t getLog()
727  {
728  return _blog;
729  }
730 
731  virtual void printLog(Option_t *option = "");
733 
734  ClassDef(URANIE::Modeler::TGPBuilder, ID_MODELER)
735  //Classe de
736  };
737  } // Fin du namespace Modeler
738 } // Fin du namespace URANIE
739 #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:322
void setVariance(Double_t var)
Set the variance of the gaussian process.
Definition: TGPBuilder.h:443
TString _sInputAttributes
input attributes names separated by ":"
Definition: TGPBuilder.h:101
std::vector< double > const & getTrendPriorMeans()
Return the array of the a priori means of the trend parameters.
Definition: TGPBuilder.h:352
TString getTrendString()
Return the deterministic trend list of formulas.
Definition: TGPBuilder.h:346
Bool_t bUse_prior
kFALSE if no prior is used, kTRUE otherwise
Definition: TGPBuilder.h:136
URANIE::DataServer::TDataServer * _tds
Pointer to the data server containing the observations.
Definition: TGPBuilder.h:99
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:328
Int_t _retOptim
optimisation results
Definition: TGPBuilder.h:139
std::vector< double > yObs
Array representation of the outputs (size: _nS)
Definition: TGPBuilder.h:122
void computeCovarianceMatrix(bool print=false)
Computes the covariance matrix of the observed data.
Double_t sigma2
Variance of the gaussian process.
Definition: TGPBuilder.h:115
Bool_t hasVarianceInFactor()
Return kTRUE if the variance of the GP is determined analytically.
Definition: TGPBuilder.h:229
Int_t getNObs()
Return the number of observations used to build the GP model.
Definition: TGPBuilder.h:271
std::vector< double > QPrior
a priori variance matrix for the deterministic trend parameters (size: _nP x _nP) ...
Definition: TGPBuilder.h:135
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:116
void setTolerance(Double_t tol)
Set the tolerance used to stop the optimisation process.
Definition: TGPBuilder.h:452
Double_t getVariance()
Return the variance of the gaussian process.
Definition: TGPBuilder.h:364
Int_t _nSeed
number of seed for pseudo-random generator
Definition: TGPBuilder.h:137
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:114
Int_t _nX
Number of input attributes.
Definition: TGPBuilder.h:103
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:241
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:134
Bool_t _blog
Boolean for edit the log.
Definition: TGPBuilder.h:94
void setCorrFunc(TCorrelationFunction *corrFunc)
Definition: TGPBuilder.h:315
Double_t getTolerance()
Return the tolerance used to stop the optimisation process.
Definition: TGPBuilder.h:370
int getReturnOptim()
Return the return of the optimisation procedure.
Definition: TGPBuilder.h:253
TCorrelationFunction * correlationFunction
The correlation function.
Definition: TGPBuilder.h:106
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:132
std::vector< double > const & getTrendMatrix()
Return the array of the deterministic trend matrix.
Definition: TGPBuilder.h:334
Bool_t bIsbuildWithDuplicateKandCandCorrFunc
kTRUE if build with duplicate K, C and CorrFunc, otherwise kFALSE
Definition: TGPBuilder.h:141
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:117
Bool_t bVariance_in_factor
if kTRUE, GP variance is calculated analytically. Otherwise, it is a parameter of the cost function...
Definition: TGPBuilder.h:108
void setAlpha(Double_t a)
Set the value of the alpha parameter.
Definition: TGPBuilder.h:496
std::vector< double > const & getObsOutput()
Return the observations outputs.
Definition: TGPBuilder.h:293
TString _sOutputAttributes
output attribute name
Definition: TGPBuilder.h:102
void setLog()
Definition: TGPBuilder.h:711
TString trendString
Character string containing the trend parameters separated by ":" or a trend type ("Const" or "Linear...
Definition: TGPBuilder.h:129
std::vector< double > const & getObsInputs()
Return the observations inputs as a flat array.
Definition: TGPBuilder.h:287
Double_t _modlengthMin
Minimum value of the correlation length if default behaviour is requested to be changed by the user...
Definition: TGPBuilder.h:126
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:109
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:376
void cleanAttributeList(TList *inputlist)
Clean list of attribute.
std::vector< double > C
The correlation matrix (dimension: _nS * _nS)
Definition: TGPBuilder.h:119
Int_t _nS
Number of observations.
Definition: TGPBuilder.h:111
Int_t getNTrendParams()
Return the number of parameters of the deterministic trend.
Definition: TGPBuilder.h:277
void setRegularisation(Double_t regularisation)
Set the regularisation coefficient.
Definition: TGPBuilder.h:466
Bool_t getLog()
Definition: TGPBuilder.h:726
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:113
Bool_t hasMeasurementError()
Return kTRUE if a measurement error should be considered to build the GP.
Definition: TGPBuilder.h:235
Double_t getRegularisation()
Return the regularisation coefficient.
Definition: TGPBuilder.h:382
std::vector< double > const & getCovarianceMatrix()
Return the array of the covariance matrix of the gaussian process.
Definition: TGPBuilder.h:310
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:118
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:125
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:121
std::vector< double > const & getTrendPriorCovMat()
Return the array of the a priori covariance matrix of the trend parameters.
Definition: TGPBuilder.h:358
void unsetLog()
Definition: TGPBuilder.h:716
TGPBuilder()
Default Constructor.
Int_t _nP
Number of coefficients of the deterministic trend.
Definition: TGPBuilder.h:112
Description of the class TGPBuilder.
Definition: TGPBuilder.h:89
Bool_t bHas_corrFunction
kFALSE if correlationFunction == NULL, kTRUE otherwise
Definition: TGPBuilder.h:107
std::vector< double > const & getObsNormParams()
Return the normalisation parameters of the inputs as a flat array.
Definition: TGPBuilder.h:304
std::vector< double > K
The covariance matrix (dimension: _nS * _nS)
Definition: TGPBuilder.h:120
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:482
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:104
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:247
Int_t getSeed()
Return the seed number.
Definition: TGPBuilder.h:389
Int_t getNInputs()
Return the number of input variables.
Definition: TGPBuilder.h:265
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:542
TCorrelationFunction * getCorrFunction()
Return a pointer to the correlation function of the gaussian process.
Definition: TGPBuilder.h:259
Double_t _lengthMaxScreeningSF
Scale factor to change the Maximum value of the correlation length during screening procedure: new le...
Definition: TGPBuilder.h:127
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:131
Bool_t hasCorrFunction()
Return kTRUE if a correlation function is defined, kFALSE otherwise.
Definition: TGPBuilder.h:223
std::vector< double > xNormParams
The matrix of normalisation parameters. Contains min and max value of each input variable: [x1Min...
Definition: TGPBuilder.h:123
void changeLog()
Definition: TGPBuilder.h:721
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()