]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJet.cxx
7c4d50389149c0fbd02752cc6903ace181d13df8
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskEmcalJet.cxx
1 // $Id$
2 //
3 // Emcal jet analysis base task.
4 //
5 // Author: S.Aiola
6
7 #include "AliAnalysisTaskEmcalJet.h"
8
9 #include <TChain.h>
10 #include <TClonesArray.h>
11 #include <TList.h>
12 #include <TObject.h>
13
14 #include "AliAnalysisManager.h"
15 #include "AliCentrality.h"
16 #include "AliEMCALGeometry.h"
17 #include "AliESDEvent.h"
18 #include "AliEmcalJet.h"
19 #include "AliLog.h"
20 #include "AliRhoParameter.h"
21 #include "AliVCluster.h"
22 #include "AliVEventHandler.h"
23 #include "AliVParticle.h"
24
25 ClassImp(AliAnalysisTaskEmcalJet)
26
27 //________________________________________________________________________
28 AliAnalysisTaskEmcalJet::AliAnalysisTaskEmcalJet() : 
29   AliAnalysisTaskEmcal("AliAnalysisTaskEmcalJet"),
30   fJetRadius(0.4),
31   fJetsName(),
32   fRhoName(),
33   fPtBiasJetTrack(0),
34   fPtBiasJetClus(0),
35   fJetPtCut(1),
36   fJetAreaCut(-1),
37   fPercAreaCut(-1),
38   fAreaEmcCut(0),
39   fJetMinEta(-0.9),
40   fJetMaxEta(0.9),
41   fJetMinPhi(-10),
42   fJetMaxPhi(10),
43   fMaxClusterPt(1000),
44   fMaxTrackPt(100),
45   fLeadingHadronType(0),
46   fNLeadingJets(1),
47   fJetBitMap(0),
48   fJets(0),
49   fRho(0),
50   fRhoVal(0)
51 {
52   // Default constructor.
53 }
54
55 //________________________________________________________________________
56 AliAnalysisTaskEmcalJet::AliAnalysisTaskEmcalJet(const char *name, Bool_t histo) : 
57   AliAnalysisTaskEmcal(name, histo),
58   fJetRadius(0.4),
59   fJetsName(),
60   fRhoName(),
61   fPtBiasJetTrack(0),
62   fPtBiasJetClus(0),
63   fJetPtCut(1),
64   fJetAreaCut(-1),
65   fPercAreaCut(-1),
66   fAreaEmcCut(0),
67   fJetMinEta(-0.9),
68   fJetMaxEta(0.9),
69   fJetMinPhi(-10),
70   fJetMaxPhi(10),
71   fMaxClusterPt(1000),
72   fMaxTrackPt(100),
73   fLeadingHadronType(0),
74   fNLeadingJets(1),
75   fJetBitMap(0),
76   fJets(0),
77   fRho(0),
78   fRhoVal(0)
79 {
80   // Standard constructor.
81 }
82
83 //________________________________________________________________________
84 AliAnalysisTaskEmcalJet::~AliAnalysisTaskEmcalJet()
85 {
86   // Destructor
87 }
88
89 //________________________________________________________________________
90 Bool_t AliAnalysisTaskEmcalJet::AcceptBiasJet(AliEmcalJet *jet) const
91
92   // Accept jet with a bias.
93
94   if (fLeadingHadronType == 0) {
95     if (jet->MaxTrackPt() < fPtBiasJetTrack) return kFALSE;
96   }
97   else if (fLeadingHadronType == 1) {
98     if (jet->MaxClusterPt() < fPtBiasJetClus) return kFALSE;
99   }
100   else {
101     if (jet->MaxTrackPt() < fPtBiasJetTrack && jet->MaxClusterPt() < fPtBiasJetClus) return kFALSE;
102   }
103
104   return kTRUE;
105 }
106
107 //________________________________________________________________________
108 Float_t* AliAnalysisTaskEmcalJet::GenerateFixedBinArray(Int_t n, Float_t min, Float_t max) const
109 {
110   Float_t *bins = new Float_t[n+1];
111
112   Float_t binWidth = (max-min)/n;
113   bins[0] = min;
114   for (Int_t i = 1; i <= n; i++) {
115     bins[i] = bins[i-1]+binWidth;
116   }
117
118   return bins;
119 }
120
121 //________________________________________________________________________
122 Double_t AliAnalysisTaskEmcalJet::GetLeadingHadronPt(AliEmcalJet *jet) const
123 {
124   if (fLeadingHadronType == 0)       // charged leading hadron
125     return jet->MaxTrackPt();
126   else if (fLeadingHadronType == 1)  // neutral leading hadron
127     return jet->MaxClusterPt();
128   else                               // charged or neutral
129     return jet->MaxPartPt();
130 }
131
132 //________________________________________________________________________
133 Bool_t AliAnalysisTaskEmcalJet::AcceptJet(AliEmcalJet *jet) const
134 {   
135   // Return true if jet is accepted.
136   if (!jet)
137     return kFALSE;
138
139   if (jet->TestBits(fJetBitMap) != (Int_t)fJetBitMap)
140     return kFALSE;
141
142   if (jet->Pt() < fJetPtCut) 
143     return kFALSE;
144  
145   if (jet->Area() < fJetAreaCut) 
146     return kFALSE;
147
148   if (jet->AreaEmc()<fAreaEmcCut)
149     return kFALSE;
150
151   if (!AcceptBiasJet(jet))
152     return kFALSE;
153   
154   if (jet->MaxTrackPt() > fMaxTrackPt || jet->MaxClusterPt() > fMaxClusterPt)
155     return kFALSE;
156
157   Double_t jetPhi = jet->Phi();
158   Double_t jetEta = jet->Eta();
159   
160   if (fJetMinPhi < 0) // if limits are given in (-pi, pi) range
161     jetPhi -= TMath::Pi() * 2;
162
163   return (Bool_t)(jetEta > fJetMinEta && jetEta < fJetMaxEta && jetPhi > fJetMinPhi && jetPhi < fJetMaxPhi);
164 }
165
166 //________________________________________________________________________
167 AliRhoParameter *AliAnalysisTaskEmcalJet::GetRhoFromEvent(const char *name)
168 {
169   // Get rho from event.
170
171   AliRhoParameter *rho = 0;
172   TString sname(name);
173   if (!sname.IsNull()) {
174     rho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(sname));
175     if (!rho) {
176       AliWarning(Form("%s: Could not retrieve rho with name %s!", GetName(), name)); 
177       return 0;
178     }
179   }
180   return rho;
181 }
182
183 //________________________________________________________________________
184 void AliAnalysisTaskEmcalJet::ExecOnce()
185 {
186   // Init the analysis.
187
188   AliAnalysisTaskEmcal::ExecOnce();
189
190   if (fPercAreaCut >= 0) {
191     if (fJetAreaCut >= 0)
192       AliInfo(Form("%s: jet area cut will be calculated as a percentage of the average area, given value will be overwritten", GetName()));
193     fJetAreaCut = fPercAreaCut * fJetRadius * fJetRadius * TMath::Pi();
194   }
195   if (fJetAreaCut < 0)
196     fJetAreaCut = 0;
197
198   if (fAnaType == kTPC) {
199     SetJetEtaLimits(-0.5, 0.5);
200     SetJetPhiLimits(-10, 10);
201   } 
202   else if (fAnaType == kEMCAL && fGeom) {
203     SetJetEtaLimits(fGeom->GetArm1EtaMin() + fJetRadius, fGeom->GetArm1EtaMax() - fJetRadius);
204     SetJetPhiLimits(fGeom->GetArm1PhiMin() * TMath::DegToRad() + fJetRadius, fGeom->GetArm1PhiMax() * TMath::DegToRad() - fJetRadius);
205   }
206
207   if (!fRhoName.IsNull() && !fRho) {
208     fRho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(fRhoName));
209     if (!fRho) {
210       AliError(Form("%s: Could not retrieve rho %s!", GetName(), fRhoName.Data()));
211       fInitialized = kFALSE;
212       return;
213     }
214   }
215
216   if (!fJetsName.IsNull() && !fJets) {
217     fJets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetsName));
218     if (!fJets) {
219       AliError(Form("%s: Could not retrieve jets %s!", GetName(), fJetsName.Data()));
220       fInitialized = kFALSE;
221       return;
222     }
223     else if (!fJets->GetClass()->GetBaseClass("AliEmcalJet")) {
224       AliError(Form("%s: Collection %s does not contain AliEmcalJet objects!", GetName(), fJetsName.Data())); 
225       fJets = 0;
226       fInitialized = kFALSE;
227       return;
228     }
229   }
230 }
231
232 //________________________________________________________________________
233 Bool_t AliAnalysisTaskEmcalJet::GetSortedArray(Int_t indexes[], TClonesArray *array, Double_t rho) const
234 {
235   // Get the leading jets.
236
237   static Float_t pt[9999] = {0};
238
239   if (!array)
240     return 0;
241
242   const Int_t n = array->GetEntriesFast();
243
244   if (n < 1)
245     return kFALSE;
246
247   if (array->GetClass()->GetBaseClass("AliEmcalJet")) {
248
249     for (Int_t i = 0; i < n; i++) {
250
251       pt[i] = -FLT_MAX;
252
253       AliEmcalJet* jet = static_cast<AliEmcalJet*>(array->At(i));
254       
255       if (!jet) {
256         AliError(Form("Could not receive jet %d", i));
257         continue;
258       }
259       
260       if (!AcceptJet(jet))
261         continue;
262       
263       pt[i] = jet->Pt() - rho * jet->Area();
264     }
265   }
266
267   else if (array->GetClass()->GetBaseClass("AliVTrack")) {
268
269     for (Int_t i = 0; i < n; i++) {
270
271       pt[i] = -FLT_MAX;
272
273       AliVTrack* track = static_cast<AliVTrack*>(array->At(i));
274       
275       if (!track) {
276         AliError(Form("Could not receive track %d", i));
277         continue;
278       }  
279       
280       if (!AcceptTrack(track))
281         continue;
282       
283       pt[i] = track->Pt();
284     }
285   }
286
287   else if (array->GetClass()->GetBaseClass("AliVCluster")) {
288
289     for (Int_t i = 0; i < n; i++) {
290
291       pt[i] = -FLT_MAX;
292
293       AliVCluster* cluster = static_cast<AliVCluster*>(array->At(i));
294       
295       if (!cluster) {
296         AliError(Form("Could not receive cluster %d", i));
297         continue;
298       }  
299       
300       if (!AcceptCluster(cluster))
301         continue;
302
303       TLorentzVector nPart;
304       cluster->GetMomentum(nPart, const_cast<Double_t*>(fVertex));
305       
306       pt[i] = nPart.Pt();
307     }
308   }
309
310   TMath::Sort(n, pt, indexes);
311
312   if (pt[indexes[0]] == -FLT_MAX) 
313     return 0;
314
315   return kTRUE;
316 }
317
318 //________________________________________________________________________
319 Bool_t AliAnalysisTaskEmcalJet::IsJetCluster(AliEmcalJet* jet, Int_t iclus, Bool_t sorted) const
320 {
321   // Return true if cluster is in jet.
322
323   for (Int_t i = 0; i < jet->GetNumberOfClusters(); ++i) {
324     Int_t ijetclus = jet->ClusterAt(i);
325     if (sorted && ijetclus > iclus)
326       return kFALSE;
327     if (ijetclus == iclus)
328       return kTRUE;
329   }
330   return kFALSE;
331 }
332
333 //________________________________________________________________________
334 Bool_t AliAnalysisTaskEmcalJet::IsJetTrack(AliEmcalJet* jet, Int_t itrack, Bool_t sorted) const
335 {
336   // Return true if track is in jet.
337
338   for (Int_t i = 0; i < jet->GetNumberOfTracks(); ++i) {
339     Int_t ijettrack = jet->TrackAt(i);
340     if (sorted && ijettrack > itrack)
341       return kFALSE;
342     if (ijettrack == itrack)
343       return kTRUE;
344   }
345   return kFALSE;
346 }
347
348 //________________________________________________________________________
349 Bool_t AliAnalysisTaskEmcalJet::RetrieveEventObjects()
350 {
351   // Retrieve objects from event.
352
353   if (!AliAnalysisTaskEmcal::RetrieveEventObjects())
354     return kFALSE;
355
356   if (fRho)
357     fRhoVal = fRho->GetVal();
358
359   return kTRUE;
360 }