]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/EMCALTasks/AliAnalysisTaskPi0V2.cxx
remove obsolete restriction on aplication of time cuts in case of AOD analysis
[u/mrichter/AliRoot.git] / PWGGA / EMCALTasks / AliAnalysisTaskPi0V2.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id: AliAnalysisTaskPi0V2.cxx 55404 2012-03-29 10:10:19Z fca $ */
17
18 /* AliAnalysisTaskPi0V2.cxx
19  *
20  * Template task producing a P_t spectrum and pseudorapidity distribution.
21  * Includes explanations of physics and primary track selections
22  *
23  * Instructions for adding histograms can be found below, starting with NEW HISTO
24  *
25  * Based on tutorial example from offline pages
26  * Edited by Arvinder Palaha
27  */
28 #include "AliAnalysisTaskPi0V2.h"
29
30 #include "Riostream.h"
31 #include "TChain.h"
32 #include "TTree.h"
33 #include "TH1F.h"
34 #include "TH2F.h"
35 #include "TH3F.h"
36 #include "TCanvas.h"
37 #include "TList.h"
38
39 #include "AliAnalysisTaskSE.h"
40 #include "AliAnalysisManager.h"
41 #include "AliStack.h"
42 #include "AliESDtrackCuts.h"
43 #include "AliESDEvent.h"
44 #include "AliESDInputHandler.h"
45 #include "AliAODEvent.h"
46 #include "AliMCEvent.h"
47
48 #include "AliEventplane.h"
49 #include "AliEMCALGeometry.h"
50 #include "THnSparse.h"
51 #include "TClonesArray.h"
52 #include "TString.h"
53
54 using std::cout;
55 using std::endl;
56
57 ClassImp(AliAnalysisTaskPi0V2)
58
59 //________________________________________________________________________
60 AliAnalysisTaskPi0V2::AliAnalysisTaskPi0V2(const char *name) // All data members should be initialised here
61    :AliAnalysisTaskSE(name),
62     fOutput(0),
63     fESD(0),
64     fTracksName("PicoTrack"),
65     fTrigClass("CVLN_|CSEMI_|CCENT|CVHN"),
66     fTracks(0),
67     fRunNumber(-999.),
68     fEvtSelect(1),
69     fVtxCut(10.),
70     fNcellCut(2), fECut(1), fEtaCut(0.65), fM02Cut(0.5), fPi0AsyCut(0),
71     fCentrality(99.),
72     fEPTPC(-999.),
73     fEPTPCreso(0.), 
74     fEPV0(-999.), fEPV0A(-999.), fEPV0C(-999.), fEPV0Ar(-999.), fEPV0Cr(-999.), fEPV0r(-999.),
75     fEPV0AR4(-999.), fEPV0AR5(-999.), fEPV0AR6(-999.), fEPV0AR7(-999.), fEPV0CR0(-999.), fEPV0CR1(-999.), fEPV0CR2(-999.), fEPV0CR3(-999.),
76     hEvtCount(0), hAllcentV0(0), hAllcentV0r(0), hAllcentV0A(0), hAllcentV0C(0), hAllcentTPC(0),
77     h2DcosV0r(0), h2DsinV0r(0), h2DcosV0A(0), h2DsinV0A(0), h2DcosV0C(0), h2DsinV0C(0), h2DcosTPC(0), h2DsinTPC(0), 
78     hEPTPC(0), hresoTPC(0),
79     hEPV0(0), hEPV0A(0), hEPV0C(0), hEPV0Ar(0), hEPV0Cr(0), hEPV0r(0), hEPV0AR4(0), hEPV0AR7(0), hEPV0CR0(0), hEPV0CR3(0),
80     hdifV0A_V0CR0(0), hdifV0A_V0CR3(0), hdifV0ACR0_V0CR3(0), hdifV0C_V0AR4(0), hdifV0C_V0AR7(0), hdifV0AR4_V0AR7(0),
81     hdifV0A_V0C(0), hdifV0A_TPC(0), hdifV0C_TPC(0), hdifV0C_V0A(0), 
82     hdifEMC_EPV0(0), hdifEMC_EPV0A(0), hdifEMC_EPV0C(0), hdifful_EPV0(0), hdifful_EPV0A(0), hdifful_EPV0C(0), 
83     hdifout_EPV0(0), hdifout_EPV0A(0), hdifout_EPV0C(0), hdifEMC_EPTPC(0), hdifful_EPTPC(0), hdifout_EPTPC(0),
84     hdifClus_EPV0(0), hdifClus_EPV0A(0), hdifClus_EPV0C(0), hdifClus_EPTPC(0),
85     fHEPV0r(0), fHEPV0A(0), fHEPV0C(0), fHEPTPC(0)
86
87 {
88     // Dummy constructor ALWAYS needed for I/O.
89     DefineInput(0, TChain::Class());
90     DefineOutput(1, TList::Class());                                            // for output list
91 }
92
93 //________________________________________________________________________
94 AliAnalysisTaskPi0V2::AliAnalysisTaskPi0V2() // All data members should be initialised here
95    :AliAnalysisTaskSE("default_name"),
96     fOutput(0),
97     fESD(0),
98     fTracksName("PicoTracks"),
99     fTrigClass("CVLN_|CSEMI_|CCENT|CVHN"),
100     fTracks(0),
101     fRunNumber(-999.),
102     fEvtSelect(1),
103     fVtxCut(10.),
104     fNcellCut(2), fECut(1), fEtaCut(0.65), fM02Cut(0.5), fPi0AsyCut(0),
105     fCentrality(99.),
106     fEPTPC(-999.),
107     fEPTPCreso(0.),
108     fEPV0(-999.), fEPV0A(-999.), fEPV0C(-999.), fEPV0Ar(-999.), fEPV0Cr(-999.), fEPV0r(-999.),
109     fEPV0AR4(-999.), fEPV0AR5(-999.), fEPV0AR6(-999.), fEPV0AR7(-999.), fEPV0CR0(-999.), fEPV0CR1(-999.), fEPV0CR2(-999.), fEPV0CR3(-999.),
110     hEvtCount(0), hAllcentV0(0), hAllcentV0r(0), hAllcentV0A(0), hAllcentV0C(0), hAllcentTPC(0),
111     h2DcosV0r(0), h2DsinV0r(0), h2DcosV0A(0), h2DsinV0A(0), h2DcosV0C(0), h2DsinV0C(0), h2DcosTPC(0), h2DsinTPC(0), 
112     hEPTPC(0), hresoTPC(0),
113     hEPV0(0), hEPV0A(0), hEPV0C(0), hEPV0Ar(0), hEPV0Cr(0), hEPV0r(0), hEPV0AR4(0), hEPV0AR7(0), hEPV0CR0(0), hEPV0CR3(0),
114     hdifV0A_V0CR0(0), hdifV0A_V0CR3(0), hdifV0ACR0_V0CR3(0), hdifV0C_V0AR4(0), hdifV0C_V0AR7(0), hdifV0AR4_V0AR7(0),
115     hdifV0A_V0C(0), hdifV0A_TPC(0), hdifV0C_TPC(0), hdifV0C_V0A(0),  
116     hdifEMC_EPV0(0), hdifEMC_EPV0A(0), hdifEMC_EPV0C(0), hdifful_EPV0(0), hdifful_EPV0A(0), hdifful_EPV0C(0), 
117     hdifout_EPV0(0), hdifout_EPV0A(0), hdifout_EPV0C(0), hdifEMC_EPTPC(0), hdifful_EPTPC(0), hdifout_EPTPC(0),
118     hdifClus_EPV0(0), hdifClus_EPV0A(0), hdifClus_EPV0C(0), hdifClus_EPTPC(0),
119     fHEPV0r(0), fHEPV0A(0), fHEPV0C(0), fHEPTPC(0)
120 {
121     // Constructor
122     // Define input and output slots here (never in the dummy constructor)
123     // Input slot #0 works with a TChain - it is connected to the default input container
124     // Output slot #1 writes into a TH1 container
125     DefineInput(0, TChain::Class());
126     DefineOutput(1, TList::Class());                                            // for output list
127 }
128
129 //________________________________________________________________________
130 AliAnalysisTaskPi0V2::~AliAnalysisTaskPi0V2()
131 {
132     // Destructor. Clean-up the output list, but not the histograms that are put inside
133     // (the list is owner and will clean-up these histograms). Protect in PROOF case.
134      delete fOutput;
135 }
136 //_____________________________________________________________________
137 Double_t AliAnalysisTaskPi0V2::GetMaxCellEnergy(const AliVCluster *cluster, Short_t &id) const
138 {
139   // Get maximum energy of attached cell.
140
141   id = -1;
142
143   AliVCaloCells *cells = 0;
144   if (fESD)
145     cells = fESD->GetEMCALCells();
146   if (!cells)
147     return 0;
148
149   Double_t maxe = 0;
150   const Int_t ncells = cluster->GetNCells();
151   for (Int_t i=0; i<ncells; i++) {
152     Double_t e = cells->GetCellAmplitude(TMath::Abs(cluster->GetCellAbsId(i)));
153     if (e>maxe) {
154       maxe = e;
155       id   = cluster->GetCellAbsId(i);
156     }
157   }
158   return maxe;
159 }
160 //_____________________________________________________________________
161 Double_t AliAnalysisTaskPi0V2::GetCrossEnergy(const AliVCluster *cluster, Short_t &idmax) const
162 {
163   // Calculate the energy of cross cells around the leading cell.
164
165   AliVCaloCells *cells = 0;
166   if (fESD)
167     cells = fESD->GetEMCALCells();
168   if (!cells)
169     return 0;
170
171   AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
172   if (!geom)
173     return 0;
174
175   Int_t iSupMod = -1;
176   Int_t iTower  = -1;
177   Int_t iIphi   = -1;
178   Int_t iIeta   = -1;
179   Int_t iphi    = -1;
180   Int_t ieta    = -1;
181   Int_t iphis   = -1;
182   Int_t ietas   = -1;
183
184   Double_t crossEnergy = 0.;
185
186   geom->GetCellIndex(idmax,iSupMod,iTower,iIphi,iIeta);
187   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphis,ietas);
188
189   Int_t ncells = cluster->GetNCells();
190   for (Int_t i=0; i<ncells; i++) {
191     Int_t cellAbsId = cluster->GetCellAbsId(i);
192     geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
193     geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
194     Int_t aphidiff = TMath::Abs(iphi-iphis);
195     if (aphidiff>1)
196       continue;
197     Int_t aetadiff = TMath::Abs(ieta-ietas);
198     if (aetadiff>1)
199       continue;
200     if ( (aphidiff==1 && aetadiff==0) ||
201         (aphidiff==0 && aetadiff==1) ) {
202       crossEnergy += cells->GetCellAmplitude(cellAbsId);
203     }
204   }
205
206   return crossEnergy;
207 }
208 //_____________________________________________________________________
209 Bool_t AliAnalysisTaskPi0V2::IsWithinFiducialVolume(Short_t id) const
210 {
211   // Check if cell is within given fiducial volume.
212
213   Double_t fNFiducial = 1;
214
215   Int_t iSupMod = -1;
216   Int_t iTower  = -1;
217   Int_t iIphi   = -1;
218   Int_t iIeta   = -1;
219   Int_t iphi    = -1;
220   Int_t ieta    = -1;
221
222   Bool_t okrow = kFALSE;
223   Bool_t okcol = kFALSE;
224
225   AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
226   if (!geom)
227     return kFALSE;
228
229   Int_t cellAbsId = id;
230   geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
231   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
232
233   // Check rows/phi
234   if (iSupMod < 10) {
235     if (iphi >= fNFiducial && iphi < 24-fNFiducial)
236       okrow = kTRUE;
237   } else {
238     if (iphi >= fNFiducial && iphi < 12-fNFiducial)
239       okrow = kTRUE;
240   }
241   // Check columns/eta
242   Bool_t noEMCALBorderAtEta0 = kTRUE;
243   if (!noEMCALBorderAtEta0) {
244     if (ieta > fNFiducial && ieta < 48-fNFiducial)
245       okcol = kTRUE;
246   } else {
247     if (iSupMod%2==0) {
248       if (ieta >= fNFiducial)
249         okcol = kTRUE;
250     } else {
251       if (ieta < 48-fNFiducial)
252         okcol = kTRUE;
253     }
254   }
255   if (okrow && okcol)
256      return kTRUE;
257
258   return kFALSE;
259 }
260 //______________________________________________________________________
261 Bool_t AliAnalysisTaskPi0V2::IsGoodCluster(const AliESDCaloCluster *c) const
262 {
263
264   if(!c)
265     return kFALSE;
266
267   if(c->GetNCells() < fNcellCut)
268    return kFALSE;
269
270   if(c->E() < fECut)
271    return kFALSE;
272
273   Short_t id = -1;
274   Double_t maxE = GetMaxCellEnergy(c, id); 
275      if((1. - double(GetCrossEnergy(c,id))/maxE) > 0.97)
276     return kFALSE;
277
278
279   Float_t pos1[3];
280   c->GetPosition(pos1);
281   TVector3 clsPos(pos1);
282   Double_t eta = clsPos.Eta();
283
284   if(TMath::Abs(eta) > fEtaCut)
285     return kFALSE;  
286
287   if (!IsWithinFiducialVolume(id))
288     return kFALSE;
289
290   if(c->GetM02() >fM02Cut)
291     return kFALSE;
292
293 //  if(c->M20 >)
294
295   return kTRUE;
296
297 }
298 //_____________________________________________________________________
299 Bool_t AliAnalysisTaskPi0V2::IsGoodPion(const TLorentzVector &p1, const TLorentzVector &p2) const
300 {
301   // Is good pion?
302
303   if(fPi0AsyCut){
304     Double_t asym = TMath::Abs(p1.E()-p2.E())/(p1.E()+p2.E());
305     if (asym>0.7)
306       return kFALSE;
307   }
308   TLorentzVector pion;
309   pion = p1 + p2;
310   Double_t eta = pion.Eta();
311   if(TMath::Abs(eta) > fEtaCut)
312     return kFALSE;
313
314   return kTRUE;
315 }
316 //_______________________________________________________________________
317 void AliAnalysisTaskPi0V2::FillPion(const TLorentzVector& p1, const TLorentzVector& p2, Double_t EPV0r, Double_t EPV0A, Double_t EPV0C, Double_t EPTPC)
318 {
319   // Fill histogram.
320
321   if (!IsGoodPion(p1,p2))
322     return;
323   TLorentzVector pion;
324   pion = p1 + p2;
325
326   Double_t mass = pion.M();
327   Double_t pt   = pion.Pt();
328   Double_t phi  = pion.Phi();
329
330   Double_t dphiV0   = phi-EPV0r;
331   Double_t dphiV0A  = phi-EPV0A;
332   Double_t dphiV0C  = phi-EPV0C;
333   Double_t dphiTPC  = phi-EPTPC;
334
335   Double_t cos2phiV0  = TMath::Cos(2.*(dphiV0));
336   Double_t cos2phiV0A = TMath::Cos(2.*(dphiV0A));
337   Double_t cos2phiV0C = TMath::Cos(2.*(dphiV0C));
338   Double_t cos2phiTPC = TMath::Cos(2.*(dphiTPC));
339
340   dphiV0  = TVector2::Phi_0_2pi(dphiV0);  if(dphiV0  >TMath::Pi())  dphiV0 -= TMath::Pi();
341   dphiV0A = TVector2::Phi_0_2pi(dphiV0A); if(dphiV0A >TMath::Pi())  dphiV0A -= TMath::Pi();
342   dphiV0C = TVector2::Phi_0_2pi(dphiV0C); if(dphiV0C >TMath::Pi())  dphiV0C -= TMath::Pi();
343   dphiTPC = TVector2::Phi_0_2pi(dphiTPC); if(dphiTPC >TMath::Pi())  dphiTPC -= TMath::Pi();
344
345   Double_t xV0[5]; // Match ndims in fH  V0 EP
346   xV0[0]       = mass;
347   xV0[1]       = pt;
348   xV0[2]       = fCentrality;
349   xV0[3]       = dphiV0;
350   xV0[4]       = cos2phiV0;
351   fHEPV0r->Fill(xV0);
352
353   Double_t xV0A[5]; // Match ndims in fH V0A EP
354   xV0A[0]       = mass;
355   xV0A[1]       = pt;
356   xV0A[2]       = fCentrality;
357   xV0A[3]       = dphiV0A;
358   xV0A[4]       = cos2phiV0A;
359   fHEPV0A->Fill(xV0A);
360
361   Double_t xV0C[5]; // Match ndims in fH V0C EP
362   xV0C[0]       = mass;
363   xV0C[1]       = pt;
364   xV0C[2]       = fCentrality;
365   xV0C[3]       = dphiV0C;
366   xV0C[4]       = cos2phiV0C;
367   fHEPV0C->Fill(xV0C);
368
369   Double_t xTPC[5]; // Match ndims in fH TPC EP
370   xTPC[0]       = mass;
371   xTPC[1]       = pt;
372   xTPC[2]       = fCentrality;
373   xTPC[3]       = dphiTPC;
374   xTPC[4]       = cos2phiTPC;
375   fHEPTPC->Fill(xTPC);
376
377
378 }
379 void AliAnalysisTaskPi0V2::FillCluster(const TLorentzVector& p1, Double_t EPV0r, Double_t EPV0A, Double_t EPV0C, Double_t EPTPC)
380 {
381   //cluster(photon) v2 method
382   Double_t Pt   = p1.Pt();
383   Double_t Phi  = p1.Phi();
384
385   Double_t difClusV0 = TVector2::Phi_0_2pi(Phi-EPV0r);   if(difClusV0 >TMath::Pi()) difClusV0  -= TMath::Pi();
386   Double_t difClusV0A = TVector2::Phi_0_2pi(Phi-EPV0A);  if(difClusV0A >TMath::Pi()) difClusV0A -= TMath::Pi();
387   Double_t difClusV0C = TVector2::Phi_0_2pi(Phi-EPV0C);  if(difClusV0C >TMath::Pi()) difClusV0C -= TMath::Pi();
388   Double_t difClusTPC = TVector2::Phi_0_2pi(Phi-EPTPC);  if(difClusTPC >TMath::Pi()) difClusTPC -= TMath::Pi();
389
390   hdifClus_EPV0->Fill(fCentrality,  difClusV0, Pt);
391   hdifClus_EPV0A->Fill(fCentrality, difClusV0A, Pt);
392   hdifClus_EPV0C->Fill(fCentrality, difClusV0C, Pt);
393   hdifClus_EPTPC->Fill(fCentrality, difClusTPC, Pt);
394
395
396 }
397 //_________________________________________________________________________________________________
398 void AliAnalysisTaskPi0V2::GetMom(TLorentzVector& p, const AliESDCaloCluster *c, Double_t *vertex)
399 {
400   // Calculate momentum.
401   Float_t posMom[3];
402   c->GetPosition(posMom);
403   TVector3 clsPos2(posMom);
404
405   Double_t e   = c->E();
406   Double_t r   = clsPos2.Perp();
407   Double_t eta = clsPos2.Eta();
408   Double_t phi = clsPos2.Phi();
409
410   TVector3 pos;
411   pos.SetPtEtaPhi(r,eta,phi);
412
413   if (vertex) { //calculate direction relative to vertex
414     pos -= vertex;
415   }
416
417   Double_t rad = pos.Mag();
418   p.SetPxPyPzE(e*pos.x()/rad, e*pos.y()/rad, e*pos.z()/rad, e);
419
420 }
421 //________________________________________________________________________
422 void AliAnalysisTaskPi0V2::UserCreateOutputObjects()
423 {
424     // Create histograms
425     // Called once (on the worker node)
426         
427     fOutput = new TList();
428     fOutput->SetOwner();  // IMPORTANT!
429
430     hEvtCount = new TH1F("hEvtCount", " Event Plane", 10, 0.5, 10.5);
431     hEvtCount->GetXaxis()->SetBinLabel(1,"All");
432     hEvtCount->GetXaxis()->SetBinLabel(2,"Evt Cut");
433     hEvtCount->GetXaxis()->SetBinLabel(3,"Trg Class");
434     hEvtCount->GetXaxis()->SetBinLabel(4,"Vtx");
435     hEvtCount->GetXaxis()->SetBinLabel(5,"Cent");
436     hEvtCount->GetXaxis()->SetBinLabel(5,"EPtask");
437     hEvtCount->GetXaxis()->SetBinLabel(7,"EPvalue");
438     hEvtCount->GetXaxis()->SetBinLabel(8,"Pass");
439     fOutput->Add(hEvtCount);
440     
441     hEPTPC   = new TH2F("hEPTPC",   "EPTPC     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
442     hresoTPC = new TH2F("hresoTPC", "TPc reso  vs cent", 100, 0., 100., 100, 0., 1.);
443     hEPV0    = new TH2F("hEPV0",    "EPV0      vs cent", 100, 0., 100., 100, 0., TMath::Pi());
444     hEPV0A   = new TH2F("hEPV0A",   "EPV0A     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
445     hEPV0C   = new TH2F("hEPV0C",   "EPV0C     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
446     hEPV0Ar  = new TH2F("hEPV0Ar",  "EPV0Ar    vs cent", 100, 0., 100., 100, 0., TMath::Pi());
447     hEPV0Cr  = new TH2F("hEPV0Cr",  "EPV0Cr    vs cent", 100, 0., 100., 100, 0., TMath::Pi());
448     hEPV0r   = new TH2F("hEPV0r",   "EPV0r     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
449     hEPV0AR4 = new TH2F("hEPV0AR4", "EPV0AR4   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
450     hEPV0AR7 = new TH2F("hEPV0AR7", "EPV0AR7   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
451     hEPV0CR0 = new TH2F("hEPV0CR0", "EPV0CR0   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
452     hEPV0CR3 = new TH2F("hEPV0CR3", "EPV0CR3   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
453     fOutput->Add(hEPTPC);
454     fOutput->Add(hresoTPC);
455     fOutput->Add(hEPV0);
456     fOutput->Add(hEPV0A);
457     fOutput->Add(hEPV0C);
458     fOutput->Add(hEPV0Ar);
459     fOutput->Add(hEPV0Cr);
460     fOutput->Add(hEPV0r);
461     fOutput->Add(hEPV0AR4);
462     fOutput->Add(hEPV0AR7);
463     fOutput->Add(hEPV0CR0);
464     fOutput->Add(hEPV0CR3);
465
466     hdifV0A_V0CR0    = new TH2F("hdifV0A_V0CR0",    "EP A-R0 ",  100, 0., 100., 100, -1., 1.);    
467     hdifV0A_V0CR3    = new TH2F("hdifV0A_V0CR3",    "EP A-R3 ",  100, 0., 100., 100, -1., 1.);    
468     hdifV0ACR0_V0CR3 = new TH2F("hdifV0ACR0_V0CR3", "EP R0-R3 ", 100, 0., 100., 100, -1., 1.);    
469     hdifV0C_V0AR4    = new TH2F("hdifV0C_V0AR4",    "EP C-R4 ",  100, 0., 100., 100, -1., 1.);    
470     hdifV0C_V0AR7    = new TH2F("hdifV0C_V0AR7",    "EP C-R7 ",  100, 0., 100., 100, -1., 1.);    
471     hdifV0AR4_V0AR7  = new TH2F("hdifV0AR4_V0AR7",  "EP R4-R7 ", 100, 0., 100., 100, -1., 1.);    
472     fOutput->Add(hdifV0A_V0CR0);
473     fOutput->Add(hdifV0A_V0CR3);
474     fOutput->Add(hdifV0ACR0_V0CR3);
475     fOutput->Add(hdifV0C_V0AR4);
476     fOutput->Add(hdifV0C_V0AR7);
477     fOutput->Add(hdifV0AR4_V0AR7);
478
479     hdifV0A_V0C = new TH2F("hdifV0A_V0C", "EP A-C  ", 100, 0., 100., 100, -1., 1.);
480     hdifV0A_TPC = new TH2F("hdifV0A_TPC", "EP A-TPC", 100, 0., 100., 100, -1., 1.);
481     hdifV0C_TPC = new TH2F("hdifV0C_TPC", "EP C-TPC", 100, 0., 100., 100, -1., 1.);
482     hdifV0C_V0A = new TH2F("hdifV0C_V0A", "EP C-A  ", 100, 0., 100., 100, -1., 1.);
483     fOutput->Add(hdifV0A_V0C);
484     fOutput->Add(hdifV0A_TPC);
485     fOutput->Add(hdifV0C_TPC);
486     fOutput->Add(hdifV0C_V0A);
487
488
489
490     hdifEMC_EPV0  = new TH3F("hdifEMC_EPV0",  "dif phi in EMC with EP",  100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
491     hdifEMC_EPV0A = new TH3F("hdifEMC_EPV0A", "dif phi in EMC with EP",  100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
492     hdifEMC_EPV0C = new TH3F("hdifEMC_EPV0C", "dif phi in EMC with EP",  100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
493     fOutput->Add(hdifEMC_EPV0);
494     fOutput->Add(hdifEMC_EPV0A);
495     fOutput->Add(hdifEMC_EPV0C);
496
497     hdifful_EPV0 = new TH3F("hdifful_EPV0",    "dif phi in full with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
498     hdifful_EPV0A = new TH3F("hdifful_EPV0A",  "dif phi in full with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
499     hdifful_EPV0C = new TH3F("hdifful_EPV0C",  "dif phi in full with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
500     fOutput->Add(hdifful_EPV0);
501     fOutput->Add(hdifful_EPV0A);
502     fOutput->Add(hdifful_EPV0C);
503
504     hdifout_EPV0  = new TH3F("hdifout_EPV0",  "dif phi NOT in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
505     hdifout_EPV0A = new TH3F("hdifout_EPV0A", "dif phi NOT in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
506     hdifout_EPV0C = new TH3F("hdifout_EPV0C", "dif phi NOT in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
507     fOutput->Add(hdifout_EPV0);
508     fOutput->Add(hdifout_EPV0A);
509     fOutput->Add(hdifout_EPV0C);
510
511     hdifEMC_EPTPC = new TH3F("hdifEMC_EPTPC", "dif phi in EMC with EP",  100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
512     hdifful_EPTPC = new TH3F("hdifful_EPTPC",  "dif phi in full with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
513     hdifout_EPTPC = new TH3F("hdifout_EPTPC", "dif phi NOT in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
514     fOutput->Add(hdifEMC_EPTPC);
515     fOutput->Add(hdifful_EPTPC);
516     fOutput->Add(hdifout_EPTPC);
517
518     hdifClus_EPV0 = new TH3F("hdifClus_EPV0", "dif phi in EMC Clus with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
519     hdifClus_EPV0A = new TH3F("hdifClus_EPV0A", "dif phi in EMC Clus with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
520     hdifClus_EPV0C = new TH3F("hdifClus_EPV0C", "dif phi in EMC Clus with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
521     hdifClus_EPTPC = new TH3F("hdifClus_EPTPC", "dif phi in EMC Clus with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
522     fOutput->Add(hdifClus_EPV0);
523     fOutput->Add(hdifClus_EPV0A);
524     fOutput->Add(hdifClus_EPV0C);
525     fOutput->Add(hdifClus_EPTPC);
526
527     hAllcentV0  = new TH1F("hAllcentV0",  "All cent EP V0",  100, 0., TMath::Pi());
528     hAllcentV0r = new TH1F("hAllcentV0r", "All cent EP V0r", 100, 0., TMath::Pi());
529     hAllcentV0A = new TH1F("hAllcentV0A", "All cent EP V0A", 100, 0., TMath::Pi());
530     hAllcentV0C = new TH1F("hAllcentV0C", "All cent EP V0C", 100, 0., TMath::Pi());
531     hAllcentTPC = new TH1F("hAllcentTPC", "All cent EP TPC", 100, 0., TMath::Pi());
532     fOutput->Add(hAllcentV0);
533     fOutput->Add(hAllcentV0r);
534     fOutput->Add(hAllcentV0A);
535     fOutput->Add(hAllcentV0C);
536     fOutput->Add(hAllcentTPC);
537
538     h2DcosV0r = new TH2F("h2DcosV0r", "cos(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
539     h2DsinV0r = new TH2F("h2DsinV0r", "sin(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
540     h2DcosV0A = new TH2F("h2DcosV0A", "cos(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
541     h2DsinV0A = new TH2F("h2DsinV0A", "sin(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
542     h2DcosV0C = new TH2F("h2DcosV0C", "cos(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
543     h2DsinV0C = new TH2F("h2DsinV0C", "sin(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
544     h2DcosTPC = new TH2F("h2DcosTPC", "cos(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
545     h2DsinTPC = new TH2F("h2DsinTPC", "sin(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
546     fOutput->Add(h2DcosV0r);
547     fOutput->Add(h2DsinV0r);
548     fOutput->Add(h2DcosV0A);
549     fOutput->Add(h2DsinV0A);
550     fOutput->Add(h2DcosV0C);
551     fOutput->Add(h2DsinV0C);
552     fOutput->Add(h2DcosTPC);
553     fOutput->Add(h2DsinTPC);
554
555     h2DsinV0C = new TH2F("h2DsinV0C", "sin(Phi) V0r vs Run NUmber", 200, 0, 200, 100, -1, 1);
556
557     const Int_t ndims = 5;
558     Int_t nMgg=500, nPt=40, nCent=20, nDeltaPhi=315,  ncos2phi=500;
559     Int_t bins[ndims] = {nMgg, nPt, nCent, nDeltaPhi, ncos2phi};
560     Double_t xmin[ndims] = { 0,   0.,  0,   0.,     -1.};
561     Double_t xmax[ndims] = { 0.5, 20., 100, 3.15,   1.};
562     fHEPV0r  = new THnSparseF("fHEPV0r",  "Flow histogram EPV0",  ndims, bins, xmin, xmax);
563     fHEPV0A = new THnSparseF("fHEPV0A",   "Flow histogram EPV0A", ndims, bins, xmin, xmax);
564     fHEPV0C = new THnSparseF("fHEPV0C",   "Flow histogram EPV0C", ndims, bins, xmin, xmax);
565     fHEPTPC = new THnSparseF("fHEPTPC",   "Flow histogram EPTPC", ndims, bins, xmin, xmax);
566     fOutput->Add(fHEPV0r);
567     fOutput->Add(fHEPV0A);
568     fOutput->Add(fHEPV0C);
569     fOutput->Add(fHEPTPC);
570
571     
572
573     PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram
574 }
575
576 //________________________________________________________________________
577 void AliAnalysisTaskPi0V2::UserExec(Option_t *) 
578 {
579     // Main loop
580     // Called for each event
581         
582     // Create pointer to reconstructed event
583    AliVEvent *event = InputEvent();
584    if (!event) { Printf("ERROR: Could not retrieve event"); return; }
585    // create pointer to event
586    fESD = dynamic_cast<AliESDEvent*>(event);
587    if (!fESD) {
588      AliError("Cannot get the ESD event");
589      return;
590    }
591    hEvtCount->Fill(1);
592     
593   Int_t AbsRunNumber = fESD->GetRunNumber();
594   fRunNumber = ConvertToInternalRunNumber(AbsRunNumber);
595
596   Bool_t isSelected =0;      
597   if(fEvtSelect == 1){  //MB+SemiCentral
598     isSelected = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & (AliVEvent::kMB | AliVEvent::kSemiCentral));
599   } else if (fEvtSelect == 2){  //MB+Central+SemiCentral
600     isSelected = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & (AliVEvent::kMB | AliVEvent::kSemiCentral | AliVEvent::kCentral));
601   } else if(fEvtSelect == 3){  //MB
602  isSelected = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & (AliVEvent::kMB ));
603   }
604   if(!isSelected )
605         return; 
606
607   hEvtCount->Fill(2);
608   if(!fTrigClass.IsNull()){
609     TString fired;
610     fired = fESD->GetFiredTriggerClasses();
611     if (!fired.Contains("-B-"))
612       return;
613     TObjArray *arr = fTrigClass.Tokenize("|");
614     if (!arr)
615       return;
616     Bool_t match = 0;
617     for (Int_t i=0;i<arr->GetEntriesFast();++i) {
618       TObject *obj = arr->At(i);
619       if (!obj)
620         continue;
621       if (fired.Contains(obj->GetName())) {
622         match = 1;
623         break;
624       }
625     }
626     delete arr;
627     if (
628         !match || //select by Trigger classes in KCentral and KSemiCentral
629         !(((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & (AliVEvent::kMB ))              // always accept MB
630         ) 
631       return; //Not match skip this event
632   }
633
634     hEvtCount->Fill(3);
635     const AliESDVertex* fvertex = fESD->GetPrimaryVertex();
636     if(TMath::Abs(fvertex->GetZ())>fVtxCut)
637       return;
638     Double_t vertex[3] = {fvertex->GetX(), fvertex->GetY(), fvertex->GetZ()};
639
640     hEvtCount->Fill(4);
641
642     if(fESD->GetCentrality()) {
643       fCentrality = 
644         fESD->GetCentrality()->GetCentralityPercentile("CL1"); //spd vertex
645     } else{
646            return;
647     }
648
649     hEvtCount->Fill(5);
650     AliEventplane *ep = fESD->GetEventplane();
651     if (ep) {
652       if (ep->GetQVector())
653         fEPTPC    = ep->GetQVector()->Phi()/2. ;
654       else
655         fEPTPC = -999.;
656       if (ep->GetQsub1()&&ep->GetQsub2())
657         fEPTPCreso  = TMath::Cos(2.*(ep->GetQsub1()->Phi()/2.-ep->GetQsub2()->Phi()/2.));
658       else
659         fEPTPCreso = -1;
660
661       fEPV0    = ep->GetEventplane("V0",  fESD);
662       fEPV0A   = ep->GetEventplane("V0A", fESD);
663       fEPV0C   = ep->GetEventplane("V0C", fESD);
664       Double_t qx=0, qy=0;
665       Double_t qxr=0, qyr=0;
666       fEPV0Ar  = ep->CalculateVZEROEventPlane(fESD, 4, 5, 2, qxr, qyr);
667       fEPV0Cr  = ep->CalculateVZEROEventPlane(fESD, 2, 3, 2, qx,  qy);
668       qxr += qx;
669       qyr += qy;
670       fEPV0r   = TMath::ATan2(qyr,qxr)/2.;
671       fEPV0AR4 = ep->CalculateVZEROEventPlane(fESD, 4, 2, qx, qy);
672       fEPV0AR5 = ep->CalculateVZEROEventPlane(fESD, 5, 2, qx, qy);
673       fEPV0AR6 = ep->CalculateVZEROEventPlane(fESD, 6, 2, qx, qy);
674       fEPV0AR7 = ep->CalculateVZEROEventPlane(fESD, 7, 2, qx, qy);
675       fEPV0CR0 = ep->CalculateVZEROEventPlane(fESD, 0, 2, qx, qy);
676       fEPV0CR1 = ep->CalculateVZEROEventPlane(fESD, 1, 2, qx, qy);
677       fEPV0CR2 = ep->CalculateVZEROEventPlane(fESD, 2, 2, qx, qy);
678       fEPV0CR3 = ep->CalculateVZEROEventPlane(fESD, 3, 2, qx, qy);
679     }
680     FillEPQA(); //Fill the EP QA
681
682     hEvtCount->Fill(6);
683
684     if( fEPV0A<-2. || fEPV0C<-2. || fEPTPC<-2. || fEPV0r<-2.) 
685       return;
686
687     hEvtCount->Fill(7);
688
689     fEPV0   = TVector2::Phi_0_2pi(fEPV0);    if(fEPV0>TMath::Pi())   fEPV0  = fEPV0  - TMath::Pi();
690     fEPV0r  = TVector2::Phi_0_2pi(fEPV0r);   if(fEPV0r>TMath::Pi())  fEPV0r = fEPV0r - TMath::Pi();
691     fEPV0A  = TVector2::Phi_0_2pi(fEPV0A);   if(fEPV0A>TMath::Pi())  fEPV0A = fEPV0A - TMath::Pi();
692     fEPV0C  = TVector2::Phi_0_2pi(fEPV0C);   if(fEPV0C>TMath::Pi())  fEPV0C = fEPV0C - TMath::Pi();
693     fEPV0Ar = TVector2::Phi_0_2pi(fEPV0Ar);  if(fEPV0Ar>TMath::Pi()) fEPV0Ar = fEPV0Ar - TMath::Pi();
694     fEPV0Cr = TVector2::Phi_0_2pi(fEPV0Cr);  if(fEPV0Cr>TMath::Pi()) fEPV0Cr = fEPV0Cr - TMath::Pi();
695     fEPV0AR4   = TVector2::Phi_0_2pi(fEPV0AR4);    if(fEPV0AR4>TMath::Pi())   fEPV0AR4  = fEPV0AR4 - TMath::Pi();
696     fEPV0AR7   = TVector2::Phi_0_2pi(fEPV0AR7);    if(fEPV0AR7>TMath::Pi())   fEPV0AR7  = fEPV0AR7 - TMath::Pi();
697     fEPV0CR0   = TVector2::Phi_0_2pi(fEPV0CR0);    if(fEPV0CR0>TMath::Pi())   fEPV0CR0  = fEPV0CR0 - TMath::Pi();
698     fEPV0CR3   = TVector2::Phi_0_2pi(fEPV0CR3);    if(fEPV0CR3>TMath::Pi())   fEPV0CR3  = fEPV0CR3 - TMath::Pi();
699
700    if(fEPTPC != -999.)
701    hEPTPC->Fill(fCentrality,  fEPTPC); 
702    if(fEPTPCreso!=-1) hresoTPC->Fill(fCentrality, fEPTPCreso);
703    hEPV0->Fill(fCentrality,   fEPV0);
704    hEPV0A->Fill(fCentrality,  fEPV0A);
705    hEPV0C->Fill(fCentrality,  fEPV0C);
706    hEPV0Ar->Fill(fCentrality, fEPV0Ar);
707    hEPV0Cr->Fill(fCentrality, fEPV0Cr);
708    hEPV0r->Fill(fCentrality,  fEPV0r);
709    hEPV0AR4->Fill(fCentrality, fEPV0AR4);
710    hEPV0AR7->Fill(fCentrality, fEPV0AR7);
711    hEPV0CR0->Fill(fCentrality, fEPV0CR0);
712    hEPV0CR3->Fill(fCentrality, fEPV0CR3);
713
714    hAllcentV0->Fill(fEPV0);
715    hAllcentV0r->Fill(fEPV0r);
716    hAllcentV0A->Fill(fEPV0A);
717    hAllcentV0C->Fill(fEPV0C);  
718    hAllcentTPC->Fill(fEPTPC);
719
720    hdifV0A_V0CR0->Fill(fCentrality, TMath::Cos(2.*(fEPV0A - fEPV0CR0)));
721    hdifV0A_V0CR3->Fill(fCentrality, TMath::Cos(2.*(fEPV0A - fEPV0CR3)));
722    hdifV0ACR0_V0CR3->Fill(fCentrality, TMath::Cos(2*(fEPV0CR0 - fEPV0CR3)));
723    hdifV0C_V0AR4->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0AR4)));
724    hdifV0C_V0AR7->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0AR7)));
725    hdifV0AR4_V0AR7->Fill(fCentrality, TMath::Cos(2*(fEPV0AR4 - fEPV0AR7)));
726         
727    hdifV0A_V0C->Fill(fCentrality, TMath::Cos(2*(fEPV0A - fEPV0C)));
728    hdifV0A_TPC->Fill(fCentrality, TMath::Cos(2*(fEPV0A - fEPTPC)));
729    hdifV0C_TPC->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPTPC)));
730    hdifV0C_V0A->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0A)));
731     // Cluster loop for reconstructed event
732
733    Int_t nCluster =  fESD->GetNumberOfCaloClusters(); 
734    for(Int_t i=0; i<nCluster; ++i){
735      AliESDCaloCluster *c1 = fESD->GetCaloCluster(i);
736      if(!c1->IsEMCAL()) continue;
737      if(!IsGoodCluster(c1)) continue;
738      TLorentzVector p1;
739      GetMom(p1, c1, vertex);
740      FillCluster(p1, fEPV0r, fEPV0A, fEPV0C, fEPTPC);
741      for(Int_t j=i+1; j<nCluster; ++j){
742        AliESDCaloCluster *c2 = fESD->GetCaloCluster(j);
743        if(!c2->IsEMCAL()) continue;
744        if(!IsGoodCluster(c2)) continue;
745        TLorentzVector p2;
746        GetMom(p2, c2, vertex);
747        FillPion(p1, p2, fEPV0r, fEPV0A, fEPV0C, fEPTPC);
748      } 
749    }
750
751    if (!fTracksName.IsNull() && !fTracks) {
752      fTracks = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTracksName));
753      if (!fTracks) {
754        AliError(Form("%s: Could not retrieve tracks %s!", GetName(), fTracksName.Data())); 
755        return;
756      }
757   }
758
759    Int_t ntracks = fTracks->GetEntries();
760    for(Int_t i=0; i<ntracks; ++i){
761      AliVTrack *track = static_cast<AliVTrack*>(fTracks->At(i));
762      Double_t tPhi = track->Phi();
763      Double_t tPt  = track->Pt();
764
765      Double_t difTrackV0  = TVector2::Phi_0_2pi(tPhi-fEPV0);   if(difTrackV0  >TMath::Pi()) difTrackV0  -= TMath::Pi();
766      Double_t difTrackV0A = TVector2::Phi_0_2pi(tPhi-fEPV0A);  if(difTrackV0A >TMath::Pi()) difTrackV0A -= TMath::Pi();
767      Double_t difTrackV0C = TVector2::Phi_0_2pi(tPhi-fEPV0C);  if(difTrackV0C >TMath::Pi()) difTrackV0C -= TMath::Pi();
768      Double_t difTrackTPC = TVector2::Phi_0_2pi(tPhi-fEPTPC);  if(difTrackTPC >TMath::Pi()) difTrackTPC -= TMath::Pi();
769      if(track->IsEMCAL()){      
770        hdifEMC_EPV0->Fill(fCentrality, difTrackV0, tPt);
771        hdifEMC_EPV0A->Fill(fCentrality, difTrackV0A, tPt);
772        hdifEMC_EPV0C->Fill(fCentrality, difTrackV0C, tPt);
773        hdifEMC_EPTPC->Fill(fCentrality, difTrackTPC, tPt);
774      }else{
775        hdifout_EPV0->Fill(fCentrality, difTrackV0, tPt);
776        hdifout_EPV0A->Fill(fCentrality, difTrackV0A, tPt);
777        hdifout_EPV0C->Fill(fCentrality, difTrackV0C, tPt);
778        hdifout_EPTPC->Fill(fCentrality, difTrackTPC, tPt);
779      }
780      hdifful_EPV0->Fill(fCentrality,    difTrackV0, tPt);
781      hdifful_EPV0A->Fill(fCentrality,   difTrackV0A, tPt);
782      hdifful_EPV0C->Fill(fCentrality,   difTrackV0C, tPt);
783      hdifful_EPTPC->Fill(fCentrality,   difTrackTPC, tPt);
784     }
785     hEvtCount->Fill(8);
786
787     // NEW HISTO should be filled before this point, as PostData puts the
788     // information for this iteration of the UserExec in the container
789     PostData(1, fOutput);
790 }
791 //____________________________________________________________________
792 Int_t AliAnalysisTaskPi0V2::ConvertToInternalRunNumber(Int_t n)
793 {
794     switch(n){
795     case  170593 : return 179 ;
796     case  170572 : return 178 ;
797     case  170556 : return 177 ;
798     case  170552 : return 176 ;
799     case  170546 : return 175 ;
800     case  170390 : return 174 ;
801     case  170389 : return 173 ;
802     case  170388 : return 172 ;
803     case  170387 : return 171 ;
804     case  170315 : return 170 ;
805     case  170313 : return 169 ;
806     case  170312 : return 168 ;
807     case  170311 : return 167 ;
808     case  170309 : return 166 ;
809     case  170308 : return 165 ;
810     case  170306 : return 164 ;
811     case  170270 : return 163 ;
812     case  170269 : return 162 ;
813     case  170268 : return 161 ;
814     case  170267 : return 160 ;
815     case  170264 : return 159 ;
816     case  170230 : return 158 ;
817     case  170228 : return 157 ;
818     case  170208 : return 156 ;
819     case  170207 : return 155 ;
820     case  170205 : return 154 ;
821     case  170204 : return 153 ;
822     case  170203 : return 152 ;
823     case  170195 : return 151 ;
824     case  170193 : return 150 ;
825     case  170163 : return 149 ;
826     case  170162 : return 148 ;
827     case  170159 : return 147 ;
828     case  170155 : return 146 ;
829     case  170152 : return 145 ;
830     case  170091 : return 144 ;
831     case  170089 : return 143 ;
832     case  170088 : return 142 ;
833     case  170085 : return 141 ;
834     case  170084 : return 140 ;
835     case  170083 : return 139 ;
836     case  170081 : return 138 ;
837     case  170040 : return 137 ;
838     case  170038 : return 136 ;
839     case  170036 : return 135 ;
840     case  170027 : return 134 ;
841     case  169981 : return 133 ;
842     case  169975 : return 132 ;
843     case  169969 : return 131 ;
844     case  169965 : return 130 ;
845     case  169961 : return 129 ;
846     case  169956 : return 128 ;
847     case  169926 : return 127 ;
848     case  169924 : return 126 ;
849     case  169923 : return 125 ;
850     case  169922 : return 124 ;
851     case  169919 : return 123 ;
852     case  169918 : return 122 ;
853     case  169914 : return 121 ;
854     case  169859 : return 120 ;
855     case  169858 : return 119 ;
856     case  169855 : return 118 ;
857     case  169846 : return 117 ;
858     case  169838 : return 116 ;
859     case  169837 : return 115 ;
860     case  169835 : return 114 ;
861     case  169683 : return 113 ;
862     case  169628 : return 112 ;
863     case  169591 : return 111 ;
864     case  169590 : return 110 ;
865     case  169588 : return 109 ;
866     case  169587 : return 108 ;
867     case  169586 : return 107 ;
868     case  169584 : return 106 ;
869     case  169557 : return 105 ;
870     case  169555 : return 104 ;
871     case  169554 : return 103 ;
872     case  169553 : return 102 ;
873     case  169550 : return 101 ;
874     case  169515 : return 100 ;
875     case  169512 : return 99 ;
876     case  169506 : return 98 ;
877     case  169504 : return 97 ;
878     case  169498 : return 96 ;
879     case  169475 : return 95 ;
880     case  169420 : return 94 ;
881     case  169419 : return 93 ;
882     case  169418 : return 92 ;
883     case  169417 : return 91 ;
884     case  169415 : return 90 ;
885     case  169411 : return 89 ;
886     case  169238 : return 88 ;
887     case  169236 : return 87 ;
888     case  169167 : return 86 ;
889     case  169160 : return 85 ;
890     case  169156 : return 84 ;
891     case  169148 : return 83 ;
892     case  169145 : return 82 ;
893     case  169144 : return 81 ;
894     case  169143 : return 80 ;
895     case  169138 : return 79 ;
896     case  169099 : return 78 ;
897     case  169094 : return 77 ;
898     case  169091 : return 76 ;
899     case  169045 : return 75 ;
900     case  169044 : return 74 ;
901     case  169040 : return 73 ;
902     case  169035 : return 72 ;
903     case  168992 : return 71 ;
904     case  168988 : return 70 ;
905     case  168984 : return 69 ;
906     case  168826 : return 68 ;
907     case  168777 : return 67 ;
908     case  168514 : return 66 ;
909     case  168512 : return 65 ;
910     case  168511 : return 64 ;
911     case  168467 : return 63 ;
912     case  168464 : return 62 ;
913     case  168461 : return 61 ;
914     case  168460 : return 60 ;
915     case  168458 : return 59 ;
916     case  168362 : return 58 ;
917     case  168361 : return 57 ;
918     case  168356 : return 56 ;
919     case  168342 : return 55 ;
920     case  168341 : return 54 ;
921     case  168325 : return 53 ;
922     case  168322 : return 52 ;
923     case  168318 : return 51 ;
924     case  168311 : return 50 ;
925     case  168310 : return 49 ;
926     case  168213 : return 48 ;
927     case  168212 : return 47 ;
928     case  168208 : return 46 ;
929     case  168207 : return 45 ;
930     case  168206 : return 44 ;
931     case  168205 : return 43 ;
932     case  168204 : return 42 ;
933     case  168203 : return 41 ;
934     case  168181 : return 40 ;
935     case  168177 : return 39 ;
936     case  168175 : return 38 ;
937     case  168173 : return 37 ;
938     case  168172 : return 36 ;
939     case  168171 : return 35 ;
940     case  168115 : return 34 ;
941     case  168108 : return 33 ;
942     case  168107 : return 32 ;
943     case  168105 : return 31 ;
944     case  168104 : return 30 ;
945     case  168103 : return 29 ;
946     case  168076 : return 28 ;
947     case  168069 : return 27 ;
948     case  168068 : return 26 ;
949     case  168066 : return 25 ;
950     case  167988 : return 24 ;
951     case  167987 : return 23 ;
952     case  167986 : return 22 ;
953     case  167985 : return 21 ;
954     case  167921 : return 20 ;
955     case  167920 : return 19 ;
956     case  167915 : return 18 ;
957     case  167909 : return 17 ;
958     case  167903 : return 16 ;
959     case  167902 : return 15 ;
960     case  167818 : return 14 ;
961     case  167814 : return 13 ;
962     case  167813 : return 12 ;
963     case  167808 : return 11 ;
964     case  167807 : return 10 ;
965     case  167806 : return 9 ;
966     case  167713 : return 8 ;
967     case  167712 : return 7 ;
968     case  167711 : return 6 ;
969     case  167706 : return 5 ;
970     case  167693 : return 4 ;
971     case  166532 : return 3 ;
972     case  166530 : return 2 ;
973     case  166529 : return 1 ;
974
975     default : return 199;
976     }
977 }
978 //_______________________________________________________________________
979 void AliAnalysisTaskPi0V2::FillEPQA()
980 {
981   
982   h2DcosV0r->Fill(fRunNumber, TMath::Cos(fEPV0r));
983   h2DsinV0r->Fill(fRunNumber, TMath::Sin(fEPV0r));
984   h2DcosV0A->Fill(fRunNumber, TMath::Cos(fEPV0A));
985   h2DsinV0A->Fill(fRunNumber, TMath::Sin(fEPV0A));
986   h2DcosV0C->Fill(fRunNumber, TMath::Cos(fEPV0C));
987   h2DsinV0C->Fill(fRunNumber, TMath::Sin(fEPV0C));
988   h2DcosTPC->Fill(fRunNumber, TMath::Cos(fEPTPC));
989   h2DsinTPC->Fill(fRunNumber, TMath::Sin(fEPTPC));
990
991
992 }
993 //________________________________________________________________________
994 void AliAnalysisTaskPi0V2::Terminate(Option_t *) 
995 {
996     // Draw result to screen, or perform fitting, normalizations
997     // Called once at the end of the query
998 //    fOutput = dynamic_cast<TList*> (GetOutputData(1));
999  //   if(!fOutput) { Printf("ERROR: could not retrieve TList fOutput"); return; }
1000         
1001     // Get the physics selection histograms with the selection statistics
1002     //AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1003     //AliESDInputHandler *inputH = dynamic_cast<AliESDInputHandler*>(mgr->GetInputEventHandler());
1004     //TH2F *histStat = (TH2F*)inputH->GetStatistics();
1005    
1006    
1007     // NEW HISTO should be retrieved from the TList container in the above way,
1008     // so it is available to draw on a canvas such as below
1009
1010 }