]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/AliFMDAnaFlowRing.h
Removed generated files@
[u/mrichter/AliRoot.git] / FMD / scripts / AliFMDAnaFlowRing.h
1 // -*- mode: C++ -*-
2 #ifndef ALIFMDANAFLOWRING_H
3 #define ALIFMDANAFLOWRING_H
4 #include <flow/AliFMDFlowBinned1D.h>
5 #include <flow/AliFMDFlowTrue.h>
6 #include <flow/AliFMDFlowBin.h>
7 #include <flow/AliFMDFlowAxis.h>
8 #include <flow/AliFMDFlowUtil.h>
9 #include <AliFMDAnaRing.h>
10 #include <TH1.h>
11 #include <TH2.h>
12 #include <TProfile.h>
13 #include <TMath.h>
14 #include <TBrowser.h>
15
16 //====================================================================
17 class AliFMDAnaFlowRing : public AliFMDAnaRing
18 {
19 public:
20   /** Constructor 
21       @param d      Detector number 
22       @param r      Ring identifier 
23       @param n      Number of bins
24       @param etamin Minimum @f$\eta@f$
25       @param etamax Maximum @f$\eta@f$
26       @param bg     Background 
27       @param c1     Lower cut 
28       @param c2     higher cut */
29   AliFMDAnaFlowRing(UShort_t d=0, Char_t r='\0', 
30                     Int_t n=0, Float_t etamin=0, Float_t etamax=0,
31                     TH2* bg=0, Float_t c1=0, Float_t c2=0)
32     : AliFMDAnaRing(d, r, bg, c1, c2), 
33       fAxis(n, etamin, etamax),
34       fV2(Name(), Name(), 2, 1, fAxis),
35       fTrue(Form("True%s", Name()), Name(), 2, fAxis),
36       fNObserv(0), 
37       fPsi(Form("%s_psi", fName), "Event planes", 
38            fAxis.N(), fAxis.Bins(), 2*fNSeq, -TMath::Pi(), TMath::Pi()),
39       fDPsi(Form("%s_dpsi", fName), "Distance to real #Psi_{R}", 
40             fAxis.N(), fAxis.Bins(), -TMath::Pi(), TMath::Pi()),
41       fMult(Form("%s_mult", fName), "Multiplicity per hit", 
42             fAxis.N(), fAxis.Bins(), 100, 0, 10)
43   {}
44   /** Constructor 
45       @param d     Detector number 
46       @param r     Ring identifier 
47       @param a     Axis 
48       @param bg    Background 
49       @param c1    Lower cut 
50       @param c2    higher cut */
51   AliFMDAnaFlowRing(UShort_t d, Char_t r, const AliFMDFlowAxis& a,
52                     TH2* bg, Float_t c1, Float_t c2)
53     : AliFMDAnaRing(d, r, bg, c1, c2), 
54       fAxis(a),
55       fV2(Name(), Name(), 2, 1, fAxis),
56       fTrue(Form("True%s", Name()), Name(), 2, fAxis),
57       fNObserv(0), 
58       fPsi(Form("%s_psi", fName), "Event planes", 
59            fAxis.N(), fAxis.Bins(), 2*fNSeq, -TMath::Pi(), TMath::Pi()),
60       fDPsi(Form("%s_dpsi", fName), "Distance to real #Psi_{R}", 
61             fAxis.N(), fAxis.Bins(), -TMath::Pi(), TMath::Pi()),
62       fMult(Form("%s_mult", fName), "Multiplicity per hit", 
63             fAxis.N(), fAxis.Bins(), 100, 0, 10)
64   {}
65   /** Initialize */
66   void Init()
67   {
68     fPhis.Set(10240);
69     fEtas.Set(10240);
70     fWeights.Set(10240);
71     fNObserv = 0;
72
73     fV2.SetLineColor(Color());
74     fV2.SetMarkerColor(Color());
75     fV2.SetFillColor(Color());
76     fV2.SetFillStyle(3003+(fRing == 'I' ? 1 : 2));
77
78     fPsi.SetXTitle("#eta");
79     fPsi.SetYTitle("#Psi");
80     fPsi.SetZTitle("Frequency");
81
82     fDPsi.SetXTitle("#eta");
83     fDPsi.SetYTitle("#Delta#Psi=|#Psi-#Psi_{R}|");
84
85     fMult.SetXTitle("#eta");
86     fMult.SetYTitle("M");
87
88     ModHist(fPsi);
89     ModHist(fDPsi);
90     ModHist(fMult);
91
92     AliFMDAnaRing::Init();
93   }
94   /** Modify histogram 
95       @param h Histogram */
96   void ModHist(TH1& h)
97   {
98     h.GetXaxis()->SetTitleFont(132);
99     h.GetXaxis()->SetLabelFont(132);
100     h.GetXaxis()->SetNdivisions(8);
101     h.GetYaxis()->SetTitleFont(132);
102     h.GetYaxis()->SetLabelFont(132);
103     h.GetYaxis()->SetNdivisions(8);
104     h.GetZaxis()->SetTitleFont(132);
105     h.GetZaxis()->SetLabelFont(132);
106     h.GetZaxis()->SetNdivisions(8);
107     h.SetLineColor(Color());
108     h.SetMarkerColor(Color());
109     h.SetFillColor(Color());
110     h.SetMarkerStyle(20);
111     h.SetDirectory(0);
112   }
113   /** Called at beginning of event */
114   void Begin() { fNObserv = 0; }
115   /** Fill in one strip 
116       @param phi   Azimuthal angle @f$ \varphi@f$ 
117       @param eta   Pseudo rapidity @f$ \eta@f$
118       @param mult  Signal */
119   void Fill(Float_t phi, Float_t eta, Float_t mult)
120   {
121     // if (mult < 0.2 || mult > 2) return;
122     // if (fBg && mult < fCut0) return;
123     if (mult <= 0.01) return;
124     Int_t imult = Int_t(mult+.1);
125     if (fNObserv + imult >= fEtas.fN) { 
126       ULong_t n = Int_t(fEtas.fN * 1.5)+imult;
127       fEtas.Set(n);
128       fPhis.Set(n);
129       fWeights.Set(n);
130     }
131     for (Int_t i = 0; i < imult; i++) { 
132       fEtas[fNObserv]    = eta;
133       fPhis[fNObserv]    = phi;
134       fWeights[fNObserv] = 1;
135       fMult.Fill(eta, mult);
136       fNObserv++;
137     }
138   }
139   /** Called at end of event */
140   void End() 
141   { 
142     fTrue.Event(fNObserv, fPhis.fArray, fEtas.fArray, 
143                 fWeights.fArray, 0); // fWeights.fArray);
144     fV2.Event(fNObserv, fPhis.fArray, fEtas.fArray, 
145               fWeights.fArray, 0); // fWeights.fArray);
146     for (UShort_t i = 0; i < fAxis.N(); i++) { 
147       AliFMDFlowBin* bin  = fV2.GetBin(i);
148       AliFMDFlowBin* tbin = fTrue.GetBin(i);
149       if (bin && bin->Counts() > 0) {
150         Double_t psir = tbin->Psi() / fV2.PsiOrder();
151         fPsi.Fill(fAxis.BinCenter(i), bin->Psi() - psir);
152         fDPsi.Fill(fAxis.BinCenter(i), bin->Psi() - psir);
153       }
154     }
155   }  
156   /** Browse this object 
157       @param b Browser */
158   void Browse(TBrowser* b)
159   {
160     AliFMDAnaRing::Browse(b);
161     b->Add(&fV2,   "V2");
162     b->Add(&fTrue, "TrueV2");
163     b->Add(&fPsi);
164     b->Add(&fDPsi);
165     b->Add(&fMult);
166   }
167 protected:
168   /** Flow axis */
169   AliFMDFlowAxis     fAxis;
170   /** Flow histogram */
171   AliFMDFlowBinned1D fV2;
172   /** Flow histogram */
173   AliFMDFlowTrue1D   fTrue;
174   /** Cache of phi */
175   TArrayD            fPhis;
176   /** Cache of eta */
177   TArrayD            fEtas;
178   /** Cache of weights */
179   TArrayD            fWeights;
180   /** Number of observables */
181   Int_t              fNObserv;
182   /** Histogram of Psi */
183   TH2D               fPsi;       // Event planes
184   /** Histogram of dPsi */
185   TProfile           fDPsi;   // Mean distance to Psi
186   /** Histogram of multiplicities */
187   TH2D               fMult;
188
189   ClassDef(AliFMDAnaFlowRing,1) // Flow ring
190 };
191
192 #endif
193 //
194 // EOF
195 //