English Français

Documentation / Manuel développeur

Modules disponibles

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

#include <TPMCABC.h>

Inheritance diagram for URANIE::Calibration::TPMCABC:
Collaboration diagram for URANIE::Calibration::TPMCABC:

Public Member Functions

Constructor and Destructor
 TPMCABC (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns=1, Option_t *option="")
 
 TPMCABC (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), const char *varexpinput, const char *varexpoutput, int ns=100, Option_t *option="")
 
 TPMCABC (URANIE::DataServer::TDataServer *tds, const char *fcn, const char *varexpinput, const char *varexpoutput, int ns=100, Option_t *option="")
 
 TPMCABC (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, int ns=1, Option_t *option="")
 
virtual ~TPMCABC ()
 Default destructor.
 
PMC ABC function
void computeParameters (Option_t *option="")
 Generate the sample.
 
void initPMCABC ()
 
void generateNewSample ()
 
double pdfGaussianMultivariate (TMatrixD X, TMatrixD MU, TMatrixD Sigma)
 
int getNSteps ()
 
URANIE::DataServer::TDataServer * getTdsStep ()
 
void setThresholds (vector< double > vectTh)
 
void setNSteps (int steps)
 
- Public Member Functions inherited from URANIE::Calibration::TABC
 TABC (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns=100, Option_t *option="")
 Default ABC constructor with the runner argument: it contains the assessor to be used.
 
 TABC (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), const char *varexpinput, const char *varexpoutput, int ns=100, Option_t *option="")
 Default ABC constructor with the function argument: it contains the assessor to be used

 
 TABC (URANIE::DataServer::TDataServer *tds, const char *fcn, const char *varexpinput, const char *varexpoutput, int ns=100, Option_t *option="")
 Default ABC constructor with the function argument: it contains the assessor to be used.
 
 TABC (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, int ns=1, Option_t *option="")
 Default ABC constructor with the code argument: it contains the assessor to be used.
 
virtual ~TABC ()
 Default destructor.
 
void initABC ()
 
void setUpABC ()
 
void computeABC (int nEval)
 
void computeABCscore (Double_t newEpsilon)
 
void posteriorToPar ()
 
void checktdsParContent ()
 
void setGaussianNoise (const char *stdname)
 
void setPercentile (Double_t eps)
 Set the percentile.
 
URANIE::DataServer::TDataServer * getTdsPosterior ()
 Get the posterior tds.
 
double getThreshold ()
 Get the posterior tds

 
int getNEval ()
 
- 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.
 
virtual void setDistanceAndReference (const char *funcName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *reference, const char *weight="")
 Set the distance function and some needed informations.
 
virtual void setDistanceAndReference (URANIE::Calibration::TDistanceFunction *distFunc, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *reference, const char *weight="")
 Set the distance function and some needed informations.
 
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.
 
virtual void drawParameters (TString sTitre, const char *variable="*", const char *select="1>0", Option_t *option="")
 
void drawResidues (TString sTitre, const char *variable="*", const char *select="1>0", Option_t *option="")
 
void setLog ()
 
void unsetLog ()
 
void changeLog ()
 
Bool_t getLog ()
 
virtual void printLog (Option_t *option="")
 
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::TDistanceFunctiongetDistanceFunction ()
 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

URANIE::DataServer::TDataServer * _tdsParStep
 Pointer of the intermediate TDS.
 
vector< double > _frequencies
 Vector of the frequencies for sampling.
 
vector< double > _indices
 Vector of possible indices for sampling.
 
vector< double > _vEpsilon
 Vector of the thresholds.
 
TMatrixD _Sigmat
 Sqrt covariance matrix.
 
int _kStep
 Iterator of steps.
 
int _nSteps
 Number of steps.
 
- Public Attributes inherited from URANIE::Calibration::TABC
vector< double > _values
 
URANIE::DataServer::TDataServer * _tdsPosterior
 TDS containing posterior members for calibration.
 
- 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
 

Additional Inherited Members

- Public Types inherited from URANIE::Calibration::TCalibration
enum  ELauncher { kCode , kFunction , kRun , kUnknown }
 
- 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.
 
virtual void parseOption (Option_t *option="")
 Read the possible options.
 
void setMethodName (const char *str)
 Set the Method name.
 
- Protected Attributes inherited from URANIE::Calibration::TABC
double _threshold
 
Int_t _kPosterior
 
Int_t _kEval
 
TRandom3 * _rand
 
Bool_t _boolEpsilon
 
Double_t _Epsilon
 
- 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

◆ TPMCABC() [1/4]

URANIE::Calibration::TPMCABC::TPMCABC ( URANIE::DataServer::TDataServer *  tds,
URANIE::Relauncher::TRun *  run,
int  ns = 1,
Option_t *  option = "" 
)

Constructor of PMC ABC where nABC is set to be the number of the produced a posteriori sample.

