]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/CaloCalib/AliAnalysisTaskEMCALPi0CalibSelection.h
Possibility to recalculate cluster position and correct for the non linearity of...
[u/mrichter/AliRoot.git] / PWG4 / CaloCalib / AliAnalysisTaskEMCALPi0CalibSelection.h
1 #ifndef ALIANALYSISTASKEMCALPI0CALIBSELECTION_H
2 #define ALIANALYSISTASKEMCALPI0CALIBSELECTION_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 //---------------------------------------------------------------------------// 
8 // Fill histograms with two-cluster invariant mass                           //
9 // using calibration coefficients of the previous iteration.                 //
10 //---------------------------------------------------------------------------//
11
12 // Root includes
13 class TH1F;
14 #include "TH2I.h"
15 #include "TObjArray.h"
16
17 // AliRoot includes
18 #include "AliAnalysisTaskSE.h"
19 class AliEMCALGeometry;
20 class AliAODCaloCluster;
21 class AliAODCaloCells;
22 //class AliEMCALCalibData ;
23 #include "AliEMCALGeoParams.h"
24 class AliEMCALRecoUtils;
25
26 class AliAnalysisTaskEMCALPi0CalibSelection : public AliAnalysisTaskSE
27 {
28 public:
29
30   AliAnalysisTaskEMCALPi0CalibSelection(const char* name);
31   virtual ~AliAnalysisTaskEMCALPi0CalibSelection();
32
33 private:
34   
35   AliAnalysisTaskEMCALPi0CalibSelection(const AliAnalysisTaskEMCALPi0CalibSelection&); 
36   AliAnalysisTaskEMCALPi0CalibSelection& operator=(const AliAnalysisTaskEMCALPi0CalibSelection&); 
37   
38 public:
39   
40   // Implementation of interface methods
41   virtual void UserCreateOutputObjects();
42   virtual void UserExec(Option_t * opt);
43   virtual void LocalInit() ;
44
45   void SetClusterMinEnergy(Float_t emin) {fEmin=emin;}
46   void SetClusterMaxEnergy(Float_t emax) {fEmax=emax;}
47   void SetClusterMinNCells(Int_t n)      {fMinNCells=n;}
48   void SetNCellsGroup(Int_t n)           {fGroupNCells=n;}
49
50   void SetLogWeight(Float_t weight) {fLogWeight=weight;}
51   //void SetCalibCorrections(AliEMCALCalibData* const cdata);
52   void CreateAODFromESD();
53   void CreateAODFromAOD();      
54
55   void CopyAOD(Bool_t copy)   { fCopyAOD = copy ; }
56   Bool_t IsAODCopied() const { return fCopyAOD ; }
57         
58   void SwitchOnSameSM()    {fSameSM = kTRUE  ; }
59   void SwitchOffSameSM()   {fSameSM = kFALSE ; }
60   
61   Int_t  GetEMCALClusters(AliVEvent* event, TRefArray *clusters) const;
62   Bool_t IsEMCALCluster(AliVCluster *clus) const;
63   void SwitchOnOldAODs()   {fOldAOD = kTRUE  ; }
64   void SwitchOffOldAODs()  {fOldAOD = kFALSE ; }  
65   
66   void SetGeometryName(TString name)   { fEMCALGeoName = name ; }
67   TString GeometryName() const { return fEMCALGeoName ; }
68  
69   //Modules fiducial region
70   Bool_t CheckCellFiducialRegion(AliVCluster* cluster, AliVCaloCells* cells) ;
71   void   SetNumberOfCellsFromEMCALBorder(Int_t n) {fNCellsFromEMCALBorder = n; }
72   Int_t  GetNumberOfCellsFromEMCALBorder() const  {return fNCellsFromEMCALBorder; }
73   
74   // Bad channels, copy from PWG4/PartCorrBase/AliCalorimeterUtils
75   Bool_t IsBadChannelsRemovalSwitchedOn()  const { return fRemoveBadChannels ; }
76   void SwitchOnBadChannelsRemoval ()  {fRemoveBadChannels = kTRUE  ; InitEMCALBadChannelStatusMap();}
77   void SwitchOffBadChannelsRemoval()  {fRemoveBadChannels = kFALSE ; }
78         
79   void InitEMCALBadChannelStatusMap() ;
80         
81   Int_t GetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow) const { 
82         if(fEMCALBadChannelMap) return (Int_t) ((TH2I*)fEMCALBadChannelMap->At(iSM))->GetBinContent(iCol,iRow); 
83         else return 0;}//Channel is ok by default
84         
85   void SetEMCALChannelStatus(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) { 
86         if(!fEMCALBadChannelMap)InitEMCALBadChannelStatusMap() ;
87         ((TH2I*)fEMCALBadChannelMap->At(iSM))->SetBinContent(iCol,iRow,c);}
88         
89   TH2I * GetEMCALChannelStatusMap(Int_t iSM) const {return (TH2I*)fEMCALBadChannelMap->At(iSM);}
90         
91   void SetEMCALChannelStatusMap(TObjArray *map) {fEMCALBadChannelMap = map;}
92         
93   Bool_t ClusterContainsBadChannel(UShort_t* cellList, Int_t nCells);
94         
95   // Recalibration
96   Bool_t IsRecalibrationOn()  const { return fRecalibration ; }
97   void SwitchOnRecalibration()    {fRecalibration = kTRUE ; InitEMCALRecalibrationFactors();}
98   void SwitchOffRecalibration()   {fRecalibration = kFALSE ; }
99         
100   void InitEMCALRecalibrationFactors() ;
101   
102   Float_t GetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow) const { 
103         if(fEMCALRecalibrationFactors) return (Float_t) ((TH2F*)fEMCALRecalibrationFactors->At(iSM))->GetBinContent(iCol,iRow); 
104         else return 1;}
105         
106   void SetEMCALChannelRecalibrationFactor(Int_t iSM , Int_t iCol, Int_t iRow, Double_t c = 1) { 
107         if(!fEMCALRecalibrationFactors) InitEMCALRecalibrationFactors();
108         ((TH2F*)fEMCALRecalibrationFactors->At(iSM))->SetBinContent(iCol,iRow,c);}
109         
110   void SetEMCALChannelRecalibrationFactors(Int_t iSM , TH2F* h) {fEMCALRecalibrationFactors->AddAt(h,iSM);}
111         
112   TH2F * GetEMCALChannelRecalibrationFactors(Int_t iSM) const {return (TH2F*)fEMCALRecalibrationFactors->At(iSM);}
113         
114   void SetEMCALChannelRecalibrationFactors(TObjArray *map) {fEMCALRecalibrationFactors = map;}
115   Float_t RecalibrateClusterEnergy(AliAODCaloCluster* cluster, AliAODCaloCells * cells);
116         
117   void SetEMCALRecoUtils(AliEMCALRecoUtils * ru) {fRecoUtils = ru;}
118   AliEMCALRecoUtils* GetEMCALRecoUtils() const {return fRecoUtils;}
119   
120   void SetInvariantMassHistoBinRange(Int_t nBins, Float_t minbin, Float_t maxbin){
121         fNbins = nBins; fMinBin = minbin; fMaxBin = maxbin; }
122           
123 private:
124   void GetMaxEnergyCellPosAndClusterPos(AliVCaloCells* cells, AliVCluster* clu, Int_t& iSM, Int_t& ieta, Int_t& iphi);
125
126 private:
127
128   AliEMCALGeometry * fEMCALGeo;  //! EMCAL geometry
129   //AliEMCALCalibData* fCalibData; // corrections to CC from the previous iteration
130         
131   Float_t fEmin;          // min. cluster energy
132   Float_t fEmax;          // max. cluster energy
133   Int_t   fMinNCells;     // min. ncells in cluster
134   Int_t   fGroupNCells;   // group n cells
135   Float_t fLogWeight;     // log weight used in cluster recalibration
136   Bool_t  fCopyAOD;       // Copy calo information only to AOD?
137   Bool_t  fSameSM;        // Combine clusters in channels on same SM
138   Bool_t  fOldAOD;        // Reading Old AODs, created before release 4.20
139   TString fEMCALGeoName;  // Name of geometry to use.
140   Int_t   fNCellsFromEMCALBorder; //  Number of cells from EMCAL border the cell with maximum amplitude has to be.
141
142   Bool_t     fRemoveBadChannels;         // Check the channel status provided and remove clusters with bad channels
143   TObjArray *fEMCALBadChannelMap;        // Array of histograms with map of bad channels, EMCAL
144   Bool_t     fRecalibration;             // Switch on or off the recalibration
145   TObjArray *fEMCALRecalibrationFactors; // Array of histograms with map of recalibration factors, EMCAL                 
146  
147   AliEMCALRecoUtils * fRecoUtils;  // Access to reconstruction utilities
148   
149   //Output histograms   
150   Int_t   fNbins;  // N       mass bins of invariant mass histograms
151   Float_t fMinBin; // Minimum mass bins of invariant mass histograms
152   Float_t fMaxBin; // Maximum mass bins of invariant mass histograms
153
154   TList*  fOutputContainer; //!histogram container
155   TH1F*   fHmpi0[AliEMCALGeoParams::fgkEMCALModules][AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows];//! two-cluster inv. mass assigned to each cell.
156
157   TH2F*   fHmgg;            //! two-cluster inv.mass vs pt of pair
158   TH2F*   fHmggDifferentSM; //! two-cluster inv.mass vs pt of pair, each cluster in different SM
159   TH2F*   fHmggSM[4];       //! two-cluster inv.mass per SM
160   TH2F*   fHmggPairSM[4];   //! two-cluster inv.mass per Pair
161   
162   TH2F*   fHOpeningAngle;            //! two-cluster opening angle vs pt of pair, with mass close to pi0
163   TH2F*   fHOpeningAngleDifferentSM; //! two-cluster opening angle vs pt of pair, each cluster in different SM, with mass close to pi0
164   TH2F*   fHOpeningAngleSM[4];       //! two-cluster opening angle vs pt per SM,with mass close to pi0
165   TH2F*   fHOpeningAnglePairSM[4];   //! two-cluster opening angle vs pt per Pair,with mass close to pi0
166
167   TH2F*   fHIncidentAngle;            //! cluster incident angle vs pt of pair, with mass close to pi0
168   TH2F*   fHIncidentAngleDifferentSM; //! cluster incident angle vs pt of pair, each cluster in different SM, with mass close to pi0
169   TH2F*   fHIncidentAngleSM[4];       //! cluster incident angle vs pt per SM,with mass close to pi0
170   TH2F*   fHIncidentAnglePairSM[4];   //! cluster incident angle vs pt per Pair,with mass close to pi0
171   
172   TH2F*   fHAsymmetry;            //! two-cluster asymmetry vs pt of pair, with mass close to pi0
173   TH2F*   fHAsymmetryDifferentSM; //! two-cluster asymmetry vs pt of pair, each cluster in different SM, with mass close to pi0
174   TH2F*   fHAsymmetrySM[4];       //! two-cluster asymmetry vs pt per SM,with mass close to pi0
175   TH2F*   fHAsymmetryPairSM[4];   //! two-cluster asymmetry vs pt per Pair,with mass close to pi0
176   
177   TH2F*   fhTowerDecayPhotonHit[4] ;       //! Cells ordered in column/row for different module, number of times a decay photon hits
178   TH2F*   fhTowerDecayPhotonEnergy[4] ;    //! Cells ordered in column/row for different module, accumulated energy in the tower by decay photons.
179   TH2F*   fhTowerDecayPhotonAsymmetry[4] ; //! Cells ordered in column/row for different module, accumulated asymmetry in the tower by decay photons.
180
181   TH1I*   fhNEvents;        //! Number of events counter histogram
182   TList * fCuts ;           //! List with analysis cuts
183
184   ClassDef(AliAnalysisTaskEMCALPi0CalibSelection,7);
185
186 };
187
188 #endif //ALIANALYSISTASKEMCALPI0CALIBSELECTION_H