1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
21 #include "AliFastGlauber.h"
30 ClassImp(AliFastGlauber)
32 TF1* AliFastGlauber::fWSb = NULL;
33 TF2* AliFastGlauber::fWSbz = NULL;
34 TF1* AliFastGlauber::fWSz = NULL;
35 TF1* AliFastGlauber::fWSta = NULL;
36 TF2* AliFastGlauber::fWStarfi = NULL;
37 TF1* AliFastGlauber::fWStaa = NULL;
38 TF1* AliFastGlauber::fWSgeo = NULL;
39 TF1* AliFastGlauber::fWSbinary = NULL;
40 TF1* AliFastGlauber::fWSN = NULL;
41 Float_t AliFastGlauber::fbMax = 0.;
43 AliFastGlauber::AliFastGlauber()
45 // Default Constructor
48 SetWoodSaxonParameters(6.624, 0.549, 0.00, 7.69e-4);
49 SetHardCrossSection();
53 void AliFastGlauber::Init(Int_t mode)
59 fWSb = new TF1("WSb", WSb, 0, fbMax, 4);
60 fWSb->SetParameter(0, fWSr0);
61 fWSb->SetParameter(1, fWSd);
62 fWSb->SetParameter(2, fWSw);
63 fWSb->SetParameter(3, fWSn);
65 fWSbz = new TF2("WSbz", WSbz, 0, fbMax, 4);
66 fWSbz->SetParameter(0, fWSr0);
67 fWSbz->SetParameter(1, fWSd);
68 fWSbz->SetParameter(2, fWSw);
69 fWSbz->SetParameter(3, fWSn);
71 fWSz = new TF1("WSz", WSz, 0, fbMax, 5);
72 fWSz->SetParameter(0, fWSr0);
73 fWSz->SetParameter(1, fWSd);
74 fWSz->SetParameter(2, fWSw);
75 fWSz->SetParameter(3, fWSn);
80 fWSta = new TF1("WSta", WSta, 0., fbMax, 0);
85 fWStarfi = new TF2("WStarfi", WStarfi, 0., fbMax, 0., TMath::Pi(), 1);
86 fWStarfi->SetParameter(0, 0.);
87 fWStarfi->SetNpx(200);
93 fWStaa = new TF1("WStaa", WStaa, 0., fbMax, 0);
96 TFile* f = new TFile("$(ALICE_ROOT)/FASTSIM/data/glauberPbPb.root");
97 fWStaa = (TF1*) f->Get("WStaa");
101 // Geometrical Cross-Section
103 fWSgeo = new TF1("WSgeo", WSgeo, 0., fbMax, 0);
106 // Hard cross section (~ binary collisions)
108 fWSbinary = new TF1("WSbinary", WSbinary, 0., fbMax, 1);
109 fWSbinary->SetParameter(0, fSigmaHard); // mb
110 fWSbinary->SetNpx(100);
112 // Hard collisions per event
114 fWSN = new TF1("WSN", WSN, 0., fbMax, 1);
118 void AliFastGlauber::DrawWSb()
121 // Draw Wood-Saxon Nuclear Density Function
123 TCanvas *c1 = new TCanvas("c1","Wood Saxon",400,10,600,700);
128 void AliFastGlauber::DrawOverlap()
131 // Draw Overlap Function
133 TCanvas *c2 = new TCanvas("c2","Overlap",400,10,600,700);
138 void AliFastGlauber::DrawThickness()
141 // Draw Thickness Function
143 TCanvas *c3 = new TCanvas("c3","Thickness",400,10,600,700);
148 void AliFastGlauber::DrawGeo()
151 // Draw Geometrical Cross-Section
153 TCanvas *c3 = new TCanvas("c3","Geometrical Cross-Section",400,10,600,700);
158 void AliFastGlauber::DrawBinary()
161 // Draw Wood-Saxon Nuclear Density Function
163 TCanvas *c4 = new TCanvas("c4","Binary Cross-Section",400,10,600,700);
168 void AliFastGlauber::DrawN()
171 // Draw Binaries per event
173 TCanvas *c5 = new TCanvas("c5","Binaries per event",400,10,600,700);
178 void AliFastGlauber::DrawKernel()
183 TCanvas *c6 = new TCanvas("c6","Kernel",400,10,600,700);
188 Double_t AliFastGlauber::WSb(Double_t* x, Double_t* par)
191 // Woods-Saxon Parameterisation
192 // as a function of radius
195 Double_t r0 = par[0];
200 Double_t y = n * (1.+w*(xx/r0)*(xx/r0))/(1.+TMath::Exp((xx-r0)/d));
204 Double_t AliFastGlauber::WSbz(Double_t* x, Double_t* par)
207 // Wood Saxon Parameterisation
208 // as a function of z and b
212 Double_t r0 = par[0];
216 Double_t xx = TMath::Sqrt(bb*bb+zz*zz);
217 Double_t y = n * (1.+w*(xx/r0)*(xx/r0))/(1.+TMath::Exp((xx-r0)/d));
221 Double_t AliFastGlauber::WSz(Double_t* x, Double_t* par)
224 // Wood Saxon Parameterisation
225 // as a function of z for fixed b
227 Double_t bb = par[4];
229 Double_t r0 = par[0];
233 Double_t xx = TMath::Sqrt(bb*bb+zz*zz);
234 Double_t y = n * (1.+w*(xx/r0)*(xx/r0))/(1.+TMath::Exp((xx-r0)/d));
235 // if (y < 1.e-6) y = 0;
239 Double_t AliFastGlauber::WSta(Double_t* x, Double_t* par)
242 // Thickness function
245 fWSz->SetParameter(4, b);
246 Double_t y = 2. * fWSz->Integral(0., fbMax);
252 Double_t AliFastGlauber::WStarfi(Double_t* x, Double_t* par)
255 // Kernel for overlap function
260 Double_t r2 = TMath::Sqrt(r1 * r1 + b * b - 2. * r1 * b * TMath::Cos(phi));
261 Double_t y = r1 * fWSta->Eval(r1) * fWSta->Eval(r2);
266 Double_t AliFastGlauber::WStaa(Double_t* x, Double_t* par)
272 fWStarfi->SetParameter(0, b);
282 Double_t y = 2. * fWStarfi->IntegralMultiple(2, al, bl, 0.001, err);
283 printf("WStaa: %f %f %f\n", b, y, err);
289 for (Int_t i = 0; i < 100000; i++)
291 Double_t phi = TMath::Pi() * gRandom->Rndm();
292 Double_t b1 = fbMax * gRandom->Rndm();
293 y += fWStarfi->Eval(b1, phi);
295 y *= 2. * 0.1 * 208. * 208. * TMath::Pi() * fbMax / 100000.;
299 Double_t AliFastGlauber::WSgeo(Double_t* x, Double_t* par)
302 // Geometrical Cross-Section
305 Double_t taa = fWStaa->Eval(b);
306 const Double_t sigma = 55.6; // mbarn
308 Double_t y = 2. * TMath::Pi() * b * (1. - TMath::Exp(- sigma * taa)); // fm
313 Double_t AliFastGlauber::WSbinary(Double_t* x, Double_t* par)
316 // Geometrical Cross-Section
319 Double_t sigma = par[0];
320 Double_t taa = fWStaa->Eval(b);
322 Double_t y = 2. * TMath::Pi() * b * sigma * taa; // fm
326 Double_t AliFastGlauber::WSN(Double_t* x, Double_t* par)
329 // Geometrical Cross-Section
332 Double_t y = fWSbinary->Eval(b)/fWSgeo->Eval(b);
336 void AliFastGlauber::SimulateTrigger(Int_t n)
341 TH1F* mbtH = new TH1F("mbtH", "MB Trigger b-Distribution", 100, 0., 20.);
342 TH1F* hdtH = new TH1F("hdtH", "Hard Trigger b-Distribution", 100, 0., 20.);
343 TH1F* mbmH = new TH1F("mbmH", "MB Trigger Multiplicity Distribution", 100, 0., 8000.);
344 TH1F* hdmH = new TH1F("hdmH", "Hard Trigger Multiplicity Distribution", 100, 0., 8000.);
346 mbtH->SetXTitle("b [fm]");
347 hdtH->SetXTitle("b [fm]");
348 mbmH->SetXTitle("Multiplicity");
349 hdmH->SetXTitle("Multiplicity");
351 TCanvas *c0 = new TCanvas("c0","Trigger Simulation",400,10,600,700);
353 TCanvas *c1 = new TCanvas("c1","Trigger Simulation",400,10,600,700);
359 for (Int_t iev = 0; iev < n; iev++)
362 GetRandom(b, p, mult);
365 mbmH->Fill(mult, 1.);
382 void AliFastGlauber::GetRandom(Float_t& b, Float_t& p, Float_t& mult)
385 // Gives back a random impact parameter, hard trigger probability and multiplicity
387 b = fWSgeo->GetRandom();
388 Float_t mu = fWSN->Eval(b);
389 p = 1.-TMath::Exp(-mu);
390 mult = 6000./fWSN->Eval(1.) * mu;
393 Float_t AliFastGlauber::GetRandomImpactParameter(Float_t bmin, Float_t bmax)
396 // Gives back a random impact parameter in the range bmin .. bmax
400 while(b < bmin || b > bmax)
401 b = fWSgeo->GetRandom();
405 Double_t AliFastGlauber::CrossSection(Double_t b1, Double_t b2)
408 // Return cross-section integrated from b1 to b2
411 return fWSgeo->Integral(b1, b2)/100.;
414 Double_t AliFastGlauber::FractionOfHardCrossSection(Double_t b1, Double_t b2)
417 // Return raction of hard cross-section integrated from b1 to b2
420 return fWSbinary->Integral(b1, b2)/fWSbinary->Integral(0., 100.);
424 Double_t AliFastGlauber::Binaries(Double_t b)
427 // Return number of binary collisions normalized to 1 at b=0
430 return fWSN->Eval(b)/fWSN->Eval(1.);