]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAnalysisTaskConversionQA.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAnalysisTaskConversionQA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Authors: Svein Lindal, Daniel Lohner                                   *
5  * Version 1.0                                                            *
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 ////////////////////////////////////////////////
17 //---------------------------------------------
18 // QA Task for V0 Reader V1
19 //---------------------------------------------
20 ////////////////////////////////////////////////
21
22 #include "AliAnalysisTaskConversionQA.h"
23 #include "TChain.h"
24 #include "AliAnalysisManager.h"
25 #include "TParticle.h"
26 #include "TVectorF.h"
27 #include "AliPIDResponse.h"
28 #include "TFile.h"
29 #include "AliESDtrackCuts.h"
30 #include "AliAODMCParticle.h"
31 #include "AliAODMCHeader.h"
32
33
34 class iostream;
35
36 using namespace std;
37
38 ClassImp(AliAnalysisTaskConversionQA)
39
40 //________________________________________________________________________
41 AliAnalysisTaskConversionQA::AliAnalysisTaskConversionQA() : AliAnalysisTaskSE(),
42    fV0Reader(NULL),
43    fConversionGammas(NULL),
44    fConversionCuts(NULL),
45    fInputEvent(NULL),
46    fNumberOfESDTracks(0),
47    fMCEvent(NULL),
48    fMCStack(NULL),
49    fTreeQA(NULL),
50    fIsHeavyIon(kFALSE),
51    ffillTree(kFALSE),
52    ffillHistograms(kFALSE),
53    fOutputList(NULL),
54    fTreeList(NULL),
55    fESDList(NULL),
56    hVertexZ(NULL),
57    hNGoodESDTracks(NULL),
58    hNV0Tracks(NULL),
59    hNContributorsVertex(NULL),
60    hITSClusterPhi(NULL),
61    hGammaPt(NULL),
62    hGammaPhi(NULL),
63    hGammaEta(NULL),
64    hGammaChi2perNDF(NULL),
65    hGammaPsiPair(NULL),
66    hGammaArmenteros(NULL),
67    hGammaCosinePointingAngle(NULL),
68    hGammaInvMass(NULL),
69    hElecPt(NULL),
70    hElecEta(NULL),
71    hElecPhi(NULL),
72    hElecNfindableClsTPC(NULL),
73    hPosiNfindableClsTPC(NULL),
74    hElecClsTPC(NULL),
75    hPosiClsTPC(NULL),
76    hElectrondEdxP(NULL),
77    hElectronITSdEdxP(NULL),
78    hElectronTOFP(NULL),
79    hElectronNSigmadEdxP(NULL),
80    hElectronNSigmadEdxEta(NULL),
81    hElectronNSigmaPiondEdxP(NULL),
82    hElectronNSigmaITSP(NULL),
83    hElectronNSigmaTOFP(NULL),
84    hPositrondEdxP(NULL),
85    hPositronITSdEdxP(NULL),
86    hPositronTOFP(NULL),
87    hPositronNSigmadEdxP(NULL),
88    hPositronNSigmadEdxEta(NULL),
89    hPositronNSigmaPiondEdxP(NULL),
90    hPositronNSigmaITSP(NULL),
91    hPositronNSigmaTOFP(NULL),
92 //    hElecAsymP(NULL),
93 //    fTrueList(NULL),
94 //    hTrueResolutionR(NULL),
95 //    hTrueResolutionZ(NULL),
96 //    hTrueResolutionPhi(NULL),
97 //    hTrueGammaPt(NULL),
98 //    hTrueGammaPhi(NULL),
99 //    hTrueGammaEta(NULL),
100 //    hTrueGammaMass(NULL),
101 //    hTrueGammaChi2perNDF(NULL),
102 //    hTrueGammaPsiPair(NULL),
103 //    hTrueGammaQt(NULL),
104 //    hTrueGammaCosinePointingAngle(NULL),
105 //    hTrueGammaXY(NULL),
106 //    hTrueGammaZR(NULL),
107 //    hTrueElecPt(NULL),
108 //    hTrueElecEta(NULL),
109 //    hTrueElecPhi(NULL),
110 //    hTrueElecNfindableClsTPC(NULL),
111 //    hTruePosiNfindableClsTPC(NULL),
112 //    hTrueElecAsymP(NULL),
113    fGammaPt(0),
114    fGammaTheta(0),
115    fGammaChi2NDF(0),
116    fGammaPhotonProp(5),
117    fGammaConvCoord(5),
118    fDaughterProp(24),
119    fKind(0),
120    fIsMC(kFALSE),
121    fnGammaCandidates(1),
122    fMCStackPos(NULL),
123    fMCStackNeg(NULL)
124 {
125
126 }
127
128 AliAnalysisTaskConversionQA::AliAnalysisTaskConversionQA(const char *name) : AliAnalysisTaskSE(name),
129    fV0Reader(NULL),
130    fConversionGammas(NULL),
131    fConversionCuts(NULL),
132    fInputEvent(NULL),
133    fNumberOfESDTracks(0),
134    fMCEvent(NULL),
135    fMCStack(NULL),
136    fTreeQA(NULL),
137    fIsHeavyIon(kFALSE),
138    ffillTree(kFALSE),
139    ffillHistograms(kFALSE),
140    fOutputList(NULL),
141    fTreeList(NULL),
142    fESDList(NULL),
143    hVertexZ(NULL),
144    hNGoodESDTracks(NULL),
145    hNV0Tracks(NULL),
146    hNContributorsVertex(NULL),
147    hITSClusterPhi(NULL),
148    hGammaPt(NULL),
149    hGammaPhi(NULL),
150    hGammaEta(NULL),
151    hGammaChi2perNDF(NULL),
152    hGammaPsiPair(NULL),
153    hGammaArmenteros(NULL),
154    hGammaCosinePointingAngle(NULL),
155    hGammaInvMass(NULL),
156    hElecPt(NULL),
157    hElecEta(NULL),
158    hElecPhi(NULL),
159    hElecNfindableClsTPC(NULL),
160    hPosiNfindableClsTPC(NULL),
161    hElecClsTPC(NULL),
162    hPosiClsTPC(NULL),
163    hElectrondEdxP(NULL),
164    hElectronITSdEdxP(NULL),
165    hElectronTOFP(NULL),
166    hElectronNSigmadEdxP(NULL),
167    hElectronNSigmadEdxEta(NULL),
168    hElectronNSigmaPiondEdxP(NULL),
169    hElectronNSigmaITSP(NULL),
170    hElectronNSigmaTOFP(NULL),
171    hPositrondEdxP(NULL),
172    hPositronITSdEdxP(NULL),
173    hPositronTOFP(NULL),
174    hPositronNSigmadEdxP(NULL),
175    hPositronNSigmadEdxEta(NULL),
176    hPositronNSigmaPiondEdxP(NULL),
177    hPositronNSigmaITSP(NULL),
178    hPositronNSigmaTOFP(NULL),
179 //    hGammaXY(NULL),
180 //    hGammaZR(NULL),
181 //    hElecAsymP(NULL),
182 //    fTrueList(NULL),
183 //    hTrueResolutionR(NULL),
184 //    hTrueResolutionZ(NULL),
185 //    hTrueResolutionPhi(NULL),
186 //    hTrueGammaPt(NULL),
187 //    hTrueGammaPhi(NULL),
188 //    hTrueGammaEta(NULL),
189 //    hTrueGammaMass(NULL),
190 //    hTrueGammaChi2perNDF(NULL),
191 //    hTrueGammaPsiPair(NULL),
192 //    hTrueGammaQt(NULL),
193 //    hTrueGammaCosinePointingAngle(NULL),
194 //    hTrueGammaXY(NULL),
195 //    hTrueGammaZR(NULL),
196 //    hTrueElecPt(NULL),
197 //    hTrueElecEta(NULL),
198 //    hTrueElecPhi(NULL),
199 //    hTrueElecNfindableClsTPC(NULL),
200 //    hTruePosiNfindableClsTPC(NULL),
201 //    hTrueElecAsymP(NULL),
202    fGammaPt(0),
203    fGammaTheta(0),
204    fGammaChi2NDF(0),
205    fGammaPhotonProp(5),
206    fGammaConvCoord(5),
207    fDaughterProp(24),
208    fKind(0),
209    fIsMC(kFALSE),
210    fnGammaCandidates(1),
211    fMCStackPos(NULL),
212    fMCStackNeg(NULL)
213 {
214    // Default constructor
215
216    DefineInput(0, TChain::Class());
217    DefineOutput(1, TList::Class());
218 }
219
220 //________________________________________________________________________
221 AliAnalysisTaskConversionQA::~AliAnalysisTaskConversionQA()
222 {
223    // default deconstructor
224   
225 }
226 //________________________________________________________________________
227 void AliAnalysisTaskConversionQA::UserCreateOutputObjects()
228 {
229    // Create User Output Objects
230
231    if(fOutputList != NULL){
232       delete fOutputList;
233       fOutputList = NULL;
234    }
235    if(fOutputList == NULL){
236       fOutputList = new TList();
237       fOutputList->SetOwner(kTRUE);
238    }
239    
240    if(ffillHistograms){
241          
242       fESDList = new TList();
243       fESDList->SetOwner(kTRUE);
244       fESDList->SetName("ESD QA");
245       fOutputList->Add(fESDList);
246
247       hVertexZ = new TH1F("Vertex_Z","Vertex_Z",300,-15,15);
248       fESDList->Add(hVertexZ);
249       hNContributorsVertex = new TH1I("ContrVertex_Z","ContrVertex_Z",3000,0,3000);
250       fESDList->Add(hNContributorsVertex);
251       if(fIsHeavyIon) hNGoodESDTracks = new TH1I("GoodESDTracks","GoodESDTracks",4000,0,4000);
252       else hNGoodESDTracks = new TH1I("GoodESDTracks","GoodESDTracks",200,0,200);
253       fESDList->Add(hNGoodESDTracks);
254       if(fIsHeavyIon) hNV0Tracks = new TH1I("V0 Multiplicity","V0 Multiplicity",30000,0,30000);
255       else hNV0Tracks = new TH1I("V0 Multiplicity","V0 Multiplicity",2000,0,2000);
256       fESDList->Add(hNV0Tracks);
257
258       hITSClusterPhi = new TH2F("ITSClusterPhi","hITSClusterPhi",72,0,2*TMath::Pi(),7,0,7);
259       fESDList->Add(hITSClusterPhi);
260       hGammaPt = new TH1F("Gamma_Pt","Gamma_Pt",250,0,25);
261       fESDList->Add(hGammaPt);
262       hGammaPhi = new TH1F("Gamma_Phi","Gamma_Phi",360,0,2*TMath::Pi());
263       fESDList->Add(hGammaPhi);
264       hGammaEta = new TH1F("Gamma_Eta","Gamma_Eta",600,-1.5,1.5);
265       fESDList->Add(hGammaEta);
266       hGammaChi2perNDF = new TH1F("Gamma_Chi2perNDF","Gamma_Chi2perNDF",500,0,100);
267       fESDList->Add(hGammaChi2perNDF);
268       hGammaPsiPair = new TH1F("Gamma_PsiPair","Gamma_PsiPair",500,0,2);
269       fESDList->Add(hGammaPsiPair);
270       hGammaArmenteros = new TH2F("Gamma_Armenteros","Gamma_Armenteros",200,-1,1,400,0,0.1);
271       fESDList->Add(hGammaArmenteros);
272       hGammaCosinePointingAngle = new TH1F("Gamma_CosinePointingAngle","Gamma_CosinePointingAngle",1000,-1.,1.);
273       fESDList->Add(hGammaCosinePointingAngle);
274       hGammaInvMass = new TH1F( "Gamma_InvMass","",200, 0, 0.2);
275       fESDList->Add(hGammaInvMass);
276       hElecPt = new TH2F("Electron_Positron_Pt","Electron_Positron_Pt",250,0,25,250,0,25);
277       fESDList->Add(hElecPt);
278       hElecEta = new TH2F("Electron_Positron_Eta","Electron_Positron_Eta",600,-1.5,1.5,600,-1.5,1.5);
279       fESDList->Add(hElecEta);
280       hElecPhi = new TH2F("Electron_Positron_Phi","Electron_Positron_Phi",360,0,2*TMath::Pi(),360,0,2*TMath::Pi());
281       fESDList->Add(hElecPhi);
282       hElecClsTPC = new TH1F("Electron_ClusterTPC","Electron_ClusterTPC",200,0,200);
283       fESDList->Add(hElecClsTPC);
284       hPosiClsTPC = new TH1F("Positron_ClusterTPC","Positron_ClusterTPC",200,0,200);
285       fESDList->Add(hPosiClsTPC);
286       
287       hElecNfindableClsTPC = new TH1F("Electron_findableClusterTPC","Electron_findableClusterTPC",100,0,1);
288       fESDList->Add(hElecNfindableClsTPC);
289       hPosiNfindableClsTPC = new TH1F("Positron_findableClusterTPC","Positron_findableClusterTPC",100,0,1);
290       fESDList->Add(hPosiNfindableClsTPC);
291       
292       hElectrondEdxP =  new TH2F("Electron_dEdx_P","Electron_dEdx_P",100, 0.05, 20, 200, 0, 200);
293       SetLogBinningXTH2(hElectrondEdxP);
294       fESDList->Add(hElectrondEdxP);
295       hPositrondEdxP =  new TH2F("Positron_dEdx_P","Positron_dEdx_P",100, 0.05, 20, 200, 0, 200);
296       SetLogBinningXTH2(hPositrondEdxP);
297       fESDList->Add(hPositrondEdxP);
298       hElectronNSigmadEdxP =  new TH2F("Electron_NSigmadEdx_P","Electron_NSigmadEdx_P",100, 0.05, 20, 200, -10, 10);  
299       SetLogBinningXTH2(hElectronNSigmadEdxP);
300       fESDList->Add(hElectronNSigmadEdxP);
301           hElectronNSigmadEdxEta =  new TH2F("Electron_NSigmadEdx_Eta","Electron_NSigmadEdx_Eta",140, -1.4, 1.4, 200, -10, 10);  
302           fESDList->Add(hElectronNSigmadEdxEta);
303       hPositronNSigmadEdxP =  new TH2F("Positron_NSigmadEdx_P","Positron_NSigmadEdx_P",100, 0.05, 20, 200, -10, 10);
304       SetLogBinningXTH2(hPositronNSigmadEdxP);
305       fESDList->Add(hPositronNSigmadEdxP);
306           hPositronNSigmadEdxEta =  new TH2F("Positron_NSigmadEdx_Eta","Positron_NSigmadEdx_Eta",140, -1.4, 1.4, 200, -10, 10);  
307           fESDList->Add(hPositronNSigmadEdxEta);
308       hElectronNSigmaPiondEdxP =  new TH2F("Electron_NSigmaPiondEdx_P","Electron_NSigmaPiondEdx_P",100, 0.05, 20, 200, -10, 10);  
309       SetLogBinningXTH2(hElectronNSigmaPiondEdxP);
310       fESDList->Add(hElectronNSigmaPiondEdxP);
311       hPositronNSigmaPiondEdxP =  new TH2F("Positron_NSigmaPiondEdx_P","Positron_NSigmaPiondEdx_P",100, 0.05, 20, 200, -10, 10);
312       SetLogBinningXTH2(hPositronNSigmaPiondEdxP);
313       fESDList->Add(hPositronNSigmaPiondEdxP);
314       
315       hElectronTOFP =  new TH2F("Electron_TOF_P","Electron_TOF_P",100, 0.05, 20, 600, -1000, 29000);
316       SetLogBinningXTH2(hElectronTOFP);
317       fESDList->Add(hElectronTOFP);
318       hPositronTOFP =  new TH2F("Positron_TOF_P","Positron_TOF_P",100, 0.05, 20, 600, -1000, 29000);
319       SetLogBinningXTH2(hPositronTOFP);
320       fESDList->Add(hPositronTOFP);
321       hElectronNSigmaTOFP =  new TH2F("Electron_NSigmaTOF_P","Electron_NSigmaTOF_P",100, 0.05, 20, 200, -10, 10);
322       SetLogBinningXTH2(hElectronNSigmaTOFP);
323       fESDList->Add(hElectronNSigmaTOFP);
324       hPositronNSigmaTOFP =  new TH2F("Positron_NSigmaTOF_P","Positron_NSigmaTOF_P",100, 0.05, 20, 200, -10, 10);
325       SetLogBinningXTH2(hPositronNSigmaTOFP);
326       fESDList->Add(hPositronNSigmaTOFP);
327       
328       hElectronITSdEdxP  =  new TH2F("Electron_ITSdEdx_P","Electron_ITSdEdx_P",100, 0.05, 20, 200, 0, 200);
329       SetLogBinningXTH2(hElectronITSdEdxP);
330       fESDList->Add(hElectronITSdEdxP);
331       hPositronITSdEdxP =  new TH2F("Positron_ITSdEdx_P","Positron_ITSdEdx_P",100, 0.05, 20, 200, 0, 200);
332       SetLogBinningXTH2(hPositronITSdEdxP);
333       fESDList->Add(hPositronITSdEdxP);
334       hElectronNSigmaITSP =  new TH2F("Electron_NSigmaITS_P","Electron_NSigmaITS_P",100, 0.05, 20, 200, -10, 10);
335       SetLogBinningXTH2(hElectronNSigmaITSP);
336       fESDList->Add(hElectronNSigmaITSP);
337       hPositronNSigmaITSP =  new TH2F("Positron_NSigmaITS_P","Positron_NSigmaITS_P",100, 0.05, 20, 200, -10, 10);
338       SetLogBinningXTH2(hPositronNSigmaITSP);
339       fESDList->Add(hPositronNSigmaITSP);
340       
341       
342       
343 //       hGammaXY = new TH2F("Gamma_ConversionPoint_XY","Gamma_ConversionPoint_XY",960,-120,120,960,-120,120);
344 //       fESDList->Add(hGammaXY);
345 //       hGammaZR= new TH2F("Gamma_ConversionPoint_ZR","Gamma_ConversionPoint_ZR",1200,-150,150,480,0,120);
346 //       fESDList->Add(hGammaZR);
347
348
349 //       hElecAsymP = new TH2F("Electron_Asym_vs_P", "Electron_Asym_vs_P",200,0.,20.,200,0.,1.); 
350 //       fESDList->Add(hElecAsymP);
351
352 //       if(fIsMC){
353 //          fTrueList = new TList();
354 //          fTrueList->SetOwner(kTRUE);
355 //          fTrueList->SetName("True QA");
356 //          fOutputList->Add(fTrueList);
357 // 
358 //          hTrueResolutionR = new TH2F("True_ConversionPointResolution_R","True_ConversionPointResolution_R",240,0,120,200,-20,20);
359 //          fTrueList->Add(hTrueResolutionR);
360 //          hTrueResolutionZ = new TH2F("True_ConversionPointResolution_Z","True_ConversionPointResolution_Z",480,-120,120,200,-20,20);
361 //          fTrueList->Add(hTrueResolutionZ);
362 //          hTrueResolutionPhi = new TH2F("True_ConversionPointResolution_Phi","True_ConversionPointResolution_Phi",360,0,2*TMath::Pi(),200,-TMath::Pi()/30., TMath::Pi()/30.);
363 //          fTrueList->Add(hTrueResolutionPhi);
364 // 
365 //          hTrueGammaPt = new TH1F("True_Gamma_Pt","True_Gamma_Pt",250,0,25);
366 //          fTrueList->Add(hTrueGammaPt);
367 //          hTrueGammaPhi = new TH1F("True_Gamma_Phi","True_Gamma_Phi",360,0,2*TMath::Pi());
368 //          fTrueList->Add(hTrueGammaPhi);
369 //          hTrueGammaEta = new TH1F("True_Gamma_Eta","True_Gamma_Eta",600,-1.5,1.5);
370 //          fTrueList->Add(hTrueGammaEta);
371 //          hTrueGammaMass = new TH1F("True_Gamma_Mass","True_Gamma_Mass",1000,0,0.3);
372 //          fTrueList->Add(hTrueGammaMass);
373 //          hTrueGammaChi2perNDF = new TH1F("True_Gamma_Chi2perNDF","True_Gamma_Chi2perNDF",500,0,100);
374 //          fTrueList->Add(hTrueGammaChi2perNDF);
375 //          hTrueGammaPsiPair = new TH1F("True_Gamma_PsiPair","True_Gamma_PsiPair",500,0,2);
376 //          fTrueList->Add(hTrueGammaPsiPair);
377 //          hTrueGammaQt = new TH1F("True_Gamma_Qt","True_Gamma_Qt",400,0,0.1);
378 //          fTrueList->Add(hTrueGammaQt);
379 //          hTrueGammaCosinePointingAngle = new TH1F("True_Gamma_CosinePointingAngle","True_Gamma_CosinePointingAngle",900,0.7,1.);
380 //          fTrueList->Add(hTrueGammaCosinePointingAngle);
381 //          hTrueGammaXY = new TH2F("True_Gamma_ConversionPoint_XY","True_Gamma_ConversionPoint_XY",960,-120,120,960,-120,120);
382 //          fTrueList->Add(hTrueGammaXY);
383 //          hTrueGammaZR= new TH2F("TrueGamma_ConversionPoint_ZR","TrueGamma_ConversionPoint_ZR",1200,-150,150,480,0,120);
384 //          fTrueList->Add(hTrueGammaZR);
385 // 
386 //          hTrueElecPt = new TH2F("True_Electron_Positron_Pt","True_Electron_Positron_Pt",250,0,25,250,0,25);
387 //          fTrueList->Add(hTrueElecPt);
388 //          hTrueElecEta = new TH2F("True_Electron_Positron_Eta","True_Electron_Positron_Eta",600,-1.5,1.5,600,-1.5,1.5);
389 //          fTrueList->Add(hTrueElecEta);
390 //          hTrueElecPhi = new TH2F("True_Electron_Positron_Phi","True_Electron_Positron_Phi",360,0,2*TMath::Pi(),360,0,2*TMath::Pi());
391 //          fTrueList->Add(hTrueElecPhi);
392 //          hTrueElecNfindableClsTPC = new TH1F("True_Electron_findableClusterTPC","True_Electron_findableClusterTPC",100,0,1);
393 //          fTrueList->Add(hTrueElecNfindableClsTPC);
394 //          hTruePosiNfindableClsTPC = new TH1F("True_Positron_findableClusterTPC","True_Positron_findableClusterTPC",100,0,1);
395 //          fTrueList->Add(hTruePosiNfindableClsTPC);
396 //                               hTrueElecAsymP = new TH2F("True_Electron_Asym_vs_P", "True_Electron_Asym_vs_P",200,0.,20.,200,0.,1.); 
397 //                               fTrueList->Add(hTrueElecAsymP);
398 //       }
399       if(fConversionCuts->GetCutHistograms()){
400          fOutputList->Add(fConversionCuts->GetCutHistograms());
401       }
402    }
403    
404    if(ffillTree){
405       fTreeList = new TList();
406       fTreeList->SetOwner(kTRUE);
407       fTreeList->SetName("TreeList");
408       fOutputList->Add(fTreeList);
409
410       fTreeQA = new TTree("PhotonQA","PhotonQA");   
411      
412       fTreeQA->Branch("daughterProp",&fDaughterProp);
413       fTreeQA->Branch("recCords",&fGammaConvCoord);
414       fTreeQA->Branch("photonProp",&fGammaPhotonProp);
415       fTreeQA->Branch("pt",&fGammaPt,"fGammaPt/F");
416       fTreeQA->Branch("theta",&fGammaTheta,"fGammaTheta/F");
417       fTreeQA->Branch("chi2ndf",&fGammaChi2NDF,"fGammaChi2NDF/F");
418       if (fIsMC) {
419          fTreeQA->Branch("kind",&fKind,"fKind/b");
420       }   
421       fTreeList->Add(fTreeQA);
422
423    }
424
425    fV0Reader=(AliV0ReaderV1*)AliAnalysisManager::GetAnalysisManager()->GetTask("V0ReaderV1");
426    
427    PostData(1, fOutputList);
428 }
429 //_____________________________________________________________________________
430 Bool_t AliAnalysisTaskConversionQA::Notify()
431 {
432    if(!fConversionCuts->GetDoEtaShift()) return kTRUE;; // No Eta Shift requested, continue
433       
434    if(fConversionCuts->GetEtaShift() == 0.0){ // Eta Shift requested but not set, get shift automatically
435       fConversionCuts->GetCorrectEtaShiftFromPeriod(fV0Reader->GetPeriodName());
436       fConversionCuts->DoEtaShift(kFALSE); // Eta Shift Set, make sure that it is called only once
437       return kTRUE;
438    }
439    else{
440       printf(" Gamma Conversion QA Task %s :: Eta Shift Manually Set to %f \n\n",
441              (fConversionCuts->GetCutNumber()).Data(),fConversionCuts->GetEtaShift());
442       fConversionCuts->DoEtaShift(kFALSE); // Eta Shift Set, make sure that it is called only once
443    }
444    
445    return kTRUE;
446 }
447 //________________________________________________________________________
448 void AliAnalysisTaskConversionQA::UserExec(Option_t *){
449
450
451
452    Int_t eventQuality = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetEventQuality();
453    if(eventQuality != 0){// Event Not Accepted
454       return;
455    }
456    fInputEvent = InputEvent();
457    if(fIsMC) fMCEvent = MCEvent();
458    if(fMCEvent && fInputEvent->IsA()==AliESDEvent::Class()){ fMCStack = fMCEvent->Stack(); }
459
460    Int_t eventNotAccepted =
461       fConversionCuts->IsEventAcceptedByConversionCut(fV0Reader->GetConversionCuts(),fInputEvent,fMCEvent,fIsHeavyIon);
462    if(eventNotAccepted) return; // Check Centrality, PileUp, SDD and V0AND --> Not Accepted => eventQuality = 1
463
464    fConversionGammas=fV0Reader->GetReconstructedGammas();
465
466    if(fMCEvent){
467       if(fConversionCuts->GetSignalRejection() != 0){
468          if(fInputEvent->IsA()==AliESDEvent::Class()){
469                fConversionCuts->GetNotRejectedParticles(fConversionCuts->GetSignalRejection(),
470                                                         fConversionCuts->GetAcceptedHeader(),
471                                                         fMCEvent);
472          }
473          else if(fInputEvent->IsA()==AliAODEvent::Class()){
474             fConversionCuts->GetNotRejectedParticles(fConversionCuts->GetSignalRejection(),
475                                                       fConversionCuts->GetAcceptedHeader(),
476                                                       fInputEvent);
477          }
478       }
479    }
480
481    if(ffillHistograms){
482       CountTracks();
483       hVertexZ->Fill(fInputEvent->GetPrimaryVertex()->GetZ());
484       hNContributorsVertex->Fill(fConversionCuts->GetNumberOfContributorsVtx(fInputEvent));
485       hNGoodESDTracks->Fill(fNumberOfESDTracks);
486       hNV0Tracks->Fill(fInputEvent->GetVZEROData()->GetMTotV0A()+fInputEvent->GetVZEROData()->GetMTotV0C());
487    }
488
489    if(fMCEvent && fInputEvent->IsA()==AliAODEvent::Class() && !(fV0Reader->AreAODsRelabeled())){
490       RelabelAODPhotonCandidates(kTRUE);    // In case of AODMC relabeling MC
491       fV0Reader->RelabelAODs(kTRUE);
492    }
493           
494           
495    for(Int_t firstGammaIndex=0;firstGammaIndex<fConversionGammas->GetEntriesFast();firstGammaIndex++){
496       AliAODConversionPhoton *gamma=dynamic_cast<AliAODConversionPhoton*>(fConversionGammas->At(firstGammaIndex));
497       if (gamma==NULL) continue;
498       if(fMCEvent && fConversionCuts->GetSignalRejection() != 0){
499          if(!fConversionCuts->IsParticleFromBGEvent(gamma->GetMCLabelPositive(), fMCStack, fInputEvent))
500             continue;
501          if(!fConversionCuts->IsParticleFromBGEvent(gamma->GetMCLabelNegative(), fMCStack, fInputEvent))
502             continue;
503       }
504       if(!fConversionCuts->PhotonIsSelected(gamma,fInputEvent)){
505          continue;
506       }
507
508       if(ffillTree) ProcessQATree(gamma);
509       if(ffillHistograms) ProcessQA(gamma);
510    }
511    if(fMCEvent && fInputEvent->IsA()==AliAODEvent::Class())
512       RelabelAODPhotonCandidates(kFALSE);    // In case of AODMC relabeling MC
513       
514    PostData(1, fOutputList);
515 }
516
517
518 ///________________________________________________________________________
519 void AliAnalysisTaskConversionQA::ProcessQATree(AliAODConversionPhoton *gamma){
520
521    // Fill Histograms for QA and MC
522    AliVEvent* event = (AliVEvent*) InputEvent();
523    
524    
525    
526    AliPIDResponse* pidResonse = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetPIDResponse();
527
528    fGammaPt = gamma->GetPhotonPt();
529    
530    fGammaTheta = gamma->Theta();
531    fGammaChi2NDF = gamma->GetChi2perNDF();
532    
533    fGammaPhotonProp(0)  = gamma->GetArmenterosQt();
534    fGammaPhotonProp(1)  = gamma->GetArmenterosAlpha();
535    fGammaPhotonProp(2)  = gamma->GetPsiPair();
536    fGammaPhotonProp(3) = fConversionCuts->GetCosineOfPointingAngle(gamma,event);
537         fGammaPhotonProp(4) = gamma->GetMass();
538    
539    fGammaConvCoord(0) = gamma->GetConversionX();
540    fGammaConvCoord(1) = gamma->GetConversionY();
541    fGammaConvCoord(2) = gamma->GetConversionZ();
542    fGammaConvCoord(3) = gamma->GetConversionRadius();
543    fGammaConvCoord(4) = gamma->GetPhotonPhi();
544    
545    AliVTrack * negTrack = fConversionCuts->GetTrack(event, gamma->GetTrackLabelNegative());
546    AliVTrack * posTrack = fConversionCuts->GetTrack(event, gamma->GetTrackLabelPositive());
547
548    
549    if(!negTrack||!posTrack)return;
550
551    fKind = 9;
552    if(fMCEvent && fInputEvent->IsA()==AliESDEvent::Class()){
553       fKind = IsTruePhotonESD(gamma);
554    } else if (fMCEvent && fInputEvent->IsA()==AliAODEvent::Class()){
555 //        cout << "entering IsTruePhotonAOD" << endl;
556       fKind = IsTruePhotonAOD(gamma);   
557    }
558
559    fDaughterProp(0) =  posTrack->Pt();
560    fDaughterProp(7) =  negTrack->Pt();
561    fDaughterProp(1) =  posTrack->Theta();
562    fDaughterProp(8) =  negTrack->Theta();
563    // dEdx TPC
564    fDaughterProp(2) =  posTrack->GetTPCsignal();
565    fDaughterProp(3) =  pidResonse->NumberOfSigmasTPC(posTrack,AliPID::kElectron);
566    fDaughterProp(22) =  pidResonse->NumberOfSigmasTPC(posTrack,AliPID::kPion);
567    fDaughterProp(9) =  negTrack->GetTPCsignal();
568    fDaughterProp(10) =  pidResonse->NumberOfSigmasTPC(negTrack,AliPID::kElectron);
569    fDaughterProp(23) =  pidResonse->NumberOfSigmasTPC(negTrack,AliPID::kPion);
570    Int_t nPosClusterITS = 0;
571    Int_t nNegClusterITS = 0;
572    for(Int_t itsLayer = 0; itsLayer<6;itsLayer++){
573       if(TESTBIT(negTrack->GetITSClusterMap(),itsLayer)){
574          nNegClusterITS++;
575       }
576       if(TESTBIT(posTrack->GetITSClusterMap(),itsLayer)){
577          nPosClusterITS++;
578       }
579    }
580    
581    // ITS signal
582    fDaughterProp(14) =  (Float_t)nPosClusterITS;
583    fDaughterProp(15) =  (Float_t)nNegClusterITS;
584    if (nPosClusterITS > 0 ){
585       fDaughterProp(16) =  posTrack->GetITSsignal();
586       fDaughterProp(20) =  pidResonse->NumberOfSigmasITS(posTrack,AliPID::kElectron);
587    } else {
588       fDaughterProp(16) =  1000;
589       fDaughterProp(20) =  20;
590    }
591    if (nNegClusterITS > 0 ){
592       fDaughterProp(17) =  negTrack->GetITSsignal();
593       fDaughterProp(21) =  pidResonse->NumberOfSigmasITS(negTrack,AliPID::kElectron);
594    } else {
595       fDaughterProp(17) =  1000;
596       fDaughterProp(21) =  20;
597    }
598
599    // TOF 
600    if((posTrack->GetStatus() & AliESDtrack::kTOFpid) && !(posTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
601       Double_t t0pos = pidResonse->GetTOFResponse().GetStartTime(posTrack->P());
602       Double_t timesPos[9];
603       posTrack->GetIntegratedTimes(timesPos);
604       Double_t TOFsignalPos =   posTrack->GetTOFsignal();
605       Double_t dTpos = TOFsignalPos - t0pos - timesPos[0];
606       fDaughterProp(4) =  dTpos;
607       fDaughterProp(5) =  pidResonse->NumberOfSigmasTOF(posTrack, AliPID::kElectron);
608    } else {
609       fDaughterProp(4) =  20000;
610       fDaughterProp(5) =  -20;
611    }
612    if((negTrack->GetStatus() & AliESDtrack::kTOFpid) && !(negTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
613       Double_t t0neg = pidResonse->GetTOFResponse().GetStartTime(negTrack->P());
614       Double_t timesNeg[9];
615       negTrack->GetIntegratedTimes(timesNeg);
616       Double_t TOFsignalNeg =   negTrack->GetTOFsignal();
617       Double_t dTneg = TOFsignalNeg - t0neg - timesNeg[0];
618       fDaughterProp(11) =  dTneg;
619       fDaughterProp(12) =  pidResonse->NumberOfSigmasTOF(negTrack, AliPID::kElectron);
620    } else {
621       fDaughterProp(11) =  20000;
622       fDaughterProp(12) =  -20;
623    }
624
625    fDaughterProp(6) =  (Float_t)posTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius()));
626    fDaughterProp(18) =  posTrack->GetNcls(1);
627    fDaughterProp(13) =  (Float_t)negTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius()));
628    fDaughterProp(19) =  negTrack->GetNcls(1);
629    
630    if (fTreeQA){
631       fTreeQA->Fill();
632    }
633 }
634
635 //_____________________________________________________________________________________________________
636 void AliAnalysisTaskConversionQA::ProcessQA(AliAODConversionPhoton *gamma){
637
638    AliPIDResponse* pidResonse = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetPIDResponse();
639    
640    // Fill Histograms for QA and MC
641
642    hGammaPt->Fill(gamma->GetPhotonPt());
643    hGammaPhi->Fill(gamma->GetPhotonPhi());
644    hGammaEta->Fill(gamma->Eta());
645    hGammaChi2perNDF->Fill(gamma->GetChi2perNDF());
646    hGammaPsiPair->Fill(gamma->GetPsiPair());
647    hGammaArmenteros->Fill(gamma->GetArmenterosAlpha(),gamma->GetArmenterosQt());
648    hGammaCosinePointingAngle->Fill(fConversionCuts->GetCosineOfPointingAngle(gamma,fInputEvent));
649    hGammaInvMass->Fill(gamma->GetMass());
650 //    hGammaXY->Fill(gamma->GetConversionX(),gamma->GetConversionY());
651 //    hGammaZR->Fill(gamma->GetConversionZ(),gamma->GetConversionRadius());
652
653    AliVTrack * negTrack = fConversionCuts->GetTrack(fInputEvent, gamma->GetTrackLabelNegative());
654    AliVTrack * posTrack = fConversionCuts->GetTrack(fInputEvent, gamma->GetTrackLabelPositive());
655    if(!negTrack||!posTrack)return;
656
657
658    hElecPt->Fill(negTrack->Pt(),posTrack->Pt());
659    hElecEta->Fill(negTrack->Eta(),posTrack->Eta());
660    hElecPhi->Fill(negTrack->Phi(),posTrack->Phi());
661
662    hElecNfindableClsTPC->Fill((Float_t)posTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius())));
663    hPosiNfindableClsTPC->Fill((Float_t)negTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius())));
664    hElecClsTPC->Fill(negTrack->GetNcls(1));
665    hPosiClsTPC->Fill(posTrack->GetNcls(1));
666    //TPC dEdx
667    hElectrondEdxP->Fill(negTrack->P() ,negTrack->GetTPCsignal());
668    hElectronNSigmadEdxP->Fill(negTrack->P() ,pidResonse->NumberOfSigmasTPC(negTrack, AliPID::kElectron));
669    hElectronNSigmadEdxEta->Fill(negTrack->Eta() ,pidResonse->NumberOfSigmasTPC(negTrack, AliPID::kElectron));
670    hElectronNSigmaPiondEdxP->Fill(negTrack->P() ,pidResonse->NumberOfSigmasTPC(negTrack, AliPID::kPion));
671    hPositrondEdxP->Fill(posTrack->P() ,posTrack->GetTPCsignal());
672    hPositronNSigmadEdxP->Fill(posTrack->P() ,pidResonse->NumberOfSigmasTPC(posTrack, AliPID::kElectron));
673    hPositronNSigmadEdxEta->Fill(posTrack->Eta() ,pidResonse->NumberOfSigmasTPC(posTrack, AliPID::kElectron));
674    hPositronNSigmaPiondEdxP->Fill(posTrack->P() ,pidResonse->NumberOfSigmasTPC(posTrack, AliPID::kPion));
675    
676    //TOF signal
677    if((negTrack->GetStatus() & AliESDtrack::kTOFpid) && !(negTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
678       Double_t t0neg = pidResonse->GetTOFResponse().GetStartTime(negTrack->P());
679       Double_t timesNeg[9];
680       negTrack->GetIntegratedTimes(timesNeg);
681       Double_t TOFsignalNeg = negTrack->GetTOFsignal();
682       Double_t dTneg = TOFsignalNeg - t0neg - timesNeg[0];
683       hElectronTOFP->Fill(negTrack->P() ,dTneg);
684       hElectronNSigmaTOFP->Fill(negTrack->P() ,pidResonse->NumberOfSigmasTOF(negTrack, AliPID::kElectron));
685    }
686    if((posTrack->GetStatus() & AliESDtrack::kTOFpid) && !(posTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
687       Double_t t0pos = pidResonse->GetTOFResponse().GetStartTime(posTrack->P());
688       Double_t timesPos[9];
689       posTrack->GetIntegratedTimes(timesPos);
690       Double_t TOFsignalPos = posTrack->GetTOFsignal();
691       Double_t dTpos = TOFsignalPos - t0pos - timesPos[0];
692       hPositronTOFP->Fill(posTrack->P() ,dTpos);
693       hPositronNSigmaTOFP->Fill(posTrack->P() ,pidResonse->NumberOfSigmasTOF(posTrack, AliPID::kElectron));
694    }
695    
696    Int_t nPosClusterITS = 0;
697    Int_t nNegClusterITS = 0;
698    for(Int_t itsLayer = 0; itsLayer<6;itsLayer++){
699       if(TESTBIT(negTrack->GetITSClusterMap(),itsLayer)){
700          nNegClusterITS++;
701       }
702       if(TESTBIT(posTrack->GetITSClusterMap(),itsLayer)){
703          nPosClusterITS++;
704       }
705    }
706    Double_t negtrackPhi = negTrack->Phi();
707    Double_t postrackPhi = posTrack->Phi();
708    hITSClusterPhi->Fill(negtrackPhi,nNegClusterITS);
709    hITSClusterPhi->Fill(postrackPhi,nPosClusterITS);
710    
711    // ITS signal
712    if (nPosClusterITS > 0 ){
713       hPositronITSdEdxP->Fill(posTrack->P() ,posTrack->GetITSsignal());
714       hPositronNSigmaITSP->Fill(posTrack->P() ,pidResonse->NumberOfSigmasITS(posTrack,AliPID::kElectron));
715    } 
716    if (nNegClusterITS > 0 ){
717       hElectronITSdEdxP->Fill(negTrack->P() ,negTrack->GetITSsignal());
718       hElectronNSigmaITSP->Fill(negTrack->P() ,pidResonse->NumberOfSigmasITS(negTrack,AliPID::kElectron));
719    }
720
721    
722 }
723
724
725 //________________________________________________________________________
726 void AliAnalysisTaskConversionQA::CountTracks(){
727
728    if(fInputEvent->IsA()==AliESDEvent::Class()){
729    // Using standard function for setting Cuts
730    Bool_t selectPrimaries=kTRUE;
731    AliESDtrackCuts *EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
732    EsdTrackCuts->SetMaxDCAToVertexZ(2);
733    EsdTrackCuts->SetEtaRange(-0.8, 0.8);
734    EsdTrackCuts->SetPtRange(0.15);
735       fNumberOfESDTracks = 0;
736       for(Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++){
737          AliESDtrack* curTrack = (AliESDtrack*) fInputEvent->GetTrack(iTracks);
738          if(!curTrack) continue;
739          // if(fMCEvent && ((AliConversionCuts*)fCutArray->At(fiCut))->GetSignalRejection() != 0){
740          //    if(!((AliConversionCuts*)fCutArray->At(fiCut))->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack)) continue;
741          // }
742          if(EsdTrackCuts->AcceptTrack(curTrack) ) fNumberOfESDTracks++;
743       }
744       delete EsdTrackCuts;
745       EsdTrackCuts=0x0;
746    }
747    else if(fInputEvent->IsA()==AliAODEvent::Class()){
748       fNumberOfESDTracks = 0;
749       for(Int_t iTracks = 0; iTracks<fInputEvent->GetNumberOfTracks(); iTracks++){
750          AliAODTrack* curTrack = (AliAODTrack*) fInputEvent->GetTrack(iTracks);
751          if(!curTrack->IsPrimaryCandidate()) continue;
752          if(abs(curTrack->Eta())>0.8) continue;
753          if(curTrack->Pt()<0.15) continue;
754          if(abs(curTrack->ZAtDCA())>2) continue;
755          fNumberOfESDTracks++;
756       }
757    }
758    
759    return;
760 }
761
762 UInt_t AliAnalysisTaskConversionQA::IsTruePhotonESD(AliAODConversionPhoton *TruePhotonCandidate)
763 {
764
765
766         UInt_t kind = 9;
767         TParticle *posDaughter = TruePhotonCandidate->GetPositiveMCDaughter(fMCStack);
768         TParticle *negDaughter = TruePhotonCandidate->GetNegativeMCDaughter(fMCStack);
769         TParticle *Photon = TruePhotonCandidate->GetMCParticle(fMCStack);
770         Int_t motherLabelPhoton; 
771         Int_t pdgCodePos = 0; 
772         Int_t pdgCodeNeg = 0; 
773         Int_t pdgCode = 0; 
774
775
776         if(posDaughter == NULL || negDaughter == NULL) {
777                 kind = 9;
778                 //              return kFALSE; // One particle does not exist
779    
780   } else if( posDaughter->GetMother(0) != negDaughter->GetMother(0)  || (posDaughter->GetMother(0) == negDaughter->GetMother(0) && posDaughter->GetMother(0) ==-1)) {
781                 //              kind = 1;
782           return 1;
783                 pdgCodePos=TMath::Abs(posDaughter->GetPdgCode());
784                 pdgCodeNeg=TMath::Abs(negDaughter->GetPdgCode());
785                 if(pdgCodePos==11 && pdgCodeNeg==11) return 10; //Electron Combinatorial
786                 if(pdgCodePos==11 && pdgCodeNeg==11 && 
787                         (posDaughter->GetMother(0) == negDaughter->GetMother(0) && posDaughter->GetMother(0) ==-1)) return 15; //direct Electron Combinatorial
788                                 
789                 if(pdgCodePos==211 && pdgCodeNeg==211) return 11; //Pion Combinatorial
790                 if((pdgCodePos==211 && pdgCodeNeg==2212) ||(pdgCodePos==2212 && pdgCodeNeg==211))       return 12; //Pion, Proton Combinatorics
791                 if((pdgCodePos==211 && pdgCodeNeg==11) ||(pdgCodePos==11 && pdgCodeNeg==211)) return 13; //Pion, Electron Combinatorics
792                 if(pdgCodePos==321 || pdgCodeNeg==321) return 14; //Kaon combinatorics
793         }else{          
794  
795                 pdgCodePos=posDaughter->GetPdgCode();
796                 pdgCodeNeg=negDaughter->GetPdgCode();
797           motherLabelPhoton= Photon->GetMother(0);
798                 if ( TruePhotonCandidate->GetMCParticle(fMCStack)->GetPdgCode()) pdgCode = TruePhotonCandidate->GetMCParticle(fMCStack)->GetPdgCode(); 
799
800                 if(TMath::Abs(pdgCodePos)!=11 || TMath::Abs(pdgCodeNeg)!=11) return 2; // true from hadronic decays
801                 else if ( !(pdgCodeNeg==pdgCodePos)){
802                         if(pdgCode == 111) return 3; // pi0 Dalitz
803                         else if (pdgCode == 221) return 4; // eta Dalitz
804                         else if (!(negDaughter->GetUniqueID() != 5 || posDaughter->GetUniqueID() !=5)){
805                                 if(pdgCode == 22 && motherLabelPhoton < fMCStack->GetNprimary()){
806                                         return 0; // primary photons
807                                 } else if (pdgCode == 22){
808                                         return 5; //secondary photons
809                                 }
810                         }
811                 }
812         }
813
814         return kind;
815 }
816
817 //________________________________________________________________________
818 UInt_t AliAnalysisTaskConversionQA::IsTruePhotonAOD(AliAODConversionPhoton *TruePhotonCandidate)
819 {   
820
821         UInt_t kind = 9;
822         TClonesArray *AODMCTrackArray = dynamic_cast<TClonesArray*>(fInputEvent->FindListObject(AliAODMCParticle::StdBranchName()));
823         AliAODMCParticle *posDaughter = (AliAODMCParticle*) AODMCTrackArray->At(TruePhotonCandidate->GetMCLabelPositive());
824         AliAODMCParticle *negDaughter = (AliAODMCParticle*) AODMCTrackArray->At(TruePhotonCandidate->GetMCLabelNegative());
825         AliAODMCParticle *Photon = (AliAODMCParticle*) AODMCTrackArray->At(posDaughter->GetMother());
826         Int_t pdgCodePos = 0; 
827         Int_t pdgCodeNeg = 0; 
828         Int_t pdgCode = 0; 
829         if(posDaughter == NULL || negDaughter == NULL) {
830                 kind = 9;
831    } else if( posDaughter->GetMother() != negDaughter->GetMother()  || (posDaughter->GetMother() == negDaughter->GetMother() && posDaughter->GetMother() ==-1)) {
832                 kind = 1;
833                 pdgCodePos=TMath::Abs(posDaughter->GetPdgCode());
834                 pdgCodeNeg=TMath::Abs(negDaughter->GetPdgCode());
835                 if(pdgCodePos==11 && pdgCodeNeg==11)    kind = 10; //Electron Combinatorial
836                 if(pdgCodePos==11 && pdgCodeNeg==11 && 
837                         (posDaughter->GetMother() == negDaughter->GetMother() && posDaughter->GetMother() ==-1))kind = 15; //direct Electron Combinatorial
838                                 
839                 if(pdgCodePos==211 && pdgCodeNeg==211) kind = 11; //Pion Combinatorial
840                 if((pdgCodePos==211 && pdgCodeNeg==2212) ||(pdgCodePos==2212 && pdgCodeNeg==211))       kind = 12; //Pion, Proton Combinatorics
841                 if((pdgCodePos==211 && pdgCodeNeg==11) ||(pdgCodePos==11 && pdgCodeNeg==211)) kind = 13; //Pion, Electron Combinatorics
842                 if(pdgCodePos==321 || pdgCodeNeg==321) kind = 14; //Kaon combinatorics
843         }else{          
844  
845                 pdgCodePos=posDaughter->GetPdgCode();
846                 pdgCodeNeg=negDaughter->GetPdgCode();
847
848                 if ( Photon->GetPdgCode()) 
849                         pdgCode = Photon->GetPdgCode(); 
850                 if(TMath::Abs(pdgCodePos)!=11 || TMath::Abs(pdgCodeNeg)!=11) kind = 2; // true from hadronic decays
851                 else if ( !(pdgCodeNeg==pdgCodePos)){
852                         if(pdgCode == 111) kind = 3; // pi0 Dalitz
853                         else if (pdgCode == 221) kind = 4; // eta Dalitz
854                         else if (!(negDaughter->GetMCProcessCode() != 5 || posDaughter->GetMCProcessCode() !=5)){
855                                 if(pdgCode == 22 && Photon->IsPrimary()){
856                                         kind = 0; // primary photons
857                                 } else if (pdgCode == 22){
858                                         kind = 5; //secondary photons
859                                 }
860                         }
861                 }
862         }
863
864         return kind;
865
866
867 }
868
869 //________________________________________________________________________
870 void AliAnalysisTaskConversionQA::RelabelAODPhotonCandidates(Bool_t mode){
871
872    // Relabeling For AOD Event
873    // ESDiD -> AODiD
874    // MCLabel -> AODMCLabel
875    
876    if(mode){
877       fMCStackPos = new Int_t[fConversionGammas->GetEntries()];
878       fMCStackNeg = new Int_t[fConversionGammas->GetEntries()];
879    }
880    
881    for(Int_t iGamma = 0;iGamma<fConversionGammas->GetEntries();iGamma++){
882       AliAODConversionPhoton* PhotonCandidate = (AliAODConversionPhoton*) fConversionGammas->At(iGamma);
883       if(!PhotonCandidate) continue;
884       if(!mode){// Back to ESD Labels
885          PhotonCandidate->SetMCLabelPositive(fMCStackPos[iGamma]);
886          PhotonCandidate->SetMCLabelNegative(fMCStackNeg[iGamma]);
887          //PhotonCandidate->IsAODMCLabel(kFALSE);
888          continue;
889       }
890       fMCStackPos[iGamma] =  PhotonCandidate->GetMCLabelPositive();
891       fMCStackNeg[iGamma] =  PhotonCandidate->GetMCLabelNegative();
892
893       Bool_t AODLabelPos = kFALSE;
894       Bool_t AODLabelNeg = kFALSE;
895
896       for(Int_t i = 0; i<fInputEvent->GetNumberOfTracks();i++){
897          AliAODTrack *tempDaughter = static_cast<AliAODTrack*>(fInputEvent->GetTrack(i));
898          if(!AODLabelPos){
899             if( tempDaughter->GetID() == PhotonCandidate->GetTrackLabelPositive() ){
900                PhotonCandidate->SetMCLabelPositive(abs(tempDaughter->GetLabel()));
901                AODLabelPos = kTRUE;
902             }
903          }
904          if(!AODLabelNeg){
905             if( tempDaughter->GetID() == PhotonCandidate->GetTrackLabelNegative()){
906                PhotonCandidate->SetMCLabelNegative(abs(tempDaughter->GetLabel()));
907                AODLabelNeg = kTRUE;
908             }
909          }
910          if(AODLabelNeg && AODLabelPos){
911             break;
912          }
913       } // Both ESD Tracks have AOD Tracks with Positive IDs
914       if(!AODLabelPos || !AODLabelNeg){
915          for(Int_t i = 0; i<fInputEvent->GetNumberOfTracks();i++){
916             AliAODTrack *tempDaughter = static_cast<AliAODTrack*>(fInputEvent->GetTrack(i));
917             if(tempDaughter->GetID()<0){
918                if(!AODLabelPos){
919                   if( (abs(tempDaughter->GetID())-1) == PhotonCandidate->GetTrackLabelPositive()){
920                      PhotonCandidate->SetMCLabelPositive(abs(tempDaughter->GetLabel()));
921                      AODLabelPos = kTRUE;
922                   }
923                }
924                if(!AODLabelNeg){
925                   if( (abs(tempDaughter->GetID())-1) == PhotonCandidate->GetTrackLabelNegative()){
926                      PhotonCandidate->SetMCLabelNegative(abs(tempDaughter->GetLabel()));
927                      AODLabelNeg = kTRUE;
928                   }
929                }
930             }
931             if(AODLabelNeg && AODLabelPos){
932                break;
933             }
934          }
935          if(!AODLabelPos || !AODLabelNeg){
936             cout<<"WARNING!!! AOD TRACKS NOT FOUND FOR"<<endl;
937          }
938       }
939    }
940    
941    if(!mode){
942       delete[] fMCStackPos;
943       delete[] fMCStackNeg;
944    }
945 }
946
947 void AliAnalysisTaskConversionQA::SetLogBinningXTH2(TH2* histoRebin){
948    TAxis *axisafter = histoRebin->GetXaxis(); 
949    Int_t bins = axisafter->GetNbins();
950    Double_t from = axisafter->GetXmin();
951    Double_t to = axisafter->GetXmax();
952    Double_t *newbins = new Double_t[bins+1];
953    newbins[0] = from;
954    Double_t factor = TMath::Power(to/from, 1./bins);
955    for(Int_t i=1; i<=bins; ++i) newbins[i] = factor * newbins[i-1];
956    axisafter->Set(bins, newbins);
957    delete [] newbins;
958
959 }
960
961 //________________________________________________________________________
962 void AliAnalysisTaskConversionQA::Terminate(Option_t *)
963 {
964
965 }