]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAnalysisTaskConversionQA.cxx
sync with GSI svn
[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
31 class iostream;
32
33 using namespace std;
34
35 ClassImp(AliAnalysisTaskConversionQA)
36
37 //________________________________________________________________________
38 AliAnalysisTaskConversionQA::AliAnalysisTaskConversionQA(const char *name) : AliAnalysisTaskSE(name),
39    fV0Reader(NULL),
40    fConversionGammas(NULL),
41    fConversionCuts(NULL),
42    fInputEvent(NULL),
43    fNumberOfESDTracks(0),
44    fMCEvent(NULL),
45    fMCStack(NULL),
46    fStreamQA(NULL),
47    fIsHeavyIon(kFALSE),
48    ffillTree(kFALSE),
49    ffillHistograms(kFALSE),
50    fOutputList(NULL),
51    fESDList(NULL),
52    hVertexZ(NULL),
53    hNGoodESDTracks(NULL),
54    hNV0Tracks(NULL),
55    hNContributorsVertex(NULL),
56    hITSClusterPhi(NULL),
57    hGammaPt(NULL),
58    hGammaPhi(NULL),
59    hGammaEta(NULL),
60    hGammaChi2perNDF(NULL),
61    hGammaPsiPair(NULL),
62    hGammaQt(NULL),
63    hGammaCosinePointingAngle(NULL),
64    hGammaXY(NULL),
65    hGammaZR(NULL),
66    hElecPt(NULL),
67    hElecEta(NULL),
68    hElecPhi(NULL),
69    hElecNfindableClsTPC(NULL),
70    hPosiNfindableClsTPC(NULL),
71    fTrueList(NULL),
72    hTrueResoulutionR(NULL),
73    hTrueResoulutionZ(NULL),
74    hTrueResoulutionPhi(NULL),
75    hTrueGammaPt(NULL),
76    hTrueGammaPhi(NULL),
77    hTrueGammaEta(NULL),
78    hTrueGammaMass(NULL),
79    hTrueGammaChi2perNDF(NULL),
80    hTrueGammaPsiPair(NULL),
81    hTrueGammaQt(NULL),
82    hTrueGammaCosinePointingAngle(NULL),
83    hTrueGammaXY(NULL),
84    hTrueGammaZR(NULL),
85    hTrueElecPt(NULL),
86    hTrueElecEta(NULL),
87    hTrueElecPhi(NULL),
88    hTrueElecNfindableClsTPC(NULL),
89    hTruePosiNfindableClsTPC(NULL)
90 {
91    // Default constructor
92
93    DefineInput(0, TChain::Class());
94    DefineOutput(1, TList::Class());
95 }
96
97 //________________________________________________________________________
98 AliAnalysisTaskConversionQA::~AliAnalysisTaskConversionQA()
99 {
100    // default deconstructor
101    if(fStreamQA){
102       delete fStreamQA;
103       fStreamQA = 0x0;
104    }
105 }
106 //________________________________________________________________________
107 void AliAnalysisTaskConversionQA::UserCreateOutputObjects()
108 {
109    // Create User Output Objects
110
111    if(fOutputList != NULL){
112       delete fOutputList;
113       fOutputList = NULL;
114    }
115    if(fOutputList == NULL){
116       fOutputList = new TList();
117       fOutputList->SetOwner(kTRUE);
118    }
119    
120    if(ffillHistograms){
121          
122       fESDList = new TList();
123       fESDList->SetOwner(kTRUE);
124       fESDList->SetName("ESD QA");
125       fOutputList->Add(fESDList);
126
127       hVertexZ = new TH1F("Vertex_Z","Vertex_Z",300,-15,15);
128       fESDList->Add(hVertexZ);
129       hNContributorsVertex = new TH1I("ContrVertex_Z","ContrVertex_Z",3000,0,3000);
130       fESDList->Add(hNContributorsVertex);
131       if(fIsHeavyIon) hNGoodESDTracks = new TH1I("GoodESDTracks","GoodESDTracks",3000,0,3000);
132       else hNGoodESDTracks = new TH1I("GoodESDTracks","GoodESDTracks",200,0,200);
133       fESDList->Add(hNGoodESDTracks);
134       if(fIsHeavyIon) hNV0Tracks = new TH1I("V0 Multiplicity","V0 Multiplicity",30000,0,30000);
135       else hNV0Tracks = new TH1I("V0 Multiplicity","V0 Multiplicity",2000,0,2000);
136       fESDList->Add(hNV0Tracks);
137
138       hITSClusterPhi = new TH2F("ITSClusterPhi","hITSClusterPhi",72,0,2*TMath::Pi(),7,0,7);
139       fESDList->Add(hITSClusterPhi);
140       hGammaPt = new TH1F("Gamma_Pt","Gamma_Pt",250,0,25);
141       fESDList->Add(hGammaPt);
142       hGammaPhi = new TH1F("Gamma_Phi","Gamma_Phi",360,0,2*TMath::Pi());
143       fESDList->Add(hGammaPhi);
144       hGammaEta = new TH1F("Gamma_Eta","Gamma_Eta",600,-1.5,1.5);
145       fESDList->Add(hGammaEta);
146       hGammaChi2perNDF = new TH1F("Gamma_Chi2perNDF","Gamma_Chi2perNDF",500,0,100);
147       fESDList->Add(hGammaChi2perNDF);
148       hGammaPsiPair = new TH1F("Gamma_PsiPair","Gamma_PsiPair",500,0,2);
149       fESDList->Add(hGammaPsiPair);
150       hGammaQt = new TH1F("Gamma_Qt","Gamma_Qt",400,0,0.1);
151       fESDList->Add(hGammaQt);
152       hGammaCosinePointingAngle = new TH1F("Gamma_CosinePointingAngle","Gamma_CosinePointingAngle",900,0.7,1.);
153       fESDList->Add(hGammaCosinePointingAngle);
154       hGammaXY = new TH2F("Gamma_ConversionPoint_XY","Gamma_ConversionPoint_XY",960,-120,120,960,-120,120);
155       fESDList->Add(hGammaXY);
156       hGammaZR= new TH2F("Gamma_ConversionPoint_ZR","Gamma_ConversionPoint_ZR",1200,-150,150,480,0,120);
157       fESDList->Add(hGammaZR);
158
159       hElecPt = new TH2F("Electron_Positron_Pt","Electron_Positron_Pt",250,0,25,250,0,25);
160       fESDList->Add(hElecPt);
161       hElecEta = new TH2F("Electron_Positron_Eta","Electron_Positron_Eta",600,-1.5,1.5,600,-1.5,1.5);
162       fESDList->Add(hElecEta);
163       hElecPhi = new TH2F("Electron_Positron_Phi","Electron_Positron_Phi",360,0,2*TMath::Pi(),360,0,2*TMath::Pi());
164       fESDList->Add(hElecPhi);
165       hElecNfindableClsTPC = new TH1F("Electron_findableClusterTPC","Electron_findableClusterTPC",100,0,1);
166       fESDList->Add(hElecNfindableClsTPC);
167       hPosiNfindableClsTPC = new TH1F("Positron_findableClusterTPC","Positron_findableClusterTPC",100,0,1);
168       fESDList->Add(hPosiNfindableClsTPC);
169
170       if(MCEvent()){
171          fTrueList = new TList();
172          fTrueList->SetOwner(kTRUE);
173          fTrueList->SetName("True QA");
174          fOutputList->Add(fTrueList);
175
176          hTrueResoulutionR = new TH2F("True_ConversionPointResolution_R","True_ConversionPointResolution_R",240,0,120,200,-20,20);
177          fTrueList->Add(hTrueResoulutionR);
178          hTrueResoulutionZ = new TH2F("True_ConversionPointResolution_Z","True_ConversionPointResolution_Z",480,-120,120,200,-20,20);
179          fTrueList->Add(hTrueResoulutionZ);
180          hTrueResoulutionPhi = new TH2F("True_ConversionPointResolution_Phi","True_ConversionPointResolution_Phi",360,0,2*TMath::Pi(),200,-20,20);
181          fTrueList->Add(hTrueResoulutionPhi);
182
183          hTrueGammaPt = new TH1F("True_Gamma_Pt","True_Gamma_Pt",250,0,25);
184          fTrueList->Add(hTrueGammaPt);
185          hTrueGammaPhi = new TH1F("True_Gamma_Phi","True_Gamma_Phi",360,0,2*TMath::Pi());
186          fTrueList->Add(hTrueGammaPhi);
187          hTrueGammaEta = new TH1F("True_Gamma_Eta","True_Gamma_Eta",600,-1.5,1.5);
188          fTrueList->Add(hTrueGammaEta);
189          hTrueGammaMass = new TH1F("True_Gamma_Mass","True_Gamma_Mass",1000,0,0.3);
190          fTrueList->Add(hTrueGammaMass);
191          hTrueGammaChi2perNDF = new TH1F("True_Gamma_Chi2perNDF","True_Gamma_Chi2perNDF",500,0,100);
192          fTrueList->Add(hTrueGammaChi2perNDF);
193          hTrueGammaPsiPair = new TH1F("True_Gamma_PsiPair","True_Gamma_PsiPair",500,0,2);
194          fTrueList->Add(hTrueGammaPsiPair);
195          hTrueGammaQt = new TH1F("True_Gamma_Qt","True_Gamma_Qt",400,0,0.1);
196          fTrueList->Add(hTrueGammaQt);
197          hTrueGammaCosinePointingAngle = new TH1F("True_Gamma_CosinePointingAngle","True_Gamma_CosinePointingAngle",900,0.7,1.);
198          fTrueList->Add(hTrueGammaCosinePointingAngle);
199          hTrueGammaXY = new TH2F("True_Gamma_ConversionPoint_XY","True_Gamma_ConversionPoint_XY",960,-120,120,960,-120,120);
200          fTrueList->Add(hTrueGammaXY);
201          hTrueGammaZR= new TH2F("TrueGamma_ConversionPoint_ZR","TrueGamma_ConversionPoint_ZR",1200,-150,150,480,0,120);
202          fTrueList->Add(hTrueGammaZR);
203
204          hTrueElecPt = new TH2F("True_Electron_Positron_Pt","True_Electron_Positron_Pt",250,0,25,250,0,25);
205          fTrueList->Add(hTrueElecPt);
206          hTrueElecEta = new TH2F("True_Electron_Positron_Eta","True_Electron_Positron_Eta",600,-1.5,1.5,600,-1.5,1.5);
207          fTrueList->Add(hTrueElecEta);
208          hTrueElecPhi = new TH2F("True_Electron_Positron_Phi","True_Electron_Positron_Phi",360,0,2*TMath::Pi(),360,0,2*TMath::Pi());
209          fTrueList->Add(hTrueElecPhi);
210          hTrueElecNfindableClsTPC = new TH1F("True_Electron_findableClusterTPC","True_Electron_findableClusterTPC",100,0,1);
211          fTrueList->Add(hTrueElecNfindableClsTPC);
212          hTruePosiNfindableClsTPC = new TH1F("True_Positron_findableClusterTPC","True_Positron_findableClusterTPC",100,0,1);
213          fTrueList->Add(hTruePosiNfindableClsTPC);
214       }
215       if(fConversionCuts->GetCutHistograms()){
216          fOutputList->Add(fConversionCuts->GetCutHistograms());
217       }
218    }
219    
220    if(ffillTree){
221       TString cutnumber = fConversionCuts->GetCutNumber();
222       fStreamQA = new TTreeSRedirector(Form("GammaConvV1_QATree_%s.root",cutnumber.Data()),"recreate");
223    }
224
225    fV0Reader=(AliV0ReaderV1*)AliAnalysisManager::GetAnalysisManager()->GetTask("V0ReaderV1");
226    
227    PostData(1, fOutputList);
228 }
229 //_____________________________________________________________________________
230 Bool_t AliAnalysisTaskConversionQA::Notify()
231 {
232    if(!fConversionCuts->GetDoEtaShift()) return kTRUE;; // No Eta Shift requested, continue
233       
234    if(fConversionCuts->GetEtaShift() == 0.0){ // Eta Shift requested but not set, get shift automatically
235       fConversionCuts->GetCorrectEtaShiftFromPeriod(fV0Reader->GetPeriodName());
236       fConversionCuts->DoEtaShift(kFALSE); // Eta Shift Set, make sure that it is called only once
237       return kTRUE;
238    }
239    else{
240       printf(" Gamma Conversion QA Task %s :: Eta Shift Manually Set to %f \n\n",
241              (fConversionCuts->GetCutNumber()).Data(),fConversionCuts->GetEtaShift());
242       fConversionCuts->DoEtaShift(kFALSE); // Eta Shift Set, make sure that it is called only once
243    }
244    
245    return kTRUE;
246 }
247 //________________________________________________________________________
248 void AliAnalysisTaskConversionQA::UserExec(Option_t *){
249
250
251
252    Int_t eventQuality = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetEventQuality();
253    if(eventQuality != 0){// Event Not Accepted
254       return;
255    }
256    fInputEvent = InputEvent();
257    fMCEvent = MCEvent();
258    if(fMCEvent) fMCStack = fMCEvent->Stack();
259
260    Int_t eventNotAccepted =
261       fConversionCuts->IsEventAcceptedByConversionCut(fV0Reader->GetConversionCuts(),fInputEvent,fMCEvent,fIsHeavyIon);
262    if(eventNotAccepted) return; // Check Centrality, PileUp, SDD and V0AND --> Not Accepted => eventQuality = 1
263
264    fConversionGammas=fV0Reader->GetReconstructedGammas();
265
266    if(fMCEvent){
267       if(fConversionCuts->GetSignalRejection() != 0){
268          fConversionCuts->GetNotRejectedParticles(fConversionCuts->GetSignalRejection(),
269                                                   fConversionCuts->GetAcceptedHeader(),
270                                                   fMCEvent);
271       }
272    }
273
274    if(ffillHistograms){
275       CountESDTracks();
276       hVertexZ->Fill(fInputEvent->GetPrimaryVertex()->GetZ());
277       hNContributorsVertex->Fill(fConversionCuts->GetNumberOfContributorsVtx(fInputEvent));
278       hNGoodESDTracks->Fill(fNumberOfESDTracks);
279       hNV0Tracks->Fill(fInputEvent->GetVZEROData()->GetMTotV0A()+fInputEvent->GetVZEROData()->GetMTotV0C());
280    }
281
282    for(Int_t firstGammaIndex=0;firstGammaIndex<fConversionGammas->GetEntriesFast();firstGammaIndex++){
283       AliAODConversionPhoton *gamma=dynamic_cast<AliAODConversionPhoton*>(fConversionGammas->At(firstGammaIndex));
284       if (gamma==NULL) continue;
285       if(fMCEvent && fConversionCuts->GetSignalRejection() != 0){
286          if(!fConversionCuts->IsParticleFromBGEvent(gamma->GetMCLabelPositive(), fMCStack))
287             continue;
288          if(!fConversionCuts->IsParticleFromBGEvent(gamma->GetMCLabelNegative(), fMCStack))
289             continue;
290       }
291       if(!fConversionCuts->PhotonIsSelected(gamma,fInputEvent)){
292          continue;
293       }
294
295       if(ffillTree) ProcessQATree(gamma);
296       if(ffillHistograms) ProcessQA(gamma);
297    }
298
299    PostData(1, fOutputList);
300 }
301
302
303 ///________________________________________________________________________
304 void AliAnalysisTaskConversionQA::ProcessQATree(AliAODConversionPhoton *gamma){
305
306    // Fill Histograms for QA and MC
307    AliESDEvent* event = (AliESDEvent*) InputEvent();
308    AliPIDResponse* pidResonse = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetPIDResponse();
309
310    Float_t gammaPt = gamma->GetPhotonPt();
311    Float_t gammaPhi = gamma->GetPhotonPhi();
312    Float_t gammaTheta = gamma->Theta();
313    Float_t gammaChi2NDF = gamma->GetChi2perNDF();
314    Float_t gammaQt = gamma->GetArmenterosQt();
315    Float_t gammaAlpha = gamma->GetArmenterosAlpha();
316    Float_t gammaPsiPair = gamma->GetPsiPair();
317    Float_t gammaCosPointing = fConversionCuts->GetCosineOfPointingAngle(gamma,event);
318    TVectorF conversionPoint(3);
319    conversionPoint(0) = gamma->GetConversionX();
320    conversionPoint(1) = gamma->GetConversionY();
321    conversionPoint(2) = gamma->GetConversionZ();
322    TVectorF daughterProp(20);
323    AliESDtrack * negTrack = fConversionCuts->GetESDTrack(event, gamma->GetTrackLabelNegative());
324    AliESDtrack * posTrack = fConversionCuts->GetESDTrack(event, gamma->GetTrackLabelPositive());
325
326    if(!negTrack||!posTrack)return;
327
328    Bool_t isTruePhoton = kFALSE;
329    if(fMCEvent){
330       if(IsTruePhoton(gamma)) isTruePhoton = kTRUE;
331    }
332
333    daughterProp(0) = posTrack->Pt();
334    daughterProp(7) = negTrack->Pt();
335    daughterProp(1) = posTrack->Theta();
336    daughterProp(8) = negTrack->Theta();
337    Double32_t signalPos[4] = {0,0,0,0};
338    Char_t nclPos[3];
339    Char_t nrowsPos[3];
340    if (posTrack->GetTPCdEdxInfo()) {
341       posTrack->GetTPCdEdxInfo()->GetTPCSignalRegionInfo(signalPos,nclPos,nrowsPos);
342       daughterProp(2) = signalPos[0];
343       daughterProp(14) = signalPos[1];
344       daughterProp(16) = signalPos[2];
345    } else {
346       daughterProp(2) = posTrack->GetTPCsignal();
347       daughterProp(14) = 0;
348       daughterProp(16) = 0;
349    }
350    Double32_t signalNeg[4] = {0,0,0,0};
351    Char_t nclNeg[3];
352    Char_t nrowsNeg[3];
353    if (negTrack->GetTPCdEdxInfo()) {
354       negTrack->GetTPCdEdxInfo()->GetTPCSignalRegionInfo(signalNeg,nclNeg,nrowsNeg);
355       daughterProp(9) = signalNeg[0];
356       daughterProp(15) = signalNeg[1];
357       daughterProp(17) = signalNeg[2];
358    } else {
359       daughterProp(9) = negTrack->GetTPCsignal();
360       daughterProp(15) = 0;
361       daughterProp(17) = 0;
362    }
363    daughterProp(3) = pidResonse->NumberOfSigmasTPC(posTrack,AliPID::kElectron);
364    daughterProp(10) = pidResonse->NumberOfSigmasTPC(negTrack,AliPID::kElectron);
365    if((posTrack->GetStatus() & AliESDtrack::kTOFpid) && !(posTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
366       Double_t t0pos = pidResonse->GetTOFResponse().GetStartTime(posTrack->P());
367       Double_t timesPos[5];
368       posTrack->GetIntegratedTimes(timesPos);
369       Double_t TOFsignalPos =   posTrack->GetTOFsignal();
370       Double_t dTpos = TOFsignalPos - t0pos - timesPos[0];
371       daughterProp(4) = dTpos;
372       daughterProp(5) = pidResonse->NumberOfSigmasTOF(posTrack, AliPID::kElectron);
373    } else {
374       daughterProp(4) = 20000;
375       daughterProp(5) = -20;
376    }
377    if((negTrack->GetStatus() & AliESDtrack::kTOFpid) && !(negTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
378       Double_t t0neg = pidResonse->GetTOFResponse().GetStartTime(negTrack->P());
379       Double_t timesNeg[5];
380       negTrack->GetIntegratedTimes(timesNeg);
381       Double_t TOFsignalNeg =   negTrack->GetTOFsignal();
382       Double_t dTneg = TOFsignalNeg - t0neg - timesNeg[0];
383       daughterProp(11) = dTneg;
384       daughterProp(12) = pidResonse->NumberOfSigmasTOF(negTrack, AliPID::kElectron);
385    } else {
386       daughterProp(11) = 20000;
387       daughterProp(12) = -20;
388    }
389
390    daughterProp(6) = (Float_t)posTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius()));
391    daughterProp(18) = posTrack->GetNcls(1);
392    daughterProp(13) = (Float_t)negTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius()));
393    daughterProp(19) = negTrack->GetNcls(1);
394
395    if (fStreamQA){
396       if(fMCEvent){
397          (*fStreamQA)<<"PhotonQA"
398                      << "pt=" << gammaPt
399                      << "phi=" << gammaPhi
400                      << "theta=" << gammaTheta
401                      << "chi2ndf=" << gammaChi2NDF
402                      << "qt="<< gammaQt
403                      << "alpha=" << gammaAlpha
404                      << "psipair=" << gammaPsiPair
405                      << "cosPoint=" << gammaCosPointing
406                      << "TruePhoton=" << isTruePhoton
407                      << "conversionPoint=" << &conversionPoint
408                      << "daugtherProp.=" << &daughterProp
409                      << "\n";
410       }
411       else{
412          (*fStreamQA)<<"PhotonQA"
413                      << "pt=" << gammaPt
414                      << "phi=" << gammaPhi
415                      << "theta=" << gammaTheta
416                      << "chi2ndf=" << gammaChi2NDF
417                      << "qt="<< gammaQt
418                      << "alpha=" << gammaAlpha
419                      << "psipair=" << gammaPsiPair
420                      << "cosPoint=" << gammaCosPointing
421                      << "conversionPoint=" << &conversionPoint
422                      << "daugtherProp.=" << &daughterProp
423                      << "\n";
424       }
425    }
426 }
427
428 //_____________________________________________________________________________________________________
429 void AliAnalysisTaskConversionQA::ProcessQA(AliAODConversionPhoton *gamma){
430
431
432    // Fill Histograms for QA and MC
433    hGammaPt->Fill(gamma->GetPhotonPt());
434    hGammaPhi->Fill(gamma->GetPhotonPhi());
435    hGammaEta->Fill(gamma->Eta());
436    hGammaChi2perNDF->Fill(gamma->GetChi2perNDF());
437    hGammaPsiPair->Fill(gamma->GetPsiPair());
438    hGammaQt->Fill(gamma->GetArmenterosQt());
439    hGammaCosinePointingAngle->Fill(fConversionCuts->GetCosineOfPointingAngle(gamma,fInputEvent));
440    hGammaXY->Fill(gamma->GetConversionX(),gamma->GetConversionY());
441    hGammaZR->Fill(gamma->GetConversionZ(),gamma->GetConversionRadius());
442
443    AliESDtrack * negTrack = fConversionCuts->GetESDTrack((AliESDEvent*)fInputEvent, gamma->GetTrackLabelNegative());
444    AliESDtrack * posTrack = fConversionCuts->GetESDTrack((AliESDEvent*)fInputEvent, gamma->GetTrackLabelPositive());
445    if(!negTrack||!posTrack)return;
446
447    Double_t negtrackPhi = negTrack->Phi();
448    Double_t postrackPhi = posTrack->Phi();
449    hITSClusterPhi->Fill(negtrackPhi,6);
450    hITSClusterPhi->Fill(postrackPhi,6);
451    for(Int_t itsLayer = 0; itsLayer<6;itsLayer++){
452       if(TESTBIT(negTrack->GetITSClusterMap(),itsLayer)){
453          hITSClusterPhi->Fill(negtrackPhi,itsLayer);
454       }
455       if(TESTBIT(posTrack->GetITSClusterMap(),itsLayer)){
456          hITSClusterPhi->Fill(postrackPhi,itsLayer);
457       }
458    }
459
460    hElecPt->Fill(negTrack->Pt(),posTrack->Pt());
461    hElecEta->Fill(negTrack->Eta(),posTrack->Eta());
462    hElecPhi->Fill(negTrack->Phi(),posTrack->Phi());
463
464    hElecNfindableClsTPC->Fill((Float_t)posTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius())));
465    hPosiNfindableClsTPC->Fill((Float_t)negTrack->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(gamma->GetConversionRadius())));
466
467    // hElecNfindableClsTPC->Fill((Float_t)negTrack->GetNcls(1)/(Float_t)negTrack->GetTPCNclsF());
468    // hPosiNfindableClsTPC->Fill((Float_t)posTrack->GetNcls(1)/(Float_t)posTrack->GetTPCNclsF());
469
470    if(fMCEvent) ProcessTrueQA(gamma,negTrack,posTrack);
471 }
472
473 //________________________________________________________________________
474 void AliAnalysisTaskConversionQA::ProcessTrueQA(AliAODConversionPhoton *TruePhotonCandidate, AliESDtrack *elec, AliESDtrack *posi)
475 {
476
477    if(!IsTruePhoton(TruePhotonCandidate)) return;
478
479    TParticle *negDaughter = TruePhotonCandidate->GetNegativeMCDaughter(fMCStack);
480    TParticle *mcPhoton = TruePhotonCandidate->GetMCParticle(fMCStack);
481
482    // True Photon
483    hTrueResoulutionR->Fill(TruePhotonCandidate->GetConversionRadius(),
484                            TruePhotonCandidate->GetConversionRadius()-negDaughter->R());
485    hTrueResoulutionZ->Fill(TruePhotonCandidate->GetConversionZ(),
486                            TruePhotonCandidate->GetConversionZ()-negDaughter->Vz());
487    hTrueResoulutionPhi->Fill(TruePhotonCandidate->Phi(),
488                              TruePhotonCandidate->Phi()-mcPhoton->Phi());
489    hTrueGammaPt->Fill(TruePhotonCandidate->Pt());
490    hTrueGammaPhi->Fill(TruePhotonCandidate->Phi());
491    hTrueGammaEta->Fill(TruePhotonCandidate->Eta());
492    hTrueGammaMass->Fill(TruePhotonCandidate->GetMass());
493    hTrueGammaChi2perNDF->Fill(TruePhotonCandidate->GetChi2perNDF());
494    hTrueGammaPsiPair->Fill(TruePhotonCandidate->GetPsiPair());
495    hTrueGammaQt->Fill(TruePhotonCandidate->GetArmenterosQt());
496    hTrueGammaCosinePointingAngle->Fill(fConversionCuts->GetCosineOfPointingAngle(TruePhotonCandidate,fInputEvent));
497    hTrueGammaXY->Fill(TruePhotonCandidate->GetConversionX(),TruePhotonCandidate->GetConversionY());
498    hTrueGammaZR->Fill(TruePhotonCandidate->GetConversionZ(),TruePhotonCandidate->GetConversionRadius());
499
500    hTrueElecPt->Fill(elec->Pt(),posi->Pt());
501    hTrueElecEta->Fill(elec->Eta(),posi->Eta());
502    hTrueElecPhi->Fill(elec->Phi(),posi->Phi());
503    hTrueElecNfindableClsTPC
504       ->Fill((Float_t)elec->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(TruePhotonCandidate->GetConversionRadius())));
505    hTruePosiNfindableClsTPC
506       ->Fill((Float_t)posi->GetTPCClusterInfo(2,0,fConversionCuts->GetFirstTPCRow(TruePhotonCandidate->GetConversionRadius())));
507    // hTrueElecNfindableClsTPC->Fill((Float_t)elec->GetNcls(1)/(Float_t)elec->GetTPCNclsF());
508    // hTruePosiNfindableClsTPC->Fill((Float_t)posi->GetNcls(1)/(Float_t)posi->GetTPCNclsF());
509
510 }
511 //________________________________________________________________________
512 void AliAnalysisTaskConversionQA::CountESDTracks(){
513
514    // Using standard function for setting Cuts
515    Bool_t selectPrimaries=kTRUE;
516    AliESDtrackCuts *EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
517    EsdTrackCuts->SetMaxDCAToVertexZ(2);
518    EsdTrackCuts->SetEtaRange(-0.8, 0.8);
519    EsdTrackCuts->SetPtRange(0.15);
520
521    fNumberOfESDTracks = 0;
522    for(Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++){
523       AliESDtrack* curTrack = (AliESDtrack*)fInputEvent->GetTrack(iTracks);
524       if(!curTrack) continue;
525       // if(fMCEvent && fConversionCuts->GetSignalRejection() != 0){
526       //    if(!fConversionCuts->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack)) continue;
527       // }
528       if(EsdTrackCuts->AcceptTrack(curTrack) ) fNumberOfESDTracks++;
529    }
530    delete EsdTrackCuts;
531    EsdTrackCuts=0x0;
532
533    return;
534 }
535 //________________________________________________________________________
536 Bool_t AliAnalysisTaskConversionQA::IsTruePhoton(AliAODConversionPhoton *TruePhotonCandidate)
537 {
538    TParticle *posDaughter = TruePhotonCandidate->GetPositiveMCDaughter(fMCStack);
539    TParticle *negDaughter = TruePhotonCandidate->GetNegativeMCDaughter(fMCStack);
540
541    if(posDaughter == NULL || negDaughter == NULL) return kFALSE; // One particle does not exist
542    Int_t pdgCode[2] = {abs(posDaughter->GetPdgCode()),abs(negDaughter->GetPdgCode())};
543    if(posDaughter->GetMother(0) != negDaughter->GetMother(0)) return kFALSE;
544    else if(posDaughter->GetMother(0) == -1) return kFALSE;
545
546    if(pdgCode[0]!=11 || pdgCode[1]!=11) return kFALSE; //One Particle is not electron
547    if(posDaughter->GetPdgCode()==negDaughter->GetPdgCode()) return kFALSE; // Same Charge
548    if(posDaughter->GetUniqueID() != 5 || negDaughter->GetUniqueID() !=5) return kFALSE;// check if the daughters come from a conversion
549
550    TParticle *Photon = TruePhotonCandidate->GetMCParticle(fMCStack);
551    if(Photon->GetPdgCode() != 22) return kFALSE; // Mother is no Photon
552
553    return kTRUE;
554 }
555
556 //________________________________________________________________________
557 void AliAnalysisTaskConversionQA::Terminate(Option_t *)
558 {
559
560 }