]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/CaloTrackCorrelations/AliAnaPi0.h
add pt dependent histograms where already exist energy dependent
[u/mrichter/AliRoot.git] / PWGGA / CaloTrackCorrelations / AliAnaPi0.h
1 #ifndef ALIANAPI0_H
2 #define ALIANAPI0_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5
6 //_________________________________________________________________________
7 // Class to fill two-photon invariant mass histograms 
8 // to be used to extract pi0 raw yield.
9 // Input is produced by AliAnaPhoton (or any other analysis producing output AliAODPWG4Particles), 
10 // it will do nothing if executed alone
11 //
12 //-- Author: Dmitri Peressounko (RRC "KI")
13 //-- Adapted to CaloTrackCorr frame by Lamia Benhabib (SUBATECH)
14 //-- and Gustavo Conesa (INFN-Frascati)
15
16 //Root
17 class TList;
18 class TH3F ;
19 class TH2F ;
20 class TObjString;
21
22 //Analysis
23 #include "AliAnaCaloTrackCorrBaseClass.h"
24 class AliAODEvent ;
25 class AliESDEvent ;
26 class AliAODPWG4Particle ;
27
28 class AliAnaPi0 : public AliAnaCaloTrackCorrBaseClass {
29   
30  public:   
31   AliAnaPi0() ; // default ctor
32   virtual ~AliAnaPi0() ;//virtual dtor
33   
34   //-------------------------------
35   // General analysis frame methods
36   //-------------------------------
37
38   TObjString * GetAnalysisCuts();
39   
40   TList      * GetCreateOutputObjects(); 
41   
42   void         Print(const Option_t * opt) const;
43   
44   void         MakeAnalysisFillHistograms();
45   
46   void         InitParameters();
47
48   //Calorimeter options
49   TString      GetCalorimeter()         const   { return fCalorimeter           ; }
50   void         SetCalorimeter(TString & det)    { fCalorimeter         = det    ; }
51   void         SetNumberOfModules(Int_t nmod)   { fNModules            = nmod   ; }
52   
53   //-------------------------------
54   // EVENT Bin Methods
55   //-------------------------------
56
57   Int_t        GetEventIndex(AliAODPWG4Particle * part, Double_t * vert)  ;  
58
59   //-------------------------------
60         //Opening angle pair selection
61   //-------------------------------
62   void         SwitchOnAngleSelection()         { fUseAngleCut         = kTRUE  ; }
63   void         SwitchOffAngleSelection()        { fUseAngleCut         = kFALSE ; }
64   
65   void         SwitchOnAngleEDepSelection()     { fUseAngleEDepCut     = kTRUE  ; }
66   void         SwitchOffAngleEDepSelection()    { fUseAngleEDepCut     = kFALSE ; }
67     
68   void         SetAngleCut(Float_t a)           { fAngleCut            = a      ; }
69   void         SetAngleMaxCut(Float_t a)        { fAngleMaxCut         = a      ; }
70
71   void         SwitchOnFillAngleHisto()         { fFillAngleHisto      = kTRUE  ; }
72   void         SwitchOffFillAngleHisto()        { fFillAngleHisto      = kFALSE ; }
73   
74   //------------------------------------------
75   //Do analysis only with clusters in same SM or different combinations of SM
76   //------------------------------------------
77   void         SwitchOnSameSM()                 { fSameSM              = kTRUE  ; }
78   void         SwitchOffSameSM()                { fSameSM              = kFALSE ; }
79   
80   void         SwitchOnSMCombinations()         { fFillSMCombinations  = kTRUE  ; }
81   void         SwitchOffSMCombinations()        { fFillSMCombinations  = kFALSE ; }
82   
83   //-------------------------------
84   //Histogram filling options off by default
85   //-------------------------------
86   void         SwitchOnInvPtWeight()            { fMakeInvPtPlots      = kTRUE  ; }
87   void         SwitchOffInvPtWeight()           { fMakeInvPtPlots      = kFALSE ; }
88   
89   void         SwitchOnFillBadDistHisto()       { fFillBadDistHisto    = kTRUE  ; }
90   void         SwitchOffFillBadDistHisto()      { fFillBadDistHisto    = kFALSE ; }
91   
92   //-------------------------------------------
93   //Cuts for multiple analysis, off by default
94   //-------------------------------------------
95   void         SwitchOnMultipleCutAnalysis()    { fMultiCutAna         = kTRUE  ; }
96   void         SwitchOffMultipleCutAnalysis()   { fMultiCutAna         = kFALSE ; }
97
98   void         SetNPtCuts   (Int_t s)           { if(s <= 10)fNPtCuts    = s    ; }
99   void         SetNAsymCuts (Int_t s)           { if(s <= 10)fNAsymCuts  = s    ; }
100   void         SetNNCellCuts(Int_t s)           { if(s <= 10)fNCellNCuts = s    ; }
101   void         SetNPIDBits  (Int_t s)           { if(s <= 10)fNPIDBits   = s    ; }
102   
103   void         SetPtCutsAt  (Int_t p,Float_t v) { if(p < 10)fPtCuts[p]   = v    ; }
104   void         SetAsymCutsAt(Int_t p,Float_t v) { if(p < 10)fAsymCuts[p] = v    ; }
105   void         SetNCellCutsAt(Int_t p,Int_t v)  { if(p < 10)fCellNCuts[p]= v    ; }
106   void         SetPIDBitsAt  (Int_t p,Int_t v)  { if(p < 10)fPIDBits[p]  = v    ; }
107   
108   void         SwitchOnFillSSCombinations()     { fFillSSCombinations  = kTRUE  ; }
109   void         SwitchOffFillSSCombinations()    { fFillSSCombinations  = kFALSE ; }
110   
111   void         SwitchOnFillAsymmetryHisto()     { fFillAsymmetryHisto  = kTRUE  ; }
112   void         SwitchOffFillAsymmetryHisto()    { fFillAsymmetryHisto  = kFALSE ; }
113
114   void         SwitchOnFillOriginHisto()        { fFillOriginHisto     = kTRUE  ; }
115   void         SwitchOffFillOriginHisto()       { fFillOriginHisto     = kFALSE ; }
116   
117   //MC analysis related methods
118     
119   void         SwitchOnConversionChecker()      { fCheckConversion     = kTRUE  ; }
120   void         SwitchOffConversionChecker()     { fCheckConversion     = kFALSE ; }  
121   
122   void         SwitchOnMultipleCutAnalysisInSimulation()  { fMultiCutAnaSim = kTRUE  ; }
123   void         SwitchOffMultipleCutAnalysisInSimulation() { fMultiCutAnaSim = kFALSE ; }
124   
125   void         FillAcceptanceHistograms();
126   void         FillMCVersusRecDataHistograms(const Int_t    index1,  const Int_t    index2,
127                                              const Float_t  pt1,     const Float_t  pt2, 
128                                              const Int_t    ncells1, const Int_t    ncells2, 
129                                              const Double_t mass,    const Double_t pt,  const Double_t asym,    
130                                              const Double_t deta,    const Double_t dphi);
131   
132   private:
133
134   TList ** fEventsList ;               //![GetNCentrBin()*GetNZvertBin()*GetNRPBin()] Containers for photons in stored events
135
136   TString  fCalorimeter ;              // Select Calorimeter for IM
137   Int_t    fNModules ;                 // Number of EMCAL/PHOS modules, set as many histogras as modules 
138   
139   Bool_t   fUseAngleCut ;              // Select pairs depending on their opening angle
140   Bool_t   fUseAngleEDepCut ;          // Select pairs depending on their opening angle
141   Float_t  fAngleCut ;                 // Select pairs with opening angle larger than a threshold
142   Float_t  fAngleMaxCut ;              // Select pairs with opening angle smaller than a threshold
143   
144   //Multiple cuts analysis
145   Bool_t   fMultiCutAna;               // Do analysis with several or fixed cut
146   Bool_t   fMultiCutAnaSim;            // Do analysis with several or fixed cut, in the simulation related part
147   Int_t    fNPtCuts;                   // Number of pt cuts
148   Float_t  fPtCuts[10];                // Array with different pt cuts
149   Int_t    fNAsymCuts;                 // Number of assymmetry cuts
150   Float_t  fAsymCuts[10];              // Array with different assymetry cuts
151   Int_t    fNCellNCuts;                // Number of cuts with number of cells in cluster
152   Int_t    fCellNCuts[10];             // Array with different cell number cluster cuts
153   Int_t    fNPIDBits ;                       // Number of possible PID bit combinations
154   Int_t    fPIDBits[10];               // Array with different PID bits
155   
156   //Switchs of different analysis options
157   Bool_t   fMakeInvPtPlots;            // D plots with inverse pt weight
158   Bool_t   fSameSM;                    // Select only pairs in same SM;
159   Bool_t   fFillSMCombinations;        // Fill histograms with different cluster pairs in SM combinations
160   Bool_t   fCheckConversion;           // Fill histograms with tagged photons as conversion
161   Bool_t   fFillBadDistHisto;          // Do plots for different distances to bad channels
162   Bool_t   fFillSSCombinations;        // Do invariant mass for different combination of shower shape clusters
163   Bool_t   fFillAngleHisto;            // Fill histograms with pair opening angle
164   Bool_t   fFillAsymmetryHisto;        // Fill histograms with asymmetry vs pt
165   Bool_t   fFillOriginHisto;           // Fill histograms depending on their origin
166
167   //Histograms
168   
169   //Event characterization
170   TH1F *   fhAverTotECluster;          //! Average number of clusters in SM
171   TH1F *   fhAverTotECell;             //! Average number of cells    in SM
172   TH2F *   fhAverTotECellvsCluster;    //! Average number of cells    in SM
173   TH1F *   fhEDensityCluster;          //! Deposited energy in event per cluster
174   TH1F *   fhEDensityCell;             //! Deposited energy in event per cell vs cluster
175   TH2F *   fhEDensityCellvsCluster;    //! Deposited energy in event per cell vs cluster
176
177   TH2F **  fhReMod ;                   //![fNModules]   REAL  two-photon invariant mass distribution for different calorimeter modules.
178   TH2F **  fhReSameSideEMCALMod ;      //![fNModules-2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
179   TH2F **  fhReSameSectorEMCALMod ;    //![fNModules/2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
180   TH2F **  fhReDiffPHOSMod ;           //![fNModules]   REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
181   TH2F **  fhMiMod ;                   //![fNModules]   MIXED two-photon invariant mass distribution for different calorimeter modules.
182   TH2F **  fhMiSameSideEMCALMod ;      //![fNModules-2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
183   TH2F **  fhMiSameSectorEMCALMod ;    //![fNModules/2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
184   TH2F **  fhMiDiffPHOSMod ;           //![fNModules-1] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
185   
186   // Pairs with at least one cluster tagged as conversion
187   TH2F *   fhReConv ;                  //! REAL  two-photon invariant mass distribution one of the pair was 2 clusters with small mass 
188   TH2F *   fhMiConv ;                  //! MIXED two-photon invariant mass distribution one of the pair was 2 clusters with small mass
189   TH2F *   fhReConv2 ;                 //! REAL  two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass 
190   TH2F *   fhMiConv2 ;                 //! MIXED two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass
191
192   TH2F **  fhRe1 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
193   TH2F **  fhMi1 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
194   TH2F **  fhRe2 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
195   TH2F **  fhMi2 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
196   TH2F **  fhRe3 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
197   TH2F **  fhMi3 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
198
199   //Histograms weighted by inverse pT
200   TH2F **  fhReInvPt1 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
201   TH2F **  fhMiInvPt1 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
202   TH2F **  fhReInvPt2 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT 
203   TH2F **  fhMiInvPt2 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
204   TH2F **  fhReInvPt3 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
205   TH2F **  fhMiInvPt3 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
206   
207   //Multiple cuts: Assymmetry, pt, n cells, PID
208   TH2F **  fhRePtNCellAsymCuts ;       //![fNPtCuts*fNAsymCuts*fNCellNCuts*] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
209   TH2F **  fhMiPtNCellAsymCuts ;       //![fNPtCuts*fNAsymCuts*fNCellNCuts] Mixed two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
210   TH2F **  fhRePtNCellAsymCutsSM[12] ; //![fNPtCuts*fNAsymCuts*fNCellNCutsfNModules] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry for each module
211  
212   TH2F **  fhRePIDBits ;               //![fNPIDBits]  REAL two-photon invariant mass distribution for different PID bits
213   TH3F **  fhRePtMult ;                //![fNAsymCuts] REAL two-photon invariant mass distribution for different track multiplicity and assymetry cuts
214   TH2F *   fhReSS[3] ;                 //! Combine clusters with 3 different cuts on shower shape
215   
216   // Asymmetry vs pt, in pi0/eta regions
217   TH2F *   fhRePtAsym    ;             //! REAL two-photon pt vs asymmetry
218   TH2F *   fhRePtAsymPi0 ;             //! REAL two-photon pt vs asymmetry, close to pi0 mass
219   TH2F *   fhRePtAsymEta ;             //! REAL two-photon pt vs asymmetry, close to eta mass
220   
221   //Centrality, Event plane bins
222   TH1I *   fhEventBin;                 //! Number of real  pairs in a particular bin (cen,vz,rp)
223   TH1I *   fhEventMixBin;              //! Number of mixed pairs in a particular bin (cen,vz,rp)
224   TH1F *   fhCentrality;               //! Histogram with centrality bins with at least one pare
225   TH1F *   fhCentralityNoPair;         //! Histogram with centrality bins with no pair
226
227   TH2F *   fhEventPlaneResolution;     //! Histogram with Event plane resolution vs centrality
228   
229   // Pair opening angle
230   TH2F *   fhRealOpeningAngle ;        //! Opening angle of pair versus pair energy
231   TH2F *   fhRealCosOpeningAngle ;     //! Cosinus of opening angle of pair version pair energy
232   TH2F *   fhMixedOpeningAngle ;       //! Opening angle of pair versus pair energy
233   TH2F *   fhMixedCosOpeningAngle ;    //! Cosinus of opening angle of pair version pair energy
234   
235   //MC analysis histograms
236   //Pi0 Acceptance
237   TH1F *   fhPrimPi0E ;                //! Spectrum of Primary
238   TH1F *   fhPrimPi0Pt ;               //! Spectrum of Primary
239   TH1F *   fhPrimPi0AccE ;             //! Spectrum of primary with accepted daughters
240   TH1F *   fhPrimPi0AccPt ;            //! Spectrum of primary with accepted daughters
241   TH2F *   fhPrimPi0Y ;                //! Rapidity distribution of primary particles  vs pT
242   TH2F *   fhPrimPi0AccY ;             //! Rapidity distribution of primary with accepted daughters  vs pT
243   TH2F *   fhPrimPi0Phi ;              //! Azimutal distribution of primary particles  vs pT
244   TH2F *   fhPrimPi0AccPhi;            //! Azimutal distribution of primary with accepted daughters  vs pT
245   TH2F *   fhPrimPi0OpeningAngle ;     //! Opening angle of pair versus pair energy, primaries
246   TH2F *   fhPrimPi0OpeningAngleAsym ; //! Opening angle of pair versus pair E asymmetry, pi0 primaries
247   TH2F *   fhPrimPi0CosOpeningAngle ;  //! Cosinus of opening angle of pair version pair energy, pi0 primaries
248   TH2F *   fhPrimPi0PtCentrality ;     //! primary pi0 reconstructed centrality  vs pT
249   TH2F *   fhPrimPi0PtEventPlane ;     //! primary pi0 reconstructed event plane vs pT
250   TH2F *   fhPrimPi0AccPtCentrality ;  //! primary pi0 with accepted daughters reconstructed centrality  vs pT
251   TH2F *   fhPrimPi0AccPtEventPlane ;  //! primary pi0 with accepted daughters reconstructed event plane vs pT
252
253   //Eta acceptance
254   TH1F *   fhPrimEtaE ;                //! Spectrum of Primary
255   TH1F *   fhPrimEtaPt ;               //! Spectrum of Primary
256   TH1F *   fhPrimEtaAccE ;             //! Spectrum of primary with accepted daughters
257   TH1F *   fhPrimEtaAccPt ;            //! Spectrum of primary with accepted daughters
258   TH2F *   fhPrimEtaY ;                //! Rapidity distribution of primary particles vs pT
259   TH2F *   fhPrimEtaAccY ;             //! Rapidity distribution of primary with accepted daughters  vs pT
260   TH2F *   fhPrimEtaPhi ;              //! Azimutal distribution of primary particles  vs pT
261   TH2F *   fhPrimEtaAccPhi;            //! Azimutal distribution of primary with accepted daughters      vs pT
262   TH2F *   fhPrimEtaOpeningAngle ;     //! Opening angle of pair versus pair energy, eta primaries
263   TH2F *   fhPrimEtaOpeningAngleAsym ; //! Opening angle of pair versus pair E asymmetry, eta primaries
264   TH2F *   fhPrimEtaCosOpeningAngle ;  //! Cosinus of opening angle of pair version pair energy, eta primaries
265   TH2F *   fhPrimEtaPtCentrality ;     //! primary eta reconstructed centrality  vs pT
266   TH2F *   fhPrimEtaPtEventPlane ;     //! primary eta reconstructed event plane vs pT
267   TH2F *   fhPrimEtaAccPtCentrality ;  //! primary eta with accepted daughters reconstructed centrality  vs pT
268   TH2F *   fhPrimEtaAccPtEventPlane ;  //! primary eta with accepted daughters reconstructed event plane vs pT
269   
270   // Primaries origin
271   TH2F *   fhPrimPi0PtOrigin ;         //! Spectrum of generated pi0 vs mother
272   TH2F *   fhPrimEtaPtOrigin ;         //! Spectrum of generated eta vs mother
273   
274   //Pair origin
275   //Array of histograms ordered as follows: 0-Photon, 1-electron, 2-pi0, 3-eta, 4-a-proton, 5-a-neutron, 6-stable particles, 
276   // 7-other decays, 8-string, 9-final parton, 10-initial parton, intermediate, 11-colliding proton, 12-unrelated
277   TH2F *   fhMCOrgMass[13];            //! Mass vs pt of real pairs, check common origin of pair
278   TH2F *   fhMCOrgAsym[13];            //! Asymmetry vs pt of real pairs, check common origin of pair
279   TH2F *   fhMCOrgDeltaEta[13];        //! Delta Eta vs pt of real pairs, check common origin of pair
280   TH2F *   fhMCOrgDeltaPhi[13];        //! Delta Phi vs pt of real pairs, check common origin of pair
281   
282   //Multiple cuts in simulation, origin pi0 or eta
283   TH2F **  fhMCPi0MassPtRec;           //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs reconstructed pt of original pair  
284   TH2F **  fhMCPi0MassPtTrue;          //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs generated pt of original pair  
285   TH2F **  fhMCPi0PtTruePtRec;         //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed pt vs generated pt of pair
286   TH2F **  fhMCEtaMassPtRec;           //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs reconstructed pt of original pair  
287   TH2F **  fhMCEtaMassPtTrue;          //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs generated pt of original pair  
288   TH2F **  fhMCEtaPtTruePtRec;         //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed pt vs generated pt of pair
289
290   TH2F *   fhMCPi0PtOrigin ;           //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
291   TH2F *   fhMCEtaPtOrigin ;           //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
292
293   TH2F *   fhReMCFromConversion ;      //! Invariant mass of 2 clusters originated in conversions
294   TH2F *   fhReMCFromNotConversion ;   //! Invariant mass of 2 clusters not originated in conversions
295   TH2F *   fhReMCFromMixConversion ;   //! Invariant mass of 2 clusters one from conversion and the other not
296
297   AliAnaPi0(              const AliAnaPi0 & api0) ; // cpy ctor
298   AliAnaPi0 & operator = (const AliAnaPi0 & api0) ; // cpy assignment
299   
300   ClassDef(AliAnaPi0,24)
301 } ;
302
303
304 #endif //ALIANAPI0_H
305
306
307