]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/EMCALTasks/AliAnalysisTaskEMCALPi0PbPb.h
remove obsolete restriction on aplication of time cuts in case of AOD analysis
[u/mrichter/AliRoot.git] / PWGGA / EMCALTasks / AliAnalysisTaskEMCALPi0PbPb.h
1 #ifndef AliAnalysisTaskEMCALPi0PbPb_h
2 #define AliAnalysisTaskEMCALPi0PbPb_h
3
4 // $Id$
5
6 class TAxis;
7 class TClonesArray;
8 class TH1;
9 class TH2;
10 class TNtuple;
11 class TObjArray;
12 class AliAODCaloCells;
13 class AliAODCaloCluster;
14 class AliAODEvent;
15 class AliAODTrack;
16 class AliAODVertex;
17 class AliEMCALGeometry;
18 class AliEMCALRecoUtils;
19 class AliESDCaloCells;
20 class AliESDCaloCluster;
21 class AliESDEvent;
22 class AliESDTrack;
23 class AliESDVertex;
24 class AliESDtrackCuts;
25 class AliMCEvent;
26 class AliMCParticle;
27 class AliStaHeader;
28 class AliStaVertex;
29
30 #include "AliAnalysisTaskSE.h"
31 #include "AliStaObjects.h"
32
33 class AliAnalysisTaskEMCALPi0PbPb : public AliAnalysisTaskSE {
34  public:
35   AliAnalysisTaskEMCALPi0PbPb();
36   AliAnalysisTaskEMCALPi0PbPb(const char *name);
37   virtual ~AliAnalysisTaskEMCALPi0PbPb(); 
38   
39   void         UserCreateOutputObjects();
40   void         UserExec(Option_t *option);
41   void         Terminate(Option_t *);
42
43   void         SetAsymMax(Double_t asymMax)                   { fAsymMax       = asymMax;   }
44   void         SetCentrality(const char *n)                   { fCentVar       = n;         }
45   void         SetCentralityRange(Double_t from, Double_t to) { fCentFrom=from; fCentTo=to; }
46   void         SetClusName(const char *n)                     { fClusName      = n;         }
47   void         SetDoAfterburner(Bool_t b)                     { fDoAfterburner = b;         }
48   void         SetDoPhysicsSelection(Bool_t b)                { fDoPSel        = b;         }
49   void         SetDoTrackMatWithGeom(Bool_t b)                { fDoTrMatGeom   = b;         }
50   void         SetEmbedMode(Bool_t b)                         { fEmbedMode     = b;         }
51   void         SetFillNtuple(Bool_t b)                        { fDoNtuple      = b;         }
52   void         SetGeoName(const char *n)                      { fGeoName       = n;         }
53   void         SetGeoUtils(AliEMCALGeometry *geo)             { fGeom          = geo;       }
54   void         SetIsoDist(Double_t d)                         { fIsoDist       = d;         }
55   void         SetL0TimeRange(Int_t l, Int_t h)               { fMinL0Time=l; fMaxL0Time=h; }
56   void         SetMarkCells(const char *n)                    { fMarkCells     = n;         }
57   void         SetMcMode(Bool_t b)                            { fMcMode        = b;         }
58   void         SetMinClusEnergy(Double_t e)                   { fMinE          = e;         }
59   void         SetMinEcc(Double_t ecc)                        { fMinEcc        = ecc;       }
60   void         SetMinErat(Double_t erat)                      { fMinErat       = erat;      }
61   void         SetMinNClustersPerTrack(Double_t m)            { fMinNClusPerTr = m;         }
62   void         SetNminCells(Int_t n)                          { fNminCells     = n;         }
63   void         SetPrimTrackCuts(AliESDtrackCuts *c)           { fPrimTrCuts    = c;         }
64   void         SetPrimTracksName(const char *n)               { fPrimTracksName = n;        }
65   void         SetRecoUtils(AliEMCALRecoUtils *reco)          { fReco          = reco;      }
66   void         SetTrClassNames(const char *n)                 { fTrClassNames  = n;         }
67   void         SetTrackCuts(AliESDtrackCuts *c)               { fTrCuts        = c;         }
68   void         SetTrainMode(Bool_t b)                         { fTrainMode     = b;         }
69   void         SetTrigName(const char *n)                     { fTrigName      = n;         }
70   void         SetUseQualFlag(Bool_t b)                       { fUseQualFlag   = b;         }
71   void         SetVertexRange(Double_t z1, Double_t z2)       { fVtxZMin=z1; fVtxZMax=z2;   }
72
73  protected:
74   virtual void CalcCaloTriggers();
75   virtual void CalcClusterProps();
76   virtual void CalcPrimTracks();
77   virtual void CalcMcInfo();
78   virtual void CalcTracks();
79   virtual void ClusterAfterburner();
80   virtual void FillCellHists();
81   virtual void FillClusHists();
82   virtual void FillNtuple();
83   virtual void FillOtherHists();
84   virtual void FillPionHists();
85   virtual void FillMcHists();
86   virtual void FillTrackHists();
87   void         FillVertex(AliStaVertex *v, const AliESDVertex *esdv);
88   void         FillVertex(AliStaVertex *v, const AliAODVertex *aodv);
89   Double_t     GetCellIsolation(Double_t cEta, Double_t cPhi, Double_t radius=0.2)                        const;
90   Double_t     GetCellIsoNxM(Double_t cEta, Double_t cPhi, Int_t N, Int_t M)                              const;
91   Double_t     GetCellEnergy(const AliVCluster *c)    const;
92   Double_t     GetMaxCellEnergy(const AliVCluster *c) const { Short_t id=-1; return GetMaxCellEnergy(c,id); }
93   Double_t     GetMaxCellEnergy(const AliVCluster *c, Short_t &id)                                        const;
94   Double_t     GetSecondMaxCellEnergy(AliVCluster *clus, Short_t &id)                                     const;
95   Int_t        GetNCells(const AliVCluster *c, Double_t emin=0.)                                          const;
96   Int_t        GetNCells(Int_t sm, Double_t emin=0.)                                                      const;
97   void         GetSigma(const AliVCluster *c, Double_t &sigmaMax, Double_t &sigmaMin)                     const;
98   void         GetSigmaEtaEta(const AliVCluster *c, Double_t &sigmaEtaEta, Double_t &sigmaPhiPhi)         const;
99   Double_t     GetTrackIsolation(Double_t cEta, Double_t cPhi, Double_t radius=0.2, Double_t pt=0.)       const;
100   Double_t     GetTrackIsoStrip(Double_t cEta, Double_t cPhi, Double_t dEta=0.015, Double_t dPhi=0.3, Double_t pt=0.)       const;
101   Bool_t       IsShared(const AliVCluster *c)                                                             const;
102   Bool_t       IsIdPartOfCluster(const AliVCluster *c, Short_t id)                                        const;
103   void         PrintDaughters(const AliVParticle *p, const TObjArray *arr, Int_t level=0)                 const;
104   void         PrintDaughters(const AliMCParticle *p, const AliMCEvent *arr, Int_t level=0)               const;
105   void         PrintTrackRefs(AliMCParticle *p)                                                           const;
106   void         ProcessDaughters(AliVParticle *p, Int_t index, const TObjArray *arr);
107   void         ProcessDaughters(AliMCParticle *p, Int_t index, const AliMCEvent *arr);
108
109     // input members
110   TString                fCentVar;                // variable for centrality determination
111   Double_t               fCentFrom;               // min centrality (def=0)
112   Double_t               fCentTo;                 // max centrality (def=100)
113   Double_t               fVtxZMin;                // min primary vertex z (def=-10cm)
114   Double_t               fVtxZMax;                // max primary vertex z (def=+10cm)
115   Bool_t                 fUseQualFlag;            // if true use quality flag for centrality
116   TString                fClusName;               // cluster branch name (def="")
117   Bool_t                 fDoNtuple;               // if true write out ntuple
118   Bool_t                 fDoAfterburner;          // if true run after burner
119   Double_t               fAsymMax;                // maximum energy asymmetry (def=1)
120   Int_t                  fNminCells;              // minimum number of cells attached to cluster (def=1)
121   Double_t               fMinE;                   // minimum cluster energy (def=0.1 GeV/c)
122   Double_t               fMinErat;                // minimum emax/ec ratio (def=0)
123   Double_t               fMinEcc;                 // minimum eccentricity (def=0)
124   TString                fGeoName;                // geometry name (def = EMCAL_FIRSTYEARV1)
125   Double_t               fMinNClusPerTr;          // minimum number of cluster per track (def=50)
126   Double_t               fIsoDist;                // isolation distance (def=0.2)
127   TString                fTrClassNames;           // trigger class names
128   AliESDtrackCuts       *fTrCuts;                 // track cuts
129   AliESDtrackCuts       *fPrimTrCuts;             // track cuts
130   TString                fPrimTracksName;         // name of track collection (if "" use branch)
131   Bool_t                 fDoTrMatGeom;            // track matching including geometry
132   Bool_t                 fTrainMode;              // train mode with minimal number of resources
133   TString                fMarkCells;              // list of mark cells to monitor
134   Int_t                  fMinL0Time;              // minimum accepted time for trigger
135   Int_t                  fMaxL0Time;              // maximum accepted time for trigger
136   Bool_t                 fMcMode;                 // monte carlo mode
137   Bool_t                 fEmbedMode;              // embedding mode
138   AliEMCALGeometry      *fGeom;                   // geometry utils
139   AliEMCALRecoUtils     *fReco;                   // reco utils
140   TString                fTrigName;               // trigger name
141   Bool_t                 fDoPSel;                 // if false then accept all events
142     // derived members (ie with ! after //)
143   Bool_t                 fIsGeoMatsSet;           //!indicate that geo matrices are set 
144   ULong64_t              fNEvs;                   //!accepted events 
145   TList                 *fOutput;                 //!container of output histograms
146   TObjArray             *fTrClassNamesArr;        //!array of trig class names  
147   AliESDEvent           *fEsdEv;                  //!pointer to input esd event
148   AliAODEvent           *fAodEv;                  //!pointer to input aod event
149   const TObjArray       *fRecPoints;              //!pointer to rec points (AliAnalysisTaskEMCALClusterizeFast)
150   const TClonesArray    *fDigits;                 //!pointer to digits     (AliAnalysisTaskEMCALClusterizeFast)
151   TObjArray             *fEsdClusters;            //!pointer to esd clusters
152   AliESDCaloCells       *fEsdCells;               //!pointer to esd cells
153   TObjArray             *fAodClusters;            //!pointer to aod clusters
154   AliAODCaloCells       *fAodCells;               //!pointer to aod cells
155   TAxis                 *fPtRanges;               //!pointer to pt ranges
156   TObjArray             *fSelTracks;              //!pointer to selected tracks
157   TObjArray             *fSelPrimTracks;          //!pointer to selected primary tracks
158     // ntuple
159   TTree                 *fNtuple;                 //!pointer to ntuple
160   AliStaHeader          *fHeader;                 //!pointer to header
161   AliStaVertex          *fPrimVert;               //!pointer to primary vertex
162   AliStaVertex          *fSpdVert;                //!pointer to SPD vertex
163   AliStaVertex          *fTpcVert;                //!pointer to TPC vertex
164   TClonesArray          *fClusters;               //!pointer to clusters
165   TClonesArray          *fTriggers;               //!pointer to triggers
166   TClonesArray          *fMcParts;                //!pointer to mc particles
167     // histograms
168   TH1                   *fHCuts;                  //!histo for cuts
169   TH1                   *fHVertexZ;               //!histo for vtxz
170   TH1                   *fHVertexZ2;              //!histo for vtxz after vtx cuts
171   TH1                   *fHCent;                  //!histo for cent
172   TH1                   *fHCentQual;              //!histo for cent after quality flag cut
173   TH1                   *fHTclsBeforeCuts;        //!histo for trigger classes before cuts
174   TH1                   *fHTclsAfterCuts;         //!histo for trigger classes after cuts
175
176     // histograms for cells
177   TH2                  **fHColuRow;               //!histo for cell column and row
178   TH2                  **fHColuRowE;              //!histo for cell column and row weight energy
179   TH1                  **fHCellMult;              //!histo for cell multiplicity in module
180   TH1                   *fHCellE;                 //!histo for cell energy
181   TH1                   *fHCellH;                 //!histo for highest cell energy
182   TH1                   *fHCellM;                 //!histo for mean cell energy (normalized to hit cells)
183   TH1                   *fHCellM2;                //!histo for mean cell energy (normalized to all cells)
184   TH1                  **fHCellFreqNoCut;         //!histo for cell frequency without cut
185   TH1                  **fHCellFreqCut100M;       //!histo for cell frequency with cut 100MeV
186   TH1                  **fHCellFreqCut300M;       //!histo for cell frequency with cut 300MeV
187   TH1                  **fHCellFreqE;             //!histo for cell frequency weighted with energy
188   TH1                  **fHCellCheckE;            //!histo for cell E distribution for given channels
189     // histograms for clusters
190   TH1                   *fHClustEccentricity;     //!histo for cluster eccentricity
191   TH2                   *fHClustEtaPhi;           //!histo for cluster eta vs. phi
192   TH2                   *fHClustEnergyPt;         //!histo for cluster energy vs. pT
193   TH2                   *fHClustEnergySigma;      //!histo for cluster energy vs. variance over long axis 
194   TH2                   *fHClustSigmaSigma;       //!histo for sigma vs. lambda_0 comparison
195   TH2                   *fHClustNCellEnergyRatio; //!histo for cluster n cells vs. energy ratio
196   TH2                   *fHClustEnergyNCell;      //!histo for cluster energy vs. cluster n cells
197     // histograms for primary tracks
198   TH1                   *fHPrimTrackPt;           //!histo for primary track pt
199   TH1                   *fHPrimTrackEta;          //!histo for primary track eta
200   TH1                   *fHPrimTrackPhi;           //!histo for primary track phi
201     // histograms for track matching
202   TH1                   *fHMatchDr;               //!histo for dR track cluster matching
203   TH1                   *fHMatchDz;               //!histo for dZ track cluster matching
204   TH1                   *fHMatchEp;               //!histo for E/p track cluster matching
205     // histograms for pion candidates
206   TH2                   *fHPionEtaPhi;            //!histo for pion eta vs. phi
207   TH2                   *fHPionMggPt;             //!histo for pion mass vs. pT
208   TH2                   *fHPionMggAsym;           //!histo for pion mass vs. asym
209   TH2                   *fHPionMggDgg;            //!histo for pion mass vs. opening angle
210   TH1                   *fHPionInvMasses[21];     //!histos for invariant mass plots 
211     // histograms for MC
212
213  private:
214   AliAnalysisTaskEMCALPi0PbPb(const AliAnalysisTaskEMCALPi0PbPb&);            // not implemented
215   AliAnalysisTaskEMCALPi0PbPb &operator=(const AliAnalysisTaskEMCALPi0PbPb&); // not implemented
216
217   ClassDef(AliAnalysisTaskEMCALPi0PbPb, 13) // Analysis task for neutral pions in Pb+Pb
218 };
219 #endif