English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Calibration: TPMCABC.h Source File
Uranie / Calibration  v4.11.0
/* @license-end */
TPMCABC.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 by
6 // 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/>.
18 // $Id$
19 // $Author$
20 // $Date$
21 // $Revision 1.2 $
22 // $State$
24 
32 #ifndef __PMC_ABC__
33 #define __PMC_ABC__
34 
35 // Root
36 #include "TTree.h"
37 #include "TMatrixD.h"
38 
39 // Uranie
40 #include "TABC.h"
41 #include "TDataServer.h"
42 
43 namespace URANIE
44 {
45 namespace Calibration
46 {
47 
49 {
50 
51 public :
52  URANIE::DataServer::TDataServer *_tdsParStep;
53  vector<double> _frequencies;
54  vector<double> _indices;
55  vector<double> _vEpsilon;
56  TMatrixD _Sigmat;
57  int _kStep;
58  int _nSteps;
59 
60  //---------------------------------------------
64 
71  TPMCABC(URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns = 1, Option_t *option = "");
72 
79  TPMCABC(URANIE::DataServer::TDataServer *tds, void (*fcn)(Double_t*,Double_t*), const char *varexpinput, const char *varexpoutput, int ns = 100, Option_t *option = "");
80 
81 
88  TPMCABC(URANIE::DataServer::TDataServer *tds, const char *fcn, const char *varexpinput, const char *varexpoutput, int ns = 100, Option_t *option = "");
89 
96  TPMCABC(URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, int ns = 1, Option_t *option = "");
97 
98 
100  virtual ~TPMCABC();
102 
103  //---------------------------------------------
120  virtual void setDistance(const char *distName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *reference, const char *weight="");
121 
130  virtual void setDistance(URANIE::Calibration::TDistanceLikelihoodFunction *distFunc, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *reference, const char *weight="");
131 
133 
139  void computeParameters(Option_t *option = "");
141 
142  void initPMCABC();
143  void generateNewSample();
144  double pdfGaussianMultivariate(TMatrixD X,TMatrixD MU,TMatrixD Sigma);
145 
146  int getNSteps()
147  {
148  return _nSteps;
149  }
150 
151 
152  URANIE::DataServer::TDataServer *getTdsStep()
153  {
154  return _tdsParStep;
155  }
156 
157  void setThresholds(vector<double> vectTh)
158  {
159  _vEpsilon = vectTh;
160  try
161  {
162  for(int i=0; i < (int)_vEpsilon.size(); i++)
163  {
164  if ((_vEpsilon[i]>1.0) or (_vEpsilon[i]<0.0))
165  {
166  throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
167  Form( "TABC::setThresholds: the percentile must be set between 0 and 1: however it is fix at %4.2f at step %d ! ",_vEpsilon[i],i));
168  }
169  }
170  }
171  catch (URANIE::Exceptions::UErrorExceptions& ue)
172  {
173  ue.printMessage();
174  throw ue;
175  }
176  _nSteps = (int)_vEpsilon.size();
177  }
178  void setNSteps(int steps)
179  {
180  try
181  {
182  if ((steps<1))
183  {
184  throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
185  Form( "TABC::setNSteps: the number of steps must be strictly upper than 0: however it is fix at %d ! ", steps));
186  }
187  }
188  catch (URANIE::Exceptions::UErrorExceptions& ue)
189  {
190  ue.printMessage();
191  throw ue;
192  }
193  _nSteps = steps;
194  }
195 
196 
197  ClassDef(URANIE::Calibration::TPMCABC, ID_CALIBRATION)
198 
199 };
200 
201 } // Fin du namespace ABC
202 } // Fin du namespace URANIE
203 #endif
Definition: TABC.cxx:45
double pdfGaussianMultivariate(TMatrixD X, TMatrixD MU, TMatrixD Sigma)
int getNSteps()
Definition: TPMCABC.h:146
Description of the class TDistanceLikelihoodFunction.
Definition: TDistanceLikelihoodFunction.h:67
int _kStep
Iterator of steps.
Definition: TPMCABC.h:57
Definition: TABC.h:53
void setNSteps(int steps)
Definition: TPMCABC.h:178
void computeParameters(Option_t *option="")
Generate the sample.
URANIE::DataServer::TDataServer * _tdsParStep
Pointer of the intermediate TDS.
Definition: TPMCABC.h:52
Interface of class URANIE::ABC::TABC.
virtual void setDistance(const char *distName, URANIE::DataServer::TDataServer *tdsRef, const char *input, const char *reference, const char *weight="")
Set the distance function and some needed informations.
int _nSteps
Number of steps.
Definition: TPMCABC.h:58
virtual ~TPMCABC()
Default destructor.
vector< double > _indices
Vector of possible indices for sampling.
Definition: TPMCABC.h:54
URANIE::DataServer::TDataServer * getTdsStep()
Definition: TPMCABC.h:152
TPMCABC(URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns=1, Option_t *option="")
TMatrixD _Sigmat
Sqrt covariance matrix.
Definition: TPMCABC.h:56
vector< double > _frequencies
Vector of the frequencies for sampling.
Definition: TPMCABC.h:53
vector< double > _vEpsilon
Vector of the thresholds.
Definition: TPMCABC.h:55
Definition: TPMCABC.h:48
void setThresholds(vector< double > vectTh)
Definition: TPMCABC.h:157