]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAnalysisTaskGammaConvV1.cxx
added new conversion task (gamma, gamma & dalitz task) + corresponding dependencies
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAnalysisTaskGammaConvV1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: Martin Wilde, 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 // Class used to do analysis on conversion pairs
19 //---------------------------------------------
20 ///////////////////////////////////////////////
21 #include "TChain.h"
22 #include "TTree.h"
23 #include "TH1F.h"
24 #include "TH2F.h"
25 #include "TH3F.h"
26 #include "THnSparse.h"
27 #include "TCanvas.h"
28 #include "TNtuple.h"
29 #include "AliAnalysisTask.h"
30 #include "AliAnalysisManager.h"
31 #include "AliESDEvent.h"
32 #include "AliESDInputHandler.h"
33 #include "AliMCEventHandler.h"
34 #include "AliMCEvent.h"
35 #include "AliMCParticle.h"
36 #include "AliCentrality.h"
37 #include "AliESDVZERO.h"
38 #include "AliESDpid.h"
39 #include "AliAnalysisTaskGammaConvV1.h"
40 #include "AliVParticle.h"
41 #include "AliESDtrackCuts.h"
42 #include "AliKFVertex.h"
43 #include "AliV0ReaderV1.h"
44 #include "AliGammaConversionAODBGHandler.h"
45
46 ClassImp(AliAnalysisTaskGammaConvV1)
47
48 //________________________________________________________________________
49 AliAnalysisTaskGammaConvV1::AliAnalysisTaskGammaConvV1(): AliAnalysisTaskSE(),
50    fV0Reader(NULL),
51    fBGHandler(NULL),
52    fESDEvent(NULL),
53    fMCEvent(NULL),
54    fMCStack(NULL),
55    fCutFolder(NULL),
56    fESDList(NULL),
57    fBackList(NULL),
58    fTrueList(NULL),
59    fMCList(NULL),
60    fOutputContainer(0),
61    fReaderGammas(NULL),
62    fGoodGammas(NULL),
63    fCutArray(NULL),
64    fConversionCuts(NULL),
65    hESDConvGammaPt(NULL),
66    hESDMotherInvMassPt(NULL),
67    sESDMotherInvMassPtZM(NULL),
68    hESDMotherBackInvMassPt(NULL),
69    sESDMotherBackInvMassPtZM(NULL),
70    hESDMotherInvMassEalpha(NULL),
71    hMCAllGammaPt(NULL),
72    hMCAllGammaPtAddedSig(NULL),
73    hMCDecayGammaPi0Pt(NULL),
74    hMCDecayGammaPi0PtAddedSig(NULL),
75    hMCDecayGammaRhoPt(NULL),
76    hMCDecayGammaRhoPtAddedSig(NULL),
77    hMCDecayGammaEtaPt(NULL),
78    hMCDecayGammaEtaPtAddedSig(NULL),
79    hMCDecayGammaOmegaPt(NULL),
80    hMCDecayGammaOmegaPtAddedSig(NULL),
81    hMCDecayGammaEtapPt(NULL),
82    hMCDecayGammaEtapPtAddedSig(NULL),
83    hMCDecayGammaPhiPt(NULL),
84    hMCDecayGammaPhiPtAddedSig(NULL),
85    hMCConvGammaPt(NULL),
86    hMCConvGammaPtAddedSig(NULL),
87    hMCPi0Pt(NULL),
88    hMCEtaPt(NULL),
89    hMCPi0InAccPt(NULL),
90    hMCEtaInAccPt(NULL),
91    hMCPi0PtAddedSig(NULL),
92    hMCEtaPtAddedSig(NULL),
93    hMCPi0InAccPtAddedSig(NULL),
94    hMCEtaInAccPtAddedSig(NULL),
95    hESDTrueMotherInvMassPt(NULL),
96    hESDTrueMotherInvMassPtAddedSig(NULL),
97    hESDTruePrimaryMotherInvMassMCPt(NULL),
98    hESDTruePrimaryMotherInvMassMCPtAddedSig(NULL),
99    hESDTruePrimaryPi0ESDPtMCPt(NULL),
100    hESDTruePrimaryPi0ESDPtMCPtAddedSig(NULL),
101    hESDTrueSecondaryMotherInvMassPt(NULL),
102    hESDTrueSecondaryMotherInvMassPtAddedSig(NULL),
103    hESDTrueSecondaryMotherFromK0sInvMassPt(NULL),
104    hESDTrueSecondaryMotherFromK0sInvMassPtAddedSig(NULL),
105    hESDTrueBckGGInvMassPt(NULL),
106    hESDTrueBckContInvMassPt(NULL),
107    hESDTrueMotherDalitzInvMassPt(NULL),
108    hESDTrueMotherDalitzInvMassPtAddedSig(NULL),
109    hESDTrueConvGammaPt(NULL),
110    hESDTrueConvGammaPtAddedSig(NULL),
111    hESDTrueElecCombPt(NULL),
112    hESDTruePionCombPt(NULL),
113    hESDTruePrimaryConvGammaPt(NULL),
114    hESDTruePrimaryConvGammaPtAddedSig(NULL),
115    hESDTruePrimaryConvGammaESDPtMCPt(NULL),
116    hESDTruePrimaryConvGammaESDPtMCPtAddedSig(NULL),
117    hESDTrueSecondaryConvGammaPt(NULL),
118    hESDTrueSecondaryConvGammaPtAddedSig(NULL),
119    hESDTrueSecondaryConvGammaFromK0sPt(NULL),
120    hESDTrueSecondaryConvGammaFromK0sPtAddedSig(NULL),
121    hESDTrueSecondaryConvGammaFromXFromK0sPt(NULL),
122    hESDTrueSecondaryConvGammaFromXFromK0sPtAddedSig(NULL),
123    hNEvents(NULL),
124    hNGoodESDTracks(NULL),
125    hNV0Tracks(NULL),
126    fRandom(0),
127    fUnsmearedPx(NULL),
128    fUnsmearedPy(NULL),
129    fUnsmearedPz(NULL),
130    fUnsmearedE(NULL),
131    fIsHeavyIon(kFALSE),
132    fDoMesonAnalysis(kTRUE)
133
134 {
135    // default Constructor
136    DefineInput(0, TChain::Class());
137    DefineOutput(1, TList::Class());
138 }
139
140 //________________________________________________________________________
141 AliAnalysisTaskGammaConvV1::AliAnalysisTaskGammaConvV1(const char *name):
142    AliAnalysisTaskSE(name),
143    fV0Reader(NULL),
144    fBGHandler(NULL),
145    fESDEvent(NULL),
146    fMCEvent(NULL),
147    fMCStack(NULL),
148    fCutFolder(NULL),
149    fESDList(NULL),
150    fBackList(NULL),
151    fTrueList(NULL),
152    fMCList(NULL),
153    fOutputContainer(0),
154    fReaderGammas(NULL),
155    fGoodGammas(NULL),
156    fCutArray(NULL),
157    fConversionCuts(NULL),
158    hESDConvGammaPt(NULL),
159    hESDMotherInvMassPt(NULL),
160    sESDMotherInvMassPtZM(NULL),
161    hESDMotherBackInvMassPt(NULL),
162    sESDMotherBackInvMassPtZM(NULL),
163    hESDMotherInvMassEalpha(NULL),
164    hMCAllGammaPt(NULL),
165    hMCAllGammaPtAddedSig(NULL),
166    hMCDecayGammaPi0Pt(NULL),
167    hMCDecayGammaPi0PtAddedSig(NULL),
168    hMCDecayGammaRhoPt(NULL),
169    hMCDecayGammaRhoPtAddedSig(NULL),
170    hMCDecayGammaEtaPt(NULL),
171    hMCDecayGammaEtaPtAddedSig(NULL),
172    hMCDecayGammaOmegaPt(NULL),
173    hMCDecayGammaOmegaPtAddedSig(NULL),
174    hMCDecayGammaEtapPt(NULL),
175    hMCDecayGammaEtapPtAddedSig(NULL),
176    hMCDecayGammaPhiPt(NULL),
177    hMCDecayGammaPhiPtAddedSig(NULL),
178    hMCConvGammaPt(NULL),
179    hMCConvGammaPtAddedSig(NULL),
180    hMCPi0Pt(NULL),
181    hMCEtaPt(NULL),
182    hMCPi0InAccPt(NULL),
183    hMCEtaInAccPt(NULL),
184    hMCPi0PtAddedSig(NULL),
185    hMCEtaPtAddedSig(NULL),
186    hMCPi0InAccPtAddedSig(NULL),
187    hMCEtaInAccPtAddedSig(NULL),
188    hESDTrueMotherInvMassPt(NULL),
189    hESDTrueMotherInvMassPtAddedSig(NULL),
190    hESDTruePrimaryMotherInvMassMCPt(NULL),
191    hESDTruePrimaryMotherInvMassMCPtAddedSig(NULL),
192    hESDTruePrimaryPi0ESDPtMCPt(NULL),
193    hESDTruePrimaryPi0ESDPtMCPtAddedSig(NULL),
194    hESDTrueSecondaryMotherInvMassPt(NULL),
195    hESDTrueSecondaryMotherInvMassPtAddedSig(NULL),
196    hESDTrueSecondaryMotherFromK0sInvMassPt(NULL),
197    hESDTrueSecondaryMotherFromK0sInvMassPtAddedSig(NULL),
198    hESDTrueBckGGInvMassPt(NULL),
199    hESDTrueBckContInvMassPt(NULL),
200    hESDTrueMotherDalitzInvMassPt(NULL),
201    hESDTrueMotherDalitzInvMassPtAddedSig(NULL),
202    hESDTrueConvGammaPt(NULL),
203    hESDTrueConvGammaPtAddedSig(NULL),
204    hESDTrueElecCombPt(NULL),
205    hESDTruePionCombPt(NULL),
206    hESDTruePrimaryConvGammaPt(NULL),
207    hESDTruePrimaryConvGammaPtAddedSig(NULL),
208    hESDTruePrimaryConvGammaESDPtMCPt(NULL),
209    hESDTruePrimaryConvGammaESDPtMCPtAddedSig(NULL),
210    hESDTrueSecondaryConvGammaPt(NULL),
211    hESDTrueSecondaryConvGammaPtAddedSig(NULL),
212    hESDTrueSecondaryConvGammaFromK0sPt(NULL),
213    hESDTrueSecondaryConvGammaFromK0sPtAddedSig(NULL),
214    hESDTrueSecondaryConvGammaFromXFromK0sPt(NULL),
215    hESDTrueSecondaryConvGammaFromXFromK0sPtAddedSig(NULL),
216    hNEvents(NULL),
217    hNGoodESDTracks(NULL),
218    hNV0Tracks(NULL),
219    fRandom(0),
220    fUnsmearedPx(NULL),
221    fUnsmearedPy(NULL),
222    fUnsmearedPz(NULL),
223    fUnsmearedE(NULL),
224    fIsHeavyIon(kFALSE),
225    fDoMesonAnalysis(kTRUE)
226
227 {
228    // Define input and output slots here
229    DefineInput(0, TChain::Class());
230    DefineOutput(1, TList::Class());
231 }
232
233 AliAnalysisTaskGammaConvV1::~AliAnalysisTaskGammaConvV1()
234 {
235    if(fGoodGammas){
236       delete fGoodGammas;
237       fGoodGammas = 0x0;
238    }
239    if(fBGHandler){
240       delete[] fBGHandler;
241       fBGHandler = 0x0;
242    }
243 }
244 //___________________________________________________________
245 void AliAnalysisTaskGammaConvV1::InitBack(){
246
247    Double_t *zBinLimitsArray = new Double_t[9];
248    zBinLimitsArray[0] = -50.00;
249    zBinLimitsArray[1] = -3.375;
250    zBinLimitsArray[2] = -1.605;
251    zBinLimitsArray[3] = -0.225;
252    zBinLimitsArray[4] = 1.065;
253    zBinLimitsArray[5] = 2.445;
254    zBinLimitsArray[6] = 4.245;
255    zBinLimitsArray[7] = 50.00;
256    zBinLimitsArray[8] = 1000.00;
257
258    Double_t *multiplicityBinLimitsArrayTracks = new Double_t[6];
259    multiplicityBinLimitsArrayTracks[0] = 0;
260    multiplicityBinLimitsArrayTracks[1] = 8.5;
261    multiplicityBinLimitsArrayTracks[2] = 16.5;
262    multiplicityBinLimitsArrayTracks[3] = 27.5;
263    multiplicityBinLimitsArrayTracks[4] = 41.5;
264    multiplicityBinLimitsArrayTracks[5] = 200.;
265    if(fIsHeavyIon){
266       multiplicityBinLimitsArrayTracks[0] = 0;
267       multiplicityBinLimitsArrayTracks[1] = 200.;
268       multiplicityBinLimitsArrayTracks[2] = 500.;
269       multiplicityBinLimitsArrayTracks[3] = 1000.;
270       multiplicityBinLimitsArrayTracks[4] = 1500.;
271       multiplicityBinLimitsArrayTracks[5] = 5000.;
272    }
273
274    Double_t *multiplicityBinLimitsArrayV0s = new Double_t[5];
275    multiplicityBinLimitsArrayV0s[0] = 2;
276    multiplicityBinLimitsArrayV0s[1] = 3;
277    multiplicityBinLimitsArrayV0s[2] = 4;
278    multiplicityBinLimitsArrayV0s[3] = 5;
279    multiplicityBinLimitsArrayV0s[4] = 9999;
280    if(fIsHeavyIon){
281       multiplicityBinLimitsArrayV0s[0] = 2;
282       multiplicityBinLimitsArrayV0s[1] = 10;
283       multiplicityBinLimitsArrayV0s[2] = 30;
284       multiplicityBinLimitsArrayV0s[3] = 50;
285       multiplicityBinLimitsArrayV0s[4] = 9999;
286    }
287
288    const Int_t nDim = 4;
289    Int_t nBins[nDim] = {1000,250,8,5};
290    Double_t xMin[nDim] = {0,0, 0,0};
291    Double_t xMax[nDim] = {1,25,8,5};
292
293    sESDMotherInvMassPtZM = new THnSparseF*[fnCuts];
294    sESDMotherBackInvMassPtZM = new THnSparseF*[fnCuts];
295
296    fBGHandler = new AliGammaConversionAODBGHandler*[fnCuts];
297    for(Int_t iCut = 0; iCut<fnCuts;iCut++){
298       TString cutstring = ((AliConversionCuts*)fCutArray->At(iCut))->GetCutNumber();
299
300       fBackList[iCut] = new TList();
301       fBackList[iCut]->SetName(Form("%s Back histograms",cutstring.Data()));
302       fBackList[iCut]->SetOwner(kTRUE);
303       fCutFolder[iCut]->Add(fBackList[iCut]);
304
305       sESDMotherInvMassPtZM[iCut] = new THnSparseF("Back_Mother_InvMass_Pt_z_m","Back_Mother_InvMass_Pt_z_m",nDim,nBins,xMin,xMax);
306       sESDMotherInvMassPtZM[iCut]->Sumw2();
307       fBackList[iCut]->Add(sESDMotherInvMassPtZM[iCut]);
308       sESDMotherBackInvMassPtZM[iCut] = new THnSparseF("Back_Back_InvMass_Pt_z_m","Back_Back_InvMass_Pt_z_m",nDim,nBins,xMin,xMax);
309       sESDMotherBackInvMassPtZM[iCut]->Sumw2();
310       fBackList[iCut]->Add(sESDMotherBackInvMassPtZM[iCut]);
311
312       if(((AliConversionCuts*)fCutArray->At(iCut))->UseTrackMultiplicity()){
313          fBGHandler[iCut] = new AliGammaConversionAODBGHandler(9,6,((AliConversionCuts*)fCutArray->At(iCut))->NumberOfRotationEvents());
314          fBGHandler[iCut]->Initialize(zBinLimitsArray, multiplicityBinLimitsArrayTracks);
315       }
316       else{
317          fBGHandler[iCut] = new AliGammaConversionAODBGHandler(9,5,((AliConversionCuts*)fCutArray->At(iCut))->NumberOfRotationEvents());
318          fBGHandler[iCut]->Initialize(zBinLimitsArray, multiplicityBinLimitsArrayV0s);
319       }
320    }
321 }
322 //________________________________________________________________________
323 void AliAnalysisTaskGammaConvV1::UserCreateOutputObjects()
324 {
325
326    // Create histograms
327    if(fOutputContainer != NULL){
328       delete fOutputContainer;
329       fOutputContainer = NULL;
330    }
331    if(fOutputContainer == NULL){
332       fOutputContainer = new TList();
333       fOutputContainer->SetOwner(kTRUE);
334    }
335
336    // Array of current cut's gammas
337    fGoodGammas = new TList();
338
339    fCutFolder = new TList*[fnCuts];
340    fESDList = new TList*[fnCuts];
341    fBackList = new TList*[fnCuts];
342    hESDConvGammaPt = new TH1F*[fnCuts];
343    hNEvents = new TH1I*[fnCuts];
344    hNGoodESDTracks = new TH1I*[fnCuts];
345         hNV0Tracks = new TH1I*[fnCuts];
346         
347    if(fDoMesonAnalysis){
348       hESDMotherInvMassPt = new TH2F*[fnCuts];
349       hESDMotherBackInvMassPt = new TH2F*[fnCuts];
350       hESDMotherInvMassEalpha = new TH2F*[fnCuts];
351    }
352    for(Int_t iCut = 0; iCut<fnCuts;iCut++){
353
354       TString cutstring = ((AliConversionCuts*)fCutArray->At(iCut))->GetCutNumber();
355
356       fCutFolder[iCut] = new TList();
357       fCutFolder[iCut]->SetName(Form("Cut Number %s",cutstring.Data()));
358       fCutFolder[iCut]->SetOwner(kTRUE);
359       fOutputContainer->Add(fCutFolder[iCut]);
360       fESDList[iCut] = new TList();
361       fESDList[iCut]->SetName(Form("%s ESD histograms",cutstring.Data()));
362       fESDList[iCut]->SetOwner(kTRUE);
363
364       hNEvents[iCut] = new TH1I("NEvents","NEvents",7,-0.5,6.5);
365       fESDList[iCut]->Add(hNEvents[iCut]);
366       if(fIsHeavyIon) hNGoodESDTracks[iCut] = new TH1I("GoodESDTracks","GoodESDTracks",3000,0,3000);
367       else hNGoodESDTracks[iCut] = new TH1I("GoodESDTracks","GoodESDTracks",200,0,200);
368       fESDList[iCut]->Add(hNGoodESDTracks[iCut]);
369                 if(fIsHeavyIon) hNV0Tracks[iCut] = new TH1I("V0 Multiplicity","V0 Multiplicity",25000,0,25000);
370       else hNV0Tracks[iCut] = new TH1I("V0 Multiplicity","V0 Multiplicity",2000,0,2000);
371       fESDList[iCut]->Add(hNV0Tracks[iCut]);
372    
373       hESDConvGammaPt[iCut] = new TH1F("ESD_ConvGamma_Pt","ESD_ConvGamma_Pt",250,0,25);
374       fESDList[iCut]->Add(hESDConvGammaPt[iCut]);
375
376       if(fDoMesonAnalysis){
377          hESDMotherInvMassPt[iCut] = new TH2F("ESD_Mother_InvMass_Pt","ESD_Mother_InvMass_Pt",1000,0,1,250,0,25);
378          fESDList[iCut]->Add(hESDMotherInvMassPt[iCut]);
379          hESDMotherBackInvMassPt[iCut] = new TH2F("ESD_Background_InvMass_Pt","ESD_Background_InvMass_Pt",1000,0,1,250,0,25);
380          fESDList[iCut]->Add(hESDMotherBackInvMassPt[iCut]);
381          hESDMotherInvMassEalpha[iCut] = new TH2F("ESD_Mother_InvMass_vs_E_alpha","ESD_Mother_InvMass_vs_E_alpha",1000,0,1,250,0,25);
382          fESDList[iCut]->Add(hESDMotherInvMassEalpha[iCut]);
383       }
384
385       fCutFolder[iCut]->Add(fESDList[iCut]);
386    }
387
388    if(fDoMesonAnalysis){
389       InitBack(); // Init Background Handler
390    }
391
392    if(MCEvent()){
393       // MC Histogramms
394       fMCList = new TList*[fnCuts];
395       // True Histogramms
396       fTrueList = new TList*[fnCuts];
397
398       hMCAllGammaPt = new TH1F*[fnCuts];
399                 hMCAllGammaPtAddedSig = new TH1F*[fnCuts];
400       hMCDecayGammaPi0Pt = new TH1F*[fnCuts];
401                 hMCDecayGammaPi0PtAddedSig = new TH1F*[fnCuts];
402       hMCDecayGammaRhoPt = new TH1F*[fnCuts];
403       hMCDecayGammaRhoPtAddedSig = new TH1F*[fnCuts];
404                 hMCDecayGammaEtaPt = new TH1F*[fnCuts];
405                 hMCDecayGammaEtaPtAddedSig = new TH1F*[fnCuts];
406       hMCDecayGammaOmegaPt = new TH1F*[fnCuts];
407                 hMCDecayGammaOmegaPtAddedSig = new TH1F*[fnCuts];
408       hMCDecayGammaEtapPt = new TH1F*[fnCuts];
409                 hMCDecayGammaEtapPtAddedSig = new TH1F*[fnCuts];
410       hMCDecayGammaPhiPt = new TH1F*[fnCuts];
411                 hMCDecayGammaPhiPtAddedSig = new TH1F*[fnCuts];
412       hMCConvGammaPt = new TH1F*[fnCuts];
413                 hMCConvGammaPtAddedSig = new TH1F*[fnCuts];
414       hESDTrueConvGammaPt = new TH1F*[fnCuts];
415                 hESDTrueConvGammaPtAddedSig = new TH1F*[fnCuts];
416       hESDTrueElecCombPt = new TH1F*[fnCuts];
417            hESDTruePionCombPt = new TH1F*[fnCuts];
418            hESDTruePrimaryConvGammaPt = new TH1F*[fnCuts];
419                 hESDTruePrimaryConvGammaPtAddedSig = new TH1F*[fnCuts];
420       hESDTruePrimaryConvGammaESDPtMCPt = new TH2F*[fnCuts];
421                 hESDTruePrimaryConvGammaESDPtMCPtAddedSig = new TH2F*[fnCuts];
422       hESDTrueSecondaryConvGammaPt = new TH1F*[fnCuts];
423                 hESDTrueSecondaryConvGammaPtAddedSig = new TH1F*[fnCuts];
424       hESDTrueSecondaryConvGammaFromK0sPt = new TH1F*[fnCuts];
425                 hESDTrueSecondaryConvGammaFromK0sPtAddedSig = new TH1F*[fnCuts];
426       hESDTrueSecondaryConvGammaFromXFromK0sPt = new TH1F*[fnCuts];
427                 hESDTrueSecondaryConvGammaFromXFromK0sPtAddedSig = new TH1F*[fnCuts];
428                 
429       if(fDoMesonAnalysis){
430          hMCPi0Pt = new TH1F*[fnCuts];
431          hMCEtaPt = new TH1F*[fnCuts];
432          hMCPi0InAccPt = new TH1F*[fnCuts];
433          hMCEtaInAccPt = new TH1F*[fnCuts];
434                         hMCPi0PtAddedSig = new TH1F*[fnCuts];
435          hMCEtaPtAddedSig = new TH1F*[fnCuts];
436          hMCPi0InAccPtAddedSig = new TH1F*[fnCuts];
437          hMCEtaInAccPtAddedSig = new TH1F*[fnCuts];
438
439          hESDTrueMotherInvMassPt = new TH2F*[fnCuts];
440                         hESDTrueMotherInvMassPtAddedSig = new TH2F*[fnCuts];
441          hESDTruePrimaryPi0ESDPtMCPt = new TH2F*[fnCuts];
442                         hESDTruePrimaryPi0ESDPtMCPtAddedSig = new TH2F*[fnCuts];
443          hESDTruePrimaryMotherInvMassMCPt = new TH2F*[fnCuts];
444                         hESDTruePrimaryMotherInvMassMCPtAddedSig = new TH2F*[fnCuts];
445          hESDTrueSecondaryMotherInvMassPt = new TH2F*[fnCuts];
446                         hESDTrueSecondaryMotherInvMassPtAddedSig = new TH2F*[fnCuts];
447          hESDTrueSecondaryMotherFromK0sInvMassPt = new TH2F*[fnCuts];
448                         hESDTrueSecondaryMotherFromK0sInvMassPtAddedSig = new TH2F*[fnCuts];
449          hESDTrueBckGGInvMassPt = new TH2F*[fnCuts];
450          hESDTrueBckContInvMassPt = new TH2F*[fnCuts];
451          hESDTrueMotherDalitzInvMassPt = new TH2F*[fnCuts];
452                         hESDTrueMotherDalitzInvMassPtAddedSig = new TH2F*[fnCuts];
453       }
454
455       for(Int_t iCut = 0; iCut<fnCuts;iCut++){
456          TString cutstring = ((AliConversionCuts*)fCutArray->At(iCut))->GetCutNumber();
457          fMCList[iCut] = new TList();
458          fMCList[iCut]->SetName(Form("%s MC histograms",cutstring.Data()));
459          fMCList[iCut]->SetOwner(kTRUE);
460          fCutFolder[iCut]->Add(fMCList[iCut]);
461
462          hMCAllGammaPt[iCut] = new TH1F("MC_AllGamma_Pt","MC_AllGamma_Pt",250,0,25);
463          fMCList[iCut]->Add(hMCAllGammaPt[iCut]);
464                         hMCAllGammaPtAddedSig[iCut] = new TH1F("MC_AllGamma_Pt_AddedSig","MC_AllGamma_Pt_AddedSig",250,0,25);
465          fMCList[iCut]->Add(hMCAllGammaPtAddedSig[iCut]);
466          hMCDecayGammaPi0Pt[iCut] = new TH1F("MC_DecayGammaPi0_Pt","MC_DecayGammaPi0_Pt",250,0,25);
467          fMCList[iCut]->Add(hMCDecayGammaPi0Pt[iCut]);
468          hMCDecayGammaPi0PtAddedSig[iCut] = new TH1F("MC_DecayGammaPi0_Pt_AddedSig","MC_DecayGammaPi0_Pt_AddedSig",250,0,25);
469          fMCList[iCut]->Add(hMCDecayGammaPi0PtAddedSig[iCut]);
470          hMCDecayGammaRhoPt[iCut] = new TH1F("MC_DecayGammaRho_Pt","MC_DecayGammaRho_Pt",250,0,25);
471          fMCList[iCut]->Add(hMCDecayGammaRhoPt[iCut]);
472          hMCDecayGammaRhoPtAddedSig[iCut] = new TH1F("MC_DecayGammaRho_Pt_AddedSig","MC_DecayGammaRho_Pt_AddedSig",250,0,25);
473          fMCList[iCut]->Add(hMCDecayGammaRhoPtAddedSig[iCut]);
474          hMCDecayGammaEtaPt[iCut] = new TH1F("MC_DecayGammaEta_Pt","MC_DecayGammaEta_Pt",250,0,25);
475          fMCList[iCut]->Add(hMCDecayGammaEtaPt[iCut]);
476          hMCDecayGammaEtaPtAddedSig[iCut] = new TH1F("MC_DecayGammaEta_Pt_AddedSig","MC_DecayGammaEta_Pt_AddedSig",250,0,25);
477          fMCList[iCut]->Add(hMCDecayGammaEtaPtAddedSig[iCut]);
478          hMCDecayGammaOmegaPt[iCut] = new TH1F("MC_DecayGammaOmega_Pt","MC_DecayGammaOmmega_Pt",250,0,25);
479          fMCList[iCut]->Add(hMCDecayGammaOmegaPt[iCut]);
480          hMCDecayGammaOmegaPtAddedSig[iCut] = new TH1F("MC_DecayGammaOmega_Pt_AddedSig","MC_DecayGammaOmmega_Pt_AddedSig",250,0,25);
481          fMCList[iCut]->Add(hMCDecayGammaOmegaPtAddedSig[iCut]);
482          hMCDecayGammaEtapPt[iCut] = new TH1F("MC_DecayGammaEtap_Pt","MC_DecayGammaEtap_Pt",250,0,25);
483          fMCList[iCut]->Add(hMCDecayGammaEtapPt[iCut]);
484          hMCDecayGammaEtapPtAddedSig[iCut] = new TH1F("MC_DecayGammaEtap_Pt_AddedSig","MC_DecayGammaEtap_Pt_AddedSig",250,0,25);
485          fMCList[iCut]->Add(hMCDecayGammaEtapPtAddedSig[iCut]);
486          hMCDecayGammaPhiPt[iCut] = new TH1F("MC_DecayGammaPhi_Pt","MC_DecayGammaPhi_Pt",250,0,25);
487          fMCList[iCut]->Add(hMCDecayGammaPhiPt[iCut]);
488          hMCDecayGammaPhiPtAddedSig[iCut] = new TH1F("MC_DecayGammaPhi_Pt_AddedSig","MC_DecayGammaPhi_Pt_AddedSig",250,0,25);
489          fMCList[iCut]->Add(hMCDecayGammaPhiPtAddedSig[iCut]);
490          hMCConvGammaPt[iCut] = new TH1F("MC_ConvGamma_Pt","MC_ConvGamma_Pt",250,0,25);
491          fMCList[iCut]->Add(hMCConvGammaPt[iCut]);
492          hMCConvGammaPtAddedSig[iCut] = new TH1F("MC_ConvGamma_Pt_AddedSig","MC_ConvGamma_Pt_AddedSig",250,0,25);
493          fMCList[iCut]->Add(hMCConvGammaPtAddedSig[iCut]);
494          if(fDoMesonAnalysis){
495             hMCPi0Pt[iCut] = new TH1F("MC_Pi0_Pt","MC_Pi0_Pt",250,0,25);
496             fMCList[iCut]->Add(hMCPi0Pt[iCut]);
497             hMCEtaPt[iCut] = new TH1F("MC_Eta_Pt","MC_Eta_Pt",250,0,25);
498             fMCList[iCut]->Add(hMCEtaPt[iCut]);
499             hMCPi0InAccPt[iCut] = new TH1F("MC_Pi0InAcc_Pt","MC_Pi0InAcc_Pt",250,0,25);
500             fMCList[iCut]->Add(hMCPi0InAccPt[iCut]);
501             hMCEtaInAccPt[iCut] = new TH1F("MC_EtaInAcc_Pt","MC_EtaInAcc_Pt",250,0,25);
502             fMCList[iCut]->Add(hMCEtaInAccPt[iCut]);
503                                 hMCPi0PtAddedSig[iCut] = new TH1F("MC_Pi0_Pt_AddedSig","MC_Pi0_Pt_AddedSig",250,0,25);
504             fMCList[iCut]->Add(hMCPi0PtAddedSig[iCut]);
505             hMCEtaPtAddedSig[iCut] = new TH1F("MC_Eta_Pt_AddedSig","MC_Eta_Pt_AddedSig",250,0,25);
506             fMCList[iCut]->Add(hMCEtaPtAddedSig[iCut]);
507             hMCPi0InAccPtAddedSig[iCut] = new TH1F("MC_Pi0InAcc_Pt_AddedSig","MC_Pi0InAcc_Pt_AddedSig",250,0,25);
508             fMCList[iCut]->Add(hMCPi0InAccPtAddedSig[iCut]);
509             hMCEtaInAccPtAddedSig[iCut] = new TH1F("MC_EtaInAcc_Pt_AddedSig","MC_EtaInAcc_Pt_AddedSig",250,0,25);
510             fMCList[iCut]->Add(hMCEtaInAccPtAddedSig[iCut]);
511          }
512          fTrueList[iCut] = new TList();
513          fTrueList[iCut]->SetName(Form("%s True histograms",cutstring.Data()));
514          fTrueList[iCut]->SetOwner(kTRUE);
515          fCutFolder[iCut]->Add(fTrueList[iCut]);
516
517          hESDTrueConvGammaPt[iCut] = new TH1F("ESD_TrueConvGamma_Pt","ESD_TrueConvGamma_Pt",250,0,25);
518          fTrueList[iCut]->Add(hESDTrueConvGammaPt[iCut]);
519                         hESDTrueConvGammaPtAddedSig[iCut] = new TH1F("ESD_TrueConvGamma_Pt_AddedSig","ESD_TrueConvGamma_Pt_AddedSig",250,0,25);
520          fTrueList[iCut]->Add(hESDTrueConvGammaPtAddedSig[iCut]);
521          hESDTrueElecCombPt[iCut] = new TH1F("ESD_TrueElecComb_Pt","ESD_TrueElecComb_Pt",250,0,25);
522          fTrueList[iCut]->Add(hESDTrueElecCombPt[iCut]);
523          hESDTruePionCombPt[iCut] = new TH1F("ESD_TruePionComb_Pt","ESD_TruePionComb_Pt",250,0,25);
524          fTrueList[iCut]->Add(hESDTruePionCombPt[iCut]);
525          hESDTruePrimaryConvGammaPt[iCut] = new TH1F("ESD_TruePrimaryConvGamma_Pt","ESD_TruePrimaryConvGamma_Pt",250,0,25);
526          fTrueList[iCut]->Add(hESDTruePrimaryConvGammaPt[iCut]);
527          hESDTruePrimaryConvGammaPtAddedSig[iCut] = new TH1F("ESD_TruePrimaryConvGamma_Pt_AddedSig","ESD_TruePrimaryConvGamma_Pt_AddedSig",250,0,25);
528          fTrueList[iCut]->Add(hESDTruePrimaryConvGammaPtAddedSig[iCut]);
529          hESDTrueSecondaryConvGammaPt[iCut] = new TH1F("ESD_TrueSecondaryConvGamma_Pt","ESD_TrueSecondaryConvGamma_Pt",250,0,25);
530          fTrueList[iCut]->Add(hESDTrueSecondaryConvGammaPt[iCut]);
531          hESDTrueSecondaryConvGammaPtAddedSig[iCut] = new TH1F("ESD_TrueSecondaryConvGamma_Pt_AddedSig", "ESD_TrueSecondaryConvGamma_Pt_AddedSig",250,0,25);
532          fTrueList[iCut]->Add(hESDTrueSecondaryConvGammaPtAddedSig[iCut]);
533          hESDTrueSecondaryConvGammaFromK0sPt[iCut] = new TH1F("ESD_TrueSecondaryConvGammaFromK0s_Pt","ESD_TrueSecondaryConvGammaFromK0s_Pt",250,0,25);
534          fTrueList[iCut]->Add(hESDTrueSecondaryConvGammaFromK0sPt[iCut]);
535          hESDTrueSecondaryConvGammaFromK0sPtAddedSig[iCut] = new TH1F("ESD_TrueSecondaryConvGammaFromK0s_Pt_AddedSig", "ESD_TrueSecondaryConvGammaFromK0s_Pt_AddedSig",250,0,25);
536          fTrueList[iCut]->Add(hESDTrueSecondaryConvGammaFromK0sPtAddedSig[iCut]);
537          hESDTrueSecondaryConvGammaFromXFromK0sPt[iCut] = new TH1F("ESD_TrueSecondaryConvGammaFromXFromK0s_Pt", "ESD_TrueSecondaryConvGammaFromXFromK0s_Pt",250,0,25);
538          fTrueList[iCut]->Add(hESDTrueSecondaryConvGammaFromXFromK0sPt[iCut]);
539          hESDTrueSecondaryConvGammaFromXFromK0sPtAddedSig[iCut] = new TH1F("ESD_TrueSecondaryConvGammaFromXFromK0s_Pt_AddedSig", "ESD_TrueSecondaryConvGammaFromXFromK0s_Pt_AddedSig",250,0,25);
540          fTrueList[iCut]->Add(hESDTrueSecondaryConvGammaFromXFromK0sPtAddedSig[iCut]);
541          hESDTruePrimaryConvGammaESDPtMCPt[iCut] = new TH2F("ESD_TruePrimaryConvGammaESD_PtMCPt", "ESD_TruePrimaryConvGammaESD_PtMCPt",250,0,25,250,0,25);
542          fTrueList[iCut]->Add(hESDTruePrimaryConvGammaESDPtMCPt[iCut]);
543                         hESDTruePrimaryConvGammaESDPtMCPtAddedSig[iCut] = new TH2F("ESD_TruePrimaryConvGammaESD_PtMCPt_AddedSig", "ESD_TruePrimaryConvGammaESD_PtMCPt_AddedSig",250,0,25,250,0,25);
544          fTrueList[iCut]->Add(hESDTruePrimaryConvGammaESDPtMCPtAddedSig[iCut]);
545
546          if(fDoMesonAnalysis){
547             hESDTrueMotherInvMassPt[iCut] = new TH2F("ESD_TrueMother_InvMass_Pt","ESD_TrueMother_InvMass_Pt",1000,0,1,250,0,25);
548             fTrueList[iCut]->Add(hESDTrueMotherInvMassPt[iCut]);
549                                 hESDTrueMotherInvMassPtAddedSig[iCut] = new TH2F("ESD_TrueMother_InvMass_Pt_AddedSig", "ESD_TrueMother_InvMass_Pt_AddedSig", 1000,0,1,250,0,25);
550             fTrueList[iCut]->Add(hESDTrueMotherInvMassPtAddedSig[iCut]);
551             hESDTruePrimaryPi0ESDPtMCPt[iCut] = new TH2F("ESD_TruePrimaryPi0_ESDPt_MCPt","ESD_TruePrimaryPi0_ESDPt_MCPt",250,0,25,250,0,25);
552             fTrueList[iCut]->Add(hESDTruePrimaryPi0ESDPtMCPt[iCut]);
553             hESDTruePrimaryPi0ESDPtMCPtAddedSig[iCut] = new TH2F("ESD_TruePrimaryPi0_ESDPt_MCPt_AddedSig", "ESD_TruePrimaryPi0_ESDPt_MCPt_AddedSig", 250,0,25,250,0,25);
554             fTrueList[iCut]->Add(hESDTruePrimaryPi0ESDPtMCPtAddedSig[iCut]);
555             hESDTruePrimaryMotherInvMassMCPt[iCut] = new TH2F("ESD_TruePrimaryMother_InvMass_MCPt", "ESD_TruePrimaryMother_InvMass_MCPt", 1000,0,1,250,0,25);
556             fTrueList[iCut]->Add(hESDTruePrimaryMotherInvMassMCPt[iCut]);
557             hESDTruePrimaryMotherInvMassMCPtAddedSig[iCut] = new TH2F("ESD_TruePrimaryMother_InvMass_MCPt_AddedSig", "ESD_TruePrimaryMother_InvMass_MCPt_AddedSig", 1000,0,1,250,0,25);
558             fTrueList[iCut]->Add(hESDTruePrimaryMotherInvMassMCPtAddedSig[iCut]);
559             hESDTrueSecondaryMotherInvMassPt[iCut] = new TH2F("ESD_TrueSecondaryMother_InvMass_Pt", "ESD_TrueSecondaryMother_InvMass_Pt", 1000,0,1,250,0,25);
560             fTrueList[iCut]->Add(hESDTrueSecondaryMotherInvMassPt[iCut]);
561             hESDTrueSecondaryMotherFromK0sInvMassPt[iCut] = new TH2F("ESD_TrueSecondaryMotherFromK0s_InvMass_Pt","ESD_TrueSecondaryMotherFromK0s_InvMass_Pt",1000,0,1,250,0,25);
562             fTrueList[iCut]->Add(hESDTrueSecondaryMotherFromK0sInvMassPt[iCut]);
563             hESDTrueSecondaryMotherInvMassPtAddedSig[iCut] = new TH2F("ESD_TrueSecondaryMother_InvMass_Pt_AddedSig","ESD_TrueSecondaryMother_InvMass_Pt_AddedSig",1000,0,1,250,0,25);
564             fTrueList[iCut]->Add(hESDTrueSecondaryMotherInvMassPtAddedSig[iCut]);
565             hESDTrueSecondaryMotherFromK0sInvMassPtAddedSig[iCut] = new TH2F("ESD_TrueSecondaryMotherFromK0s_InvMass_Pt_AddedSig","ESD_TrueSecondaryMotherFromK0s_InvMass_Pt_AddedSig",1000,0,1,250,0,25);
566             fTrueList[iCut]->Add(hESDTrueSecondaryMotherFromK0sInvMassPtAddedSig[iCut]);
567                                 hESDTrueBckGGInvMassPt[iCut] = new TH2F("ESD_TrueBckGG_InvMass_Pt","ESD_TrueBckGG_InvMass_Pt",1000,0,1,250,0,25);
568             fTrueList[iCut]->Add(hESDTrueBckGGInvMassPt[iCut]);
569                  hESDTrueBckContInvMassPt[iCut] = new TH2F("ESD_TrueBckCont_InvMass_Pt","ESD_TrueBckCont_InvMass_Pt",1000,0,1,250,0,25);
570             fTrueList[iCut]->Add(hESDTrueBckContInvMassPt[iCut]);
571             hESDTrueMotherDalitzInvMassPt[iCut] = new TH2F("ESD_TrueDalitz_InvMass_Pt","ESD_TrueDalitz_InvMass_Pt",1000,0,1,250,0,25);
572             fTrueList[iCut]->Add(hESDTrueMotherDalitzInvMassPt[iCut]);
573             hESDTrueMotherDalitzInvMassPtAddedSig[iCut] = new TH2F("ESD_TrueDalitz_InvMass_Pt_AddedSig", "ESD_TrueDalitz_InvMass_Pt_AddedSig", 1000,0,1,250,0,25);
574             fTrueList[iCut]->Add(hESDTrueMotherDalitzInvMassPtAddedSig[iCut]);
575
576          }
577       }
578    }
579
580    PostData(1, fOutputContainer);
581 }
582
583 //_____________________________________________________________________________
584 void AliAnalysisTaskGammaConvV1::UserExec(Option_t *)
585 {
586    //
587    // Called for each event
588    //
589    fV0Reader=(AliV0ReaderV1*)AliAnalysisManager::GetAnalysisManager()->GetTask("V0ReaderV1");
590    if(!fV0Reader){printf("Error: No V0 Reader");return;} // GetV0Reader
591
592    Int_t eventQuality = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetEventQuality();
593    if(eventQuality != 0){// Event Not Accepted
594       for(Int_t iCut = 0; iCut<fnCuts; iCut++){
595          hNEvents[iCut]->Fill(eventQuality);
596       }
597       return;
598    }
599
600    fMCEvent = MCEvent();
601    fESDEvent = (AliESDEvent*) InputEvent();
602    fReaderGammas = fV0Reader->GetReconstructedGammas(); // Gammas from default Cut
603    CountESDTracks(); // Estimate Event Multiplicity
604
605    for(Int_t iCut = 0; iCut<fnCuts; iCut++){
606       fiCut = iCut;
607       if(fIsHeavyIon && !((AliConversionCuts*)fCutArray->At(iCut))->IsCentralitySelected(fESDEvent)){
608          hNEvents[iCut]->Fill(1); // Check Centrality --> Not Accepted => eventQuality = 1
609          continue;
610       }
611       hNEvents[iCut]->Fill(eventQuality);
612
613       hNGoodESDTracks[iCut]->Fill(fNumberOfESDTracks);
614                 hNV0Tracks[iCut]->Fill(fESDEvent->GetVZEROData()->GetMTotV0A()+fESDEvent->GetVZEROData()->GetMTotV0C());
615       if(fMCEvent){ // Process MC Particle
616          fMCStack = fMCEvent->Stack();
617          ProcessMCParticles();
618       }
619
620       ProcessPhotonCandidates(); // Process this cuts gammas
621
622       if(fDoMesonAnalysis){ // Meson Analysis
623          if(((AliConversionCuts*)fCutArray->At(iCut))->UseMCPSmearing() && fMCEvent){
624             fUnsmearedPx = new Double_t[fGoodGammas->GetEntries()]; // Store unsmeared Momenta
625             fUnsmearedPy = new Double_t[fGoodGammas->GetEntries()];
626             fUnsmearedPz = new Double_t[fGoodGammas->GetEntries()];
627             fUnsmearedE =  new Double_t[fGoodGammas->GetEntries()];
628
629             for(Int_t gamma=0;gamma<fGoodGammas->GetEntries();gamma++){ // Smear the AODPhotons in MC
630                fUnsmearedPx[gamma] = ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->Px();
631                fUnsmearedPy[gamma] = ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->Py();
632                fUnsmearedPz[gamma] = ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->Pz();
633                fUnsmearedE[gamma] =  ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->E();
634                ((AliConversionCuts*)fCutArray->At(iCut))->SmearParticle(dynamic_cast<AliAODConversionPhoton*>(fGoodGammas->At(gamma)));
635             }
636          }
637          CalculatePi0Candidates(); // Combine Gammas
638          CalculateBackground(); // Combinatorial Background
639          UpdateEventByEventData(); // Store Event for mixed Events
640
641          if(((AliConversionCuts*)fCutArray->At(iCut))->UseMCPSmearing() && fMCEvent){
642             for(Int_t gamma=0;gamma<fGoodGammas->GetEntries();gamma++){ // Smear the AODPhotons in MC
643                ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->SetPx(fUnsmearedPx[gamma]); // Reset Unsmeared Momenta
644                ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->SetPy(fUnsmearedPy[gamma]);
645                ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->SetPz(fUnsmearedPz[gamma]);
646                ((AliAODConversionPhoton*)fGoodGammas->At(gamma))->SetE(fUnsmearedE[gamma]);
647             }
648             delete[] fUnsmearedPx; fUnsmearedPx = 0x0;
649             delete[] fUnsmearedPy; fUnsmearedPy = 0x0;
650             delete[] fUnsmearedPz; fUnsmearedPz = 0x0;
651             delete[] fUnsmearedE;  fUnsmearedE  = 0x0;
652          }
653       }
654       fGoodGammas->Clear(); // delete this cuts good gammas
655    }
656
657    PostData(1, fOutputContainer);
658 }
659 //________________________________________________________________________
660 void AliAnalysisTaskGammaConvV1::ProcessPhotonCandidates()
661 {
662    Int_t nV0 = 0;
663    TList *GoodGammasStepOne = new TList();
664    TList *GoodGammasStepTwo = new TList();
665    // Loop over Photon Candidates allocated by ReaderV1
666    for(Int_t i = 0; i < fReaderGammas->GetEntriesFast(); i++){
667       AliAODConversionPhoton* PhotonCandidate = (AliAODConversionPhoton*) fReaderGammas->At(i);
668       if(!PhotonCandidate) continue;
669       if(!((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelected(PhotonCandidate,fESDEvent)) continue;
670       
671       if(!((AliConversionCuts*)fCutArray->At(fiCut))->UseElecSharingCut() && 
672          !((AliConversionCuts*)fCutArray->At(fiCut))->UseToCloseV0sCut()){ // if no post reader loop is required add to events good gammas
673          fGoodGammas->Add(PhotonCandidate);
674          hESDConvGammaPt[fiCut]->Fill(PhotonCandidate->Pt());
675          if(fMCEvent){
676             ProcessTruePhotonCandidates(PhotonCandidate);
677          }
678       }
679       else if(((AliConversionCuts*)fCutArray->At(fiCut))->UseElecSharingCut()){ // if Shared Electron cut is enabled, Fill array, add to step one
680          ((AliConversionCuts*)fCutArray->At(fiCut))->FillElectonLabelArray(PhotonCandidate,nV0);
681          nV0++;
682          GoodGammasStepOne->Add(PhotonCandidate);
683       }
684       else if(!((AliConversionCuts*)fCutArray->At(fiCut))->UseElecSharingCut() && 
685               ((AliConversionCuts*)fCutArray->At(fiCut))->UseToCloseV0sCut()){ // shared electron is disabled, step one not needed -> step two
686          GoodGammasStepTwo->Add(PhotonCandidate);
687       }
688    }
689    if(((AliConversionCuts*)fCutArray->At(fiCut))->UseElecSharingCut()){
690       for(Int_t i = 0;i<GoodGammasStepOne->GetEntries();i++){
691          AliAODConversionPhoton *PhotonCandidate= (AliAODConversionPhoton*) GoodGammasStepOne->At(i);
692          if(!PhotonCandidate) continue;
693          if(!((AliConversionCuts*)fCutArray->At(fiCut))->RejectSharedElectronV0s(PhotonCandidate,i,GoodGammasStepOne->GetEntries())) continue;
694          if(!((AliConversionCuts*)fCutArray->At(fiCut))->UseToCloseV0sCut()){ // To Colse v0s cut diabled, step two not needed 
695             fGoodGammas->Add(PhotonCandidate);
696             hESDConvGammaPt[fiCut]->Fill(PhotonCandidate->Pt());
697             if(fMCEvent){
698                ProcessTruePhotonCandidates(PhotonCandidate);
699             }
700          }
701          else GoodGammasStepTwo->Add(PhotonCandidate); // Close v0s cut enabled -> add to list two
702       }
703    }
704    if(((AliConversionCuts*)fCutArray->At(fiCut))->UseToCloseV0sCut()){
705       for(Int_t i = 0;i<GoodGammasStepTwo->GetEntries();i++){
706          AliAODConversionPhoton* PhotonCandidate = (AliAODConversionPhoton*) GoodGammasStepTwo->At(i);
707          if(!PhotonCandidate) continue;
708          if(!((AliConversionCuts*)fCutArray->At(fiCut))->RejectToCloseV0s(PhotonCandidate,GoodGammasStepTwo,i)) continue;
709          fGoodGammas->Add(PhotonCandidate); // Add gamma to current cut TList
710          hESDConvGammaPt[fiCut]->Fill(PhotonCandidate->Pt()); // Differences to old V0Reader in p_t due to conversion KF->TLorentzVector
711          if(fMCEvent){
712             ProcessTruePhotonCandidates(PhotonCandidate);
713          }
714       }
715    }
716
717    delete GoodGammasStepOne;
718    GoodGammasStepOne = 0x0;
719    delete GoodGammasStepTwo;
720    GoodGammasStepTwo = 0x0;
721
722 }
723 //________________________________________________________________________
724 void AliAnalysisTaskGammaConvV1::ProcessTruePhotonCandidates(AliAODConversionPhoton *TruePhotonCandidate)
725 {
726    // Process True Photons
727    AliStack *MCStack = fMCEvent->Stack();
728    TParticle *posDaughter = TruePhotonCandidate->GetPositiveMCDaughter(MCStack);
729    TParticle *negDaughter = TruePhotonCandidate->GetNegativeMCDaughter(MCStack);
730
731    if(posDaughter == NULL || negDaughter == NULL) return; // One particle does not exist
732    if(posDaughter->GetMother(0) != negDaughter->GetMother(0)){  // Not Same Mother == Combinatorial Bck
733       if(TMath::Abs(posDaughter->GetPdgCode())==11 && TMath::Abs(negDaughter->GetPdgCode())==11)
734          hESDTrueElecCombPt[fiCut]->Fill(TruePhotonCandidate->Pt()); //Electron Combinatorial
735       if(TMath::Abs(posDaughter->GetPdgCode())==211 || TMath::Abs(negDaughter->GetPdgCode())==211)
736          hESDTruePionCombPt[fiCut]->Fill(TruePhotonCandidate->Pt()); //At least on Pion Combinatorial
737       return;
738    }
739    if(TMath::Abs(posDaughter->GetPdgCode())!=11 || TMath::Abs(negDaughter->GetPdgCode())!=11) return; //One Particle is not electron
740    if(posDaughter->GetPdgCode()==negDaughter->GetPdgCode()) return; // Same Charge
741    if(posDaughter->GetUniqueID() != 5 || negDaughter->GetUniqueID() !=5) return;// check if the daughters come from a conversion
742
743    TParticle *Photon = TruePhotonCandidate->GetMCParticle(MCStack);
744    if(Photon->GetPdgCode() != 22) return; // Mother is no Photon
745
746         if (!fMCEvent->IsFromBGEvent(TruePhotonCandidate->GetMCParticleLabel(fMCStack)) && fIsHeavyIon){
747       // True Photon
748                 hESDTrueConvGammaPtAddedSig[fiCut]->Fill(TruePhotonCandidate->Pt());
749
750                 if(posDaughter->GetMother(0) <= MCStack->GetNprimary()){
751                         // Count just primary MC Gammas as true --> For Ratio esdtruegamma / mcconvgamma
752                         hESDTruePrimaryConvGammaPtAddedSig[fiCut]->Fill(TruePhotonCandidate->Pt());
753                         hESDTruePrimaryConvGammaESDPtMCPtAddedSig[fiCut]->Fill(TruePhotonCandidate->Pt(),Photon->Pt());
754                 }
755                 else{
756                         hESDTrueSecondaryConvGammaPtAddedSig[fiCut]->Fill(TruePhotonCandidate->Pt());
757                         if(MCStack->Particle(Photon->GetMother(0))->GetPdgCode() == 310){
758                                 hESDTrueSecondaryConvGammaFromK0sPtAddedSig[fiCut]->Fill(TruePhotonCandidate->Pt());
759                         }
760                         if(MCStack->Particle(Photon->GetMother(0))->GetMother(0) > -1 &&
761                                 MCStack->Particle(MCStack->Particle(Photon->GetMother(0))->GetMother(0))->GetPdgCode() == 310){
762                                 hESDTrueSecondaryConvGammaFromXFromK0sPtAddedSig[fiCut]->Fill(TruePhotonCandidate->Pt());
763                         }
764                 }
765         } else {
766       // True Photon
767                 hESDTrueConvGammaPt[fiCut]->Fill(TruePhotonCandidate->Pt());
768
769                 if(posDaughter->GetMother(0) <= MCStack->GetNprimary()){
770                         // Count just primary MC Gammas as true --> For Ratio esdtruegamma / mcconvgamma
771                         hESDTruePrimaryConvGammaPt[fiCut]->Fill(TruePhotonCandidate->Pt());
772                         hESDTruePrimaryConvGammaESDPtMCPt[fiCut]->Fill(TruePhotonCandidate->Pt(),Photon->Pt());
773                 }
774                 else{
775                         hESDTrueSecondaryConvGammaPt[fiCut]->Fill(TruePhotonCandidate->Pt());
776                         if(MCStack->Particle(Photon->GetMother(0))->GetPdgCode() == 310){
777                                 hESDTrueSecondaryConvGammaFromK0sPt[fiCut]->Fill(TruePhotonCandidate->Pt());
778                         }
779                         if(MCStack->Particle(Photon->GetMother(0))->GetMother(0) > -1 &&
780                                 MCStack->Particle(MCStack->Particle(Photon->GetMother(0))->GetMother(0))->GetPdgCode() == 310){
781                                 hESDTrueSecondaryConvGammaFromXFromK0sPt[fiCut]->Fill(TruePhotonCandidate->Pt());
782                         }
783                 }               
784         }
785 }
786 //________________________________________________________________________
787 void AliAnalysisTaskGammaConvV1::ProcessMCParticles()
788 {
789
790    // Loop over all primary MC particle
791    for(Int_t i = 0; i < fMCStack->GetNprimary(); i++) {
792       TParticle* particle = (TParticle *)fMCStack->Particle(i);
793       if (!particle) continue;
794
795                 if (!fMCEvent->IsFromBGEvent(i) && fIsHeavyIon){
796                         if(((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(particle,fMCStack,kFALSE)){
797                                 hMCAllGammaPtAddedSig[fiCut]->Fill(particle->Pt()); // All MC Gamma
798                                 if(particle->GetMother(0) >-1){ // Meson Decay Gamma
799                                         switch(fMCStack->Particle(particle->GetMother(0))->GetPdgCode()){
800                                         case 111: // Pi0
801                                                 hMCDecayGammaPi0PtAddedSig[fiCut]->Fill(particle->Pt());
802                                                 break;
803                                         case 113: // Rho0
804                                                 hMCDecayGammaRhoPtAddedSig[fiCut]->Fill(particle->Pt());
805                                                 break;
806                                         case 221: // Eta
807                                                 hMCDecayGammaEtaPtAddedSig[fiCut]->Fill(particle->Pt());
808                                                 break;
809                                         case 223: // Omega
810                                                 hMCDecayGammaOmegaPtAddedSig[fiCut]->Fill(particle->Pt());
811                                                 break;
812                                         case 331: // Eta'
813                                                 hMCDecayGammaEtapPtAddedSig[fiCut]->Fill(particle->Pt());
814                                                 break;
815                                         case 333: // Phi
816                                                 hMCDecayGammaPhiPtAddedSig[fiCut]->Fill(particle->Pt());
817                                                 break;
818                                         }
819                                 }
820                         }
821                         if(((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(particle,fMCStack,kTRUE)){
822                                 hMCConvGammaPtAddedSig[fiCut]->Fill(particle->Pt());
823                         } // Converted MC Gamma
824                         if(fDoMesonAnalysis){
825                                 if(((AliConversionCuts*)fCutArray->At(fiCut))->MesonIsSelectedMC(particle,fMCStack,kFALSE)){
826                                         if(particle->GetPdgCode() == 111)hMCPi0PtAddedSig[fiCut]->Fill(particle->Pt()); // All MC Pi0
827                                         if(particle->GetPdgCode() == 221)hMCEtaPtAddedSig[fiCut]->Fill(particle->Pt()); // All MC Eta
828                                         // Check the acceptance for both gammas
829                                         if(particle->GetNDaughters() == 2){
830                                                 TParticle* daughter0 = (TParticle*)fMCStack->Particle(particle->GetFirstDaughter());
831                                                 TParticle* daughter1 = (TParticle*)fMCStack->Particle(particle->GetLastDaughter());
832                                                 if(((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(daughter0,fMCStack,kFALSE) &&
833                                                         ((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(daughter1,fMCStack,kFALSE) ){
834                                                         if(particle->GetPdgCode() == 111)hMCPi0InAccPtAddedSig[fiCut]->Fill(particle->Pt()); // MC Pi0 with gamma in acc
835                                                         if(particle->GetPdgCode() == 221)hMCEtaInAccPtAddedSig[fiCut]->Fill(particle->Pt()); // MC Eta with gamma in acc
836                                                 }
837                                         }
838                                 }
839                         }
840                 } else {
841                         if(((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(particle,fMCStack,kFALSE)){
842                                 hMCAllGammaPt[fiCut]->Fill(particle->Pt()); // All MC Gamma
843                                 if(particle->GetMother(0) >-1){ // Meson Decay Gamma
844                                         switch(fMCStack->Particle(particle->GetMother(0))->GetPdgCode()){
845                                         case 111: // Pi0
846                                                 hMCDecayGammaPi0Pt[fiCut]->Fill(particle->Pt());
847                                                 break;
848                                         case 113: // Rho0
849                                                 hMCDecayGammaRhoPt[fiCut]->Fill(particle->Pt());
850                                                 break;
851                                         case 221: // Eta
852                                                 hMCDecayGammaEtaPt[fiCut]->Fill(particle->Pt());
853                                                 break;
854                                         case 223: // Omega
855                                                 hMCDecayGammaOmegaPt[fiCut]->Fill(particle->Pt());
856                                                 break;
857                                         case 331: // Eta'
858                                                 hMCDecayGammaEtapPt[fiCut]->Fill(particle->Pt());
859                                                 break;
860                                         case 333: // Phi
861                                                 hMCDecayGammaPhiPt[fiCut]->Fill(particle->Pt());
862                                                 break;
863                                         }
864                                 }
865                         }
866                         if(((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(particle,fMCStack,kTRUE)){
867                                 hMCConvGammaPt[fiCut]->Fill(particle->Pt());
868                         } // Converted MC Gamma
869                         if(fDoMesonAnalysis){
870                                 if(((AliConversionCuts*)fCutArray->At(fiCut))->MesonIsSelectedMC(particle,fMCStack,kFALSE)){
871                                         if(particle->GetPdgCode() == 111)hMCPi0Pt[fiCut]->Fill(particle->Pt()); // All MC Pi0
872                                         if(particle->GetPdgCode() == 221)hMCEtaPt[fiCut]->Fill(particle->Pt()); // All MC Eta
873                                         // Check the acceptance for both gammas
874                                         if(particle->GetNDaughters() == 2){
875                                                 TParticle* daughter0 = (TParticle*)fMCStack->Particle(particle->GetFirstDaughter());
876                                                 TParticle* daughter1 = (TParticle*)fMCStack->Particle(particle->GetLastDaughter());
877                                                 if(((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(daughter0,fMCStack,kFALSE) &&
878                                                         ((AliConversionCuts*)fCutArray->At(fiCut))->PhotonIsSelectedMC(daughter1,fMCStack,kFALSE) ){
879                                                         if(particle->GetPdgCode() == 111)hMCPi0InAccPt[fiCut]->Fill(particle->Pt()); // MC Pi0 with gamma in acc
880                                                         if(particle->GetPdgCode() == 221)hMCEtaInAccPt[fiCut]->Fill(particle->Pt()); // MC Eta with gamma in acc
881                                                 }
882                                         }
883                                 }
884                         }
885                 }
886    }
887 }
888 //________________________________________________________________________
889 void AliAnalysisTaskGammaConvV1::CalculatePi0Candidates(){
890
891    // Conversion Gammas
892    if(fGoodGammas->GetEntries()>1){
893       for(Int_t firstGammaIndex=0;firstGammaIndex<fGoodGammas->GetEntries()-1;firstGammaIndex++){
894          AliAODConversionPhoton *gamma0=dynamic_cast<AliAODConversionPhoton*>(fGoodGammas->At(firstGammaIndex));
895          for(Int_t secondGammaIndex=firstGammaIndex+1;secondGammaIndex<fGoodGammas->GetEntries();secondGammaIndex++){
896             AliAODConversionPhoton *gamma1=dynamic_cast<AliAODConversionPhoton*>(fGoodGammas->At(secondGammaIndex));
897             //Check for same Electron ID
898             if(gamma0->GetTrackLabelPositive() == gamma1->GetTrackLabelPositive() ||
899                gamma0->GetTrackLabelNegative() == gamma1->GetTrackLabelNegative() ||
900                gamma0->GetTrackLabelNegative() == gamma1->GetTrackLabelPositive() ||
901                gamma0->GetTrackLabelPositive() == gamma1->GetTrackLabelNegative() ) continue;
902
903             AliAODConversionMother *pi0cand = new AliAODConversionMother(gamma0,gamma1);
904             pi0cand->SetLabels(firstGammaIndex,secondGammaIndex);
905
906             if((((AliConversionCuts*)fCutArray->At(fiCut))->MesonIsSelected(pi0cand,kTRUE))){
907                hESDMotherInvMassPt[fiCut]->Fill(pi0cand->M(),pi0cand->Pt());
908                if(pi0cand->GetAlpha()<0.1){
909                   hESDMotherInvMassEalpha[fiCut]->Fill(pi0cand->M(),pi0cand->E());
910                }
911                Int_t zbin= fBGHandler[fiCut]->GetZBinIndex(fESDEvent->GetPrimaryVertex()->GetZ());
912                Int_t mbin = 0;
913                if(((AliConversionCuts*)fCutArray->At(fiCut))->UseTrackMultiplicity()){
914                   mbin = fBGHandler[fiCut]->GetMultiplicityBinIndex(fNumberOfESDTracks);
915                } else {
916                   mbin = fBGHandler[fiCut]->GetMultiplicityBinIndex(fGoodGammas->GetEntries());
917                }
918                Double_t sparesFill[4] = {pi0cand->M(),pi0cand->Pt(),zbin,mbin};
919                sESDMotherInvMassPtZM[fiCut]->Fill(sparesFill,1);
920                if(fMCEvent){
921                   ProcessTrueMesonCandidates(pi0cand,gamma0,gamma1);
922                }
923             }
924             delete pi0cand;
925             pi0cand=0x0;
926          }
927       }
928    }
929 }
930 //______________________________________________________________________
931 void AliAnalysisTaskGammaConvV1::ProcessTrueMesonCandidates(AliAODConversionMother *Pi0Candidate, AliAODConversionPhoton *TrueGammaCandidate0, AliAODConversionPhoton *TrueGammaCandidate1)
932 {
933    // Process True Mesons
934    AliStack *MCStack = fMCEvent->Stack();
935
936    if(TrueGammaCandidate0->GetV0Index()<fESDEvent->GetNumberOfV0s()){
937       Bool_t isTruePi0 = kFALSE;
938       Bool_t isTrueEta = kFALSE;
939       Int_t gamma0MCLabel = TrueGammaCandidate0->GetMCParticleLabel(MCStack);
940       Int_t gamma0MotherLabel = -1;
941       if(gamma0MCLabel != -1){ // Gamma is Combinatorial; MC Particles don't belong to the same Mother
942          // Daughters Gamma 0
943          TParticle * negativeMC = (TParticle*)TrueGammaCandidate0->GetNegativeMCDaughter(MCStack);
944          TParticle * positiveMC = (TParticle*)TrueGammaCandidate0->GetPositiveMCDaughter(MCStack);
945          TParticle * gammaMC0 = (TParticle*)MCStack->Particle(gamma0MCLabel);
946          if(TMath::Abs(negativeMC->GetPdgCode())==11 && TMath::Abs(positiveMC->GetPdgCode())==11){  // Electrons ...
947             if(negativeMC->GetUniqueID() == 5 && positiveMC->GetUniqueID() ==5){ // ... From Conversion ...
948                if(gammaMC0->GetPdgCode() == 22){ // ... with Gamma Mother
949                   gamma0MotherLabel=gammaMC0->GetFirstMother();
950                }
951             }
952             if(gammaMC0->GetPdgCode() ==111){ // Conversion but Pi0 Mother
953                gamma0MotherLabel=-111;
954             }
955             if(gammaMC0->GetPdgCode() ==221){ // Conversion but Eta Mother
956                gamma0MotherLabel=-221;
957             }
958          }
959       }
960       if(TrueGammaCandidate1->GetV0Index()<fESDEvent->GetNumberOfV0s()){
961          Int_t gamma1MCLabel = TrueGammaCandidate1->GetMCParticleLabel(MCStack);
962          Int_t gamma1MotherLabel = -1;
963          if(gamma1MCLabel != -1){ // Gamma is Combinatorial; MC Particles don't belong to the same Mother
964             // Daughters Gamma 1
965             TParticle * negativeMC = (TParticle*)TrueGammaCandidate1->GetNegativeMCDaughter(MCStack);
966             TParticle * positiveMC = (TParticle*)TrueGammaCandidate1->GetPositiveMCDaughter(MCStack);
967             TParticle * gammaMC1 = (TParticle*)MCStack->Particle(gamma1MCLabel);
968             if(TMath::Abs(negativeMC->GetPdgCode())==11 && TMath::Abs(positiveMC->GetPdgCode())==11){  // Electrons ...
969                if(negativeMC->GetUniqueID() == 5 && positiveMC->GetUniqueID() ==5){ // ... From Conversion ...
970                   if(gammaMC1->GetPdgCode() == 22){ // ... with Gamma Mother
971                      gamma1MotherLabel=gammaMC1->GetFirstMother();
972                   }
973                }
974                if(gammaMC1->GetPdgCode() ==111){ // Conversion but Pi0 Mother
975                   gamma1MotherLabel=-111;
976                }
977                if(gammaMC1->GetPdgCode() ==221){ // Conversion but Eta Mother
978                   gamma1MotherLabel=-221;
979                }
980             }
981          }
982          if(gamma0MotherLabel>=0 && gamma0MotherLabel==gamma1MotherLabel){
983             if(((TParticle*)MCStack->Particle(gamma1MotherLabel))->GetPdgCode() == 111){
984                isTruePi0=kTRUE;
985             }
986             if(((TParticle*)MCStack->Particle(gamma1MotherLabel))->GetPdgCode() == 221){
987                isTrueEta=kTRUE;
988             }
989          }
990          if(isTruePi0 || isTrueEta){ // True Pion or Eta
991                                 if (!fMCEvent->IsFromBGEvent(gamma1MotherLabel) && fIsHeavyIon){
992                                         hESDTrueMotherInvMassPtAddedSig[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
993                                         if(gamma0MotherLabel > MCStack->GetNprimary()){ // Secondary Meson
994                                                 hESDTrueSecondaryMotherInvMassPtAddedSig[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
995                                                 if (((TParticle*)MCStack->Particle(gamma1MotherLabel))->GetMother(0) >-1){
996                                                         if(MCStack->Particle(((TParticle*)MCStack->Particle(gamma1MotherLabel))->GetMother(0))->GetPdgCode()==kK0Short){
997                                                                 hESDTrueSecondaryMotherFromK0sInvMassPtAddedSig[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
998                                                         }
999                                                 }
1000                                         }                       
1001                                         if(gamma0MotherLabel <= MCStack->GetNprimary()){ // Only primary pi0 for efficiency calculation
1002                                                 hESDTruePrimaryMotherInvMassMCPtAddedSig[fiCut]->Fill(Pi0Candidate->M(),((TParticle*)MCStack->Particle(gamma1MotherLabel))->Pt());
1003                                                 if(isTruePi0){ // Only primaries for unfolding
1004                                                         hESDTruePrimaryPi0ESDPtMCPtAddedSig[fiCut]->Fill(Pi0Candidate->Pt(),((TParticle*)MCStack->Particle(gamma1MotherLabel))->Pt());
1005                                                 }
1006                                         }
1007                                 } else {
1008                                         hESDTrueMotherInvMassPt[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1009                                         if(gamma0MotherLabel > MCStack->GetNprimary()){ // Secondary Meson
1010                                                 hESDTrueSecondaryMotherInvMassPt[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1011                                                 if (((TParticle*)MCStack->Particle(gamma1MotherLabel))->GetMother(0) >-1){
1012                                                         if(MCStack->Particle(((TParticle*)MCStack->Particle(gamma1MotherLabel))->GetMother(0))->GetPdgCode()==kK0Short){
1013                                                                 hESDTrueSecondaryMotherFromK0sInvMassPt[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1014                                                         }
1015                                                 }
1016                                         }
1017                                         if(gamma0MotherLabel <= MCStack->GetNprimary()){ // Only primary pi0 for efficiency calculation
1018                                                 hESDTruePrimaryMotherInvMassMCPt[fiCut]->Fill(Pi0Candidate->M(),((TParticle*)MCStack->Particle(gamma1MotherLabel))->Pt());
1019                                                 if(isTruePi0){ // Only primaries for unfolding
1020                                                         hESDTruePrimaryPi0ESDPtMCPt[fiCut]->Fill(Pi0Candidate->Pt(),((TParticle*)MCStack->Particle(gamma1MotherLabel))->Pt());
1021                                                 }
1022                                         }
1023                                 }
1024          }
1025          if(!isTruePi0 && !isTrueEta){ // Background
1026             if(gamma0MotherLabel>-1 && gamma1MotherLabel>-1){ // Both Tracks are Photons and have a mother but not Pi0 or Eta
1027                hESDTrueBckGGInvMassPt[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1028             } else { // No photon or without mother
1029                hESDTrueBckContInvMassPt[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1030             }
1031             if((gamma0MotherLabel==-111 || gamma1MotherLabel==-111 || gamma0MotherLabel==-221 || gamma1MotherLabel==-221) && (!fMCEvent->IsFromBGEvent(gamma1MotherLabel) || !fMCEvent->IsFromBGEvent(gamma0MotherLabel)) && fIsHeavyIon){
1032                // Dalitz
1033                hESDTrueMotherDalitzInvMassPtAddedSig[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1034             } else if (gamma0MotherLabel==-111 || gamma1MotherLabel==-111 || gamma0MotherLabel==-221 || gamma1MotherLabel==-221){
1035                                         // Dalitz
1036                hESDTrueMotherDalitzInvMassPt[fiCut]->Fill(Pi0Candidate->M(),Pi0Candidate->Pt());
1037             }
1038          }
1039       }
1040    }
1041 }
1042
1043 //________________________________________________________________________
1044 void AliAnalysisTaskGammaConvV1::CalculateBackground(){
1045
1046    Int_t zbin= fBGHandler[fiCut]->GetZBinIndex(fESDEvent->GetPrimaryVertex()->GetZ());
1047    Int_t mbin = 0;
1048
1049    if(((AliConversionCuts*)fCutArray->At(fiCut))->UseTrackMultiplicity()){
1050       mbin = fBGHandler[fiCut]->GetMultiplicityBinIndex(fNumberOfESDTracks);
1051    } else {
1052       mbin = fBGHandler[fiCut]->GetMultiplicityBinIndex(fGoodGammas->GetEntries());
1053    }
1054
1055    if(((AliConversionCuts*)fCutArray->At(fiCut))->UseRotationMethod()){
1056
1057       for(Int_t iCurrent=0;iCurrent<fGoodGammas->GetEntries();iCurrent++){
1058          AliAODConversionPhoton currentEventGoodV0 = *(AliAODConversionPhoton*)(fGoodGammas->At(iCurrent));
1059          for(Int_t iCurrent2=iCurrent+1;iCurrent2<fGoodGammas->GetEntries();iCurrent2++){
1060             for(Int_t nRandom=0;nRandom<((AliConversionCuts*)fCutArray->At(fiCut))->NumberOfRotationEvents();nRandom++){
1061                AliAODConversionPhoton currentEventGoodV02 = *(AliAODConversionPhoton*)(fGoodGammas->At(iCurrent2));
1062
1063                if(((AliConversionCuts*)fV0Reader->GetConversionCuts())->DoBGProbability()){
1064                   AliAODConversionMother *backgroundCandidateProb = new AliAODConversionMother(&currentEventGoodV0,&currentEventGoodV02);
1065                   Double_t massBGprob = backgroundCandidateProb->M();
1066                   if(massBGprob>0.1 && massBGprob<0.14){
1067                      if(fRandom.Rndm()>fBGHandler[fiCut]->GetBGProb(zbin,mbin)){
1068                         delete backgroundCandidateProb;
1069                         continue;
1070                      }
1071                   }
1072                   delete backgroundCandidateProb;
1073                   backgroundCandidateProb = 0x0;
1074                }
1075
1076                RotateParticle(&currentEventGoodV02);
1077                AliAODConversionMother *backgroundCandidate = new AliAODConversionMother(&currentEventGoodV0,&currentEventGoodV02);
1078                if((((AliConversionCuts*)fCutArray->At(fiCut))->MesonIsSelected(backgroundCandidate,kFALSE))){
1079                   hESDMotherBackInvMassPt[fiCut]->Fill(backgroundCandidate->M(),backgroundCandidate->Pt());
1080                   Double_t sparesFill[4] = {backgroundCandidate->M(),backgroundCandidate->Pt(),zbin,mbin};
1081                   sESDMotherBackInvMassPtZM[fiCut]->Fill(sparesFill,1);
1082                }
1083                delete backgroundCandidate;
1084                backgroundCandidate = 0x0;
1085             }
1086          }
1087       }
1088    }
1089    else{
1090       AliGammaConversionAODBGHandler::GammaConversionVertex *bgEventVertex = NULL;
1091
1092       if(((AliConversionCuts*)fCutArray->At(fiCut))->UseTrackMultiplicity()){
1093          for(Int_t nEventsInBG=0;nEventsInBG<fBGHandler[fiCut]->GetNBGEvents();nEventsInBG++){
1094             AliGammaConversionAODVector *previousEventV0s = fBGHandler[fiCut]->GetBGGoodV0s(zbin,mbin,nEventsInBG);
1095             if(fMoveParticleAccordingToVertex == kTRUE){
1096                bgEventVertex = fBGHandler[fiCut]->GetBGEventVertex(zbin,mbin,nEventsInBG);
1097             }
1098
1099             for(Int_t iCurrent=0;iCurrent<fGoodGammas->GetEntries();iCurrent++){
1100                AliAODConversionPhoton currentEventGoodV0 = *(AliAODConversionPhoton*)(fGoodGammas->At(iCurrent));
1101                for(UInt_t iPrevious=0;iPrevious<previousEventV0s->size();iPrevious++){
1102                   AliAODConversionPhoton previousGoodV0 = (AliAODConversionPhoton)(*(previousEventV0s->at(iPrevious)));
1103                   if(fMoveParticleAccordingToVertex == kTRUE){
1104                      MoveParticleAccordingToVertex(&previousGoodV0,bgEventVertex);
1105                   }
1106
1107                   AliAODConversionMother *backgroundCandidate = new AliAODConversionMother(&currentEventGoodV0,&previousGoodV0);
1108                   if((((AliConversionCuts*)fCutArray->At(fiCut))->MesonIsSelected(backgroundCandidate,kFALSE))){
1109                      hESDMotherBackInvMassPt[fiCut]->Fill(backgroundCandidate->M(),backgroundCandidate->Pt());
1110                      Double_t sparesFill[4] = {backgroundCandidate->M(),backgroundCandidate->Pt(),zbin,mbin};
1111                      sESDMotherBackInvMassPtZM[fiCut]->Fill(sparesFill,1);
1112                   }
1113                   delete backgroundCandidate;
1114                   backgroundCandidate = 0x0;
1115                }
1116             }
1117          }
1118       }
1119       else{
1120          for(Int_t nEventsInBG=0;nEventsInBG <fBGHandler[fiCut]->GetNBGEvents();nEventsInBG++){
1121             AliGammaConversionAODVector *previousEventV0s = fBGHandler[fiCut]->GetBGGoodV0s(zbin,mbin,nEventsInBG);
1122             if(previousEventV0s){
1123                if(fMoveParticleAccordingToVertex == kTRUE){
1124                   bgEventVertex = fBGHandler[fiCut]->GetBGEventVertex(zbin,mbin,nEventsInBG);
1125                }
1126                for(Int_t iCurrent=0;iCurrent<fGoodGammas->GetEntries();iCurrent++){
1127                   AliAODConversionPhoton currentEventGoodV0 = *(AliAODConversionPhoton*)(fGoodGammas->At(iCurrent));
1128                   for(UInt_t iPrevious=0;iPrevious<previousEventV0s->size();iPrevious++){
1129
1130                      AliAODConversionPhoton previousGoodV0 = (AliAODConversionPhoton)(*(previousEventV0s->at(iPrevious)));
1131
1132                      if(fMoveParticleAccordingToVertex == kTRUE){
1133                         MoveParticleAccordingToVertex(&previousGoodV0,bgEventVertex);
1134                      }
1135
1136                      AliAODConversionMother *backgroundCandidate = new AliAODConversionMother(&currentEventGoodV0,&previousGoodV0);
1137
1138                      if((((AliConversionCuts*)fCutArray->At(fiCut))->MesonIsSelected(backgroundCandidate,kFALSE))){
1139                         hESDMotherBackInvMassPt[fiCut]->Fill(backgroundCandidate->M(),backgroundCandidate->Pt());
1140                         Double_t sparesFill[4] = {backgroundCandidate->M(),backgroundCandidate->Pt(),zbin,mbin};
1141                         sESDMotherBackInvMassPtZM[fiCut]->Fill(sparesFill,1);
1142                      }
1143                      delete backgroundCandidate;
1144                      backgroundCandidate = 0x0;
1145                   }
1146                }
1147             }
1148          }
1149       }
1150    }
1151 }
1152 //________________________________________________________________________
1153 void AliAnalysisTaskGammaConvV1::RotateParticle(AliAODConversionPhoton *gamma){
1154    Int_t fNDegreesPMBackground= ((AliConversionCuts*)fV0Reader->GetConversionCuts())->NDegreesRotation();
1155    Double_t nRadiansPM = fNDegreesPMBackground*TMath::Pi()/180;
1156    Double_t rotationValue = fRandom.Rndm()*2*nRadiansPM + TMath::Pi()-nRadiansPM;
1157    gamma->RotateZ(rotationValue);
1158 }
1159 //________________________________________________________________________
1160 void AliAnalysisTaskGammaConvV1::MoveParticleAccordingToVertex(AliAODConversionPhoton* particle,const AliGammaConversionAODBGHandler::GammaConversionVertex *vertex){
1161    //see header file for documentation
1162
1163    Double_t dx = vertex->fX - fESDEvent->GetPrimaryVertex()->GetX();
1164    Double_t dy = vertex->fY - fESDEvent->GetPrimaryVertex()->GetY();
1165    Double_t dz = vertex->fZ - fESDEvent->GetPrimaryVertex()->GetZ();
1166
1167    Double_t movedPlace[3] = {particle->GetConversionX() - dx,particle->GetConversionY() - dy,particle->GetConversionZ() - dz};
1168    particle->SetConversionPoint(movedPlace);
1169 }
1170 //________________________________________________________________________
1171 void AliAnalysisTaskGammaConvV1::UpdateEventByEventData(){
1172    //see header file for documentation
1173    if(fGoodGammas->GetEntries() >0 ){
1174       if(((AliConversionCuts*)fCutArray->At(fiCut))->UseTrackMultiplicity()){
1175          fBGHandler[fiCut]->AddEvent(fGoodGammas,fESDEvent->GetPrimaryVertex()->GetX(),fESDEvent->GetPrimaryVertex()->GetY(),fESDEvent->GetPrimaryVertex()->GetZ(),fNumberOfESDTracks);
1176       }
1177       else{ // means we use #V0s for multiplicity
1178          fBGHandler[fiCut]->AddEvent(fGoodGammas,fESDEvent->GetPrimaryVertex()->GetX(),fESDEvent->GetPrimaryVertex()->GetY(),fESDEvent->GetPrimaryVertex()->GetZ(),fGoodGammas->GetEntries());
1179       }
1180    }
1181 }
1182
1183 //________________________________________________________________________
1184 void AliAnalysisTaskGammaConvV1::CountESDTracks(){
1185
1186    AliESDtrackCuts *EsdTrackCuts = new AliESDtrackCuts("AliESDtrackCuts");
1187    // Using standard function for setting Cuts
1188    Bool_t selectPrimaries=kTRUE;
1189    EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
1190    EsdTrackCuts->SetMaxDCAToVertexZ(2);
1191    EsdTrackCuts->SetEtaRange(-0.8, 0.8);
1192    EsdTrackCuts->SetPtRange(0.15);
1193
1194    fNumberOfESDTracks = 0;
1195    for(Int_t iTracks = 0; iTracks < fESDEvent->GetNumberOfTracks(); iTracks++){
1196       AliESDtrack* curTrack = fESDEvent->GetTrack(iTracks);
1197       if(!curTrack) continue;
1198       if(EsdTrackCuts->AcceptTrack(curTrack) ) fNumberOfESDTracks++;
1199    }
1200    delete EsdTrackCuts;
1201    EsdTrackCuts=0x0;
1202
1203    return;
1204 }
1205 //________________________________________________________________________
1206 void AliAnalysisTaskGammaConvV1::Terminate(const Option_t *)
1207 {
1208    fOutputContainer->Add(((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetCutHistograms());
1209
1210    for(Int_t iCut = 0; iCut<fnCuts;iCut++){
1211       if(((AliConversionCuts*)fCutArray->At(iCut))->GetCutHistograms()){
1212          fCutFolder[iCut]->Add(((AliConversionCuts*)fCutArray->At(iCut))->GetCutHistograms());
1213       }
1214    }
1215    fOutputContainer->Print();
1216 }
1217