]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/PIDFluctuation/task/AliEbyEPidRatioHelper.h
Minor fix pr task: sjena
[u/mrichter/AliRoot.git] / PWGCF / EBYE / PIDFluctuation / task / AliEbyEPidRatioHelper.h
1 #ifndef ALIEBYEPIDRATIOHELPER_H
2 #define ALIEBYEPIDRATIOHELPER_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6 //=========================================================================//
7 //             AliEbyE Analysis for Particle Ratio Fluctuation             //
8 //                   Deepika Rathee  | Satyajit Jena                       //
9 //                   drathee@cern.ch | sjena@cern.ch                       //
10 //                  Date: Wed Jul  9 18:38:30 CEST 2014                    // 
11 //          New approch to find particle ratio to reduce memory            //
12 //                             (Test Only)                                 //
13 //=========================================================================//
14
15 #include "THnBase.h"
16 #include "THn.h"
17 #include "TH1F.h"
18 #include "TF1.h"
19 #include "TProfile2D.h"
20 #include "TRandom3.h"
21
22 class AliESDtrack;
23 class AliMCEvent;
24 class AliStack;
25 class AliPIDResponse;
26 class AliESDtrackCuts;
27 class AliInputEventHandler;
28 class AliESDInputHandler;
29 class AliAODInputHandler;
30 class AliAODEvent;
31 class AliAODTrack;
32 class AliAODMCParticle;
33 class AliMCParticle;
34
35 class AliEbyEPidRatioHelper : public TNamed {
36
37  public:
38
39   AliEbyEPidRatioHelper();
40   virtual ~AliEbyEPidRatioHelper();
41
42   void SetCentralityBinMax(Int_t d)            {fCentralityBinMax    = d;}
43   void SetVertexZMax(Float_t f)                {fVertexZMax          = f;}
44   void SetRapidityMax(Float_t f)               {fRapidityMax         = f;}
45   void SetMinTrackLengthMC(Float_t f)          {fMinTrackLengthMC    = f;}
46   void SetNSigmaMaxCdd(Float_t f)              {fNSigmaMaxCdd        = f;}
47   void SetNSigmaMaxCzz(Float_t f)              {fNSigmaMaxCzz        = f;}
48   void SetPIDStrategy(Int_t i)                 {fPIDStrategy         = i;}
49   void SetNSigmaMaxITS(Float_t f)              {fNSigmaMaxITS        = f;}
50   void SetNSigmaMaxTPC(Float_t f)              {fNSigmaMaxTPC        = f;}
51   void SetNSigmaMaxTPClow(Float_t f)           {fNSigmaMaxTPClow     = f;}
52   void SetNSigmaMaxTOF(Float_t f)              {fNSigmaMaxTOF        = f;}
53   void SetMinPtForTOFRequired(Float_t f)       {fMinPtForTOFRequired = f;}
54   void SetMaxPtForTPClow(Float_t f)            {fMaxPtForTPClow      = f;}
55
56   
57   TH1F*    GetHEventStat0()                    {return fHEventStat0;     }
58   TH1F*    GetHEventStat1()                    {return fHEventStat1;     }
59   TH1F*    GetHTriggerStat()                   {return fHTriggerStat;    }
60   TH1F*    GetHCentralityStat()                {return fHCentralityStat; }
61   TH1F*    GetHCentralityPercentile()          {return fHCentralityPer; }
62   TH1F*    GetHCentralityPercentileAll()       {return fHCentralityPerAll; }
63   Int_t    GetCentralityBin()                  {return fCentralityBin;   }
64   Float_t  GetMaxPtForTPClow()                 {return fMaxPtForTPClow;  }
65   Float_t  GetRapidityMax()                    {return fRapidityMax;     }
66   Float_t  GetPhiMin()                         {return fPhiMin;          }
67   Float_t  GetPhiMax()                         {return fPhiMax;          }
68   AliESDtrackCuts* GetESDTrackCuts()           {return fESDTrackCuts;    }
69   Bool_t           GetIsMC()                   {return fIsMC;            }
70   Bool_t           GetIsRatio()                {return fIsRatio;         }
71   Bool_t           GetIsPtBin()                {return fIsPtBin;         }
72   Bool_t           GetDetWise()                {return fIsDetectorWise ; }
73
74   Int_t            GetAODtrackCutBit()         {return fAODtrackCutBit;  }
75   AliMCEvent*           GetMCEvent()           {return fMCEvent;         }
76   AliInputEventHandler* GetInputEventHandler() {return fInputEventHandler;}
77   void SetPhiRange(Float_t f1, Float_t f2);
78   Float_t  GetCentralityPercentile()           {return fCentralityPercentile;}
79   Float_t  GetMinPtForTOFRequired()            {return fMinPtForTOFRequired;}
80   /** Initialize Helper */
81   Int_t Initialize(AliESDtrackCuts *cuts, Bool_t isMC, Bool_t isRatio, Bool_t isPtBin, Bool_t isDetWise, Int_t trackCutBit, Int_t modeDistCreation);
82
83   /** Setup Event */
84   Int_t SetupEvent(AliESDInputHandler *esdHandler, AliAODInputHandler *aodHandler, AliMCEvent *mcEvent);
85
86   /** Check if event is triggred */
87   Bool_t IsEventTriggered();
88
89   /** Fill event cut statistics */
90   Bool_t IsEventRejected();
91
92   /** Check if charged MC particle is accepted for basic parameters */
93   Bool_t IsParticleAcceptedBasicCharged(AliVParticle *particle, Int_t idxMC);
94
95   /** Check if neutral MC particle is accepted for basic parameters */
96   Bool_t IsParticleAcceptedBasicNeutral(AliVParticle *particle, Int_t idxMC);
97  
98   /** Check if MC particle is accepted for Rapidity */
99   Bool_t IsParticleAcceptedRapidity(AliVParticle *particle, Double_t &yP, Int_t gCurPid);
100
101   /** Check if MC particle is accepted for Phi */
102   Bool_t IsParticleAcceptedPhi(AliVParticle *particle);
103
104   /** Check if MC particle is findable tracks */
105   Bool_t IsParticleFindable(Int_t label);
106   
107   /** Check if track is accepted for basic parameters */
108   Bool_t IsTrackAcceptedBasicCharged(AliVTrack *track);
109   
110   /** Check if track is accepted for Rapidity */
111   Bool_t IsTrackAcceptedRapidity(AliVTrack *track, Double_t &yP, Int_t gCurPid);
112
113   /** Check if track is accepted for DCA */
114   Bool_t IsTrackAcceptedDCA(AliVTrack *track);
115
116   /** Check if track is accepted for PID */
117   Bool_t IsTrackAcceptedPID(AliVTrack *track, Double_t *pid, AliPID::EParticleType gCurPid);
118
119   /** Check if trackis  accepted for Phi */
120   Bool_t IsTrackAcceptedPhi(AliVTrack *track);
121
122   /** Method for the correct logarithmic binning of histograms 
123    *  and Update MinPtForTOFRequired, using the pT log-scale 
124    */
125   void BinLogAxis(const THnBase *h, Int_t axisNumber, AliESDtrackCuts* cuts = NULL);
126
127   // void SetIsRatio(Bool_t b)                    {fIsRatio             = b;}  
128   // void SetIsPtBin(Bool_t   b)                  {fIsPtBin             = b;}  
129   
130   static const Float_t fgkfHistBinWitdthRap;   // Histogram std bin width for rapidity/eta
131   static const Float_t fgkfHistBinWitdthPt;    // Histogram std bin width for pt
132
133   static const Float_t fgkfHistRangeCent[];    // Histogram range for centrality
134   static const Int_t   fgkfHistNBinsCent;      // Histogram N bins for centrality
135   static const Float_t fgkfHistRangeEta[];     // Histogram range for eta
136   static const Int_t   fgkfHistNBinsEta;       // Histogram N bins for eta
137   static const Float_t fgkfHistRangeRap[];     // Histogram range for rapidity
138   static const Int_t   fgkfHistNBinsRap;       // Histogram N bins for rapidity
139   static const Float_t fgkfHistRangePhi[];     // Histogram range for phi
140   static const Int_t   fgkfHistNBinsPhi;       // Histogram N bins for phi
141   static const Float_t fgkfHistRangePt[];      // Histogram range for pt
142   static const Int_t   fgkfHistNBinsPt;        // Histogram N bins for pt
143   static const Float_t fgkfHistRangeSign[];    // Histogram range for sign
144   static const Int_t   fgkfHistNBinsSign;      // Histogram N bins for sign
145
146   static const Char_t* fgkEventNames[];         // Event names 
147   static const Char_t* fgkCentralityMaxNames[]; // Centrality names 
148   static const Char_t* fgkTriggerNames[];       // Trigger names 
149   static const Char_t* fgkCentralityNames[];    // Centrality names  
150   static const Char_t* fgkPidName[4];
151   static const Char_t* fgkPidLatex[4][2];
152   static const Char_t* fgkPidTitles[4][2];
153  
154  private:
155
156   AliEbyEPidRatioHelper(const AliEbyEPidRatioHelper&); // not implemented
157   AliEbyEPidRatioHelper& operator=(const AliEbyEPidRatioHelper&); // not implemented
158
159   void InitializeEventStats();
160   void InitializeTriggerStats();
161   void InitializeCentralityStats();
162   Bool_t FillEventStats(Int_t *aEventCuts);
163
164   Int_t                 fModeDistCreation;         //  Dist creation mode       : 1 = on | 0 = off 
165   AliInputEventHandler *fInputEventHandler;        //! Ptr to input event handler (ESD or AOD)
166   AliPIDResponse       *fPIDResponse;              //! Ptr to PID response Object
167   AliESDEvent          *fESD;                      //! Ptr to ESD event
168   AliESDtrackCuts      *fESDTrackCuts;             //! Ptr to ESD cuts  
169   AliAODEvent          *fAOD;                      //! Ptr to AOD event
170   Int_t                 fAODtrackCutBit;           //  Track filter bit for AOD tracks
171   Bool_t                fIsMC;                     //  Is MC event
172   AliMCEvent           *fMCEvent;                  //! Ptr to MC event
173   AliStack             *fStack;                    //! Ptr to stack
174   Int_t                 fCentralityBin;            //  Centrality bin of current event within max centrality bin
175   Float_t               fCentralityPercentile;     //  Centrality percentile of current event
176   Int_t                 fCentralityBinMax;         //  Max centrality bin to be used
177   Float_t               fVertexZMax;               //  VertexZ cut
178   Float_t               fRapidityMax;              //  Rapidity cut
179   Float_t               fPhiMin;                   //  Phi min cut
180   Float_t               fPhiMax;                   //  Phi max cut
181   Float_t               fMinTrackLengthMC;         //  Min track length for MC tracks
182   Float_t               fNSigmaMaxCdd;             //  N Sigma for dcar / sqrt(cdd) - turn off with 0.
183   Float_t               fNSigmaMaxCzz;             //  N Sigma for dcaz / sqrt(czz) - turn off with 0.
184  
185   Int_t                 fPIDStrategy;              //  PID Strategy to be used
186   Float_t               fNSigmaMaxITS;             //  N Sigma for ITS PID
187   Float_t               fNSigmaMaxTPC;             //  N Sigma for TPC PID
188   Float_t               fNSigmaMaxTPClow;          //  N Sigma for TPC PID lower part
189   Float_t               fNSigmaMaxTOF;             //  N Sigma for TOF PID
190   Float_t               fMinPtForTOFRequired;      //  Min pt from where TOF is required
191   Float_t               fMaxPtForTPClow;           //  Max pt until TPClow is used
192   TH1F                 *fHEventStat0;              //  Event cut statistics
193   TH1F                 *fHEventStat1;              //  Event cut statistics - incremental
194   Int_t                 fHEventStatMax;            //  Max N cuts to be included in HEventStat
195   TH1F                 *fHTriggerStat;             //  Trigger statistics
196   Int_t                 fNTriggers;                //  N triggers used
197   TH1F                 *fHCentralityStat;          //  Centrality statistics
198   TH1F                 *fHCentralityPer;           //  Centrality Percentile
199   TH1F                 *fHCentralityPerAll;        //  Centrality Percentile
200   Int_t                 fNCentralityBins;          //  N centrality bins used
201   TRandom3             *fRandom;                   //  Random generator
202   Bool_t                fIsRatio;                  //  Is ratio
203   Bool_t                fIsPtBin;                  //  Is Pt Bin
204   Bool_t                fIsDetectorWise;           // is detector wise
205
206   ClassDef(AliEbyEPidRatioHelper, 1);
207 };
208
209 #endif