Logo ROOT   6.16/01
Reference Guide
foam_kanwa.C File Reference

Detailed Description

View in nbviewer Open in SWAN This program can be execute from the command line as folows:

root -l foam_kanwa.C
auto * l
Definition: textangle.C:4
#include "Riostream.h"
#include "TFoam.h"
#include "TCanvas.h"
#include "TH2.h"
#include "TMath.h"
#include "TFoamIntegrand.h"
#include "TRandom3.h"
//_____________________________________________________________________________
return x*x;
}
//_____________________________________________________________________________
Double_t Camel2(Int_t nDim, Double_t *Xarg){
// 2-dimensional distribution for Foam, normalized to one (within 1e-5)
Double_t x=Xarg[0];
Double_t y=Xarg[1];
Double_t GamSq= sqr(0.100e0);
Dist +=exp(-(sqr(x-1./3) +sqr(y-1./3))/GamSq)/GamSq/TMath::Pi();
Dist +=exp(-(sqr(x-2./3) +sqr(y-2./3))/GamSq)/GamSq/TMath::Pi();
return 0.5*Dist;
}
//_____________________________________________________________________________
Int_t foam_kanwa(){
cout<<"--- kanwa started ---"<<endl;
TH2D *hst_xy = new TH2D("hst_xy" , "x-y plot", 50,0,1.0, 50,0,1.0);
Double_t *MCvect =new Double_t[2]; // 2-dim vector generated in the MC run
//
TRandom *PseRan = new TRandom3(); // Create random number generator
PseRan->SetSeed(4357);
TFoam *FoamX = new TFoam("FoamX"); // Create Simulator
FoamX->SetkDim(2); // No. of dimensions, obligatory!
FoamX->SetnCells(500); // Optionally No. of cells, default=2000
FoamX->SetRhoInt(Camel2); // Set 2-dim distribution, included below
FoamX->SetPseRan(PseRan); // Set random number generator
FoamX->Initialize(); // Initialize simulator, may take time...
//
// visualising generated distribution
TCanvas *cKanwa = new TCanvas("cKanwa","Canvas for plotting",600,600);
cKanwa->cd();
// From now on FoamX is ready to generate events
int nshow=5000;
for(long loop=0; loop<100000; loop++){
FoamX->MakeEvent(); // generate MC event
FoamX->GetMCvect( MCvect); // get generated vector (x,y)
Double_t x=MCvect[0];
Double_t y=MCvect[1];
if(loop<10) cout<<"(x,y) = ( "<< x <<", "<< y <<" )"<<endl;
hst_xy->Fill(x,y);
// live plot
if(loop == nshow){
nshow += 5000;
hst_xy->Draw("lego2");
cKanwa->Update();
}
}// loop
//
hst_xy->Draw("lego2"); // final plot
cKanwa->Update();
//
Double_t MCresult, MCerror;
FoamX->GetIntegMC( MCresult, MCerror); // get MC integral, should be one
cout << " MCresult= " << MCresult << " +- " << MCerror <<endl;
cout<<"--- kanwa ended ---"<<endl;
return 0;
}
int Int_t
Definition: RtypesCore.h:41
double Double_t
Definition: RtypesCore.h:55
double exp(double)
The Canvas class.
Definition: TCanvas.h:31
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2286
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:693
Definition: TFoam.h:27
virtual void GetMCvect(Double_t *)
User may get generated MC point/vector with help of this method.
Definition: TFoam.cxx:1202
virtual void MakeEvent()
User subprogram.
Definition: TFoam.cxx:1152
virtual void Initialize()
Basic initialization of FOAM invoked by the user.
Definition: TFoam.cxx:357
virtual void GetIntegMC(Double_t &, Double_t &)
User subprogram.
Definition: TFoam.cxx:1237
virtual void SetnCells(Long_t nCells)
Definition: TFoam.h:125
virtual void SetRhoInt(Double_t(*fun)(Int_t, Double_t *))
User may use this method to set the distribution object as a global function pointer (and not as an i...
Definition: TFoam.cxx:1060
virtual void SetPseRan(TRandom *PseRan)
Definition: TFoam.h:121
virtual void SetkDim(Int_t kDim)
Definition: TFoam.h:124
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:291
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
Random number generator class based on M.
Definition: TRandom3.h:27
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:27
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:589
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
double Dist(void *xp, void *yp)
VecExpr< UnaryOp< Sqr< T >, VecExpr< A, T, D >, T >, T, D > sqr(const VecExpr< A, T, D > &rhs)
constexpr Double_t Pi()
Definition: TMath.h:38
Author
Stascek Jadach

Definition in file foam_kanwa.C.