]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/EMCALTasks/AliAnalysisTaskPi0V2.cxx
0a3053679cfcb97c610a45bc2af2a66e7d62f573
[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
52 using std::cout;
53 using std::endl;
54
55 ClassImp(AliAnalysisTaskPi0V2)
56
57 //________________________________________________________________________
58 AliAnalysisTaskPi0V2::AliAnalysisTaskPi0V2() // All data members should be initialised here
59  //  :AliAnalysisTaskSE(),
60    :AliAnalysisTaskSE(),
61     fOutput(0),
62     fESD(0),
63     fcheckEP2sub(1),
64     fCentrality(99.),
65     fEPTPC(-999.),
66     fEPTPCreso(0.), 
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     hEvtCount(0), 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)
76
77 {
78     // Dummy constructor ALWAYS needed for I/O.
79     DefineInput(0, TChain::Class());
80     DefineOutput(1, TList::Class());                                            // for output list
81 }
82
83 //________________________________________________________________________
84 AliAnalysisTaskPi0V2::AliAnalysisTaskPi0V2(const char *name) // All data members should be initialised here
85    :AliAnalysisTaskSE(name),
86     fOutput(0),
87     fESD(0),
88     fcheckEP2sub(1),
89     fCentrality(99.),
90     fEPTPC(-999.),
91     fEPTPCreso(0.),
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     hEvtCount(0), 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)
101 {
102     // Constructor
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
108 }
109
110 //________________________________________________________________________
111 AliAnalysisTaskPi0V2::~AliAnalysisTaskPi0V2()
112 {
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()) {
116         delete fOutput;
117     }
118 }
119 //_____________________________________________________________________
120 Double_t AliAnalysisTaskPi0V2::GetMaxCellEnergy(const AliVCluster *cluster, Short_t &id) const
121 {
122   // Get maximum energy of attached cell.
123
124   id = -1;
125
126   AliVCaloCells *cells = 0;
127   if (fESD)
128     cells = fESD->GetEMCALCells();
129 //  else
130 //    cells = fAOD->GetEMCALCells();
131   if (!cells)
132     return 0;
133
134   Double_t maxe = 0;
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)));
138     if (e>maxe) {
139       maxe = e;
140       id   = cluster->GetCellAbsId(i);
141     }
142   }
143   return maxe;
144 }
145 //_____________________________________________________________________
146 Double_t AliAnalysisTaskPi0V2::GetCrossEnergy(const AliVCluster *cluster, Short_t &idmax) const
147 {
148   // Calculate the energy of cross cells around the leading cell.
149
150   AliVCaloCells *cells = 0;
151   if (fESD)
152     cells = fESD->GetEMCALCells();
153 //  else
154 //    cells = fAOD->GetEMCALCells();
155   if (!cells)
156     return 0;
157
158   AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
159   if (!geom)
160     return 0;
161
162   Int_t iSupMod = -1;
163   Int_t iTower  = -1;
164   Int_t iIphi   = -1;
165   Int_t iIeta   = -1;
166   Int_t iphi    = -1;
167   Int_t ieta    = -1;
168   Int_t iphis   = -1;
169   Int_t ietas   = -1;
170
171   Double_t crossEnergy = 0;
172
173   geom->GetCellIndex(idmax,iSupMod,iTower,iIphi,iIeta);
174   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphis,ietas);
175
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);
182     if (aphidiff>1)
183       continue;
184     Int_t aetadiff = TMath::Abs(ieta-ietas);
185     if (aetadiff>1)
186       continue;
187     if ( (aphidiff==1 && aetadiff==0) ||
188         (aphidiff==0 && aetadiff==1) ) {
189       crossEnergy += cells->GetCellAmplitude(cellAbsId);
190     }
191   }
192
193   return crossEnergy;
194 }
195 //_____________________________________________________________________
196 Bool_t AliAnalysisTaskPi0V2::IsWithinFiducialVolume(Short_t id) const
197 {
198   // Check if cell is within given fiducial volume.
199
200   Double_t fNFiducial = 1;
201
202   Int_t iSupMod = -1;
203   Int_t iTower  = -1;
204   Int_t iIphi   = -1;
205   Int_t iIeta   = -1;
206   Int_t iphi    = -1;
207   Int_t ieta    = -1;
208
209   Bool_t okrow = kFALSE;
210   Bool_t okcol = kFALSE;
211
212   AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
213   if (!geom)
214     return kFALSE;
215
216   Int_t cellAbsId = id;
217   geom->GetCellIndex(cellAbsId,iSupMod,iTower,iIphi,iIeta);
218   geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,iIphi, iIeta,iphi,ieta);
219
220   // Check rows/phi
221   if (iSupMod < 10) {
222     if (iphi >= fNFiducial && iphi < 24-fNFiducial)
223       okrow = kTRUE;
224   } else {
225     if (iphi >= fNFiducial && iphi < 12-fNFiducial)
226       okrow = kTRUE;
227   }
228   // Check columns/eta
229   Bool_t noEMCALBorderAtEta0 = kTRUE;
230   if (!noEMCALBorderAtEta0) {
231     if (ieta > fNFiducial && ieta < 48-fNFiducial)
232       okcol = kTRUE;
233   } else {
234     if (iSupMod%2==0) {
235       if (ieta >= fNFiducial)
236         okcol = kTRUE;
237     } else {
238       if (ieta < 48-fNFiducial)
239         okcol = kTRUE;
240     }
241   }
242   if (okrow && okcol)
243      return kTRUE;
244
245   return kFALSE;
246 }
247 //______________________________________________________________________
248 Bool_t AliAnalysisTaskPi0V2::IsGoodCluster(const AliESDCaloCluster *c) const
249 {
250
251   if(!c)
252     return kFALSE;
253
254   if(c->GetNCells() < 2)
255    return kFALSE;
256
257   if(c->E() < 1.)
258    return kFALSE;
259
260   Short_t id = -1;
261   Double_t maxE = GetMaxCellEnergy(c, id); 
262   if((1. - GetCrossEnergy(c,id) / maxE) > 0.97)
263     return kFALSE;
264
265   Float_t pos1[3];
266   c->GetPosition(pos1);
267   TVector3 clsPos(pos1);
268   Double_t eta = clsPos.Eta();
269
270   if(eta > 0.65 || eta < -0.65)
271     return kFALSE;  
272
273   if (!IsWithinFiducialVolume(id))
274     return kFALSE;
275
276   if(c->GetM02() >0.5 )
277     return kFALSE;
278
279 //  if(c->M20 >)
280
281   return kTRUE;
282
283 }
284 //_____________________________________________________________________
285 Bool_t AliAnalysisTaskPi0V2::IsGoodPion(const TLorentzVector &p1, const TLorentzVector &p2) const
286 {
287   // Is good pion?
288
289
290 //  Double_t asym = TMath::Abs(p1.E()-p2.E())/(p1.E()+p2.E());
291 //  if (asym>0.7)
292 //    return kFALSE;
293
294 //    if (TMath::Abs(p1.Eta()-p2.Eta())>0.5)
295 //      return kFALSE;
296
297   TLorentzVector pion;
298   pion = p1 + p2;
299   Double_t eta = pion.Eta();
300   if (eta<-0.65)
301     return kFALSE;
302   if (eta>0.65)
303     return kFALSE;
304
305   return kTRUE;
306 }
307 //_______________________________________________________________________
308 void AliAnalysisTaskPi0V2::FillPion(const TLorentzVector& p1, const TLorentzVector& p2, Double_t EPV0r, Double_t EPV0A, Double_t EPV0C, Double_t EPTPC)
309 {
310   // Fill histogram.
311
312   if (!IsGoodPion(p1,p2))
313     return;
314   TLorentzVector pion;
315   pion = p1 + p2;
316
317   Double_t mass = pion.M();
318   Double_t pt   = pion.Pt();
319   Double_t phi  = pion.Phi();
320
321   Double_t dphiV0   = phi-EPV0r;
322   Double_t dphiV0A  = phi-EPV0A;
323   Double_t dphiV0C  = phi-EPV0C;
324   Double_t dphiTPC  = phi-EPTPC;
325
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));
330
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();
335
336   //cout<<"cos2V0: "<<cos2phiV0<<"  cos2V0A: "<<cos2phiV0A<<"  cos2V0C: "<<cos2phiV0C<<endl;
337   //cout<<mass<<"  "<<pt<<"  "<<phi<<"  "<<endl;
338   //cout<<" dphiV0: "<<dphiV0<<"    dphiV0A: "<<dphiV0A<<"  dphiV0C: "<<dphiV0C<<"+++++++"<<endl;
339
340   Double_t xV0[5]; // Match ndims in fH  V0 EP
341   xV0[0]       = mass;
342   xV0[1]       = pt;
343   xV0[2]       = fCentrality;
344   xV0[3]       = dphiV0;
345   xV0[4]       = cos2phiV0;
346   fHEPV0r->Fill(xV0);
347
348   Double_t xV0A[5]; // Match ndims in fH V0A EP
349   xV0A[0]       = mass;
350   xV0A[1]       = pt;
351   xV0A[2]       = fCentrality;
352   xV0A[3]       = dphiV0A;
353   xV0A[4]       = cos2phiV0A;
354   fHEPV0A->Fill(xV0A);
355
356   Double_t xV0C[5]; // Match ndims in fH V0C EP
357   xV0C[0]       = mass;
358   xV0C[1]       = pt;
359   xV0C[2]       = fCentrality;
360   xV0C[3]       = dphiV0C;
361   xV0C[4]       = cos2phiV0C;
362   fHEPV0C->Fill(xV0C);
363
364   Double_t xTPC[5]; // Match ndims in fH TPC EP
365   xTPC[0]       = mass;
366   xTPC[1]       = pt;
367   xTPC[2]       = fCentrality;
368   xTPC[3]       = dphiTPC;
369   xTPC[4]       = cos2phiTPC;
370   fHEPTPC->Fill(xTPC);
371
372
373 }
374 //_________________________________________________________________________________________________
375 void AliAnalysisTaskPi0V2::GetMom(TLorentzVector& p, const AliESDCaloCluster *c, Double_t *vertex)
376 {
377   // Calculate momentum.
378   Float_t posMom[3];
379   c->GetPosition(posMom);
380   TVector3 clsPos2(posMom);
381
382   Double_t e   = c->E();
383   Double_t r   = clsPos2.Perp();
384   Double_t eta = clsPos2.Eta();
385   Double_t phi = clsPos2.Phi();
386
387   TVector3 pos;
388   pos.SetPtEtaPhi(r,eta,phi);
389
390   if (vertex) { //calculate direction relative to vertex
391     pos -= vertex;
392   }
393
394   Double_t rad = pos.Mag();
395   p.SetPxPyPzE(e*pos.x()/rad, e*pos.y()/rad, e*pos.z()/rad, e);
396
397 }
398 //________________________________________________________________________
399 void AliAnalysisTaskPi0V2::UserCreateOutputObjects()
400 {
401     // Create histograms
402     // Called once (on the worker node)
403         
404     fOutput = new TList();
405     fOutput->SetOwner();  // IMPORTANT!
406
407     hEvtCount = new TH1F("hEvtCount", " Event Plane", 10, 0.5, 10.5);
408     hEvtCount->GetXaxis()->SetBinLabel(1,"All");
409     hEvtCount->GetXaxis()->SetBinLabel(2,"vert");
410     hEvtCount->GetXaxis()->SetBinLabel(3,"cent");
411     hEvtCount->GetXaxis()->SetBinLabel(4,"EPtask");
412     hEvtCount->GetXaxis()->SetBinLabel(5,"EPvalue");
413     hEvtCount->GetXaxis()->SetBinLabel(6,"Pass");
414     fOutput->Add(hEvtCount);
415     
416     hEPTPC   = new TH2F("hEPTPC",   "EPTPC     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
417     hresoTPC = new TH2F("hresoTPC", "TPc reso  vs cent", 100, 0., 100., 100, 0., 1.);
418     hEPV0    = new TH2F("hEPV0",    "EPV0      vs cent", 100, 0., 100., 100, 0., TMath::Pi());
419     hEPV0A   = new TH2F("hEPV0A",   "EPV0A     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
420     hEPV0C   = new TH2F("hEPV0C",   "EPV0C     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
421     hEPV0Ar  = new TH2F("hEPV0Ar",  "EPV0Ar    vs cent", 100, 0., 100., 100, 0., TMath::Pi());
422     hEPV0Cr  = new TH2F("hEPV0Cr",  "EPV0Cr    vs cent", 100, 0., 100., 100, 0., TMath::Pi());
423     hEPV0r   = new TH2F("hEPV0r",   "EPV0r     vs cent", 100, 0., 100., 100, 0., TMath::Pi());
424     hEPV0AR4 = new TH2F("hEPV0AR4", "EPV0AR4   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
425     hEPV0AR7 = new TH2F("hEPV0AR7", "EPV0AR7   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
426     hEPV0CR0 = new TH2F("hEPV0CR0", "EPV0CR0   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
427     hEPV0CR3 = new TH2F("hEPV0CR3", "EPV0CR3   vs cent", 100, 0., 100., 100, 0., TMath::Pi());
428     fOutput->Add(hEPTPC);
429     fOutput->Add(hresoTPC);
430     fOutput->Add(hEPV0);
431     fOutput->Add(hEPV0A);
432     fOutput->Add(hEPV0C);
433     fOutput->Add(hEPV0Ar);
434     fOutput->Add(hEPV0Cr);
435     fOutput->Add(hEPV0r);
436     fOutput->Add(hEPV0AR4);
437     fOutput->Add(hEPV0AR7);
438     fOutput->Add(hEPV0CR0);
439     fOutput->Add(hEPV0CR3);
440
441     hdifV0A_V0CR0    = new TH2F("hdifV0A_V0CR0",    "EP A-R0 ",  100, 0., 100., 100, -1., 1.);    
442     hdifV0A_V0CR3    = new TH2F("hdifV0A_V0CR3",    "EP A-R3 ",  100, 0., 100., 100, -1., 1.);    
443     hdifV0ACR0_V0CR3 = new TH2F("hdifV0ACR0_V0CR3", "EP R0-R3 ", 100, 0., 100., 100, -1., 1.);    
444     hdifV0C_V0AR4    = new TH2F("hdifV0C_V0AR4",    "EP C-R4 ",  100, 0., 100., 100, -1., 1.);    
445     hdifV0C_V0AR7    = new TH2F("hdifV0C_V0AR7",    "EP C-R7 ",  100, 0., 100., 100, -1., 1.);    
446     hdifV0AR4_V0AR7  = new TH2F("hdifV0AR4_V0AR7",  "EP R4-R7 ", 100, 0., 100., 100, -1., 1.);    
447     fOutput->Add(hdifV0A_V0CR0);
448     fOutput->Add(hdifV0A_V0CR3);
449     fOutput->Add(hdifV0ACR0_V0CR3);
450     fOutput->Add(hdifV0C_V0AR4);
451     fOutput->Add(hdifV0C_V0AR7);
452     fOutput->Add(hdifV0AR4_V0AR7);
453
454     hdifV0A_V0C = new TH2F("hdifV0A_V0C", "EP A-C  ", 100, 0., 100., 100, -1., 1.);
455     hdifV0A_TPC = new TH2F("hdifV0A_TPC", "EP A-TPC", 100, 0., 100., 100, -1., 1.);
456     hdifV0C_TPC = new TH2F("hdifV0C_TPC", "EP C-TPC", 100, 0., 100., 100, -1., 1.);
457     hdifV0C_V0A = new TH2F("hdifV0C_V0A", "EP C-A  ", 100, 0., 100., 100, -1., 1.);
458     fOutput->Add(hdifV0A_V0C);
459     fOutput->Add(hdifV0A_TPC);
460     fOutput->Add(hdifV0C_TPC);
461     fOutput->Add(hdifV0C_V0A);
462
463
464
465     hdifEMC_EP = new TH3F("hdifEMC_EP", "dif phi in EMC with EP",  100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
466     hdifful_EP = new TH3F("hdifful_EP", "dif phi in full with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
467     hdifout_EP = new TH3F("hdifout_EP", "dif phi NOT in EMC with EP", 100, 0., 100., 100, 0., TMath::Pi(), 15, 0., 15.);
468     fOutput->Add(hdifEMC_EP);
469     fOutput->Add(hdifful_EP);
470     fOutput->Add(hdifout_EP);
471
472     hAllcentV0  = new TH1F("hAllcentV0",  "All cent EP V0",  100, 0., TMath::Pi());
473     hAllcentV0r = new TH1F("hAllcentV0r", "All cent EP V0r", 100, 0., TMath::Pi());
474     hAllcentV0A = new TH1F("hAllcentV0A", "All cent EP V0A", 100, 0., TMath::Pi());
475     hAllcentV0C = new TH1F("hAllcentV0C", "All cent EP V0C", 100, 0., TMath::Pi());
476     hAllcentTPC = new TH1F("hAllcentTPC", "All cent EP TPC", 100, 0., TMath::Pi());
477     fOutput->Add(hAllcentV0);
478     fOutput->Add(hAllcentV0r);
479     fOutput->Add(hAllcentV0A);
480     fOutput->Add(hAllcentV0C);
481     fOutput->Add(hAllcentTPC);
482
483     const Int_t ndims = 5;
484     Int_t nMgg=500, nPt=40, nCent=20, nDeltaPhi=315,  ncos2phi=500;
485     Int_t bins[ndims] = {nMgg, nPt, nCent, nDeltaPhi, ncos2phi};
486     Double_t xmin[ndims] = { 0,   0.,  0,   0.,     -1.};
487     Double_t xmax[ndims] = { 0.5, 20., 100, 3.15,   1.};
488     fHEPV0r  = new THnSparseF("fHEPV0r",  "Flow histogram EPV0",  ndims, bins, xmin, xmax);
489     fHEPV0A = new THnSparseF("fHEPV0A",   "Flow histogram EPV0A", ndims, bins, xmin, xmax);
490     fHEPV0C = new THnSparseF("fHEPV0C",   "Flow histogram EPV0C", ndims, bins, xmin, xmax);
491     fHEPTPC = new THnSparseF("fHEPTPC",   "Flow histogram EPTPC", ndims, bins, xmin, xmax);
492     fOutput->Add(fHEPV0r);
493     fOutput->Add(fHEPV0A);
494     fOutput->Add(fHEPV0C);
495     fOutput->Add(fHEPTPC);
496
497     
498
499     PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram
500 }
501
502 //________________________________________________________________________
503 void AliAnalysisTaskPi0V2::UserExec(Option_t *) 
504 {
505     // Main loop
506     // Called for each event
507         
508         
509     // Create pointer to reconstructed event
510    AliVEvent *event = InputEvent();
511    if (!event) { Printf("ERROR: Could not retrieve event"); return; }
512
513     // create pointer to event
514     fESD = dynamic_cast<AliESDEvent*>(event);
515     if (!fESD) {
516         AliError("Cannot get the ESD event");
517         return;
518     }  
519
520     hEvtCount->Fill(1);
521     
522     const AliESDVertex* fvertex = fESD->GetPrimaryVertex();
523     if(TMath::Abs(fvertex->GetZ())>10.)
524       return;
525     Double_t vertex[3] = {fvertex->GetX(), fvertex->GetY(), fvertex->GetZ()};
526
527     hEvtCount->Fill(2);
528
529     if(fESD->GetCentrality()) {
530       fCentrality = 
531         fESD->GetCentrality()->GetCentralityPercentile("V0M");
532     } else{
533            return;
534     }
535
536     hEvtCount->Fill(3);
537
538     AliEventplane *ep = fESD->GetEventplane();
539       if (ep) {
540       if (ep->GetQVector())
541         fEPTPC    = ep->GetQVector()->Phi()/2. ;
542       else
543         fEPTPC = -999.;
544       if (ep->GetQsub1()&&ep->GetQsub2())
545         fEPTPCreso  = TMath::Cos(2.*(ep->GetQsub1()->Phi()/2.-ep->GetQsub2()->Phi()/2.));
546       else
547         fEPTPCreso = -1;
548
549       fEPV0    = ep->GetEventplane("V0",  fESD);
550       fEPV0A   = ep->GetEventplane("V0A", fESD);
551       fEPV0C   = ep->GetEventplane("V0C", fESD);
552       Double_t qx=0, qy=0;
553       Double_t qxr=0, qyr=0;
554       fEPV0Ar  = ep->CalculateVZEROEventPlane(fESD, 4, 5, 2, qxr, qyr);
555       fEPV0Cr  = ep->CalculateVZEROEventPlane(fESD, 2, 3, 2, qx,  qy);
556       qxr += qx;
557       qyr += qy;
558       fEPV0r   = TMath::ATan2(qyr,qxr)/2.;
559       fEPV0AR4 = ep->CalculateVZEROEventPlane(fESD, 4, 2, qx, qy);
560       fEPV0AR5 = ep->CalculateVZEROEventPlane(fESD, 5, 2, qx, qy);
561       fEPV0AR6 = ep->CalculateVZEROEventPlane(fESD, 6, 2, qx, qy);
562       fEPV0AR7 = ep->CalculateVZEROEventPlane(fESD, 7, 2, qx, qy);
563       fEPV0CR0 = ep->CalculateVZEROEventPlane(fESD, 0, 2, qx, qy);
564       fEPV0CR1 = ep->CalculateVZEROEventPlane(fESD, 1, 2, qx, qy);
565       fEPV0CR2 = ep->CalculateVZEROEventPlane(fESD, 2, 2, qx, qy);
566       fEPV0CR3 = ep->CalculateVZEROEventPlane(fESD, 3, 2, qx, qy);
567
568     }
569 //cout<<" fEPV0:"<<fEPV0<<" fEPV0A:"<<fEPV0A<<" fEPV0C:"<<fEPV0C<<" fEPV0Ar:"<<fEPV0Ar<<" fEPV0Cr:"<<fEPV0Cr<<" fEPV0r:"<<fEPV0AR4<<" fEPV0AR7:"<<fEPV0AR7<<" fEPV0CR0:"<<fEPV0CR0<<" fEPV0CR3:"<<fEPV0CR3<<"--------------------------------------------"<<endl;
570     
571     hEvtCount->Fill(4);
572
573     if(fcheckEP2sub){
574       if(fEPV0r<-2. || fEPV0Ar<-2. || fEPV0Cr<-2.) return; 
575     }
576     if(!fcheckEP2sub){
577       if(fEPV0A<-2. || fEPV0C<-2. || fEPV0AR4<-2. || fEPV0AR7<-2. || fEPV0CR0<-2. || fEPV0CR3<-2. || fEPTPC<-2.) return;
578     }
579
580     hEvtCount->Fill(5);
581
582     fEPV0   = TVector2::Phi_0_2pi(fEPV0);    if(fEPV0>TMath::Pi())   fEPV0  = fEPV0  - TMath::Pi();
583     fEPV0r  = TVector2::Phi_0_2pi(fEPV0r);   if(fEPV0r>TMath::Pi())  fEPV0r = fEPV0r - TMath::Pi();
584     fEPV0A  = TVector2::Phi_0_2pi(fEPV0A);   if(fEPV0A>TMath::Pi())  fEPV0A = fEPV0A - TMath::Pi();
585     fEPV0C  = TVector2::Phi_0_2pi(fEPV0C);   if(fEPV0C>TMath::Pi())  fEPV0C = fEPV0C - TMath::Pi();
586     fEPV0Ar = TVector2::Phi_0_2pi(fEPV0Ar);  if(fEPV0Ar>TMath::Pi()) fEPV0Ar = fEPV0Ar - TMath::Pi();
587     fEPV0Cr = TVector2::Phi_0_2pi(fEPV0Cr);  if(fEPV0Cr>TMath::Pi()) fEPV0Cr = fEPV0Cr - TMath::Pi();
588     fEPV0AR4   = TVector2::Phi_0_2pi(fEPV0AR4);    if(fEPV0AR4>TMath::Pi())   fEPV0AR4  = fEPV0AR4 - TMath::Pi();
589     fEPV0AR7   = TVector2::Phi_0_2pi(fEPV0AR7);    if(fEPV0AR7>TMath::Pi())   fEPV0AR7  = fEPV0AR7 - TMath::Pi();
590     fEPV0CR0   = TVector2::Phi_0_2pi(fEPV0CR0);    if(fEPV0CR0>TMath::Pi())   fEPV0CR0  = fEPV0CR0 - TMath::Pi();
591     fEPV0CR3   = TVector2::Phi_0_2pi(fEPV0CR3);    if(fEPV0CR3>TMath::Pi())   fEPV0CR3  = fEPV0CR3 - TMath::Pi();
592
593 //cout<<" EPTPC: "<<fEPTPC<<" reso: "<<fEPTPCreso<<" -------------------"<<endl;
594 //cout<<" cent: "<<fCentrality<<" fEPV0:"<<fEPV0<<" fEPV0A:"<<fEPV0A<<" fEPV0C:"<<fEPV0C<<" fEPV0Ar:"<<fEPV0Ar<<" fEPV0Cr:"<<fEPV0Cr<<" fEPV0r:"<<fEPV0AR4<<" fEPV0AR7:"<<fEPV0AR7<<" fEPV0CR0:"<<fEPV0CR0<<" fEPV0CR3:"<<fEPV0CR3<<"--------------------------------------------"<<endl;
595
596    hEPTPC->Fill(fCentrality,  fEPTPC); 
597    if(fEPTPCreso!=-1) hresoTPC->Fill(fCentrality, fEPTPCreso);
598    hEPV0->Fill(fCentrality,   fEPV0);
599    hEPV0A->Fill(fCentrality,  fEPV0A);
600    hEPV0C->Fill(fCentrality,  fEPV0C);
601    hEPV0Ar->Fill(fCentrality, fEPV0Ar);
602    hEPV0Cr->Fill(fCentrality, fEPV0Cr);
603    hEPV0r->Fill(fCentrality,  fEPV0r);
604    hEPV0AR4->Fill(fCentrality, fEPV0AR4);
605    hEPV0AR7->Fill(fCentrality, fEPV0AR7);
606    hEPV0CR0->Fill(fCentrality, fEPV0CR0);
607    hEPV0CR3->Fill(fCentrality, fEPV0CR3);
608
609    hAllcentV0->Fill(fEPV0);
610    hAllcentV0r->Fill(fEPV0r);
611    hAllcentV0A->Fill(fEPV0A);
612    hAllcentV0C->Fill(fEPV0C);  
613    hAllcentTPC->Fill(fEPTPC);
614
615    hdifV0A_V0CR0->Fill(fCentrality, TMath::Cos(2.*(fEPV0A - fEPV0CR0)));
616    hdifV0A_V0CR3->Fill(fCentrality, TMath::Cos(2.*(fEPV0A - fEPV0CR3)));
617    hdifV0ACR0_V0CR3->Fill(fCentrality, TMath::Cos(2*(fEPV0CR0 - fEPV0CR3)));
618    hdifV0C_V0AR4->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0AR4)));
619    hdifV0C_V0AR7->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0AR7)));
620    hdifV0AR4_V0AR7->Fill(fCentrality, TMath::Cos(2*(fEPV0AR4 - fEPV0AR7)));
621         
622    hdifV0A_V0C->Fill(fCentrality, TMath::Cos(2*(fEPV0A - fEPV0C)));
623    hdifV0A_TPC->Fill(fCentrality, TMath::Cos(2*(fEPV0A - fEPTPC)));
624    hdifV0C_TPC->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPTPC)));
625    hdifV0C_V0A->Fill(fCentrality, TMath::Cos(2*(fEPV0C - fEPV0A)));
626     // Cluster loop for reconstructed event
627    
628
629     Int_t nCluster =  fESD->GetNumberOfCaloClusters(); 
630     for(Int_t i=0; i<nCluster; ++i){
631       AliESDCaloCluster *c1 = fESD->GetCaloCluster(i);
632       if(!IsGoodCluster(c1)) continue;
633       for(Int_t j=i+1; j<nCluster; ++j){
634         AliESDCaloCluster *c2 = fESD->GetCaloCluster(j);
635         if(!IsGoodCluster(c2)) continue;
636         TLorentzVector p1;
637         GetMom(p1, c1, vertex);
638         TLorentzVector p2;
639         GetMom(p2, c2, vertex);
640         FillPion(p1, p2, fEPV0r, fEPV0A, fEPV0C, fEPTPC);
641       } 
642     }
643
644     Int_t nTrack = fESD->GetNumberOfTracks();
645     for(Int_t i=0; i<nTrack; ++i){
646         AliESDtrack* esdtrack = fESD->GetTrack(i); // pointer to reconstructed to track          
647         if(!esdtrack) {
648             AliError(Form("ERROR: Could not retrieve esdtrack %d",i));
649             continue;
650         }
651         Double_t tPhi = esdtrack->Phi();
652         Double_t tPt  = esdtrack->Pt();
653
654         Double_t difTrack = TVector2::Phi_0_2pi(tPhi-fEPV0);  if(difTrack >TMath::Pi()) difTrack -= TMath::Pi();
655         if(esdtrack->IsEMCAL()){        
656           hdifEMC_EP->Fill(fCentrality, difTrack, tPt);
657         }else{
658           hdifout_EP->Fill(fCentrality, difTrack, tPt);
659         }
660         hdifful_EP->Fill(fCentrality,   difTrack, tPt);
661     }
662     hEvtCount->Fill(6);
663
664     // NEW HISTO should be filled before this point, as PostData puts the
665     // information for this iteration of the UserExec in the container
666     PostData(1, fOutput);
667 }
668
669
670 //________________________________________________________________________
671 void AliAnalysisTaskPi0V2::Terminate(Option_t *) 
672 {
673     // Draw result to screen, or perform fitting, normalizations
674     // Called once at the end of the query
675 //    fOutput = dynamic_cast<TList*> (GetOutputData(1));
676  //   if(!fOutput) { Printf("ERROR: could not retrieve TList fOutput"); return; }
677         
678     // Get the physics selection histograms with the selection statistics
679     //AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
680     //AliESDInputHandler *inputH = dynamic_cast<AliESDInputHandler*>(mgr->GetInputEventHandler());
681     //TH2F *histStat = (TH2F*)inputH->GetStatistics();
682    
683    
684     // NEW HISTO should be retrieved from the TList container in the above way,
685     // so it is available to draw on a canvas such as below
686
687 }