English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Modeler: URANIE::Modeler::TGPBuilder Class Reference
Uranie / Modeler v4.9.0
/* @license-end */
URANIE::Modeler::TGPBuilder Class Reference

Description of the class TGPBuilder. More...

#include <TGPBuilder.h>

Collaboration diagram for URANIE::Modeler::TGPBuilder:

Public Member Functions

Constructor and Destructor
 TGPBuilder ()
 Default Constructor.
 
 TGPBuilder (URANIE::DataServer::TDataServer *tdsObs, const char *varexpinput, const char *varexpoutput, TCorrelationFunction *corrFunc, const char *trend="", Option_t *option="")
 Construction from a data server, a list of input and output variables and a correlation function object.
 
 TGPBuilder (URANIE::DataServer::TDataServer *tdsObs, const char *varexpinput, const char *varexpoutput, const char *corrFunc, const char *trend="", Option_t *option="")
 Construction from a data server, a list of input and output variables and a correlation function name.
 
 TGPBuilder (URANIE::DataServer::TDataServer *tdsModelData)
 Construction from a data server containing model parameters.
 
TGPBuilderbuildWithDuplicateKandCandCorrFunc ()
 
void initVariables ()
 Initialisation of the GPBuilder variables.
 
void initMatrices ()
 Initialisation of the matrices.
 
void createTrend ()
 Create the matrices for the deterministic trend.
 
virtual ~TGPBuilder ()
 Default destructor.
 
Get methods
Bool_t hasCorrFunction ()
 Return kTRUE if a correlation function is defined, kFALSE otherwise.
 
Bool_t hasVarianceInFactor ()
 Return kTRUE if the variance of the GP is determined analytically.
 
Bool_t hasMeasurementError ()
 Return kTRUE if a measurement error should be considered to build the GP.
 
Bool_t hasTrend ()
 Return kTRUE if a deterministic trend has been defined.
 
Bool_t usePrior ()
 Return kTRUE if prior knowledge on the deterministic trend is used (bayesian approach).
 
int getReturnOptim ()
 Return the return of the optimisation procedure.
 
TCorrelationFunctiongetCorrFunction ()
 Return a pointer to the correlation function of the gaussian process.
 
Int_t getNInputs ()
 Return the number of input variables.
 
Int_t getNObs ()
 Return the number of observations used to build the GP model.
 
Int_t getNTrendParams ()
 Return the number of parameters of the deterministic trend.
 
Double_t * getObsInputs ()
 Return the observations inputs as a flat array.
 
Double_t * getObsOutput ()
 Return the observations outputs.
 
Double_t * getObsNormParams ()
 Return the normalisation parameters of the inputs as a flat array.
 
Double_t * getCovarianceMatrix ()
 Return the array of the covariance matrix of the gaussian process.
 
void setKandCandCorrFunc (Double_t *inK, Double_t *inC, TCorrelationFunction *corrFunc)
 
Double_t * getCorrelationMatrix ()
 Return the array of the correlation matrix of the gaussian process.
 
Double_t * getMeasurementErrorCovMatrix ()
 Return the array of the covariance matrix of the measurement errors.
 
Double_t * getTrendMatrix ()
 Return the array of the deterministic trend matrix.
 
TObjArray * getTrendFormula ()
 Return the deterministic trend list of formulas.
 
TString getTrendString ()
 Return the deterministic trend character string.
 
Double_t * getTrendPriorMeans ()
 Return the array of the a priori means of the trend parameters.
 
Double_t * getTrendPriorCovMat ()
 Return the array of the a priori covariance matrix of the trend parameters.
 
Double_t getVariance ()
 Return the variance of the gaussian process.
 
Double_t getTolerance ()
 Return the tolerance used to stop the optimisation process

 
Double_t getAlpha ()
 Return the alpha parameter, quotient of the variance of the measurement error and the variance of the gaussian process.
 
Double_t getRegularisation ()
 Return the regularisation coefficient.
 
Int_t getSeed ()
 Return the seed number.
 
Set methods
void setLengthRange (double themin, double themax)
 Set the Minimun and maximum boundaries for the correlation length.
 
void setLengthMaxScreeningSF (double sf)
 Set the scale factor to reduce the maximum value of the correlation length during screening.
 
void setCorrFunction (TCorrelationFunction *corrFunc)
 Set the correlation function of the gaussian process.
 
void setCorrFunction (const char *corrFuncType)
 Set the correlation function of the gaussian process from a character string.
 
void setVariance (Double_t var)
 Set the variance of the gaussian process.
 
void setTolerance (Double_t tol)
 Set the tolerance used to stop the optimisation process.
 
void setRegularisation (Double_t regularisation)
 Set the regularisation coefficient.
 
void setHasMeasurementError (Bool_t has_measurement_error=kTRUE)
 Indicate if a measurement error should be taken into account or not.
 
void setAlpha (Double_t a)
 Set the value of the alpha parameter.
 
void setMeasurementErrorCovMatrix (Double_t *covMes)
 Set the covariance matrix of the measurement errors.
 
void setMeasurementErrorVariance (Double_t varMes)
 Set the variance of the measurement errors.
 
