Documentation / Developer's manual
Available modules
Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,  ![]() |
Uranie / Calibration
v4.11.0
|
TMCMC.h
Go to the documentation of this file.
67 TMCMC(URANIE::DataServer::TDataServer* tds, URANIE::Relauncher::TRun* run, int ns=100, Option_t* option = "");
78 TMCMC(URANIE::DataServer::TDataServer* tds, void (*fcn)(Double_t*,Double_t*), const char* varexpinput, const char* varexpoutput, int ns = 100, Option_t* option = "");
88 TMCMC(URANIE::DataServer::TDataServer* tds, const char* fcn, const char* varexpinput, const char* varexpoutput, int ns = 100, Option_t* option = "");
96 TMCMC(URANIE::DataServer::TDataServer* tds, URANIE::Launcher::TCode* fcode, int ns = 100, Option_t* option = "");
168 void setLikelihood(const char *likelihoodName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="");
178 void setLikelihood(URANIE::Calibration::TDistanceLikelihoodFunction *likelihoodFunc, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="");
234 void MH_multiD(TRandom3* rand,vector<double> dBase,vector<int> accept,vector<int> reject, int Nstart, int Nend);
245 void MH_1D(TRandom3* rand,vector<double> dBase,vector<int> accept,vector<int> reject, int Nstart, int Nend);
Definition: TABC.cxx:45
void setLikelihood(const char *likelihoodName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="")
Set the likelihood function and some needed informations.
std::string _simuName
name of the simulation (will be the name of the results repository)
Definition: TMCMC.h:46
void logPriorPdf(double &ret)
Compute the logarithm of the prior.
string getDefaultCut()
Print the actual cut and lag values.
Definition: TMCMC.h:352
void drawParameters(TString sTitre, const char *variable="*", Option_t *option="")
Draw samples of the posterior distribution (if converged)
Description of the class TDistanceLikelihoodFunction.
Definition: TDistanceLikelihoodFunction.h:67
URANIE::DataServer::TDataServer * _tdsPar
TDS containing parameters properties (parameters that should be calibrated)
Definition: TCalibration.h:79
void setAcceptationRatioRange(double lower, double higher)
Define the range of the acceptation ratio (the proposal will be modified so that the acceptation rati...
void export_chain_MCMC(const char *fileName)
Save the current state of the Markov Chain in a file.
double diagAutoCorrelation(int lag, int nPoints, Double_t *Data)
Compute the autocorrelation.
std::unordered_map< string, double > diagGelmanRubin()
Compute Gelman-Rubin statistic (multistart only)
void continueCalculation(int new_Ns)
Continue the MCMC computation for new_Ns iterations.
TMCMC(URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns=100, Option_t *option="")
void drawTrace(TString sTitre, const char *variable="*", Option_t *option="")
Draw the evolution of the parameters as a function of the iterator.
int _nbDump
Frequency to which the algo dump a line.
Definition: TMCMC.h:49
std::string _algoMCMC
MCMC algo to use between ("MH_1D" and "MH_multiD")
Definition: TMCMC.h:53
void setProposalStd(int ichain, vector< double > standDev)
Initialise the proposal std from vector.
void setStartingPoints(int ichain, vector< double > values)
Initialise the parameters values.
Definition: TMCMC.h:31
std::unordered_map< string, vector< int > > diagESS()
Compute Effective Sample Size.
void draw2DTrace(TString sTitre, const char *variable, Option_t *option="")
Draw the evolution of 2 parameters on a 2D plot.
void computeParameters(Option_t *option="")
Compute ns iterations (with ns fixed in the constructor TMCMC)
double _lowAccRange
lower acceptation ratio bound
Definition: TMCMC.h:51
double _higAccRange
higher acceptation ratio bound
Definition: TMCMC.h:52
void MH_multiD(TRandom3 *rand, vector< double > dBase, vector< int > accept, vector< int > reject, int Nstart, int Nend)
Run Metropolis-Hastings algorithm with candidates drawn in multiple directions (default method) ...
Interface of class URANIE::Calibration::TCalibration.
Description of the class TCalibration.
Definition: TCalibration.h:64
int _multiStart
The number of chains initialised.
Definition: TMCMC.h:50
vector< double > _values
Vector of values to be tested.
Definition: TMCMC.h:37
void clearDefaultCut()
Reinitialisation of the cut and lag values.
Definition: TMCMC.h:368
void MH_1D(TRandom3 *rand, vector< double > dBase, vector< int > accept, vector< int > reject, int Nstart, int Nend)
Run Metropolis-Hastings algorithm with candidates drawn one direction at a time.
void setAlgo(const char *algoMCMC)
Define the seed of the random generator for one or several chains.
vector< int > _rejectSaved
Vector containing the number of rejected candidates for each parameter.
Definition: TMCMC.h:41
void read_chain_MCMC(const char *fileName)
Read a saved Markov Chain.
void drawAcceptationRatio(TString sTitre, const char *variable="*", Option_t *option="")
Draw the evolution of the acceptation ratio as a function of the iterator.
int _burnin
The number of iterations considered as warm-up or burn-in.
Definition: TMCMC.h:47
bool _bcleaningAtt
Do not store the underlying att.
Definition: TMCMC.h:54
vector< int > _acceptSaved
Vector containing the number of accepted candidates for each parameter.
Definition: TMCMC.h:40
vector< double > _vstd
Vector of standard deviation.
Definition: TMCMC.h:38
TRandom3 _randSaved
Random generator saved.
Definition: TMCMC.h:44
void setMultistart(int nb_multistart)
Initialise several chains.
void setSeed(int ichain, int seed)
Define the seed of the random generator for one or several chains.
void checktdsParContent()
Initialise the TMCMC object and the strcuture of file to save the result.
vector< double > _dBaseSaved
Database containing the values, acceptation rate, ...
Definition: TMCMC.h:39
