1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /* $Id: AliAnalysisTaskPi0V2.cxx 55404 2012-03-29 10:10:19Z fca $ */
18 /* AliAnalysisTaskPi0V2.cxx
20 * Template task producing a P_t spectrum and pseudorapidity distribution.
21 * Includes explanations of physics and primary track selections
23 * Instructions for adding histograms can be found below, starting with NEW HISTO
25 * Based on tutorial example from offline pages
26 * Edited by Arvinder Palaha
28 #include "AliAnalysisTaskPi0V2.h"
30 #include "Riostream.h"
39 #include "AliAnalysisTaskSE.h"
40 #include "AliAnalysisManager.h"
42 #include "AliESDtrackCuts.h"
43 #include "AliESDEvent.h"
44 #include "AliESDInputHandler.h"
45 #include "AliAODEvent.h"
46 #include "AliMCEvent.h"
48 #include "AliEventplane.h"
49 #include "AliEMCALGeometry.h"
50 #include "THnSparse.h"
55 ClassImp(AliAnalysisTaskPi0V2)
57 //________________________________________________________________________
58 AliAnalysisTaskPi0V2::AliAnalysisTaskPi0V2() // All data members should be initialised here
59 // :AliAnalysisTaskSE(),
67 fEPV0(-999.), fEPV0A(-999.), fEPV0C(-999.), fEPV0Ar(-999.), fEPV0Cr(-999.), fEPV0r(-999.),
68 fEPV0AR4(-999.), fEPV0AR5(-999.), fEPV0AR6(-999.), fEPV0AR7(-999.), fEPV0CR0(-999.), fEPV0CR1(-999.), fEPV0CR2(-999.), fEPV0CR3(-999.),
69 hAllcentV0(0), hAllcentV0r(0), hAllcentV0A(0), hAllcentV0C(0), hAllcentTPC(0),
70 hEPTPC(0), hresoTPC(0),
71 hEPV0(0), hEPV0A(0), hEPV0C(0), hEPV0Ar(0), hEPV0Cr(0), hEPV0r(0), hEPV0AR4(0), hEPV0AR7(0), hEPV0CR0(0), hEPV0CR3(0),
72 hdifV0A_V0CR0(0), hdifV0A_V0CR3(0), hdifV0ACR0_V0CR3(0), hdifV0C_V0AR4(0), hdifV0C_V0AR7(0), hdifV0AR4_V0AR7(0),
73 hdifV0A_V0C(0), hdifV0A_TPC(0), hdifV0C_TPC(0), hdifV0C_V0A(0),
74 hdifEMC_EP(0), hdifful_EP(0), hdifout_EP(0),
75 fHEPV0r(0), fHEPV0A(0), fHEPV0C(0), fHEPTPC(0)
78 // Dummy constructor ALWAYS needed for I/O.
79 DefineInput(0, TChain::Class());
80 DefineOutput(1, TList::Class()); // for output list
83 //________________________________________________________________________
84 AliAnalysisTaskPi0V2::AliAnalysisTaskPi0V2(const char *name) // All data members should be initialised here
85 :AliAnalysisTaskSE(name),
92 fEPV0(-999.), fEPV0A(-999.), fEPV0C(-999.), fEPV0Ar(-999.), fEPV0Cr(-999.), fEPV0r(-999.),
93 fEPV0AR4(-999.), fEPV0AR5(-999.), fEPV0AR6(-999.), fEPV0AR7(-999.), fEPV0CR0(-999.), fEPV0CR1(-999.), fEPV0CR2(-999.), fEPV0CR3(-999.),
94 hAllcentV0(0), hAllcentV0r(0), hAllcentV0A(0), hAllcentV0C(0), hAllcentTPC(0),
95 hEPTPC(0), hresoTPC(0),
96 hEPV0(0), hEPV0A(0), hEPV0C(0), hEPV0Ar(0), hEPV0Cr(0), hEPV0r(0), hEPV0AR4(0), hEPV0AR7(0), hEPV0CR0(0), hEPV0CR3(0),
97 hdifV0A_V0CR0(0), hdifV0A_V0CR3(0), hdifV0ACR0_V0CR3(0), hdifV0C_V0AR4(0), hdifV0C_V0AR7(0), hdifV0AR4_V0AR7(0),
98 hdifV0A_V0C(0), hdifV0A_TPC(0), hdifV0C_TPC(0), hdifV0C_V0A(0),
99 hdifEMC_EP(0), hdifful_EP(0), hdifout_EP(0),
100 fHEPV0r(0), fHEPV0A(0), fHEPV0C(0), fHEPTPC(0)
103 // Define input and output slots here (never in the dummy constructor)
104 // Input slot #0 works with a TChain - it is connected to the default input container
105 // Output slot #1 writes into a TH1 container
106 DefineInput(0, TChain::Class());
107 DefineOutput(1, TList::Class()); // for output list
110 //________________________________________________________________________
111 AliAnalysisTaskPi0V2::~AliAnalysisTaskPi0V2()
113 // Destructor. Clean-up the output list, but not the histograms that are put inside
114 // (the list is owner and will clean-up these histograms). Protect in PROOF case.
115 if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
119 //_____________________________________________________________________
120 Double_t AliAnalysisTaskPi0V2::GetMaxCellEnergy(const AliVCluster *cluster, Short_t &id) const
122 // Get maximum energy of attached cell.
126 AliVCaloCells *cells = 0;
128 cells = fESD->GetEMCALCells();
130 // cells = fAOD->GetEMCALCells();
135 const Int_t ncells = cluster->GetNCells();
136 for (Int_t i=0; i<ncells; i++) {
137 Double_t e = cells->GetCellAmplitude(TMath::Abs(cluster->GetCellAbsId(i)));
140 id = cluster->GetCellAbsId(i);
145 //_____________________________________________________________________
146 Double_t AliAnalysisTaskPi0V2::GetCrossEnergy(const AliVCluster *cluster, Short_t &idmax) const
148 // Calculate the energy of cross cells around the leading cell.
150 AliVCaloCells *cells = 0;
152 cells = fESD->GetEMCALCells();
154 // cells = fAOD->GetEMCALCells();
158 AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
171 Double_t crossEnergy = 0;
173 geom->GetCellIndex(idmax,iSupMod,iTower,iIphi,iIeta);
174 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphis,ietas);
176 Int_t ncells = cluster->GetNCells();
177 for (Int_t i=0; i<ncells; i++) {
178 Int_t cellAbsId = cluster->GetCellAbsId(i);
179 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
180 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
181 Int_t aphidiff = TMath::Abs(iphi-iphis);
184 Int_t aetadiff = TMath::Abs(ieta-ietas);
187 if ( (aphidiff==1 && aetadiff==0) ||
188 (aphidiff==0 && aetadiff==1) ) {
189 crossEnergy += cells->GetCellAmplitude(cellAbsId);
195 //_____________________________________________________________________
196 Bool_t AliAnalysisTaskPi0V2::IsWithinFiducialVolume(Short_t id) const
198 // Check if cell is within given fiducial volume.
200 Double_t fNFiducial = 1;
209 Bool_t okrow = kFALSE;
210 Bool_t okcol = kFALSE;
212 AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
216 Int_t cellAbsId = id;
217 geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
218 geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
222 if (iphi >= fNFiducial && iphi < 24-fNFiducial)
225 if (iphi >= fNFiducial && iphi < 12-fNFiducial)
229 Bool_t noEMCALBorderAtEta0 = kTRUE;
230 if (!noEMCALBorderAtEta0) {
231 if (ieta > fNFiducial && ieta < 48-fNFiducial)
235 if (ieta >= fNFiducial)
238 if (ieta < 48-fNFiducial)
247 //______________________________________________________________________
248 Bool_t AliAnalysisTaskPi0V2::IsGoodCluster(const AliESDCaloCluster *c) const
254 if(c->GetNCells() < 2)
261 Double_t maxE = GetMaxCellEnergy(c, id);
262 if((1. - GetCrossEnergy(c,id) / maxE) > 0.97)
266 c->GetPosition(pos1);
267 TVector3 clsPos(pos1);
268 Double_t eta = clsPos.Eta();
270 if(eta > 0.65 || eta < -0.65)
273 if (!IsWithinFiducialVolume(id))
276 if(c->GetM02() >0.5 )
284 //_____________________________________________________________________
285 Bool_t AliAnalysisTaskPi0V2::IsGoodPion(const TLorentzVector &p1, const TLorentzVector &p2) const
290 // Double_t asym = TMath::Abs(p1.E()-p2.E())/(p1.E()+p2.E());
294 // if (TMath::Abs(p1.Eta()-p2.Eta())>0.5)
299 Double_t eta = pion.Eta();
307 //_______________________________________________________________________
308 void AliAnalysisTaskPi0V2::FillPion(const TLorentzVector& p1, const TLorentzVector& p2, Double_t EPV0r, Double_t EPV0A, Double_t EPV0C, Double_t EPTPC)
312 if (!IsGoodPion(p1,p2))
317 Double_t mass = pion.M();
318 Double_t pt = pion.Pt();
319 Double_t phi = pion.Phi();
321 Double_t dphiV0 = phi-EPV0r;
322 Double_t dphiV0A = phi-EPV0A;
323 Double_t dphiV0C = phi-EPV0C;
324 Double_t dphiTPC = phi-EPTPC;
326 Double_t cos2phiV0 = TMath::Cos(2.*(dphiV0));
327 Double_t cos2phiV0A = TMath::Cos(2.*(dphiV0A));
328 Double_t cos2phiV0C = TMath::Cos(2.*(dphiV0C));
329 Double_t cos2phiTPC = TMath::Cos(2.*(dphiTPC));
331 dphiV0 = TVector2::Phi_0_2pi(dphiV0); if(dphiV0 >TMath::Pi()) dphiV0 -= TMath::Pi();
332 dphiV0A = TVector2::Phi_0_2pi(dphiV0A); if(dphiV0A >TMath::Pi()) dphiV0A -= TMath::Pi();
333 dphiV0C = TVector2::Phi_0_2pi(dphiV0C); if(dphiV0C >TMath::Pi()) dphiV0C -= TMath::Pi();
334 dphiTPC = TVector2::Phi_0_2pi(dphiTPC); if(dphiTPC >TMath::Pi()) dphiTPC -= TMath::Pi();
336 //cout<<"cos2V0: "<<cos2phiV0<<" cos2V0A: "<<cos2phiV0A<<" cos2V0C: "<<cos2phiV0C<<endl;
337 //cout<<mass<<" "<<pt<<" "<<phi<<" "<<endl;
338 //cout<<" dphiV0: "<<dphiV0<<" dphiV0A: "<<dphiV0A<<" dphiV0C: "<<dphiV0C<<"+++++++"<<endl;
340 Double_t xV0[5]; // Match ndims in fH V0 EP
343 xV0[2] = fCentrality;
348 Double_t xV0A[5]; // Match ndims in fH V0A EP
351 xV0A[2] = fCentrality;
353 xV0A[4] = cos2phiV0A;
356 Double_t xV0C[5]; // Match ndims in fH V0C EP
359 xV0C[2] = fCentrality;
361 xV0C[4] = cos2phiV0C;
364 Double_t xTPC[5]; // Match ndims in fH TPC EP
367 xTPC[2] = fCentrality;
369 xTPC[4] = cos2phiTPC;
374 //_________________________________________________________________________________________________
375 void AliAnalysisTaskPi0V2::GetMom(TLorentzVector& p, const AliESDCaloCluster *c, Double_t *vertex)
377 // Calculate momentum.
379 c->GetPosition(posMom);
380 TVector3 clsPos2(posMom);
383 Double_t r = clsPos2.Perp();
384 Double_t eta = clsPos2.Eta();
385 Double_t phi = clsPos2.Phi();
388 pos.SetPtEtaPhi(r,eta,phi);
390 if (vertex) { //calculate direction relative to vertex
394 Double_t rad = pos.Mag();
395 p.SetPxPyPzE(e*pos.x()/rad, e*pos.y()/rad, e*pos.z()/rad, e);
398 //________________________________________________________________________
399 void AliAnalysisTaskPi0V2::UserCreateOutputObjects()
402 // Called once (on the worker node)
404 fOutput = new TList();
405 fOutput->SetOwner(); // IMPORTANT!
407 hEPTPC = new TH2F("hEPTPC", "EPTPC vs cent", 100, 0., 100., 100, 0., TMath::Pi());
408 hresoTPC = new TH2F("hresoTPC", "TPc reso vs cent", 100, 0., 100., 100, 0., 1.);
409 hEPV0 = new TH2F("hEPV0", "EPV0 vs cent", 100, 0., 100., 100, 0., TMath::Pi());
410 hEPV0A = new TH2F("hEPV0A", "EPV0A vs cent", 100, 0., 100., 100, 0., TMath::Pi());
411 hEPV0C = new TH2F("hEPV0C", "EPV0C vs cent", 100, 0., 100., 100, 0., TMath::Pi());
412 hEPV0Ar = new TH2F("hEPV0Ar", "EPV0Ar vs cent", 100, 0., 100., 100, 0., TMath::Pi());
413 hEPV0Cr = new TH2F("hEPV0Cr", "EPV0Cr vs cent", 100, 0., 100., 100, 0., TMath::Pi());
414 hEPV0r = new TH2F("hEPV0r", "EPV0r vs cent", 100, 0., 100., 100, 0., TMath::Pi());
415 hEPV0AR4 = new TH2F("hEPV0AR4", "EPV0AR4 vs cent", 100, 0., 100., 100, 0., TMath::Pi());
416 hEPV0AR7 = new TH2F("hEPV0AR7", "EPV0AR7 vs cent", 100, 0., 100., 100, 0., TMath::Pi());
417 hEPV0CR0 = new TH2F("hEPV0CR0", "EPV0CR0 vs cent", 100, 0., 100., 100, 0., TMath::Pi());
418 hEPV0CR3 = new TH2F("hEPV0CR3", "EPV0CR3 vs cent", 100, 0., 100., 100, 0., TMath::Pi());
419 fOutput->Add(hEPTPC);
420 fOutput->Add(hresoTPC);
422 fOutput->Add(hEPV0A);
423 fOutput->Add(hEPV0C);
424 fOutput->Add(hEPV0Ar);
425 fOutput->Add(hEPV0Cr);
426 fOutput->Add(hEPV0r);
427 fOutput->Add(hEPV0AR4);
428 fOutput->Add(hEPV0AR7);
429 fOutput->Add(hEPV0CR0);
430 fOutput->Add(hEPV0CR3);
432 hdifV0A_V0CR0 = new TH2F("hdifV0A_V0CR0", "EP A-R0 ", 100, 0., 100., 100, -1., 1.);
433 hdifV0A_V0CR3 = new TH2F("hdifV0A_V0CR3", "EP A-R3 ", 100, 0., 100., 100, -1., 1.);
434 hdifV0ACR0_V0CR3 = new TH2F("hdifV0ACR0_V0CR3", "EP R0-R3 ", 100, 0., 100., 100, -1., 1.);
435 hdifV0C_V0AR4 = new TH2F("hdifV0C_V0AR4", "EP C-R4 ", 100, 0., 100., 100, -1., 1.);
436 hdifV0C_V0AR7 = new TH2F("hdifV0C_V0AR7", "EP C-R7 ", 100, 0., 100., 100, -1., 1.);
437 hdifV0AR4_V0AR7 = new TH2F("hdifV0AR4_V0AR7", "EP R4-R7 ", 100, 0., 100., 100, -1., 1.);
438 fOutput->Add(hdifV0A_V0CR0);
439 fOutput->Add(hdifV0A_V0CR3);
440 fOutput->Add(hdifV0ACR0_V0CR3);
441 fOutput->Add(hdifV0C_V0AR4);
442 fOutput->Add(hdifV0C_V0AR7);
443 fOutput->Add(hdifV0AR4_V0AR7);
445 hdifV0A_V0C = new TH2F("hdifV0A_V0C", "EP A-C ", 100, 0., 100., 100, -1., 1.);
446 hdifV0A_TPC = new TH2F("hdifV0A_TPC", "EP A-TPC", 100, 0., 100., 100, -1., 1.);
447 hdifV0C_TPC = new TH2F("hdifV0C_TPC", "EP C-TPC", 100, 0., 100., 100, -1., 1.);
448 hdifV0C_V0A = new TH2F("hdifV0C_V0A", "EP C-A ", 100, 0., 100., 100, -1., 1.);
449 fOutput->Add(hdifV0A_V0C);
450 fOutput->Add(hdifV0A_TPC);
451 fOutput->Add(hdifV0C_TPC);
452 fOutput->Add(hdifV0C_V0A);
456 hdifEMC_EP = new TH3F("hdifEMC_EP", "dif phi in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
457 hdifful_EP = new TH3F("hdifful_EP", "dif phi in full with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
458 hdifout_EP = new TH3F("hdifout_EP", "dif phi NOT in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
459 fOutput->Add(hdifEMC_EP);
460 fOutput->Add(hdifful_EP);
461 fOutput->Add(hdifout_EP);
463 hAllcentV0 = new TH1F("hAllcentV0", "All cent EP V0", 100, 0., TMath::Pi());
464 hAllcentV0r = new TH1F("hAllcentV0r", "All cent EP V0r", 100, 0., TMath::Pi());
465 hAllcentV0A = new TH1F("hAllcentV0A", "All cent EP V0A", 100, 0., TMath::Pi());
466 hAllcentV0C = new TH1F("hAllcentV0C", "All cent EP V0C", 100, 0., TMath::Pi());
467 hAllcentTPC = new TH1F("hAllcentTPC", "All cent EP TPC", 100, 0., TMath::Pi());
468 fOutput->Add(hAllcentV0);
469 fOutput->Add(hAllcentV0r);
470 fOutput->Add(hAllcentV0A);
471 fOutput->Add(hAllcentV0C);
472 fOutput->Add(hAllcentTPC);
474 const Int_t ndims = 5;
475 Int_t nMgg=500, nPt=40, nCent=20, nDeltaPhi=315, ncos2phi=500;
476 Int_t bins[ndims] = {nMgg, nPt, nCent, nDeltaPhi, ncos2phi};
477 Double_t xmin[ndims] = { 0, 0., 0, 0., -1.};
478 Double_t xmax[ndims] = { 0.5, 20., 100, 3.15, 1.};
479 fHEPV0r = new THnSparseF("fHEPV0r", "Flow histogram EPV0", ndims, bins, xmin, xmax);
480 fHEPV0A = new THnSparseF("fHEPV0A", "Flow histogram EPV0A", ndims, bins, xmin, xmax);
481 fHEPV0C = new THnSparseF("fHEPV0C", "Flow histogram EPV0C", ndims, bins, xmin, xmax);
482 fHEPTPC = new THnSparseF("fHEPTPC", "Flow histogram EPTPC", ndims, bins, xmin, xmax);
483 fOutput->Add(fHEPV0r);
484 fOutput->Add(fHEPV0A);
485 fOutput->Add(fHEPV0C);
486 fOutput->Add(fHEPTPC);
490 PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram
493 //________________________________________________________________________
494 void AliAnalysisTaskPi0V2::UserExec(Option_t *)
497 // Called for each event
500 // Create pointer to reconstructed event
501 AliVEvent *event = InputEvent();
502 if (!event) { Printf("ERROR: Could not retrieve event"); return; }
504 // create pointer to event
505 fESD = dynamic_cast<AliESDEvent*>(event);
507 AliError("Cannot get the ESD event");
510 const AliESDVertex* fvertex = fESD->GetPrimaryVertex();
512 if(!(fvertex->GetStatus())) return;
513 // if vertex is from spd vertexZ, require more stringent cut
514 if (fvertex->IsFromVertexerZ()) {
515 if (fvertex->GetDispersion()>0.02 || fvertex->GetZRes()>0.25 ) return; // bad vertex from VertexerZ
517 Double_t vertex[3] = {fvertex->GetX(), fvertex->GetY(), fvertex->GetZ()};
519 if(fESD->GetCentrality()) {
521 fESD->GetCentrality()->GetCentralityPercentile("V0M");
524 AliEventplane *ep = fESD->GetEventplane();
526 if (ep->GetQVector())
527 fEPTPC = ep->GetQVector()->Phi()/2. ;
530 if (ep->GetQsub1()&&ep->GetQsub2())
531 fEPTPCreso = TMath::Cos(2.*(ep->GetQsub1()->Phi()/2.-ep->GetQsub2()->Phi()/2.));
535 fEPV0 = ep->GetEventplane("V0", fESD);
536 fEPV0A = ep->GetEventplane("V0A", fESD);
537 fEPV0C = ep->GetEventplane("V0C", fESD);
539 Double_t qxr=0, qyr=0;
540 fEPV0Ar = ep->CalculateVZEROEventPlane(fESD, 4, 5, 2, qxr, qyr);
541 fEPV0Cr = ep->CalculateVZEROEventPlane(fESD, 2, 3, 2, qx, qy);
544 fEPV0r = TMath::ATan2(qyr,qxr)/2.;
545 fEPV0AR4 = ep->CalculateVZEROEventPlane(fESD, 4, 2, qx, qy);
546 fEPV0AR5 = ep->CalculateVZEROEventPlane(fESD, 5, 2, qx, qy);
547 fEPV0AR6 = ep->CalculateVZEROEventPlane(fESD, 6, 2, qx, qy);
548 fEPV0AR7 = ep->CalculateVZEROEventPlane(fESD, 7, 2, qx, qy);
549 fEPV0CR0 = ep->CalculateVZEROEventPlane(fESD, 0, 2, qx, qy);
550 fEPV0CR1 = ep->CalculateVZEROEventPlane(fESD, 1, 2, qx, qy);
551 fEPV0CR2 = ep->CalculateVZEROEventPlane(fESD, 2, 2, qx, qy);
552 fEPV0CR3 = ep->CalculateVZEROEventPlane(fESD, 3, 2, qx, qy);
555 //cout<<" fEPV0:"<<fEPV0<<" fEPV0A:"<<fEPV0A<<" fEPV0C:"<<fEPV0C<<" fEPV0Ar:"<<fEPV0Ar<<" fEPV0Cr:"<<fEPV0Cr<<" fEPV0r:"<<fEPV0AR4<<" fEPV0AR7:"<<fEPV0AR7<<" fEPV0CR0:"<<fEPV0CR0<<" fEPV0CR3:"<<fEPV0CR3<<"--------------------------------------------"<<endl;
557 if(fEPV0r<-2. || fEPV0Ar<-2. || fEPV0Cr<-2.) return;
560 if(fEPV0A<-2. || fEPV0C<-2. || fEPV0AR4<-2. || fEPV0AR7<-2. || fEPV0CR0<-2. || fEPV0CR3<-2.) return;
562 fEPV0 = TVector2::Phi_0_2pi(fEPV0); if(fEPV0>TMath::Pi()) fEPV0 = fEPV0 - TMath::Pi();
563 fEPV0r = TVector2::Phi_0_2pi(fEPV0r); if(fEPV0r>TMath::Pi()) fEPV0r = fEPV0r - TMath::Pi();
564 fEPV0A = TVector2::Phi_0_2pi(fEPV0A); if(fEPV0A>TMath::Pi()) fEPV0A = fEPV0A - TMath::Pi();
565 fEPV0C = TVector2::Phi_0_2pi(fEPV0C); if(fEPV0C>TMath::Pi()) fEPV0C = fEPV0C - TMath::Pi();
566 fEPV0Ar = TVector2::Phi_0_2pi(fEPV0Ar); if(fEPV0Ar>TMath::Pi()) fEPV0Ar = fEPV0Ar - TMath::Pi();
567 fEPV0Cr = TVector2::Phi_0_2pi(fEPV0Cr); if(fEPV0Cr>TMath::Pi()) fEPV0Cr = fEPV0Cr - TMath::Pi();
568 fEPV0AR4 = TVector2::Phi_0_2pi(fEPV0AR4); if(fEPV0AR4>TMath::Pi()) fEPV0AR4 = fEPV0AR4 - TMath::Pi();
569 fEPV0AR7 = TVector2::Phi_0_2pi(fEPV0AR7); if(fEPV0AR7>TMath::Pi()) fEPV0AR7 = fEPV0AR7 - TMath::Pi();
570 fEPV0CR0 = TVector2::Phi_0_2pi(fEPV0CR0); if(fEPV0CR0>TMath::Pi()) fEPV0CR0 = fEPV0CR0 - TMath::Pi();
571 fEPV0CR3 = TVector2::Phi_0_2pi(fEPV0CR3); if(fEPV0CR3>TMath::Pi()) fEPV0CR3 = fEPV0CR3 - TMath::Pi();
573 //cout<<" EPTPC: "<<fEPTPC<<" reso: "<<fEPTPCreso<<" -------------------"<<endl;
574 //cout<<" cent: "<<fCentrality<<" fEPV0:"<<fEPV0<<" fEPV0A:"<<fEPV0A<<" fEPV0C:"<<fEPV0C<<" fEPV0Ar:"<<fEPV0Ar<<" fEPV0Cr:"<<fEPV0Cr<<" fEPV0r:"<<fEPV0AR4<<" fEPV0AR7:"<<fEPV0AR7<<" fEPV0CR0:"<<fEPV0CR0<<" fEPV0CR3:"<<fEPV0CR3<<"--------------------------------------------"<<endl;
576 hEPTPC->Fill(fCentrality, fEPTPC);
577 if(fEPTPCreso!=-1) hresoTPC->Fill(fCentrality, fEPTPCreso);
578 hEPV0->Fill(fCentrality, fEPV0);
579 hEPV0A->Fill(fCentrality, fEPV0A);
580 hEPV0C->Fill(fCentrality, fEPV0C);
581 hEPV0Ar->Fill(fCentrality, fEPV0Ar);
582 hEPV0Cr->Fill(fCentrality, fEPV0Cr);
583 hEPV0r->Fill(fCentrality, fEPV0r);
584 hEPV0AR4->Fill(fCentrality, fEPV0AR4);
585 hEPV0AR7->Fill(fCentrality, fEPV0AR7);
586 hEPV0CR0->Fill(fCentrality, fEPV0CR0);
587 hEPV0CR3->Fill(fCentrality, fEPV0CR3);
589 hAllcentV0->Fill(fEPV0);
590 hAllcentV0r->Fill(fEPV0r);
591 hAllcentV0A->Fill(fEPV0A);
592 hAllcentV0C->Fill(fEPV0C);
593 hAllcentTPC->Fill(fEPTPC);
595 hdifV0A_V0CR0->Fill(fCentrality, TMath::Cos(2.*(fEPV0A - fEPV0CR0)));
596 hdifV0A_V0CR3->Fill(fCentrality, TMath::Cos(2.*(fEPV0A - fEPV0CR3)));
597 hdifV0ACR0_V0CR3->Fill(fCentrality, TMath::Cos(2*(fEPV0CR0 - fEPV0CR3)));
598 hdifV0C_V0AR4->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0AR4)));
599 hdifV0C_V0AR7->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0AR7)));
600 hdifV0AR4_V0AR7->Fill(fCentrality, TMath::Cos(2*(fEPV0AR4 - fEPV0AR7)));
602 hdifV0A_V0C->Fill(fCentrality, TMath::Cos(2*(fEPV0A - fEPV0C)));
603 hdifV0A_TPC->Fill(fCentrality, TMath::Cos(2*(fEPV0A - fEPTPC)));
604 hdifV0C_TPC->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPTPC)));
605 hdifV0C_V0A->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0A)));
606 // Cluster loop for reconstructed event
608 Int_t nCluster = fESD->GetNumberOfCaloClusters();
609 for(Int_t i=0; i<nCluster; ++i){
610 AliESDCaloCluster *c1 = fESD->GetCaloCluster(i);
611 if(!IsGoodCluster(c1)) continue;
612 for(Int_t j=i+1; j<nCluster; ++j){
613 AliESDCaloCluster *c2 = fESD->GetCaloCluster(j);
614 if(!IsGoodCluster(c2)) continue;
616 GetMom(p1, c1, vertex);
618 GetMom(p2, c2, vertex);
619 FillPion(p1, p2, fEPV0r, fEPV0A, fEPV0C, fEPTPC);
624 Int_t nTrack = fESD->GetNumberOfTracks();
625 for(Int_t i=0; i<nTrack; ++i){
626 AliESDtrack* esdtrack = fESD->GetTrack(i); // pointer to reconstructed to track
628 AliError(Form("ERROR: Could not retrieve esdtrack %d",i));
631 Double_t tPhi = esdtrack->Phi();
632 Double_t tPt = esdtrack->Pt();
634 Double_t difTrack = TVector2::Phi_0_2pi(tPhi-fEPV0); if(difTrack >TMath::Pi()) difTrack -= TMath::Pi();
635 if(esdtrack->IsEMCAL()){
636 hdifEMC_EP->Fill(fCentrality, difTrack, tPt);
638 hdifout_EP->Fill(fCentrality, difTrack, tPt);
640 hdifful_EP->Fill(fCentrality, difTrack, tPt);
643 // NEW HISTO should be filled before this point, as PostData puts the
644 // information for this iteration of the UserExec in the container
645 PostData(1, fOutput);
649 //________________________________________________________________________
650 void AliAnalysisTaskPi0V2::Terminate(Option_t *)
652 // Draw result to screen, or perform fitting, normalizations
653 // Called once at the end of the query
654 // fOutput = dynamic_cast<TList*> (GetOutputData(1));
655 // if(!fOutput) { Printf("ERROR: could not retrieve TList fOutput"); return; }
657 // Get the physics selection histograms with the selection statistics
658 //AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
659 //AliESDInputHandler *inputH = dynamic_cast<AliESDInputHandler*>(mgr->GetInputEventHandler());
660 //TH2F *histStat = (TH2F*)inputH->GetStatistics();
663 // NEW HISTO should be retrieved from the TList container in the above way,
664 // so it is available to draw on a canvas such as below