void setTrendString (const char *trend)
 Set the deterministic trend character string.
 
void setPriorData (Double_t *betaPrior, Double_t *qPrior)
 Set the a priori mean and variance of the trend coefficients.
 
void setUsePrior (Bool_t use_prior=kTRUE)
 Set wether or not to use a priori knowledge on the deterministic trend parameters.
 
void computeCovarianceMatrix (bool print=false)
 Computes the covariance matrix of the observed data.
 
void setSeed (Int_t nSeed=0)
 set a seed number
 
Construction of the Gaussian Process
void updateObservations ()
 Update the GP builder matrices.
 
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.
 
URANIE::Modeler::TKrigingbuildGP (Bool_t computeLooErrors=kTRUE)
 Build the gaussian process model.
 
void cleanAttributeList (TList *inputlist)
 Clean list of attribute.
 
Printing Log
void setLog ()
 
void unsetLog ()
 
void changeLog ()
 
Bool_t getLog ()
 
virtual void printLog (Option_t *option="")
 

Public Attributes

URANIE::DataServer::TDataServer * _tds
 Pointer to the data server containing the observations.
 
TString _sInputAttributes
 input attributes names separated by ":"
 
TString _sOutputAttributes
 output attribute name
 
Int_t _nX
 Number of input attributes.
 
Int_t _nY
 Number of output attributes (should be equal to 1)
 
TCorrelationFunctioncorrelationFunction
 The correlation function.
 
Bool_t bHas_corrFunction
 kFALSE if correlationFunction == NULL, kTRUE otherwise
 
Bool_t bVariance_in_factor
 if kTRUE, GP variance is calculated analytically. Otherwise, it is a parameter of the cost function.
 
Bool_t bHas_measurement_error
 if kTRUE, the construction of the GP will take into account the existence of a measurement error variance.
 
Int_t _nS
 Number of observations.
 
Int_t _nP
 Number of coefficients of the deterministic trend.
 
Double_t reg
 regularisation parameter, used to improve numerical stability of matrix operations
 
Double_t tolerance
 Tolerance value used to stop the optimisation process.
 
Double_t sigma2
 Variance of the gaussian process.
 
Double_t alpha
 Quotient of the variance of the measurment error and the variance of the gaussian process.
 
Double_t * Vmes
 The covariance matrix of the measurement errors (dimension: _nS * _nS)
 
Double_t * H
 Deterministic trend matrix (dimension: _nS * _nP)
 
Double_t * C
 The correlation matrix (dimension: _nS * _nS)
 
Double_t * K
 The covariance matrix (dimension: _nS * _nS)
 
Double_t * xObs
 The matrix of normalised inputs (dimension: _nX * _nS)

 
Double_t * yObs
 Array representation of the outputs (size: _nS)
 
Double_t * xNormParams
 The matrix of normalisation parameters. Contains min and max value of each input variable: [x1Min, x1Max, x2Min, x2Max,...] (size: _nX * 2)

 
Double_t _modlengthMax
 Maximum value of the correlation length if default behaviour is requested to be changed by the user.
 
Double_t _modlengthMin
 Minimum value of the correlation length if default behaviour is requested to be changed by the user.
 
Double_t _lengthMaxScreeningSF
 Scale factor to change the Maximum value of the correlation length during screening procedure: new length is SF * lengthMax (0 < SF <=1)
 
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)
 
Bool_t bHas_trend
 kFALSE if no deterministic trend is defined, kTRUE otherwise.
 
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.
 
Double_t * BetaPrior
 a priori mean vector for the deterministic trend parameters (size: _nP)
 
Double_t * QPrior
 a priori variance matrix for the deterministic trend parameters (size: _nP x _nP)
 
Bool_t bUse_prior
 kFALSE if no prior is used, kTRUE otherwise
 
Int_t _nSeed
 number of seed for pseudo-random generator
 
Int_t _retOptim
 optimisation results
 
Bool_t bIsbuildWithDuplicateKandCandCorrFunc
 kTRUE if build with duplicate K, C and CorrFunc, otherwise kFALSE
 

Private Member Functions

void LSort (Double_t *costValues, Int_t *minIndexes, Int_t nbthreads, Int_t n)
 

Private Attributes

Bool_t _blog
 Boolean for edit the log.
 

Detailed Description

Description of the class TGPBuilder.

This class allows to build a gaussian processus (kriging) model from a set of observations and a correlation function. It provides tools to search for the optimal parameters of the correlation function and compute the parameters of the model. It wraps the gpModx functions of the "gaussian process Library" (gpLib), by Jean-Marc Martinez.

Remark about normalisation: the observed input variable are automatically normalised such that they belong to [0,1]. The normalisation coefficients are exported into the TKriging model and new input coordinates are thus automatically normalised. As a consequence, the correlation lengths of the correlation function are normalised as well. A special case happens when the deterministic trend is not a constant or a linear combination of the inputs. In that case, the values of the trend coefficients are calculated in the "original" space.

Constructor & Destructor Documentation

◆ TGPBuilder() [1/4]

URANIE::Modeler::TGPBuilder::TGPBuilder ( )

Default Constructor.

Referenced by ClassImp().

