]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/CaloTrackCorrelations/AliAnaCalorimeterQA.h
clean include dependencies, rename histogram
[u/mrichter/AliRoot.git] / PWGGA / CaloTrackCorrelations / AliAnaCalorimeterQA.h
1 #ifndef ALIANACALORIMETERQA_H
2 #define ALIANACALORIMETERQA_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 check results from simulations or reconstructed real data. 
8 // Fill few histograms and do some checking plots
9 //
10 //-- Author: Gustavo Conesa (INFN-LNF)
11
12 // --- Root system ---
13 class TH3F;
14 class TH2F;
15 class TH1F;
16 class TObjString;
17 class TObjArray;
18
19 // --- Analysis system --- 
20 class AliVCaloCells;
21 class AliVCaloCluster;
22 class AliVTrack;
23
24 #include "AliAnaCaloTrackCorrBaseClass.h"
25  
26 class AliAnaCalorimeterQA : public AliAnaCaloTrackCorrBaseClass {
27   
28 public: 
29   AliAnaCalorimeterQA() ; // default ctor       
30   virtual ~AliAnaCalorimeterQA() {;} //virtual dtor
31     
32   // General methods
33   
34   TObjString * GetAnalysisCuts();
35
36   TList      * GetCreateOutputObjects();
37   
38   void         Init();
39   
40   void         InitParameters();
41     
42   void         MakeAnalysisFillHistograms() ;
43   
44   void         Print(const Option_t * opt) const;
45     
46   // Main methods
47   
48   void         BadClusterHistograms(AliVCluster* clus,    const TObjArray *caloClusters,  AliVCaloCells * cells, 
49                                     const Int_t absIdMax, const Double_t maxCellFraction, const Float_t eCrossFrac,
50                                     const Double_t tmax, Double_t timeAverages[2]);  
51     
52   void         CalculateAverageTime(AliVCluster *clus, AliVCaloCells *cells, Double_t timeAverages[2]);
53   
54   void         CellHistograms(AliVCaloCells * cells);
55
56   void         CellInClusterPositionHistograms(AliVCluster* cluster);
57     
58   void         ClusterAsymmetryHistograms(AliVCluster* clus, const Int_t absIdMax, const Bool_t goodCluster );
59   
60   void         ClusterHistograms(AliVCluster* cluster, const TObjArray *caloClusters,  AliVCaloCells * cells, 
61                                  const Int_t absIdMax, const Double_t maxCellFraction, const Float_t eCrossFrac,
62                                  const Double_t tmax, Double_t timeAverages[2]);
63   
64   void         ClusterLoopHistograms(const TObjArray * clusters, AliVCaloCells * cells);
65   
66   Bool_t       ClusterMCHistograms(const TLorentzVector mom,const Bool_t matched,
67                                    const Int_t * labels, const Int_t nLabels, Int_t & pdg );
68
69   void         ClusterMatchedWithTrackHistograms(AliVCluster* clus, TLorentzVector mom, 
70                                                  const Bool_t mcOK, const Int_t pdg);
71
72   void         Correlate();
73   
74   Float_t      GetECross(const Int_t absId, AliVCaloCells* cells);
75   
76   void         InvariantMassHistograms(const Int_t iclus, const TLorentzVector mom, const Int_t nModule,
77                                        const TObjArray* caloClusters, AliVCaloCells * cells);
78
79   Bool_t       IsGoodCluster(const Int_t absIdMax, AliVCaloCells *cells);
80   
81   void         MCHistograms();  
82   
83   void         MCHistograms(const TLorentzVector mom, const Int_t pdg);
84     
85   void         WeightHistograms(AliVCluster *clus, AliVCaloCells* cells);
86
87   // Setters and Getters
88
89   
90   Float_t      GetEMCALCellAmpMin()      const  { return fEMCALCellAmpMin    ; }
91   void         SetEMCALCellAmpMin(Float_t amp)  { fEMCALCellAmpMin = amp     ; }
92   
93   Float_t      GetPHOSCellAmpMin()       const  { return fPHOSCellAmpMin     ; }
94   void         SetPHOSCellAmpMin (Float_t amp)  { fPHOSCellAmpMin  = amp     ; }
95     
96   TString      GetCalorimeter()          const  { return fCalorimeter        ; }
97   void         SetCalorimeter(TString calo)     { fCalorimeter = calo        ; }
98     
99   void         SetNumberOfModules(Int_t nmod)   { fNModules   = nmod         ; }
100   
101   Double_t     GetTimeCutMin()           const  { return fTimeCutMin         ; }
102   Double_t     GetTimeCutMax()           const  { return fTimeCutMax         ; }
103   void         SetTimeCut(Double_t min, Double_t max) {
104                           fTimeCutMin = min ; fTimeCutMax = max              ; }
105     
106   // Histogram Switchs
107   
108   void SwitchOnFillAllPositionHistogram()       { fFillAllPosHisto  = kTRUE  ; }
109   void SwitchOffFillAllPositionHistogram()      { fFillAllPosHisto  = kFALSE ; }
110   
111   void SwitchOnFillAllPositionHistogram2()      { fFillAllPosHisto2 = kTRUE  ; }
112   void SwitchOffFillAllPositionHistogram2()     { fFillAllPosHisto2 = kFALSE ; }
113   
114   void SwitchOnFillAllTH12Histogram()           { fFillAllTH12      = kTRUE  ; }
115   void SwitchOffFillAllTH12Histogram()          { fFillAllTH12      = kFALSE ; }
116   
117   void SwitchOnFillAllTH3Histogram()            { fFillAllTH3       = kTRUE  ; }
118   void SwitchOffFillAllTH3Histogram()           { fFillAllTH3       = kFALSE ; }
119   
120   void SwitchOnFillAllTrackMatchingHistogram()  { fFillAllTMHisto   = kTRUE  ; }
121   void SwitchOffFillAllTrackMatchingHistogram() { fFillAllTMHisto   = kFALSE ; }
122   
123   void SwitchOnFillAllPi0Histogram()            { fFillAllPi0Histo  = kTRUE  ; }
124   void SwitchOffFillAllPi0Histogram()           { fFillAllPi0Histo  = kFALSE ; }
125
126   void SwitchOnCorrelation()                    { fCorrelate        = kTRUE  ; }
127   void SwitchOffCorrelation()                   { fCorrelate        = kFALSE ; }
128
129   void SwitchOnStudyBadClusters()               { fStudyBadClusters = kTRUE  ; }
130   void SwitchOffStudyBadClusters()              { fStudyBadClusters = kFALSE ; }
131   
132   void SwitchOnStudyClustersAsymmetry()         { fStudyClustersAsymmetry = kTRUE  ; }
133   void SwitchOffStudyClustersAsymmetry()        { fStudyClustersAsymmetry = kFALSE ; }
134
135   void SwitchOnStudyWeight()                    { fStudyWeight      = kTRUE  ; }
136   void SwitchOffStudyWeight()                   { fStudyWeight      = kFALSE ; }
137
138   
139  private:
140   
141   TString  fCalorimeter ;                     // Calorimeter selection
142   
143   //Switches
144   Bool_t   fFillAllPosHisto;                  // Fill all the position related histograms 
145   Bool_t   fFillAllPosHisto2;                 // Fill all the position related histograms 2
146   Bool_t   fFillAllTH12 ;                     // Fill simple histograms which information is already in TH3 histograms
147   Bool_t   fFillAllTH3 ;                      // Fill TH3 histograms
148   Bool_t   fFillAllTMHisto ;                  // Fill track matching histograms
149   Bool_t   fFillAllPi0Histo ;                 // Fill invariant mass histograms
150   Bool_t   fCorrelate   ;                     // Correlate PHOS/EMCAL cells/clusters, also with V0 and track multiplicity
151   Bool_t   fStudyBadClusters;                 // Study bad clusters
152   Bool_t   fStudyClustersAsymmetry;           // Study asymmetry of clusters
153   Bool_t   fStudyWeight;                      // Study the energy weight used in different cluster calculations
154   
155   // Parameters
156   Int_t    fNModules    ;                     // Number of EMCAL/PHOS modules
157   Int_t    fNRCU        ;                     // Number of EMCAL/PHOS RCU 
158   Int_t    fNMaxCols    ;                     // Number of EMCAL/PHOS rows 
159   Int_t    fNMaxRows    ;                     // Number of EMCAL/PHOS columns
160   
161   //Cuts
162   Double_t fTimeCutMin  ;                     // Remove clusters/cells with time smaller than this value, in ns
163   Double_t fTimeCutMax  ;                     // Remove clusters/cells with time larger than this value, in ns
164   Float_t  fEMCALCellAmpMin;                  // amplitude Threshold on emcal cells
165   Float_t  fPHOSCellAmpMin ;                  // amplitude Threshold on phos cells
166   
167   //CaloClusters 
168   TH1F *   fhE  ;                             //! E distribution, Reco
169   TH1F *   fhPt ;                             //! pT distribution, Reco
170   TH1F *   fhPhi;                             //! phi distribution, Reco 
171   TH1F *   fhEta;                             //! eta distribution, Reco 
172   TH3F *   fhEtaPhiE  ;                       //! eta vs phi vs E, Reco
173   TH1F *   fhECharged  ;                      //! E distribution, Reco, matched with track
174   TH1F *   fhPtCharged ;                      //! pT distribution, Reco, matched with track
175   TH1F *   fhPhiCharged;                      //! phi distribution, Reco, matched with track 
176   TH1F *   fhEtaCharged;                      //! eta distribution, Reco, matched with track 
177   TH3F *   fhEtaPhiECharged;                  //! eta vs phi vs E, Reco, matched with track 
178     
179   TH2F *   fhIM;                              //! cluster pairs invariant mass
180   TH2F *   fhAsym;                            //! cluster pairs invariant mass  
181   
182   TH2F *   fhNCellsPerCluster;                //! N cells per cluster vs cluster energy vs eta of cluster       
183   TH2F *   fhNCellsPerClusterNoCut;           //! N cells per cluster vs cluster energy vs eta of cluster       
184
185   TH1F *   fhNClusters;                       //! Number of clusters
186
187   TH2F *   fhClusterTimeEnergy;               //! Cluster Time vs Energy 
188   TH2F *   fhCellTimeSpreadRespectToCellMax;  //! Difference of the time of cell with maximum dep energy and the rest of cells
189   TH1F *   fhCellIdCellLargeTimeSpread;       //! Cells with large time respect to max (diff > 100 ns)
190   TH2F *   fhClusterPairDiffTimeE;            //! Pair of clusters time difference vs E
191   
192   TH2F *   fhClusterMaxCellCloseCellRatio;    //! Ratio between max cell energy and cell energy of the same cluster 
193   TH2F *   fhClusterMaxCellCloseCellDiff;     //! Difference between max cell energy and cell energy of the same cluster   
194   TH2F *   fhClusterMaxCellDiff;              //! Difference between cluster energy and energy of cell with more energy, good clusters only
195   TH2F *   fhClusterMaxCellDiffNoCut;         //! Difference between cluster energy and energy of cell with more energy, no bad cluster rejection
196   
197   TH2F *   fhClusterMaxCellDiffAverageTime;      //! Difference between cluster average time and time of cell with more energy
198   TH2F *   fhClusterMaxCellDiffWeightedTime;     //! Difference between cluster weighted time and time of cell with more energy
199   TH2F *   fhClusterMaxCellECross;               //! 1 - Energy in cross around max energy cell / max energy cell vs cluster energy, good clusters
200   
201   TH2F *   fhLambda0;                         //! cluster Lambda0    vs Energy
202   TH2F *   fhLambda1;                         //! cluster Lambda1    vs Energy
203   TH2F *   fhDispersion;                      //! cluster Dispersion vs Energy
204   
205   // Bad clusters histograms
206   TH1F *   fhBadClusterEnergy;                //! energy of bad cluster
207   TH2F *   fhBadClusterTimeEnergy;            //! Time Max cell of bad cluster
208   TH2F *   fhBadClusterPairDiffTimeE;         //! Pair of clusters time difference vs E, bad cluster
209   TH2F *   fhBadCellTimeSpreadRespectToCellMax; //! Difference of the time of cell with maximum dep energy and the rest of cells for bad clusters
210   
211   TH2F *   fhBadClusterMaxCellCloseCellRatio; //! Ratio between max cell energy and cell energy of the same cluster for bad clusters 
212   TH2F *   fhBadClusterMaxCellCloseCellDiff ; //! Difference between max cell energy and cell energy of the same cluster for bad clusters 
213   TH2F *   fhBadClusterMaxCellDiff;           //! Difference between cluster energy and energy of cell with more energy
214   
215   TH2F *   fhBadClusterMaxCellDiffAverageTime;      //! Difference between cluster average time and time of cell with more energy
216   TH2F *   fhBadClusterMaxCellDiffWeightedTime;     //! Difference between cluster weighted time and time of cell with more energy
217   TH2F *   fhBadClusterMaxCellECross;               //! 1 - Energy in cross around max energy cell / max energy cell vs cluster energy, bad clusters
218
219   // Cluster cell size
220   TH2F *   fhDeltaIEtaDeltaIPhiE0[2];         //! Difference between max cell index and farthest cell, eta vs phi, E < 2 GeV, with and without matching; 
221   TH2F *   fhDeltaIEtaDeltaIPhiE2[2];         //! Difference between max cell index and farthest cell, eta vs phi, 2 < E < 6 GeV, with and without matching; 
222   TH2F *   fhDeltaIEtaDeltaIPhiE6[2];         //! Difference between max cell index and farthest cell, eta vs phi, E > 6 GeV, with and without matching; 
223   TH2F *   fhDeltaIA[2];                      //! Cluster "asymmetry" in cell terms vs E, with and without matching
224   TH2F *   fhDeltaIAL0[2];                    //! Cluster "asymmetry" in cell units vs Lambda0    for E > 0.5 GeV, n cells in cluster > 3, with and without matching
225   TH2F *   fhDeltaIAL1[2];                    //! Cluster "asymmetry" in cell units vs Lambda1    for E > 0.5 GeV, n cells in cluster > 3, with and without matching
226   TH2F *   fhDeltaIANCells[2] ;               //! Cluster "asymmetry" in cell units vs number of cells in cluster for E > 0.5, with and without matching
227   TH2F *   fhDeltaIAMC[4];                    //! Cluster "asymmetry" in cell terms vs E, from MC photon, electron, conversion or hadron.
228   TH2F *   fhBadClusterDeltaIEtaDeltaIPhiE0;  //! Difference between max cell index and farthest cell, eta vs phi, E < 2 GeV, with and without matching; bad clusters. 
229   TH2F *   fhBadClusterDeltaIEtaDeltaIPhiE2;  //! Difference between max cell index and farthest cell, eta vs phi, 2 < E < 6 GeV, with and without matching; bad clusters.
230   TH2F *   fhBadClusterDeltaIEtaDeltaIPhiE6;  //! Difference between max cell index and farthest cell, eta vs phi, E > 6 GeV, with and without matching; bad clusters.
231   TH2F *   fhBadClusterDeltaIA;               //! Cluster "asymmetry" in cell terms vs E, with and without matching; bad clusters.
232   
233   //Cluster/cell Position
234   TH2F *   fhRNCells ;                        //! R=sqrt(x^2+y^2) (cm) cluster distribution vs N cells in cluster
235   TH2F *   fhXNCells ;                        //! X (cm) cluster distribution vs N cells in cluster
236   TH2F *   fhYNCells ;                        //! Y (cm) cluster distribution vs N cells in cluster
237   TH2F *   fhZNCells ;                        //! Z (cm) cluster distribution vs N cells in cluster
238         
239   TH2F *   fhRE ;                             //! R=sqrt(x^2+y^2) (cm) cluster distribution vs cluster energy
240   TH2F *   fhXE ;                             //! X (cm) cluster distribution vs cluster energy
241   TH2F *   fhYE ;                             //! Y (cm) cluster distribution vs cluster energy
242   TH2F *   fhZE ;                             //! Z (cm) cluster distribution vs cluster energy
243   TH3F *   fhXYZ;                             //! cluster X vs Y vs Z (cm)
244         
245   TH2F *   fhRCellE ;                         //! R=sqrt(x^2+y^2) (cm) cell distribution vs cell energy
246   TH2F *   fhXCellE ;                         //! X (cm) cell distribution vs cell energy
247   TH2F *   fhYCellE ;                         //! Y (cm) cell distribution vs cell energy
248   TH2F *   fhZCellE ;                         //! Z (cm) cell distribution vs cell energy
249   TH3F *   fhXYZCell;                         //! cell X vs Y vs Z (cm)
250   
251   TH2F *   fhDeltaCellClusterRNCells ;        //! R cluster - R cell distribution (cm) vs N cells in cluster
252   TH2F *   fhDeltaCellClusterXNCells ;        //! X cluster - X cell distribution (cm) vs N cells in cluster
253   TH2F *   fhDeltaCellClusterYNCells ;        //! Y cluster - Y cell distribution (cm) vs N cells in cluster
254   TH2F *   fhDeltaCellClusterZNCells ;        //! Z cluster - Z cell distribution (cm) vs N cells in cluster
255         
256   TH2F *   fhDeltaCellClusterRE ;             //! R cluster - R cell distribution (cm) vs cluster energy
257   TH2F *   fhDeltaCellClusterXE ;             //! X cluster - X cell distribution (cm) vs cluster energy
258   TH2F *   fhDeltaCellClusterYE ;             //! Y cluster - Y cell distribution (cm) vs cluster energy
259   TH2F *   fhDeltaCellClusterZE ;             //! Z cluster - Z cell distribution (cm) vs cluster energy
260         
261   //Calo Cells
262   TH1F *   fhNCells;                          //! Number of towers/crystals with signal
263   TH1F *   fhAmplitude;                       //! Amplitude measured in towers/crystals
264   TH2F *   fhAmpId;                           //! Amplitude measured in towers/crystals vs id of tower. 
265   TH3F *   fhEtaPhiAmp;                       //! eta vs phi vs amplitude, cells
266    
267   TH1F *   fhTime;                            //! Time measured in towers/crystals
268   TH2F *   fhTimeVz;                          //! Time measured in towers/crystals vs vertex z component, for E > 0.5
269   TH2F *   fhTimeId;                          //! Time vs Absolute cell Id
270   TH2F *   fhTimeAmp;                         //! Time vs Amplitude 
271   
272   TH2F *   fhCellECross;                      //! 1 - Energy in cross around cell /  cell energy 
273   
274   //Calorimeters Correlation
275   TH2F *   fhCaloCorrNClusters;               //! EMCAL vs PHOS, number of clusters     
276   TH2F *   fhCaloCorrEClusters;               //! EMCAL vs PHOS, total measured cluster energy
277   TH2F *   fhCaloCorrNCells;                  //! EMCAL vs PHOS, number of cells
278   TH2F *   fhCaloCorrECells;                  //! EMCAL vs PHOS,  total measured cell energy
279         
280   //V0 Correlation
281   TH2F *   fhCaloV0SCorrNClusters;            //! Calo vs V0 signal , number of clusters        
282   TH2F *   fhCaloV0SCorrEClusters;            //! Calo vs V0 signal, total measured cluster energy
283   TH2F *   fhCaloV0SCorrNCells;               //! Calo vs V0 signal, number of cells
284   TH2F *   fhCaloV0SCorrECells;               //! Calo vs V0 signal,  total measured cell energy
285   TH2F *   fhCaloV0MCorrNClusters;            //! Calo vs V0 multiplicity , number of clusters  
286   TH2F *   fhCaloV0MCorrEClusters;            //! Calo vs V0 multiplicity, total measured cluster energy
287   TH2F *   fhCaloV0MCorrNCells;               //! Calo vs V0 multiplicity, number of cells
288   TH2F *   fhCaloV0MCorrECells;               //! Calo vs V0 multiplicity,  total measured cell energy
289   
290   //Track Correlation
291   TH2F *   fhCaloTrackMCorrNClusters;         //! Calo vs Track Multiplicity, number of clusters        
292   TH2F *   fhCaloTrackMCorrEClusters;         //! Calo vs Track Multiplicity, total measured cluster energy
293   TH2F *   fhCaloTrackMCorrNCells;            //! Calo vs V0 Track Multiplicity, number of cells
294   TH2F *   fhCaloTrackMCorrECells;            //! Calo vs V0 Track Multipliticy,  total measured cell energy
295   
296   //Module histograms
297   TH2F *   fhEMod  ;                          //! cluster E distribution for different module, Reco
298   TH2F *   fhAmpMod ;                         //! cell amplitude distribution for different module, Reco
299   TH2F *   fhTimeMod ;                        //! cell time distribution for different module, Reco
300   TH2F *   fhNClustersMod ;                   //! Number of clusters for different module, Reco
301   TH2F *   fhNCellsMod ;                      //! Number of towers/crystals with signal different module, Reco
302   TH2F **  fhNCellsPerClusterMod ;            //! N cells per clusters different module, Reco
303   TH2F **  fhNCellsPerClusterModNoCut ;       //! N cells per clusters different module, Reco, No cut
304   TH2F *   fhGridCells ;                      //! Cells ordered in column/row for different module, Reco
305   TH2F *   fhGridCellsE ;                     //! Cells ordered in column/row for different module, weighted with energy, Reco
306   TH2F *   fhGridCellsTime ;                  //! Cells ordered in column/row for different module, weighted with time, Reco
307   TH2F **  fhTimeAmpPerRCU;                   //! Time vs Amplitude measured in towers/crystals different RCU
308   TH2F **  fhIMMod;                           //! cluster pairs invariant mass, different module,
309         
310   // Weight studies
311   
312   TH2F* fhECellClusterRatio;                  //! e cell / e cluster vs e cluster
313   TH2F* fhECellClusterLogRatio;               //! log (e cell / e cluster)  vs e cluster
314   TH2F* fhEMaxCellClusterRatio;               //! e max cell / e cluster vs e cluster
315   TH2F* fhEMaxCellClusterLogRatio;            //! log (e max cell / e cluster) vs e cluster
316   
317   TH2F* fhLambda0ForW0[14];                    //! L0 for 7 defined w0= 3, 3.5 ... 6
318   //TH2F* fhLambda1ForW0[7];                    //! L1 for 7 defined w0= 3, 3.5 ... 6
319
320   TH2F* fhLambda0ForW0MC[14][5];               //! L0 for 7 defined w0= 3, 3.5 ... 6, depending on the particle of origin
321   //TH2F* fhLambda1ForW0MC[7][5];               //! L1 for 7 defined w0= 3, 3.5 ... 6, depending on the particle of origin
322   
323   //Pure MC
324
325   enum mcTypes {kmcPhoton = 0, kmcPi0 = 1, kmcEta = 2, kmcElectron = 3, kmcNeHadron = 4, kmcChHadron = 5 };
326   
327   TH2F *   fhRecoMCE[6][2]  ;                 //! E   generated particle vs reconstructed E
328   TH2F *   fhRecoMCPhi[6][2] ;                //! phi generated particle vs reconstructed phi
329   TH2F *   fhRecoMCEta[6][2] ;                //! eta generated particle vs reconstructed Eta
330   TH2F *   fhRecoMCDeltaE[6][2]  ;            //! Gen-Reco E    generated particle vs reconstructed E
331   TH2F *   fhRecoMCRatioE[6][2]  ;            //! Reco/Gen E    generated particle vs reconstructed E
332   TH2F *   fhRecoMCDeltaPhi[6][2];            //! Gen-Reco phi  generated particle vs reconstructed E
333   TH2F *   fhRecoMCDeltaEta[6][2];            //! Gen-Reco eta  generated particle vs reconstructed E
334   
335   TH1F *   fhGenMCE[4]     ;                  //! pt of primary particle
336   TH2F *   fhGenMCEtaPhi[4] ;                 //! eta vs phi of primary particle
337   TH1F *   fhGenMCAccE[4]     ;               //! pt of primary particle, in acceptance
338   TH2F *   fhGenMCAccEtaPhi[4] ;              //! eta vs phi of primary particle, in acceptance
339   
340   TH2F *   fhEMVxyz    ;                      //! Electromagnetic particle production vertex
341   TH2F *   fhEMR       ;                      //! Electromagnetic distance to vertex vs rec energy  
342   TH2F *   fhHaVxyz    ;                      //! Hadron production vertex
343   TH2F *   fhHaR       ;                      //! Hadron distance to vertex vs rec energy  
344         
345   //Histograms for MC track-matching
346   TH2F *   fh1EOverP;                         //! p/E for track-cluster matches
347   TH2F *   fh2dR;                             //! distance between projected track and cluster (eta-phi units)
348   TH2F *   fh2EledEdx;                        //! dE/dx vs. momentum for electron candidates
349   TH2F *   fh2MatchdEdx;                      //! dE/dx vs. momentum for all matches
350         
351   TH2F *   fhMCEle1EOverP;                    //! p/E for track-cluster matches, MC electrons
352   TH1F *   fhMCEle1dR;                        //! distance between projected track and cluster, MC electrons
353   TH2F *   fhMCEle2MatchdEdx;                 //! dE/dx vs. momentum for all matches, MC electrons      
354         
355   TH2F *   fhMCChHad1EOverP;                  //! p/E for track-cluster matches, MC charged hadrons
356   TH1F *   fhMCChHad1dR;                      //! distance between projected track and cluster, MC charged hadrons
357   TH2F *   fhMCChHad2MatchdEdx;               //! dE/dx vs. momentum for all matches, MC charged
358         
359   TH2F *   fhMCNeutral1EOverP;                //! p/E for track-cluster matches, MC neutral
360   TH1F *   fhMCNeutral1dR;                    //! distance between projected track and cluster, MC neutral
361   TH2F *   fhMCNeutral2MatchdEdx;             //! dE/dx vs. momentum for all matches, MC neutral        
362         
363   TH2F *   fh1EOverPR02;                      //! p/E for track-cluster matches, dR > 0.2       
364   TH2F *   fhMCEle1EOverPR02;                 //! p/E for track-cluster matches, dR > 0.2, MC electrons
365   TH2F *   fhMCChHad1EOverPR02;               //! p/E for track-cluster matches, dR > 0.2, MC charged hadrons
366   TH2F *   fhMCNeutral1EOverPR02;             //! p/E for track-cluster matches, dR > 0.2, MC neutral
367         
368   AliAnaCalorimeterQA & operator = (const AliAnaCalorimeterQA & qa) ;//cpy assignment
369   AliAnaCalorimeterQA(              const AliAnaCalorimeterQA & qa) ; // cpy ctor
370   
371   ClassDef(AliAnaCalorimeterQA,23)
372 } ;
373
374
375 #endif //ALIANACALORIMETERQA_H
376
377
378