]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/CaloCalib/AliAnalysisTaskEMCALPi0PbPb.cxx
first histograms filled
[u/mrichter/AliRoot.git] / PWG4 / CaloCalib / AliAnalysisTaskEMCALPi0PbPb.cxx
CommitLineData
ea3fd2d5 1// $Id$
2
ea3fd2d5 3#include "AliAnalysisTaskEMCALPi0PbPb.h"
fa443410 4#include <TAxis.h>
717fe7de 5#include <TChain.h>
6#include <TClonesArray.h>
7#include <TH1F.h>
8#include <TH2F.h>
9#include <TList.h>
10#include <TLorentzVector.h>
11#include "AliAODEvent.h"
12#include "AliAODVertex.h"
13#include "AliAnalysisManager.h"
14#include "AliAnalysisTaskEMCALClusterizeFast.h"
ea3fd2d5 15#include "AliCentrality.h"
ea3fd2d5 16#include "AliEMCALGeoUtils.h"
717fe7de 17#include "AliESDEvent.h"
18#include "AliESDVertex.h"
43cfaa06 19#include "AliLog.h"
ea3fd2d5 20
21ClassImp(AliAnalysisTaskEMCALPi0PbPb)
717fe7de 22
23//________________________________________________________________________
24AliAnalysisTaskEMCALPi0PbPb::AliAnalysisTaskEMCALPi0PbPb()
25 : AliAnalysisTaskSE(),
26 fCentVar(),
27 fCentFrom(0),
28 fCentTo(100),
29 fVtxZMin(-7),
30 fVtxZMax(+7),
43cfaa06 31 fUseQualFlag(1),
717fe7de 32 fClusName(),
33 fOutput(0),
34 fEsdEv(0),
35 fAodEv(0),
43cfaa06 36 fRecPoints(0),
717fe7de 37 fEsdClusters(0),
38 fEsdCells(0),
39 fAodClusters(0),
286b47a5 40 fAodCells(0),
fa443410 41 fPtRanges(0),
42 fHCuts(0x0),
43 fHVertexZ(0x0),
44 fHCent(0x0),
45 fHCentQual(0x0),
46 fHClustEtaPhi(0x0),
47 fHClustEnergyPt(0x0),
48 fHClustEnergySigma(0x0),
49 fHClustNTowEnergyRatio(0x0),
50 fHPionEtaPhi(0x0),
51 fHPionMggPt(0x0),
52 fHPionMggAsym(0x0)
717fe7de 53{
54 // ROOT constructor.
55}
56
ea3fd2d5 57//________________________________________________________________________
58AliAnalysisTaskEMCALPi0PbPb::AliAnalysisTaskEMCALPi0PbPb(const char *name)
59 : AliAnalysisTaskSE(name),
717fe7de 60 fCentVar("V0M"),
61 fCentFrom(0),
62 fCentTo(100),
63 fVtxZMin(-7),
64 fVtxZMax(+7),
43cfaa06 65 fUseQualFlag(1),
717fe7de 66 fClusName(),
67 fOutput(0),
68 fEsdEv(0),
69 fAodEv(0),
43cfaa06 70 fRecPoints(0),
717fe7de 71 fEsdClusters(0),
72 fEsdCells(0),
73 fAodClusters(0),
286b47a5 74 fAodCells(0),
fa443410 75 fPtRanges(0),
76 fHCuts(0x0),
77 fHVertexZ(0x0),
78 fHCent(0x0),
79 fHCentQual(0x0),
80 fHClustEtaPhi(0x0),
81 fHClustEnergyPt(0x0),
82 fHClustEnergySigma(0x0),
83 fHClustNTowEnergyRatio(0x0),
84 fHPionEtaPhi(0x0),
85 fHPionMggPt(0x0),
86 fHPionMggAsym(0x0)
ea3fd2d5 87{
717fe7de 88 // Constructor.
ea3fd2d5 89
ea3fd2d5 90 DefineInput(0, TChain::Class());
ea3fd2d5 91 DefineOutput(1, TList::Class());
fa443410 92 fBranchNames="ESD:AliESDRun.,AliESDHeader.,PrimaryVertex.,SPDVertex.,TPCVertex.,EMCALCells. "
29b496d5 93 "AOD:header,vertices,emcalCells";
ea3fd2d5 94}
ea3fd2d5 95
717fe7de 96//________________________________________________________________________
97AliAnalysisTaskEMCALPi0PbPb::~AliAnalysisTaskEMCALPi0PbPb()
98{
99 // Destructor.
ea3fd2d5 100
717fe7de 101 delete fOutput; fOutput = 0;
fa443410 102 delete fPtRanges; fPtRanges = 0;
ea3fd2d5 103}
717fe7de 104
ea3fd2d5 105//________________________________________________________________________
106void AliAnalysisTaskEMCALPi0PbPb::UserCreateOutputObjects()
107{
717fe7de 108 // Create user objects here.
ea3fd2d5 109
717fe7de 110 fOutput = new TList();
111 fOutput->SetOwner();
ea3fd2d5 112
fa443410 113 fHCuts = new TH1F("hEventCuts","",4,0.5,4.5);
114 fHCuts->GetXaxis()->SetBinLabel(1,"All (PS)");
115 fHCuts->GetXaxis()->SetBinLabel(2,Form("%s: %.0f-%.0f",fCentVar.Data(),fCentFrom,fCentTo));
116 fHCuts->GetXaxis()->SetBinLabel(3,"QFlag");
117 fHCuts->GetXaxis()->SetBinLabel(4,Form("zvtx: %.0f-%.0f",fVtxZMin,fVtxZMax));
118 fOutput->Add(fHCuts);
119 fHVertexZ = new TH1F("hVertexZBeforeCut","",100,-25,25);
120 fHVertexZ->SetXTitle("z [cm]");
121 fOutput->Add(fHVertexZ);
122 fHCent = new TH1F("hCentBeforeCut","",101,-1,100);
123 fHCent->SetXTitle(fCentVar.Data());
124 fHCentQual = new TH1F("hCentAfterCut","",101,-1,100);
125 fHCentQual->SetXTitle(fCentVar.Data());
126 fOutput->Add(fHCentQual);
127
128 fHClustEtaPhi = new TH2F("hClustEtaPhi","",100,-0.8,0.8,100,1.2,2.2);
129 fHClustEtaPhi->SetXTitle("#eta");
130 fHClustEtaPhi->SetYTitle("#varphi");
131 fOutput->Add(fHClustEtaPhi);
132 fHClustEnergyPt = new TH2F("hClustEnergyPt","",250,0,50,250,0,50);
133 fHClustEnergyPt->SetXTitle("E [GeV]");
134 fHClustEnergyPt->SetYTitle("p_{T}^{} [GeV/c]");
135 fOutput->Add(fHClustEnergyPt);
136 fHClustEnergySigma = new TH2F("hClustEnergySigma","",100,0,100,100,0,100);
137 fHClustEnergySigma->SetXTitle("E*#sigma_{max}^{} [GeV*cm]");
138 fHClustEnergySigma->SetYTitle("#sigma_{max}^{} [cm]");
139 //fOutput->Add(fHClustEnergySigma);
140 fHClustNTowEnergyRatio = new TH2F("hClustNTowEnergyRatio","",15,-0.5,14.5,101,-0.05,1.05);
141 fHClustNTowEnergyRatio->SetXTitle("N towers");
142 fHClustNTowEnergyRatio->SetYTitle("Energy ratio");
143 //fOutput->Add(fHClustNTowEnergyRatio);
144
145 fHPionEtaPhi = new TH2F("hPionEtaPhi","",100,-0.8,0.8,100,1.2,2.2);
146 fHPionEtaPhi->SetXTitle("#eta^{#gamma#gamma}");
147 fHPionEtaPhi->SetYTitle("#varphi^{#gamma#gamma}");
148 fOutput->Add(fHPionEtaPhi);
149 fHPionMggPt = new TH2F("hPionMggPt","",100,0,2,100,0,20.0);
150 fHPionMggPt->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
151 fHPionMggPt->SetYTitle("p_{T}^{#gamma#gamma} [GeV/c]");
152 fOutput->Add(fHPionMggPt);
153 fHPionMggAsym = new TH2F("hPionMggAsym","",100,0,2,100,0,1);
154 fHPionMggAsym->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
155 fHPionMggAsym->SetYTitle("Z_{#gamma#gamma} [GeV]");
156 fOutput->Add(fHPionMggAsym);
157 const Int_t nbins = 19;
158 Double_t xbins[nbins] = {0.5,1,1.5,2,2.5,3,3.5,4,4.5,6,7,8,9,10,12.5,15,20,25,50};
159 fPtRanges = new TAxis(nbins-1,xbins);
160 for (Int_t i = 0; i<=nbins; ++i) {
161 fHPionInvMasses[i] = new TH1F(Form("HPionInvMass%d",i),"",100,0,2);
162 fHPionInvMasses[i]->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
163 if (i==0)
164 fHPionInvMasses[i]->SetTitle(Form("0 < p_{T}^{#gamma#gamma} <%.1f",xbins[0]));
165 else if (i==nbins)
166 fHPionInvMasses[i]->SetTitle(Form("p_{T}^{#gamma#gamma} > 50"));
167 else
168 fHPionInvMasses[i]->SetTitle(Form("%.1f < p_{T}^{#gamma#gamma} <%.1f",xbins[i-1],xbins[i]));
169 //cout << i << ": " << fHPionInvMasses[i]->GetTitle() << endl;
170 fOutput->Add(fHPionInvMasses[i]);
171 }
286b47a5 172
717fe7de 173 PostData(1, fOutput);
ea3fd2d5 174}
175
176//________________________________________________________________________
177void AliAnalysisTaskEMCALPi0PbPb::UserExec(Option_t *)
178{
717fe7de 179 // Called for each event.
180
181 if (!InputEvent())
182 return;
183
43cfaa06 184 AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager();
185 fEsdEv = dynamic_cast<AliESDEvent*>(InputEvent());
186 if (fEsdEv) {
187 am->LoadBranch("AliESDRun.");
188 am->LoadBranch("AliESDHeader.");
189 } else {
190 fAodEv = dynamic_cast<AliAODEvent*>(InputEvent());
191 am->LoadBranch("header");
192 }
193
286b47a5 194 Int_t cut = 1;
fa443410 195 fHCuts->Fill(cut++);
286b47a5 196
717fe7de 197 const AliCentrality *centP = InputEvent()->GetCentrality();
198 Double_t cent = centP->GetCentralityPercentileUnchecked(fCentVar);
fa443410 199 fHCent->Fill(cent);
717fe7de 200 if (cent<fCentFrom||cent>fCentTo)
286b47a5 201 return;
43cfaa06 202
fa443410 203 fHCuts->Fill(cut++);
286b47a5 204
43cfaa06 205 if (fUseQualFlag) {
206 if (centP->GetQuality()>0)
207 return;
208 }
286b47a5 209
fa443410 210 fHCentQual->Fill(cent);
211 fHCuts->Fill(cut++);
717fe7de 212
43cfaa06 213 if (fEsdEv) {
fa443410 214 am->LoadBranch("PrimaryVertex.");
215 am->LoadBranch("SPDVertex.");
216 am->LoadBranch("TPCVertex.");
43cfaa06 217 } else {
218 fAodEv = dynamic_cast<AliAODEvent*>(InputEvent());
219 am->LoadBranch("vertices");
220 }
221
222 const AliVVertex *vertex = InputEvent()->GetPrimaryVertex();
717fe7de 223 if (!vertex)
224 return;
225
fa443410 226 fHVertexZ->Fill(vertex->GetZ());
286b47a5 227
717fe7de 228 if(vertex->GetZ()<fVtxZMin||vertex->GetZ()>fVtxZMax)
229 return;
286b47a5 230
fa443410 231 fHCuts->Fill(cut++);
717fe7de 232
43cfaa06 233 fRecPoints = 0; // will be set if fClusName is given and AliAnalysisTaskEMCALClusterizeFast is used
234 fEsdClusters = 0; // will be set if ESD input used and if fRecPoints are not set of if clusters are attached
235 fEsdCells = 0; // will be set if ESD input used
236 fAodClusters = 0; // will be set if AOD input used and if fRecPoints are not set of if clusters are attached
237 // or if fClusName is given and AliAnalysisTaskEMCALClusterizeFast in AOD output mode
238 fAodCells = 0; // will be set if AOD input used
717fe7de 239
43cfaa06 240 // deal with special output from AliAnalysisTaskEMCALClusterizeFast first
241 Bool_t clusattached = 0;
242 Bool_t recalibrated = 0;
243 if (1 && !fClusName.IsNull()) {
244 AliAnalysisTaskEMCALClusterizeFast *cltask = 0;
245 TObjArray *ts = am->GetTasks();
246 cltask = dynamic_cast<AliAnalysisTaskEMCALClusterizeFast*>(ts->FindObject(fClusName));
247 if (cltask && cltask->GetClusters()) {
248 fRecPoints = cltask->GetClusters();
249 clusattached = cltask->GetAttachClusters();
250 if (cltask->GetCalibData()!=0)
251 recalibrated = kTRUE;
252 }
253 }
254 if (1 && AODEvent() && !fClusName.IsNull()) {
717fe7de 255 TList *l = AODEvent()->GetList();
256 TClonesArray *clus = 0;
257 if (l) {
258 clus = dynamic_cast<TClonesArray*>(l->FindObject(fClusName));
259 fAodClusters = clus;
260 }
ea3fd2d5 261 }
717fe7de 262
263 if (fEsdEv) { // ESD input mode
43cfaa06 264 if (1 && (!fRecPoints||clusattached)) {
265 if (!clusattached)
266 am->LoadBranch("CaloClusters");
717fe7de 267 TList *l = fEsdEv->GetList();
268 TClonesArray *clus = 0;
269 if (l) {
270 clus = dynamic_cast<TClonesArray*>(l->FindObject("CaloClusters"));
271 fEsdClusters = clus;
ea3fd2d5 272 }
273 }
43cfaa06 274 if (1) {
275 if (!recalibrated)
276 am->LoadBranch("EMCALCells.");
717fe7de 277 fEsdCells = fEsdEv->GetEMCALCells();
278 }
279 } else if (fAodEv) { // AOD input mode
43cfaa06 280 if (1 && (!fAodClusters || clusattached)) {
281 if (!clusattached)
282 am->LoadBranch("caloClusters");
717fe7de 283 TList *l = fAodEv->GetList();
284 TClonesArray *clus = 0;
285 if (l) {
286 clus = dynamic_cast<TClonesArray*>(l->FindObject("caloClusters"));
287 fAodClusters = clus;
ea3fd2d5 288 }
289 }
43cfaa06 290 if (1) {
291 if (!recalibrated)
292 am->LoadBranch("emcalCells");
717fe7de 293 fAodCells = fAodEv->GetEMCALCells();
294 }
295 } else {
296 AliFatal("Impossible to not have either pointer to ESD or AOD event");
297 }
ea3fd2d5 298
43cfaa06 299 if (1) {
300 AliDebug(2,Form("fRecPoints set: %p", fRecPoints));
301 AliDebug(2,Form("fEsdClusters set: %p", fEsdClusters));
302 AliDebug(2,Form("fEsdCells set: %p", fEsdCells));
303 AliDebug(2,Form("fAodClusters set: %p", fAodClusters));
304 AliDebug(2,Form("fAodCells set: %p", fAodCells));
305 }
306
286b47a5 307 FillCellHists();
308 FillClusHists();
309 FillPionHists();
ea3fd2d5 310
717fe7de 311 PostData(1, fOutput);
ea3fd2d5 312}
313
314//________________________________________________________________________
315void AliAnalysisTaskEMCALPi0PbPb::Terminate(Option_t *)
316{
717fe7de 317 // Terminate called at the end of analysis.
ea3fd2d5 318}
ea3fd2d5 319
286b47a5 320//________________________________________________________________________
321void AliAnalysisTaskEMCALPi0PbPb::FillCellHists()
322{
323 // Fill histograms related to cell properties.
324}
325
326//________________________________________________________________________
327void AliAnalysisTaskEMCALPi0PbPb::FillClusHists()
328{
fa443410 329 // Fill histograms related to clusters.
330
331 Double_t vertex[3] = {0,0,0};
332 InputEvent()->GetPrimaryVertex()->GetXYZ(vertex);
333
286b47a5 334 // Fill histograms related to cluster properties.
fa443410 335 TLorentzVector clusterVec;
336 if (fEsdClusters) {
337 Int_t nclus = fEsdClusters->GetEntries();
338 for(Int_t i = 0; i<nclus; ++i) {
339 AliESDCaloCluster *clus = static_cast<AliESDCaloCluster*>(fEsdClusters->At(i));
340 if (!clus)
341 continue;
342 if (!clus->IsEMCAL())
343 continue;
344 clus->GetMomentum(clusterVec,vertex);
345
346 fHClustEtaPhi->Fill(clusterVec.Eta(),clusterVec.Phi());
347 fHClustEnergyPt->Fill(clusterVec.E(),clusterVec.Pt());
348 //fHClustEnergySigma->Fill();
349 //fHClustNTowEnergyRatio->Fill();
350 }
351 } else if (fAodClusters) {
352 Int_t nclus = fAodClusters->GetEntries();
353 for (Int_t i = 0; i<nclus; ++i) {
354 AliAODCaloCluster *clus = static_cast<AliAODCaloCluster*>(fAodClusters->At(i));
355 if (!clus)
356 continue;
357 if (!clus->IsEMCAL())
358 continue;
359 clus->GetMomentum(clusterVec,vertex);
360
361 fHClustEtaPhi->Fill(clusterVec.Eta(),clusterVec.Phi());
362 fHClustEnergyPt->Fill(clusterVec.E(),clusterVec.Pt());
363 //fHClustEnergySigma->Fill();
364 //fHClustNTowEnergyRatio->Fill();
365 }
366 }
286b47a5 367}
368
369//________________________________________________________________________
370void AliAnalysisTaskEMCALPi0PbPb::FillPionHists()
371{
372 // Fill histograms related to pions.
286b47a5 373
fa443410 374 Double_t vertex[3] = {0,0,0};
375 InputEvent()->GetPrimaryVertex()->GetXYZ(vertex);
376
377 TLorentzVector clusterVec1;
378 TLorentzVector clusterVec2;
379 TLorentzVector pionVec;
380
381 if (fEsdClusters) {
382 Int_t nclus = fEsdClusters->GetEntries();
383 for (Int_t i = 0; i<nclus; ++i) {
384 AliESDCaloCluster *clus1 = static_cast<AliESDCaloCluster*>(fEsdClusters->At(i));
385 if (!clus1)
386 continue;
387 if (!clus1->IsEMCAL())
388 continue;
389 clus1->GetMomentum(clusterVec1,vertex);
390 for (Int_t j = i+1; j<nclus; ++j) {
391 AliESDCaloCluster *clus2 = static_cast<AliESDCaloCluster*>(fEsdClusters->At(j));
392 if (!clus2)
393 continue;
394 if (!clus2->IsEMCAL())
395 continue;
396 clus2->GetMomentum(clusterVec2,vertex);
397 pionVec = clusterVec1 + clusterVec2;
398 Double_t pionZgg = TMath::Abs(clusterVec1.E()-clusterVec2.E())/pionVec.E();
399 fHPionEtaPhi->Fill(pionVec.Eta(),pionVec.Phi());
400 fHPionMggPt->Fill(pionVec.M(),pionVec.Pt());
401 fHPionMggAsym->Fill(pionVec.M(),pionZgg);
402 Int_t bin = fPtRanges->FindBin(pionVec.Pt());
403 fHPionInvMasses[bin]->Fill(pionVec.M());
404 }
405 }
406 } else if (fAodClusters) {
407 Int_t nclus = fAodClusters->GetEntries();
408 for (Int_t i = 0; i<nclus; ++i) {
409 AliAODCaloCluster *clus1 = static_cast<AliAODCaloCluster*>(fAodClusters->At(i));
410 if (!clus1)
411 continue;
412 if (!clus1->IsEMCAL())
413 continue;
414 clus1->GetMomentum(clusterVec1,vertex);
415 for (Int_t j = i+1; j<nclus; ++j) {
416 AliAODCaloCluster *clus2 = static_cast<AliAODCaloCluster*>(fAodClusters->At(j));
417 if (!clus2)
418 continue;
419 if (!clus2->IsEMCAL())
420 continue;
421 clus2->GetMomentum(clusterVec2,vertex);
422 pionVec = clusterVec1 + clusterVec2;
423 Double_t pionZgg = TMath::Abs(clusterVec1.E()-clusterVec2.E())/pionVec.E();
424 fHPionEtaPhi->Fill(pionVec.Eta(),pionVec.Phi());
425 fHPionMggPt->Fill(pionVec.M(),pionVec.Pt());
426 fHPionMggAsym->Fill(pionVec.M(),pionZgg);
427 Int_t bin = fPtRanges->FindBin(pionVec.Pt());
428 fHPionInvMasses[bin]->Fill(pionVec.M());
429 }
430 }
431 }
432}