]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/Base/AliUEHistograms.h
pt,a < pt,t condition optional
[u/mrichter/AliRoot.git] / PWGCF / Correlations / Base / AliUEHistograms.h
1 #ifndef AliUEHistograms_H
2 #define AliUEHistograms_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 /* $Id: AliUEHistograms.h 20164 2007-08-14 15:31:50Z morsch $ */
8
9 // encapsulates several AliUEHist objects for a full UE analysis plus additional control histograms
10
11 #include "TNamed.h"
12 #include "AliUEHist.h"
13 #include "TMath.h"
14 #include "THn.h" // in cxx file causes .../THn.h:257: error: conflicting declaration ‘typedef class THnT<float> THnF’
15
16 class AliVParticle;
17
18 class TList;
19 class TSeqCollection;
20 class TObjArray;
21 class TH1F;
22 class TH2F;
23 class TH3F;
24
25 class AliUEHistograms : public TNamed
26 {
27  public:
28   AliUEHistograms(const char* name = "AliUEHistograms", const char* histograms = "", const char* binning = 0);
29   virtual ~AliUEHistograms();
30   
31   void Fill(Int_t eventType, Float_t zVtx, AliUEHist::CFStep step, AliVParticle* leading, TList* toward, TList* away, TList* min, TList* max);
32   void FillCorrelations(Double_t centrality, Float_t zVtx, AliUEHist::CFStep step, TObjArray* particles, TObjArray* mixed = 0, Float_t weight = 1, Bool_t firstTime = kTRUE, Bool_t twoTrackEfficiencyCut = kFALSE, Float_t bSign = 0, Float_t twoTrackEfficiencyCutValue = 0.02, Bool_t applyEfficiency = kFALSE);
33   void Fill(AliVParticle* leadingMC, AliVParticle* leadingReco);
34   void FillEvent(Int_t eventType, Int_t step);
35   void FillEvent(Double_t centrality, Int_t step);
36   void FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, TObjArray* recoPrimPID, TObjArray* recoAllPID, TObjArray* fake, Int_t particleType, Double_t centrality = 0, Double_t zVtx = 0);
37   void FillFakePt(TObjArray* fake, Double_t centrality);
38  
39   void CopyReconstructedData(AliUEHistograms* from);
40   void DeepCopy(AliUEHistograms* from);
41   
42   AliUEHist* GetUEHist(Int_t id);
43   
44   AliUEHist* GetNumberDensitypT() { return fNumberDensitypT; }
45   AliUEHist* GetSumpT() { return fSumpT; }
46   AliUEHist* GetNumberDensityPhi() { return fNumberDensityPhi; }
47   
48   void SetNumberDensitypT(AliUEHist* obj) { fNumberDensitypT = obj; }
49   void SetSumpT(AliUEHist* obj) { fSumpT = obj; }
50   void SetNumberDensityPhi(AliUEHist* obj) { fNumberDensityPhi = obj; }
51   
52   void SetRunNumber(Long64_t runNumber) { fRunNumber = runNumber; }
53   
54   void SetEfficiencyCorrectionTriggers(THnF* hist)   { fEfficiencyCorrectionTriggers = hist;   }
55   void SetEfficiencyCorrectionAssociated(THnF* hist) { fEfficiencyCorrectionAssociated = hist; }
56   
57   TH2F* GetCorrelationpT()  { return fCorrelationpT; }
58   TH2F* GetCorrelationEta() { return fCorrelationEta; }
59   TH2F* GetCorrelationPhi() { return fCorrelationPhi; }
60   TH2F* GetCorrelationR()   { return fCorrelationR; }
61   TH2F* GetCorrelationLeading2Phi() { return fCorrelationLeading2Phi; }
62   TH2F* GetCorrelationMultiplicity() { return fCorrelationMultiplicity; }
63   TH3F* GetYield() { return fYields; }
64   
65   TH2F* GetEventCount()     { return fEventCount; }
66   TH3F* GetEventCountDifferential() { return fEventCountDifferential; }
67   TH1F* GetVertexContributors() { return fVertexContributors; }
68   TH1F* GetCentralityDistribution() { return fCentralityDistribution; }
69   TH2F* GetCentralityCorrelation() { return fCentralityCorrelation; }
70   Long64_t GetRunNumber() { return fRunNumber; }
71   Int_t GetMergeCount() { return fMergeCount; }
72   TH3F* GetTwoTrackDistance(Int_t i) { return fTwoTrackDistancePt[i]; }
73   Bool_t GetWeightPerEvent() { return fWeightPerEvent; }
74   
75   void Correct(AliUEHistograms* corrections);
76   
77   void SetEtaRange(Float_t etaMin, Float_t etaMax);
78   void SetPtRange(Float_t ptMin, Float_t ptMax);
79   void SetPartSpecies(Int_t species);
80   void SetZVtxRange(Float_t min, Float_t max);
81   void SetContaminationEnhancement(TH1F* hist);
82   void SetCombineMinMax(Bool_t flag);
83   void SetTrackEtaCut(Float_t value);
84   void SetWeightPerEvent(Bool_t flag);
85   void SetSelectCharge(Int_t selectCharge) { fSelectCharge = selectCharge; }
86   void SetSelectTriggerCharge(Int_t selectCharge) { fTriggerSelectCharge = selectCharge; }
87   void SetSelectAssociatedCharge(Int_t selectCharge) { fAssociatedSelectCharge = selectCharge; }
88   void SetTriggerRestrictEta(Float_t eta) { fTriggerRestrictEta = eta; }
89   void SetEtaOrdering(Bool_t flag) { fEtaOrdering = flag; }
90   void SetPairCuts(Bool_t conversions, Bool_t resonances) { fCutConversions = conversions; fCutResonances = resonances; }
91   void SetOnlyOneEtaSide(Int_t flag)    { fOnlyOneEtaSide = flag; }
92   void SetPtOrder(Bool_t flag) { fPtOrder = flag; }
93   
94   void ExtendTrackingEfficiency(Bool_t verbose = kFALSE);
95   void Reset();
96
97   AliUEHistograms(const AliUEHistograms &c);
98   AliUEHistograms& operator=(const AliUEHistograms& c);
99   virtual void Copy(TObject& c) const;
100
101   virtual Long64_t Merge(TCollection* list);
102   void Scale(Double_t factor);
103   
104 protected:
105   void FillRegion(AliUEHist::Region region, Float_t zVtx, AliUEHist::CFStep step, AliVParticle* leading, TList* list, Int_t multiplicity);
106   Int_t CountParticles(TList* list, Float_t ptMin);
107   void DeleteContainers();
108   inline Float_t GetInvMassSquared(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2);
109   inline Float_t GetInvMassSquaredCheap(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2);
110   inline Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1, Float_t phi2, Float_t pt2, Float_t charge2, Float_t radius, Float_t bSign);
111   
112   static const Int_t fgkUEHists; // number of histograms
113
114   AliUEHist* fNumberDensitypT;   // d^2N/dphideta vs pT,lead
115   AliUEHist* fSumpT;             // d^2 sum(pT)/dphideta vs pT,lead
116   AliUEHist* fNumberDensityPhi;  // d^2N/dphideta vs delta phi,lead (in pT,lead bins)
117   
118   TH2F* fCorrelationpT;         // pT,lead: true vs reco
119   TH2F* fCorrelationEta;        // #eta,lead; true vs reco
120   TH2F* fCorrelationPhi;        // #phi,lead; true vs reco
121   TH2F* fCorrelationR;          // R = sqrt(delta eta^2 + delta phi^2) (true vs reco) vs pT,lead,MC
122   TH2F* fCorrelationLeading2Phi;// delta phi (true vs reco) vs pT,lead,MC
123   TH2F* fCorrelationMultiplicity; // number of mc particls vs reco particles (for pT > 0.5 GeV/c)
124   TH3F* fYields;                // centrality vs pT vs eta
125   
126   TH2F* fEventCount;            // event count as function of step, (for pp: event type (plus additional step -1 for all events without vertex range even in MC)) (for PbPb: centrality)
127   TH3F* fEventCountDifferential;// event count as function of leading pT, step, event type
128   
129   TH1F* fVertexContributors;    // number of contributors to the vertex
130   TH1F* fCentralityDistribution; // distribution of the variable used for centrality selection
131   TH2F* fCentralityCorrelation;  // centrality vs multiplicity
132   
133   TH3F* fITSClusterMap;          // its cluster map vs centrality vs pT
134   
135   TH3F* fTwoTrackDistancePt[2];    // control histograms for two-track efficiency study: dphi*_min vs deta (0 = before cut, 1 = after cut)
136   TH2F* fControlConvResoncances; // control histograms for cuts on conversions and resonances
137   
138   THnF* fEfficiencyCorrectionTriggers;   // if non-0 this efficiency correction is applied on the fly to the filling for trigger particles. The factor is multiplicative, i.e. should contain 1/efficiency
139   THnF* fEfficiencyCorrectionAssociated;   // if non-0 this efficiency correction is applied on the fly to the filling for associated particles. The factor is multiplicative, i.e. should contain 1/efficiency
140   
141   Int_t fSelectCharge;           // (un)like sign selection when building correlations: 0: no selection; 1: unlike sign; 2: like sign
142   Int_t fTriggerSelectCharge;    // select charge of trigger particle
143   Int_t fAssociatedSelectCharge; // select charge of associated particle
144   Float_t fTriggerRestrictEta;   // restrict eta range for trigger particle (default: -1 [off])
145   Bool_t fEtaOrdering;           // activate eta ordering to prevent shape distortions. see FillCorrelation for the details
146   Bool_t fCutConversions;        // cut on conversions (inv mass)
147   Bool_t fCutResonances;         // cut on resonances (inv mass)
148   Int_t fOnlyOneEtaSide;       // decides that only trigger particle from one eta side are considered (0 = all; -1 = negative, 1 = positive)
149   Bool_t fWeightPerEvent;       // weight with the number of trigger particles per event
150   Bool_t fPtOrder;              // apply pT,a < pT,t condition
151   
152   Long64_t fRunNumber;           // run number that has been processed
153   
154   Int_t fMergeCount;            // counts how many objects have been merged together
155   
156   ClassDef(AliUEHistograms, 25)  // underlying event histogram container
157 };
158
159 Float_t AliUEHistograms::GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1, Float_t phi2, Float_t pt2, Float_t charge2, Float_t radius, Float_t bSign)
160
161   //
162   // calculates dphistar
163   //
164   
165   Float_t dphistar = phi1 - phi2 - charge1 * bSign * TMath::ASin(0.075 * radius / pt1) + charge2 * bSign * TMath::ASin(0.075 * radius / pt2);
166   
167   static const Double_t kPi = TMath::Pi();
168   
169   // circularity
170 //   if (dphistar > 2 * kPi)
171 //     dphistar -= 2 * kPi;
172 //   if (dphistar < -2 * kPi)
173 //     dphistar += 2 * kPi;
174   
175   if (dphistar > kPi)
176     dphistar = kPi * 2 - dphistar;
177   if (dphistar < -kPi)
178     dphistar = -kPi * 2 - dphistar;
179   if (dphistar > kPi) // might look funny but is needed
180     dphistar = kPi * 2 - dphistar;
181   
182   return dphistar;
183 }
184
185 Float_t AliUEHistograms::GetInvMassSquared(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2)
186 {
187   // calculate inv mass squared
188   // same can be achieved, but with more computing time with
189   /*TLorentzVector photon, p1, p2;
190   p1.SetPtEtaPhiM(triggerParticle->Pt(), triggerEta, triggerParticle->Phi(), 0.510e-3);
191   p2.SetPtEtaPhiM(particle->Pt(), eta[j], particle->Phi(), 0.510e-3);
192   photon = p1+p2;
193   photon.M()*/
194   
195   Float_t tantheta1 = 1e10;
196   
197   if (eta1 < -1e-10 || eta1 > 1e-10)
198   {
199     Float_t expTmp = TMath::Exp(-eta1);
200     tantheta1 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
201   }
202   
203   Float_t tantheta2 = 1e10;
204   if (eta2 < -1e-10 || eta2 > 1e-10)
205   {
206     Float_t expTmp = TMath::Exp(-eta2);
207     tantheta2 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
208   }
209   
210   Float_t e1squ = m0_1 * m0_1 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
211   Float_t e2squ = m0_2 * m0_2 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
212   
213   Float_t mass2 = m0_1 * m0_1 + m0_2 * m0_2 + 2 * ( TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( TMath::Cos(phi1 - phi2) + 1.0 / tantheta1 / tantheta2 ) ) );
214   
215 //   Printf(Form("%f %f %f %f %f %f %f %f %f", pt1, eta1, phi1, pt2, eta2, phi2, m0_1, m0_2, mass2));
216   
217   return mass2;
218 }
219
220 Float_t AliUEHistograms::GetInvMassSquaredCheap(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2)
221 {
222   // calculate inv mass squared approximately
223   
224   Float_t tantheta1 = 1e10;
225   
226   if (eta1 < -1e-10 || eta1 > 1e-10)
227   {
228     Float_t expTmp = 1.0-eta1+eta1*eta1/2-eta1*eta1*eta1/6+eta1*eta1*eta1*eta1/24;
229     tantheta1 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
230   }
231   
232   Float_t tantheta2 = 1e10;
233   if (eta2 < -1e-10 || eta2 > 1e-10)
234   {
235     Float_t expTmp = 1.0-eta2+eta2*eta2/2-eta2*eta2*eta2/6+eta2*eta2*eta2*eta2/24;
236     tantheta2 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
237   }
238   
239   Float_t e1squ = m0_1 * m0_1 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
240   Float_t e2squ = m0_2 * m0_2 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
241   
242   // fold onto 0...pi
243   Float_t deltaPhi = TMath::Abs(phi1 - phi2);
244   while (deltaPhi > TMath::TwoPi())
245     deltaPhi -= TMath::TwoPi();
246   if (deltaPhi > TMath::Pi())
247     deltaPhi = TMath::TwoPi() - deltaPhi;
248   
249   Float_t cosDeltaPhi = 0;
250   if (deltaPhi < TMath::Pi()/3)
251     cosDeltaPhi = 1.0 - deltaPhi*deltaPhi/2 + deltaPhi*deltaPhi*deltaPhi*deltaPhi/24;
252   else if (deltaPhi < 2*TMath::Pi()/3)
253     cosDeltaPhi = -(deltaPhi - TMath::Pi()/2) + 1.0/6 * TMath::Power((deltaPhi - TMath::Pi()/2), 3);
254   else
255     cosDeltaPhi = -1.0 + 1.0/2.0*(deltaPhi - TMath::Pi())*(deltaPhi - TMath::Pi()) - 1.0/24.0 * TMath::Power(deltaPhi - TMath::Pi(), 4);
256   
257   Float_t mass2 = m0_1 * m0_1 + m0_2 * m0_2 + 2 * ( TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( cosDeltaPhi + 1.0 / tantheta1 / tantheta2 ) ) );
258   
259 //   Printf(Form("%f %f %f %f %f %f %f %f %f", pt1, eta1, phi1, pt2, eta2, phi2, m0_1, m0_2, mass2));
260   
261   return mass2;
262 }
263
264 #endif