Parameters
tdsthe pointer of the TDataServer. It must contains attribuutes that inherit from TStochasticAttribute
runthe runner that would be used to produce evaluations
nABCsize of a posteriori sample [100 by default], by default also size of the sample evaluated at each bloc,

Referenced by ClassImp().

◆ TPMCABC() [2/4]

URANIE::Calibration::TPMCABC::TPMCABC ( URANIE::DataServer::TDataServer *  tds,
void(*)(Double_t *, Double_t *)  fcn,
const char *  varexpinput,
const char *  varexpoutput,
int  ns = 100,
Option_t *  option = "" 
)

Constructor of PMC ABC where nABC is set to be the number of the produced a posteriori sample.

Parameters
tdsthe pointer of the TDataServer. It must contains attribuutes that inherit from TStochasticAttribute
fcnthe function that would be used to produce evaluations
nABCsize of a posteriori sample [100 by default], by default also size of the sample evaluated at each bloc,

◆ TPMCABC() [3/4]

URANIE::Calibration::TPMCABC::TPMCABC ( URANIE::DataServer::TDataServer *  tds,
const char *  fcn,
const char *  varexpinput,
const char *  varexpoutput,
int  ns = 100,
Option_t *  option = "" 
)

Constructor of PMC ABC where nABC is set to be the number of the produced a posteriori sample.

Parameters
tdsthe pointer of the TDataServer. It must contains attribuutes that inherit from TStochasticAttribute
fcnthe function that would be used to produce evaluations
nABCsize of a posteriori sample [100 by default], by default also size of the sample evaluated at each bloc,

◆ TPMCABC() [4/4]

URANIE::Calibration::TPMCABC::TPMCABC ( URANIE::DataServer::TDataServer *  tds,
URANIE::Launcher::TCode *  fcode,
int  ns = 1,
Option_t *  option = "" 
)

Constructor of PMC ABC where nABC is set to be the number of the produced a posteriori sample.

Parameters
tdsthe pointer of the TDataServer. It must contains attribuutes that inherit from TStochasticAttribute
codethe code that would be used to produce evaluations
nABCsize of a posteriori sample [100 by default], by default also size of the sample evaluated at each bloc,

◆ ~TPMCABC()

virtual URANIE::Calibration::TPMCABC::~TPMCABC ( )
virtual

Default destructor.

Referenced by ClassImp().

Member Function Documentation

◆ computeParameters()

void URANIE::Calibration::TPMCABC::computeParameters ( Option_t *  option = "")
virtual

Generate the sample.

Parameters
optionpossble 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().

◆ generateNewSample()

void URANIE::Calibration::TPMCABC::generateNewSample ( )

Referenced by ClassImp().

◆ getNSteps()

int URANIE::Calibration::TPMCABC::getNSteps ( )
inline

References _nSteps.

◆ getTdsStep()

URANIE::DataServer::TDataServer * URANIE::Calibration::TPMCABC::getTdsStep ( )
inline

References _tdsParStep.

◆ initPMCABC()

void URANIE::Calibration::TPMCABC::initPMCABC ( )

Referenced by ClassImp().

◆ pdfGaussianMultivariate()

double URANIE::Calibration::TPMCABC::pdfGaussianMultivariate ( TMatrixD  X,
TMatrixD  MU,
TMatrixD  Sigma 
)

Referenced by ClassImp().

◆ setNSteps()

void URANIE::Calibration::TPMCABC::setNSteps ( int  steps)
inline

References _nSteps.

◆ setThresholds()

void URANIE::Calibration::TPMCABC::setThresholds ( vector< double >  vectTh)
inline

References _nSteps, and _vEpsilon.

Member Data Documentation

◆ _frequencies

vector<double> URANIE::Calibration::TPMCABC::_frequencies

Vector of the frequencies for sampling.

Referenced by ClassImp().

◆ _indices

vector<double> URANIE::Calibration::TPMCABC::_indices

Vector of possible indices for sampling.

Referenced by ClassImp().

◆ _kStep

int URANIE::Calibration::TPMCABC::_kStep

Iterator of steps.

Referenced by ClassImp().

◆ _nSteps

int URANIE::Calibration::TPMCABC::_nSteps

Number of steps.

Referenced by ClassImp(), getNSteps(), setNSteps(), and setThresholds().

◆ _Sigmat

TMatrixD URANIE::Calibration::TPMCABC::_Sigmat

Sqrt covariance matrix.

Referenced by ClassImp().

◆ _tdsParStep

URANIE::DataServer::TDataServer* URANIE::Calibration::TPMCABC::_tdsParStep

Pointer of the intermediate TDS.

Referenced by ClassImp(), and getTdsStep().

◆ _vEpsilon

vector<double> URANIE::Calibration::TPMCABC::_vEpsilon

Vector of the thresholds.

Referenced by ClassImp(), and setThresholds().