◆ TGPBuilder() [2/4]

URANIE::Modeler::TGPBuilder::TGPBuilder ( URANIE::DataServer::TDataServer *  tdsObs,
const char *  varexpinput,
const char *  varexpoutput,
TCorrelationFunction corrFunc,
const char *  trend = "",
Option_t *  option = "" 
)

Construction from a data server, a list of input and output variables and a correlation function object.

Parameters
tdsObsa pointer to the data server containing the observations.
varexpinputa character string containing the names of the input variables, separated by colons (":").
varexpoutputa character string containing the name of the output variable.
corrFunca pointer to the correlation function.
trenda string containing the list of deterministic trend coefficients separated by colons ":" (default = ""). (cf. TGPBuilder::createTrend)
optiona character string containing optional commands. None is avaiblable at this time (default = "").
Exceptions
UErrorExceptionif varexpoutput contains more than one variable name.

◆ TGPBuilder() [3/4]

URANIE::Modeler::TGPBuilder::TGPBuilder ( URANIE::DataServer::TDataServer *  tdsObs,
const char *  varexpinput,
const char *  varexpoutput,
const char *  corrFunc,
const char *  trend = "",
Option_t *  option = "" 
)

Construction from a data server, a list of input and output variables and a correlation function name.

Parameters
tdsObsa pointer to the data server containing the observations.
varexpinputa character string containing the names of the input variables, separated by colons (":").
varexpoutputa character string containing the name of the output variable.
corrFuncname of the correlation function. Available names are: "gauss", "isogauss", "exponential", "maternI", "maternII", "maternIII", "matern3/2", "matern5/2", "matern7/2". The class name of a correlation function is also valid (e.g. "TExponentialCorrFunction", "TMatern32CorrFunction", etc.)
trenda string containing the list of deterministic trend coefficients separated by colons ":" (default = ""). (cf. TGPBuilder::createTrend)
optiona character string containing optional commands. None is available at this time (default = "").
Exceptions
UErrorExceptionif varexpoutput contains more than one variable name.

◆ TGPBuilder() [4/4]

URANIE::Modeler::TGPBuilder::TGPBuilder ( URANIE::DataServer::TDataServer *  tdsModelData)

Construction from a data server containing model parameters.

This constructor works ONLY with a data server generated by the function TGPBuilder::exportGPData. A set of parameters is written in the #TITLE: field of the file header. These parameters are extracted and used to initialise the GPBuilder. The rest of the file contains the observations (inputs and output). After the construction of the TGPBuilder object, the user can call the function TGPBuilder::buildGP to generate the TKriging object.

Parameters
tdsObs(URANIE::DataServer::TDataServer*): a pointer to the data server containing the observations.
Note
a message is displayed if the user needs to define the measurement error variance or the a priori information on the trend parameters.

◆ ~TGPBuilder()

virtual URANIE::Modeler::TGPBuilder::~TGPBuilder ( )
virtual

Default destructor.

Referenced by ClassImp().

Member Function Documentation

◆ buildGP()

URANIE::Modeler::TKriging * URANIE::Modeler::TGPBuilder::buildGP ( Bool_t  computeLooErrors = kTRUE)

Build the gaussian process model.

This function computes the parameters of the gaussian process, necessary to perform new predictions. It then creates a TKriging model and return its pointer.

Note
Calls gpMod0 or gpModh function from gpL
Warning
The user must delete the TKriging object by himself.
Parameters
computeLooErrors(Bool_t): if kTRUE, the Leave One Out errors of the kriging model are computed and stored in the model. (default = kTRUE).

Referenced by ClassImp().

◆ buildWithDuplicateKandCandCorrFunc()

TGPBuilder * URANIE::Modeler::TGPBuilder::buildWithDuplicateKandCandCorrFunc ( )

Referenced by ClassImp(), and ClassImp().

◆ changeLog()

void URANIE::Modeler::TGPBuilder::changeLog ( )
inline

References _blog.

◆ cleanAttributeList()

void URANIE::Modeler::TGPBuilder::cleanAttributeList ( TList *  inputlist)

Clean list of attribute.

Referenced by ClassImp().

◆ computeCovarianceMatrix()

void URANIE::Modeler::TGPBuilder::computeCovarianceMatrix ( bool  print = false)

Computes the covariance matrix of the observed data.

Referenced by ClassImp(), ClassImp(), ClassImp(), and ClassImp().

◆ createTrend()

void URANIE::Modeler::TGPBuilder::createTrend ( )

Create the matrices for the deterministic trend.

Referenced by ClassImp().

◆ exportGPData()

void URANIE::Modeler::TGPBuilder::exportGPData ( const char *  fileName)

Export the minimum information necessary to build a new TKriging object.

This function allows to save into an ASCII file (Salome table format) a set of information and data that can later be read by a TGPBuilder object to recreate the full TKriging object. The saved information is:

  • the type of the correlation function
  • the correlation function parameters (not normalised)
  • the gaussian process variance
  • if needed, the deterministic trend character string
  • if needed, the variance of the measurement error (except if it was user defined)
  • if needed, the regularisation coefficient
  • the observations inputs and outputs.
    Parameters
    fileName(const char*): name of the output file.

