]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJet.cxx
protections against failures in deleting event content
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskEmcalJet.cxx
1 // $Id: AliAnalysisTaskEmcalJet.cxx 56756 2012-05-30 05:03:02Z loizides $
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   fPtBiasJetTrack(5),
33   fPtBiasJetClus(5),
34   fJetPtCut(1),
35   fJetAreaCut(-1),
36   fPercAreaCut(0.8),
37   fMinEta(-0.9),
38   fMaxEta(0.9),
39   fMinPhi(-10),
40   fMaxPhi(10),
41   fMaxClusterPt(100),
42   fMaxTrackPt(100),
43   fJets(0)
44 {
45   // Default constructor.
46 }
47
48 //________________________________________________________________________
49 AliAnalysisTaskEmcalJet::AliAnalysisTaskEmcalJet(const char *name, Bool_t histo) : 
50   AliAnalysisTaskEmcal(name, histo),
51   fJetRadius(0.4),
52   fJetsName(),
53   fPtBiasJetTrack(5),
54   fPtBiasJetClus(5),
55   fJetPtCut(1),
56   fJetAreaCut(-1),
57   fPercAreaCut(0.8),
58   fMinEta(-0.9),
59   fMaxEta(0.9),
60   fMinPhi(-10),
61   fMaxPhi(10),
62   fMaxClusterPt(100),
63   fMaxTrackPt(100),
64   fJets(0)
65 {
66   // Standard constructor.
67 }
68
69 //________________________________________________________________________
70 AliAnalysisTaskEmcalJet::~AliAnalysisTaskEmcalJet()
71 {
72   // Destructor
73 }
74
75 //________________________________________________________________________
76 Bool_t AliAnalysisTaskEmcalJet::AcceptBiasJet(AliEmcalJet *jet) const
77
78   // Accept jet with a bias.
79
80   if (jet->MaxTrackPt() < fPtBiasJetTrack && (fAnaType == kTPC || jet->MaxClusterPt() < fPtBiasJetClus))
81     return kFALSE;
82   else
83     return kTRUE;
84 }
85
86 //________________________________________________________________________
87 Bool_t AliAnalysisTaskEmcalJet::AcceptJet(AliEmcalJet *jet, Bool_t bias, Bool_t upCut) const
88 {   
89   // Return true if jet is accepted.
90
91   if (jet->Pt() <= fJetPtCut)
92     return kFALSE;
93   if (jet->Area() <= fJetAreaCut)
94     return kFALSE;
95   if (bias && !AcceptBiasJet(jet))
96     return kFALSE;
97   if (upCut && (jet->MaxTrackPt() > fMaxTrackPt || jet->MaxClusterPt() > fMaxClusterPt))
98     return kFALSE;
99
100   return (Bool_t)(jet->Eta() > fMinEta && jet->Eta() < fMaxEta && jet->Phi() > fMinPhi && jet->Phi() < fMaxPhi);
101 }
102
103 //________________________________________________________________________
104 AliRhoParameter *AliAnalysisTaskEmcalJet::GetRhoFromEvent(const char *name)
105 {
106   // Get rho from event.
107
108   AliRhoParameter *rho = 0;
109   TString sname(name);
110   if (!sname.IsNull()) {
111     rho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(sname));
112     if (!rho) {
113       AliWarning(Form("%s: Could not retrieve rho with name %s!", GetName(), name)); 
114       return 0;
115     }
116   }
117   return rho;
118 }
119
120 //________________________________________________________________________
121 void AliAnalysisTaskEmcalJet::ExecOnce()
122 {
123   // Init the analysis.
124
125   if (fPercAreaCut >= 0) {
126     if (fJetAreaCut >= 0)
127       AliInfo(Form("%s: jet area cut will be calculated as a percentage of the average area, given value will be overwritten", GetName()));
128     fJetAreaCut = fPercAreaCut * fJetRadius * fJetRadius * TMath::Pi();
129   }
130
131   if (fAnaType == kTPC) {
132     SetEtaLimits(-0.9 + fJetRadius, 0.9 - fJetRadius);
133     SetPhiLimits(-10, 10);
134   } else if (fAnaType == kEMCAL || fAnaType == kTPCSmall || fAnaType == kEMCALOnly) {
135     AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
136     if (geom) {
137       SetEtaLimits(geom->GetArm1EtaMin() + fJetRadius, geom->GetArm1EtaMax() - fJetRadius);
138       SetPhiLimits(geom->GetArm1PhiMin() * TMath::DegToRad() + fJetRadius, geom->GetArm1PhiMax() * TMath::DegToRad() - fJetRadius);
139     }
140     else {
141       AliWarning(Form("%s: Can not create geometry", GetName()));
142     }
143   } else {
144     AliWarning(Form("%s: Analysis type not recognized! Assuming kTPC!", GetName()));
145     SetAnaType(kTPC);
146     ExecOnce();
147     return;
148   }
149
150   if (fAnaType == kTPCSmall)
151     fAnaType = kTPC;
152
153   AliAnalysisTaskEmcal::ExecOnce();
154 }
155
156 //________________________________________________________________________
157 Bool_t AliAnalysisTaskEmcalJet::IsJetCluster(AliEmcalJet* jet, Int_t iclus, Bool_t sorted) const
158 {
159   // Return true if cluster is in jet.
160
161   for (Int_t i = 0; i < jet->GetNumberOfClusters(); ++i) {
162     Int_t ijetclus = jet->ClusterAt(i);
163     if (sorted && ijetclus > iclus)
164       return kFALSE;
165     if (ijetclus == iclus)
166       return kTRUE;
167   }
168   return kFALSE;
169 }
170
171 //________________________________________________________________________
172 Bool_t AliAnalysisTaskEmcalJet::IsJetTrack(AliEmcalJet* jet, Int_t itrack, Bool_t sorted) const
173 {
174   // Return true if track is in jet.
175
176   for (Int_t i = 0; i < jet->GetNumberOfTracks(); ++i) {
177     Int_t ijettrack = jet->TrackAt(i);
178     if (sorted && ijettrack > itrack)
179       return kFALSE;
180     if (ijettrack == itrack)
181       return kTRUE;
182   }
183   return kFALSE;
184 }
185
186 //________________________________________________________________________
187 Bool_t AliAnalysisTaskEmcalJet::RetrieveEventObjects()
188 {
189   // Retrieve objects from event.
190
191   if (!AliAnalysisTaskEmcal::RetrieveEventObjects())
192     return kFALSE;
193
194   if (!fJetsName.IsNull() && !fJets) {
195     fJets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetsName));
196     if (!fJets) {
197       AliError(Form("%s: Could not retrieve jets %s!", GetName(), fJetsName.Data()));
198       return kFALSE;
199     }
200     else if (!fJets->GetClass()->GetBaseClass("AliEmcalJet")) {
201       AliError(Form("%s: Collection %s does not contain AliEmcalJet objects!", GetName(), fJetsName.Data())); 
202       fJets = 0;
203       return kFALSE;
204     }
205   }
206
207   return kTRUE;
208 }