]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliConversionCuts.cxx
78a34026fcc943e9e2e89933722efc64482e6f63
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliConversionCuts.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 // Class handling all kinds of selection cuts for
19 // Gamma Conversion analysis
20 //---------------------------------------------
21 ////////////////////////////////////////////////
22
23 #include "AliConversionCuts.h"
24
25 #include "AliKFVertex.h"
26 #include "AliAODTrack.h"
27 #include "AliESDtrack.h"
28 #include "AliAnalysisManager.h"
29 #include "AliInputEventHandler.h"
30 #include "AliMCEventHandler.h"
31 #include "AliAODHandler.h"
32 #include "AliPIDResponse.h"
33 #include "TH1.h"
34 #include "TH2.h"
35 #include "AliStack.h"
36 #include "AliAODConversionMother.h"
37 #include "TObjString.h"
38 #include "AliAODEvent.h"
39 #include "AliESDEvent.h"
40 #include "AliCentrality.h"
41 #include "TList.h"
42 #include "AliLog.h"
43
44 class iostream;
45
46 using namespace std;
47
48 ClassImp(AliConversionCuts)
49
50
51 const char* AliConversionCuts::fgkCutNames[AliConversionCuts::kNCuts] = {
52 "GoodId",
53 "V0FinderType",
54 "eProbCut",
55 "ededxSigmaCut",
56 "pidedxSigmaCut",
57 "piMomdedxSigmaCut",
58 "Chi2GammaCut",
59 "SinglePtCut",
60 "ClsTPCCut",
61 "EtaCut",
62 "Chi2MesonCut",
63 "LowPRejectionSigmaCut",
64 "QtMaxCut",
65 "piMaxMomdedxSigmaCut",
66 "AlphaMesonCut",
67 "MinRCut",
68 "RapidityMesonCut",
69 "BackgroundScheme",
70 "DegreesForRotationMethod",
71 "NumberOfRotations",
72 "RemovePileUp",
73 "SelectV0AND",
74 "MultiplicityMethod",
75 "HeavyIon",
76 "CentralityMin",
77 "CentralityMax",
78 "TOFelectronPID",
79 "UseMCPSmearing",
80 "DoPhotonAsymmetryCut",
81 "PsiPair",
82 "CosinePointingAngle",
83 "SharedElectronCuts",
84 "RejectToCloseV0s",
85 };
86
87
88 //________________________________________________________________________
89 AliConversionCuts::AliConversionCuts(const char *name,const char *title) : AliAnalysisCuts(name,title),
90     fHistograms(NULL),
91     fPIDResponse(NULL),
92     fEventQuality(-1),
93     fMaxR(200),
94     fMinR(0),
95     fEtaCut(0.9),
96     fEtaCutMin(-0.1),
97     fPtCut(0.02),
98     fSinglePtCut(0),
99     fMaxZ(1000),
100     fMinClsTPC(0.),
101     fMinClsTPCToF(0.),
102     fLineCutZRSlope(0.),
103     fLineCutZValue(0),
104     fLineCutZRSlopeMin(0.),
105     fLineCutZValueMin(0),
106     fChi2CutConversion(1000),
107     fChi2CutMeson(1000),
108     fPIDProbabilityCutNegativeParticle(0),
109     fPIDProbabilityCutPositiveParticle(0),
110     fDodEdxSigmaCut(kTRUE),
111     fDoTOFsigmaCut(kFALSE), // RRnewTOF
112     fPIDTRDEfficiency(1),
113     fDoTRDPID(kFALSE),
114     fPIDnSigmaAboveElectronLine(100),
115     fPIDnSigmaBelowElectronLine(-100),
116     fTofPIDnSigmaAboveElectronLine(100), // RRnewTOF
117     fTofPIDnSigmaBelowElectronLine(-100), // RRnewTOF
118     fPIDnSigmaAbovePionLine(0),
119     fPIDnSigmaAbovePionLineHighPt(-100),
120     fPIDMinPnSigmaAbovePionLine(0),
121     fPIDMaxPnSigmaAbovePionLine(0),
122     fDoKaonRejectionLowP(kFALSE),
123     fDoProtonRejectionLowP(kFALSE),
124     fDoPionRejectionLowP(kFALSE),
125     fPIDnSigmaAtLowPAroundKaonLine(0),
126     fPIDnSigmaAtLowPAroundProtonLine(0),
127     fPIDnSigmaAtLowPAroundPionLine(0),
128     fPIDMinPKaonRejectionLowP(0),
129     fPIDMinPProtonRejectionLowP(0),
130     fPIDMinPPionRejectionLowP(0),
131     fDoQtGammaSelection(kTRUE),
132     fDoHighPtQtGammaSelection(kFALSE), // RRnew
133     fQtMax(100),
134     fHighPtQtMax(0.), // RRnew
135     fPtBorderForQt(0), // RRnew
136     fXVertexCut(0.),
137     fYVertexCut(0.),
138     fZVertexCut(0.),
139     fNSigmaMass(0.),
140     fUseEtaMinCut(kFALSE),
141     fUseOnFlyV0Finder(kTRUE),
142     fDoPhotonAsymmetryCut(kTRUE),
143     fMinPPhotonAsymmetryCut(100.),
144     fMinPhotonAsymmetry(0.),
145     fIsHeavyIon(kFALSE),
146     fDetectorCentrality(0),
147     fModCentralityClass(0),
148     fMaxVertexZ(10),
149     fCentralityMin(0),
150     fCentralityMax(0),
151     fUseCorrectedTPCClsInfo(kFALSE),
152     fUseTOFpid(kFALSE),
153     fAlphaMinCutMeson(0),
154     fAlphaCutMeson(1),
155     fRapidityCutMeson(1),
156     fUseRotationMethodInBG(kFALSE),
157     fdoBGProbability(kFALSE),
158     fUseTrackMultiplicityForBG(kFALSE),
159     fnDegreeRotationPMForBG(0),
160     fnumberOfRotationEventsForBG(0),
161     fUseMCPSmearing(kFALSE),
162     fPBremSmearing(0),
163     fPSigSmearing(0),
164     fPSigSmearingCte(0),
165     fBrem(NULL),
166     fMultiplicityMethod(0),
167     fSelectV0AND(kFALSE),
168     fRemovePileUp(kFALSE),
169     fOpeningAngle(0.005),
170     fPsiPairCut(10000),
171     fCosPAngleCut(10000),
172     fDoToCloseV0sCut(kFALSE),
173     fminV0Dist(200.),
174     fDoSharedElecCut(kFALSE),
175     fOfflineTriggerMask(0),
176     fRandom(0),
177     fElectronLabelArray(NULL),
178     fCutString(NULL),
179     hdEdxCuts(NULL),
180     hTPCdEdxbefore(NULL),
181     hTPCdEdxafter(NULL),
182     hTOFbefore(NULL),
183     hTOFafter(NULL),
184     hTrackCuts(NULL),
185     hPhotonCuts(NULL),
186     hInvMassbefore(NULL),
187     hArmenterosbefore(NULL),
188     hInvMassafter(NULL),
189     hArmenterosafter(NULL),
190     hAcceptanceCuts(NULL),
191     hCutIndex(NULL),
192     hV0EventCuts(NULL),
193     hCentrality(NULL),
194     hVertexZ(NULL),
195     hTriggerClass(NULL),
196     hMesonCuts(NULL),
197     hMesonBGCuts(NULL)
198 {
199     InitPIDResponse();
200     for(Int_t jj=0;jj<kNCuts;jj++){fCuts[jj]=0;}
201     fCutString=new TObjString((GetCutNumber()).Data());
202
203     if (fBrem == NULL){
204        fBrem = new TF1("fBrem","pow(-log(x),[0]/log(2.0)-1.0)/TMath::Gamma([0]/log(2.0))",0.00001,0.999999999);
205        // tests done with 1.0e-14
206        fBrem->SetParameter(0,fPBremSmearing);
207        fBrem->SetNpx(100000);
208     }
209
210     fElectronLabelArray = new Int_t[500];
211 }
212
213 //________________________________________________________________________
214 AliConversionCuts::~AliConversionCuts() {
215     // Destructor
216   //Deleting fHistograms leads to seg fault it it's added to output collection of a task
217   // if(fHistograms)
218   //    delete fHistograms;
219   // fHistograms = NULL;
220    if(fCutString != NULL){
221       delete fCutString;
222       fCutString = NULL;
223    }
224    if(fBrem){
225       delete fBrem;
226       fBrem = NULL;
227    }
228    if(fElectronLabelArray){
229       delete fElectronLabelArray;
230       fElectronLabelArray = NULL;
231    }
232
233 }
234
235 //________________________________________________________________________
236 void AliConversionCuts::InitCutHistograms(TString name, Bool_t preCut){
237
238     // Initialize Cut Histograms for QA (only initialized and filled if function is called)
239
240     if(fHistograms != NULL){
241         delete fHistograms;
242         fHistograms=NULL;
243     }
244     if(fHistograms==NULL){
245         fHistograms=new TList();
246         if(name=="")fHistograms->SetName(Form("ConvCuts_%s",GetCutNumber().Data()));
247         else fHistograms->SetName(Form("%s_%s",name.Data(),GetCutNumber().Data()));
248     }
249
250     // IsPhotonSelected
251     hCutIndex=new TH1F(Form("IsPhotonSelected %s",GetCutNumber().Data()),"IsPhotonSelected",10,-0.5,9.5);
252     hCutIndex->GetXaxis()->SetBinLabel(kPhotonIn+1,"in");
253     hCutIndex->GetXaxis()->SetBinLabel(kOnFly+1,"onfly");
254     hCutIndex->GetXaxis()->SetBinLabel(kNoTracks+1,"no tracks");
255     hCutIndex->GetXaxis()->SetBinLabel(kdEdxCuts+1,"dEdx");
256     hCutIndex->GetXaxis()->SetBinLabel(kTrackCuts+1,"Track cuts");
257     hCutIndex->GetXaxis()->SetBinLabel(kConvPointFail+1,"ConvPoint fail");
258     hCutIndex->GetXaxis()->SetBinLabel(kPhotonCuts+1,"PhotonCuts");
259     hCutIndex->GetXaxis()->SetBinLabel(kPhotonOut+1,"out");
260     fHistograms->Add(hCutIndex);
261
262     // Track Cuts
263     hTrackCuts=new TH1F(Form("TrackCuts %s",GetCutNumber().Data()),"TrackCuts",10,-0.5,9.5);
264     hTrackCuts->GetXaxis()->SetBinLabel(1,"in");
265     hTrackCuts->GetXaxis()->SetBinLabel(2,"likesign");
266     hTrackCuts->GetXaxis()->SetBinLabel(3,"ntpccl");
267     hTrackCuts->GetXaxis()->SetBinLabel(4,"acceptance");
268     hTrackCuts->GetXaxis()->SetBinLabel(5,"singlept");
269     hTrackCuts->GetXaxis()->SetBinLabel(6,"TPCrefit");
270     hTrackCuts->GetXaxis()->SetBinLabel(7,"kink");
271     hTrackCuts->GetXaxis()->SetBinLabel(8,"out");
272     fHistograms->Add(hTrackCuts);
273
274     // Photon Cuts
275     hPhotonCuts=new TH1F(Form("PhotonCuts %s",GetCutNumber().Data()),"PhotonCuts",12,-0.5,11.5);
276     hPhotonCuts->GetXaxis()->SetBinLabel(1,"in");
277     hPhotonCuts->GetXaxis()->SetBinLabel(2,"qtcut");
278     hPhotonCuts->GetXaxis()->SetBinLabel(3,"chi2");
279     hPhotonCuts->GetXaxis()->SetBinLabel(4,"acceptance");
280     hPhotonCuts->GetXaxis()->SetBinLabel(5,"asymmetry");
281     hPhotonCuts->GetXaxis()->SetBinLabel(6,"pidprob");
282     hPhotonCuts->GetXaxis()->SetBinLabel(7,"cortpcclinfo");
283     hPhotonCuts->GetXaxis()->SetBinLabel(8,"PsiPair");
284     hPhotonCuts->GetXaxis()->SetBinLabel(9,"CosPAngle");
285     hPhotonCuts->GetXaxis()->SetBinLabel(10,"out");
286     fHistograms->Add(hPhotonCuts);
287
288     if(preCut){
289        hInvMassbefore=new TH1F(Form("InvMass_before %s",GetCutNumber().Data()),"InvMass_before",100,0,0.3);
290        fHistograms->Add(hInvMassbefore);
291        hArmenterosbefore=new TH2F(Form("Armenteros_before %s",GetCutNumber().Data()),"Armenteros_before",200,-1,1,250,0,0.25);
292        fHistograms->Add(hArmenterosbefore);
293     }
294     hInvMassafter=new TH1F(Form("InvMass_after %s",GetCutNumber().Data()),"InvMass_after",100,0,0.3);
295     fHistograms->Add(hInvMassafter);
296     hArmenterosafter=new TH2F(Form("Armenteros_after %s",GetCutNumber().Data()),"Armenteros_after",200,-1,1,250,0,0.25);
297     fHistograms->Add(hArmenterosafter);
298
299     hAcceptanceCuts=new TH1F(Form("PhotonAcceptanceCuts %s",GetCutNumber().Data()),"PhotonAcceptanceCuts",10,-0.5,9.5);
300     hAcceptanceCuts->GetXaxis()->SetBinLabel(1,"in");
301     hAcceptanceCuts->GetXaxis()->SetBinLabel(2,"maxR");
302     hAcceptanceCuts->GetXaxis()->SetBinLabel(3,"minR");
303     hAcceptanceCuts->GetXaxis()->SetBinLabel(4,"line");
304     hAcceptanceCuts->GetXaxis()->SetBinLabel(5,"maxZ");
305     hAcceptanceCuts->GetXaxis()->SetBinLabel(6,"eta");
306     hAcceptanceCuts->GetXaxis()->SetBinLabel(7,"minpt");
307     hAcceptanceCuts->GetXaxis()->SetBinLabel(8,"out");
308     fHistograms->Add(hAcceptanceCuts);
309
310     // dEdx Cuts
311     hdEdxCuts=new TH1F(Form("dEdxCuts %s",GetCutNumber().Data()),"dEdxCuts",10,-0.5,9.5);
312     hdEdxCuts->GetXaxis()->SetBinLabel(1,"in");
313     hdEdxCuts->GetXaxis()->SetBinLabel(2,"TPCelectron");
314     hdEdxCuts->GetXaxis()->SetBinLabel(3,"TPCpion");
315     hdEdxCuts->GetXaxis()->SetBinLabel(4,"TPCpionhighp");
316     hdEdxCuts->GetXaxis()->SetBinLabel(5,"TPCkaonlowprej");
317     hdEdxCuts->GetXaxis()->SetBinLabel(6,"TPCprotonlowprej");
318     hdEdxCuts->GetXaxis()->SetBinLabel(7,"TPCpionlowprej");
319     hdEdxCuts->GetXaxis()->SetBinLabel(8,"TOFelectron");
320     hdEdxCuts->GetXaxis()->SetBinLabel(9,"TRDelectron");
321     hdEdxCuts->GetXaxis()->SetBinLabel(10,"out");
322     fHistograms->Add(hdEdxCuts);
323     
324     TAxis *AxisBeforedEdx = NULL;
325     TAxis *AxisBeforeTOF = NULL;
326     if(preCut){
327        hTPCdEdxbefore=new TH2F(Form("Gamma_dEdx_before %s",GetCutNumber().Data()),"dEdx Gamma before" ,150,0.05,20,400,-10,10);
328        fHistograms->Add(hTPCdEdxbefore);
329        AxisBeforedEdx = hTPCdEdxbefore->GetXaxis(); 
330
331        hTOFbefore=new TH2F(Form("Gamma_TOF_before %s",GetCutNumber().Data()),"TOF Gamma before" ,150,0.05,20,400,-6,10);
332        fHistograms->Add(hTOFbefore);
333        AxisBeforeTOF = hTOFbefore->GetXaxis(); 
334     }
335     hTPCdEdxafter=new TH2F(Form("Gamma_dEdx_after %s",GetCutNumber().Data()),"dEdx Gamma after" ,150,0.05,20,400, -10,10);
336     fHistograms->Add(hTPCdEdxafter);
337
338     hTOFafter=new TH2F(Form("Gamma_TOF_after %s",GetCutNumber().Data()),"TOF Gamma after" ,150,0.05,20,400,-6,10);
339     fHistograms->Add(hTOFafter);
340
341     TAxis *AxisAfter = hTPCdEdxafter->GetXaxis(); 
342     Int_t bins = AxisAfter->GetNbins();
343     Double_t from = AxisAfter->GetXmin();
344     Double_t to = AxisAfter->GetXmax();
345     Double_t *newBins = new Double_t[bins+1];
346     newBins[0] = from;
347     Double_t factor = TMath::Power(to/from, 1./bins);
348     for(Int_t i=1; i<=bins; ++i) newBins[i] = factor * newBins[i-1];
349     AxisAfter->Set(bins, newBins);
350     AxisAfter = hTOFafter->GetXaxis(); 
351     AxisAfter->Set(bins, newBins);
352     if(preCut){
353        AxisBeforedEdx->Set(bins, newBins);
354        AxisBeforeTOF->Set(bins, newBins);
355     }
356     delete [] newBins;
357         
358     // Event Cuts and Info
359     if(preCut){
360        hV0EventCuts=new TH1F(Form("ESD_EventCuts %s",GetCutNumber().Data()),"Event Cuts",10,-0.5,9.5);
361        hV0EventCuts->GetXaxis()->SetBinLabel(1,"in");
362        hV0EventCuts->GetXaxis()->SetBinLabel(2,"OfflineTrigger");
363        hV0EventCuts->GetXaxis()->SetBinLabel(3,"VertexZ");
364        hV0EventCuts->GetXaxis()->SetBinLabel(4,"nvtxcontr");
365        hV0EventCuts->GetXaxis()->SetBinLabel(5,"pileup");
366        hV0EventCuts->GetXaxis()->SetBinLabel(6,"centrsel");
367        hV0EventCuts->GetXaxis()->SetBinLabel(7,"out");
368        fHistograms->Add(hV0EventCuts);
369        
370        hCentrality=new TH1F(Form("Centrality %s",GetCutNumber().Data()),"Centrality",100,0,100);
371        fHistograms->Add(hCentrality);
372        hVertexZ=new TH1F(Form("VertexZ %s",GetCutNumber().Data()),"VertexZ",1000,-50,50);
373        fHistograms->Add(hVertexZ);
374        
375        hTriggerClass= new TH1F(Form("OfflineTrigger %s",GetCutNumber().Data()),"OfflineTrigger",4,-0.5,3.5);
376        hTriggerClass->GetXaxis()->SetBinLabel(1,"kAny");
377        hTriggerClass->GetXaxis()->SetBinLabel(2,"kMB");
378        hTriggerClass->GetXaxis()->SetBinLabel(3,"kCentral");
379        hTriggerClass->GetXaxis()->SetBinLabel(4,"kSemiCentral");
380        fHistograms->Add(hTriggerClass);
381     }
382     
383     // Meson Cuts
384     hMesonCuts=new TH1F(Form("MesonCuts %s",GetCutNumber().Data()),"MesonCuts",10,-0.5,9.5);
385     hMesonCuts->GetXaxis()->SetBinLabel(1,"in");
386     hMesonCuts->GetXaxis()->SetBinLabel(2,"undef rapidity");
387     hMesonCuts->GetXaxis()->SetBinLabel(3,"rapidity cut");
388     hMesonCuts->GetXaxis()->SetBinLabel(4,"opening angle");
389     hMesonCuts->GetXaxis()->SetBinLabel(5,"alpha max");
390     hMesonCuts->GetXaxis()->SetBinLabel(6,"alpha min");
391     hMesonCuts->GetXaxis()->SetBinLabel(7,"out");
392     fHistograms->Add(hMesonCuts);
393
394     hMesonBGCuts=new TH1F(Form("MesonBGCuts %s",GetCutNumber().Data()),"MesonBGCuts",10,-0.5,9.5);
395     hMesonBGCuts->GetXaxis()->SetBinLabel(1,"in");
396     hMesonBGCuts->GetXaxis()->SetBinLabel(2,"undef rapidity");
397     hMesonBGCuts->GetXaxis()->SetBinLabel(3,"rapidity cut");
398     hMesonBGCuts->GetXaxis()->SetBinLabel(4,"opening angle");
399     hMesonBGCuts->GetXaxis()->SetBinLabel(5,"alpha max");
400     hMesonBGCuts->GetXaxis()->SetBinLabel(6,"alpha min");
401     hMesonBGCuts->GetXaxis()->SetBinLabel(7,"out");
402     fHistograms->Add(hMesonBGCuts);
403 }
404
405 //________________________________________________________________________
406 Bool_t AliConversionCuts::InitPIDResponse(){
407     // Set Pointer to AliPIDResponse
408
409   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
410   if(man) { 
411         AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
412     fPIDResponse = (AliPIDResponse*)inputHandler->GetPIDResponse();
413     if(fPIDResponse)return kTRUE;
414     
415   }
416   
417   return kFALSE;
418 }
419 ///________________________________________________________________________
420 Bool_t AliConversionCuts::EventIsSelected(AliVEvent *fInputEvent, AliVEvent *fMCEvent){
421     // Process Event Selection
422
423     Int_t cutindex=0;
424     if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
425     cutindex++;
426     
427     // Check for MC event
428     if(fMCEvent){
429        // Check if MC event is correctly loaded
430        AliMCEventHandler* mcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
431        if (!mcHandler){
432           fEventQuality = 2;
433           return kFALSE;
434        }
435        if (!mcHandler->InitOk() ){
436           fEventQuality = 2;
437           return kFALSE;
438        }
439        if (!mcHandler->TreeK() ){
440           fEventQuality = 2;
441           return kFALSE;
442        }
443        if (!mcHandler->TreeTR() ) {
444           fEventQuality = 2;
445           return kFALSE;
446        }
447     }
448     
449     // Event Trigger
450
451     if(!IsTriggerSelected()){
452         if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
453         fEventQuality = 3;
454         return kFALSE;
455     }
456     cutindex++;
457
458     // Z Vertex Position Cut
459     if(!VertexZCut(fInputEvent)){
460         if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
461         fEventQuality = 4;
462         return kFALSE;
463     }
464     cutindex++;
465
466     // Number of Contributors Cut
467     if(GetNumberOfContributorsVtx(fInputEvent)<=0) {
468         if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
469         fEventQuality = 5;
470         return kFALSE;
471     }
472     cutindex++;
473
474     // Pile Up Rejection
475
476     if(fRemovePileUp){
477         if(fInputEvent->IsPileupFromSPD(3,0.8,3.,2.,5.)){
478             if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
479             fEventQuality = 6;
480             return kFALSE;
481         }
482     }
483     cutindex++;
484
485     // Centrality Selection
486     if(!IsCentralitySelected(fInputEvent)){
487         if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
488         fEventQuality = 1;
489         return kFALSE;
490     }
491     cutindex++;
492
493     // Fill Event Histograms
494     if(hV0EventCuts)hV0EventCuts->Fill(cutindex);
495     if(hVertexZ)hVertexZ->Fill(fInputEvent->GetPrimaryVertex()->GetZ());
496     if(hCentrality)hCentrality->Fill(GetCentrality(fInputEvent));
497
498     fEventQuality = 0;
499     return kTRUE;
500 }
501
502 ///________________________________________________________________________
503 Bool_t AliConversionCuts::PhotonIsSelectedMC(TParticle *particle,AliStack *fMCStack,Bool_t checkForConvertedGamma){
504     // MonteCarlo Photon Selection
505
506     if(!fMCStack)return kFALSE;
507
508     if (particle->GetPdgCode() == 22){
509               
510         if(particle->R() > fMaxR)       return kFALSE;
511         if(TMath::Abs(particle->Eta())> fEtaCut || TMath::Abs(particle->Eta())< fEtaCutMin)     return kFALSE;
512
513         if(particle->GetMother(0) >-1 && fMCStack->Particle(particle->GetMother(0))->GetPdgCode() == 22){
514             return kFALSE; // no photon as mothers!
515         }
516
517         if(particle->GetMother(0) >= fMCStack->GetNprimary()){
518             return kFALSE; // the gamma has a mother, and it is not a primary particle
519         }
520
521         if(!checkForConvertedGamma) return kTRUE; // return in case of accepted gamma
522
523         // looking for conversion gammas (electron + positron from pairbuilding (= 5) )
524         TParticle* ePos = NULL;
525         TParticle* eNeg = NULL;
526
527         if(particle->GetNDaughters() >= 2){
528             for(Int_t daughterIndex=particle->GetFirstDaughter();daughterIndex<=particle->GetLastDaughter();daughterIndex++){
529                 TParticle *tmpDaughter = fMCStack->Particle(daughterIndex);
530                 if(tmpDaughter->GetUniqueID() == 5){
531                     if(tmpDaughter->GetPdgCode() == 11){
532                         eNeg = tmpDaughter;
533                     } else if(tmpDaughter->GetPdgCode() == -11){
534                         ePos = tmpDaughter;
535                     }
536                 }
537             }
538         }
539
540         if(ePos == NULL || eNeg == NULL){ // means we do not have two daughters from pair production
541             return kFALSE;
542         }
543         
544         if(ePos->Pt()<fSinglePtCut || eNeg->Pt()<fSinglePtCut){
545            return kFALSE; // no reconstruction below the Pt cut
546         }
547         
548         if( TMath::Abs(ePos->Eta())> fEtaCut || TMath::Abs(ePos->Eta())< fEtaCutMin || 
549             TMath::Abs(eNeg->Eta())> fEtaCut || TMath::Abs(eNeg->Eta())< fEtaCutMin ) {
550            return kFALSE;
551         }
552         
553         if(ePos->R()>fMaxR){
554            return kFALSE; // cuts on distance from collision point
555         }
556
557         if(TMath::Abs(ePos->Vz()) > fMaxZ){
558            return kFALSE;        // outside material
559         }
560         if(TMath::Abs(eNeg->Vz()) > fMaxZ){
561            return kFALSE;        // outside material
562         }
563
564         if( ePos->R() <= ((TMath::Abs(ePos->Vz()) * fLineCutZRSlope) - fLineCutZValue)){
565            return kFALSE;  // line cut to exclude regions where we do not reconstruct
566         } else if ( fEtaCutMin != -0.1 &&   ePos->R() >= ((TMath::Abs(ePos->Vz()) * fLineCutZRSlopeMin) - fLineCutZValueMin)){
567            return kFALSE;
568         }
569         
570         if( eNeg->R() <= ((TMath::Abs(eNeg->Vz()) * fLineCutZRSlope) - fLineCutZValue)){
571            return kFALSE; // line cut to exclude regions where we do not reconstruct
572         } else if ( fEtaCutMin != -0.1 &&   eNeg->R() >= ((TMath::Abs(eNeg->Vz()) * fLineCutZRSlopeMin) - fLineCutZValueMin)){
573            return kFALSE;
574         }
575         
576         return kTRUE;
577         //if(AcceptanceCut(particle,ePos,eNeg))return kTRUE;
578     }
579     return kFALSE;
580 }
581
582 //________________________________________________________________________
583 Bool_t AliConversionCuts::MesonIsSelectedMC(TParticle *fMCMother,AliStack *fMCStack,Bool_t bMCDaughtersInAcceptance){
584     // Returns true for all pions within acceptance cuts for decay into 2 photons
585     // If bMCDaughtersInAcceptance is selected, it requires in addition that both daughter photons are within acceptance cuts
586
587     if(!fMCStack)return kFALSE;
588     
589     if(fMCMother->GetPdgCode()==111 || fMCMother->GetPdgCode()==221){
590        
591        if(fMCMother->R()>fMaxR) return kFALSE; // cuts on distance from collision point
592
593        Double_t rapidity = 10.;
594        if(fMCMother->Energy() - fMCMother->Pz() == 0 || fMCMother->Energy() + fMCMother->Pz() == 0){
595           rapidity=8.;
596        }
597        else{
598           rapidity = 0.5*(TMath::Log((fMCMother->Energy()+fMCMother->Pz()) / (fMCMother->Energy()-fMCMother->Pz())));
599        }        
600        
601        // Rapidity Cut
602        if(TMath::Abs(rapidity)>fRapidityCutMeson)return kFALSE;
603
604         // Select only -> 2y decay channel
605         if(fMCMother->GetNDaughters()!=2)return kFALSE;
606
607         for(Int_t i=0;i<2;i++){
608             TParticle *MDaughter=fMCStack->Particle(fMCMother->GetDaughter(i));
609
610             // Is Daughter a Photon?
611             if(MDaughter->GetPdgCode()!=22)return kFALSE;
612             // Is Photon in Acceptance?
613             if(bMCDaughtersInAcceptance){
614                 if(!PhotonIsSelectedMC(MDaughter,fMCStack)){return kFALSE;}
615             }
616         }
617         return kTRUE;
618     }
619     return kFALSE;
620 }
621
622 ///________________________________________________________________________
623 Bool_t AliConversionCuts::PhotonCuts(AliConversionPhotonBase *photon,AliVEvent *event)
624 {   // Specific Photon Cuts
625
626     Int_t cutIndex = 0;
627     if(hPhotonCuts)hPhotonCuts->Fill(cutIndex);
628     cutIndex++;
629
630     // Fill Histos before Cuts
631     if(hInvMassbefore)hInvMassbefore->Fill(photon->GetMass());
632     if(hArmenterosbefore)hArmenterosbefore->Fill(photon->GetArmenterosAlpha(),photon->GetArmenterosQt());
633
634     // Gamma selection based on QT from Armenteros
635     if(fDoQtGammaSelection == kTRUE){
636         if(!ArmenterosQtCut(photon)){
637             if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //1
638             return kFALSE;
639         }
640     }
641     cutIndex++; //2
642
643     // Chi Cut
644     if(photon->GetChi2perNDF() > fChi2CutConversion || photon->GetChi2perNDF() <=0){
645         {
646             if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //2
647             return kFALSE;
648         }
649     }
650     cutIndex++;//3
651
652     // Reconstruction Acceptance Cuts
653     if(!AcceptanceCuts(photon)){
654         if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //3
655         return kFALSE;
656     }
657
658     cutIndex++; //4
659     // Asymmetry Cut
660     if(fDoPhotonAsymmetryCut == kTRUE){
661         if(!AsymmetryCut(photon,event)){
662             if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //4
663             return kFALSE;
664         }
665     }
666
667     //Check the pid probability
668     cutIndex++; //5
669     if(!PIDProbabilityCut(photon, event)) {
670         if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //5
671         return kFALSE;
672     }
673
674     cutIndex++; //6
675     if(!CorrectedTPCClusterCut(photon, event)) {
676         if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //6
677         return kFALSE;
678     }
679
680
681     cutIndex++; //7
682     if(!PsiPairCut(photon)) {
683           if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //7
684           return kFALSE;
685     }
686
687     cutIndex++; //8
688     if(!CosinePAngleCut(photon, event)) {
689           if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //8
690           return kFALSE;
691     }
692
693     cutIndex++; //9
694     if(hPhotonCuts)hPhotonCuts->Fill(cutIndex); //9
695
696     // Histos after Cuts
697     if(hInvMassafter)hInvMassafter->Fill(photon->GetMass());
698     if(hArmenterosafter)hArmenterosafter->Fill(photon->GetArmenterosAlpha(),photon->GetArmenterosQt());
699
700
701     return kTRUE;
702
703 }
704
705 ///________________________________________________________________________
706 Bool_t AliConversionCuts::CorrectedTPCClusterCut(AliConversionPhotonBase *photon, AliVEvent * event)
707 {   //Cut on corrected TPC Cluster Info
708
709     AliVTrack * negTrack = GetTrack(event, photon->GetTrackLabelNegative());
710     AliVTrack * posTrack = GetTrack(event, photon->GetTrackLabelPositive());
711
712     if(!negTrack||!posTrack)return kFALSE;
713
714     Double_t negclsToF=0;
715
716     if (!fUseCorrectedTPCClsInfo ){
717         if(negTrack->GetTPCNclsF()!=0){
718             negclsToF = (Double_t)negTrack->GetNcls(1)/(Double_t)negTrack->GetTPCNclsF();}// Ncluster/Nfindablecluster
719     }
720     else {
721         negclsToF = negTrack->GetTPCClusterInfo(2,0,GetFirstTPCRow(photon->GetConversionRadius()));
722     }
723
724     Double_t posclsToF = 0.;
725     if (!fUseCorrectedTPCClsInfo ){
726         if(posTrack->GetTPCNclsF()!=0   ){
727             posclsToF = (Double_t)posTrack->GetNcls(1)/(Double_t)posTrack->GetTPCNclsF();
728         }
729     }else{
730         posclsToF = posTrack->GetTPCClusterInfo(2,0,GetFirstTPCRow(photon->GetConversionRadius()));
731     }
732
733     if( negclsToF < fMinClsTPCToF ||    posclsToF < fMinClsTPCToF ){
734     return kFALSE;
735     }
736
737     return kTRUE;
738 }
739
740 ///________________________________________________________________________
741 Bool_t AliConversionCuts::PhotonIsSelected(AliConversionPhotonBase *photon, AliVEvent * event)
742 {
743     //Selection of Reconstructed Photons
744
745     FillPhotonCutIndex(kPhotonIn);
746
747     // Get Tracks
748     AliVTrack * negTrack = GetTrack(event, photon->GetTrackLabelNegative());
749     AliVTrack * posTrack = GetTrack(event, photon->GetTrackLabelPositive());
750
751     if(!negTrack || !posTrack) {
752         FillPhotonCutIndex(kNoTracks);
753         return kFALSE;
754     }
755
756     // dEdx Cuts
757     if(!dEdxCuts(negTrack) || !dEdxCuts(posTrack)) {
758         FillPhotonCutIndex(kdEdxCuts);
759         return kFALSE;
760     }
761
762     // Track Cuts
763     if(!TracksAreSelected(negTrack, posTrack)){
764         FillPhotonCutIndex(kTrackCuts);
765         return kFALSE;
766     }
767
768     // Photon Cuts
769     if(!PhotonCuts(photon,event)){
770         FillPhotonCutIndex(kPhotonCuts);
771         return kFALSE;
772     }
773
774     // Photon passed cuts
775     FillPhotonCutIndex(kPhotonOut);
776     return kTRUE;
777 }
778
779 ///________________________________________________________________________
780 Bool_t AliConversionCuts::MesonIsSelected(AliAODConversionMother *pi0,Bool_t IsSignal)
781 {
782     // Selection of reconstructed Meson candidates
783     // Use flag IsSignal in order to fill Fill different
784     // histograms for Signal and Background
785     TH1 *hist=0x0;
786
787     if(IsSignal){hist=hMesonCuts;
788     }
789     else{hist=hMesonBGCuts;}
790
791     Int_t cutIndex=0;
792     if(hist)hist->Fill(cutIndex);
793     cutIndex++;
794
795     // Undefined Rapidity -> Floating Point exception
796     if((pi0->E()+pi0->Pz())/(pi0->E()-pi0->Pz())<=0){
797         if(hist)hist->Fill(cutIndex);
798         cutIndex++;
799         return kFALSE;
800     }
801     else{
802         // PseudoRapidity Cut --> But we cut on Rapidity !!!
803         cutIndex++;
804         if(TMath::Abs(pi0->Rapidity())>fRapidityCutMeson){
805         //if(TMath::Abs(pi0->PseudoRapidity())>fRapidityCutMeson){
806             if(hist)hist->Fill(cutIndex);
807             return kFALSE;
808         }
809     }
810     cutIndex++;
811
812     // Opening Angle Cut
813     //fOpeningAngle=2*TMath::ATan(0.134/pi0->P());// physical minimum opening angle
814     if(pi0->GetOpeningAngle()<fOpeningAngle){
815         if(hist)hist->Fill(cutIndex);
816         return kFALSE;
817     }
818     cutIndex++;
819
820     // Alpha Max Cut
821     if(pi0->GetAlpha()>fAlphaCutMeson){
822         if(hist)hist->Fill(cutIndex);
823         return kFALSE;
824     }
825     cutIndex++;
826
827     // Alpha Min Cut
828     if(pi0->GetAlpha()<fAlphaMinCutMeson){
829         if(hist)hist->Fill(cutIndex);
830         return kFALSE;
831     }
832     cutIndex++;
833
834     if(hist)hist->Fill(cutIndex);
835     return kTRUE;
836 }
837
838
839 ///________________________________________________________________________
840 Bool_t AliConversionCuts::ArmenterosQtCut(AliConversionPhotonBase *photon)
841 {   // Armenteros Qt Cut
842
843   if(fDoHighPtQtGammaSelection){
844         if(photon->GetPhotonPt() < fPtBorderForQt){
845           if(photon->GetArmenterosQt()>fQtMax){
846                 return kFALSE;
847           }
848         } else {
849           if(photon->GetArmenterosQt()>fHighPtQtMax){
850                 return kFALSE;
851           }
852         }
853   } else {
854
855         if(photon->GetArmenterosQt()>fQtMax){
856           return kFALSE;
857         }
858   }
859   return kTRUE;
860 }
861
862
863 ///________________________________________________________________________
864 Bool_t AliConversionCuts::AcceptanceCuts(AliConversionPhotonBase *photon) {
865     // Exclude certain areas for photon reconstruction
866
867     Int_t cutIndex=0;
868     if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
869     cutIndex++;
870
871     if(photon->GetConversionRadius()>fMaxR){ // cuts on distance from collision point
872         if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
873         return kFALSE;
874     }
875     cutIndex++;
876
877     if(photon->GetConversionRadius()<fMinR){ // cuts on distance from collision point
878         if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
879         return kFALSE;
880     }
881     cutIndex++;
882
883     if(photon->GetConversionRadius() <= ((TMath::Abs(photon->GetConversionZ())*fLineCutZRSlope)-fLineCutZValue)){
884       if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
885       return kFALSE;
886     }
887     else if (fUseEtaMinCut &&  photon->GetConversionRadius() >= ((TMath::Abs(photon->GetConversionZ())*fLineCutZRSlopeMin)-fLineCutZValueMin )){
888         if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
889         return kFALSE;
890     }
891     cutIndex++;
892
893   if(TMath::Abs(photon->GetConversionZ()) > fMaxZ ){ // cuts out regions where we do not reconstruct
894       if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
895       return kFALSE;
896   }
897     cutIndex++;
898
899
900   if(TMath::Abs(photon->GetPhotonEta())> fEtaCut || TMath::Abs(photon->GetPhotonEta())< fEtaCutMin){
901       if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
902       return kFALSE;
903   }
904     cutIndex++;
905
906
907   if(photon->GetPhotonPt()<fPtCut){
908       if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
909       return kFALSE;
910   }
911     cutIndex++;
912
913   if(hAcceptanceCuts)hAcceptanceCuts->Fill(cutIndex);
914  
915   return kTRUE;
916 }
917
918
919 ///________________________________________________________________________
920 Bool_t AliConversionCuts::SpecificTrackCuts(AliAODTrack * negTrack, AliAODTrack * posTrack,Int_t &cutIndex) {
921     // Track Cuts which require AOD/ESD specific implementation
922
923   if( !negTrack->IsOn(AliESDtrack::kTPCrefit)  || !posTrack->IsOn(AliESDtrack::kTPCrefit)   )  {
924       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
925       return kFALSE;
926   }
927   cutIndex++;
928
929   AliAODVertex * NegVtxType=negTrack->GetProdVertex();
930   AliAODVertex * PosVtxType=posTrack->GetProdVertex();
931   if((NegVtxType->GetType())==AliAODVertex::kKink  || (PosVtxType->GetType())==AliAODVertex::kKink) {
932       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
933       return kFALSE;
934   }
935   return kTRUE;
936
937 }
938
939
940 ///________________________________________________________________________
941 Bool_t AliConversionCuts::SpecificTrackCuts(AliESDtrack * negTrack, AliESDtrack * posTrack,Int_t &cutIndex) {
942     // Track Cuts which require AOD/ESD specific implementation
943
944    if( !negTrack->IsOn(AliESDtrack::kTPCrefit)  || !posTrack->IsOn(AliESDtrack::kTPCrefit)   )  {
945       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
946       return kFALSE;
947   }
948   cutIndex++;
949
950   if(negTrack->GetKinkIndex(0) > 0  || posTrack->GetKinkIndex(0) > 0 ) {
951       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
952       return kFALSE;
953   }
954   return kTRUE;
955 }
956
957
958
959 ///________________________________________________________________________
960 Bool_t AliConversionCuts::TracksAreSelected(AliVTrack * negTrack, AliVTrack * posTrack) {
961     // Track Selection for Photon Reconstruction
962
963     Int_t cutIndex=0;
964     if(hTrackCuts)hTrackCuts->Fill(cutIndex);
965     cutIndex++;
966
967   // avoid like sign
968   if(negTrack->Charge() == posTrack->Charge()) {
969        if(hTrackCuts)hTrackCuts->Fill(cutIndex);
970         return kFALSE;
971   }
972   cutIndex++;
973
974   // Number of TPC Clusters
975   if( negTrack->GetNcls(1) < fMinClsTPC || posTrack->GetNcls(1) < fMinClsTPC ) {
976       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
977         return kFALSE;
978   }
979   cutIndex++;
980
981   // Acceptance
982
983   if(TMath::Abs(negTrack->Eta()) > fEtaCut || TMath::Abs(negTrack->Eta()) < fEtaCutMin ||
984      TMath::Abs(posTrack->Eta())> fEtaCut || TMath::Abs(posTrack->Eta())< fEtaCutMin) {
985        if(hTrackCuts)hTrackCuts->Fill(cutIndex);
986       return kFALSE;
987   }
988   cutIndex++;
989
990   // Single Pt Cut
991   if( negTrack->Pt()< fSinglePtCut ||   posTrack->Pt()< fSinglePtCut){
992       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
993       return kFALSE;
994   }
995   cutIndex++;
996
997   // AOD ESD specific cuts
998   Bool_t passCuts = kTRUE;
999
1000   if(negTrack->IsA()==AliAODTrack::Class()) {
1001         passCuts = passCuts * SpecificTrackCuts(static_cast<AliAODTrack*>(negTrack), static_cast<AliAODTrack*>(posTrack),cutIndex);
1002   } else { 
1003         passCuts = passCuts * SpecificTrackCuts(static_cast<AliESDtrack*>(negTrack), static_cast<AliESDtrack*>(posTrack),cutIndex);
1004   }     
1005
1006   if(!passCuts){
1007       if(hTrackCuts)hTrackCuts->Fill(cutIndex);
1008       return kFALSE;
1009   }
1010   cutIndex++;
1011
1012   if(hTrackCuts)hTrackCuts->Fill(cutIndex);
1013
1014   return kTRUE;
1015                     
1016 }
1017
1018 ///________________________________________________________________________
1019 Bool_t AliConversionCuts::dEdxCuts(AliVTrack *fCurrentTrack){
1020     // Electron Identification Cuts for Photon reconstruction
1021
1022     if(!fPIDResponse){InitPIDResponse();}// Try to reinitialize PID Response
1023     if(!fPIDResponse){AliError("No PID Response"); return kTRUE;}// if still missing fatal error
1024
1025     Int_t cutIndex=0;
1026     if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1027     if(hTPCdEdxbefore)hTPCdEdxbefore->Fill(fCurrentTrack->P(),fPIDResponse->NumberOfSigmasTPC(fCurrentTrack, AliPID::kElectron));
1028     cutIndex++;
1029     
1030
1031   if(fDodEdxSigmaCut == kTRUE){
1032       // TPC Electron Line
1033       if( fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kElectron)<fPIDnSigmaBelowElectronLine ||
1034                 fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kElectron)>fPIDnSigmaAboveElectronLine){
1035
1036           if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1037           return kFALSE;
1038       }
1039       cutIndex++;
1040
1041       // TPC Pion Line
1042         if( fCurrentTrack->P()>fPIDMinPnSigmaAbovePionLine && fCurrentTrack->P()<fPIDMaxPnSigmaAbovePionLine ){
1043           if(fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kElectron)>fPIDnSigmaBelowElectronLine &&
1044                  fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kElectron)<fPIDnSigmaAboveElectronLine&&
1045                  fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kPion)<fPIDnSigmaAbovePionLine){
1046
1047               if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1048               return kFALSE;
1049           }
1050         }
1051         cutIndex++;
1052    
1053         // High Pt Pion rej
1054         if( fCurrentTrack->P()>fPIDMaxPnSigmaAbovePionLine ){
1055           if(fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kElectron)>fPIDnSigmaBelowElectronLine &&
1056                  fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kElectron)<fPIDnSigmaAboveElectronLine&&
1057                  fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kPion)<fPIDnSigmaAbovePionLineHighPt){
1058
1059                 if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1060                 return kFALSE;
1061           }
1062         }
1063         cutIndex++;
1064   }
1065   else{cutIndex+=3;}
1066
1067   if(fDoKaonRejectionLowP == kTRUE){
1068         if(fCurrentTrack->P()<fPIDMinPKaonRejectionLowP ){
1069           if( TMath::Abs(fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kKaon))<fPIDnSigmaAtLowPAroundKaonLine){
1070
1071               if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1072               return kFALSE;
1073           }
1074         }
1075   }
1076   cutIndex++;
1077    
1078   if(fDoProtonRejectionLowP == kTRUE){
1079         if( fCurrentTrack->P()<fPIDMinPProtonRejectionLowP ){
1080           if( TMath::Abs(fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kProton))<fPIDnSigmaAtLowPAroundProtonLine){
1081
1082               if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1083               return kFALSE;
1084           }
1085         }
1086   }
1087    cutIndex++;
1088    
1089   if(fDoPionRejectionLowP == kTRUE){
1090         if( fCurrentTrack->P()<fPIDMinPPionRejectionLowP ){
1091           if( TMath::Abs(fPIDResponse->NumberOfSigmasTPC(fCurrentTrack,AliPID::kPion))<fPIDnSigmaAtLowPAroundPionLine){
1092
1093               if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1094                 return kFALSE;
1095           }
1096         }
1097   }
1098   cutIndex++;
1099    
1100   if((fCurrentTrack->GetStatus() & AliESDtrack::kTOFpid) && !(fCurrentTrack->GetStatus() & AliESDtrack::kTOFmismatch)){
1101      if(hTOFbefore) hTOFbefore->Fill(fCurrentTrack->P(),fPIDResponse->NumberOfSigmasTOF(fCurrentTrack, AliPID::kElectron));
1102      if(fUseTOFpid){
1103         if(fPIDResponse->NumberOfSigmasTOF(fCurrentTrack, AliPID::kElectron)>fTofPIDnSigmaAboveElectronLine ||
1104            fPIDResponse->NumberOfSigmasTOF(fCurrentTrack, AliPID::kElectron)<fTofPIDnSigmaBelowElectronLine ){
1105            if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1106            return kFALSE;
1107         }
1108      }
1109      if(hTOFafter)hTOFafter->Fill(fCurrentTrack->P(),fPIDResponse->NumberOfSigmasTOF(fCurrentTrack, AliPID::kElectron));
1110   }
1111      cutIndex++;
1112    
1113     // Apply TRD PID
1114   if(fDoTRDPID){
1115         if(!fPIDResponse->IdentifiedAsElectronTRD(fCurrentTrack,fPIDTRDEfficiency)){
1116
1117             if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1118             return kFALSE;
1119         }
1120   }
1121   cutIndex++;
1122
1123   if(hdEdxCuts)hdEdxCuts->Fill(cutIndex);
1124   if(hTPCdEdxafter)hTPCdEdxafter->Fill(fCurrentTrack->P(),fPIDResponse->NumberOfSigmasTPC(fCurrentTrack, AliPID::kElectron));
1125
1126   return kTRUE;
1127 }
1128
1129 ///________________________________________________________________________
1130 Bool_t AliConversionCuts::AsymmetryCut(AliConversionPhotonBase * photon,AliVEvent *event) {
1131     // Cut on Energy Assymetry
1132
1133     for(Int_t ii=0;ii<2;ii++){
1134
1135         AliVTrack *track=GetTrack(event,photon->GetTrackLabel(ii));
1136
1137         if( track->P() > fMinPPhotonAsymmetryCut ){
1138             Double_t trackNegAsy=0;
1139             if (photon->GetPhotonP()!=0.){
1140                 trackNegAsy= track->P()/photon->GetPhotonP();
1141             }
1142
1143             if( trackNegAsy<fMinPhotonAsymmetry ||trackNegAsy>(1.- fMinPhotonAsymmetry)){
1144                 return kFALSE;
1145             }
1146         }
1147
1148     }
1149     return kTRUE;
1150 }
1151
1152 ///________________________________________________________________________
1153 AliVTrack *AliConversionCuts::GetTrack(AliVEvent * event, Int_t label){
1154     //Returns pointer to the track with given ESD label
1155     //(Important for AOD implementation, since Track array in AOD data is different
1156     //from ESD array, but ESD tracklabels are stored in AOD Tracks)
1157
1158   AliESDEvent * esdEvent = dynamic_cast<AliESDEvent*>(event);
1159   if(esdEvent) {
1160         if(label > event->GetNumberOfTracks() ) return NULL;
1161         AliESDtrack * track = esdEvent->GetTrack(label);
1162         return track;
1163         
1164   } else { 
1165         for(Int_t ii=0; ii<event->GetNumberOfTracks(); ii++) {
1166           AliVTrack * track = dynamic_cast<AliVTrack*>(event->GetTrack(ii));
1167           
1168           if(track) { 
1169                 if(track->GetID() == label) {
1170                   return track;
1171                 }
1172           }
1173         }
1174   }
1175   
1176   //AliDebug(5,(Form("track not found %d %d",label,event->GetNumberOfTracks()));
1177   return NULL;
1178 }
1179
1180
1181 ///________________________________________________________________________
1182 Bool_t AliConversionCuts::PIDProbabilityCut(AliConversionPhotonBase *photon, AliVEvent * event){
1183     // Cut on Electron Probability for Photon Reconstruction
1184
1185   AliESDEvent * esdEvent = dynamic_cast<AliESDEvent*>(event);
1186
1187   if(esdEvent){
1188         
1189         Bool_t iResult=kFALSE;
1190         
1191         Double_t *posProbArray = new Double_t[AliPID::kSPECIES];
1192         Double_t *negProbArray = new Double_t[AliPID::kSPECIES];
1193         
1194         AliESDtrack* negTrack   = esdEvent->GetTrack(photon->GetTrackLabelNegative());
1195         AliESDtrack* posTrack   = esdEvent->GetTrack(photon->GetTrackLabelPositive());
1196         
1197         if(negProbArray && posProbArray){
1198           
1199           negTrack->GetTPCpid(negProbArray);
1200           posTrack->GetTPCpid(posProbArray);
1201
1202           if(negProbArray[AliPID::kElectron]>=fPIDProbabilityCutNegativeParticle && posProbArray[AliPID::kElectron]>=fPIDProbabilityCutPositiveParticle){
1203                 iResult=kTRUE;
1204           }
1205         }
1206         
1207         delete [] posProbArray;
1208         delete [] negProbArray;
1209         return iResult;
1210
1211   } else {
1212       ///Not possible for AODs
1213       return kTRUE;
1214   }
1215
1216
1217
1218   
1219 }
1220
1221
1222 ///________________________________________________________________________
1223 Bool_t AliConversionCuts::AcceptanceCut(TParticle *particle, TParticle * ePos,TParticle* eNeg){
1224     // MC Acceptance Cuts
1225     //(Certain areas were excluded for photon reconstruction)
1226
1227   if(particle->R()>fMaxR){
1228         return kFALSE;}
1229
1230   if(ePos->R()>fMaxR){
1231         return kFALSE;
1232   }
1233
1234   if(ePos->R()<fMinR){
1235         return kFALSE;
1236   }
1237
1238   if( ePos->R() <= ((TMath::Abs(ePos->Vz())*fLineCutZRSlope)-fLineCutZValue)){
1239         return kFALSE;
1240   }
1241   else if (fUseEtaMinCut &&  ePos->R() >= ((TMath::Abs(ePos->Vz())*fLineCutZRSlopeMin)-fLineCutZValueMin )){
1242       return kFALSE;
1243   }
1244
1245   if(TMath::Abs(eNeg->Vz()) > fMaxZ ){ // cuts out regions where we do not reconstruct
1246         return kFALSE;
1247   }
1248
1249   if(eNeg->Vz()!=ePos->Vz()||eNeg->R()!=ePos->R()){
1250         return kFALSE;
1251   }
1252
1253   if(TMath::Abs(ePos->Vz()) > fMaxZ ){ // cuts out regions where we do not reconstruct
1254         return kFALSE;
1255   }
1256
1257   if(TMath::Abs(particle->Eta())> fEtaCut || TMath::Abs(particle->Eta())< fEtaCutMin){
1258         return kFALSE;
1259   }
1260
1261   if(TMath::Abs(ePos->Eta())> fEtaCut || TMath::Abs(ePos->Eta())< fEtaCutMin){
1262         return kFALSE;
1263   }
1264
1265   if(TMath::Abs(eNeg->Eta())> fEtaCut || TMath::Abs(eNeg->Eta())< fEtaCutMin){
1266         return kFALSE;
1267   }
1268
1269   if( ePos->Pt()< fSinglePtCut ||  eNeg->Pt()< fSinglePtCut){
1270         return kFALSE;
1271   }
1272
1273   if(particle->Pt()<fPtCut){
1274         return kFALSE;
1275   }
1276
1277   return kTRUE;
1278 }
1279 ///________________________________________________________________________
1280 Bool_t AliConversionCuts::UpdateCutString(cutIds cutID, Int_t value) {
1281 ///Update the cut string (if it has been created yet)
1282
1283   if(fCutString && fCutString->GetString().Length() == kNCuts) {
1284         fCutString->SetString(GetCutNumber());
1285   } else {
1286         return kFALSE;
1287   }
1288   return kTRUE;
1289 }
1290
1291 ///________________________________________________________________________
1292 Bool_t AliConversionCuts::InitializeCutsFromCutString(const TString analysisCutSelection ) {
1293    // Initialize Cuts from a given Cut string
1294
1295     AliInfo(Form("Set Cut Number: %s",analysisCutSelection.Data()));
1296   if(analysisCutSelection.Length()!=kNCuts) {
1297         AliError(Form("Cut selection has the wrong length! size is %d, number of cuts is %d", analysisCutSelection.Length(), kNCuts));
1298         return kFALSE;
1299   }
1300   if(!analysisCutSelection.IsDigit()){
1301         AliError("Cut selection contains characters");
1302         return kFALSE;
1303   }
1304   
1305   const char *cutSelection = analysisCutSelection.Data();
1306   #define ASSIGNARRAY(i)        fCuts[i] = cutSelection[i] - '0'
1307   for(Int_t ii=0;ii<kNCuts;ii++){
1308       ASSIGNARRAY(ii);
1309   }
1310
1311   // TestFlag
1312   if(fCuts[0] !=9){
1313     AliError("Analysis Cut Selection does not start with 9");
1314         PrintCuts();
1315     return kFALSE;
1316   }
1317
1318   // Set Individual Cuts
1319   for(Int_t ii=0;ii<kNCuts;ii++){
1320       if(!SetCut(cutIds(ii),fCuts[ii]))return kFALSE;
1321   }
1322
1323   //PrintCuts();
1324
1325   // Set StandardTriggers
1326
1327   if(fIsHeavyIon)SelectCollisionCandidates(AliVEvent::kMB | AliVEvent::kCentral | AliVEvent::kSemiCentral);
1328   else SelectCollisionCandidates(AliVEvent::kMB);
1329
1330   return kTRUE;
1331 }
1332 ///________________________________________________________________________
1333 Bool_t AliConversionCuts::SetCut(cutIds cutID, const Int_t value) {
1334   ///Set individual cut ID
1335
1336
1337   switch (cutID) {
1338   case kgoodId:
1339         fCuts[kgoodId] = value;
1340         if(value != 9) {
1341           AliError("First value of cut string is wrong, aborting!!");
1342           return kFALSE;
1343         } else {
1344           return kTRUE;
1345         }
1346
1347   case kv0FinderType:
1348         if( SetV0Finder(value)) {
1349           fCuts[kv0FinderType] = value;
1350           UpdateCutString(cutID, value);
1351           return kTRUE;
1352         } else return kFALSE;
1353
1354   case keProbCut:
1355         if( SetElectronProbCut(value)) {
1356           fCuts[keProbCut] = value;
1357           UpdateCutString(cutID, value);
1358           return kTRUE;
1359         } else return kFALSE;
1360
1361   case kededxSigmaCut:
1362         if( SetTPCdEdxCutElectronLine(value)) {
1363           fCuts[kededxSigmaCut] = value;
1364           UpdateCutString(cutID, value);
1365           return kTRUE;
1366         } else return kFALSE;
1367
1368   case kpidedxSigmaCut:
1369         if( SetTPCdEdxCutPionLine(value)) {
1370           fCuts[kpidedxSigmaCut] = value;
1371           UpdateCutString(cutID, value);
1372           return kTRUE;
1373         } else return kFALSE;
1374
1375   case kpiMomdedxSigmaCut:
1376         if( SetMinMomPiondEdxCut(value)) {
1377           fCuts[kpiMomdedxSigmaCut] = value;
1378           UpdateCutString(cutID, value);
1379           return kTRUE;
1380         } else return kFALSE;
1381
1382   case kchi2GammaCut:
1383         if( SetChi2GammaCut(value)) {
1384           fCuts[kchi2GammaCut] = value;
1385           UpdateCutString(cutID, value);
1386           return kTRUE;
1387         } else return kFALSE;
1388
1389   case ksinglePtCut:
1390         if( SetSinglePtCut(value)) {
1391           fCuts[ksinglePtCut] = value;
1392           UpdateCutString(cutID, value);
1393           return kTRUE;
1394         } else return kFALSE;
1395
1396   case kclsTPCCut:
1397         if( SetTPCClusterCut(value)) {
1398           fCuts[kclsTPCCut] = value;
1399           UpdateCutString(cutID, value);
1400           return kTRUE;
1401         } else return kFALSE;
1402
1403   case ketaCut:
1404         if( SetEtaCut(value)) {
1405           fCuts[ketaCut] = value;
1406           UpdateCutString(cutID, value);
1407           return kTRUE;
1408         } else return kFALSE;
1409
1410   case kchi2MesonCut:
1411         if( SetChi2MesonCut(value)) {
1412           fCuts[kchi2MesonCut] = value;
1413           UpdateCutString(cutID, value);
1414           return kTRUE;
1415         } else return kFALSE;
1416
1417   case kLowPRejectionSigmaCut:
1418         if( SetLowPRejectionCuts(value)) {
1419           fCuts[kLowPRejectionSigmaCut] = value;
1420           UpdateCutString(cutID, value);
1421           return kTRUE;
1422         } else return kFALSE;
1423
1424   case kQtMaxCut:
1425         if( SetQtMaxCut(value)) {
1426           fCuts[kQtMaxCut] = value;
1427           UpdateCutString(cutID, value);
1428           return kTRUE;
1429         } else return kFALSE;
1430
1431   case kpiMaxMomdedxSigmaCut:
1432         if( SetMaxMomPiondEdxCut(value)) {
1433           fCuts[kpiMaxMomdedxSigmaCut] = value;
1434           UpdateCutString(cutID, value);
1435           return kTRUE;
1436         } else return kFALSE;
1437
1438   case kalphaMesonCut:
1439         if( SetAlphaMesonCut(value)) {
1440           fCuts[kalphaMesonCut] = value;
1441           UpdateCutString(cutID, value);
1442           return kTRUE;
1443         } else return kFALSE;
1444
1445   case kminRCut:
1446         if( SetRCut(value)) {
1447           fCuts[kminRCut] = value;
1448           UpdateCutString(cutID, value);
1449           return kTRUE;
1450         } else return kFALSE;
1451
1452   case kRapidityMesonCut:
1453         if( SetRapidityMesonCut(value)) {
1454           fCuts[kRapidityMesonCut] = value;
1455           UpdateCutString(cutID, value);
1456           return kTRUE;
1457         } else return kFALSE;
1458
1459   case kBackgroundScheme:
1460         if( SetBackgroundScheme(value)) {
1461           fCuts[kBackgroundScheme] = value;
1462           UpdateCutString(cutID, value);
1463           return kTRUE;
1464         } else return kFALSE;
1465
1466   case kDegreesForRotationMethod:
1467         if( SetNDegreesForRotationMethod(value)) {
1468           fCuts[kDegreesForRotationMethod] = value;
1469           UpdateCutString(cutID, value);
1470           return kTRUE;
1471         } else return kFALSE;
1472
1473   case kNumberOfRotations:
1474         if( SetNumberOfRotations(value)) {
1475           fCuts[kNumberOfRotations] = value;
1476           UpdateCutString(cutID, value);
1477           return kTRUE;
1478         } else return kFALSE;
1479
1480   case kremovePileUp:
1481         if( SetRemovePileUp(value)) {
1482           fCuts[kremovePileUp] = value;
1483           UpdateCutString(cutID, value);
1484           return kTRUE;
1485         } else return kFALSE;
1486
1487   case kselectV0AND:
1488         if( SetSelectV0AND(value)) {
1489           fCuts[kselectV0AND] = value;
1490           UpdateCutString(cutID, value);
1491           return kTRUE;
1492         } else return kFALSE;
1493
1494   case kmultiplicityMethod:
1495         if( SetMultiplicityMethod(value)) {
1496           fCuts[kmultiplicityMethod] = value;
1497           UpdateCutString(cutID, value);
1498           return kTRUE;
1499         } else return kFALSE;
1500
1501   case kisHeavyIon:
1502         if( SetIsHeavyIon(value)) {
1503           fCuts[kisHeavyIon] = value;
1504           UpdateCutString(cutID, value);
1505           return kTRUE;
1506         } else return kFALSE;
1507
1508   case kCentralityMin:
1509         if( SetCentralityMin(value)) {
1510           fCuts[kCentralityMin] = value;
1511           UpdateCutString(cutID, value);
1512           return kTRUE;
1513         } else return kFALSE;
1514
1515   case kCentralityMax:
1516         if( SetCentralityMax(value)) {
1517           fCuts[kCentralityMax] = value;
1518           UpdateCutString(cutID, value);
1519           return kTRUE;
1520         } else return kFALSE;
1521
1522   case kTOFelectronPID:
1523         if( SetTOFElectronPIDCut(value)) {
1524           fCuts[kTOFelectronPID] = value;
1525           UpdateCutString(cutID, value);
1526           return kTRUE;
1527         } else return kFALSE;
1528
1529   case kuseMCPSmearing:
1530         if( SetMCPSmearing(value)) {
1531           fCuts[kuseMCPSmearing] = value;
1532           UpdateCutString(cutID, value);
1533           return kTRUE;
1534         } else return kFALSE;
1535
1536   case kdoPhotonAsymmetryCut:
1537         if( SetPhotonAsymmetryCut(value)) {
1538           fCuts[kdoPhotonAsymmetryCut] = value;
1539           UpdateCutString(cutID, value);
1540           return kTRUE;
1541         } else return kFALSE;
1542
1543   case kPsiPair:
1544         if( SetPsiPairCut(value)) {
1545           fCuts[kPsiPair] = value;
1546           UpdateCutString(cutID, value);
1547           return kTRUE;
1548         } else return kFALSE;
1549
1550   case kCosPAngle:
1551         if( SetCosPAngleCut(value)) {
1552           fCuts[kCosPAngle] = value;
1553           UpdateCutString(cutID, value);
1554           return kTRUE;
1555         } else return kFALSE;
1556
1557
1558   case kElecShare:
1559         if( SetSharedElectronCut(value)) {
1560           fCuts[kElecShare] = value;
1561           UpdateCutString(cutID, value);
1562           return kTRUE;
1563         } else return kFALSE;
1564
1565
1566   case kToCloseV0s:
1567         if( SetToCloseV0sCut(value)) {
1568           fCuts[kToCloseV0s] = value;
1569           UpdateCutString(cutID, value);
1570           return kTRUE;
1571         } else return kFALSE;
1572
1573   case kNCuts:
1574       AliError("Cut id out of range");
1575         return kFALSE;
1576   }
1577
1578   AliError("Cut id %d not recognized");
1579   return kFALSE;
1580
1581 }
1582 ///________________________________________________________________________
1583 Bool_t AliConversionCuts::SetRemovePileUp(Int_t removePileUp)
1584 {// Set Cut
1585     switch(removePileUp){
1586     case 0:
1587         fRemovePileUp=kFALSE;
1588         break;
1589     case 1:
1590         fRemovePileUp=kTRUE;
1591         break;
1592     default:
1593         AliError("RemovePileUpCut not defined");
1594         return kFALSE;
1595     }
1596     return kTRUE;
1597 }
1598
1599 ///________________________________________________________________________
1600 Bool_t AliConversionCuts::SetSelectV0AND(Int_t selectV0AND)
1601 {// Set Cut
1602   switch(selectV0AND){
1603   case 0:
1604       fSelectV0AND=kFALSE;
1605       break;
1606   case 1:
1607       fSelectV0AND=kTRUE;
1608       break;
1609   default:
1610       AliError("Warning: V0ANDCut not defined");
1611       return kFALSE;
1612   }
1613   return kTRUE;
1614 }
1615
1616 ///________________________________________________________________________
1617 Bool_t AliConversionCuts::SetMultiplicityMethod(Int_t multiplicityMethod)
1618 {
1619     // Set Cut
1620     fMultiplicityMethod=multiplicityMethod;
1621
1622     // 0 Photon Multiplicity
1623     // 1 TPC Track multiplicity
1624     // 2 V0 Mult
1625     // 3 SPD Mult
1626
1627     return kTRUE;
1628 }
1629
1630 ///________________________________________________________________________
1631 Bool_t AliConversionCuts::SetMCPSmearing(Int_t useMCPSmearing)
1632 {// Set Cut
1633     switch(useMCPSmearing){
1634     case 0:
1635         fUseMCPSmearing=0;
1636         fPBremSmearing=1.;
1637         fPSigSmearing=0.;
1638         fPSigSmearingCte=0.;
1639         break;
1640     case 1:
1641         fUseMCPSmearing=1;
1642         fPBremSmearing=1.0e-14;
1643         fPSigSmearing=0.;
1644         fPSigSmearingCte=0.;
1645         break;
1646     case 2:
1647         fUseMCPSmearing=1;
1648         fPBremSmearing=1.0e-15;
1649         fPSigSmearing=0.0;
1650         fPSigSmearingCte=0.;
1651         break;
1652     case 3:
1653         fUseMCPSmearing=1;
1654         fPBremSmearing=1.;
1655         fPSigSmearing=0.003;
1656         fPSigSmearingCte=0.002;
1657         break;
1658     case 4:
1659         fUseMCPSmearing=1;
1660         fPBremSmearing=1.;
1661         fPSigSmearing=0.003;
1662         fPSigSmearingCte=0.007;
1663         break;
1664     case 5:
1665         fUseMCPSmearing=1;
1666         fPBremSmearing=1.;
1667         fPSigSmearing=0.003;
1668         fPSigSmearingCte=0.016;
1669         break;
1670     case 6:
1671         fUseMCPSmearing=1;
1672         fPBremSmearing=1.;
1673         fPSigSmearing=0.007;
1674         fPSigSmearingCte=0.016;
1675         break;
1676     case 7:
1677         fUseMCPSmearing=1;
1678         fPBremSmearing=1.0e-16;
1679         fPSigSmearing=0.0;
1680         fPSigSmearingCte=0.;
1681         break;
1682     case 8:
1683         fUseMCPSmearing=1;
1684         fPBremSmearing=1.;
1685         fPSigSmearing=0.007;
1686         fPSigSmearingCte=0.014;
1687         break;
1688     case 9:
1689         fUseMCPSmearing=1;
1690         fPBremSmearing=1.;
1691         fPSigSmearing=0.007;
1692         fPSigSmearingCte=0.011;
1693         break;
1694
1695     default:
1696         AliError("Warning: UseMCPSmearing not defined");
1697         return kFALSE;
1698     }
1699     return kTRUE;
1700 }
1701
1702 ///________________________________________________________________________
1703 void AliConversionCuts::PrintCuts() {
1704     // Print out current Cut Selection
1705   for(Int_t ic = 0; ic < kNCuts; ic++) {
1706         printf("%-30s : %d \n", fgkCutNames[ic], fCuts[ic]);
1707   }
1708
1709 }
1710
1711 ///________________________________________________________________________
1712 Bool_t AliConversionCuts::SetRCut(Int_t RCut){
1713     // Set Cut
1714     switch(RCut){
1715     case 0:
1716         fMinR=0;
1717         fMaxR = 180.;
1718         break;
1719     case 1:
1720         fMinR=2.8;
1721         fMaxR = 180.;
1722         break;
1723     case 2:
1724         fMinR=5.;
1725         fMaxR = 180.;
1726         break;
1727     case 3:
1728         fMaxR = 70.;
1729         fMinR = 10.;
1730         break;
1731     case 4:
1732         fMaxR = 70.;
1733         fMinR = 5.;
1734         break;
1735     case 5:
1736         fMaxR = 180.;
1737         fMinR = 10.;
1738         break;
1739     // High purity cuts for PbPb
1740     case 9:
1741         fMaxR = 180.;
1742         fMinR = 60.;
1743         break;
1744
1745     default:
1746         AliError("RCut not defined");
1747         return kFALSE;
1748     }
1749     return kTRUE;
1750 }
1751
1752 ///________________________________________________________________________
1753 Bool_t AliConversionCuts::SetTPCdEdxCutElectronLine(Int_t ededxSigmaCut)
1754 {   // Set Cut
1755     switch(ededxSigmaCut){
1756     case 0: // -10,10
1757         fPIDnSigmaBelowElectronLine=-10;
1758         fPIDnSigmaAboveElectronLine=10;
1759         break;
1760     case 1: // -5,5
1761         fPIDnSigmaBelowElectronLine=-5;
1762         fPIDnSigmaAboveElectronLine=5;
1763         break;
1764     case 2: // -3,5
1765         fPIDnSigmaBelowElectronLine=-3;
1766         fPIDnSigmaAboveElectronLine=5;
1767         break;
1768     case 3: // -4,5
1769         fPIDnSigmaBelowElectronLine=-4;
1770         fPIDnSigmaAboveElectronLine=5;
1771         break;
1772     case 4: // -6,7
1773        fPIDnSigmaBelowElectronLine=-6;
1774        fPIDnSigmaAboveElectronLine=7;
1775        break;
1776     case 5: // -4,4
1777        fPIDnSigmaBelowElectronLine=-4;
1778        fPIDnSigmaAboveElectronLine=4;
1779        break;
1780     case 6: // -2.5,4
1781        fPIDnSigmaBelowElectronLine=-2.5;
1782        fPIDnSigmaAboveElectronLine=4;
1783        break;
1784     case 7: // -2,3.5
1785        fPIDnSigmaBelowElectronLine=-2;
1786        fPIDnSigmaAboveElectronLine=3.5;
1787        break;
1788     default:
1789         AliError("TPCdEdxCutElectronLine not defined");
1790         return kFALSE;
1791         
1792     }
1793     return kTRUE;
1794 }
1795
1796 ///________________________________________________________________________
1797 Bool_t AliConversionCuts::SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut)
1798 {   // Set Cut
1799
1800     switch(pidedxSigmaCut){
1801     case 0:  // -10
1802         fPIDnSigmaAbovePionLine=-10;
1803         fPIDnSigmaAbovePionLineHighPt=-10;
1804         break;
1805     case 1:   // 0
1806         fPIDnSigmaAbovePionLine=0;
1807         fPIDnSigmaAbovePionLineHighPt=-10;
1808         break;
1809     case 2:  // 1
1810         fPIDnSigmaAbovePionLine=1;
1811         fPIDnSigmaAbovePionLineHighPt=-10;
1812         break;
1813     case 3:  // 1
1814         fPIDnSigmaAbovePionLine=-1;
1815         fPIDnSigmaAbovePionLineHighPt=-10;
1816         break;
1817     case 4:  // 1
1818         fPIDnSigmaAbovePionLine=2.5;
1819         fPIDnSigmaAbovePionLineHighPt=-10;
1820         break;
1821     case 5:  // 1
1822         fPIDnSigmaAbovePionLine=2.;
1823         fPIDnSigmaAbovePionLineHighPt=-10;
1824         break;
1825     case 6:  // 1
1826         fPIDnSigmaAbovePionLine=2.;
1827         fPIDnSigmaAbovePionLineHighPt=0.5;
1828         break;
1829     case 7:  // 1
1830         fPIDnSigmaAbovePionLine=3.5;
1831         fPIDnSigmaAbovePionLineHighPt=-10;
1832         break;
1833     case 8:  // 1
1834         fPIDnSigmaAbovePionLine=2.;
1835         fPIDnSigmaAbovePionLineHighPt=1.;
1836         break;
1837     case 9:
1838         fPIDnSigmaAbovePionLine=3.0; // We need a bit less tight cut on dE/dx
1839         fPIDnSigmaAbovePionLineHighPt=-10;
1840         break;
1841     default:
1842         AliError(Form("Warning: pidedxSigmaCut not defined %d",pidedxSigmaCut));
1843         return kFALSE;
1844     }
1845     return kTRUE;
1846 }
1847
1848 ///________________________________________________________________________
1849 Bool_t AliConversionCuts::SetMinMomPiondEdxCut(Int_t piMomdedxSigmaCut)
1850 {   // Set Cut
1851     switch(piMomdedxSigmaCut){
1852     case 0:  // 0.5 GeV
1853         fPIDMinPnSigmaAbovePionLine=0.5;
1854         break;
1855     case 1:  // 1. GeV
1856         fPIDMinPnSigmaAbovePionLine=1.;
1857         break;
1858     case 2:  // 1.5 GeV
1859         fPIDMinPnSigmaAbovePionLine=1.5;
1860         break;
1861     case 3:  // 20.0 GeV
1862         fPIDMinPnSigmaAbovePionLine=20.;
1863         break;
1864     case 4:  // 50.0 GeV
1865         fPIDMinPnSigmaAbovePionLine=50.;
1866         break;
1867     case 5:  // 0.3 GeV
1868         fPIDMinPnSigmaAbovePionLine=0.3;
1869         break;
1870     case 6:  // 0.25 GeV
1871         fPIDMinPnSigmaAbovePionLine=0.25;
1872         break;
1873     case 7:  // 0.4 GeV
1874         fPIDMinPnSigmaAbovePionLine=0.4;
1875         break;
1876     default:
1877         AliError(Form("piMomdedxSigmaCut not defined %d",piMomdedxSigmaCut));
1878         return kFALSE;
1879     }
1880     return kTRUE;
1881 }
1882
1883 ///________________________________________________________________________
1884 Bool_t AliConversionCuts::SetChi2GammaCut(Int_t chi2GammaCut)
1885 {   // Set Cut
1886
1887     switch(chi2GammaCut){
1888     case 0: // 100
1889         fChi2CutConversion = 100.;
1890         break;
1891     case 1:  // 50
1892         fChi2CutConversion = 50.;
1893         break;
1894     case 2:  // 30
1895         fChi2CutConversion = 30.;
1896         break;
1897     case 3:
1898         fChi2CutConversion = 200.;
1899         break;
1900     case 4:
1901         fChi2CutConversion = 500.;
1902         break;
1903     case 5:
1904         fChi2CutConversion = 100000.;
1905         break;
1906     case 6:
1907         fChi2CutConversion = 5.;
1908         break;
1909     case 7:
1910         fChi2CutConversion = 10.;
1911         break;
1912     case 8:
1913         fChi2CutConversion = 20.;
1914         break;
1915     case 9:
1916         fChi2CutConversion = 15.;
1917         break;
1918     default:
1919         AliError(Form("Warning: Chi2GammaCut not defined %d",chi2GammaCut));
1920         return kFALSE;
1921     }
1922     return kTRUE;
1923 }
1924
1925 ///________________________________________________________________________
1926 Bool_t AliConversionCuts::SetV0Finder(Int_t v0FinderType)
1927 {   // Set Cut
1928     switch (v0FinderType){
1929     case 0:  // on fly V0 finder
1930         fUseOnFlyV0Finder=kTRUE;
1931         break;
1932     case 1:  // offline V0 finder
1933         fUseOnFlyV0Finder=kFALSE;
1934         break;
1935     default:
1936         AliError(Form(" v0FinderType not defined %d",v0FinderType));
1937         return kFALSE;
1938     }
1939     return kTRUE;
1940 }
1941
1942 ///________________________________________________________________________
1943 Bool_t AliConversionCuts::SetElectronProbCut(Int_t eProbCut)
1944 {   // Set Cut
1945
1946     switch(eProbCut){
1947     case 0:
1948         fPIDProbabilityCutNegativeParticle=0;
1949         fPIDProbabilityCutPositiveParticle=0;
1950         break;
1951     case 1:
1952         fPIDProbabilityCutNegativeParticle=0.1;
1953         fPIDProbabilityCutPositiveParticle=0.1;
1954         break;
1955     case 2:
1956         fPIDProbabilityCutNegativeParticle=0.5;
1957         fPIDProbabilityCutPositiveParticle=0.5;
1958         break;
1959     case 3:
1960         fPIDProbabilityCutNegativeParticle=0.7;
1961         fPIDProbabilityCutPositiveParticle=0.7;
1962         break;
1963     default:
1964         AliError(Form("Warning: eProbCut not defined %d",eProbCut));
1965         return kFALSE;
1966     }
1967     return kTRUE;
1968 }
1969
1970 ///________________________________________________________________________
1971 Bool_t AliConversionCuts::SetSinglePtCut(Int_t singlePtCut)
1972 {   // Set Cut
1973     switch(singlePtCut){
1974     case 0: // 0.050 GeV
1975         fSinglePtCut = 0.050;
1976         break;
1977     case 1:  // 0.100 GeV
1978         fSinglePtCut = 0.100;
1979         break;
1980     case 2:  // 0.150 GeV
1981         fSinglePtCut = 0.150;
1982         break;
1983     case 3:  // 0.200 GeV
1984         fSinglePtCut = 0.200;
1985         break;
1986     case 4:  // 0.075 GeV
1987         fSinglePtCut = 0.075;
1988         break;
1989     case 5:  // 0.125 GeV
1990         fSinglePtCut = 0.125;
1991         break;
1992     default:
1993         AliError(Form("singlePtCut not defined %d",singlePtCut));
1994         return kFALSE;
1995     }
1996     return kTRUE;
1997 }
1998 ///________________________________________________________________________
1999 Bool_t AliConversionCuts::SetTPCClusterCut(Int_t clsTPCCut)
2000 {   // Set Cut
2001     switch(clsTPCCut){
2002     case 0: // 0
2003         fMinClsTPC= 0.;
2004         break;
2005     case 1:  // 70
2006         fMinClsTPC= 70.;
2007         break;
2008     case 2:  // 80
2009         fMinClsTPC= 80.;
2010         break;
2011     case 3:  // 100
2012         fMinClsTPC= 100.;
2013         break;
2014     case 4:  // 60% of findable clusters
2015         fMinClsTPCToF= 0.6;
2016         fUseCorrectedTPCClsInfo=0;
2017         break;
2018     case 5:  // 0% of findable clusters
2019         fMinClsTPCToF= 0.0;
2020         fUseCorrectedTPCClsInfo=1;
2021         break;
2022     case 6:  // 0% of findable clusters
2023         fMinClsTPCToF= 0.7;
2024         fUseCorrectedTPCClsInfo=1;
2025         break;
2026     case 7:  // 0% of findable clusters
2027         fMinClsTPCToF= 0.35;
2028         fUseCorrectedTPCClsInfo=0;
2029         break;
2030     case 8:
2031         fMinClsTPCToF= 0.35;
2032         fUseCorrectedTPCClsInfo=1;
2033         break;
2034     case 9:
2035         fMinClsTPCToF= 0.6;
2036         fUseCorrectedTPCClsInfo=1;
2037         break;
2038     default:
2039         AliError(Form("Warning: clsTPCCut not defined %d",clsTPCCut));
2040         return kFALSE;
2041     }
2042     return kTRUE;
2043 }
2044
2045 ///________________________________________________________________________
2046 Bool_t AliConversionCuts::SetEtaCut(Int_t etaCut)
2047 {   // Set Cut
2048
2049    //Set Standard LineCutZValues
2050    fLineCutZValueMin = -2;
2051    fLineCutZValue = 7.;
2052    
2053     switch(etaCut){
2054     case 0: // 0.9
2055         fEtaCut         = 0.9;
2056         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2057         fEtaCutMin              = -0.1;
2058         fLineCutZRSlopeMin = 0.;
2059         break;
2060     case 1:     // 1.2
2061         fEtaCut         = 1.2;
2062         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2063         fEtaCutMin              = -0.1;
2064         fLineCutZRSlopeMin = 0.;
2065         break;
2066     case 2:     // 1.4
2067         fEtaCut         = 1.4;
2068         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2069         fEtaCutMin              = -0.1;
2070         fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCut)));
2071         break;
2072     case 3: // 0.8
2073         fEtaCut         = 0.8;
2074         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2075         fEtaCutMin              = -0.1;
2076         fLineCutZRSlopeMin = 0.;
2077         break;
2078     case 4: // 0.75
2079         fEtaCut         = 0.75;
2080         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2081         fEtaCutMin              = -0.1;
2082         fLineCutZRSlopeMin = 0.;
2083         break;
2084     case 5: // 0.9 - 1.4
2085         fEtaCut         = 1.4;
2086         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2087         fEtaCutMin              = 0.9;
2088         fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCutMin)));
2089         break;
2090     case 6: // 0.9 - 1.2
2091         fEtaCut         = 1.2;
2092         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2093         fEtaCutMin              = 0.9;
2094         fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCutMin)));
2095         break;
2096     case 7: // 0.1 - 0.8
2097         fEtaCut         = 0.8;
2098         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2099         fEtaCutMin              = 0.1;
2100         fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCutMin)));
2101         break;
2102     case 8: // 0.1 - 0.8
2103         fEtaCut         = 0.9;
2104         fLineCutZRSlope = tan(2*atan(exp(-fEtaCut)));
2105         fEtaCutMin              = 0.1;
2106         fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCutMin)));
2107         break;
2108     default:
2109         AliError(Form(" EtaCut not defined %d",etaCut));
2110         return kFALSE;
2111     }
2112     return kTRUE;
2113 }
2114
2115 ///________________________________________________________________________
2116 Bool_t AliConversionCuts::SetChi2MesonCut(Int_t chi2MesonCut)
2117 {   // Set Cut
2118     switch(chi2MesonCut){
2119     case 0:  // 100.
2120         fChi2CutMeson = 100.;
2121         break;
2122     case 1:  // 50.
2123         fChi2CutMeson = 50.;
2124         break;
2125     case 2:  // 30.
2126         fChi2CutMeson = 30.;
2127         break;
2128     case 3:
2129         fChi2CutMeson = 200.;
2130         break;
2131     case 4:
2132         fChi2CutMeson = 500.;
2133         break;
2134     case 5:
2135         fChi2CutMeson = 1000.;
2136         break;
2137     default:
2138         AliError(Form("Chi2MesonCut not defined %d",chi2MesonCut));
2139         return kFALSE;
2140     }
2141     return kTRUE;
2142 }
2143
2144 ///________________________________________________________________________
2145 Bool_t AliConversionCuts::SetMaxMomPiondEdxCut(Int_t piMaxMomdedxSigmaCut)
2146 {   // Set Cut
2147     switch(piMaxMomdedxSigmaCut){
2148     case 0:  // 100. GeV
2149         fPIDMaxPnSigmaAbovePionLine=100.;
2150         break;
2151     case 1:  // 5. GeV
2152         fPIDMaxPnSigmaAbovePionLine=5.;
2153         break;
2154     case 2:  // 4. GeV
2155         fPIDMaxPnSigmaAbovePionLine=4.;
2156         break;
2157     case 3:  // 3.5 GeV
2158         fPIDMaxPnSigmaAbovePionLine=3.5;
2159         break;
2160     case 4:  // 3. GeV
2161         fPIDMaxPnSigmaAbovePionLine=3.;
2162         break;
2163     default:
2164         AliError(Form("piMaxMomdedxSigmaCut not defined %d",piMaxMomdedxSigmaCut));
2165         return kFALSE;
2166     }
2167     return kTRUE;
2168 }
2169
2170 ///________________________________________________________________________
2171 Bool_t AliConversionCuts::SetIsHeavyIon(Int_t isHeavyIon)
2172 {   // Set Cut
2173    switch(isHeavyIon){
2174                 case 0:
2175                         fIsHeavyIon=0;
2176                         break;
2177                 case 1:
2178                         fIsHeavyIon=1;
2179                         fDetectorCentrality=0;
2180                         break;
2181                 case 2:
2182                         fIsHeavyIon=1;
2183                         fDetectorCentrality=1;
2184                         break;
2185                 case 3: //allows to select centrality 0-45% in steps of 5% for V0 Multiplicity
2186                         fIsHeavyIon=1;
2187                         fDetectorCentrality=0;    
2188                         fModCentralityClass=1;
2189                         break;
2190                 case 4: //allows to select centrality 45-90% in steps of 5% for V0 Multiplicity
2191                         fIsHeavyIon=1;
2192                         fDetectorCentrality=0;
2193                         fModCentralityClass=2;
2194                         break;
2195                 default:
2196                         AliError(Form("SetHeavyIon not defined %d",isHeavyIon));
2197                         return kFALSE;
2198    }
2199    return kTRUE;
2200 }
2201
2202 ///________________________________________________________________________
2203 Bool_t AliConversionCuts::SetAlphaMesonCut(Int_t alphaMesonCut)
2204 {   // Set Cut
2205     switch(alphaMesonCut){
2206     case 0:     // 0- 0.7
2207         fAlphaMinCutMeson        = 0.0;
2208         fAlphaCutMeson   = 0.7;
2209         break;
2210     case 1:     // 0-0.5
2211         fAlphaMinCutMeson        = 0.0;
2212         fAlphaCutMeson   = 0.5;
2213         break;
2214     case 2:     // 0.5-1
2215         fAlphaMinCutMeson        = 0.5;
2216         fAlphaCutMeson   = 1.;
2217         break;
2218     case 3:     // 0.0-1
2219         fAlphaMinCutMeson        = 0.0;
2220         fAlphaCutMeson   = 1.;
2221         break;
2222     case 4:     // 0-0.65
2223         fAlphaMinCutMeson        = 0.0;
2224         fAlphaCutMeson   = 0.65;
2225         break;
2226     case 5:     // 0-0.75
2227         fAlphaMinCutMeson        = 0.0;
2228         fAlphaCutMeson   = 0.75;
2229         break;
2230     case 6:     // 0-0.8
2231         fAlphaMinCutMeson        = 0.0;
2232         fAlphaCutMeson   = 0.8;
2233         break;
2234     case 7:     // 0.0-0.85
2235         fAlphaMinCutMeson        = 0.0;
2236         fAlphaCutMeson   = 0.85;
2237         break;
2238     case 8:     // 0.0-0.6
2239         fAlphaMinCutMeson        = 0.0;
2240         fAlphaCutMeson   = 0.6;
2241         break;
2242     default:
2243         AliError(Form("AlphaMesonCut not defined %d",alphaMesonCut));
2244         return kFALSE;
2245     }
2246     return kTRUE;
2247 }
2248
2249 ///________________________________________________________________________
2250 Bool_t AliConversionCuts::SetRapidityMesonCut(Int_t RapidityMesonCut)
2251 {   // Set Cut
2252     switch(RapidityMesonCut){
2253     case 0:  //
2254         fRapidityCutMeson   = 0.9;
2255         break;
2256     case 1:  //
2257         fRapidityCutMeson   = 0.8;
2258         break;
2259     case 2:  //
2260         fRapidityCutMeson   = 0.7;
2261         break;
2262
2263     default:
2264         AliError(Form("RapidityMesonCut not defined %d",RapidityMesonCut));
2265         return kFALSE;
2266     }
2267     return kTRUE;
2268 }
2269
2270 ///________________________________________________________________________
2271 Bool_t AliConversionCuts::SetLowPRejectionCuts(Int_t LowPRejectionSigmaCut)
2272 {   // Set Cut
2273     switch(LowPRejectionSigmaCut){
2274         case 0:  //
2275         fPIDnSigmaAtLowPAroundKaonLine=0;
2276         fPIDnSigmaAtLowPAroundProtonLine=0;
2277         fPIDnSigmaAtLowPAroundPionLine=0;
2278         break;
2279     case 1:  //
2280         fPIDnSigmaAtLowPAroundKaonLine=0.5;
2281         fPIDnSigmaAtLowPAroundProtonLine=0.5;
2282         fPIDnSigmaAtLowPAroundPionLine=0.5;
2283         break;
2284     case 2:  //
2285         fPIDnSigmaAtLowPAroundKaonLine=1;
2286         fPIDnSigmaAtLowPAroundProtonLine=1;
2287         fPIDnSigmaAtLowPAroundPionLine=1;
2288         break;
2289     case 3:  //
2290         fPIDnSigmaAtLowPAroundKaonLine=2.;
2291         fPIDnSigmaAtLowPAroundProtonLine=2.;
2292         fPIDnSigmaAtLowPAroundPionLine=2.;
2293         break;
2294     case 4:  //
2295         fPIDnSigmaAtLowPAroundKaonLine=0.;
2296         fPIDnSigmaAtLowPAroundProtonLine=0.;
2297         fPIDnSigmaAtLowPAroundPionLine=1;
2298         break;
2299     case 5:  //
2300         fPIDnSigmaAtLowPAroundKaonLine=0.;
2301         fPIDnSigmaAtLowPAroundProtonLine=0.;
2302         fPIDnSigmaAtLowPAroundPionLine=1.5;
2303         break;
2304     case 6:  //
2305         fPIDnSigmaAtLowPAroundKaonLine=0.;
2306         fPIDnSigmaAtLowPAroundProtonLine=0.;
2307         fPIDnSigmaAtLowPAroundPionLine=2.;
2308         break;
2309     default:
2310         AliError(Form("LowPRejectionSigmaCut not defined %d",LowPRejectionSigmaCut));
2311         return kFALSE;
2312     }
2313     return kTRUE;
2314 }
2315
2316 ///________________________________________________________________________
2317 Bool_t AliConversionCuts::SetTOFElectronPIDCut(Int_t TOFelectronPID){
2318     // Set Cut
2319     switch(TOFelectronPID){ // RRnewTOF start //////////////////////////////////////////////////////////////////////////
2320     case 0: // no cut
2321         fUseTOFpid = kFALSE;
2322         fTofPIDnSigmaBelowElectronLine=-100;
2323         fTofPIDnSigmaAboveElectronLine=100;
2324         break;
2325     case 1: // -7,7
2326         fUseTOFpid = kTRUE;
2327         fTofPIDnSigmaBelowElectronLine=-7;
2328         fTofPIDnSigmaAboveElectronLine=7;
2329         break;
2330     case 2: // -5,5
2331         fUseTOFpid = kTRUE;
2332         fTofPIDnSigmaBelowElectronLine=-5;
2333         fTofPIDnSigmaAboveElectronLine=5;
2334         break;
2335     case 3: // -3,5
2336         fUseTOFpid = kTRUE;
2337         fTofPIDnSigmaBelowElectronLine=-3;
2338         fTofPIDnSigmaAboveElectronLine=5;
2339         break;
2340     case 4: // -2,3
2341         fUseTOFpid = kTRUE;
2342         fTofPIDnSigmaBelowElectronLine=-2;
2343         fTofPIDnSigmaAboveElectronLine=3;
2344         break;
2345     default:
2346         AliError(Form("TOFElectronCut not defined %d",TOFelectronPID));
2347         return kFALSE;
2348     } //////////////////////// RRnewTOF end //////////////////////////////////////////////////////////////////////////
2349     return kTRUE;
2350 }
2351
2352 ///________________________________________________________________________
2353 Bool_t AliConversionCuts::SetTRDElectronCut(Int_t TRDElectronCut)
2354 {   // Set Cut
2355     switch(TRDElectronCut){
2356     case 0:
2357         fDoTRDPID=kFALSE;
2358         break;
2359     case 1:
2360         fDoTRDPID=kTRUE;
2361         fPIDTRDEfficiency=0.1;
2362         break;
2363     case 8:
2364         fDoTRDPID=kTRUE;
2365         fPIDTRDEfficiency=0.8;
2366         break;
2367     case 9:
2368         fDoTRDPID=kTRUE;
2369         fPIDTRDEfficiency=0.9;
2370         break;
2371     default:
2372         AliError(Form("TRDElectronCut not defined %d",TRDElectronCut));
2373         return kFALSE;
2374     }
2375
2376     return kTRUE;
2377 }
2378 ///________________________________________________________________________
2379 Bool_t AliConversionCuts::SetQtMaxCut(Int_t QtMaxCut)
2380 {   // Set Cut
2381     switch(QtMaxCut){
2382     case 0: //
2383         fQtMax=1.;
2384         fDoQtGammaSelection=kFALSE;      //No Qt selection (true by default)
2385         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2386         fHighPtQtMax=100.;              // RRnew
2387         fPtBorderForQt=100.;            // RRnew
2388         break;
2389     case 1:
2390         fQtMax=0.1;
2391         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2392         fHighPtQtMax=100.;              // RRnew
2393         fPtBorderForQt=100.;            // RRnew
2394         break;
2395     case 2:
2396         fQtMax=0.07;
2397         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2398         fHighPtQtMax=100.;              // RRnew
2399         fPtBorderForQt=100.;            // RRnew
2400         break;
2401     case 3:
2402         fQtMax=0.05;
2403         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2404         fHighPtQtMax=100.;              // RRnew
2405         fPtBorderForQt=100.;            // RRnew
2406         break;
2407     case 4:
2408         fQtMax=0.03;
2409         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2410         fHighPtQtMax=100.;              // RRnew
2411         fPtBorderForQt=100.;            // RRnew
2412         break;
2413     case 5: // RR try to improve (get rid of) low InvMass peak in PbPb
2414         fQtMax=0.02;
2415         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2416         fHighPtQtMax=100.;              // RRnew
2417         fPtBorderForQt=100.;            // RRnew
2418         break; // end RR ///////////////////////////////////////////////
2419     case 6:  // RRnew start: pT dependent qT cut
2420         fQtMax=0.02;
2421         fDoHighPtQtGammaSelection=kTRUE;
2422         fHighPtQtMax=0.06;
2423         fPtBorderForQt=2.5;
2424         break; // RRnew end ////////////////////////////////////////////
2425     case 7:
2426         fQtMax=0.15;
2427         fDoHighPtQtGammaSelection=kFALSE; // RRnew
2428         fHighPtQtMax=100.;              // RRnew
2429         fPtBorderForQt=100.;            // RRnew
2430         break;
2431     default:
2432         AliError(Form("Warning: QtMaxCut not defined %d",QtMaxCut));
2433         return kFALSE;
2434     }
2435     return kTRUE;
2436 }
2437
2438 //-------------------------------------------------------------
2439 Double_t AliConversionCuts::GetCentrality(AliVEvent *event)
2440 {   // Get Event Centrality
2441
2442     AliESDEvent *esdEvent=dynamic_cast<AliESDEvent*>(event);
2443     if(esdEvent){
2444         AliCentrality *fESDCentrality=(AliCentrality*)esdEvent->GetCentrality();
2445
2446         if(fDetectorCentrality==0){
2447             return fESDCentrality->GetCentralityPercentile("V0M"); // default
2448         }
2449         if(fDetectorCentrality==1){
2450             return fESDCentrality->GetCentralityPercentile("CL1");
2451         }
2452     }
2453
2454     AliAODEvent *aodEvent=dynamic_cast<AliAODEvent*>(event);
2455     if(aodEvent){
2456         if(aodEvent->GetHeader()){return aodEvent->GetHeader()->GetCentrality();}
2457     }
2458
2459     return -1;
2460 }
2461
2462 //-------------------------------------------------------------
2463 Bool_t AliConversionCuts::IsCentralitySelected(AliVEvent *event)
2464 {   // Centrality Selection
2465    if(!fIsHeavyIon)return kTRUE;
2466     
2467         if(fCentralityMin == 0 && fCentralityMax == 0) return kTRUE;//0-100%
2468         if(fCentralityMin >= fCentralityMax) return kTRUE;//0-100%
2469    
2470         Double_t centrality=GetCentrality(event);
2471         if(centrality<0)return kFALSE;
2472         
2473         Int_t centralityC=0;
2474         if (fModCentralityClass == 0){
2475                 centralityC= Int_t(centrality/10);
2476                 if(centralityC >= fCentralityMin && centralityC < fCentralityMax)
2477                         return kTRUE;
2478                 else 
2479                         return kFALSE;
2480         }       else if (fModCentralityClass ==1){
2481                 centralityC= Int_t(centrality);
2482                 if(centralityC >= fCentralityMin*5 && centralityC < fCentralityMax*5){
2483                         return kTRUE;
2484                 } else { 
2485                         return kFALSE;
2486                 }
2487         } else if (fModCentralityClass ==2){
2488                 centralityC= Int_t(centrality+1);
2489                 if(centralityC >= (fCentralityMin*5+45) && centralityC < (fCentralityMax*5+45))
2490                         return kTRUE;
2491                 else 
2492                         return kFALSE;
2493         }
2494         return kFALSE;
2495 }
2496
2497 //-------------------------------------------------------------
2498 Bool_t AliConversionCuts::SetCentralityMin(Int_t minCentrality)
2499 {
2500     // Set Cut
2501     if(minCentrality<0||minCentrality>9){
2502         AliError(Form("minCentrality not defined %d",minCentrality));
2503         return kFALSE;
2504     }
2505
2506     fCentralityMin=minCentrality;
2507     return kTRUE;
2508 }
2509 //-------------------------------------------------------------
2510 Bool_t AliConversionCuts::SetCentralityMax(Int_t maxCentrality)
2511 {
2512     // Set Cut
2513     if(maxCentrality<0||maxCentrality>9){
2514         AliError(Form("maxCentrality not defined %d",maxCentrality));
2515         return kFALSE;
2516     }
2517
2518     fCentralityMax=maxCentrality;
2519     return kTRUE;
2520 }
2521
2522 ///________________________________________________________________________
2523 Bool_t AliConversionCuts::SetPhotonAsymmetryCut(Int_t doPhotonAsymmetryCut){
2524     // Set Cut
2525     switch(doPhotonAsymmetryCut){
2526     case 0:
2527         fDoPhotonAsymmetryCut=0;
2528         fMinPPhotonAsymmetryCut=100.;
2529         fMinPhotonAsymmetry=0.;
2530         break;
2531     case 1:
2532         fDoPhotonAsymmetryCut=1;
2533         fMinPPhotonAsymmetryCut=3.5;
2534         fMinPhotonAsymmetry=0.04;
2535         break;
2536     case 2:
2537         fDoPhotonAsymmetryCut=1;
2538         fMinPPhotonAsymmetryCut=3.5;
2539         fMinPhotonAsymmetry=0.06;
2540         break;
2541     default:
2542         AliError(Form("PhotonAsymmetryCut not defined %d",doPhotonAsymmetryCut));
2543         return kFALSE;
2544     }
2545     fCuts[kdoPhotonAsymmetryCut]=doPhotonAsymmetryCut;
2546     return kTRUE;
2547 }
2548
2549 ///________________________________________________________________________
2550 Bool_t AliConversionCuts::SetBackgroundScheme(Int_t BackgroundScheme){
2551     // Set Cut
2552     switch(BackgroundScheme){
2553     case 0: //Rotation
2554         fUseRotationMethodInBG=kTRUE;
2555         fdoBGProbability=kFALSE;
2556         break;
2557     case 1: // mixed event with V0 multiplicity
2558         fUseRotationMethodInBG=kFALSE;
2559         fUseTrackMultiplicityForBG=kFALSE;
2560         fdoBGProbability=kFALSE;
2561         break;
2562     case 2: // mixed event with track multiplicity
2563         fUseRotationMethodInBG=kFALSE;
2564         fUseTrackMultiplicityForBG=kTRUE;
2565         fdoBGProbability=kFALSE;
2566         break;
2567     case 3: //Rotation
2568         fUseRotationMethodInBG=kTRUE;
2569         fdoBGProbability=kTRUE;
2570         break;
2571     default:
2572         AliError(Form("BackgroundScheme not defined %d",BackgroundScheme));
2573         return kFALSE;
2574     }
2575     return kTRUE;
2576 }
2577
2578 ///________________________________________________________________________
2579 Bool_t AliConversionCuts::SetNDegreesForRotationMethod(Int_t DegreesForRotationMethod){
2580     // Set Cut
2581     switch(DegreesForRotationMethod){
2582     case 0:
2583         fnDegreeRotationPMForBG = 5;
2584         break;
2585     case 1:
2586         fnDegreeRotationPMForBG = 10;
2587         break;
2588     case 2:
2589         fnDegreeRotationPMForBG = 15;
2590         break;
2591     case 3:
2592         fnDegreeRotationPMForBG = 20;
2593         break;
2594     default:
2595         AliError(Form("DegreesForRotationMethod not defined %d",DegreesForRotationMethod));
2596         return kFALSE;
2597     }
2598     fCuts[kDegreesForRotationMethod]=DegreesForRotationMethod;
2599     return kTRUE;
2600 }
2601
2602 ///________________________________________________________________________
2603 Bool_t AliConversionCuts::SetNumberOfRotations(Int_t NumberOfRotations)
2604 {   // Set Cut
2605     switch(NumberOfRotations){
2606     case 0:
2607         fnumberOfRotationEventsForBG = 5;
2608         break;
2609     case 1:
2610         fnumberOfRotationEventsForBG = 10;
2611         break;
2612     case 2:
2613         fnumberOfRotationEventsForBG = 15;
2614         break;
2615     case 3:
2616         fnumberOfRotationEventsForBG = 20;
2617         break;
2618     case 4:
2619         fnumberOfRotationEventsForBG = 2;
2620         break;
2621     case 5:
2622         fnumberOfRotationEventsForBG = 50;
2623         break;
2624     case 6:
2625         fnumberOfRotationEventsForBG = 80;
2626         break;
2627     case 7:
2628         fnumberOfRotationEventsForBG = 100;
2629         break;
2630     default:
2631         AliError(Form("NumberOfRotations not defined %d",NumberOfRotations));
2632         return kFALSE;
2633     }
2634     return kTRUE;
2635 }
2636
2637 ///________________________________________________________________________
2638 Bool_t AliConversionCuts::SetPsiPairCut(Int_t psiCut) {
2639   
2640
2641   switch(psiCut) {
2642   case 0:
2643         fPsiPairCut = 10000; // 
2644         break;
2645   case 1:
2646         fPsiPairCut = 0.1; // 
2647         break;
2648   case 2:
2649         fPsiPairCut = 0.05; // Standard
2650         break;
2651   case 3:
2652         fPsiPairCut = 0.035; // 
2653         break;
2654   case 4:
2655         fPsiPairCut = 0.15; // 
2656         break;
2657   case 5:
2658         fPsiPairCut = 0.2; // 
2659         break;
2660   case 6:
2661         fPsiPairCut = 0.03; // 
2662         break;
2663   case 7:
2664         fPsiPairCut = 0.025; // 
2665         break;
2666   case 8:
2667         fPsiPairCut = 0.01; // 
2668         break;
2669   default:
2670       AliError(Form("PsiPairCut not defined %d",psiCut));
2671       return kFALSE;
2672   }
2673
2674   return kTRUE;
2675 }
2676
2677 ///________________________________________________________________________
2678 Bool_t AliConversionCuts::SetCosPAngleCut(Int_t cosCut) {
2679
2680     switch(cosCut){
2681     case 0:
2682         fCosPAngleCut = TMath::Pi(); //
2683         break;
2684     case 1:
2685         fCosPAngleCut = 0.1; //
2686         break;
2687     case 2:
2688         fCosPAngleCut = 0.05; //
2689         break;
2690     case 3:
2691         fCosPAngleCut = 0.025; // Standard
2692         break;
2693     case 4:
2694         fCosPAngleCut = 0.01; //
2695         break;
2696     default:
2697         AliError(Form("Cosine Pointing Angle cut not defined %d",cosCut));
2698         return kFALSE;
2699     }
2700         
2701     return kTRUE;
2702 }
2703
2704
2705 ///________________________________________________________________________
2706 Bool_t AliConversionCuts::VertexZCut(AliVEvent *event){
2707     // Cut on z position of primary vertex
2708      Double_t fVertexZ=event->GetPrimaryVertex()->GetZ();
2709
2710      if(TMath::Abs(fVertexZ)>fMaxVertexZ)return kFALSE;
2711      return kTRUE;
2712 }
2713
2714 ///________________________________________________________________________
2715 Bool_t AliConversionCuts::SetSharedElectronCut(Int_t sharedElec) {
2716
2717    switch(sharedElec){
2718     case 0:
2719        fDoSharedElecCut = kFALSE;
2720         break;
2721     case 1:
2722        fDoSharedElecCut = kTRUE;
2723         break;
2724     default:
2725         AliError(Form("Shared Electron Cut not defined %d",sharedElec));
2726         return kFALSE;
2727     }
2728         
2729     return kTRUE;
2730 }
2731
2732 ///________________________________________________________________________
2733 Bool_t AliConversionCuts::SetToCloseV0sCut(Int_t toClose) {
2734
2735    switch(toClose){
2736    case 0:
2737       fDoToCloseV0sCut = kFALSE;
2738       fminV0Dist = 250;
2739       break;
2740    case 1:
2741       fDoToCloseV0sCut = kTRUE;
2742       fminV0Dist = 1;
2743       break;
2744    case 2:
2745       fDoToCloseV0sCut = kTRUE;
2746       fminV0Dist = 2;
2747       break;
2748    case 3:
2749       fDoToCloseV0sCut = kTRUE;
2750       fminV0Dist = 3;
2751       break;
2752    default:
2753        AliError(Form("Shared Electron Cut not defined %d",toClose));
2754         return kFALSE;
2755    }
2756    return kTRUE;
2757 }
2758
2759
2760
2761 ///________________________________________________________________________
2762
2763 Int_t AliConversionCuts::GetNumberOfContributorsVtx(AliVEvent *event){
2764     // returns number of contributors to the vertex
2765     
2766     AliESDEvent *fESDEvent=dynamic_cast<AliESDEvent*>(event);
2767     if(fESDEvent){
2768         if(fESDEvent->GetPrimaryVertexTracks()->GetNContributors()>0) {
2769             return fESDEvent->GetPrimaryVertexTracks()->GetNContributors();
2770         }
2771      
2772         if(fESDEvent->GetPrimaryVertexTracks()->GetNContributors()<1) {
2773             //          return 0;
2774             //-AM test pi0s without SPD only vertex
2775             if(fESDEvent->GetPrimaryVertexSPD()->GetNContributors()>0) {
2776                 return fESDEvent->GetPrimaryVertexSPD()->GetNContributors();
2777
2778             }
2779             if(fESDEvent->GetPrimaryVertexSPD()->GetNContributors()<1) {
2780                 return 0;
2781             }
2782         }
2783     }
2784
2785     AliAODEvent *fAODEvent=dynamic_cast<AliAODEvent*>(event);
2786     if(fAODEvent){
2787         if(fAODEvent->GetPrimaryVertex()->GetNContributors()>0) {
2788             return fAODEvent->GetPrimaryVertex()->GetNContributors();
2789         }
2790         if(fAODEvent->GetPrimaryVertex()->GetNContributors()<1) {
2791             if(fAODEvent->GetPrimaryVertexSPD()->GetNContributors()>0) {
2792                 return fAODEvent->GetPrimaryVertexSPD()->GetNContributors();
2793             }
2794             if(fAODEvent->GetPrimaryVertexSPD()->GetNContributors()<1) {
2795                 AliWarning(Form("Number of contributors from bad vertex type:: %s",fAODEvent->GetPrimaryVertex()->GetName()));
2796                 return 0;
2797             }
2798         }
2799     }
2800
2801
2802   return 0;
2803 }
2804
2805 ///________________________________________________________________________
2806
2807 Bool_t AliConversionCuts::IsTriggerSelected()
2808 {
2809     AliInputEventHandler *fInputHandler=(AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
2810
2811     UInt_t isSelected = AliVEvent::kAny;
2812     if( fInputHandler && fInputHandler->GetEventSelection()) {
2813         // Get the actual offline trigger mask for the event and AND it with the
2814         // requested mask. If no mask requested select by default the event.
2815         if (fOfflineTriggerMask)
2816             isSelected = fOfflineTriggerMask & fInputHandler->IsEventSelected();
2817     }
2818    
2819     if(!isSelected)return kFALSE;
2820
2821     // Fill Histogram
2822     if(hTriggerClass){
2823         if (fInputHandler->IsEventSelected() & AliVEvent::kAny)hTriggerClass->Fill(0);
2824         if (fInputHandler->IsEventSelected() & AliVEvent::kMB)hTriggerClass->Fill(1);
2825         if (fInputHandler->IsEventSelected() & AliVEvent::kCentral)hTriggerClass->Fill(2);
2826         if (fInputHandler->IsEventSelected() & AliVEvent::kSemiCentral)hTriggerClass->Fill(3);
2827     }
2828
2829     
2830     return kTRUE;
2831
2832 }
2833
2834 ///________________________________________________________________________
2835 Int_t AliConversionCuts::GetFirstTPCRow(Double_t radius){
2836     // Get first TPC row
2837         Int_t firstTPCRow=0;
2838         Double_t radiusI        =       84.8;
2839         Double_t radiusO        = 134.6;
2840         Double_t radiusOB = 198.;
2841         Double_t rSizeI  = 0.75;
2842         Double_t rSizeO  = 1.;
2843         Double_t rSizeOB        = 1.5;
2844         Int_t nClsI=63;
2845         Int_t nClsIO=127;
2846
2847         if(radius <= radiusI){
2848                 return firstTPCRow;
2849         }
2850         if(radius>radiusI && radius<=radiusO){
2851                 firstTPCRow = (Int_t)((radius-radiusI)/rSizeI);
2852         }
2853         if(radius>radiusO && radius<=radiusOB){
2854                 firstTPCRow = (Int_t)(nClsI+(radius-radiusO)/rSizeO);
2855         }
2856
2857         if(radius>radiusOB){
2858                 firstTPCRow =(Int_t)(nClsIO+(radius-radiusOB)/rSizeOB);
2859         }
2860
2861         return firstTPCRow;
2862 }
2863
2864 Bool_t AliConversionCuts::CosinePAngleCut(const AliConversionPhotonBase * photon, AliVEvent * event) const {
2865   ///Check if passes cosine of pointing angle cut
2866    if(GetCosineOfPointingAngle(photon, event) < (TMath::Cos(fCosPAngleCut))){
2867       return kFALSE;
2868    }
2869    return kTRUE;
2870 }
2871
2872 Double_t AliConversionCuts::GetCosineOfPointingAngle( const AliConversionPhotonBase * photon, AliVEvent * event) const{
2873    // calculates the pointing angle of the recalculated V0 
2874
2875    Double_t momV0[3] = {0,0,0};
2876    if(event->IsA()==AliESDEvent::Class()){
2877       AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*>(event);
2878       AliESDv0 *v0 = esdEvent->GetV0(photon->GetV0Index());
2879       v0->GetPxPyPz(momV0[0],momV0[1],momV0[2]);
2880    }
2881    if(event->IsA()==AliAODEvent::Class()){
2882       momV0[0] = photon->GetPx();
2883       momV0[1] = photon->GetPy();
2884       momV0[2] = photon->GetPz();
2885    }
2886    
2887    //Double_t momV0[3] = { photon->GetPx(), photon->GetPy(), photon->GetPz() }; //momentum of the V0
2888    Double_t PosV0[3] = { photon->GetConversionX() - event->GetPrimaryVertex()->GetX(), 
2889                          photon->GetConversionY() - event->GetPrimaryVertex()->GetY(), 
2890                          photon->GetConversionZ() - event->GetPrimaryVertex()->GetZ() }; //Recalculated V0 Position vector
2891    
2892    Double_t momV02 = momV0[0]*momV0[0] + momV0[1]*momV0[1] + momV0[2]*momV0[2];
2893    Double_t PosV02 = PosV0[0]*PosV0[0] + PosV0[1]*PosV0[1] + PosV0[2]*PosV0[2];
2894
2895    Double_t cosinePointingAngle = (PosV0[0]*momV0[0] +  PosV0[1]*momV0[1] + PosV0[2]*momV0[2] ) / TMath::Sqrt(momV02 * PosV02);
2896   
2897    return cosinePointingAngle;
2898 }
2899
2900
2901 Bool_t AliConversionCuts::PsiPairCut(const AliConversionPhotonBase * photon) const {
2902
2903     if(photon->GetPsiPair() > fPsiPairCut){
2904         return kFALSE;}
2905     else{return kTRUE;}
2906
2907 }
2908
2909 ///________________________________________________________________________
2910 TString AliConversionCuts::GetCutNumber(){
2911     // returns TString with current cut number
2912   TString a(kNCuts);
2913   for(Int_t ii=0;ii<kNCuts;ii++){
2914       a.Append(Form("%d",fCuts[ii]));
2915   }
2916   return a;
2917 }
2918
2919 ///________________________________________________________________________
2920 void AliConversionCuts::SmearParticle(AliAODConversionPhoton* photon)
2921 {
2922    Double_t facPBrem = 1.;
2923    Double_t facPSig = 0.;
2924
2925    Double_t phi=0.;
2926    Double_t theta=0.;
2927    Double_t P=0.;
2928
2929    
2930    P=photon->P();
2931    phi=photon->Phi();
2932    if( photon->P()!=0){
2933       theta=acos( photon->Pz()/ photon->P());
2934    }
2935
2936    if( fPSigSmearing != 0. || fPSigSmearingCte!=0. ){ 
2937       facPSig = TMath::Sqrt(fPSigSmearingCte*fPSigSmearingCte+fPSigSmearing*fPSigSmearing*P*P)*fRandom.Gaus(0.,1.);
2938    }
2939         
2940    if( fPBremSmearing != 1.){
2941       if(fBrem!=NULL){
2942          facPBrem = fBrem->GetRandom();
2943       }
2944    }
2945
2946    photon->SetPx(facPBrem* (1+facPSig)* P*sin(theta)*cos(phi)) ;
2947    photon->SetPy(facPBrem* (1+facPSig)* P*sin(theta)*sin(phi)) ;
2948    photon->SetPz(facPBrem* (1+facPSig)* P*cos(theta)) ;
2949    photon->SetE(photon->P());
2950 }
2951 ///________________________________________________________________________
2952 void AliConversionCuts::FillElectonLabelArray(AliAODConversionPhoton* photon, Int_t nV0){
2953    
2954    Int_t posLabel = photon->GetTrackLabelPositive();
2955    Int_t negLabel = photon->GetTrackLabelNegative();
2956    
2957    fElectronLabelArray[nV0*2] = posLabel;
2958    fElectronLabelArray[(nV0*2)+1] = negLabel;
2959 }
2960 ///________________________________________________________________________
2961 Bool_t AliConversionCuts::RejectSharedElectronV0s(AliAODConversionPhoton* photon, Int_t nV0, Int_t nV0s){
2962
2963    Int_t posLabel = photon->GetTrackLabelPositive();
2964    Int_t negLabel = photon->GetTrackLabelNegative();
2965    
2966    for(Int_t i = 0; i<nV0s*2;i++){
2967       if(i==nV0*2)     continue;
2968       if(i==(nV0*2)+1) continue;
2969       if(fElectronLabelArray[i] == posLabel){
2970          return kFALSE;}
2971       if(fElectronLabelArray[i] == negLabel){
2972          return kFALSE;}
2973    }
2974
2975    return kTRUE;
2976 }
2977 ///________________________________________________________________________
2978 Bool_t AliConversionCuts::RejectToCloseV0s(AliAODConversionPhoton* photon, TList *photons, Int_t nV0){
2979
2980
2981    Double_t posX = photon->GetConversionX();
2982    Double_t posY = photon->GetConversionY();
2983    Double_t posZ = photon->GetConversionZ();
2984
2985    for(Int_t i = 0;i<photons->GetEntries();i++){
2986       if(nV0 == i) continue;
2987       AliAODConversionPhoton *photonComp = (AliAODConversionPhoton*) photons->At(i);
2988       Double_t posCompX = photonComp->GetConversionX();
2989       Double_t posCompY = photonComp->GetConversionY();
2990       Double_t posCompZ = photonComp->GetConversionZ();
2991       
2992       Double_t dist = pow((posX - posCompX),2)+pow((posY - posCompY),2)+pow((posZ - posCompZ),2);
2993
2994       if(dist < fminV0Dist*fminV0Dist){
2995          if(photon->GetChi2perNDF() < photonComp->GetChi2perNDF()) return kTRUE;
2996          else {
2997             return kFALSE;}
2998       }
2999       
3000    }
3001    return kTRUE;
3002 }