]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliJetTriggerSelectionTask.cxx
Add option to embed centrality
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliJetTriggerSelectionTask.cxx
1 // $Id$
2 //
3 // Jet trigger selection task.
4 //
5 // Author: S.Aiola
6
7 #include <TLorentzVector.h>
8 #include <TMath.h>
9 #include <TF1.h>
10
11 #include "AliVEvent.h"
12 #include "AliVCluster.h"
13 #include "AliEmcalJet.h"
14 #include "AliLog.h"
15 #include "AliVVZERO.h"
16 #include "AliESDUtils.h"
17
18 #include "AliJetTriggerSelectionTask.h"
19
20 ClassImp(AliJetTriggerSelectionTask)
21
22 //________________________________________________________________________
23 AliJetTriggerSelectionTask::AliJetTriggerSelectionTask() : 
24   AliAnalysisTaskEmcalJet("AliJetTriggerSelectionTask", kFALSE),
25   fEnergyThreshold(0),
26   fMaxDistance2(0.0225),
27   fTriggerBits(AliVEvent::kEMCEJE),
28   fTaskSettingsOk(kFALSE),
29   fNTriggers(0),
30   fVZERO(0),
31   fV0ATotMult(0),
32   fV0CTotMult(0)
33 {
34   // Default constructor.
35   
36   for (Int_t i = 0; i < 999; i++) {
37     fTrigPos[i][0] = -999;
38     fTrigPos[i][1] = -999;
39   }
40 }
41
42 //________________________________________________________________________
43 AliJetTriggerSelectionTask::AliJetTriggerSelectionTask(const char *name) : 
44   AliAnalysisTaskEmcalJet(name, kFALSE),
45   fEnergyThreshold(0),
46   fMaxDistance2(0.0225),
47   fTriggerBits(AliVEvent::kEMCEJE),
48   fTaskSettingsOk(kFALSE),
49   fNTriggers(0),
50   fVZERO(0),
51   fV0ATotMult(0),
52   fV0CTotMult(0)
53 {
54   // Standard constructor.
55
56   for (Int_t i = 0; i < 999; i++) {
57     fTrigPos[i][0] = -999;
58     fTrigPos[i][1] = -999;
59   }
60 }
61
62 //________________________________________________________________________
63 void AliJetTriggerSelectionTask::ExecOnce()
64 {
65   // Initialize the task.
66
67   AliAnalysisTaskEmcalJet::ExecOnce();
68
69   fTaskSettingsOk = kTRUE;
70
71   fVZERO = InputEvent()->GetVZEROData();
72   if (!fVZERO) {
73     AliError(Form("%s: AliVVZERO not available, task will not be executed!",GetName()));
74     fTaskSettingsOk = kFALSE;
75   }
76   
77   if (GetClusterArray() == 0) {
78     AliError(Form("%s: No cluster collection provided, task will not be executed!",GetName()));
79     fTaskSettingsOk = kFALSE;
80   }
81
82   if (GetJetArray() == 0) {
83     AliError(Form("%s: No jet collection provided, task will not be executed!",GetName()));
84     fTaskSettingsOk = kFALSE;
85   }
86  
87   if (!fEnergyThreshold) {
88     AliError(Form("%s: No threshold function provided, task will not be executed!",GetName()));
89     fTaskSettingsOk = kFALSE;
90   }
91 }
92
93 //________________________________________________________________________
94 Bool_t AliJetTriggerSelectionTask::RetrieveEventObjects()
95 {
96   // Retrieve event objects.
97
98   if (!AliAnalysisTaskEmcalJet::RetrieveEventObjects())
99     return kFALSE;
100
101   if (fVZERO) {
102     fV0ATotMult = fVZERO->GetMTotV0A();
103     fV0CTotMult = fVZERO->GetMTotV0C();
104   }
105
106   return kTRUE;
107 }
108
109 //________________________________________________________________________
110 Bool_t AliJetTriggerSelectionTask::Run()
111 {
112   // Run the analysis.
113
114   if (!fTaskSettingsOk) return kFALSE;
115
116   FindTriggers();
117   SelectJets();
118   
119   return kTRUE;
120 }
121
122 //________________________________________________________________________
123 void AliJetTriggerSelectionTask::FindTriggers()
124 {
125   fNTriggers = 0;
126
127   Int_t nclusters = GetNClusters();
128
129   Double_t th = fEnergyThreshold->Eval(fV0ATotMult+fV0CTotMult);
130
131   for (Int_t i = 0; i < nclusters; i++) {
132     if (fNTriggers >= 999) {
133       AliError("More than 999 triggers found!");
134       break;
135     }
136
137     AliVCluster *cluster = GetAcceptClusterFromArray(i);
138     if (!cluster)
139       continue;
140
141     //Printf("Cluster energy=%.3f, th=%.3f",cluster->E(), th);
142
143     if (cluster->E() > th) {
144       TLorentzVector vect;
145       cluster->GetMomentum(vect,fVertex);
146       fTrigPos[fNTriggers][0] = vect.Eta();
147       fTrigPos[fNTriggers][1] = vect.Phi();
148       fNTriggers++;
149     }
150   }
151
152   AliDebug(2,Form("%s: %d triggers found among %d candidates (cent=%.1f, mult=%.1f, th=%.2f)!",GetName(),fNTriggers,nclusters,fCent,fV0ATotMult+fV0CTotMult,th));
153 }
154
155 //________________________________________________________________________
156 void AliJetTriggerSelectionTask::SelectJets()
157 {
158   for (Int_t c = 0; c < fJetCollArray.GetEntriesFast(); c++) {
159
160     Int_t njets = GetNJets(c);
161
162     for (Int_t i = 0; i < njets; i++) {
163       AliEmcalJet *jet = GetAcceptJetFromArray(i,c);
164       if (IsTriggerJet(jet)) jet->AddTrigger(fTriggerBits);
165     }
166   }
167 }
168
169 //________________________________________________________________________
170 Bool_t AliJetTriggerSelectionTask::IsTriggerJet(AliEmcalJet *jet)
171 {
172   if (!jet) return kFALSE;
173
174   for (Int_t i = 0; i < fNTriggers; i++) {
175     Double_t deta = jet->Eta() - fTrigPos[i][0];
176     Double_t dphi = jet->Phi() - fTrigPos[i][1];
177     
178     Double_t d2 = deta * deta + dphi * dphi;
179
180     if (d2 < fMaxDistance2)
181       return kTRUE;
182   }
183   
184   return kFALSE;
185 }