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 **************************************************************************/
18 Revision 1.2 2003/04/14 14:23:44 morsch
19 Correction in Binaries().
21 Revision 1.1 2003/04/10 08:48:13 morsch
27 #include "AliFastGlauber.h"
36 ClassImp(AliFastGlauber)
38 TF1* AliFastGlauber::fWSb = NULL;
39 TF2* AliFastGlauber::fWSbz = NULL;
40 TF1* AliFastGlauber::fWSz = NULL;
41 TF1* AliFastGlauber::fWSta = NULL;
42 TF2* AliFastGlauber::fWStarfi = NULL;
43 TF1* AliFastGlauber::fWStaa = NULL;
44 TF1* AliFastGlauber::fWSgeo = NULL;
45 TF1* AliFastGlauber::fWSbinary = NULL;
46 TF1* AliFastGlauber::fWSN = NULL;
47 Float_t AliFastGlauber::fbMax = 0.;
49 AliFastGlauber::AliFastGlauber()
51 // Default Constructor
54 SetWoodSaxonParameters(6.624, 0.549, 0.00, 7.69e-4);
55 SetHardCrossSection();
59 void AliFastGlauber::Init(Int_t mode)
65 fWSb = new TF1("WSb", WSb, 0, fbMax, 4);
66 fWSb->SetParameter(0, fWSr0);
67 fWSb->SetParameter(1, fWSd);
68 fWSb->SetParameter(2, fWSw);
69 fWSb->SetParameter(3, fWSn);
71 fWSbz = new TF2("WSbz", WSbz, 0, fbMax, 4);
72 fWSbz->SetParameter(0, fWSr0);
73 fWSbz->SetParameter(1, fWSd);
74 fWSbz->SetParameter(2, fWSw);
75 fWSbz->SetParameter(3, fWSn);
77 fWSz = new TF1("WSz", WSz, 0, fbMax, 5);
78 fWSz->SetParameter(0, fWSr0);
79 fWSz->SetParameter(1, fWSd);
80 fWSz->SetParameter(2, fWSw);
81 fWSz->SetParameter(3, fWSn);
86 fWSta = new TF1("WSta", WSta, 0., fbMax, 0);
91 fWStarfi = new TF2("WStarfi", WStarfi, 0., fbMax, 0., TMath::Pi(), 1);
92 fWStarfi->SetParameter(0, 0.);
93 fWStarfi->SetNpx(200);
99 fWStaa = new TF1("WStaa", WStaa, 0., fbMax, 0);
102 TFile* f = new TFile("$(ALICE_ROOT)/FASTSIM/data/glauberPbPb.root");
103 fWStaa = (TF1*) f->Get("WStaa");
107 // Geometrical Cross-Section
109 fWSgeo = new TF1("WSgeo", WSgeo, 0., fbMax, 0);
112 // Hard cross section (~ binary collisions)
114 fWSbinary = new TF1("WSbinary", WSbinary, 0., fbMax, 1);
115 fWSbinary->SetParameter(0, fSigmaHard); // mb
116 fWSbinary->SetNpx(100);
118 // Hard collisions per event
120 fWSN = new TF1("WSN", WSN, 0., fbMax, 1);
124 void AliFastGlauber::DrawWSb()
127 // Draw Wood-Saxon Nuclear Density Function
129 TCanvas *c1 = new TCanvas("c1","Wood Saxon",400,10,600,700);
134 void AliFastGlauber::DrawOverlap()
137 // Draw Overlap Function
139 TCanvas *c2 = new TCanvas("c2","Overlap",400,10,600,700);
144 void AliFastGlauber::DrawThickness()
147 // Draw Thickness Function
149 TCanvas *c3 = new TCanvas("c3","Thickness",400,10,600,700);
154 void AliFastGlauber::DrawGeo()
157 // Draw Geometrical Cross-Section
159 TCanvas *c3 = new TCanvas("c3","Geometrical Cross-Section",400,10,600,700);
164 void AliFastGlauber::DrawBinary()
167 // Draw Wood-Saxon Nuclear Density Function
169 TCanvas *c4 = new TCanvas("c4","Binary Cross-Section",400,10,600,700);
174 void AliFastGlauber::DrawN()
177 // Draw Binaries per event
179 TCanvas *c5 = new TCanvas("c5","Binaries per event",400,10,600,700);
184 void AliFastGlauber::DrawKernel()
189 TCanvas *c6 = new TCanvas("c6","Kernel",400,10,600,700);
194 Double_t AliFastGlauber::WSb(Double_t* x, Double_t* par)
197 // Woods-Saxon Parameterisation
198 // as a function of radius
201 Double_t r0 = par[0];
206 Double_t y = n * (1.+w*(xx/r0)*(xx/r0))/(1.+TMath::Exp((xx-r0)/d));
210 Double_t AliFastGlauber::WSbz(Double_t* x, Double_t* par)
213 // Wood Saxon Parameterisation
214 // as a function of z and b
218 Double_t r0 = par[0];
222 Double_t xx = TMath::Sqrt(bb*bb+zz*zz);
223 Double_t y = n * (1.+w*(xx/r0)*(xx/r0))/(1.+TMath::Exp((xx-r0)/d));
227 Double_t AliFastGlauber::WSz(Double_t* x, Double_t* par)
230 // Wood Saxon Parameterisation
231 // as a function of z for fixed b
233 Double_t bb = par[4];
235 Double_t r0 = par[0];
239 Double_t xx = TMath::Sqrt(bb*bb+zz*zz);
240 Double_t y = n * (1.+w*(xx/r0)*(xx/r0))/(1.+TMath::Exp((xx-r0)/d));
241 // if (y < 1.e-6) y = 0;
245 Double_t AliFastGlauber::WSta(Double_t* x, Double_t* par)
248 // Thickness function
251 fWSz->SetParameter(4, b);
252 Double_t y = 2. * fWSz->Integral(0., fbMax);
258 Double_t AliFastGlauber::WStarfi(Double_t* x, Double_t* par)
261 // Kernel for overlap function
266 Double_t r2 = TMath::Sqrt(r1 * r1 + b * b - 2. * r1 * b * TMath::Cos(phi));
267 Double_t y = r1 * fWSta->Eval(r1) * fWSta->Eval(r2);
272 Double_t AliFastGlauber::WStaa(Double_t* x, Double_t* par)
278 fWStarfi->SetParameter(0, b);
288 Double_t y = 2. * fWStarfi->IntegralMultiple(2, al, bl, 0.001, err);
289 printf("WStaa: %f %f %f\n", b, y, err);
295 for (Int_t i = 0; i < 100000; i++)
297 Double_t phi = TMath::Pi() * gRandom->Rndm();
298 Double_t b1 = fbMax * gRandom->Rndm();
299 y += fWStarfi->Eval(b1, phi);
301 y *= 2. * 0.1 * 208. * 208. * TMath::Pi() * fbMax / 100000.;
305 Double_t AliFastGlauber::WSgeo(Double_t* x, Double_t* par)
308 // Geometrical Cross-Section
311 Double_t taa = fWStaa->Eval(b);
312 const Double_t sigma = 55.6; // mbarn
314 Double_t y = 2. * TMath::Pi() * b * (1. - TMath::Exp(- sigma * taa)); // fm
319 Double_t AliFastGlauber::WSbinary(Double_t* x, Double_t* par)
322 // Geometrical Cross-Section
325 Double_t sigma = par[0];
326 Double_t taa = fWStaa->Eval(b);
328 Double_t y = 2. * TMath::Pi() * b * sigma * taa; // fm
332 Double_t AliFastGlauber::WSN(Double_t* x, Double_t* par)
335 // Geometrical Cross-Section
340 y = fWSbinary->Eval(b)/fWSgeo->Eval(b);
342 y = fWSbinary->Eval(1.e-4)/fWSgeo->Eval(1.e-4);
347 void AliFastGlauber::SimulateTrigger(Int_t n)
352 TH1F* mbtH = new TH1F("mbtH", "MB Trigger b-Distribution", 100, 0., 20.);
353 TH1F* hdtH = new TH1F("hdtH", "Hard Trigger b-Distribution", 100, 0., 20.);
354 TH1F* mbmH = new TH1F("mbmH", "MB Trigger Multiplicity Distribution", 100, 0., 8000.);
355 TH1F* hdmH = new TH1F("hdmH", "Hard Trigger Multiplicity Distribution", 100, 0., 8000.);
357 mbtH->SetXTitle("b [fm]");
358 hdtH->SetXTitle("b [fm]");
359 mbmH->SetXTitle("Multiplicity");
360 hdmH->SetXTitle("Multiplicity");
362 TCanvas *c0 = new TCanvas("c0","Trigger Simulation",400,10,600,700);
364 TCanvas *c1 = new TCanvas("c1","Trigger Simulation",400,10,600,700);
370 for (Int_t iev = 0; iev < n; iev++)
373 GetRandom(b, p, mult);
376 mbmH->Fill(mult, 1.);
393 void AliFastGlauber::GetRandom(Float_t& b, Float_t& p, Float_t& mult)
396 // Gives back a random impact parameter, hard trigger probability and multiplicity
398 b = fWSgeo->GetRandom();
399 Float_t mu = fWSN->Eval(b);
400 p = 1.-TMath::Exp(-mu);
401 mult = 6000./fWSN->Eval(0.) * mu;
404 Float_t AliFastGlauber::GetRandomImpactParameter(Float_t bmin, Float_t bmax)
407 // Gives back a random impact parameter in the range bmin .. bmax
411 while(b < bmin || b > bmax)
412 b = fWSgeo->GetRandom();
416 Double_t AliFastGlauber::CrossSection(Double_t b1, Double_t b2)
419 // Return cross-section integrated from b1 to b2
422 return fWSgeo->Integral(b1, b2)/100.;
425 Double_t AliFastGlauber::FractionOfHardCrossSection(Double_t b1, Double_t b2)
428 // Return raction of hard cross-section integrated from b1 to b2
431 return fWSbinary->Integral(b1, b2)/fWSbinary->Integral(0., 100.);
435 Double_t AliFastGlauber::Binaries(Double_t b)
438 // Return number of binary collisions normalized to 1 at b=0
441 return fWSN->Eval(b)/fWSN->Eval(0.);