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