]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FASTSIM/AliFastGlauber.h
32c4a550a0a46972a3f05bb0259319438e69ce58
[u/mrichter/AliRoot.git] / FASTSIM / AliFastGlauber.h
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                               */
5
6 /* $Id$ */
7 //
8 // Utility class to make simple Glauber type calculations for collision geometries:
9 // Impact parameter, production points, reaction plane dependence
10 //
11 // Author: andreas.morsch@cern.ch
12
13 #include <TObject.h>
14 #include <TString.h>
15 #include <TF2.h>
16 class TF1;
17
18 class AliFastGlauber : public TObject {
19  public:
20     AliFastGlauber();
21     virtual ~AliFastGlauber();
22     void Init(Int_t mode = 0);
23
24     void SetWoodSaxonParameters(Double_t r0, Double_t d, Double_t w, Double_t n)
25         {fWSr0 = r0; fWSd = d; fWSw = w; fWSn = n;}
26     void SetWoodSaxonParametersAu()
27         {fWSr0 = 6.38; fWSd = 0.535; fWSw = 0.; fWSn = 8.59e-4;}
28     void SetWoodSaxonParametersPb()
29         {fWSr0 = 6.624; fWSd = 0.549; fWSw = 0.; fWSn = 7.69e-4;}
30     void SetMaxImpact(Float_t bmax = 20.) {fgBMax = bmax;};
31     void SetHardCrossSection(Float_t xs = 1.0) {fSigmaHard = xs;}
32     void SetNNCrossSection  (Float_t xs = 55.6) {fSigmaNN = xs;}
33     void SetNucleus(Int_t n=208) {fA=n;}
34     void SetAuAuRhic();
35     void SetPbPbLHC();
36     void SetFileName(TString &fn){fName=fn;}
37     void SetFileName(const char *fn="$(ALICE_ROOT)/FASTSIM/data/glauberPbPb.root"){fName=fn;}
38
39     static Double_t WSb            (Double_t *xx, Double_t *par);
40     static Double_t WSbz           (Double_t *xx, Double_t *par);
41     static Double_t WSz            (Double_t *xx, Double_t *par);
42     static Double_t WSta           (Double_t *xx, Double_t *par);
43     static Double_t WStarfi        (Double_t *xx, Double_t *par);
44     static Double_t WStaa          (Double_t *xx, Double_t *par);
45     static Double_t WKParticipants (Double_t *xx, Double_t *par);
46     static Double_t WParticipants  (Double_t *xx, Double_t *par);    
47     static Double_t WSgeo          (Double_t *xx, Double_t *par);
48     static Double_t WSbinary       (Double_t *xx, Double_t *par);
49     static Double_t WSN            (Double_t *xx, Double_t *par);
50     static Double_t WAlmond        (Double_t *xx, Double_t *par);
51     static Double_t WPathLength0   (Double_t *xx, Double_t *par);
52     static Double_t WPathLength    (Double_t *xx, Double_t *par);
53     static Double_t WIntRadius     (Double_t *xx, Double_t *par);
54     static Double_t WEnergyDensity (Double_t *xx, Double_t *par);
55
56     const TF1* GetWSB()            const {return fgWSb;}
57     const TF2* GetWSbz()           const {return fgWSbz;}
58     const TF1* GetWSz()            const {return fgWSz;} 
59     const TF1* GetWSta()           const {return fgWSta;}
60     const TF2* Kernel()            const {return fgWStarfi;}
61     const TF2* GetWStarfi()        const {return fgWStarfi;}
62     const TF2* GetWKParticipants() const {return fgWKParticipants;}
63     const TF1* GetWParticipants()  const {return fgWParticipants;} 
64     const TF1* Overlap()           const {return fgWStaa;}
65     const TF1* GetWStaa()          const {return fgWStaa;} 
66     const TF2* GetWAlmond()        const {return fgWAlmond;}
67     const TF1* GetWPathLength0()   const {return fgWPathLength0;} 
68     const TF1* GetWPathLength()    const {return fgWPathLength;}
69     const TF1* GetWIntRadius()     const {return fgWIntRadius;}
70     const TF1* GetWSgeo()          const {return fgWSgeo;}
71     const TF1* GetWSbinary()       const {return fgWSbinary;}
72     const TF1* GetWSN()            const {return fgWSN;}     
73     const TF1* GetWEnergyDensity() const {return fgWEnergyDensity;} 
74     const TF2* GetWAlmondFixedB(Int_t i) const {return &fgWAlmondFixedB[i];}
75     
76     Float_t GetWr0() const {return fWSr0;}
77     Float_t GetWSd() const {return fWSd;}
78     Float_t GetWSw() const {return fWSw;}
79     Float_t GetWSn() const {return fWSn;}
80     Float_t GetSigmaHard()       const {return fSigmaHard;}
81     Float_t GetSigmaNN()         const {return fSigmaNN;}
82     Int_t GetA()                 const {return fA;}
83     const TString* GetFileName() const {return &fName;}
84     Float_t GetBmin() const {return fBmin;}
85     Float_t GetBmax() const {return fBmax;}
86
87     void DrawWSb()          const;
88     void DrawThickness()    const;
89     void DrawOverlap()      const;
90     void DrawParticipants() const;
91     void DrawGeo()          const;
92     void DrawBinary()       const;
93     void DrawN()            const;    
94     void DrawKernel(Double_t b = 0.) const;
95     void DrawAlmond(Double_t b = 0.) const;
96     void DrawPathLength0(Double_t b = 0., Int_t iopt = 0)            const;
97     void DrawPathLength(Double_t b, Int_t ni = 1000, Int_t iopt = 0) const;
98     void DrawIntRadius(Double_t b = 0.) const;
99     void DrawEnergyDensity()            const;
100     
101     Double_t CrossSection(Double_t b1, Double_t b2)               const;
102     Double_t HardCrossSection(Double_t b1, Double_t b2)           const;
103     Double_t NHard(Double_t b1, Double_t b2)                      const;
104     Double_t FractionOfHardCrossSection(Double_t b1, Double_t b2) const;
105     Double_t Binaries(Double_t b)                 const;
106     Double_t GetNumberOfBinaries(Double_t b)      const;
107     Double_t Participants(Double_t b)             const;
108     Double_t GetNumberOfParticipants(Double_t  b) const;
109     Double_t GetNumberOfCollisions(Double_t  b)   const;
110     Double_t GetNumberOfCollisionsPerEvent(Double_t  b) const;
111     void SimulateTrigger(Int_t n);
112     void GetRandom(Float_t& b, Float_t& p, Float_t& mult);
113     void GetRandom(Int_t& bin, Bool_t& hard);
114     Double_t GetRandomImpactParameter(Double_t bmin, Double_t bmax);
115
116     void StoreFunctions() const;
117     void StoreAlmonds()   const;
118
119     void SetLengthDefinition(Int_t def=1) {fEllDef=def;}
120     Int_t GetLengthDef() const {return fEllDef;}
121     void SetCentralityClass(Double_t xsecFrLow=0.0,Double_t xsecFrUp=0.1);    
122     void GetRandomBHard(Double_t& b);
123     void GetRandomXY(Double_t& x,Double_t& y);
124     void GetRandomPhi(Double_t& phi);
125     Double_t CalculateLength(Double_t b=0.,Double_t x0=0.,Double_t y0=0.,
126                              Double_t phi0=0.);
127     void GetLengthAndPhi(Double_t& ell,Double_t &phi,Double_t b=-1.);
128     void GetLength(Double_t& ell,Double_t b=-1.);
129     void GetLengthsBackToBackAndPhi(Double_t& ell1,Double_t& ell2,
130                                     Double_t &phi,
131                                     Double_t b=-1.);
132     void GetLengthsBackToBack(Double_t& ell1,Double_t& ell2,
133                               Double_t b=-1.);
134     void GetLengthsForPythia(Int_t n,Double_t* phi,Double_t* ell,
135                              Double_t b=-1.);
136     void PlotBDistr(Int_t n=1000);
137     void PlotLengthDistr(Int_t n=1000,Bool_t save=kFALSE,
138                          const char *fname="length.root");
139     void PlotLengthB2BDistr(Int_t n=1000,Bool_t save=kFALSE,
140                             const char *fname="lengthB2B.root");
141     void CalculateI0I1(Double_t& integral0,Double_t& integral1,
142                        Double_t b=0.,
143                        Double_t x0=0.,Double_t y0=0.,Double_t phi0=0.,
144                        Double_t ellCut=20.) const;
145     void GetI0I1AndPhi(Double_t& integral0,Double_t& integral1,Double_t &phi,
146                  Double_t ellCut=20.,Double_t b=-1.);
147     void GetI0I1(Double_t& integral0,Double_t& integral1,
148                  Double_t ellCut=20.,Double_t b=-1.);
149     void GetI0I1BackToBackAndPhi(Double_t& integral01,Double_t& integral11,
150                                  Double_t& integral02,Double_t& integral12,
151                                  Double_t& phi,
152                                  Double_t ellCut=20.,Double_t b=-1.);
153     void GetI0I1BackToBackAndPhiAndXY(Double_t& integral01,Double_t& integral11,
154                                       Double_t& integral02,Double_t& integral12,
155                                       Double_t& phi,Double_t& x,Double_t&y,
156                                       Double_t ellCut=20.,Double_t b=-1.);
157     void GetI0I1BackToBack(Double_t& integral01,Double_t& integral11,
158                            Double_t& integral02,Double_t& integral12,
159                            Double_t ellCut=20.,Double_t b=-1.);
160     void GetI0I1ForPythia(Int_t n,Double_t* phi,
161                           Double_t* integral0,Double_t* integral1,
162                           Double_t ellCut=20.,Double_t b=-1.);
163     void GetI0I1ForPythiaAndXY(Int_t n,Double_t* phi,
164                                Double_t* integral0,Double_t* integral1,
165                                Double_t&x, Double_t &y,
166                                Double_t ellCut=20.,Double_t b=-1.);
167     void PlotI0I1Distr(Int_t n=1000,Double_t ellCut=20.,Bool_t save=kFALSE,
168                        const char *fname="i0i1.root");
169     void PlotI0I1B2BDistr(Int_t n=1000,Double_t ellCut=20.,Bool_t save=kFALSE,
170                           const char *fname="i0i1B2B.root");
171     void PlotAlmonds() const;
172  protected:
173     void Reset();
174
175     static Float_t fgBMax;           // Maximum Impact Parameter
176     static Int_t fgCounter;          // Counter to protect double instantiation
177     static const Int_t fgkMCInts;    // Number of MC integrations
178
179     static TF1*    fgWSb;            // Wood-Saxon Function (b)
180     static TF2*    fgWSbz;           // Wood-Saxon Function (b, z)
181     static TF1*    fgWSz;            // Wood-Saxon Function (b = b0, z)
182     static TF1*    fgWSta;           // Thickness Function
183     static TF2*    fgWStarfi;        // Kernel for Overlap Function
184     static TF2*    fgWKParticipants; // Kernel for number of participants
185     static TF1*    fgWParticipants;  // Number of participants
186     static TF1*    fgWStaa;          // Overlap Function
187     static TF2*    fgWAlmond;        // Interaction Almond
188     static TF1*    fgWPathLength0;   // Path Length as a function of phi
189     static TF1*    fgWPathLength;    // Path Length as a function of phi
190     static TF1*    fgWIntRadius;     // Interaction Radius
191     static TF1*    fgWSgeo;          // dSigma/db geometric
192     static TF1*    fgWSbinary;       // dSigma/db binary
193     static TF1*    fgWSN;            // dN/db binary
194     static TF1*    fgWEnergyDensity; // Energy density as a function of impact parameter
195     static TF2     fgWAlmondFixedB[40]; // Interaction Almonds read from file
196     static TF2*    fgWAlmondCurrent;    // Interaction Almond used for length
197     
198     Float_t fWSr0;      // Wood-Saxon Parameter r0
199     Float_t fWSd;       // Wood-Saxon Parameter d
200     Float_t fWSw;       // Wood-Saxon Parameter w
201     Float_t fWSn;       // Wood-Saxon Parameter n
202                         // (chosen such that integral is one)
203     Float_t fSigmaHard; // Hard Cross Section [mbarn]
204     Float_t fSigmaNN;   // NN Cross Section [mbarn]   
205     Int_t fA;           // Nucleon number of nucleus A
206
207     Float_t fBmin;      // Minimum b (set through centrality selection)
208     Float_t fBmax;      // Coresponding maximum b
209
210     Int_t fEllDef;      // definition of length (see CalculateLength())
211     TString fName;     // filename of stored distributions
212     ClassDef(AliFastGlauber,1) // Event geometry simulation in the Glauber Model
213 };
214
215 #endif