Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
unuranFoamTest.C File Reference

Detailed Description

This program must be compiled and executed with Aclic as follows.

it is an extension of tutorials foam_kanwa.C to compare generation of a 2D distribution with unuran and Foam

#include "TH2.h"
#include "TF2.h"
#include "TSystem.h"
#include "TCanvas.h"
#include "TMath.h"
#include "TRandom3.h"
#include "TFoam.h"
#include "TFoamIntegrand.h"
#include "TStopwatch.h"
#include "TROOT.h"
#include "TUnuran.h"
#include <iostream>
//_____________________________________________________________________________
Double_t sqr(Double_t x){return x*x;};
//_____________________________________________________________________________
//_____________________________________________________________________________
// 2-dimensional distribution for Foam, normalized to one (within 1e-5)
Double_t GamSq= sqr(0.100e0);
Double_t Dist= 0;
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;
}// Camel2
class FoamFunction : public TFoamIntegrand {
public:
virtual ~FoamFunction() {}
double Density(int nDim, double * x) {
return Camel2(nDim,x);
}
};
cout<<"--- kanwa started ---"<<endl;
gSystem->Load("libFoam.so");
TH2D *hst_xy = new TH2D("foam_hst_xy" , "FOAM x-y plot", 50,0,1.0, 50,0,1.0);
Double_t MCvect[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->SetRho(new FoamFunction() ); // Set 2-dim distribution, included below
FoamX->SetPseRan(PseRan); // Set random number generator
//
// From now on FoamX is ready to generate events
// test first the time
w.Start();
FoamX->Initialize(); // Initialize simulator, may take time...
//int nshow=5000;
int nshow=nev;
for(long loop=0; loop<nev; loop++){
FoamX->MakeEvent(); // generate MC event
FoamX->GetMCvect( MCvect); // get generated vector (x,y)
//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
w.Stop();
double time = w.CpuTime()*1.E9/nev;
cout << "Time using FOAM \t\t " << " \t=\t " << time << "\tns/call" << endl;
//
hst_xy->Draw("lego2"); // final plot
//
FoamX->GetIntegMC( MCresult, MCerror); // get MC integral, should be one
cout << " MCresult= " << MCresult << " +- " << MCerror <<endl;
cout<<"--- kanwa ended ---"<<endl;
return 0;
}//kanwa
double UCamel2(double * x, double *) {
return Camel2(2,x);
}
int run_unuran(int nev, std::string method = "hitro") {
// use unuran
std::cout << "run unuran " << std::endl;
gSystem->Load("libUnuran.so");
TH2D *h1 = new TH2D("unr_hst_xy" , "UNURAN x-y plot", 50,0,1.0, 50,0,1.0);
hUnr= h1;
TF2 * f = new TF2("f",UCamel2,0,1,0,1,0);
TUnuran unr(&r,2); // 2 is debug level
// test first the time
w.Start();
// init unuran
bool ret = unr.Init(dist,method);
if (!ret) {
std::cerr << "Error initializing unuran with method " << unr.MethodName() << endl;
return -1;
}
double x[2];
for (int i = 0; i < nev; ++i) {
unr.SampleMulti(x);
h1->Fill(x[0],x[1]);
// if (method == "gibbs" && i < 100)
// std::cout << x[0] << " , " << x[1] << std::endl;
}
w.Stop();
double time = w.CpuTime()*1.E9/nev;
cout << "Time using Unuran " << unr.MethodName() << " \t=\t " << time << "\tns/call" << endl;
h1->Draw("lego2");
return 0;
}
// visualising generated distribution
TCanvas *cKanwa = new TCanvas("cKanwa","Canvas for plotting",600,1000);
cKanwa->Divide(1,2);
cKanwa->cd(1);
int n = 100000;
run_unuran(n,"hitro");
cKanwa->Update();
cKanwa->cd(2);
cKanwa->Update();
std::cout <<"\nChi2 Test Results (UNURAN-FOAM):\t";
// test chi2
hFoam->Chi2Test(hUnr,"UUP");
return 0;
}
#define f(i)
Definition RSha256.hxx:104
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
#define ClassDef(name, id)
Definition Rtypes.h:342
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
R__EXTERN TSystem * gSystem
Definition TSystem.h:561
The Canvas class.
Definition TCanvas.h:23
A 2-Dim function with parameters.
Definition TF2.h:29
Abstract class representing n-dimensional real positive integrand function.
TFoam is the main class of the multi-dimensional general purpose Monte Carlo event generator (integra...
Definition TFoam.h:21
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3315
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3037
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:351
Service class for 2-D histogram classes.
Definition TH2.h:30
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
Stopwatch class.
Definition TStopwatch.h:28
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition TSystem.cxx:1857
TUnuranMultiContDist class describing multi dimensional continuous distributions.
TUnuran class.
Definition TUnuran.h:79
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5
double dist(Rotation3D const &r1, Rotation3D const &r2)
constexpr Double_t Pi()
Definition TMath.h:37
Author
Lorenzo Moneta

Definition in file unuranFoamTest.C.