English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Calibration: TABC.h Source File
Uranie / Calibration v4.9.0
/* @license-end */
TABC.h
Go to the documentation of this file.
1
2// 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/>.
17// TABC
19// $Id$
20// $Author$
21// $Date$
22// $Revision 1.2 $
23// $State$
25
35#ifndef TABC_H
36#define TABC_H
37
39#include "TTree.h"
40#include "TMatrixD.h"
41#include "TRandom3.h"
42#include "TMath.h"
43
44// Uranie
45#include "TCalibration.h"
46#include "TDataServer.h"
47
48namespace URANIE
49{
50namespace Calibration
51{
52
54{
55public :
56 vector<double> _values;
57 URANIE::DataServer::TDataServer *_tdsPosterior;
58
59protected:
60
61 double _threshold;
63 Int_t _kEval;
64
65 TRandom3 *_rand;
66
68 Double_t _Epsilon;
69
70
72public :
73
74 //---------------------------------------------
78
85 TABC(URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, int ns = 100, Option_t *option = "");
86
95 TABC(URANIE::DataServer::TDataServer *tds, void (*fcn)(Double_t*,Double_t*), const char *varexpinput, const char *varexpoutput, int ns = 100, Option_t *option = "");
96
105 TABC(URANIE::DataServer::TDataServer *tds, const char *fcn, const char *varexpinput, const char *varexpoutput, int ns = 100, Option_t *option = "");
106
113 TABC(URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, int ns = 1, Option_t *option = "");
114
116 virtual ~TABC();
118
119 //---------------------------------------------
123 void initABC();
124
126
127 //---------------------------------------------
132 //---------------------------------------------
135 void setUpABC();
136// void evalABC();
137 void computeABC(int nEval);
138 void computeABCscore(Double_t newEpsilon);
139
140 void posteriorToPar();
141
146 void checktdsParContent();
147
149 //---------------------------------------------
156 void setGaussianNoise(const char *stdname);
157
159 void setPercentile(Double_t eps)
160 {
161 try
162 {
163 if ((eps>1.0) || (eps<0.0))
164 {
165 throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
166 Form( "TABC::setPercentile: the percentile must be set between 0 and 1: however it is fix at %4.2f ! ", eps));
167 }
168 }
169 catch (URANIE::Exceptions::UErrorExceptions& ue)
170 {
171 ue.printMessage();
172 throw ue;
173 }
174
175 _boolEpsilon = true;
176 _Epsilon = eps;
177 }
178
180 URANIE::DataServer::TDataServer *getTdsPosterior()
181 {
182 return _tdsPosterior;
183 }
186 {
187 try
188 {
189 if (_threshold==12345678.9)
190 {
191 throw URANIE::Exceptions::UErrorExceptions(__FILE__, __LINE__,
192 "TABC::getThreshold: the threshold is not yet computed");
193 }
194 }
195 catch (URANIE::Exceptions::UErrorExceptions& ue)
196 {
197 ue.printMessage();
198 throw ue;
199 }
200 return _threshold;
201 }
203 {
204 return _kEval;
205 }
206
207
208
209
211
212 ClassDef(URANIE::Calibration::TABC, ID_CALIBRATION)
213
214};
215
216} // End of namespace ABC
217} // End of namespace URANIE
218#endif
Interface of class URANIE::Calibration::TCalibration.
Definition TABC.h:54
void checktdsParContent()
Definition TABC.cxx:100
vector< double > _values
Definition TABC.h:56
void setPercentile(Double_t eps)
Set the percentile.
Definition TABC.h:159
Int_t _kEval
Definition TABC.h:63
Bool_t _boolEpsilon
Definition TABC.h:67
virtual ~TABC()
Default destructor.
Definition TABC.cxx:147
double getThreshold()
Get the posterior tds
Definition TABC.h:185
void setUpABC()
Definition TABC.cxx:181
TRandom3 * _rand
Definition TABC.h:65
void computeABC(int nEval)
Definition TABC.cxx:198
int getNEval()
Definition TABC.h:202
void computeABCscore(Double_t newEpsilon)
Definition TABC.cxx:258
void setGaussianNoise(const char *stdname)
Definition TABC.cxx:292
URANIE::DataServer::TDataServer * _tdsPosterior
TDS containing posterior members for calibration.
Definition TABC.h:57
void posteriorToPar()
Definition TABC.cxx:225
Double_t _Epsilon
Definition TABC.h:68
URANIE::DataServer::TDataServer * getTdsPosterior()
Get the posterior tds.
Definition TABC.h:180
void initABC()
Definition TABC.cxx:168
Int_t _kPosterior
Definition TABC.h:62
double _threshold
Definition TABC.h:61
Description of the class TCalibration.
Definition TCalibration.h:64
Definition TABC.cxx:46