Documentation / Developer's manual
Available modules
Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,  
Uranie / Sampler v4.9.0
|
#include <TMCMC.h>
Public Member Functions | |
TMCMC () | |
virtual | ~TMCMC () |
void | setSize (Int_t n) |
Int_t | getSize () |
TNtupleD * | varTimeMetropolisSampling (TNormale *f, TVectorD X) |
Method of simulation of a gaussian vector using the algorithm variable-at-a-time of MH. The components of the gaussian vector are updated one by one in each step, the instrumental law being a real gaussian centered at the value of the component to be updated at the previous step, and the variance is the same as the one of the law to be simulated. | |
TNtupleD * | randomWalkMetropolisSampling (TNormale *f, TVectorD X) |
Method of simulation of a gaussian vector using the algorithm random walk of MH. The components are updated at the same time at each step, and the instrumental law is a gaussian vector with mean being the one of the vector simulated at the previous step and with variance being the one of the law to be simulated multiplied by , where is the dimension of the law to be simulated. | |
TNtupleD * | gibbsSampling (TNormale *f, TVectorD X) |
Method implementing the algorithm of Gibbs, used to simulate a gaussian vector. | |
TVectorD | getObs (TNormale *f, TVectorD X) |
TNtupleD * | gibbsSampling (TMelange *mel, TVectorD X) |
Method implementing the algorithm of Gibbs, in order to simulate a multidimensionnal gaussian mixing. | |
TNtupleD * | randomWalkMetropolisSampling (TMelange *mel, TVectorD X) |
Method implementing the random Walk de MH of a a multidimensionnal gaussian mixing. The instrumental law is the same as observed in the case of a gaussian vector, except that the variance used is the mean value of the variances of each element of the mixing, multiplied by , where is the dimension of the law to be simulated. | |
TNtupleD * | varTimeMetropolisSampling (TMelange *mel, TVectorD X) |
Method implementing the "variable-at-a-time" of MH for a gaussian mixing. The transformation of the instrumentale law. | |
void | Draw1 (TMelange *mel, TNtupleD *ntd, Option_t *option) |
void | Draw (TMelange *mel, TNtupleD *ntd, Option_t *option) |
Method allowing to visualize the results of the MCMC methods programmed in this class, when the law has to be simulated is 2D dimension. Four visualizations are provided, allowing to follow the evolution of the simulation versus time, to get to know whether the values generated are distributed according to an expected law and after a certain number of iterations. | |
Double_t | getScore () |
Int_t | selectComponentState () |
Method used in the NKC, to determine the vector which will be used as mean of the instrumental law. | |
TNtupleD * | NKCSampling (TMelange *mel, TVectorD *X) |
Method allowing to deal with the "Normal Kernel Coupler", a MMC methode MCMC efficient for the simulation of mixing. | |
TNtupleD * | NKCSampling2 (TMelange *mel, TVectorD *X) |
Method implementing the NKC using 100 state vectors. | |
void | SetCanvas (TCanvas *canvas) |
void | SetTuple (TNtupleD *tuple) |
void | setVarDraw (TString str) |
Public Attributes | |
Int_t | _size |
Sample size. | |
TNormale * | q |
Variable representing a normal law which is used in all intermediary calculations. | |
Double_t | _score |
TRandom3 * | _rdm |
Generator of random number. | |
TCanvas * | cvisu |
TNtupleD * | _dataorg |
TString | _sdraw |
Generator of random number. | |
Constructor & Destructor Documentation
◆ TMCMC()
TMCMC::TMCMC | ( | ) |
◆ ~TMCMC()
|
virtual |
Member Function Documentation
◆ Draw()
void TMCMC::Draw | ( | TMelange * | mel, |
TNtupleD * | ntd, | ||
Option_t * | option | ||
) |
Method allowing to visualize the results of the MCMC methods programmed in this class, when the law has to be simulated is 2D dimension. Four visualizations are provided, allowing to follow the evolution of the simulation versus time, to get to know whether the values generated are distributed according to an expected law and after a certain number of iterations.
References _dataorg, _sdraw, cvisu, TMelange::getDimMelange(), TMelange::getMatEcartTypeCorrelation(), TMelange::getMean(), TMelange::getSigma(), TMelange::getTailleMelange(), and m.
Referenced by gibbsSampling(), randomWalkMetropolisSampling(), and varTimeMetropolisSampling().
◆ Draw1()
void TMCMC::Draw1 | ( | TMelange * | mel, |
TNtupleD * | ntd, | ||
Option_t * | option | ||
) |
References cvisu.
Referenced by gibbsSampling().
◆ getObs()
TVectorD TMCMC::getObs | ( | TNormale * | f, |
TVectorD | X | ||
) |
◆ getScore()
|
inline |
References _score.
◆ getSize()
|
inline |
References _size.
◆ gibbsSampling() [1/2]
TNtupleD * TMCMC::gibbsSampling | ( | TMelange * | mel, |
TVectorD | X | ||
) |
Method implementing the algorithm of Gibbs, in order to simulate a multidimensionnal gaussian mixing.
Since the conditional densities of such a law can not be simulated, a completion step is retained : A variable is introduced which will enable to obtain conditional densities that can be simulated.
References _rdm, _size, TMelange::discrimination(), Draw(), Draw1(), TMelange::getDimMelange(), TMelange::getMatEcartTypeCorrelation(), TMelange::getMean(), TNormale::getObs(), TNormale::getObsMultiDim(), TMelange::getTailleMelange(), TMelange::getTuple(), TMelange::LoiDeZ(), q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ gibbsSampling() [2/2]
TNtupleD * TMCMC::gibbsSampling | ( | TNormale * | f, |
TVectorD | X | ||
) |
Method implementing the algorithm of Gibbs, used to simulate a gaussian vector.
This is a particular case of the method "variable-at-a-time" of MH, where the instrumental laws aimed at updating the components, are their laws of each of them, conditioned under the others.
fast convergence
References _size, TDistribution::getDim(), TDistribution::getMatCovariance(), TDistribution::getMean(), TNormale::getObs(), TDistribution::getTuple(), q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ NKCSampling()
TNtupleD * TMCMC::NKCSampling | ( | TMelange * | mel, |
TVectorD * | X | ||
) |
Method allowing to deal with the "Normal Kernel Coupler", a MMC methode MCMC efficient for the simulation of mixing.
References _size, TMelange::addLoi(), TMelange::eval(), TMelange::getDimMelange(), TMelange::getMatEcartTypeCorrelation(), TNormale::getObsMultiDim(), TMelange::getTuple(), m, q, TMelange::reset(), ResizeTo(), selectComponentState(), TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ NKCSampling2()
TNtupleD * TMCMC::NKCSampling2 | ( | TMelange * | mel, |
TVectorD * | X | ||
) |
Method implementing the NKC using 100 state vectors.
Theorically this method is supposed to work, but we have to face to implementation problems.
References _rdm, _size, TNormale::eval(), TMelange::getDimMelange(), TMelange::getMatEcartTypeCorrelation(), TMelange::getMean(), TNormale::getObsMultiDim(), TMelange::getProb(), TMelange::getTailleMelange(), TMelange::getTuple(), m, q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ randomWalkMetropolisSampling() [1/2]
TNtupleD * TMCMC::randomWalkMetropolisSampling | ( | TMelange * | mel, |
TVectorD | X | ||
) |
Method implementing the random Walk de MH of a a multidimensionnal gaussian mixing. The instrumental law is the same as observed in the case of a gaussian vector, except that the variance used is the mean value of the variances of each element of the mixing, multiplied by , where is the dimension of the law to be simulated.
References _size, TMelange::discrimination(), Draw(), TMelange::eval(), TMelange::getDimMelange(), TMelange::getMatEcartTypeCorrelation(), TNormale::getObsMultiDim(), TDistribution::getRandom(), TMelange::getTailleMelange(), TMelange::getTuple(), q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ randomWalkMetropolisSampling() [2/2]
TNtupleD * TMCMC::randomWalkMetropolisSampling | ( | TNormale * | f, |
TVectorD | X | ||
) |
Method of simulation of a gaussian vector using the algorithm random walk of MH.
The components are updated at the same time at each step, and the instrumental law is a gaussian vector with mean being the one of the vector simulated at the previous step and with variance being the one of the law to be simulated multiplied by , where is the dimension of the law to be simulated.
relatively fast convergence
References _size, TNormale::eval(), TDistribution::getDim(), TDistribution::getMatEcartTypeCorrelation(), TNormale::getObsMultiDim(), TDistribution::getRandom(), TDistribution::getTuple(), m, q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ selectComponentState()
Int_t TMCMC::selectComponentState | ( | ) |
Method used in the NKC, to determine the vector which will be used as mean of the instrumental law.
References _rdm.
Referenced by NKCSampling().
◆ SetCanvas()
void TMCMC::SetCanvas | ( | TCanvas * | canvas | ) |
◆ setSize()
|
inline |
References _size.
◆ SetTuple()
void TMCMC::SetTuple | ( | TNtupleD * | tuple | ) |
No descriptions
References _dataorg.
◆ setVarDraw()
|
inline |
References _sdraw.
◆ varTimeMetropolisSampling() [1/2]
TNtupleD * TMCMC::varTimeMetropolisSampling | ( | TMelange * | mel, |
TVectorD | X | ||
) |
Method implementing the "variable-at-a-time" of MH for a gaussian mixing. The transformation of the instrumentale law.
References _size, TMelange::discrimination(), Draw(), TMelange::eval(), TNormale::eval(), TMelange::getDimMelange(), TNormale::getObs(), TDistribution::getRandom(), TMelange::getSigma(), TMelange::getTailleMelange(), TMelange::getTuple(), q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
◆ varTimeMetropolisSampling() [2/2]
TNtupleD * TMCMC::varTimeMetropolisSampling | ( | TNormale * | f, |
TVectorD | X | ||
) |
Method of simulation of a gaussian vector using the algorithm variable-at-a-time of MH.
The components of the gaussian vector are updated one by one in each step, the instrumental law being a real gaussian centered at the value of the component to be updated at the previous step, and the variance is the same as the one of the law to be simulated.
slow convergence
References _size, TNormale::eval(), TDistribution::getDim(), TDistribution::getMean(), TNormale::getObs(), TDistribution::getRandom(), TDistribution::getSigma(), TDistribution::getTuple(), q, TDistribution::setMatEcartTypeCorrelation(), and TDistribution::setMean().
Member Data Documentation
◆ _dataorg
TNtupleD* TMCMC::_dataorg |
Tuple de base ayant servi � mod�liser
Referenced by Draw(), and SetTuple().
◆ _rdm
TRandom3* TMCMC::_rdm |
Generator of random number.
Referenced by gibbsSampling(), NKCSampling2(), and selectComponentState().
◆ _score
Double_t TMCMC::_score |
Referenced by getScore().
◆ _sdraw
TString TMCMC::_sdraw |
Generator of random number.
Referenced by Draw(), and setVarDraw().
◆ _size
Int_t TMCMC::_size |
Sample size.
Referenced by getSize(), gibbsSampling(), gibbsSampling(), NKCSampling(), NKCSampling2(), randomWalkMetropolisSampling(), randomWalkMetropolisSampling(), setSize(), varTimeMetropolisSampling(), and varTimeMetropolisSampling().
◆ cvisu
TCanvas* TMCMC::cvisu |
Association d'un canvas pour la visu
Referenced by Draw(), Draw1(), and SetCanvas().
◆ q
TNormale* TMCMC::q |
Variable representing a normal law which is used in all intermediary calculations.
Referenced by getObs(), gibbsSampling(), gibbsSampling(), NKCSampling(), NKCSampling2(), randomWalkMetropolisSampling(), randomWalkMetropolisSampling(), varTimeMetropolisSampling(), and varTimeMetropolisSampling().