]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/EMCALTasks/AliAnalysisTaskEMCALPi0CalibSelection.h
bug-fix and increase in the range of the axis for in cone mc pt sum 2D histos
[u/mrichter/AliRoot.git] / PWGGA / EMCALTasks / AliAnalysisTaskEMCALPi0CalibSelection.h
1 #ifndef ALIANALYSISTASKEMCALPI0CALIBSELECTION_H
2 #define ALIANALYSISTASKEMCALPI0CALIBSELECTION_H
3
4 // $Id$
5
6 // Root includes
7 class TH1F;
8 #include "TH2I.h"
9 #include "TObjArray.h"
10
11 // AliRoot includes
12 #include "AliAnalysisTaskSE.h"
13 class AliEMCALGeometry;
14 #include "AliEMCALGeoParams.h"
15 class AliEMCALRecoUtils;
16
17 class AliAnalysisTaskEMCALPi0CalibSelection : public AliAnalysisTaskSE
18 {
19 public:
20
21   AliAnalysisTaskEMCALPi0CalibSelection(const char* name);
22   virtual ~AliAnalysisTaskEMCALPi0CalibSelection();
23     
24   void    CorrectClusters();
25   
26   void    FillHistograms();
27   
28   void    InitGeometryMatrices();
29   
30   void    InitTemperatureCorrections();
31   
32   void    UserCreateOutputObjects();
33   
34   void    UserExec(Option_t * opt);
35     
36   void    PrintInfo();
37
38   void    Terminate(Option_t* opt);
39   
40   void    GetMaxEnergyCellPosAndClusterPos(AliVCaloCells* cells, AliVCluster* clu, Int_t& iSM, Int_t& ieta, Int_t& iphi);
41   
42   // Analysis parameter setting
43   
44   void    SetPairDTimeCut(Float_t t)                     { fDTimeCut    = t          ; }
45   void    SetClusterMinTime(Float_t tmin)                { fTimeMin     = tmin       ; }
46   void    SetClusterMaxTime(Float_t tmax)                { fTimeMax     = tmax       ; }
47
48   void    SetAsymmetryCut(Float_t asy)                   { fAsyCut      = asy        ; }
49   
50   void    SetClusterMinEnergy(Float_t emin)              { fEmin        = emin       ; }
51   void    SetClusterMaxEnergy(Float_t emax)              { fEmax        = emax       ; }
52   
53   void    SetClusterLambda0Cuts(Float_t min, Float_t max){ fL0max       = max        ;
54                                                            fL0min       = min        ; }
55   void    SetClusterMinNCells(Int_t n)                   { fMinNCells   = n          ; }
56   void    SetNCellsGroup(Int_t n)                        { fGroupNCells = n          ; }
57   void    SetLogWeight(Float_t w)                        { fLogWeight   = w          ; }
58           
59   void    SetPairMinMassCut(Float_t min)                 { fInvMassCutMin = min      ; }
60   void    SetPairMaxMassCut(Float_t max)                 { fInvMassCutMax = max      ; }
61   
62   void    SwitchOnSameSM()                               { fSameSM = kTRUE           ; }
63   void    SwitchOffSameSM()                              { fSameSM = kFALSE          ; }
64   
65   void    UseFilteredEventAsInput()                      { fFilteredInput = kTRUE    ; }
66   void    UseNormalEventAsInput()                        { fFilteredInput = kFALSE   ; }
67   
68   void    SetTriggerName(TString name)                   { fTriggerName = name       ; }
69
70   //Geometry setters
71   
72   void    SetGeometryName(TString name)                  { fEMCALGeoName = name      ; }
73   TString GeometryName() const                           { return fEMCALGeoName      ; }
74   void    SwitchOnLoadOwnGeometryMatrices()              { fLoadMatrices = kTRUE     ; }
75   void    SwitchOffLoadOwnGeometryMatrices()             { fLoadMatrices = kFALSE    ; }
76   void    SetGeometryMatrixInSM(TGeoHMatrix* m, Int_t i) { fMatrix[i]    = m         ; }
77
78   void    SetOADBFilePath(TString path)                  { fOADBFilePath      = path ; }
79   
80   // Cluster recalculation
81   
82   void    SwitchOnClusterCorrection()                    { fCorrectClusters = kTRUE  ; }
83   void    SwitchOffClusterCorrection()                   { fCorrectClusters = kFALSE ; }
84   void    SetEMCALRecoUtils(AliEMCALRecoUtils * ru)      { fRecoUtils = ru           ; }
85   AliEMCALRecoUtils* GetEMCALRecoUtils()    const        { return fRecoUtils         ; }
86   
87   void    SetInvariantMassHistoBinRange(Int_t nBins, Float_t minbin, Float_t maxbin){
88                                         fNbins     = nBins ; fMinBin     = minbin ; fMaxBin     = maxbin ; }
89
90   void    SetTimeHistoBinRange         (Int_t nBins, Float_t minbin, Float_t maxbin){
91                                         fNTimeBins = nBins ; fMinTimeBin = minbin ; fMaxTimeBin = maxbin ; }
92
93   
94   // Mask clusters
95   void    SetNMaskCellColumns(Int_t n) {
96     if(n > fNMaskCellColumns){ delete [] fMaskCellColumns ; fMaskCellColumns = new Int_t[n] ; }
97     fNMaskCellColumns = n ; }
98   
99   void    SetMaskCellColumn(Int_t ipos, Int_t icol) { if(ipos < fNMaskCellColumns) fMaskCellColumns[ipos] = icol            ;
100                                                       else printf("Not set, position larger than allocated set size first") ; }
101   
102   Bool_t  MaskFrameCluster(const Int_t iSM, const Int_t ieta) const;
103   
104 private:
105
106   AliEMCALGeometry  * fEMCALGeo;         //! EMCAL geometry
107   TGeoHMatrix       * fMatrix[AliEMCALGeoParams::fgkEMCALModules]; // Geometry matrices with alignments
108   Bool_t              fLoadMatrices;    //  Matrices set from configuration, not get from geometry.root or from ESDs/AODs
109   
110   TString             fEMCALGeoName;    //  Name of geometry to use.
111   TString             fTriggerName;     //  Trigger name must contain this name
112     
113   AliEMCALRecoUtils * fRecoUtils;       //  Access to reconstruction utilities
114   TString             fOADBFilePath ;   //  Default path $ALICE_ROOT/OADB/EMCAL, if needed change
115   Bool_t              fCorrectClusters; //  Correct clusters energy, position etc.
116
117
118   TRefArray         * fCaloClustersArr; //! list of clusters
119   AliVCaloCells     * fEMCALCells;      //! list of cells
120   
121   TList             * fCuts ;           //! List with analysis cuts
122   TList             * fOutputContainer; //! histogram container
123   Double_t            fVertex[3];       //! primary vertex
124   Bool_t              fFilteredInput;   //  Read input produced with filter.
125
126   // analysis cuts
127   
128   Float_t   fEmin;               // min. cluster energy (GeV)
129   Float_t   fEmax;               // max. cluster energy (GeV)
130   Float_t   fL0min;              // min. cluster L0
131   Float_t   fL0max;              // max. cluster L0
132
133   Float_t   fDTimeCut;           // Maximum difference between time of cluster pairs (ns)
134   Float_t   fTimeMax;            // Maximum cluster time (ns)
135   Float_t   fTimeMin;            // Minimum cluster time (ns)
136
137   Float_t   fAsyCut;             // Asymmetry cut
138   Int_t     fMinNCells;          // min. ncells in cluster
139   Int_t     fGroupNCells;        // group n cells
140   Float_t   fLogWeight;          // log weight used in cluster recalibration
141   Bool_t    fSameSM;             // Combine clusters in channels on same SM
142
143   Int_t     fNMaskCellColumns;   // Number of masked columns
144   Int_t*    fMaskCellColumns;    //[fNMaskCellColumns] list of masked cell collumn
145     
146   Float_t   fInvMassCutMin;      // Min mass cut for clusters to fill time or other histograms
147   Float_t   fInvMassCutMax;      // Mas mass cut for clusters to fill time or other histograms
148   
149   //Output histograms   and settings
150
151   Int_t     fNbins;              // N       mass bins of invariant mass histograms
152   Float_t   fMinBin;             // Minimum mass bins of invariant mass histograms
153   Float_t   fMaxBin;             // Maximum mass bins of invariant mass histograms
154   
155   Int_t     fNTimeBins;          // N       time bins of invariant mass histograms
156   Float_t   fMinTimeBin;         // Minimum time bins of invariant mass histograms
157   Float_t   fMaxTimeBin;         // Maximum time bins of invariant mass histograms
158   
159   TH1F*     fHmpi0[AliEMCALGeoParams::fgkEMCALModules][AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows];//! two-cluster inv. mass assigned to each cell.
160
161   TH2F*     fHmgg;                                                                 //! two-cluster inv.mass vs pt of pair
162   TH2F*     fHmggDifferentSM;                                                      //! two-cluster inv.mass vs pt of pair, each cluster in different SM
163   TH2F*     fHmggSM[AliEMCALGeoParams::fgkEMCALModules];                           //! two-cluster inv.mass per SM
164   TH2F*     fHmggPairSameSectorSM[AliEMCALGeoParams::fgkEMCALModules/2];           //! two-cluster inv.mass per Pair
165   TH2F*     fHmggPairSameSideSM  [AliEMCALGeoParams::fgkEMCALModules-2];           //! two-cluster inv.mass per Pair
166   
167   TH2F*     fHmggMaskFrame;                                                        //! two-cluster inv.mass vs pt of pair, mask clusters facing frames
168   TH2F*     fHmggDifferentSMMaskFrame;                                             //! two-cluster inv.mass vs pt of pair, each cluster in different SM,mask clusters facing frames
169   TH2F*     fHmggSMMaskFrame[AliEMCALGeoParams::fgkEMCALModules];                  //! two-cluster inv.mass per SM, mask clusters facing frames
170   TH2F*     fHmggPairSameSectorSMMaskFrame[AliEMCALGeoParams::fgkEMCALModules/2];  //! two-cluster inv.mass per Pair, mask clusters facing frames
171   TH2F*     fHmggPairSameSideSMMaskFrame  [AliEMCALGeoParams::fgkEMCALModules-2];  //! two-cluster inv.mass per Pair, mask clusters facing frames
172
173   TH2F*     fHOpeningAngle;                                                        //! two-cluster opening angle vs pt of pair, with mass close to pi0
174   TH2F*     fHOpeningAngleDifferentSM;                                             //! two-cluster opening angle vs pt of pair, each cluster in different SM, with mass close to pi0
175   TH2F*     fHOpeningAngleSM[AliEMCALGeoParams::fgkEMCALModules];                  //! two-cluster opening angle vs pt per SM,with mass close to pi0
176   TH2F*     fHOpeningAnglePairSM[AliEMCALGeoParams::fgkEMCALModules];              //! two-cluster opening angle vs pt per Pair,with mass close to pi0
177
178   TH2F*     fHAsymmetry;                                                           //! two-cluster asymmetry vs pt of pair, with mass close to pi0
179   TH2F*     fHAsymmetryDifferentSM;                                                //! two-cluster asymmetry vs pt of pair, each cluster in different SM, with mass close to pi0
180   TH2F*     fHAsymmetrySM[AliEMCALGeoParams::fgkEMCALModules];                     //! two-cluster asymmetry vs pt per SM,with mass close to pi0
181   TH2F*     fHAsymmetryPairSM[AliEMCALGeoParams::fgkEMCALModules];                 //! two-cluster asymmetry vs pt per Pair,with mass close to pi0
182   
183   TH2F*     fhTowerDecayPhotonHit[AliEMCALGeoParams::fgkEMCALModules] ;            //! Cells ordered in column/row for different module, number of times a decay photon hits
184   TH2F*     fhTowerDecayPhotonEnergy[AliEMCALGeoParams::fgkEMCALModules] ;         //! Cells ordered in column/row for different module, accumulated energy in the tower by decay photons.
185   TH2F*     fhTowerDecayPhotonAsymmetry[AliEMCALGeoParams::fgkEMCALModules] ;      //! Cells ordered in column/row for different module, accumulated asymmetry in the tower by decay photons.
186   TH2F*     fhTowerDecayPhotonHitMaskFrame[AliEMCALGeoParams::fgkEMCALModules] ;   //! Cells ordered in column/row for different module, number of times a decay photon hits
187
188   TH1I*     fhNEvents;                                                             //! Number of events counter histogram
189  
190   //Time
191   TH2F*     fHTpi0[4];                                                             //! Time of cell under pi0 mass, for 4 bunch crossings
192   TH2F*     fhClusterTime ;                                                        //! Timing of clusters vs energy
193   TH2F*     fhClusterTimeSM[AliEMCALGeoParams::fgkEMCALModules] ;                  //! Timing of clusters vs energy per SM
194   TH2F*     fhClusterPairDiffTime;                                                 //! Diference in time of clusters
195   TH2F*     fhClusterPairDiffTimeSameSM[AliEMCALGeoParams::fgkEMCALModules];       //! Diference in time of clusters same SM
196   TH2F*     fhClusterPairDiffTimeSameSector[AliEMCALGeoParams::fgkEMCALModules/2]; //! Diference in time of clusters same sector
197   TH2F*     fhClusterPairDiffTimeSameSide[AliEMCALGeoParams::fgkEMCALModules-2];   //! Diference in time of clusters same side
198
199   AliAnalysisTaskEMCALPi0CalibSelection(           const AliAnalysisTaskEMCALPi0CalibSelection&); 
200   AliAnalysisTaskEMCALPi0CalibSelection& operator=(const AliAnalysisTaskEMCALPi0CalibSelection&); 
201   
202   ClassDef(AliAnalysisTaskEMCALPi0CalibSelection,20);
203
204 };
205
206 #endif //ALIANALYSISTASKEMCALPI0CALIBSELECTION_H