]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalJetHadEPpid.h
Merge branch 'master' into TPCdev
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskEmcalJetHadEPpid.h
1 #ifndef AliAnalysisTaskEmcalJetHadEPpid_h
2 #define AliAnalysisTaskEmcalJetHadEPpid_h
3
4 // root classes
5 class TClonesArray;
6 class TH1F;
7 class TH2F;
8 class TH3F;
9 class THnSparse;
10 class TList;
11 class TLorentzVector;
12 class TGraph;
13
14 // AliROOT classes
15 class AliEventPoolManager;
16 class AliLocalRhoParameter;
17 class AliEMCALTrack;
18 class AliMagF;
19 class AliESDEvent;
20 class AliAODEvent;
21 class AliEMCALGeometry;
22 class AliEMCALRecoUtils;
23 class AliESDtrack;
24
25 // container classes
26 class AliJetContainer;
27 class AliParticleContainer;
28 class AliClusterContainer;
29
30 // includes
31 #include <AliAnalysisTaskEmcalJet.h>
32 #include <AliEmcalJet.h>
33 #include <AliVEvent.h>
34 #include <AliVTrack.h>
35 #include <AliVCluster.h>
36 #include <TClonesArray.h>
37 #include <TMath.h>
38 #include <TRandom3.h>
39 #include <AliLog.h>
40
41 // Local Rho includes
42 #include "AliAnalysisTaskLocalRho.h"
43 #include "AliLocalRhoParameter.h"
44
45 // PID includes
46 #include "AliPIDResponse.h"
47
48 #include "AliAnalysisFilter.h"
49
50 class AliAnalysisTaskEmcalJetHadEPpid : public AliAnalysisTaskEmcalJet {
51  public:
52   AliAnalysisTaskEmcalJetHadEPpid();
53   AliAnalysisTaskEmcalJetHadEPpid(const char *name);
54   //virtual ~AliAnalysisTaskEmcalJetHadEPpid() {}
55   virtual ~AliAnalysisTaskEmcalJetHadEPpid();
56
57   virtual void            UserCreateOutputObjects();
58   virtual THnSparse*      NewTHnSparseF(const char* name, UInt_t entries);
59   virtual void            GetDimParams(Int_t iEntry,TString &label, Int_t &nbins, Double_t &xmin, Double_t &xmax);
60   virtual THnSparse*      NewTHnSparseFPID(const char* name, UInt_t entries);
61   virtual void            GetDimParamsPID(Int_t iEntry,TString &label, Int_t &nbins, Double_t &xmin, Double_t &xmax);
62   // set a bun of histogram switches up
63   void                    SetPlotGlobalRho(Bool_t g)            { doPlotGlobalRho = g; } // plot global rho switch
64   void                    SetVariableBinning(Bool_t v)          { doVariableBinning = v; } // do variable binning switch
65   void                        SetvarbinTHnSparse(Bool_t vb)         { dovarbinTHnSparse = vb; } // variable THnSparse bin switch
66   void                                    SetmakeQAhistos(Bool_t QAhist)        { makeQAhistos = QAhist; } // make QA histos  
67   void                                    SetmakeBIAShistos(Bool_t BIAShist)    { makeBIAShistos = BIAShist; } // make bias histos
68   void                            SetmakeextraCORRhistos(Bool_t Xhist)  { makeextraCORRhistos = Xhist; } // make extra correlations histos
69
70   // set data, detectors type, and PID and PID w bias switches
71   void                        SetDataType(Bool_t data)                  { useAOD = data; }    // data type switch
72   void                                    SetcutType(TString cut)                               { fcutType = cut; }    // EMCAL / TPC acceptance cut
73   void                    SetdoPID(Bool_t p)                    { doPID = p; }   // do PID switch
74   void                                    SetdoPIDtrackBIAS(Bool_t PIDbias)     { doPIDtrackBIAS = PIDbias; } // do PID track bias switch
75
76   // give comments setter
77   void                                    SetdoComments(Bool_t comm)                    { doComments = comm; } // give comment switch
78
79   // getters
80   TString                 GetLocalRhoName() const               {return fLocalRhoName; }
81
82   // set names of some objects
83   virtual void            SetLocalRhoName(const char *ln)       { fLocalRhoName = ln; }
84   virtual void            SetTracksName(const char *tn)         { fTracksName = tn; }
85   virtual void            SetJetsName(const char *jn)           { fJetsName = jn; }
86
87   // bias and cuts - setters
88   virtual void            SetAreaCut(Double_t a)                { fAreacut    = a; }
89   virtual void            SetTrkBias(Double_t b)                { fTrkBias    = b; }  //require a track with pt > b in jet
90   virtual void            SetClusBias(Double_t b)               { fClusBias   = b; }  //require a cluster with pt > b in jet
91   virtual void            SetTrkEta(Double_t e)                 { fTrkEta   = e; }  //eta range of the associated tracks
92   virtual void            SetJetPtcut(Double_t jpt)             { fJetPtcut = jpt; } // jet pt cut
93   virtual void                    SetJetRad(Double_t jrad)                              { fJetRad = jrad; } // jet radius 
94   virtual void                    SetConstituentCut(Double_t constCut)   { fConstituentCut = constCut; } // constituent Cut
95
96   // eta and phi limits of jets - setters
97   virtual void            SetJetEta(Double_t emin, Double_t emax)  { fEtamin = emin; fEtamax = emax; }
98   virtual void            SetJetPhi(Double_t pmin, Double_t pmax)  { fPhimin = pmin; fPhimax = pmax; }
99
100   // event mixing - setters
101   virtual void            SetEventMixing(Int_t yesno)              { fDoEventMixing=yesno; }
102   virtual void            SetMixingTracks(Int_t tracks)            { fMixingTracks = tracks; }
103
104   // jet container - setters
105   void SetContainerAllJets(Int_t c)         { fContainerAllJets      = c;}
106   void SetContainerPIDJets(Int_t c)         { fContainerPIDJets      = c;}
107
108 protected:
109   // functions 
110   void                                   ExecOnce();
111   Bool_t                         Run();
112   virtual void           Terminate(Option_t *); 
113   virtual Int_t          AcceptMyJet(AliEmcalJet *jet);   // applies basic jet tests/cuts before accepting
114   virtual Int_t          GetCentBin(Double_t cent) const; // centrality bin of event
115   Double_t                RelativePhi(Double_t mphi,Double_t vphi) const; // relative jet track angle
116   Double_t                RelativeEPJET(Double_t jetAng, Double_t EPAng) const;  // relative jet event plane angle
117   virtual Int_t          GetEtaBin(Double_t eta) const;      // eta bins
118   virtual Int_t          GetpTjetBin(Double_t pt) const;     // jet pt bins
119   virtual Int_t          GetpTtrackBin(Double_t pt) const;   // track pt bins
120   virtual Int_t          GetzVertexBin(Double_t zVtx) const; // zVertex bin
121   void                   SetfHistPIDcounterLabels(TH1* fHistPID) const;  // PID counter
122
123   // parameters of detector to cut on for event
124   Double_t               fPhimin;                  // phi min
125   Double_t               fPhimax;                  // phi max
126   Double_t               fEtamin;                  // eta min
127   Double_t               fEtamax;                  // eta max
128   Double_t               fAreacut;                 // area cut
129   Double_t               fTrkBias;                 // track bias
130   Double_t               fClusBias;                // cluster bias
131   Double_t               fTrkEta;                  // eta min/max of tracks
132   Double_t                   fJetPtcut;                    // jet pt to cut on for correlations
133   Double_t                               fJetRad;                                  // jet radius
134   Double_t                               fConstituentCut;          // jet constituent cut
135
136   // event mixing
137   Int_t                  fDoEventMixing;
138   Int_t                  fMixingTracks;
139
140   // switches for plots
141   Bool_t                 doPlotGlobalRho;
142   Bool_t                 doVariableBinning;
143   Bool_t         dovarbinTHnSparse;
144   Bool_t                 makeQAhistos;
145   Bool_t                 makeBIAShistos;
146   Bool_t                 makeextraCORRhistos; 
147
148   // data type switch
149   Bool_t             useAOD;
150
151   // Cut type (EMCAL/TPC acceptance)
152   TString        fcutType;
153
154   // switches for PID
155   Bool_t                 doPID;
156   Bool_t                 doPIDtrackBIAS;
157
158   // do comment switch
159   Bool_t                 doComments;
160
161   // local rho value
162   Double_t               fLocalRhoVal;
163
164   // object names
165   TString                fTracksName;
166   TString                fJetsName;
167
168   // event counter
169   Int_t                  event;
170
171   // boolean functions for PID
172   Bool_t                 isPItpc, isKtpc, isPtpc;
173   Bool_t                 isPIits, isKits, isPits;
174   Bool_t                 isPItof, isKtof, isPtof;
175
176   // event pool
177   TObjArray*                CloneAndReduceTrackList(TObjArray* tracks);
178   AliEventPoolManager   *fPoolMgr;  // event pool Manager object
179
180   // PID
181   AliPIDResponse        *fPIDResponse;   // PID response object
182   AliTPCPIDResponse     *fTPCResponse;   // TPC pid response object
183
184  private:
185   // needed for PID, track objects
186   AliESDEvent       *fESD;          // ESD object
187   AliAODEvent       *fAOD;                // AOD object
188
189   TH2F                  *fHistTPCdEdX;
190   TH2F                  *fHistITSsignal;
191 //  TH2F                    *fHistTOFsignal;
192
193   TH2F                  *fHistRhovsCent; //!
194   TH2F                  *fHistNjetvsCent;//! number of jets versus Centrality
195   TH2F                  *fHistJetPtvsTrackPt[6];//!
196   TH2F                  *fHistRawJetPtvsTrackPt[6];//!
197   TH1F                  *fHistTrackPt[6];//!
198   TH1F                  *fHistEP0[6];//!
199   TH1F                  *fHistEP0A[6];//!
200   TH1F                  *fHistEP0C[6];//!
201   TH2F                  *fHistEPAvsC[6];//!
202   TH1F                                  *fHistJetPtcorrGlRho[6];//!
203   TH2F                  *fHistJetPtvsdEP[6];//!
204   TH2F                  *fHistJetPtvsdEPBias[6];//!
205   TH2F                  *fHistRhovsdEP[6]; //!
206   TH3F                  *fHistJetEtaPhiPt[6];//!
207   TH3F                  *fHistJetEtaPhiPtBias[6];//!
208   TH2F                  *fHistJetPtArea[6];//!
209   TH2F                  *fHistJetPtAreaBias[6];//!
210   TH2F                  *fHistJetPtNcon[6]; //!
211   TH2F                  *fHistJetPtNconBias[6]; //!
212   TH2F                  *fHistJetPtNconCh[6]; //!
213   TH2F                  *fHistJetPtNconBiasCh[6]; //!
214   TH2F                  *fHistJetPtNconEm[6]; //!
215   TH2F                  *fHistJetPtNconBiasEm[6]; //!
216   TH1F                      *fHistJetHaddPhiINcent[6];
217   TH1F                          *fHistJetHaddPhiOUTcent[6];
218   TH1F                          *fHistJetHaddPhiMIDcent[6];
219
220   TH1                   *fHistCentrality;
221   TH1                   *fHistZvtx;
222   TH1                   *fHistMult;
223   TH1                           *fHistJetPhi;
224   TH1                       *fHistTrackPhi;
225   TH1                       *fHistJetHaddPhiIN;
226   TH1                           *fHistJetHaddPhiOUT;
227   TH1                           *fHistJetHaddPhiMID;
228   TH1                                   *fHistJetHaddPhiBias;
229   TH1                                   *fHistJetHaddPhiINBias;
230   TH1                                   *fHistJetHaddPhiOUTBias;
231   TH1                                   *fHistJetHaddPhiMIDBias;
232
233   TH1                   *fHistMEdPHI; // phi distrubtion of mixed events
234   TH1                                   *fHistTrackPtallcent;
235
236   TH2                   *fHistJetEtaPhi;  
237   TH2                   *fHistTrackEtaPhi[4][7];
238   TH1                           *fHistJetHadbindPhi[9]; 
239   TH1                                   *fHistJetHadbindPhiIN[9]; 
240   TH1                                   *fHistJetHadbindPhiMID[9]; 
241   TH1                           *fHistJetHadbindPhiOUT[9]; 
242   TH2                   *fHistJetHEtaPhi;
243
244   TH1                   *fHistJetPt[6];
245   TH1                   *fHistJetPtBias[6];
246   TH1                   *fHistJetPtTT[6];
247   TH2                   *fHistAreavsRawPt[6];
248   TH2                   *fHistJetH[2][5][3];
249   TH2                   *fHistJetHBias[2][5][3];
250   TH2                   *fHistJetHTT[2][5][3];
251   TH1F                                  *fHistJetHdPHI[11];
252   TH2F                                  *fHistJetHdETAdPHI[11];
253   TH2F                  *fHistSEphieta; // single events phi-eta distributions
254   TH2F                  *fHistMEphieta; // mixed events phi-eta distributions
255   TH1F                                  *fHistJetHaddPHI;
256
257   // PID status histo's
258   TH1                                   *fHistPID;
259
260   // THn Sparse's
261   THnSparse             *fhnPID;          // PID sparse
262   THnSparse             *fhnMixedEvents;  // mixed events matrix
263   THnSparse             *fhnJH;           // jet hadron events matrix
264
265   // container objects
266   AliJetContainer            *fJetsCont;                   //!Jets
267   AliParticleContainer       *fTracksCont;                 //!Tracks
268   AliClusterContainer        *fCaloClustersCont;           //!Clusters
269
270   // container specifier
271   Int_t                                 fContainerAllJets;  // number of container with all full jets
272   Int_t                                 fContainerPIDJets;  // number of container with full jets meeting Pt cut (for PID)
273
274 // ***********************************************************
275    
276   //Declare it private to avoid compilation warning
277   AliAnalysisTaskEmcalJetHadEPpid(const AliAnalysisTaskEmcalJetHadEPpid & g) ; // cpy ctor
278
279   AliAnalysisTaskEmcalJetHadEPpid& operator=(const AliAnalysisTaskEmcalJetHadEPpid&); // not implemented
280   ClassDef(AliAnalysisTaskEmcalJetHadEPpid, 4); // Emcal jet hadron PID - Event plane dependence
281 };
282 #endif