]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrDep/AliAnaPi0.h
1bcbcf088c620720e9dac8ee87512169ede2ffe9
[u/mrichter/AliRoot.git] / PWG4 / PartCorrDep / 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 /* $Id: $ */
6
7 //_________________________________________________________________________
8 // Class to fill two-photon invariant mass histograms 
9 // to be used to extract pi0 raw yield.
10 // Input is produced by AliAnaPhoton (or any other analysis producing output AliAODPWG4Particles), 
11 // it will do nothing if executed alone
12 //
13 //-- Author: Dmitri Peressounko (RRC "KI")
14 //-- Adapted to PartCorr frame by Lamia Benhabib (SUBATECH)
15 //-- and Gustavo Conesa (INFN-Frascati)
16
17 //Root
18 class TList;
19 class TH3D ;
20 class TH2D ;
21 class TObjString;
22
23 //Analysis
24 #include "AliAnaPartCorrBaseClass.h"
25 class AliAODEvent ;
26 class AliESDEvent ;
27 class AliAODPWG4Particle ;
28
29 class AliAnaPi0 : public AliAnaPartCorrBaseClass {
30   
31  public:   
32   AliAnaPi0() ; // default ctor
33   virtual ~AliAnaPi0() ;//virtual dtor
34  private:
35   AliAnaPi0(const AliAnaPi0 & g) ; // cpy ctor
36   AliAnaPi0 & operator = (const AliAnaPi0 & api0) ;//cpy assignment
37   
38  public:
39
40   //-------------------------------
41   // General analysis frame methods
42   //-------------------------------
43
44   TObjString * GetAnalysisCuts();
45   TList      * GetCreateOutputObjects(); 
46   void Terminate(TList* outputList);
47   void ReadHistograms(TList * outputList); //Fill histograms with histograms in ouput list, needed in Terminate.
48   void Print(const Option_t * opt) const;  
49   //void MakeAnalysisFillAOD() {;} //Not needed
50   void MakeAnalysisFillHistograms();
51   //void Init();
52   void InitParameters();
53
54   //Calorimeter options
55   TString GetCalorimeter()        const { return fCalorimeter; }
56   void SetCalorimeter(TString & det)    { fCalorimeter = det ; }
57   void SetNumberOfModules(Int_t nmod)   { fNModules    = nmod; }
58   
59   //-------------------------------
60   // EVENT Bin Methods
61   //-------------------------------
62
63   virtual Int_t GetEventIndex(AliAODPWG4Particle * part, Double_t * vert)  ;
64
65   void CountAndGetAverages(Int_t &nClus,Int_t &nCell, Float_t &eClusTot,Float_t &eCellTot, Float_t &eDenClus,Float_t &eDenCell) ;
66   
67   //Setters for parameters of event buffers
68   void SetNCentrBin(Int_t n=5)    {fNCentrBin=n ;} //number of bins in centrality 
69 //void SetNZvertBin(Int_t n=5)    {fNZvertBin=n ;} //number of bins for vertex position
70 //void SetNRPBin(Int_t n=6)       {fNrpBin=n ;}    //number of bins in reaction plain
71   void SetNMaxEvMix(Int_t n=20)   {fNmaxMixEv=n ;} //Maximal number of events for mixing
72   
73   //Switchs for event multiplicity bin option, by default, centrality
74   void SwitchOnTrackMultBins()    {fUseTrackMultBins = kTRUE  ; }
75   void SwitchOffTrackMultBins()   {fUseTrackMultBins = kFALSE ; }
76   
77   void SwitchOnPhotonMultBins()   {fUsePhotonMultBins = kTRUE  ; }
78   void SwitchOffPhotonMultBins()  {fUsePhotonMultBins = kFALSE ; }
79   
80   void SwitchOnClusterEBins()     {fUseAverClusterEBins = kTRUE  ; }
81   void SwitchOffClusterEBins()    {fUseAverClusterEBins = kFALSE ; }
82   
83   void SwitchOnCellEBins()        {fUseAverCellEBins = kTRUE  ; }
84   void SwitchOffCellEBins()       {fUseAverCellEBins = kFALSE ; }
85
86   void SwitchOnClusterEDenBins()  {fUseAverClusterEDenBins = kTRUE  ; }
87   void SwitchOffClusterEDenBins() {fUseAverClusterEDenBins = kFALSE ; }
88   
89 //  void SwitchOnClusterPairRBins()     {fUseAverClusterPairRBins = kTRUE  ; }
90 //  void SwitchOffClusterPairRBins()    {fUseAverClusterPairRBins = kFALSE ; }
91 //  
92 //  void SwitchOnClusterPairRWeightBins() {fUseAverClusterPairRWeightBins = kTRUE  ; }
93 //  void SwitchOffClusterPairRWeightBins(){fUseAverClusterPairRWeightBins = kFALSE ; }
94
95 //  void SwitchOnEMaxBins()             {fUseEMaxBins = kTRUE  ; }
96 //  void SwitchOffEMaxBins()            {fUseEMaxBins = kFALSE ; }
97
98   //-------------------------------
99         //Opening angle pair selection
100   //-------------------------------
101   void SwitchOnAngleSelection()      {fUseAngleCut = kTRUE      ; }
102   void SwitchOffAngleSelection()     {fUseAngleCut = kFALSE     ; }
103   void SwitchOnAngleEDepSelection()  {fUseAngleEDepCut = kTRUE  ; }
104   void SwitchOffAngleEDepSelection() {fUseAngleEDepCut = kFALSE ; }
105   void SetAngleCut(Float_t a)        {fAngleCut    = a          ; }
106   void SetAngleMaxCut(Float_t a)     {fAngleMaxCut = a          ; }
107
108   //-------------------------------
109   // Use mixing code of this class
110   //-------------------------------
111   void SwitchOnOwnMix()              {fDoOwnMix = kTRUE  ; }
112   void SwitchOffOwnMix()             {fDoOwnMix = kFALSE ; }
113
114   //------------------------------------------
115   //Do analysis only with clusters in same SM or different combinations of SM
116   //------------------------------------------
117   void SwitchOnSameSM()              {fSameSM = kTRUE    ; }
118   void SwitchOffSameSM()             {fSameSM = kFALSE   ; }
119   
120   void SwitchOnSMCombinations()   {fFillSMCombinations = kTRUE  ; }
121   void SwitchOffSMCombinations()  {fFillSMCombinations = kFALSE ; }
122   
123   //-------------------------------
124   //Histogram filling options off by default
125   //-------------------------------
126   void SwitchOnInvPtWeight()         {fMakeInvPtPlots = kTRUE  ; }
127   void SwitchOffInvPtWeight()        {fMakeInvPtPlots = kFALSE ; }
128   
129   void SwitchOnFillBadDistHisto()    {fFillBadDistHisto    = kTRUE;}
130   void SwitchOffFillBadDistHisto()   {fFillBadDistHisto    = kFALSE;}
131   
132   //-------------------------------------------
133   //Cuts for multiple analysis, off by default
134   //-------------------------------------------
135   void SwitchOnMultipleCutAnalysis()          {fMultiCutAna    = kTRUE ;}
136   void SwitchOffMultipleCutAnalysis()         {fMultiCutAna    = kFALSE;}
137
138   void SetNPtCuts   (Int_t size)              {if(size <= 10)fNPtCuts    = size; }
139   void SetNAsymCuts (Int_t size)              {if(size <= 10)fNAsymCuts  = size; }
140   void SetNNCellCuts(Int_t size)              {if(size <= 10)fNCellNCuts = size; }
141   void SetNPIDBits  (Int_t size)              {if(size <= 10)fNPIDBits   = size; }
142   
143   void SetPtCutsAt   (Int_t pos,Float_t val)  {if(pos < 10)fPtCuts[pos]    = val;}
144   void SetAsymCutsAt (Int_t pos,Float_t val)  {if(pos < 10)fAsymCuts[pos]  = val;}
145   void SetNCellCutsAt(Int_t pos,Int_t val)    {if(pos < 10)fCellNCuts[pos] = val;}
146   void SetPIDBitsAt  (Int_t pos,Int_t val)    {if(pos < 10)fPIDBits[pos]   = val;}
147   
148   //MC analysis related methods
149   void FillAcceptanceHistograms();
150   void FillMCVersusRecDataHistograms(const Int_t    index1,  const Int_t    index2,
151                                      const Float_t  pt1,     const Float_t  pt2, 
152                                      const Int_t    ncells1, const Int_t    ncells2, 
153                                      const Double_t mass,    const Double_t pt,  const Double_t asym,    
154                                      const Double_t deta,    const Double_t dphi);
155   
156   void SwitchOnMultipleCutAnalysisInSimulation()   {fMultiCutAnaSim = kTRUE;}
157   void SwitchOffMultipleCutAnalysisInSimulation()  {fMultiCutAnaSim = kFALSE;}
158
159   void     SwitchOnConversionChecker()             { fCheckConversion = kTRUE  ; }
160   void     SwitchOffConversionChecker()            { fCheckConversion = kFALSE ; }  
161   
162   Double_t WeightPi0(Int_t pi0Id);
163   Int_t    GetMotherPi0Index(Int_t label);
164   
165   private:
166   Bool_t   fDoOwnMix;            // Do combinatorial background not the one provided by the frame
167   Int_t    fNCentrBin ;          // Number of bins in event container for centrality
168 //Int_t    fNZvertBin ;          // Number of bins in event container for vertex position
169 //Int_t    fNrpBin ;               // Number of bins in event container for reaction plain
170   Int_t    fNmaxMixEv ;          // Maximal number of events stored in buffer for mixing
171   TString  fCalorimeter ;        // Select Calorimeter for IM
172   Int_t    fNModules ;           // Number of EMCAL/PHOS modules, set as many histogras as modules 
173   Bool_t   fUseAngleCut ;        // Select pairs depending on their opening angle
174   Bool_t   fUseAngleEDepCut ;    // Select pairs depending on their opening angle
175   Float_t  fAngleCut ;           // Select pairs with opening angle larger than a threshold
176   Float_t  fAngleMaxCut ;        // Select pairs with opening angle smaller than a threshold
177   TList ** fEventsList ;         //![fNCentrBin*GetNZvertBin()*GetNRPBin()] Containers for photons in stored events
178   
179   //Multiple cuts analysis
180   Bool_t   fMultiCutAna;         // Do analysis with several or fixed cut
181   Bool_t   fMultiCutAnaSim;      // Do analysis with several or fixed cut, in the simulation related part
182   Int_t    fNPtCuts;             // Number of pt cuts
183   Float_t  fPtCuts[10];          // Array with different pt cuts
184   Int_t    fNAsymCuts;           // Number of assymmetry cuts
185   Float_t  fAsymCuts[10];        // Array with different assymetry cuts
186   Int_t    fNCellNCuts;          // Number of cuts with number of cells in cluster
187   Int_t    fCellNCuts[10];       // Array with different cell number cluster cuts
188   Int_t    fNPIDBits ;                 // Number of possible PID bit combinations
189   Int_t    fPIDBits[10];         // Array with different PID bits
190   
191   //Switchs of different analysis options
192   Bool_t   fMakeInvPtPlots;      // D plots with inverse pt weight
193   Bool_t   fSameSM;              // Select only pairs in same SM;
194   Bool_t   fFillSMCombinations;  // Fill histograms with different cluster pairs in SM combinations
195   Bool_t   fCheckConversion;     // Fill histograms with tagged photons as conversion
196   Bool_t   fUseTrackMultBins;    // Use track multiplicity and not centrality bins
197   Bool_t   fUsePhotonMultBins;   // Use photon multiplicity and not centrality bins
198   Bool_t   fUseAverCellEBins;    // Use cell average energy and not centrality bins
199   Bool_t   fUseAverClusterEBins; // Use cluster average energy and not centrality bins
200   Bool_t   fUseAverClusterEDenBins; // Use cluster average energy density and not centrality bins
201 //  Bool_t   fUseAverClusterPairRBins; // Use cluster average energy and not centrality bins
202 //  Bool_t   fUseAverClusterPairRWeightBins; // Use cluster average energy and not centrality bins
203 //  Bool_t   fUseEMaxBins;         // Use Emax bins
204   Bool_t   fFillBadDistHisto;    // Do plots for different distances to bad channels
205   
206   //Histograms
207   
208   //Event characterization
209   TH1F* fhAverTotECluster;          //! Average number of clusters in SM
210   TH1F* fhAverTotECell;             //! Average number of cells    in SM
211   TH2F* fhAverTotECellvsCluster;    //! Average number of cells    in SM
212   TH1F* fhEDensityCluster;          //! Deposited energy in event per cluster
213   TH1F* fhEDensityCell;             //! Deposited energy in event per cell vs cluster
214   TH2F* fhEDensityCellvsCluster;    //! Deposited energy in event per cell vs cluster
215 //  TH1F* fhClusterPairDist;          //! Distance between clusters
216 //  TH1F* fhClusterPairDistWeight;    //! Distance between clusters weighted by pair energy
217 //  TH1F* fhAverClusterPairDist;      //! Average distance between cluster pairs
218 //  TH1F* fhAverClusterPairDistWeight;//! Average distance between cluster pairs weighted with pair energy
219 //  TH2F* fhAverClusterPairDistvsAverE;        //! Average distance between cluster pairs vs average cluster energy
220 //  TH2F* fhAverClusterPairDistWeightvsAverE;  //! Average distance between cluster pairs weighted with pair energy vs average cluster energy
221 //  TH2F* fhAverClusterPairDistvsN;            //! Average distance between cluster pairs vs number of clusters
222 //  TH2F* fhAverClusterPairDistWeightvsN;      //! Average distance between cluster pairs weighted with pair energy vs number of clusters
223 //  TH2F* fhMaxEvsClustMult;          //!
224 //  TH2F* fhMaxEvsClustEDen;          //!
225
226   
227   TH2D ** fhReMod ;                 //![fNModules]   REAL  two-photon invariant mass distribution for different calorimeter modules.
228   TH2D ** fhReSameSideEMCALMod ;    //![fNModules-2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
229   TH2D ** fhReSameSectorEMCALMod ;  //![fNModules/2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
230   TH2D ** fhReDiffPHOSMod ;         //![fNModules]   REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
231   TH2D ** fhMiMod ;                 //![fNModules]   MIXED two-photon invariant mass distribution for different calorimeter modules.
232   TH2D ** fhMiSameSideEMCALMod ;    //![fNModules-2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
233   TH2D ** fhMiSameSectorEMCALMod ;  //![fNModules/2] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
234   TH2D ** fhMiDiffPHOSMod ;         //![fNModules-1] REAL  two-photon invariant mass distribution for different clusters in different calorimeter modules.
235   
236   // Pairs with at least one cluster tagged as conversion
237   TH2D * fhReConv ;                 //! REAL  two-photon invariant mass distribution one of the pair was 2 clusters with small mass 
238   TH2D * fhMiConv ;                 //! MIXED two-photon invariant mass distribution one of the pair was 2 clusters with small mass
239   TH2D * fhReConv2 ;                //! REAL  two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass 
240   TH2D * fhMiConv2 ;                //! MIXED two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass
241
242   TH2D ** fhRe1 ;                   //![fNCentrBin*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
243   TH2D ** fhMi1 ;                   //![fNCentrBin*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
244   TH2D ** fhRe2 ;                   //![fNCentrBin*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
245   TH2D ** fhMi2 ;                   //![fNCentrBin*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
246   TH2D ** fhRe3 ;                   //![fNCentrBin*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry 
247   TH2D ** fhMi3 ;                   //![fNCentrBin*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
248   
249   //Histograms weighted by inverse pT
250   TH2D ** fhReInvPt1 ;              //![fNCentrBin*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
251   TH2D ** fhMiInvPt1 ;              //![fNCentrBin*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
252   TH2D ** fhReInvPt2 ;              //![fNCentrBin*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT 
253   TH2D ** fhMiInvPt2 ;              //![fNCentrBin*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
254   TH2D ** fhReInvPt3 ;              //![fNCentrBin*fNPIDBits*fNAsymCuts] REAL  two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
255   TH2D ** fhMiInvPt3 ;              //![fNCentrBin*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
256   
257   //Multiple cuts: Assymmetry, pt, n cells, PID
258   TH2D ** fhRePtNCellAsymCuts ;     //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
259   TH2D ** fhRePtNCellAsymCutsSM0 ;  //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
260   TH2D ** fhRePtNCellAsymCutsSM1 ;  //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
261   TH2D ** fhRePtNCellAsymCutsSM2 ;  //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
262   TH2D ** fhRePtNCellAsymCutsSM3 ;  //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
263   TH2D ** fhMiPtNCellAsymCuts ;     //![fNPtCuts*fNAsymCuts*fNCellNCuts] Mixed two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
264   TH2D ** fhRePIDBits ;             //![fNPIDBits]  REAL two-photon invariant mass distribution for different PID bits
265   TH3D ** fhRePtMult ;              //![fNAsymCuts] REAL two-photon invariant mass distribution for different track multiplicity and assymetry cuts
266   
267   // Asymmetry vs pt, in pi0/eta regions
268   TH2D *  fhRePtAsym    ;           //! REAL two-photon pt vs asymmetry
269   TH2D *  fhRePtAsymPi0 ;           //! REAL two-photon pt vs asymmetry, close to pi0 mass
270   TH2D *  fhRePtAsymEta ;           //! REAL two-photon pt vs asymmetry, close to eta mass
271
272   TH3D * fhEvents;                  //! Number of events per centrality, RP, zbin
273   TH1D * fhCentrality;              //! Histogram with centrality bins with at least one pare
274   TH1D * fhCentralityNoPair;        //! Histogram with centrality bins with no pair
275
276   // Pair opening angle
277   TH2D * fhRealOpeningAngle ;       //! Opening angle of pair versus pair energy
278   TH2D * fhRealCosOpeningAngle ;    //! Cosinus of opening angle of pair version pair energy
279   TH2D * fhMixedOpeningAngle ;      //! Opening angle of pair versus pair energy
280   TH2D * fhMixedCosOpeningAngle ;   //! Cosinus of opening angle of pair version pair energy
281   
282   //MC analysis histograms
283   //Pi0 Acceptance
284   TH1D * fhPrimPi0Pt ;              //! Spectrum of Primary 
285   TH1D * fhPrimPi0AccPt ;           //! Spectrum of primary with accepted daughters 
286   TH2D * fhPrimPi0Y ;               //! Rapidity distribution of primary particles  vs pT
287   TH2D * fhPrimPi0AccY ;            //! Rapidity distribution of primary with accepted daughters  vs pT
288   TH2D * fhPrimPi0Phi ;             //! Azimutal distribution of primary particles  vs pT
289   TH2D * fhPrimPi0AccPhi;           //! Azimutal distribution of primary with accepted daughters  vs pT
290   TH2D * fhPrimPi0OpeningAngle ;    //! Opening angle of pair versus pair energy, primaries
291   TH2D * fhPrimPi0CosOpeningAngle ; //! Cosinus of opening angle of pair version pair energy, primaries
292   //Eta acceptance
293   TH1D * fhPrimEtaPt ;              //! Spectrum of Primary 
294   TH1D * fhPrimEtaAccPt ;           //! Spectrum of primary with accepted daughters 
295   TH2D * fhPrimEtaY ;               //! Rapidity distribution of primary particles vs pT
296   TH2D * fhPrimEtaAccY ;            //! Rapidity distribution of primary with accepted daughters  vs pT
297   TH2D * fhPrimEtaPhi ;             //! Azimutal distribution of primary particles  vs pT
298   TH2D * fhPrimEtaAccPhi;           //! Azimutal distribution of primary with accepted daughters         vs pT
299   
300   // Primaries origin
301   TH2D * fhPrimPi0PtOrigin ;        //! Spectrum of generated pi0 vs mother
302   TH2D * fhPrimEtaPtOrigin ;        //! Spectrum of generated eta vs mother
303   
304   //Pair origin
305   //Array of histograms ordered as follows: 0-Photon, 1-electron, 2-pi0, 3-eta, 4-a-proton, 5-a-neutron, 6-stable particles, 
306   // 7-other decays, 8-string, 9-final parton, 10-initial parton, intermediate, 11-colliding proton, 12-unrelated
307   TH2D * fhMCOrgMass[13];           //! Mass vs pt of real pairs, check common origin of pair
308   TH2D * fhMCOrgAsym[13];           //! Asymmetry vs pt of real pairs, check common origin of pair
309   TH2D * fhMCOrgDeltaEta[13];       //! Delta Eta vs pt of real pairs, check common origin of pair
310   TH2D * fhMCOrgDeltaPhi[13];       //! Delta Phi vs pt of real pairs, check common origin of pair
311   
312   //Multiple cuts in simulation, origin pi0 or eta
313   TH2D ** fhMCPi0MassPtRec;         //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs reconstructed pt of original pair  
314   TH2D ** fhMCPi0MassPtTrue;        //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs generated pt of original pair  
315   TH2D ** fhMCPi0PtTruePtRec;       //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed pt vs generated pt of pair
316   TH2D ** fhMCEtaMassPtRec;         //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs reconstructed pt of original pair  
317   TH2D ** fhMCEtaMassPtTrue;        //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs generated pt of original pair  
318   TH2D ** fhMCEtaPtTruePtRec;       //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed pt vs generated pt of pair
319   
320   TH2D *  fhMCPi0PtOrigin ;         //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
321   TH2D *  fhMCEtaPtOrigin ;         //! Mass of reoconstructed pi0 pairs  in calorimeter vs mother
322
323   Int_t fFirstPi0;
324   Int_t fLastPi0;
325
326   ClassDef(AliAnaPi0,18)
327 } ;
328
329
330 #endif //ALIANAPI0_H
331
332
333