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