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