English Français

Documentation / Manuel développeur

Modules disponibles

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Sampler: TMCMC.h Source File
Uranie / Sampler v4.9.0
/* @license-end */
TMCMC.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
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#ifndef TMCMC_H
18#define TMCMC_H
19
20#ifndef ROOT_TNtupleD
21#include "TNtupleD.h"
22#endif
23
24#ifndef ROOT_TNtuple
25#include "TNtuple.h"
26#endif
27
28#ifndef ROOT_TFile
29#include "TFile.h"
30#endif
31
32#ifndef ROOT_TVectorD
33#include "TVectorD.h"
34#endif
35
36#ifndef ROOT_TMatrixD
37#include "TMatrixD.h"
38#endif
39
40#ifndef ROOT_TRandom3
41#include "TRandom3.h"
42#endif
43
44#ifndef ROOT_TMath
45#include "TMath.h"
46#endif
47
48#ifndef ROOT_TCanvas
49#include "TCanvas.h"
50#endif
51
52#ifndef ROOT_TString
53#include "TString.h"
54#endif
55
56#ifndef ROOT_TEllipse
57#include "TEllipse.h"
58#endif
59
60#ifndef ROOT_TPaveText
61#include "TPaveText.h"
62#endif
63
64#ifndef ROOT_TMelange
65#include "TMelange.h"
66#endif
67
68#include "TDistribution.h"
69#include <iostream>
70#include "Rtypes.h"
71#include "TNormale.h"
72#include "TSampler.h"
73
74class TList;
75class TMCMC
76{
77 // Associations
78 // Attributes
79public:
81 Int_t _size;
84 Double_t _score;
86 TRandom3* _rdm;
88 TCanvas* cvisu;
90 TNtupleD *_dataorg;
93 TString _sdraw;
94 // TNtupleD* _vect;
95 // TFile* _rootFile;
96 //TNtuple* _tuple;
97 // Operations
98
99public:
101 virtual ~TMCMC();
102 void setSize(Int_t n)
103 {
104 _size = n;
105 }
106 Int_t getSize()
107 {
108 return _size;
109 }
110 TNtupleD* varTimeMetropolisSampling(TNormale* f, TVectorD X);
111 TNtupleD* randomWalkMetropolisSampling(TNormale* f, TVectorD X);
112 TNtupleD* gibbsSampling(TNormale* f, TVectorD X);
113 TVectorD getObs(TNormale* f, TVectorD X);
114 TNtupleD* gibbsSampling(TMelange* mel, TVectorD X);
115 TNtupleD* randomWalkMetropolisSampling(TMelange* mel, TVectorD X);
116 TNtupleD* varTimeMetropolisSampling(TMelange* mel, TVectorD X);
117 void Draw1(TMelange* mel, TNtupleD* ntd, Option_t *option);
118 void Draw(TMelange* mel, TNtupleD* ntd, Option_t *option);
119 Double_t getScore()
120 {
121 return _score;
122 }
123 Int_t selectComponentState();
124 TNtupleD* NKCSampling(TMelange* mel, TVectorD* X);
125 TNtupleD* NKCSampling2(TMelange* mel, TVectorD* X);
127 void SetCanvas(TCanvas * canvas);
129 void SetTuple(TNtupleD *tuple);
130
131 void setVarDraw(TString str)
132 {
133 _sdraw = str;
134 }
135
136 ClassDef(TMCMC, ID_SAMPLER)
137 //Classe
138};
139
140#endif
Creation of the abstract class TDistribution.
Creation of the abstract class TSampler.
Definition TMCMC.h:76
TVectorD getObs(TNormale *f, TVectorD X)
Definition TMCMC.cxx:199
TNtupleD * _dataorg
Definition TMCMC.h:90
Double_t _score
Definition TMCMC.h:84
Int_t getSize()
Definition TMCMC.h:106
void Draw1(TMelange *mel, TNtupleD *ntd, Option_t *option)
Definition TMCMC.cxx:445
TRandom3 * _rdm
Generator of random number.
Definition TMCMC.h:86
TNtupleD * gibbsSampling(TNormale *f, TVectorD X)
Method implementing the algorithm of Gibbs, used to simulate a gaussian vector.
Definition TMCMC.cxx:164
Int_t _size
Sample size.
Definition TMCMC.h:81
void SetCanvas(TCanvas *canvas)
Definition TMCMC.cxx:893
TNtupleD * varTimeMetropolisSampling(TNormale *f, TVectorD X)
Method of simulation of a gaussian vector using the algorithm variable-at-a-time of MH....
Definition TMCMC.cxx:62
void setSize(Int_t n)
Definition TMCMC.h:102
TNormale * q
Variable representing a normal law which is used in all intermediary calculations.
Definition TMCMC.h:83
TNtupleD * randomWalkMetropolisSampling(TNormale *f, TVectorD X)
Method of simulation of a gaussian vector using the algorithm random walk of MH. The components are...
Definition TMCMC.cxx:119
Double_t getScore()
Definition TMCMC.h:119
void SetTuple(TNtupleD *tuple)
Definition TMCMC.cxx:898
TNtupleD * NKCSampling(TMelange *mel, TVectorD *X)
Method allowing to deal with the "Normal Kernel Coupler", a MMC methode MCMC efficient for the simula...
Definition TMCMC.cxx:707
TNtupleD * NKCSampling2(TMelange *mel, TVectorD *X)
Method implementing the NKC using 100 state vectors.
Definition TMCMC.cxx:811
Int_t selectComponentState()
Method used in the NKC, to determine the vector which will be used as mean of the instrumental law.
Definition TMCMC.cxx:685
void setVarDraw(TString str)
Definition TMCMC.h:131
TString _sdraw
Generator of random number.
Definition TMCMC.h:93
void Draw(TMelange *mel, TNtupleD *ntd, Option_t *option)
Method allowing to visualize the results of the MCMC methods programmed in this class,...
Definition TMCMC.cxx:473
TCanvas * cvisu
Definition TMCMC.h:88
virtual ~TMCMC()
Definition TMCMC.cxx:52
Definition TMelange.h:42
Definition TNormale.h:48