Example macro to show unuran capabilities The results are compared with what is obtained using TRandom or TF1::GetRandom The macro is divided in 3 parts:
#include <iostream>
#include <cassert>
using std::cout;
using std::endl;
#define NGEN 1000000
int izone = 0;
void testStringAPI() {
TH1D *
h1 =
new TH1D(
"h1G",
"gaussian distribution from Unuran",100,-10,10);
TH1D * h2 =
new TH1D(
"h2G",
"gaussian distribution from TRandom",100,-10,10);
cout << "\nTest using UNURAN string API \n\n";
if (! unr.
Init(
"normal()",
"method=arou") ) {
cout << "Error initializing unuran" << endl;
return;
}
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using Unuran method " << unr.
MethodName() <<
"\t=\t " << w.
CpuTime() << endl;
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using TRandom::Gaus \t=\t " << w.
CpuTime() << endl;
}
double distr(
double *
x,
double *p) {
}
double cdf(
double *
x,
double *p) {
}
void testDistr1D() {
cout << "\nTest 1D Continous distributions\n\n";
TH1D *
h1 =
new TH1D(
"h1BW",
"Breit-Wigner distribution from Unuran",100,-10,10);
TH1D * h2 =
new TH1D(
"h2BW",
"Breit-Wigner distribution from GetRandom",100,-10,10);
TF1 *
f =
new TF1(
"distrFunc",distr,-10,10,2);
double par[2] = {1,0};
TF1 *
fc =
new TF1(
"cdfFunc",cdf,-10,10,2);
int logLevel = 2;
std::string method = "tdr";
if (!unr.
Init(dist,method) ) {
cout << "Error initializing unuran" << endl;
return;
}
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using Unuran method " << unr.
MethodName() <<
"\t=\t " << w.
CpuTime() << endl;
for (
int i = 0; i <
n; ++i) {
double x =
f->GetRandom();
}
cout <<
"Time using TF1::GetRandom() \t=\t " << w.
CpuTime() << endl;
std::cout << " chi2 test of UNURAN vs GetRandom generated histograms: " << std::endl;
}
double gaus3d(
double *
x,
double *p) {
double sigma_x = p[0];
double sigma_y = p[1];
double sigma_z = p[2];
double rho = p[2];
double u =
x[0] / sigma_x ;
double v =
x[1] / sigma_y ;
double w =
x[2] / sigma_z ;
double result = (1 / (2 *
TMath::Pi() * sigma_x * sigma_y * sigma_z *
sqrt(
c)))
*
exp (-(u * u - 2 * rho * u *
v +
v *
v + w*w) / (2 *
c));
return result;
}
void testDistrMultiDim() {
cout << "\nTest Multidimensional distributions\n\n";
TH3D *
h1 =
new TH3D(
"h13D",
"gaussian 3D distribution from Unuran",50,-10,10,50,-10,10,50,-10,10);
TH3D * h2 =
new TH3D(
"h23D",
"gaussian 3D distribution from GetRandom",50,-10,10,50,-10,10,50,-10,10);
TF3 *
f =
new TF3(
"g3d",gaus3d,-10,10,-10,10,-10,10,3);
double par[3] = {2,2,0.5};
std::string method = "hitro";
if ( ! unr.
Init(dist,method) ) {
cout << "Error initializing unuran" << endl;
return;
}
for (int i = 0; i < NGEN; ++i) {
}
cout <<
"Time using Unuran method " << unr.
MethodName() <<
"\t=\t\t " << w.
CpuTime() << endl;
int np = 200;
for (int i = 0; i < NGEN; ++i) {
f->GetRandom3(
x[0],
x[1],
x[2]);
}
cout <<
"Time using TF1::GetRandom \t\t=\t " << w.
CpuTime() << endl;
std::cout << " chi2 test of UNURAN vs GetRandom generated histograms: " << std::endl;
}
double poisson(
double *
x,
double * p) {
}
void testDiscDistr() {
cout << "\nTest Discrete distributions\n\n";
TH1D *
h1 =
new TH1D(
"h1PS",
"Unuran Poisson prob",20,0,20);
TH1D * h2 =
new TH1D(
"h2PS",
"Poisson dist from TRandom",20,0,20);
double mu = 5;
TF1 *
f =
new TF1(
"fps",poisson,1,0,1);
bool ret = unr.
Init(dist2,
"dari");
if (!ret) return;
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using Unuran method " << unr.
MethodName() <<
"\t=\t\t " << w.
CpuTime() << endl;
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using TRandom::Poisson " <<
"\t=\t\t " << w.
CpuTime() << endl;
std::cout << " chi2 test of UNURAN vs TRandom generated histograms: " << std::endl;
}
void testEmpDistr() {
cout << "\nTest Empirical distributions using smoothing\n\n";
const int Ndata = 1000;
for (int i = 0; i < Ndata; ++i) {
if (i < 0.5*Ndata )
else
}
TH1D * h0 =
new TH1D(
"h0Ref",
"Starting data",100,-10,10);
TH1D *
h1 =
new TH1D(
"h1Unr",
"Unuran unbin Generated data",100,-10,10);
TH1D * h1b =
new TH1D(
"h1bUnr",
"Unuran bin Generated data",100,-10,10);
TH1D * h2 =
new TH1D(
"h2GR",
"Data from TH1::GetRandom",100,-10,10);
if (!unr.
Init(dist))
return;
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using Unuran unbin " << unr.
MethodName() <<
"\t=\t\t " << w.
CpuTime() << endl;
if (!unr.
Init(binDist))
return;
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using Unuran bin " << unr.
MethodName() <<
"\t=\t\t " << w.
CpuTime() << endl;
for (
int i = 0; i <
n; ++i) {
}
cout <<
"Time using TH1::GetRandom " <<
"\t=\t\t " << w.
CpuTime() << endl;
}
void unuranDemo() {
c1 =
new TCanvas(
"c1_unuranMulti",
"Multidimensional distribution",10,10,1000,1000);
testStringAPI();
testDistr1D();
testDistrMultiDim();
testDiscDistr();
testEmpDistr();
}
R__EXTERN TRandom * gRandom
R__EXTERN TStyle * gStyle
static struct mg_connection * fc(struct mg_context *ctx)
virtual void SetLineColor(Color_t lcolor)
Set the line color.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
A 3-Dim function with parameters.
1-D histogram with a double per channel (see TH1 documentation)}
virtual Double_t Chi2Test(const TH1 *h2, Option_t *option="UU", Double_t *res=0) const
test for comparing weighted and unweighted histograms
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
virtual void Draw(Option_t *option="")
Draw this histogram with options.
virtual Double_t GetRandom(TRandom *rng=nullptr) const
Return a random number distributed according the histogram bin contents.
virtual void FillN(Int_t ntimes, const Double_t *x, const Double_t *w, Int_t stride=1)
Fill this histogram with an array x and weights w.
3-D histogram with a double per channel (see TH1 documentation)}
Int_t Fill(Double_t)
Invalid Fill method.
Random number generator class based on the maximally quidistributed combined Tausworthe generator by ...
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
virtual Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
void Stop()
Stop the stopwatch.
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
TUnuranContDist class describing one dimensional continuous distribution.
TUnuranDiscrDist class for one dimensional discrete distribution.
void SetMode(int mode)
set the mode of the distribution (location of maximum probability)
void SetProbSum(double sum)
set the value of the sum of the probabilities in the given domain
TUnuranEmpDist class for describing empiral distributions.
TUnuranMultiContDist class describing multi dimensional continuous distributions.
int SampleDiscr()
Sample discrete distributions User is responsible for having previously correctly initialized with TU...
const std::string & MethodName() const
used Unuran method
bool SampleMulti(double *x)
Sample multidimensional distributions User is responsible for having previously correctly initialized...
bool Init(const std::string &distr, const std::string &method)
initialize with Unuran string interface
double Sample()
Sample 1D distribution User is responsible for having previously correctly initialized with TUnuran::...
double breitwigner_pdf(double x, double gamma, double x0=0)
Probability density function of Breit-Wigner distribution, which is similar, just a different definit...
double poisson_pdf(unsigned int n, double mu)
Probability density function of the Poisson distribution.
double breitwigner_cdf(double x, double gamma, double x0=0)
Cumulative distribution function (lower tail) of the Breit_Wigner distribution and it is similar (jus...
double dist(Rotation3D const &r1, Rotation3D const &r2)