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