Referenced by ClassImp().

◆ findOptimalParameters()

void URANIE::Modeler::TGPBuilder::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.

This function allows to search for the parameters of the correlation function which minimise a given criterion. The function works in two steps:

  1. a screening of the parameter space is performed to find a good starting point for the optimisation;
  2. an optimisation process is performed to refine the results.

To skip the screening step, simply set the screening size to 0. The initial point of the optimisation will then be chosen as the set of default value of the parameters to optimise (see below), or their current values if "reset" if set to kFALSE.

To skip the optimisation step, simply set the optimisation algorithm to "" (empty string). The optimal parameters will then be the one corresponding to the optimal output of the screening.

This function is set to be "user-friendly", and takes care of most of the tedious part of the process. This implies some limitations:

  • hyper-parameters ranges can't be set by the user;
  • most of the optimisation parameters can't be set by the user;
  • currently supported correlation functions are: Gaussian, Isotropic Gaussian, Exponential, MaternI, MaternII, MaternIII, Matern3/2, Matern5/2 and Matern7/2.

The consequence of these limitations is that the "optimal" parameters might not be that optimal. Thus, if the resulting model isn't good enough, using the "not so user-friendly" method might be necessary.

Here is the list of the current "hard coded" ranges and default values for the parameters to optimise:

  • correlation lengths: [1e-3, 10.0], 0.01;
  • exponent: [0.5, 2.0], 1.0;
  • Matern's regularity (\nu): [0.5, 7.5], 1.5;
  • alpha: [1e-6, 10.0], 0.00316;
  • GP variance coefficient: [1e-6, 1e2], 0.01; (this value is multiplied by the observation output's variance).
Warning
This function throws an exception if the NLOpt optimisation library isn't available.
Parameters
criterion(TString): criterion to optimize. Can be "ML" (maximum likelyhood), "ReML" (restrained maximum likelyhood), or "LOO" (leave one out cross-validation). "ML" can't be used in bayesian mode, and ReML can't be used without a deterministic trend. (default = "ML").
screeningSize(Int_t): number of parameters combinaisons to evaluate. An LHS sampling of this size is performed on the parameter space, and evaluated by the cost function. (default = 100).
optimAlgo(TString): name of the optimisation algorithm. Can be: "NelderMead", "Subplexe", "Cobyla", "Bobyqa", "Praxis", "BFGS", "Newton", "MMA", "SLSQP", "VariableMetric". For details about these algorithms, cf. NLopt website: http://ab-initio.mit.edu/wiki/index.php/NLopt_Algorithms. (default = "BFGS").
nbMaxOptimSteps(Int_t): maximum number of optimisation steps. This can be used when the cost function is long to evaluate and the optimisation algorithm is slow to converge. (default = 1000).
reset(Bool_t): if kTRUE, the parameters to optimise are reset to their default values. If kFALSE, the current state of the parameters (correlation lengths, alpha, GP variance, etc.) is used. (default = kTRUE).

Referenced by ClassImp().

◆ getAlpha()

Double_t URANIE::Modeler::TGPBuilder::getAlpha ( )
inline

Return the alpha parameter, quotient of the variance of the measurement error and the variance of the gaussian process.

References alpha.

◆ getCorrelationMatrix()

Double_t * URANIE::Modeler::TGPBuilder::getCorrelationMatrix ( )
inline

Return the array of the correlation matrix of the gaussian process.

References C.

◆ getCorrFunction()

TCorrelationFunction * URANIE::Modeler::TGPBuilder::getCorrFunction ( )
inline

Return a pointer to the correlation function of the gaussian process.

References correlationFunction.

Referenced by ClassImp(), ClassImp(), and ClassImp().

◆ getCovarianceMatrix()

Double_t * URANIE::Modeler::TGPBuilder::getCovarianceMatrix ( )
inline

Return the array of the covariance matrix of the gaussian process.

References K.

◆ getLog()

Bool_t URANIE::Modeler::TGPBuilder::getLog ( )
inline

References _blog.

◆ getMeasurementErrorCovMatrix()

Double_t * URANIE::Modeler::TGPBuilder::getMeasurementErrorCovMatrix ( )
inline

Return the array of the covariance matrix of the measurement errors.

References Vmes.

◆ getNInputs()

Int_t URANIE::Modeler::TGPBuilder::getNInputs ( )
inline

Return the number of input variables.

References _nX.

◆ getNObs()

Int_t URANIE::Modeler::TGPBuilder::getNObs ( )
inline

Return the number of observations used to build the GP model.

References _nS.

◆ getNTrendParams()

Int_t URANIE::Modeler::TGPBuilder::getNTrendParams ( )
inline

Return the number of parameters of the deterministic trend.

References _nP.

◆ getObsInputs()

Double_t * URANIE::Modeler::TGPBuilder::getObsInputs ( )
inline

Return the observations inputs as a flat array.

For M observations of inputs a, b and c, the returned array is organised as follows: [a1, b1, c1, a2, b2, c2, ... , aM, bM, cM]

References xObs.

◆ getObsNormParams()

Double_t * URANIE::Modeler::TGPBuilder::getObsNormParams ( )
inline

Return the normalisation parameters of the inputs as a flat array.

Normalisation parameters are the minimum and maximum values of the input variables. For input variables x1, x2, ... xN, they are stored as follows: [x1min, x1max, x2min, x2max, ... , xNmin, xNmax].

References xNormParams.

◆ getObsOutput()

Double_t * URANIE::Modeler::TGPBuilder::getObsOutput ( )
inline

Return the observations outputs.

References yObs.

◆ getRegularisation()

Double_t URANIE::Modeler::TGPBuilder::getRegularisation ( )
inline

Return the regularisation coefficient.

References reg.

◆ getReturnOptim()

int URANIE::Modeler::TGPBuilder::getReturnOptim ( )
inline

Return the return of the optimisation procedure.

References _retOptim.

◆ getSeed()

Int_t URANIE::Modeler::TGPBuilder::getSeed ( )
inline

Return the seed number.

References _nSeed.

◆ getTolerance()

Double_t URANIE::Modeler::TGPBuilder::getTolerance ( )
inline

Return the tolerance used to stop the optimisation process

References tolerance.

◆ getTrendFormula()

TObjArray * URANIE::Modeler::TGPBuilder::getTrendFormula ( )
inline

Return the deterministic trend list of formulas.

References trendCoefList.

◆ getTrendMatrix()

Double_t * URANIE::Modeler::TGPBuilder::getTrendMatrix ( )
inline

Return the array of the deterministic trend matrix.

References H.

◆ getTrendPriorCovMat()

Double_t * URANIE::Modeler::TGPBuilder::getTrendPriorCovMat ( )
inline

Return the array of the a priori covariance matrix of the trend parameters.

References QPrior.

◆ getTrendPriorMeans()

Double_t * URANIE::Modeler::TGPBuilder::getTrendPriorMeans ( )
inline

Return the array of the a priori means of the trend parameters.

References BetaPrior.

◆ getTrendString()

TString URANIE::Modeler::TGPBuilder::getTrendString ( )
inline

Return the deterministic trend character string.

References trendString.

◆ getVariance()

Double_t URANIE::Modeler::TGPBuilder::getVariance ( )
inline

Return the variance of the gaussian process.

References sigma2.

◆ hasCorrFunction()

Bool_t URANIE::Modeler::TGPBuilder::hasCorrFunction ( )
inline

Return kTRUE if a correlation function is defined, kFALSE otherwise.

References bHas_corrFunction.

Referenced by ClassImp(), ClassImp(), and ClassImp().

◆ hasMeasurementError()

Bool_t URANIE::Modeler::TGPBuilder::hasMeasurementError ( )
inline

Return kTRUE if a measurement error should be considered to build the GP.

References bHas_measurement_error.

Referenced by ClassImp(), ClassImp(), and ClassImp().

◆ hasTrend()

Bool_t URANIE::Modeler::TGPBuilder::hasTrend ( )
inline

Return kTRUE if a deterministic trend has been defined.

References bHas_trend.

Referenced by ClassImp(), and ClassImp().

◆ hasVarianceInFactor()

Bool_t URANIE::Modeler::TGPBuilder::hasVarianceInFactor ( )
inline

Return kTRUE if the variance of the GP is determined analytically.

References bVariance_in_factor.

Referenced by ClassImp(), ClassImp(), and ClassImp().

◆ initMatrices()

void URANIE::Modeler::TGPBuilder::initMatrices ( )

Initialisation of the matrices.

Referenced by ClassImp().

◆ initVariables()

void URANIE::Modeler::TGPBuilder::initVariables ( )

Initialisation of the GPBuilder variables.

Referenced by ClassImp().

◆ LSort()

void URANIE::Modeler::TGPBuilder::LSort ( Double_t *  costValues,
Int_t *  minIndexes,
Int_t  nbthreads,
Int_t  n 
)
private

Referenced by ClassImp().

◆ printLog()

virtual void URANIE::Modeler::TGPBuilder::printLog ( Option_t *  option = "")
virtual

Referenced by ClassImp().

◆ setAlpha()

void URANIE::Modeler::TGPBuilder::setAlpha ( Double_t  a)
inline

Set the value of the alpha parameter.

This parameter is used when the variance of the measurement error can be expressed as sigma2*alpha, where sigma2 is the variance of the GP.

Note
calling this function automatically sets the flags bHas_measurement_error and bVariance_in_factor to kTRUE.
Parameters
athe value of the parameter alpha

References alpha, bHas_measurement_error, and bVariance_in_factor.

Referenced by ClassImp(), ClassImp(), ClassImp(), and ClassImp().

◆ setCorrFunction() [1/2]

void URANIE::Modeler::TGPBuilder::setCorrFunction ( const char *  corrFuncType)

Set the correlation function of the gaussian process from a character string.

The name of the correlation function can either be its class name (e.g. "TExponentialCorrFunction", "TMatern32CorrFunction", etc.), or one of these (case insensitive) keywords: "gauss", "isogauss", "exponential", "maternI", "maternII", "maternIII", "matern3/2", "matern5/2", "matern7/2".

Parameters
corrFuncname of the correlation function.

◆ setCorrFunction() [2/2]

void URANIE::Modeler::TGPBuilder::setCorrFunction ( TCorrelationFunction corrFunc)

Set the correlation function of the gaussian process.

The internal correlation function is cloned from the correlation function provided as parameter. This means that the later isn't modified by the GPBuilder.

Parameters
corrFunca pointer to the correlation function of the gaussian process.

Referenced by ClassImp().

◆ setHasMeasurementError()

void URANIE::Modeler::TGPBuilder::setHasMeasurementError ( Bool_t  has_measurement_error = kTRUE)
inline

Indicate if a measurement error should be taken into account or not.

This function sets the bHas_measurement_error flag which indicates wether a measurement error should be considered to build the Gaussian Process. If kTRUE, and the user doesn't provide a covariance matrix (using setMeasurementErrorCovMatrix or setMeasurementErrorVariance), the parameter alpha will be a parameter of the cost function, and the variance of the measurment error will be sigma2*alpha.

Parameters
has_measurement_errorflag to tell if measurement error is taken into account. (default = kTRUE)

References bHas_measurement_error.

Referenced by ClassImp().

◆ setKandCandCorrFunc()

void URANIE::Modeler::TGPBuilder::setKandCandCorrFunc ( Double_t *  inK,
Double_t *  inC,
TCorrelationFunction corrFunc 
)
inline

◆ setLengthMaxScreeningSF()

void URANIE::Modeler::TGPBuilder::setLengthMaxScreeningSF ( double  sf)

Set the scale factor to reduce the maximum value of the correlation length during screening.

This function changes the value of the scale factor that can be used to reduce the maximum value of the correlation length during the screening procedure .

Parameters
sfa double (default is 1) that is the scale factor to define a new lengthMax specifically for the screening, as lengthMaxScr = sf * lengthMax. Exception if sf < 0 or sf > 1

Referenced by ClassImp().

◆ setLengthRange()

void URANIE::Modeler::TGPBuilder::setLengthRange ( double  themin,
double  themax 
)

Set the Minimun and maximum boundaries for the correlation length.

This function changes the default value of the Minimun and maximum boundaries used to define the space definition of the correlation length.

Parameters
themina double that is the min value allowed for the correlation length (default being 1e-3).
themaxa double that is the max value allowed for the correlation length (default being 10).

Referenced by ClassImp().

◆ setLog()

void URANIE::Modeler::TGPBuilder::setLog ( )
inline

References _blog.

◆ setMeasurementErrorCovMatrix()

void URANIE::Modeler::TGPBuilder::setMeasurementErrorCovMatrix ( Double_t *  covMes)

Set the covariance matrix of the measurement errors.

The data from covMes are copied into Vmes. covMes MUST be a valid covariance matrix.

Note
Calling this function automatically sets the flag bHas_measurement_error to kTRUE and the flag bVariance_in_factor to kFALSE, meaning that the GP variance is an input parameter of the cost function.
Parameters
covMes"flat" array of nS*nS doubles.

Referenced by ClassImp().

◆ setMeasurementErrorVariance()

void URANIE::Modeler::TGPBuilder::setMeasurementErrorVariance ( Double_t  varMes)

Set the variance of the measurement errors.

The diagonal of the matrix Vmes is filled with varMes.

Note
Calling this function automatically sets the flag bHas_measurement_error to kTRUE and the flag bVariance_in_factor to kFALSE, meaning that the GP variance is an input parameter of the cost function.
Parameters
varMesvariance of the measurement error.

Referenced by ClassImp().

◆ setPriorData()

void URANIE::Modeler::TGPBuilder::setPriorData ( Double_t *  betaPrior,
Double_t *  qPrior 
)

Set the a priori mean and variance of the trend coefficients.

This function should not be used if no trend was defined. It allows to set the a priori mean and variance of the trend coefficients.

Warning
the size of the array is no checked !
Note
Calling this function automatically sets the flag bUse_prior to kTRUE.
Parameters
betaPriora priori mean vector for the deterministic trend parameters (size: _nP)
qPriora priori covariance matrix for the deterministic trend parameters (size: _nP x _nP)
Exceptions
UranieErrorif no trend was defined.

Referenced by ClassImp().

◆ setRegularisation()

void URANIE::Modeler::TGPBuilder::setRegularisation ( Double_t  regularisation)
inline

Set the regularisation coefficient.

The regularisation coefficient is a small (~ 1e-12 to 1e-6) value that is added to the diagonal of the correlation matrix when it is poorly conditionned. It should only be set when matrix decomposition problems arise.

Parameters
regularisationvalue of the regularisation coefficient. It is typically between 1e-12 and 1e-6.

References reg.

Referenced by ClassImp().

◆ setSeed()

void URANIE::Modeler::TGPBuilder::setSeed ( Int_t  nSeed = 0)

set a seed number

This function set the value of _nSeed which indicates the seed number to use by the random generator

Parameters
nSeedthe seed number to use by the random generator

Referenced by ClassImp().

◆ setTolerance()

void URANIE::Modeler::TGPBuilder::setTolerance ( Double_t  tol)
inline

Set the tolerance used to stop the optimisation process.

Parameters
varvalue of the tolerance.

References tolerance.

◆ setTrendString()

void URANIE::Modeler::TGPBuilder::setTrendString ( const char *  trend)
inline

Set the deterministic trend character string.

This function creates the deterministic trend matrix based on a list of trend coefficients. This list is created from a character string. It can be either a keyword describing the trend type ("const" for a constant value or "linear" for a linear combination of the inputs) or a set of ROOT style formulas of the data server attributes, separated by ":" (cf. ROOT documentation of TFormula::Analyze). For example, if we have attributes x1 and x2, the following strings are valid trends:

  • "const" will create a constant trend of the form: a0
  • "linear" will create a trend of the form: a0 + a1*X1 + a2*X2, where X1 et X2 are the normalised values of x1 and x2
  • "1:x1:x1*x1:log(x2)" will create a trend of the form: a0 + a1*x1 + a2*x1^2 + a3*log(x2)
    Warning
    if a coefficient formula doesn't contain any reference to an attribute of the data server, then it MUST be a constant expression (like "1", "log(2)", "pi", etc.). Otherwise, an error will occure.
    Parameters
    trendthe character string containing the trend type ("const" or "linear"), or a the list of coefficients separated by colons ":".

References trendString.

◆ setUsePrior()

void URANIE::Modeler::TGPBuilder::setUsePrior ( Bool_t  use_prior = kTRUE)

Set wether or not to use a priori knowledge on the deterministic trend parameters.

This function should not be used if no trend and no prior data were defined.

Parameters
use_priorif kTRUE, the a priori knowledge is used (default = kTRUE).
Exceptions
UranieErrorif no trend was defined.

Referenced by ClassImp().

◆ setVariance()

void URANIE::Modeler::TGPBuilder::setVariance ( Double_t  var)
inline

Set the variance of the gaussian process.

Parameters
varvalue of the variance.

References sigma2.

Referenced by ClassImp(), ClassImp(), ClassImp(), and ClassImp().

◆ unsetLog()

void URANIE::Modeler::TGPBuilder::unsetLog ( )
inline

References _blog.

◆ updateObservations()

void URANIE::Modeler::TGPBuilder::updateObservations ( )

Update the GP builder matrices.

This function shall be called when the data server containing the observations is modified. This is updating the matrices but doesn't change the GP hyper-parameters (correlation lengths, variance, etc.)

Referenced by ClassImp().

◆ usePrior()

Bool_t URANIE::Modeler::TGPBuilder::usePrior ( )
inline

Return kTRUE if prior knowledge on the deterministic trend is used (bayesian approach).

References bUse_prior.

Referenced by ClassImp().

Member Data Documentation

◆ _blog

Bool_t URANIE::Modeler::TGPBuilder::_blog
private

Boolean for edit the log.

Referenced by changeLog(), ClassImp(), getLog(), setLog(), and unsetLog().

◆ _lengthMaxScreeningSF

Double_t URANIE::Modeler::TGPBuilder::_lengthMaxScreeningSF

Scale factor to change the Maximum value of the correlation length during screening procedure: new length is SF * lengthMax (0 < SF <=1)

Referenced by ClassImp().

◆ _modlengthMax

Double_t URANIE::Modeler::TGPBuilder::_modlengthMax

Maximum value of the correlation length if default behaviour is requested to be changed by the user.

Referenced by ClassImp().

◆ _modlengthMin

Double_t URANIE::Modeler::TGPBuilder::_modlengthMin

Minimum value of the correlation length if default behaviour is requested to be changed by the user.

Referenced by ClassImp().

◆ _nP

Int_t URANIE::Modeler::TGPBuilder::_nP

Number of coefficients of the deterministic trend.

Referenced by ClassImp(), and getNTrendParams().

◆ _nS

Int_t URANIE::Modeler::TGPBuilder::_nS

Number of observations.

Referenced by ClassImp(), and getNObs().

◆ _nSeed

Int_t URANIE::Modeler::TGPBuilder::_nSeed

number of seed for pseudo-random generator

Referenced by ClassImp(), and getSeed().

◆ _nX

Int_t URANIE::Modeler::TGPBuilder::_nX

Number of input attributes.

Referenced by ClassImp(), and getNInputs().

◆ _nY

Int_t URANIE::Modeler::TGPBuilder::_nY

Number of output attributes (should be equal to 1)

Referenced by ClassImp().

◆ _retOptim

Int_t URANIE::Modeler::TGPBuilder::_retOptim

optimisation results

Referenced by ClassImp(), and getReturnOptim().

◆ _sInputAttributes

TString URANIE::Modeler::TGPBuilder::_sInputAttributes

input attributes names separated by ":"

Referenced by ClassImp().

◆ _sOutputAttributes

TString URANIE::Modeler::TGPBuilder::_sOutputAttributes

output attribute name

Referenced by ClassImp().

◆ _tds

URANIE::DataServer::TDataServer* URANIE::Modeler::TGPBuilder::_tds

Pointer to the data server containing the observations.

Referenced by ClassImp().

◆ alpha

Double_t URANIE::Modeler::TGPBuilder::alpha

Quotient of the variance of the measurment error and the variance of the gaussian process.

Referenced by ClassImp(), getAlpha(), and setAlpha().

◆ BetaPrior

Double_t* URANIE::Modeler::TGPBuilder::BetaPrior

a priori mean vector for the deterministic trend parameters (size: _nP)

Referenced by ClassImp(), and getTrendPriorMeans().

◆ bHas_corrFunction

Bool_t URANIE::Modeler::TGPBuilder::bHas_corrFunction

kFALSE if correlationFunction == NULL, kTRUE otherwise

Referenced by ClassImp(), and hasCorrFunction().

◆ bHas_measurement_error

Bool_t URANIE::Modeler::TGPBuilder::bHas_measurement_error

if kTRUE, the construction of the GP will take into account the existence of a measurement error variance.

Referenced by ClassImp(), hasMeasurementError(), setAlpha(), and setHasMeasurementError().

◆ bHas_trend

Bool_t URANIE::Modeler::TGPBuilder::bHas_trend

kFALSE if no deterministic trend is defined, kTRUE otherwise.

Referenced by ClassImp(), and hasTrend().

◆ bIsbuildWithDuplicateKandCandCorrFunc

Bool_t URANIE::Modeler::TGPBuilder::bIsbuildWithDuplicateKandCandCorrFunc

kTRUE if build with duplicate K, C and CorrFunc, otherwise kFALSE

Referenced by ClassImp(), and setKandCandCorrFunc().

◆ bUse_normalisation

Bool_t URANIE::Modeler::TGPBuilder::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().

◆ bUse_prior

Bool_t URANIE::Modeler::TGPBuilder::bUse_prior

kFALSE if no prior is used, kTRUE otherwise

Referenced by ClassImp(), and usePrior().

◆ bVariance_in_factor

Bool_t URANIE::Modeler::TGPBuilder::bVariance_in_factor

if kTRUE, GP variance is calculated analytically. Otherwise, it is a parameter of the cost function.

Referenced by ClassImp(), hasVarianceInFactor(), and setAlpha().

◆ C

Double_t* URANIE::Modeler::TGPBuilder::C

The correlation matrix (dimension: _nS * _nS)

Referenced by ClassImp(), ClassImp(), getCorrelationMatrix(), and setKandCandCorrFunc().

◆ correlationFunction

TCorrelationFunction* URANIE::Modeler::TGPBuilder::correlationFunction

The correlation function.

Referenced by ClassImp(), getCorrFunction(), and setKandCandCorrFunc().

◆ H

Double_t* URANIE::Modeler::TGPBuilder::H

Deterministic trend matrix (dimension: _nS * _nP)

Referenced by ClassImp(), and getTrendMatrix().

◆ K

Double_t* URANIE::Modeler::TGPBuilder::K

The covariance matrix (dimension: _nS * _nS)

Referenced by ClassImp(), ClassImp(), getCovarianceMatrix(), and setKandCandCorrFunc().

◆ QPrior

Double_t* URANIE::Modeler::TGPBuilder::QPrior

a priori variance matrix for the deterministic trend parameters (size: _nP x _nP)

Referenced by ClassImp(), and getTrendPriorCovMat().

◆ reg

Double_t URANIE::Modeler::TGPBuilder::reg

regularisation parameter, used to improve numerical stability of matrix operations

Referenced by ClassImp(), getRegularisation(), and setRegularisation().

◆ sigma2

Double_t URANIE::Modeler::TGPBuilder::sigma2

Variance of the gaussian process.

Referenced by ClassImp(), getVariance(), and setVariance().

◆ tolerance

Double_t URANIE::Modeler::TGPBuilder::tolerance

Tolerance value used to stop the optimisation process.

Referenced by ClassImp(), getTolerance(), and setTolerance().

◆ trendCoefList

TObjArray* URANIE::Modeler::TGPBuilder::trendCoefList

List of the formulas of the trend coefficients (_nP elements)

Referenced by ClassImp(), and getTrendFormula().

◆ trendString

TString URANIE::Modeler::TGPBuilder::trendString

Character string containing the trend parameters separated by ":" or a trend type ("Const" or "Linear")

Referenced by ClassImp(), getTrendString(), and setTrendString().

◆ Vmes

Double_t* URANIE::Modeler::TGPBuilder::Vmes

The covariance matrix of the measurement errors (dimension: _nS * _nS)

Referenced by ClassImp(), and getMeasurementErrorCovMatrix().

◆ xNormParams

Double_t* URANIE::Modeler::TGPBuilder::xNormParams

The matrix of normalisation parameters. Contains min and max value of each input variable: [x1Min, x1Max, x2Min, x2Max,...] (size: _nX * 2)

Referenced by ClassImp(), and getObsNormParams().

◆ xObs

Double_t* URANIE::Modeler::TGPBuilder::xObs

The matrix of normalised inputs (dimension: _nX * _nS)

Referenced by ClassImp(), and getObsInputs().

◆ yObs

Double_t* URANIE::Modeler::TGPBuilder::yObs

Array representation of the outputs (size: _nS)

Referenced by ClassImp(), and getObsOutput().