]>
Commit | Line | Data |
---|---|---|
1101e468 | 1 | // To make jet hadron correlations |
2 | // Author: Megan Connors | |
3 | ||
4 | #include "TChain.h" | |
5 | #include "TTree.h" | |
6 | #include "TList.h" | |
7 | #include "TH1F.h" | |
8 | #include "TH2F.h" | |
9 | #include "THnSparse.h" | |
10 | #include "TCanvas.h" | |
11 | #include <TClonesArray.h> | |
12 | #include <TParticle.h> | |
13 | #include "AliVTrack.h" | |
14 | #include "TParameter.h" | |
15 | ||
16 | #include "AliAODEvent.h" | |
17 | #include "AliAnalysisTask.h" | |
18 | #include "AliAnalysisManager.h" | |
19 | ||
20 | #include "AliESDEvent.h" | |
21 | #include "AliESDInputHandler.h" | |
22 | #include "AliESDCaloCluster.h" | |
23 | #include "AliESDVertex.h" | |
24 | #include "AliCentrality.h" | |
25 | #include "AliAODJet.h" | |
26 | #include "AliEmcalJet.h" | |
27 | #include "AliESDtrackCuts.h" | |
28 | #include "AliAnalysisTaskEmcalJetHMEC.h" | |
29 | #include "TVector3.h" | |
30 | ||
31 | ClassImp(AliAnalysisTaskEmcalJetHMEC) | |
32 | ||
33 | //________________________________________________________________________ | |
34 | AliAnalysisTaskEmcalJetHMEC::AliAnalysisTaskEmcalJetHMEC() : | |
35 | AliAnalysisTaskSE(), | |
36 | fTracksName("tracks"), | |
37 | fJetsName("jets"), | |
38 | fPhimin(-10), | |
39 | fPhimax(10), | |
40 | fEtamin(-0.9), | |
41 | fEtamax(0.9), | |
42 | fAreacut(0.0), | |
43 | fESD(0), | |
44 | fOutputList(0), | |
45 | fHistTrackPt(0), | |
46 | fHistCentrality(0), | |
47 | fHistJetEtaPhi(0), | |
48 | fHistTrackEtaPhi(0), | |
49 | fHistJetHEtaPhi(0) | |
50 | { | |
51 | // Default Constructor | |
52 | ||
53 | for(Int_t icent = 0; icent<6; ++icent){ | |
54 | fHistJetPt[icent]=0; | |
55 | fHistJetPt[icent]=0; | |
56 | for(Int_t iptjet = 0; iptjet<3; ++iptjet){ | |
57 | for(Int_t ieta = 0; ieta<3; ++ieta){ | |
58 | fHistJetH[icent][iptjet][ieta]=0; | |
59 | fHistJetHBias[icent][iptjet][ieta]=0; | |
60 | } | |
61 | } | |
62 | } | |
63 | ||
64 | } | |
65 | //________________________________________________________________________ | |
66 | AliAnalysisTaskEmcalJetHMEC::AliAnalysisTaskEmcalJetHMEC(const char *name) : | |
67 | AliAnalysisTaskSE(name), | |
68 | fTracksName("tracks"), | |
69 | fJetsName("jets"), | |
70 | fPhimin(-10), | |
71 | fPhimax(10), | |
72 | fEtamin(-0.9), | |
73 | fEtamax(0.9), | |
74 | fAreacut(0.0), | |
75 | fESD(0), | |
76 | fOutputList(0), | |
77 | fHistTrackPt(0), | |
78 | fHistCentrality(0), | |
79 | fHistJetEtaPhi(0), | |
80 | fHistTrackEtaPhi(0), | |
81 | fHistJetHEtaPhi(0) | |
82 | { | |
83 | // Constructor | |
84 | for(Int_t icent = 0; icent<6; ++icent){ | |
85 | fHistJetPt[icent]=0; | |
86 | fHistJetPtBias[icent]=0; | |
87 | for(Int_t iptjet = 0; iptjet<3; ++iptjet){ | |
88 | for(Int_t ieta = 0; ieta<3; ++ieta){ | |
89 | fHistJetH[icent][iptjet][ieta]=0; | |
90 | fHistJetHBias[icent][iptjet][ieta]=0; | |
91 | } | |
92 | } | |
93 | } | |
94 | ||
95 | DefineInput(0, TChain::Class()); | |
96 | DefineOutput(1, TList::Class()); | |
97 | ||
98 | } | |
99 | ||
100 | //________________________________________________________________________ | |
101 | void AliAnalysisTaskEmcalJetHMEC::UserCreateOutputObjects() | |
102 | { | |
103 | // Called once | |
104 | ||
105 | ||
106 | AliVEventHandler* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler(); | |
107 | if (!handler) { | |
108 | AliError("Input handler not available!"); | |
109 | return; | |
110 | } | |
111 | ||
112 | OpenFile(1); | |
113 | fOutputList = new TList(); | |
114 | fOutputList->SetOwner(); | |
115 | ||
116 | // Create histograms | |
117 | ||
118 | fHistTrackPt = new TH1F("fHistTrackPt", "P_{T} distribution", 1000, 0.0, 100.0); | |
119 | ||
120 | ||
121 | ||
122 | fHistCentrality = new TH1F("fHistCentrality","centrality",100,0,100); | |
123 | ||
124 | fHistJetEtaPhi = new TH2F("fHistJetEtaPhi","Jet eta-phi",900,-1.8,1.8,640,-3.2,3.2); | |
125 | fHistTrackEtaPhi = new TH2F("fHistTrackEtaPhi","Track eta-phi",900,-1.8,1.8,640,-3.2,3.2); | |
126 | fHistJetHEtaPhi = new TH2F("fHistJetHEtaPhi","Jet-Hadron deta-dphi",900,-1.8,1.8,640,-1.6,4.8); | |
127 | ||
128 | char name[200]; | |
129 | ||
130 | ||
131 | for(Int_t icent = 0; icent<6; ++icent){ | |
132 | sprintf(name,"fHistJetPt_%i",icent); | |
133 | fHistJetPt[icent] = new TH1F(name,name,200,0,200); | |
134 | fOutputList->Add(fHistJetPt[icent]); | |
135 | ||
136 | sprintf(name,"fHistJetPtBias_%i",icent); | |
137 | fHistJetPtBias[icent] = new TH1F(name,name,200,0,200); | |
138 | fOutputList->Add(fHistJetPtBias[icent]); | |
139 | ||
140 | for(Int_t iptjet = 0; iptjet<3; ++iptjet){ | |
141 | for(Int_t ieta = 0; ieta<3; ++ieta){ | |
142 | sprintf(name,"fHistJetH_%i_%i_%i",icent,iptjet,ieta); | |
143 | fHistJetH[icent][iptjet][ieta]=new TH2F(name,name,64,-0.5*TMath::Pi(),1.5*TMath::Pi(),300,0,30); | |
144 | fOutputList->Add(fHistJetH[icent][iptjet][ieta]); | |
145 | ||
146 | sprintf(name,"fHistJetHBias_%i_%i_%i",icent,iptjet,ieta); | |
147 | fHistJetHBias[icent][iptjet][ieta]=new TH2F(name,name,64,-0.5*TMath::Pi(),1.5*TMath::Pi(),300,0,30); | |
148 | fOutputList->Add(fHistJetHBias[icent][iptjet][ieta]); | |
149 | ||
150 | } | |
151 | } | |
152 | } | |
153 | ||
154 | ||
155 | ||
156 | fOutputList->Add(fHistTrackPt); | |
157 | fOutputList->Add(fHistCentrality); | |
158 | fOutputList->Add(fHistJetEtaPhi); | |
159 | fOutputList->Add(fHistTrackEtaPhi); | |
160 | fOutputList->Add(fHistJetHEtaPhi); | |
161 | ||
162 | ||
163 | PostData(1, fOutputList); | |
164 | ||
165 | } | |
166 | ||
167 | //________________________________________________________________________ | |
168 | ||
169 | Double_t AliAnalysisTaskEmcalJetHMEC:: RelativePhi(Double_t mphi,Double_t vphi) { | |
170 | ||
171 | ||
172 | if (vphi < -1*TMath::Pi()) vphi += (2*TMath::Pi()); | |
173 | else if (vphi > TMath::Pi()) vphi -= (2*TMath::Pi()); | |
174 | if (mphi < -1*TMath::Pi()) mphi += (2*TMath::Pi()); | |
175 | else if (mphi > TMath::Pi()) mphi -= (2*TMath::Pi()); | |
176 | Double_t dphi = vphi-mphi; | |
177 | if (dphi < -1*TMath::Pi()) dphi += (2*TMath::Pi()); | |
178 | else if (dphi > TMath::Pi()) dphi -= (2*TMath::Pi()); | |
179 | ||
180 | return dphi;//dphi in [-Pi, Pi] | |
181 | } | |
182 | ||
183 | ||
184 | //________________________________________________________________________ | |
185 | Int_t AliAnalysisTaskEmcalJetHMEC::GetCentBin(Double_t cent) const | |
186 | { | |
187 | // Get centrality bin. | |
188 | ||
189 | Int_t centbin = -1; | |
190 | if (cent>=0 && cent<10) | |
191 | centbin = 0; | |
192 | else if (cent>=10 && cent<20) | |
193 | centbin = 1; | |
194 | else if (cent>=20 && cent<30) | |
195 | centbin = 2; | |
196 | else if (cent>=30 && cent<40) | |
197 | centbin = 3; | |
198 | else if (cent>=40 && cent<50) | |
199 | centbin = 4; | |
200 | else if (cent>=50 && cent<90) | |
201 | centbin = 5; | |
202 | return centbin; | |
203 | } | |
204 | ||
205 | ||
206 | //________________________________________________________________________ | |
207 | Int_t AliAnalysisTaskEmcalJetHMEC::GetEtaBin(Double_t eta) const | |
208 | { | |
209 | // Get eta bin for histos. | |
210 | ||
211 | Int_t etabin = -1; | |
212 | if (TMath::Abs(eta)<=0.4) | |
213 | etabin = 0; | |
214 | else if (TMath::Abs(eta)>0.4 && TMath::Abs(eta)<0.8) | |
215 | etabin = 1; | |
216 | else if (TMath::Abs(eta)>=0.8) | |
217 | etabin = 2; | |
218 | return etabin; | |
219 | } | |
220 | //________________________________________________________________________ | |
221 | Int_t AliAnalysisTaskEmcalJetHMEC::GetpTjetBin(Double_t jetpt) const | |
222 | { | |
223 | // Get eta bin for histos. | |
224 | ||
225 | Int_t jetptbin = -1; | |
226 | if (jetpt>20) | |
227 | jetptbin = 0; | |
228 | else if (jetpt>30) | |
229 | jetptbin = 1; | |
230 | else if (jetpt>45) | |
231 | jetptbin = 2; | |
232 | ||
233 | return jetptbin; | |
234 | } | |
235 | ||
236 | ||
237 | //________________________________________________________________________ | |
238 | void AliAnalysisTaskEmcalJetHMEC::UserExec(Option_t *) | |
239 | { | |
240 | ||
241 | ||
242 | // Main loop called for each event | |
243 | // esd or aod mode | |
244 | Bool_t esdMode = kTRUE; | |
245 | if (dynamic_cast<AliAODEvent*>(InputEvent())) | |
246 | esdMode = kFALSE; | |
247 | ||
248 | ||
249 | if (esdMode) { | |
250 | // optimization in case autobranch loading is off | |
251 | AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager(); | |
252 | if (fTracksName == "Tracks") | |
253 | am->LoadBranch("Tracks"); | |
254 | } | |
255 | ||
256 | ||
257 | //get centrality | |
258 | TList *list = InputEvent()->GetList(); | |
259 | AliCentrality *centrality = InputEvent()->GetCentrality() ; | |
260 | Double_t fcent=-1; | |
261 | if(centrality) | |
262 | fcent = centrality->GetCentralityPercentile("V0M"); | |
263 | else | |
264 | fcent=99;//probably pp data | |
265 | ||
266 | if (fcent<0) { | |
267 | AliError(Form("Centrality negative: %f", fcent)); | |
268 | return; | |
269 | } | |
270 | ||
271 | ||
272 | fHistCentrality->Fill(fcent); | |
273 | Int_t centbin = GetCentBin(fcent); | |
274 | ||
275 | TClonesArray *jets = 0; | |
276 | TClonesArray *tracks = 0; | |
277 | ||
278 | tracks = dynamic_cast<TClonesArray*>(list->FindObject(fTracksName)); | |
279 | if (!tracks) { | |
280 | AliError(Form("Pointer to tracks %s == 0", fTracksName.Data() )); | |
281 | return; | |
282 | } | |
283 | const Int_t Ntracks=tracks->GetEntries(); | |
284 | ||
285 | jets= dynamic_cast<TClonesArray*>(list->FindObject(fJetsName)); | |
286 | const Int_t Njets = jets->GetEntries(); | |
287 | ||
288 | ||
289 | ||
290 | ||
291 | for (Int_t ijet = 0; ijet < Njets; ijet++) | |
292 | { | |
293 | AliEmcalJet *jet = static_cast<AliEmcalJet*>(jets->At(ijet)); | |
294 | ||
295 | if (!jet) | |
296 | continue; | |
297 | ||
298 | //pt,eta,phi,centrality | |
299 | float jetphi = jet->Phi(); | |
300 | if (jetphi>TMath::Pi()) | |
301 | jetphi = jetphi-2*TMath::Pi(); | |
302 | ||
303 | if ((jet->Phi()<fPhimin)||(jet->Phi()>fPhimax)) | |
304 | continue; | |
305 | if ((jet->Eta()<fEtamin)||(jet->Eta()>fEtamax)) | |
306 | continue; | |
307 | //fHistAreavsRawPt[centbin]->Fill(jet->Pt(),jet->Area()); | |
308 | if (jet->Area()<fAreacut) | |
309 | continue; | |
310 | //prevents 0 area jets from sneaking by when area cut == 0 | |
311 | if (jet->Area()==0) | |
312 | continue; | |
313 | ||
314 | ||
315 | Double_t jetPt = jet->Pt(); | |
316 | ||
317 | ||
318 | fHistJetPt[centbin]->Fill(jet->Pt()); | |
319 | ||
320 | fHistJetEtaPhi->Fill(jet->Eta(),jetphi); | |
321 | ||
322 | //fHistDeltaPtvsArea->Fill(jetPt,jet->Area()); | |
323 | ||
324 | if (jetPt<20) | |
325 | continue; | |
326 | ||
327 | for (Int_t iTracks = 0; iTracks < Ntracks; iTracks++) | |
328 | { | |
329 | AliVTrack* track = static_cast<AliVTrack*>(tracks->At(iTracks)); | |
330 | if (!track) { | |
331 | printf("ERROR: Could not receive track %d\n", iTracks); | |
332 | continue; | |
333 | } | |
334 | ||
335 | if(TMath::Abs(track->Eta())>0.9) continue; | |
336 | ||
337 | fHistTrackPt->Fill(track->Pt()); | |
338 | ||
339 | if (track->Pt()<0.5) | |
340 | continue; | |
341 | ||
342 | Double_t trackphi = track->Phi(); | |
343 | if (trackphi > TMath::Pi()) | |
344 | trackphi = trackphi-2*TMath::Pi(); | |
345 | ||
346 | Double_t tracketa=track->Eta(); | |
347 | Double_t jeteta=jet->Eta(); | |
348 | ||
349 | Double_t deta=tracketa-jeteta; | |
350 | Int_t ieta=GetEtaBin(deta); | |
351 | ||
352 | //Jet pt, track pt, dPhi,deta,fcent | |
353 | Double_t dphijh = RelativePhi(jetphi,trackphi); | |
354 | if (dphijh < -0.5*TMath::Pi()) | |
355 | dphijh= dphijh+ 2*TMath::Pi(); | |
356 | ||
357 | Int_t iptjet=-1; | |
358 | iptjet=GetpTjetBin(jetPt); | |
359 | ||
360 | fHistJetH[centbin][iptjet][ieta]->Fill(dphijh,track->Pt()); | |
361 | fHistJetHEtaPhi->Fill(deta,dphijh); | |
362 | fHistTrackEtaPhi->Fill(tracketa,trackphi); | |
363 | if ((jet->MaxTrackPt()>6) || (jet->MaxClusterPt()>6)){ | |
364 | fHistJetHBias[centbin][iptjet][ieta]->Fill(dphijh,track->Pt()); | |
365 | fHistJetPtBias[centbin]->Fill(jet->Pt()); | |
366 | ||
367 | } | |
368 | ||
369 | } //track loop | |
370 | ||
371 | }//jet loop | |
372 | ||
373 | ||
374 | PostData(1, fOutputList); | |
375 | } | |
376 | ||
377 | //________________________________________________________________________ | |
378 | void AliAnalysisTaskEmcalJetHMEC::Terminate(Option_t *) | |
379 | { | |
380 | //just terminate | |
381 | ||
382 | } | |
383 | ||
384 |