1 #ifndef ALIFASTGLAUBER_H
2 #define ALIFASTGLAUBER_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
8 // Utility class to make simple Glauber type calculations for collision geometries:
9 // Impact parameter, production points, reaction plane dependence
10 // Author: Andreas Morsch
11 // andreas.morsch@cern.ch
19 class AliFastGlauber : public TObject {
21 static AliFastGlauber* Instance();
22 virtual ~AliFastGlauber();
23 void Init(Int_t mode = 0);
25 void SetWoodSaxonParameters(Double_t r0, Double_t d, Double_t w, Double_t n)
26 {fWSr0 = r0; fWSd = d; fWSw = w; fWSn = n;}
27 void SetWoodSaxonParametersAu()
28 {fWSr0 = 6.38; fWSd = 0.535; fWSw = 0.; fWSn = 8.59e-4;}
29 void SetWoodSaxonParametersPb()
30 {fWSr0 = 6.78; fWSd = 0.54; fWSw = 0.; fWSn = 7.14e-4;}
31 void SetMaxImpact(Float_t bmax = 20.) {fgBMax = bmax;};
32 void SetHardCrossSection(Float_t xs = 1.0) {fSigmaHard = xs;}
33 void SetNNCrossSection (Float_t xs = 55.6) {fSigmaNN = xs;}
34 void SetNucleus(Int_t n=208) {fA=n;}
37 void SetFileName(const TString &fn){fName=fn;}
38 void SetFileName(const char *fn="$(ALICE_ROOT)/FASTSIM/data/glauberPbPb.root"){fName=fn;}
40 const TF1* GetWSB() const {return fgWSb;}
41 const TF1* GetRWSB() const {return fgRWSb;}
42 const TF2* GetWSbz() const {return fgWSbz;}
43 const TF1* GetWSz() const {return fgWSz;}
44 const TF1* GetWSta() const {return fgWSta;}
45 const TF2* Kernel() const {return fgWStarfi;}
46 const TF2* GetWStarfi() const {return fgWStarfi;}
47 const TF2* GetWKParticipants() const {return fgWKParticipants;}
48 const TF1* GetWParticipants() const {return fgWParticipants;}
49 const TF1* Overlap() const {return fgWStaa;}
50 const TF1* GetWStaa() const {return fgWStaa;}
51 const TF2* GetWAlmond() const {return fgWAlmond;}
52 const TF1* GetWPathLength0() const {return fgWPathLength0;}
53 const TF1* GetWPathLength() const {return fgWPathLength;}
54 const TF1* GetWIntRadius() const {return fgWIntRadius;}
55 const TF1* GetWSgeo() const {return fgWSgeo;}
56 const TF1* GetWSbinary() const {return fgWSbinary;}
57 const TF1* GetWSN() const {return fgWSN;}
58 const TF1* GetWEnergyDensity() const {return fgWEnergyDensity;}
59 const TF2* GetWAlmondFixedB(Int_t i) const {return fgWAlmondFixedB[i];}
61 Float_t GetWr0() const {return fWSr0;}
62 Float_t GetWSd() const {return fWSd;}
63 Float_t GetWSw() const {return fWSw;}
64 Float_t GetWSn() const {return fWSn;}
65 Float_t GetSigmaHard() const {return fSigmaHard;}
66 Float_t GetSigmaNN() const {return fSigmaNN;}
67 Int_t GetA() const {return fA;}
68 const TString* GetFileName() const {return &fName;}
69 Float_t GetBmin() const {return fBmin;}
70 Float_t GetBmax() const {return fBmax;}
73 void DrawThickness() const;
74 void DrawOverlap() const;
75 void DrawParticipants() const;
77 void DrawBinary() const;
79 void DrawKernel(Double_t b = 0.) const;
80 void DrawAlmond(Double_t b = 0.) const;
81 void DrawPathLength0(Double_t b = 0., Int_t iopt = 0) const;
82 void DrawPathLength(Double_t b, Int_t ni = 1000, Int_t iopt = 0) const;
83 void DrawIntRadius(Double_t b = 0.) const;
84 void DrawEnergyDensity() const;
86 Double_t CrossSection(Double_t b1, Double_t b2) const;
87 Double_t HardCrossSection(Double_t b1, Double_t b2) const;
88 Double_t NHard(Double_t b1, Double_t b2) const;
89 Double_t FractionOfHardCrossSection(Double_t b1, Double_t b2) const;
90 Double_t Binaries(Double_t b) const;
91 Double_t GetNumberOfBinaries(Double_t b) const;
92 Double_t Participants(Double_t b) const;
93 Double_t GetNumberOfParticipants(Double_t b) const;
94 Double_t GetNumberOfCollisions(Double_t b) const;
95 Double_t GetNumberOfCollisionsPerEvent(Double_t b) const;
96 Double_t MeanOverlap(Double_t b1, Double_t b2);
97 Double_t MeanNumberOfCollisionsPerEvent(Double_t b1, Double_t b2);
98 void SimulateTrigger(Int_t n);
99 void GetRandom(Float_t& b, Float_t& p, Float_t& mult);
100 void GetRandom(Int_t& bin, Bool_t& hard);
101 Double_t GetRandomImpactParameter(Double_t bmin, Double_t bmax);
103 void StoreFunctions() const;
104 void StoreAlmonds() const;
106 void SetLengthDefinition(Int_t def=1) {fEllDef=def;}
107 Int_t GetLengthDef() const {return fEllDef;}
108 void SetCentralityClass(Double_t xsecFrLow=0.0,Double_t xsecFrUp=0.1);
109 void GetRandomBHard(Double_t& b);
110 void GetRandomXY(Double_t& x,Double_t& y);
111 void GetSavedXY(Double_t xy[2]) const {xy[0] = fXY[0]; xy[1] = fXY[1];}
112 void GetSavedI0I1(Double_t i0i1[2]) const {i0i1[0] = fI0I1[0]; i0i1[1] = fI0I1[1];}
113 void SaveXY(Double_t x, Double_t y) {fXY[0] = x; fXY[1] = y;}
114 void SaveI0I1(Double_t i0, Double_t i1) {fI0I1[0] = i0; fI0I1[1] = i1;}
116 void GetRandomPhi(Double_t& phi);
117 Double_t CalculateLength(Double_t b=0.,Double_t x0=0.,Double_t y0=0.,
119 void GetLengthAndPhi(Double_t& ell,Double_t &phi,Double_t b=-1.);
120 void GetLength(Double_t& ell,Double_t b=-1.);
121 void GetLengthsBackToBackAndPhi(Double_t& ell1,Double_t& ell2,
124 void GetLengthsBackToBack(Double_t& ell1,Double_t& ell2,
126 void GetLengthsForPythia(Int_t n,Double_t* const phi,Double_t* ell,
128 void PlotBDistr(Int_t n=1000);
129 void PlotLengthDistr(Int_t n=1000,Bool_t save=kFALSE,
130 const char *fname="length.root");
131 void PlotLengthB2BDistr(Int_t n=1000,Bool_t save=kFALSE,
132 const char *fname="lengthB2B.root");
133 void CalculateI0I1(Double_t& integral0,Double_t& integral1,
135 Double_t x0=0.,Double_t y0=0.,Double_t phi0=0.,
136 Double_t ellCut=20.) const;
137 void GetI0I1AndPhi(Double_t& integral0,Double_t& integral1,Double_t &phi,
138 Double_t ellCut=20.,Double_t b=-1.);
139 void GetI0I1(Double_t& integral0,Double_t& integral1,
140 Double_t ellCut=20.,Double_t b=-1.);
141 void GetI0I1BackToBackAndPhi(Double_t& integral01,Double_t& integral11,
142 Double_t& integral02,Double_t& integral12,
144 Double_t ellCut=20.,Double_t b=-1.);
145 void GetI0I1BackToBackAndPhiAndXY(Double_t& integral01,Double_t& integral11,
146 Double_t& integral02,Double_t& integral12,
147 Double_t& phi,Double_t& x,Double_t&y,
148 Double_t ellCut=20.,Double_t b=-1.);
149 void GetI0I1BackToBack(Double_t& integral01,Double_t& integral11,
150 Double_t& integral02,Double_t& integral12,
151 Double_t ellCut=20.,Double_t b=-1.);
152 void GetI0I1ForPythia(Int_t n,Double_t* phi,
153 Double_t* integral0,Double_t* integral1,
154 Double_t ellCut=20.,Double_t b=-1.);
155 void GetI0I1ForPythiaAndXY(Int_t n,Double_t* phi,
156 Double_t* integral0,Double_t* integral1,
157 Double_t&x, Double_t &y,
158 Double_t ellCut=20.,Double_t b=-1.);
159 void PlotI0I1Distr(Int_t n=1000,Double_t ellCut=20.,Bool_t save=kFALSE,
160 const char *fname="i0i1.root");
161 void PlotI0I1B2BDistr(Int_t n=1000,Double_t ellCut=20.,Bool_t save=kFALSE,
162 const char *fname="i0i1B2B.root");
163 void PlotAlmonds() const;
165 AliFastGlauber& operator=(const AliFastGlauber & rhs);
166 void Copy(TObject&) const;
168 static Double_t RWSb (const Double_t *xx, const Double_t *par);
169 static Double_t WSb (const Double_t *xx, const Double_t *par);
170 static Double_t WSbz (const Double_t *xx, const Double_t *par);
171 static Double_t WSz (const Double_t *xx, const Double_t *par);
172 static Double_t WSta (const Double_t *xx, const Double_t *par);
173 static Double_t WStarfi (const Double_t *xx, const Double_t *par);
174 static Double_t WStaa (const Double_t *xx, const Double_t *par);
175 static Double_t WKParticipants (const Double_t *xx, const Double_t *par);
176 static Double_t WParticipants (const Double_t *xx, const Double_t *par);
177 static Double_t WSgeo (const Double_t *xx, const Double_t *par);
178 static Double_t WSbinary (const Double_t *xx, const Double_t *par);
179 static Double_t WSN (const Double_t *xx, const Double_t *par);
180 static Double_t WAlmond (const Double_t *xx, const Double_t *par);
181 static Double_t WPathLength0 (const Double_t *xx, const Double_t *par);
182 static Double_t WPathLength (const Double_t *xx, const Double_t *par);
183 static Double_t WIntRadius (const Double_t *xx, const Double_t *par);
184 static Double_t WEnergyDensity (const Double_t *xx, const Double_t *par);
189 AliFastGlauber(const AliFastGlauber& glauber);
191 static Float_t fgBMax; // Maximum Impact Parameter
192 static const Int_t fgkMCInts; // Number of MC integrations
193 static AliFastGlauber* fgGlauber; // Singleton instance
196 static TF1* fgWSb; // Wood-Saxon Function (b)
197 static TF1* fgRWSb; // Wood-Saxon Function (b) with phase space factor
198 static TF2* fgWSbz; // Wood-Saxon Function (b, z)
199 static TF1* fgWSz; // Wood-Saxon Function (b = b0, z)
200 static TF1* fgWSta; // Thickness Function
201 static TF2* fgWStarfi; // Kernel for Overlap Function
202 static TF2* fgWKParticipants; // Kernel for number of participants
203 static TF1* fgWParticipants; // Number of participants
204 static TF1* fgWStaa; // Overlap Function
205 static TF2* fgWAlmond; // Interaction Almond
206 static TF1* fgWPathLength0; // Path Length as a function of phi
207 static TF1* fgWPathLength; // Path Length as a function of phi
208 static TF1* fgWIntRadius; // Interaction Radius
209 static TF1* fgWSgeo; // dSigma/db geometric
210 static TF1* fgWSbinary; // dSigma/db binary
211 static TF1* fgWSN; // dN/db binary
212 static TF1* fgWEnergyDensity; // Energy density as a function of impact parameter
213 static TF2* fgWAlmondFixedB[40]; // Interaction Almonds read from file
214 static TF2* fgWAlmondCurrent; // Interaction Almond used for length
216 Float_t fWSr0; // Wood-Saxon Parameter r0
217 Float_t fWSd; // Wood-Saxon Parameter d
218 Float_t fWSw; // Wood-Saxon Parameter w
219 Float_t fWSn; // Wood-Saxon Parameter n
220 // (chosen such that integral is one)
221 Float_t fSigmaHard; // Hard Cross Section [mbarn]
222 Float_t fSigmaNN; // NN Cross Section [mbarn]
223 Int_t fA; // Nucleon number of nucleus A
225 Float_t fBmin; // Minimum b (set through centrality selection)
226 Float_t fBmax; // Coresponding maximum b
227 Double_t fXY[2]; // Current generated production point
228 Double_t fI0I1[2]; // Current integrals I0 and I1
229 Int_t fEllDef; // definition of length (see CalculateLength())
230 TString fName; // filename of stored distributions
231 ClassDef(AliFastGlauber,2) // Event geometry simulation in the Glauber Model