]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/CaloTrackCorrelations/AliAnaPi0.h
Merge branch 'feature-movesplit'
[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   //-------------------------------
49   // EVENT Bin Methods
50   //-------------------------------
51
52   Int_t        GetEventIndex(AliAODPWG4Particle * part, Double_t * vert)  ;  
53
54   //-------------------------------
55         //Opening angle pair selection
56   //-------------------------------
57   void         SwitchOnAngleSelection()         { fUseAngleCut         = kTRUE  ; }
58   void         SwitchOffAngleSelection()        { fUseAngleCut         = kFALSE ; }
59   
60   void         SwitchOnAngleEDepSelection()     { fUseAngleEDepCut     = kTRUE  ; }
61   void         SwitchOffAngleEDepSelection()    { fUseAngleEDepCut     = kFALSE ; }
62     
63   void         SetAngleCut(Float_t a)           { fAngleCut            = a      ; }
64   void         SetAngleMaxCut(Float_t a)        { fAngleMaxCut         = a      ; }
65
66   void         SwitchOnFillAngleHisto()         { fFillAngleHisto      = kTRUE  ; }
67   void         SwitchOffFillAngleHisto()        { fFillAngleHisto      = kFALSE ; }
68   
69   //------------------------------------------
70   //Do analysis only with clusters in same SM or different combinations of SM
71   //------------------------------------------
72   void         SwitchOnSameSM()                 { fSameSM              = kTRUE  ; }
73   void         SwitchOffSameSM()                { fSameSM              = kFALSE ; }
74   
75   void         SwitchOnSMCombinations()         { fFillSMCombinations  = kTRUE  ; }
76   void         SwitchOffSMCombinations()        { fFillSMCombinations  = kFALSE ; }
77   
78   //-------------------------------
79   //Histogram filling options off by default
80   //-------------------------------
81   void         SwitchOnInvPtWeight()            { fMakeInvPtPlots      = kTRUE  ; }
82   void         SwitchOffInvPtWeight()           { fMakeInvPtPlots      = kFALSE ; }
83   
84   void         SwitchOnFillBadDistHisto()       { fFillBadDistHisto    = kTRUE  ; }
85   void         SwitchOffFillBadDistHisto()      { fFillBadDistHisto    = kFALSE ; }
86   
87   //-------------------------------------------
88   //Cuts for multiple analysis, off by default
89   //-------------------------------------------
90   void         SwitchOnMultipleCutAnalysis()    { fMultiCutAna         = kTRUE  ; }
91   void         SwitchOffMultipleCutAnalysis()   { fMultiCutAna         = kFALSE ; }
92
93   void         SetNPtCuts   (Int_t s)           { if(s <= 10)fNPtCuts    = s    ; }
94   void         SetNAsymCuts (Int_t s)           { if(s <= 10)fNAsymCuts  = s    ; }
95   void         SetNNCellCuts(Int_t s)           { if(s <= 10)fNCellNCuts = s    ; }
96   void         SetNPIDBits  (Int_t s)           { if(s <= 10)fNPIDBits   = s    ; }
97   
98   void         SetPtCutsAt  (Int_t p,Float_t v) { if(p < 10)fPtCuts[p]   = v    ; }
99   void         SetAsymCutsAt(Int_t p,Float_t v) { if(p < 10)fAsymCuts[p] = v    ; }
100   void         SetNCellCutsAt(Int_t p,Int_t v)  { if(p < 10)fCellNCuts[p]= v    ; }
101   void         SetPIDBitsAt  (Int_t p,Int_t v)  { if(p < 10)fPIDBits[p]  = v    ; }
102   
103   void         SwitchOnFillSSCombinations()     { fFillSSCombinations  = kTRUE  ; }
104   void         SwitchOffFillSSCombinations()    { fFillSSCombinations  = kFALSE ; }
105   
106   void         SwitchOnFillAsymmetryHisto()     { fFillAsymmetryHisto  = kTRUE  ; }
107   void         SwitchOffFillAsymmetryHisto()    { fFillAsymmetryHisto  = kFALSE ; }
108
109   void         SwitchOnFillOriginHisto()        { fFillOriginHisto     = kTRUE  ; }
110   void         SwitchOffFillOriginHisto()       { fFillOriginHisto     = kFALSE ; }
111
112   void         SwitchOnFillArmenterosThetaStarHisto()  { fFillArmenterosThetaStar = kTRUE  ; }
113   void         SwitchOffFillArmenterosThetaStarHisto() { fFillArmenterosThetaStar = kFALSE ; }
114   
115   //MC analysis related methods
116     
117   void         SwitchOnConversionChecker()      { fCheckConversion     = kTRUE  ; }
118   void         SwitchOffConversionChecker()     { fCheckConversion     = kFALSE ; }  
119   
120   void         SwitchOnMultipleCutAnalysisInSimulation()  { fMultiCutAnaSim = kTRUE  ; }
121   void         SwitchOffMultipleCutAnalysisInSimulation() { fMultiCutAnaSim = kFALSE ; }
122   
123   void         SwitchOnCheckAcceptanceInSector() { fCheckAccInSector   = kTRUE  ; }
124   void         SwitchOffCheckAcceptanceInSector(){ fCheckAccInSector   = kFALSE ; }
125   
126   void         FillAcceptanceHistograms();
127   void         FillMCVersusRecDataHistograms(Int_t    index1,  Int_t    index2,
128                                              Float_t  pt1,     Float_t  pt2,
129                                              Int_t    ncells1, Int_t    ncells2,
130                                              Double_t mass,    Double_t pt,     Double_t asym,
131                                              Double_t deta,    Double_t dphi);
132   
133   void         FillArmenterosThetaStar(Int_t pdg);
134
135   
136   private:
137
138   TList ** fEventsList ;               //![GetNCentrBin()*GetNZvertBin()*GetNRPBin()] Containers for photons in stored events
139
140   Int_t    fNModules ;                 // Number of EMCAL/PHOS modules, set as many histogras as modules 
141   
142   Bool_t   fUseAngleCut ;              // Select pairs depending on their opening angle
143   Bool_t   fUseAngleEDepCut ;          // Select pairs depending on their opening angle
144   Float_t  fAngleCut ;                 // Select pairs with opening angle larger than a threshold
145   Float_t  fAngleMaxCut ;              // Select pairs with opening angle smaller than a threshold
146   
147   //Multiple cuts analysis
148   Bool_t   fMultiCutAna;               // Do analysis with several or fixed cut
149   Bool_t   fMultiCutAnaSim;            // Do analysis with several or fixed cut, in the simulation related part
150   Int_t    fNPtCuts;                   // Number of pt cuts
151   Float_t  fPtCuts[10];                // Array with different pt cuts
152   Int_t    fNAsymCuts;                 // Number of assymmetry cuts
153   Float_t  fAsymCuts[10];              // Array with different assymetry cuts
154   Int_t    fNCellNCuts;                // Number of cuts with number of cells in cluster
155   Int_t    fCellNCuts[10];             // Array with different cell number cluster cuts
156   Int_t    fNPIDBits ;                       // Number of possible PID bit combinations
157   Int_t    fPIDBits[10];               // Array with different PID bits
158   
159   //Switchs of different analysis options
160   Bool_t   fMakeInvPtPlots;            // D plots with inverse pt weight
161   Bool_t   fSameSM;                    // Select only pairs in same SM;
162   Bool_t   fFillSMCombinations;        // Fill histograms with different cluster pairs in SM combinations
163   Bool_t   fCheckConversion;           // Fill histograms with tagged photons as conversion
164   Bool_t   fFillBadDistHisto;          // Do plots for different distances to bad channels
165   Bool_t   fFillSSCombinations;        // Do invariant mass for different combination of shower shape clusters
166   Bool_t   fFillAngleHisto;            // Fill histograms with pair opening angle
167   Bool_t   fFillAsymmetryHisto;        // Fill histograms with asymmetry vs pt
168   Bool_t   fFillOriginHisto;           // Fill histograms depending on their origin
169   Bool_t   fFillArmenterosThetaStar;   // Fill armenteros histograms
170   
171   Bool_t   fCheckAccInSector;          // Check that the decay pi0 falls in the same SM or sector
172   
173   TLorentzVector fPhotonMom1;          //! photon cluster momentum
174   TLorentzVector fPhotonMom1Boost;     //! photon cluster momentum
175   TLorentzVector fPhotonMom2;          //! photon cluster momentum
176   TLorentzVector fPi0Mom;              //! pi0 cluster momentum
177   TVector3       fProdVertex;          //! production vertex
178   
179   //Histograms
180   
181   //Event characterization
182   TH1F *   fhAverTotECluster;          //! Average number of clusters in SM
183   TH1F *   fhAverTotECell;             //! Average number of cells    in SM
184   TH2F *   fhAverTotECellvsCluster;    //! Average number of cells    in SM
185   TH1F *   fhEDensityCluster;          //! Deposited energy in event per cluster
186   TH1F *   fhEDensityCell;             //! Deposited energy in event per cell vs cluster
187   TH2F *   fhEDensityCellvsCluster;    //! Deposited energy in event per cell vs cluster
188
189   TH2F **  fhReMod ;                   //![fNModules]   REAL  two-photon invariant mass distribution for different calorimeter modules.
190   TH2F **  fhReSameSideEMCALMod ;      //![fNModules-2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
191   TH2F **  fhReSameSectorEMCALMod ;    //![fNModules/2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
192   TH2F **  fhReDiffPHOSMod ;           //![fNModules]   REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
193   TH2F **  fhMiMod ;                   //![fNModules]   MIXED two-photon invariant mass distribution for different calorimeter modules.
194   TH2F **  fhMiSameSideEMCALMod ;      //![fNModules-2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
195   TH2F **  fhMiSameSectorEMCALMod ;    //![fNModules/2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
196   TH2F **  fhMiDiffPHOSMod ;           //![fNModules-1] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
197   
198   // Pairs with at least one cluster tagged as conversion
199   TH2F *   fhReConv ;                  //! REAL  two-photon invariant mass distribution one of the pair was 2 clusters with small mass 
200   TH2F *   fhMiConv ;                  //! MIXED two-photon invariant mass distribution one of the pair was 2 clusters with small mass
201   TH2F *   fhReConv2 ;                 //! REAL  two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass 
202   TH2F *   fhMiConv2 ;                 //! MIXED two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass
203
204   TH2F **  fhRe1 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
205   TH2F **  fhMi1 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
206   TH2F **  fhRe2 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
207   TH2F **  fhMi2 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
208   TH2F **  fhRe3 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
209   TH2F **  fhMi3 ;                     //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
210
211   //Histograms weighted by inverse pT
212   TH2F **  fhReInvPt1 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
213   TH2F **  fhMiInvPt1 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
214   TH2F **  fhReInvPt2 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT 
215   TH2F **  fhMiInvPt2 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
216   TH2F **  fhReInvPt3 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
217   TH2F **  fhMiInvPt3 ;                //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
218   
219   //Multiple cuts: Assymmetry, pt, n cells, PID
220   TH2F **  fhRePtNCellAsymCuts ;       //![fNPtCuts*fNAsymCuts*fNCellNCuts*] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
221   TH2F **  fhMiPtNCellAsymCuts ;       //![fNPtCuts*fNAsymCuts*fNCellNCuts] Mixed two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
222   TH2F **  fhRePtNCellAsymCutsSM[12] ; //![fNPtCuts*fNAsymCuts*fNCellNCutsfNModules] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry for each module
223  
224   TH2F **  fhRePIDBits ;               //![fNPIDBits]  REAL two-photon invariant mass distribution for different PID bits
225   TH3F **  fhRePtMult ;                //![fNAsymCuts] REAL two-photon invariant mass distribution for different track multiplicity and assymetry cuts
226   TH2F *   fhReSS[3] ;                 //! Combine clusters with 3 different cuts on shower shape
227   
228   // Asymmetry vs pt, in pi0/eta regions
229   TH2F *   fhRePtAsym    ;             //! REAL two-photon pt vs asymmetry
230   TH2F *   fhRePtAsymPi0 ;             //! REAL two-photon pt vs asymmetry, close to pi0 mass
231   TH2F *   fhRePtAsymEta ;             //! REAL two-photon pt vs asymmetry, close to eta mass
232   
233   //Centrality, Event plane bins
234   TH1I *   fhEventBin;                 //! Number of real  pairs in a particular bin (cen,vz,rp)
235   TH1I *   fhEventMixBin;              //! Number of mixed pairs in a particular bin (cen,vz,rp)
236   TH1F *   fhCentrality;               //! Histogram with centrality bins with at least one pare
237   TH1F *   fhCentralityNoPair;         //! Histogram with centrality bins with no pair
238
239   TH2F *   fhEventPlaneResolution;     //! Histogram with Event plane resolution vs centrality
240   
241   // Pair opening angle
242   TH2F *   fhRealOpeningAngle ;        //! Opening angle of pair versus pair energy
243   TH2F *   fhRealCosOpeningAngle ;     //! Cosinus of opening angle of pair version pair energy
244   TH2F *   fhMixedOpeningAngle ;       //! Opening angle of pair versus pair energy
245   TH2F *   fhMixedCosOpeningAngle ;    //! Cosinus of opening angle of pair version pair energy
246   
247   //MC analysis histograms
248   //Pi0 Acceptance
249   TH1F *   fhPrimPi0E ;                //! Spectrum of Primary
250   TH1F *   fhPrimPi0Pt ;               //! Spectrum of Primary
251   TH1F *   fhPrimPi0AccE ;             //! Spectrum of primary with accepted daughters
252   TH1F *   fhPrimPi0AccPt ;            //! Spectrum of primary with accepted daughters
253   TH2F *   fhPrimPi0Y ;                //! Rapidity distribution of primary particles  vs pT
254   TH2F *   fhPrimPi0AccY ;             //! Rapidity distribution of primary with accepted daughters  vs pT
255   TH2F *   fhPrimPi0Yeta ;             //! PseudoRapidity distribution of primary particles  vs pT
256   TH2F *   fhPrimPi0YetaYcut ;         //! PseudoRapidity distribution of primary particles  vs pT, Y<1
257   TH2F *   fhPrimPi0AccYeta ;          //! PseudoRapidity distribution of primary with accepted daughters  vs pT
258   TH2F *   fhPrimPi0Phi ;              //! Azimutal distribution of primary particles  vs pT
259   TH2F *   fhPrimPi0AccPhi;            //! Azimutal distribution of primary with accepted daughters  vs pT
260   TH2F *   fhPrimPi0OpeningAngle ;     //! Opening angle of pair versus pair energy, primaries
261   TH2F *   fhPrimPi0OpeningAngleAsym ; //! Opening angle of pair versus pair E asymmetry, pi0 primaries
262   TH2F *   fhPrimPi0CosOpeningAngle ;  //! Cosinus of opening angle of pair version pair energy, pi0 primaries
263   TH2F *   fhPrimPi0PtCentrality ;     //! primary pi0 reconstructed centrality  vs pT
264   TH2F *   fhPrimPi0PtEventPlane ;     //! primary pi0 reconstructed event plane vs pT
265   TH2F *   fhPrimPi0AccPtCentrality ;  //! primary pi0 with accepted daughters reconstructed centrality  vs pT
266   TH2F *   fhPrimPi0AccPtEventPlane ;  //! primary pi0 with accepted daughters reconstructed event plane vs pT
267
268   //Eta acceptance
269   TH1F *   fhPrimEtaE ;                //! Spectrum of Primary
270   TH1F *   fhPrimEtaPt ;               //! Spectrum of Primary
271   TH1F *   fhPrimEtaAccE ;             //! Spectrum of primary with accepted daughters
272   TH1F *   fhPrimEtaAccPt ;            //! Spectrum of primary with accepted daughters
273   TH2F *   fhPrimEtaY ;                //! Rapidity distribution of primary particles vs pT
274   TH2F *   fhPrimEtaAccY ;             //! Rapidity distribution of primary with accepted daughters  vs pT
275   TH2F *   fhPrimEtaYeta ;             //! PseudoRapidity distribution of primary particles vs pT
276   TH2F *   fhPrimEtaYetaYcut ;         //! PseudoRapidity distribution of primary particles vs pT, Y<1
277   TH2F *   fhPrimEtaAccYeta ;          //! PseudoRapidity distribution of primary with accepted daughters  vs pT
278   TH2F *   fhPrimEtaPhi ;              //! Azimutal distribution of primary particles  vs pT
279   TH2F *   fhPrimEtaAccPhi;            //! Azimutal distribution of primary with accepted daughters      vs pT
280   TH2F *   fhPrimEtaOpeningAngle ;     //! Opening angle of pair versus pair energy, eta primaries
281   TH2F *   fhPrimEtaOpeningAngleAsym ; //! Opening angle of pair versus pair E asymmetry, eta primaries
282   TH2F *   fhPrimEtaCosOpeningAngle ;  //! Cosinus of opening angle of pair version pair energy, eta primaries
283   TH2F *   fhPrimEtaPtCentrality ;     //! primary eta reconstructed centrality  vs pT
284   TH2F *   fhPrimEtaPtEventPlane ;     //! primary eta reconstructed event plane vs pT
285   TH2F *   fhPrimEtaAccPtCentrality ;  //! primary eta with accepted daughters reconstructed centrality  vs pT
286   TH2F *   fhPrimEtaAccPtEventPlane ;  //! primary eta with accepted daughters reconstructed event plane vs pT
287   
288   // Primaries origin
289   TH2F *   fhPrimPi0PtOrigin ;         //! Spectrum of generated pi0 vs mother
290   TH2F *   fhPrimEtaPtOrigin ;         //! Spectrum of generated eta vs mother
291   
292   //Pair origin
293   //Array of histograms ordered as follows: 0-Photon, 1-electron, 2-pi0, 3-eta, 4-a-proton, 5-a-neutron, 6-stable particles, 
294   // 7-other decays, 8-string, 9-final parton, 10-initial parton, intermediate, 11-colliding proton, 12-unrelated
295   TH2F *   fhMCOrgMass[13];            //! Mass vs pt of real pairs, check common origin of pair
296   TH2F *   fhMCOrgAsym[13];            //! Asymmetry vs pt of real pairs, check common origin of pair
297   TH2F *   fhMCOrgDeltaEta[13];        //! Delta Eta vs pt of real pairs, check common origin of pair
298   TH2F *   fhMCOrgDeltaPhi[13];        //! Delta Phi vs pt of real pairs, check common origin of pair
299   
300   //Multiple cuts in simulation, origin pi0 or eta
301   TH2F **  fhMCPi0MassPtRec;           //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs reconstructed pt of original pair  
302   TH2F **  fhMCPi0MassPtTrue;          //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs generated pt of original pair  
303   TH2F **  fhMCPi0PtTruePtRec;         //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed pt vs generated pt of pair
304   TH2F **  fhMCEtaMassPtRec;           //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs reconstructed pt of original pair  
305   TH2F **  fhMCEtaMassPtTrue;          //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs generated pt of original pair  
306   TH2F **  fhMCEtaPtTruePtRec;         //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed pt vs generated pt of pair
307
308   TH2F *   fhMCPi0PtOrigin ;           //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
309   TH2F *   fhMCEtaPtOrigin ;           //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
310
311   TH2F *   fhMCPi0ProdVertex;          //! Spectrum of selected pi0 vs production vertex
312   TH2F *   fhMCEtaProdVertex;          //! Spectrum of selected eta vs production vertex
313   TH2F *   fhPrimPi0ProdVertex;        //! Spectrum of primary pi0 vs production vertex
314   TH2F *   fhPrimEtaProdVertex;        //! Spectrum of primary eta vs production vertex
315   
316   TH2F *   fhReMCFromConversion ;      //! Invariant mass of 2 clusters originated in conversions
317   TH2F *   fhReMCFromNotConversion ;   //! Invariant mass of 2 clusters not originated in conversions
318   TH2F *   fhReMCFromMixConversion ;   //! Invariant mass of 2 clusters one from conversion and the other not
319
320   TH2F *   fhArmPrimPi0[4];            //! Armenteros plots for primary pi0 in 6 energy bins
321   TH2F *   fhArmPrimEta[4];            //! Armenteros plots for primary eta in 6 energy bins
322   TH2F *   fhCosThStarPrimPi0;         //! cos(theta*) plots vs E for primary pi0, same as asymmetry ...
323   TH2F *   fhCosThStarPrimEta;         //! cos(theta*) plots vs E for primary eta, same as asymmetry ...
324   
325   TH2F *   fhEPairDiffTime;            //! E pair vs Pair of clusters time difference vs E
326   
327   AliAnaPi0(              const AliAnaPi0 & api0) ; // cpy ctor
328   AliAnaPi0 & operator = (const AliAnaPi0 & api0) ; // cpy assignment
329   
330   ClassDef(AliAnaPi0,29)
331 } ;
332
333
334 #endif //ALIANAPI0_H
335
336
337