Documentation / Manuel développeur
Modules disponibles
Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,  
Uranie / Calibration v4.9.0
|
#include <TMetropHasting.h>
Public Member Functions | |
Constructor and Destructor | |
TMetropHasting (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int nS=100, Option_t *option="lhs") | |
TMetropHasting (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), const char *varexpinput, const char *varexpoutput, int ns=100, Option_t *option="") | |
Default Calibration constructor with the function argument: it contains the assessor to be used. | |
TMetropHasting (URANIE::DataServer::TDataServer *tds, const char *fcn, const char *varexpinput, const char *varexpoutput, int ns=100, Option_t *option="") | |
Default Calibration constructor with the function argument: it contains the assessor to be used. | |
TMetropHasting (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, int ns=100, Option_t *option="") | |
Default Calibration constructor with the code argument: it contains the assessor to be used. | |
virtual | ~TMetropHasting () |
Default destructor. | |
evaluate the events | |
void | computeParameters (Option_t *option="") |
Generate the sample. | |
void | checktdsParContent () |
Set reference information and initialisation methods | |
void | setBurnin (int burn) |
void | setLag (int lag) |
void | setNbDump (int nbDump) |
void | setAcceptationRatioRange (double lower, double higher) |
string | getDefaultCut () |
void | clearDefaultCut () |
void | setInitialisation (int n, double *values, double *standDev) |
Initialise the parameters. | |
void | setInitialisation (vector< double > values, vector< double > standDev) |
Initialise the parameters. | |
void | setDistanceAndReference (const char *funcName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="") |
Set the distance function and some needed informations. | |
void | setDistanceAndReference (URANIE::Calibration::TDistanceFunction *distFunc, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="") |
Set the distance function and some needed informations. | |
Visualisation and diagnostic | |
void | drawTrace (TString sTitre, const char *variable="*", const char *select="1>0", Option_t *option="") |
Draws the evolution of parameters as a function of the iterator. | |
void | drawAcceptationRatio (TString sTitre, const char *variable="*", const char *select="1>0", Option_t *option="") |
Draws the evolution of acceptation ratio as a function of the iterator. | |
void | drawParameters (TString sTitre, const char *variable="*", const char *select="1>0", Option_t *option="") |
Draws the parameters as distributions The estimateParamters method has computed the parameters, so this one plots the parameters as distributions. | |
void | getAutoCorrelation (vector< int > l, vector< double > *out, int cut=0) |
Compute the autocorrelation. | |
Printing Log | |
virtual void | printLog (Option_t *option="") |
Prints the log. | |
Public Member Functions inherited from URANIE::Calibration::TCalibration | |
TCalibration (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns=1, Option_t *option="") | |
Default constructor with the runner argument: it contains the assessor to be used. | |
TCalibration (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), const char *varexpinput, const char *varexpoutput, int ns=1, Option_t *option="") | |
Default Calibration constructor with the function argument: it contains the assessor to be used. | |
TCalibration (URANIE::DataServer::TDataServer *tds, const char *fcn, const char *varexpinput, const char *varexpoutput, int ns=1, Option_t *option="") | |
Default Calibration constructor with the function argument: it contains the assessor to be used. | |
TCalibration (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, int ns=1, const char *option="") | |
Default Calibration constructor with the code argument: it contains the assessor to be used. | |
virtual | ~TCalibration () |
Default destructor. | |
void | estimateParameters (Option_t *option="") |
void | estimateCustomResidues (string resName, int n_theta, double *theta) |
void | checkReference (URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight) |
Check the consistency of the formation provided. | |
void | drawResidues (TString sTitre, const char *variable="*", const char *select="1>0", Option_t *option="") |
void | setLog () |
void | unsetLog () |
void | changeLog () |
Bool_t | getLog () |
Int_t | getID () |
void | setSeed (UInt_t nval) |
Set the seed of the random generator if one is used. | |
UInt_t | getSeed () |
Get the seed of the random generator if one is used. | |
const char * | getMethodName () |
Get the method name. | |
void | setObservationCovarianceMatrix (TMatrixD &mat) |
Set the observatiton covariance matrix. | |
URANIE::Calibration::TDistanceFunction * | getDistanceFunction () |
Return the distance function. | |
Int_t | getNPar () |
Get the number of parameters to be calibrated. | |
URANIE::DataServer::TDataServer * | getEvaluationTDS () |
Get the tds in which evaluation will be performed. | |
Public Attributes | |
vector< double > | _values |
Vector of values to be tested. | |
vector< double > | _vstd |
Vector of standard deviation. | |
int | _burnin |
The warm-up or burn-in. | |
int | _lag |
The lag value. | |
int | _nbDump |
Frequency to which the algo dump a line. | |
double | _lowAccRange |
loweracceptation ratio bound to increase _vstd | |
double | _higAccRange |
higheracceptation ratio bound to decrease _vstd | |
bool | _mbGDVersion |
Use the Guillaume Damblin. | |
bool | _bcleaningAtt |
Do not store the underlying att. | |
Public Attributes inherited from URANIE::Calibration::TCalibration | |
URANIE::DataServer::TDataServer * | _tdsPar |
TDS containing parameters properties (parameters that should be calibrated) | |
URANIE::DataServer::TDataServer * | _tdsObs |
TDS containing observations used for calibration. | |
URANIE::DataServer::TDataServer * | _tdsEval |
TDS containing a priori / a posteriori evaluations. | |
ELauncher | _nLauncher |
The type of launcher. | |
TString | _sFunctionName |
The Name of the evaluatuor. | |
TString | _sEI |
The Name of input. | |
TString | _sEO |
The Name of output. | |
URANIE::Launcher::TCode * | _code |
The tcode. | |
URANIE::Relauncher::TRun * | _run |
Pointer to the runner to be used. | |
void(* | _pFunction )(double *, double *) |
Function pointer. | |
vector< URANIE::DataServer::TStochasticAttribute * > | _vatt |
internal vector of stochastic attribute for some methods | |
Protected Member Functions | |
private methods | |
void | logPriorPdf (double &ret) |
Logarithm of the prior. | |
void | parseOption (Option_t *option="") |
Read the possible options. | |
Protected Member Functions inherited from URANIE::Calibration::TCalibration | |
void | checkCanvasCreation (bool newcan) |
Create a canvas if needed. | |
void | initInputs () |
Initialise some common inputs. | |
void | initResults (vector< string > *ParsedLines) |
Initialise some common inputs. | |
void | computeAPosterioriForDistribution () |
Compute the a posteriori residual for many-solutions method. | |
void | setMethodName (const char *str) |
Set the Method name. | |
Additional Inherited Members | |
Public Types inherited from URANIE::Calibration::TCalibration | |
enum | ELauncher { kCode , kFunction , kRun , kUnknown } |
Protected Attributes inherited from URANIE::Calibration::TCalibration | |
URANIE::Calibration::TDistanceFunction * | _dFunc |
Pointer to chosen distance function. | |
URANIE::DataServer::TDSNtupleD * | _evTuple |
Pointer to the eval ntuple. | |
TList * | _listOfParameters |
List of the parameters to be calibrated. | |
TCanvas * | _canvas |
Canvas object to deal with. | |
TObjArray * | _drawingGarbageCollector |
Garbage collector for prints. | |
int | _nSam |
The number of sample in a posteriori distributions. | |
int | _nObs |
The number of observations in the reference database. | |
int | _nIterMax |
The maximum number of iteration allowed (meaning total number of code estimation is _nIterMax * _nObs) | |
int | _nPar |
Dimension of the parameters. | |
int | _nVar |
Dimension of the output and references to be compared with. | |
int | _nSeed |
The seed of the random generator. | |
TString | _sMethodName |
The method name. | |
TString | _referenceName |
The reference name. | |
TString | _outputName |
The output name. | |
vector< string > | _vrefName |
The reference names. | |
vector< string > | _voutName |
The output names. | |
TString | _weightName |
The weight name. | |
TMatrixD | _mObsCovMat |
Observation Covariance matrix. | |
bool | _buseMatrix |
Use matrix instead of vectors in the Distance Function. | |
bool | _bsaveAll |
Whether all evaluations should be saved, not only a priori and a posteriori. | |
bool | _bdontKeepAgreement |
Remove the agreement attribute from the tdsPar object. | |
bool | _buseMode |
Use Mode instead of Mean. | |
Bool_t | _blog |
Boolean for edit the log. | |
Constructor & Destructor Documentation
◆ TMetropHasting() [1/4]
URANIE::Calibration::TMetropHasting::TMetropHasting | ( | URANIE::DataServer::TDataServer * | tds, |
URANIE::Relauncher::TRun * | run, | ||
int | nS = 100 , |
||
Option_t * | option = "lhs" |
||
) |
Constructor usual where nS is set to be the number of sample produced by the algorithm.
- Parameters
-
tds the pointer of the TDataServer. It must contains attribuutes that inherit from TStochasticAttribute run the runner that would be used to produce evaluations nS size of the sample to be generated [100 by default]
Referenced by ClassImp().
◆ TMetropHasting() [2/4]
URANIE::Calibration::TMetropHasting::TMetropHasting | ( | URANIE::DataServer::TDataServer * | tds, |
void(*)(Double_t *, Double_t *) | fcn, | ||
const char * | varexpinput, | ||
const char * | varexpoutput, | ||
int | ns = 100 , |
||
Option_t * | option = "" |
||
) |
Default Calibration constructor with the function argument: it contains the assessor to be used.
- Parameters
-
tds : the dataserver that contains no data but only one attribute per parameter to be calibrated fcn : the pointer to a function ns : number of sample to be generated (depending on the method) varexpinput : the input variable for the function in the correct order (both input and parameters) varexpoutput : the output of the function in the correct order
◆ TMetropHasting() [3/4]
URANIE::Calibration::TMetropHasting::TMetropHasting | ( | URANIE::DataServer::TDataServer * | tds, |
const char * | fcn, | ||
const char * | varexpinput, | ||
const char * | varexpoutput, | ||
int | ns = 100 , |
||
Option_t * | option = "" |
||
) |
Default Calibration constructor with the function argument: it contains the assessor to be used.
- Parameters
-
tds : the dataserver that contains no data but only one attribute per parameter to be calibrated fcn : the name of the function ns : number of sample to be generated (depending on the method) varexpinput : the input variable for the function in the correct order (both input and parameters) varexpoutput : the output of the function in the correct order
◆ TMetropHasting() [4/4]
URANIE::Calibration::TMetropHasting::TMetropHasting | ( | URANIE::DataServer::TDataServer * | tds, |
URANIE::Launcher::TCode * | fcode, | ||
int | ns = 100 , |
||
Option_t * | option = "" |
||
) |
Default Calibration constructor with the code argument: it contains the assessor to be used.
- Parameters
-
tds : the dataserver that contains no data but only one attribute per parameter to be calibrated code : the code object that will be runned ns : number of sample to be generated (depending on the method)
◆ ~TMetropHasting()
|
virtual |
Default destructor.
Referenced by ClassImp().
Member Function Documentation
◆ checktdsParContent()
|
virtual |
Check the content of the input dataserver
Implements URANIE::Calibration::TCalibration.
Referenced by ClassImp().
◆ clearDefaultCut()
|
inline |
References setBurnin(), and setLag().
◆ computeParameters()
|
virtual |
Generate the sample.
- Parameters
-
option possble options are - "clean": remove internal branche from the tdsPar object at the end
- "gd": another implementation of this algorithm
Implements URANIE::Calibration::TCalibration.
Referenced by ClassImp().
◆ drawAcceptationRatio()
void URANIE::Calibration::TMetropHasting::drawAcceptationRatio | ( | TString | sTitre, |
const char * | variable = "*" , |
||
const char * | select = "1>0" , |
||
Option_t * | option = "" |
||
) |
Draws the evolution of acceptation ratio as a function of the iterator.
- Parameters
-
sTitre title of the plot if specified variable list of variables that should be plotted select select that can be applied on the parameters or the iterator (burnin, selection for autocorrelation) option possible option for this plot - "nonewcanvas": don't create a new canvas, but use the one active
Referenced by ClassImp().
◆ drawParameters()
|
virtual |
Draws the parameters as distributions The estimateParamters method has computed the parameters, so this one plots the parameters as distributions.
- Parameters
-
sTitre title of the plots if specified variable list of variables that should be plotted select select that can be applied on the parameters or the iterator (burnin, selection for autocorrelation) option possible option for this plot - "nonewcanvas": don't create a new canvas, but use the one active
- "vertical" : if there are several parameters to be shown put them on top of another instead of side by side
Reimplemented from URANIE::Calibration::TCalibration.
Referenced by ClassImp().
◆ drawTrace()
void URANIE::Calibration::TMetropHasting::drawTrace | ( | TString | sTitre, |
const char * | variable = "*" , |
||
const char * | select = "1>0" , |
||
Option_t * | option = "" |
||
) |
Draws the evolution of parameters as a function of the iterator.
- Parameters
-
sTitre title of the plot if specified variable list of variables that should be plotted select select that can be applied on the parameters or the iterator (burnin, selection for autocorrelation) option possible option for this plot - "nonewcanvas": don't create a new canvas, but use the one active
Referenced by ClassImp().
◆ getAutoCorrelation()
void URANIE::Calibration::TMetropHasting::getAutoCorrelation | ( | vector< int > | l, |
vector< double > * | out, | ||
int | cut = 0 |
||
) |
Compute the autocorrelation.
- Parameters
-
l vector of steps (reading 1 events out of l) out vector of estimations (one per step) cut burn-up value (events ignored from 0 to cut)
Referenced by ClassImp().
◆ getDefaultCut()
|
inline |
References _burnin, _lag, and URANIE::Calibration::TCalibration::_tdsPar.
Referenced by ClassImp().
◆ logPriorPdf()
|
protected |
Logarithm of the prior.
Referenced by ClassImp().
◆ parseOption()
|
protectedvirtual |
Read the possible options.
On top of the other, one can call
- "savealleval": keep all estimation in the _tdsEval (Extremely dangerous as it will keep _nIter branches with _nObs estimations. This option might slow down the process.
- "usematrix" : option use to specify that computation is done with matrices instead of vector<vector<double> >. Shoud be used only if one is using an home-made TDistanceFunction
Reimplemented from URANIE::Calibration::TCalibration.
Referenced by ClassImp().
◆ printLog()
|
virtual |
◆ setAcceptationRatioRange()
void URANIE::Calibration::TMetropHasting::setAcceptationRatioRange | ( | double | lower, |
double | higher | ||
) |
Referenced by ClassImp().
◆ setBurnin()
|
inline |
References _burnin.
Referenced by clearDefaultCut().
◆ setDistanceAndReference() [1/2]
|
virtual |
Set the distance function and some needed informations.
- Parameters
-
funcName name of the distance function chosen among the already implemented ones: - "LS" for least square
- "weighterLS" for weighted least square
- "relativeLS" for relative least square
- "L1" for L1 norm
- "Mahalanobis" for Mahalanobis distance function
tdsRef : the dataserver that contains all needed information detailled below - input: the input attributes used of the assessors in run
- reference: the reference attributes that will be used to compared the newly done estimations for every iterations
- weight: a weight to reweight every event, one-by-one (an experimental uncertainty for instance)
Reimplemented from URANIE::Calibration::TCalibration.
Referenced by ClassImp().
◆ setDistanceAndReference() [2/2]
|
virtual |
Set the distance function and some needed informations.
- Parameters
-
distFunc a pointer to a TDistanceFunction object (not usual, recommended only when dealing with a home-made distance unction compiled on the spot) tdsRef : the dataserver that contains all needed information detailled below - input: the input attributes used of the assessors in run
- reference: the reference attributes that will be used to compared the newly done estimations for every iterations
- weight: a weight to reweight every event, one-by-one (an experimental uncertainty for instance)
Reimplemented from URANIE::Calibration::TCalibration.
◆ setInitialisation() [1/2]
void URANIE::Calibration::TMetropHasting::setInitialisation | ( | int | n, |
double * | values, | ||
double * | standDev | ||
) |
Initialise the parameters.
- Parameters
-
n number of parameters values array of values (size must equal to n) standDev array of standard deviation (size must equal to n)
Referenced by ClassImp().
◆ setInitialisation() [2/2]
void URANIE::Calibration::TMetropHasting::setInitialisation | ( | vector< double > | values, |
vector< double > | standDev | ||
) |
Initialise the parameters.
- Parameters
-
values vector of values (size must equal to n) standDev vector of standard deviation (size must equal to n)
◆ setLag()
|
inline |
References _lag.
Referenced by clearDefaultCut().
◆ setNbDump()
|
inline |
References _nbDump.
Member Data Documentation
◆ _bcleaningAtt
bool URANIE::Calibration::TMetropHasting::_bcleaningAtt |
Do not store the underlying att.
Referenced by ClassImp().
◆ _burnin
int URANIE::Calibration::TMetropHasting::_burnin |
The warm-up or burn-in.
Referenced by getDefaultCut(), and setBurnin().
◆ _higAccRange
double URANIE::Calibration::TMetropHasting::_higAccRange |
higheracceptation ratio bound to decrease _vstd
Referenced by ClassImp().
◆ _lag
int URANIE::Calibration::TMetropHasting::_lag |
The lag value.
Referenced by getDefaultCut(), and setLag().
◆ _lowAccRange
double URANIE::Calibration::TMetropHasting::_lowAccRange |
loweracceptation ratio bound to increase _vstd
Referenced by ClassImp().
◆ _mbGDVersion
bool URANIE::Calibration::TMetropHasting::_mbGDVersion |
Use the Guillaume Damblin.
Referenced by ClassImp().
◆ _nbDump
int URANIE::Calibration::TMetropHasting::_nbDump |
Frequency to which the algo dump a line.
Referenced by ClassImp(), and setNbDump().
◆ _values
vector<double> URANIE::Calibration::TMetropHasting::_values |
Vector of values to be tested.
Referenced by ClassImp().
◆ _vstd
vector<double> URANIE::Calibration::TMetropHasting::_vstd |
Vector of standard deviation.
Referenced by ClassImp().