]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJetDev.cxx
5f54aa50f856123e7cc114ed1792963b2adb8192
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskEmcalJetDev.cxx
1 // $Id$
2 //
3 // Emcal jet analysis base task.
4 //
5 // Author: S.Aiola
6
7 #include "AliAnalysisTaskEmcalJetDev.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 #include "AliJetContainer.h"
25
26 ClassImp(AliAnalysisTaskEmcalJetDev)
27
28 //________________________________________________________________________
29 AliAnalysisTaskEmcalJetDev::AliAnalysisTaskEmcalJetDev() : 
30   AliAnalysisTaskEmcalDev("AliAnalysisTaskEmcalJetDev"),
31   fJetsName(),
32   fRhoName(),
33   fJets(0),
34   fRho(0),
35   fRhoVal(0),
36   fJetCollArray()
37 {
38   // Default constructor.
39
40   fJetCollArray.SetOwner(kTRUE);
41 }
42
43 //________________________________________________________________________
44 AliAnalysisTaskEmcalJetDev::AliAnalysisTaskEmcalJetDev(const char *name, Bool_t histo) : 
45   AliAnalysisTaskEmcalDev(name, histo),
46   fJetsName(),
47   fRhoName(),
48   fJets(0),
49   fRho(0),
50   fRhoVal(0),
51   fJetCollArray()
52 {
53   // Standard constructor.
54
55   fJetCollArray.SetOwner(kTRUE);
56 }
57
58 //________________________________________________________________________
59 AliAnalysisTaskEmcalJetDev::~AliAnalysisTaskEmcalJetDev()
60 {
61   // Destructor
62 }
63
64 //________________________________________________________________________
65 Bool_t AliAnalysisTaskEmcalJetDev::AcceptBiasJet(AliEmcalJet *jet, Int_t c)
66
67   // Accept jet with a bias.
68
69   AliJetContainer *cont = GetJetContainer(c);
70   if(!cont) {
71     AliError(Form("%s:Container %d not found",GetName(),c));
72     return 0;
73   }
74
75   return cont->AcceptBiasJet(jet);
76
77 }
78
79 //________________________________________________________________________
80 Float_t* AliAnalysisTaskEmcalJetDev::GenerateFixedBinArray(Int_t n, Float_t min, Float_t max) const
81 {
82   Float_t *bins = new Float_t[n+1];
83
84   Float_t binWidth = (max-min)/n;
85   bins[0] = min;
86   for (Int_t i = 1; i <= n; i++) {
87     bins[i] = bins[i-1]+binWidth;
88   }
89
90   return bins;
91 }
92
93 //________________________________________________________________________
94 Double_t AliAnalysisTaskEmcalJetDev::GetLeadingHadronPt(AliEmcalJet *jet, Int_t c)
95 {
96
97   AliJetContainer *cont = GetJetContainer(c);
98   if(!cont) {
99     AliError(Form("%s:Container %d not found",GetName(),c));
100     return 0;
101   }
102
103   return cont->GetLeadingHadronPt(jet);
104
105 }
106
107 //________________________________________________________________________
108 Bool_t AliAnalysisTaskEmcalJetDev::AcceptJet(AliEmcalJet *jet, Int_t c) 
109 {   
110   // Return true if jet is accepted.
111   if (!jet)
112     return kFALSE;
113
114   AliJetContainer *cont = GetJetContainer(c);
115   if(!cont) {
116     AliError(Form("%s:Container %d not found",GetName(),c));
117     return 0;
118   }
119
120   return cont->AcceptJet(jet);
121
122 }
123
124 //________________________________________________________________________
125 AliRhoParameter *AliAnalysisTaskEmcalJetDev::GetRhoFromEvent(const char *name)
126 {
127   // Get rho from event.
128
129   AliRhoParameter *rho = 0;
130   TString sname(name);
131   if (!sname.IsNull()) {
132     rho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(sname));
133     if (!rho) {
134       AliWarning(Form("%s: Could not retrieve rho with name %s!", GetName(), name)); 
135       return 0;
136     }
137   }
138   return rho;
139 }
140
141 //________________________________________________________________________
142 void AliAnalysisTaskEmcalJetDev::ExecOnce()
143 {
144   // Init the analysis.
145
146   AliAnalysisTaskEmcalDev::ExecOnce();
147
148   //Load all requested track branches - each container knows name already
149   if(fJetCollArray.GetEntriesFast()==0) {
150     AliWarning("There are no jet collections");
151     return;
152   }
153
154   for(Int_t i =0; i<fJetCollArray.GetEntriesFast(); i++) {
155     AliJetContainer *cont = static_cast<AliJetContainer*>(fJetCollArray.At(i));
156     cont->SetJetArray(InputEvent());
157     cont->SetEMCALGeometry();
158     cont->LoadRho(InputEvent());
159   }
160
161   //Get Jets, cuts and rho for first jet container
162   AliJetContainer *cont = GetJetContainer(0);
163   // fJetsName = cont->GetArrayName();
164   if (fAnaType == kTPC) {
165     cont->SetJetAcceptanceType(AliJetContainer::kTPC);
166     cont->SetJetEtaPhiTPC();
167   }
168   else if (fAnaType == kEMCAL) {
169     cont->SetJetAcceptanceType(AliJetContainer::kEMCAL);
170     cont->SetJetEtaPhiEMCAL();
171   }
172   
173   fJets = GetJetArray(0);
174   if(!fJets && fJetCollArray.GetEntriesFast()>0) {
175     AliError(Form("%s: Could not retrieve first jet branch!", GetName()));
176     fInitialized = kFALSE;
177     return;
178   }
179
180   fRhoName = cont->GetRhoName();
181   if(!fRhoName.IsNull()) {
182     fRho = cont->GetRhoParameter();
183     if(!fRho) {
184       AliError(Form("%s: Could not retrieve rho of first jet branch!", GetName()));
185       fInitialized = kFALSE;
186       return;
187     }
188   }
189
190 }
191
192 //________________________________________________________________________
193 Bool_t AliAnalysisTaskEmcalJetDev::GetSortedArray(Int_t indexes[], TClonesArray *array, Double_t rho, Int_t c)
194 {
195   // Get the leading jets.
196
197   static Float_t pt[9999] = {0};
198
199   if (!array)
200     return 0;
201
202   const Int_t n = array->GetEntriesFast();
203
204   if (n < 1)
205     return kFALSE;
206   
207   if (array->GetClass()->GetBaseClass("AliEmcalJet")) {
208
209     for (Int_t i = 0; i < n; i++) {
210
211       pt[i] = -FLT_MAX;
212
213       AliEmcalJet* jet = static_cast<AliEmcalJet*>(array->At(i));
214       
215       if (!jet) {
216         AliError(Form("Could not receive jet %d", i));
217         continue;
218       }
219       
220       if (!AcceptJet(jet,c))
221         continue;
222       
223       pt[i] = jet->Pt() - rho * jet->Area();
224     }
225   }
226
227   else if (array->GetClass()->GetBaseClass("AliVTrack")) {
228
229     for (Int_t i = 0; i < n; i++) {
230
231       pt[i] = -FLT_MAX;
232
233       AliVTrack* track = static_cast<AliVTrack*>(array->At(i));
234       
235       if (!track) {
236         AliError(Form("Could not receive track %d", i));
237         continue;
238       }  
239       
240       if (!AcceptTrack(track, c))
241         continue;
242       
243       pt[i] = track->Pt();
244     }
245   }
246
247   else if (array->GetClass()->GetBaseClass("AliVCluster")) {
248
249     for (Int_t i = 0; i < n; i++) {
250
251       pt[i] = -FLT_MAX;
252
253       AliVCluster* cluster = static_cast<AliVCluster*>(array->At(i));
254       
255       if (!cluster) {
256         AliError(Form("Could not receive cluster %d", i));
257         continue;
258       }  
259       
260       if (!AcceptCluster(cluster, c))
261         continue;
262
263       TLorentzVector nPart;
264       cluster->GetMomentum(nPart, const_cast<Double_t*>(fVertex));
265       
266       pt[i] = nPart.Pt();
267     }
268   }
269
270   TMath::Sort(n, pt, indexes);
271
272   if (pt[indexes[0]] == -FLT_MAX) 
273     return 0;
274
275   return kTRUE;
276 }
277
278 //________________________________________________________________________
279 Bool_t AliAnalysisTaskEmcalJetDev::IsJetCluster(AliEmcalJet* jet, Int_t iclus, Bool_t sorted) const
280 {
281   // Return true if cluster is in jet.
282
283   for (Int_t i = 0; i < jet->GetNumberOfClusters(); ++i) {
284     Int_t ijetclus = jet->ClusterAt(i);
285     if (sorted && ijetclus > iclus)
286       return kFALSE;
287     if (ijetclus == iclus)
288       return kTRUE;
289   }
290   return kFALSE;
291 }
292
293 //________________________________________________________________________
294 Bool_t AliAnalysisTaskEmcalJetDev::IsJetTrack(AliEmcalJet* jet, Int_t itrack, Bool_t sorted) const
295 {
296   // Return true if track is in jet.
297
298   for (Int_t i = 0; i < jet->GetNumberOfTracks(); ++i) {
299     Int_t ijettrack = jet->TrackAt(i);
300     if (sorted && ijettrack > itrack)
301       return kFALSE;
302     if (ijettrack == itrack)
303       return kTRUE;
304   }
305   return kFALSE;
306 }
307
308 //________________________________________________________________________
309 Bool_t AliAnalysisTaskEmcalJetDev::RetrieveEventObjects()
310 {
311   // Retrieve objects from event.
312
313   if (!AliAnalysisTaskEmcalDev::RetrieveEventObjects())
314     return kFALSE;
315
316   if (fRho)
317     fRhoVal = fRho->GetVal();
318
319   return kTRUE;
320 }
321
322 //________________________________________________________________________
323 void AliAnalysisTaskEmcalJetDev::AddJetContainer(const char *n, TString defaultCutType) {
324
325   // Add particle container
326   // will be called in AddTask macro
327
328   AliJetContainer *cont = 0x0;
329   cont = new AliJetContainer();
330   cont->SetArrayName(n);
331
332   if(!defaultCutType.IsNull()) {
333     if(defaultCutType.EqualTo("TPC"))
334      cont->SetJetAcceptanceType(AliJetContainer::kTPC);
335     else if(defaultCutType.EqualTo("EMCAL"))
336       cont->SetJetAcceptanceType(AliJetContainer::kEMCAL);
337     else
338       AliWarning(Form("%s: default cut type %s not recognized. Not setting cuts.",GetName(),defaultCutType.Data()));
339   } else
340     cont->SetJetAcceptanceType(AliJetContainer::kUser);
341  
342   fJetCollArray.Add(cont);
343
344 }
345
346 //________________________________________________________________________
347 AliJetContainer* AliAnalysisTaskEmcalJetDev::GetJetContainer(Int_t i) const{
348   // Get i^th jet container
349
350   if(i<0 || i>fJetCollArray.GetEntriesFast()) return 0;
351   AliJetContainer *cont = static_cast<AliJetContainer*>(fJetCollArray.At(i));
352   return cont;
353 }
354
355 //________________________________________________________________________
356 void AliAnalysisTaskEmcalJetDev::SetRhoName(const char *n, Int_t c)
357 {
358   AliJetContainer *cont = GetJetContainer(c);
359   cont->SetRhoName(n);
360 }
361
362 //________________________________________________________________________
363 void AliAnalysisTaskEmcalJetDev::SetJetEtaLimits(Float_t min, Float_t max, Int_t c)
364 {
365   AliJetContainer *cont = GetJetContainer(c);
366   cont->SetJetEtaLimits(min,max);
367 }
368
369 //________________________________________________________________________
370 void AliAnalysisTaskEmcalJetDev::SetJetPhiLimits(Float_t min, Float_t max, Int_t c)
371 {
372   AliJetContainer *cont = GetJetContainer(c);
373   cont->SetJetPhiLimits(min,max);
374 }
375
376 //________________________________________________________________________
377 void AliAnalysisTaskEmcalJetDev::SetJetAreaCut(Float_t cut, Int_t c)
378 {
379   AliJetContainer *cont = GetJetContainer(c);
380   cont->SetJetAreaCut(cut);
381 }
382
383 //________________________________________________________________________
384 void AliAnalysisTaskEmcalJetDev::SetPercAreaCut(Float_t p, Int_t c)
385 {
386   AliJetContainer *cont = GetJetContainer(c);
387   cont->SetPercAreaCut(p);
388 }
389
390 //________________________________________________________________________
391 void AliAnalysisTaskEmcalJetDev::SetAreaEmcCut(Double_t a, Int_t c)
392 {
393   AliJetContainer *cont = GetJetContainer(c);
394   cont->SetAreaEmcCut(a);
395 }
396
397 //________________________________________________________________________
398 void AliAnalysisTaskEmcalJetDev::SetJetPtCut(Float_t cut, Int_t c)
399 {
400   AliJetContainer *cont = GetJetContainer(c);
401   cont->SetJetPtCut(cut);
402 }
403
404 //________________________________________________________________________
405 void AliAnalysisTaskEmcalJetDev::SetJetRadius(Float_t r, Int_t c)
406 {
407   AliJetContainer *cont = GetJetContainer(c);
408   cont->SetJetRadius(r);
409 }
410
411 //________________________________________________________________________
412 void AliAnalysisTaskEmcalJetDev::SetMaxClusterPt(Float_t cut, Int_t c)
413 {
414   AliJetContainer *cont = GetJetContainer(c);
415   cont->SetMaxClusterPt(cut);
416 }
417
418 //________________________________________________________________________
419 void AliAnalysisTaskEmcalJetDev::SetMaxTrackPt(Float_t cut, Int_t c)
420 {
421   AliJetContainer *cont = GetJetContainer(c);
422   cont->SetMaxTrackPt(cut);
423 }
424
425 //________________________________________________________________________
426 void AliAnalysisTaskEmcalJetDev::SetPtBiasJetClus(Float_t cut, Int_t c)
427 {
428   AliJetContainer *cont = GetJetContainer(c);
429   cont->SetPtBiasJetClus(cut);
430 }
431
432 //________________________________________________________________________
433 void AliAnalysisTaskEmcalJetDev::SetPtBiasJetTrack(Float_t cut, Int_t c)
434 {
435   AliJetContainer *cont = GetJetContainer(c);
436   cont->SetPtBiasJetTrack(cut);
437 }
438
439 //________________________________________________________________________
440 void AliAnalysisTaskEmcalJetDev::SetLeadingHadronType(Int_t t, Int_t c)
441 {
442   AliJetContainer *cont = GetJetContainer(c);
443   cont->SetLeadingHadronType(t);
444 }
445
446 //________________________________________________________________________
447 void AliAnalysisTaskEmcalJetDev::SetNLeadingJets(Int_t t, Int_t c)
448 {
449   AliJetContainer *cont = GetJetContainer(c);
450   cont->SetNLeadingJets(t);
451 }
452
453 //________________________________________________________________________
454 void AliAnalysisTaskEmcalJetDev::SetJetBitMap(UInt_t m, Int_t c)
455 {
456   AliJetContainer *cont = GetJetContainer(c);
457   cont->SetJetBitMap(m);
458 }
459
460 //________________________________________________________________________
461 TClonesArray* AliAnalysisTaskEmcalJetDev::GetJetArray(Int_t i) const {
462   // Get i^th TClonesArray with AliEmcalJet
463
464   AliJetContainer *cont = GetJetContainer(i);
465   if(!cont) {
466     AliError(Form("%s:Container %d not found",GetName(),i));
467     return 0;
468   }
469   return cont->GetArray();
470
471 }
472
473 //________________________________________________________________________
474 AliEmcalJet* AliAnalysisTaskEmcalJetDev::GetAcceptJetFromArray(Int_t j, Int_t c) const {
475   // Get jet j if accepted from  container c
476   // If jet not accepted return 0
477
478   AliJetContainer *cont = GetJetContainer(c);
479   if(!cont) {
480     AliError(Form("%s:Container %d not found",GetName(),c));
481     return 0;
482   }
483   AliEmcalJet *jet = cont->GetAcceptJet(j);
484
485   return jet;
486   
487 }
488
489 //________________________________________________________________________
490 Int_t AliAnalysisTaskEmcalJetDev::GetNJets(Int_t i) const {
491   // Get number of entries in jet array i
492
493   AliJetContainer *cont = GetJetContainer(i);
494   if(!cont) {
495     AliError(Form("%s:Container %d not found",GetName(),i));
496     return 0;
497   }
498   return cont->GetNJets();
499
500 }
501
502 //________________________________________________________________________
503 Double_t AliAnalysisTaskEmcalJetDev::GetRhoVal(Int_t i) const {
504   // Get rho corresponding to jet array i
505
506   AliJetContainer *cont = GetJetContainer(i);
507   if(!cont) {
508     AliError(Form("%s:Container %d not found",GetName(),i));
509     return 0;
510   }
511   return cont->GetRhoVal();
512
513 }
514
515
516