English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Calibration: TMCMC.h Source File
Uranie / Calibration  v4.11.0
/* @license-end */
TMCMC.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/>.
17 
18 #ifndef TMCMC_H
19 #define TMCMC_H
20 
21 // Uranie
22 #include "TCalibration.h"
23 
24 #include <unordered_map>
25 
26 namespace URANIE
27 {
28 namespace Calibration
29 {
30 
32 {
33 
34 public:
35 
36  // Attributes
37  vector<double> _values;
38  vector<double> _vstd;
39  vector<double> _dBaseSaved;
40  vector<int> _acceptSaved;
41  vector<int> _rejectSaved;
42 
43  // Random generator
44  TRandom3 _randSaved;
45 
46  std::string _simuName;
47  int _burnin;
48  int _lag;
49  int _nbDump;
51  double _lowAccRange;
52  double _higAccRange;
53  std::string _algoMCMC;
55 
56  //---------------------------------------------
60 
67  TMCMC(URANIE::DataServer::TDataServer* tds, URANIE::Relauncher::TRun* run, int ns=100, Option_t* option = "");
68 
78  TMCMC(URANIE::DataServer::TDataServer* tds, void (*fcn)(Double_t*,Double_t*), const char* varexpinput, const char* varexpoutput, int ns = 100, Option_t* option = "");
79 
88  TMCMC(URANIE::DataServer::TDataServer* tds, const char* fcn, const char* varexpinput, const char* varexpoutput, int ns = 100, Option_t* option = "");
89 
96  TMCMC(URANIE::DataServer::TDataServer* tds, URANIE::Launcher::TCode* fcode, int ns = 100, Option_t* option = "");
97 
101  void checktdsParContent();
103 
104  //---------------------------------------------
111  virtual ~TMCMC();
113 
114  //---------------------------------------------
118 
125  void setAlgo(const char* algoMCMC);
126 
132  void setSeed(int ichain, int seed);
133 
138  void setBurnin( int burn );
139 
144  void setLag( int lag );
145 
150  void setNbDump( int nbDump );
151 
157  void setAcceptationRatioRange( double lower, double higher);
158 
168  void setLikelihood(const char *likelihoodName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="");
169 
178  void setLikelihood(URANIE::Calibration::TDistanceLikelihoodFunction *likelihoodFunc, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *output, const char *weight="");
179 
184  void setMultistart(int nb_multistart);
185 
191  void setStartingPoints(int ichain, vector<double> values);
192 
198  void setProposalStd(int ichain, vector<double> standDev);
200 
201  //---------------------------------------------
205 
211  void computeParameters(Option_t *option = "");
212 
217  void continueCalculation(int new_Ns);
219 
220  //---------------------------------------------
224 
234  void MH_multiD(TRandom3* rand,vector<double> dBase,vector<int> accept,vector<int> reject, int Nstart, int Nend);
235 
245  void MH_1D(TRandom3* rand,vector<double> dBase,vector<int> accept,vector<int> reject, int Nstart, int Nend);
247 
248  //---------------------------------------------
252 
257  void export_chain_MCMC(const char *fileName);
258 
263  void read_chain_MCMC(const char *fileName);
265 
266  //---------------------------------------------
270 
278  void drawTrace(TString sTitre, const char *variable="*", Option_t *option="");
279 
287  void draw2DTrace(TString sTitre, const char *variable, Option_t *option="");
288 
296  void drawAcceptationRatio(TString sTitre, const char *variable="*", Option_t *option="");
297 
306  void drawParameters(TString sTitre, const char *variable = "*", Option_t * option = "");
308 
309  //---------------------------------------------
313 
319  double diagAutoCorrelation(int lag, int nPoints, Double_t* Data);
320 
324  std::unordered_map<string, double> diagGelmanRubin();
326 
330  std::unordered_map<string, vector<int>> diagESS();
332 
333  //---------------------------------------------
337 
341  virtual void printLog(Option_t *option = "");
343 
344  //---------------------------------------------
348 
352  string getDefaultCut()
353  {
354  string cut = "";
355  if( _burnin != 0 )
356  cut += string(Form("%s > %d",_tdsPar->getIteratorName(), _burnin));
357  if( _lag != 1 )
358  cut += ((cut=="")? "": " && " ) + string(Form("%s %% %d",_tdsPar->getIteratorName(), _lag));
359 
360  if(cut!="")
361  cout<<"TMCMC drawing method apply default cut ["<<cut<<"]"<<endl;
362  return cut;
363  }
364 
368  void clearDefaultCut(){ setBurnin(0); setLag(1); setNbDump(1000);}
369 
371 
372 protected:
373  //---------------------------------------------
377 
381  void logPriorPdf(double &ret);
382 
390  void parseOption(Option_t *option = "");
392 
393 
394  ClassDef(URANIE::Calibration::TMCMC, ID_CALIBRATION)
395  //Classe
396 };
397 }
398 }
399 
400 
401 #endif
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.
void setNbDump(int nbDump)
Define nbDump.
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)
virtual ~TMCMC()
Default destructor.
void continueCalculation(int new_Ns)
Continue the MCMC computation for new_Ns iterations.
void setBurnin(int burn)
Define the burn-in size.
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
virtual void printLog(Option_t *option="")
Prints the log.
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
void parseOption(Option_t *option="")
Possible options:
int _lag
The lag value.
Definition: TMCMC.h:48
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
void setLag(int lag)
Define the lag.