]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/GammaConv/AliV0Reader.h
Fixed warnings regarding par file compilation
[u/mrichter/AliRoot.git] / PWG4 / GammaConv / AliV0Reader.h
1 #ifndef ALIV0READER_H
2 #define ALIV0READER_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5
6 ////////////////////////////////////////////////
7 //--------------------------------------------- 
8 // Class used to do analysis on conversion pairs
9 //---------------------------------------------
10 ////////////////////////////////////////////////
11
12 // --- ROOT system ---
13 #include "TObject.h" 
14 #include "AliMCEvent.h"   // for CF
15 #include "AliESDv0.h"
16 #include "AliESDEvent.h"
17 #include "AliKFParticle.h"
18 #include "TParticle.h"
19 #include "AliGammaConversionHistograms.h"
20 #include <vector>
21 #include "AliCFManager.h"
22 #include "AliGammaConversionBGHandler.h"
23 #include "AliESDpid.h"
24
25 class TClonesArray; 
26 class TFormula;
27 class Riostream;
28 class TChain;
29
30 //--- AliRoot system ---
31
32 class AliStack;
33 class AliMCEvent;       // for CF
34 class AliESDEvent; 
35 class AliMCEventHandler;
36 class AliESDInputHandler;
37 class AliESDVertex;
38 class AliLog;
39 class TChain;
40 class TChain;
41 class AliCFManager;   // for CF
42 class AliCFContainer;  // for CF
43 //class AliESDpid; // for dEdx cut based on nSigma to particle lines 
44 class AliESDtrackCuts; 
45
46 class AliV0Reader : public TObject {
47         
48  public: 
49         
50         
51   // for CF
52   enum{
53     kStepGenerated = 0,
54     kStepReconstructable = 1, 
55     kStepGetOnFly = 2,
56     kStepLikeSign = 3,
57     kStepTPCRefit = 4,
58     kStepKinks = 5,
59     kStepdEdx_electronselection = 6,
60     kStepdEdx_pionrejection = 7,
61     kStepNContributors = 8,
62     kStepTPCPID = 9,
63     kStepR = 10,
64     kStepLine = 11,
65     kStepZ = 12,
66     kStepMinClsTPC = 13,
67     kStepSinglePt= 14,  
68     kStepNDF = 15,
69     kStepChi2 = 16,
70     kStepEta = 17,
71     kStepPt = 18,
72     kStepTrueGamma = 19
73   };
74         
75   AliV0Reader();                                        //constructor
76   AliV0Reader(const AliV0Reader & g);                   //copy constructor
77   AliV0Reader & operator = (const AliV0Reader & g);     //assignment operator
78   //  virtual ~AliV0Reader() {;}                            //virtual destructor
79   virtual ~AliV0Reader();                            //virtual destructor
80   /*
81    *Initialize the reader
82    */
83   void Initialize();
84   void SetInputAndMCEvent(AliVEvent* esd, AliMCEvent* mc) ;
85
86
87   virtual void SetInputEvent(AliVEvent* const input)  {fESDEvent  = dynamic_cast<AliESDEvent*>(input);}
88   virtual void SetMC(AliMCEvent* const mc)            {fMCEvent = mc;}
89
90         
91   // for CF
92   void SetCFManager(AliCFManager * const io){fCFManager = io;};
93   AliCFManager *GetCFManager() const {return fCFManager;}
94         
95         
96         
97         
98   /*
99    * Returns  AliESDEvent
100    */                   
101   AliESDEvent* GetESDEvent() const{return fESDEvent;}   
102         
103   /*
104    *Returns the number of v0s in the event, no cuts applied.
105    */
106   Int_t GetNumberOfV0s() const{return fESDEvent->GetNumberOfV0s();}
107         
108   /*
109    *Returns the number of contributors to the vertex
110    */
111   Int_t GetNumberOfContributorsVtx() const{return fESDEvent->GetPrimaryVertex()->GetNContributors();}
112   
113   /*
114    * Check if there are any more good v0s left in the v0 stack
115    * if so, fCurrent v0 is set to this v0 and can be retrieved
116    * by GetCurrentV0 function.
117    * returns kFALSE if there is no more good v0s in the v0 stack
118    */
119   Bool_t NextV0();
120         
121   /*
122    * Returns the v0 at the given index, no checks are done on the v0. 
123    */
124   AliESDv0* GetV0(Int_t index);
125         
126   /*
127    * Returns the current v0
128    */
129   AliESDv0* GetCurrentV0() const{return fCurrentV0;}
130         
131   /*
132    * Returns the negative ESD track which belongs to fCurrentV0
133    */
134   //  AliESDtrack* GetNegativeESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetNindex());}
135   AliESDtrack* GetNegativeESDTrack(){return fCurrentNegativeESDTrack;}
136         
137   /*
138    * Returns the positive ESD track which belongs to fCurrentV0
139    */
140   //  AliESDtrack* GetPositiveESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetPindex());}
141   AliESDtrack* GetPositiveESDTrack(){return fCurrentPositiveESDTrack;}
142         
143   /*
144    * Returns the negative KF particle which belongs to fCurrentV0
145    */
146   AliKFParticle* GetNegativeKFParticle() const{return fCurrentNegativeKFParticle;}
147         
148   /*
149    * Returns the positive KF particle which belongs to fCurrentV0
150    */
151   AliKFParticle* GetPositiveKFParticle() const{return fCurrentPositiveKFParticle;}
152         
153   /*
154    * Returns the KFParticle object of the 2 tracks.
155    */
156   AliKFParticle* GetMotherCandidateKFCombination() const{return fCurrentMotherKFCandidate;}
157         
158   /*
159    * Checks the probablity that the PID of the particle is what we want it to be.
160    */
161   Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);
162         
163   /*
164    * Checks if the PID of the two particles are within our cuts.
165    */
166   void GetPIDProbability(Double_t &negPIDProb, Double_t &posPIDProb);
167
168   /*
169    * Checks if the PID of the two particles are within our cuts.
170    */
171   void GetPIDProbabilityMuonPion(Double_t &negPIDProb, Double_t &posPIDProb);
172         
173   /*
174    *Get the negative MC TParticle from the stack 
175    */
176   TParticle * GetNegativeMCParticle() const{return fNegativeMCParticle;}
177         
178   /*
179    *Get the positive MC TParticle from the stack 
180    */
181   TParticle * GetPositiveMCParticle() const{return fPositiveMCParticle;}
182         
183   /*
184    *Get the mother MC TParticle from the stack 
185    */
186   TParticle * GetMotherMCParticle() const{return fMotherMCParticle;}
187         
188   /*
189    * Flag to see if the v0 particles share the same mother
190    */
191   Bool_t HasSameMCMother();
192         
193         
194   /*
195    *Get the PID of the MC mother particle
196    */
197   Int_t GetMotherMCParticlePDGCode() const{if(fMotherMCParticle != NULL){ cout<<"MCParticle exists"<<endl;} return fMotherMCParticle->GetPdgCode();}
198         
199   /*
200    *Get the MC stack 
201    */
202   AliStack* GetMCStack() const{return fMCStack;}
203         
204         
205   /*
206    * Setup  AliMCEventHandler
207    */                   
208   //  AliMCEventHandler* GetMCTruth() const{return fMCTruth;}   // for CF
209         
210         
211   /*
212    *Get the MC stack 
213    */
214   AliMCEvent* GetMCEvent() const{return fMCEvent;}   // for CF
215         
216         
217   /*
218    *Get the magnetic field from the ESD event 
219    */
220   Double_t GetMagneticField() const{return fESDEvent->GetMagneticField();}
221         
222   /*
223    *Get the primary vertex from the esd event
224    */
225   const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}
226         
227   /*
228    * Set the PID of the negative track
229    */
230   void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}
231         
232   /*
233    * Set the PID of the positive track
234    */
235   void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}
236         
237   /*
238    * Set the flag to use the kfparticle class. Will also disable the use of esd tracks
239    */
240   void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}
241         
242   /*
243    *  Set the flag to use the esd track class. Will also disable the use of kf particles
244    */
245   void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}
246         
247   /*
248    *  Set the flag to use improved vertex or not
249    */
250   void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}
251         
252   /*
253    * Return the number in the species array belonging to the negative or positive track pid.
254    */
255   Int_t GetSpeciesIndex(Int_t chargeOfTrack);
256         
257   /*
258    * Return the x coordinate of the v0
259    */
260   Double_t GetX() const{return fCurrentXValue;}
261         
262   /*
263    * Return the y coordinate of the v0
264    */
265   Double_t GetY() const{return fCurrentYValue;}
266         
267   /*
268    * Return the Z coordinate of the v0
269    */
270   Double_t GetZ() const{return fCurrentZValue;}
271         
272   /*
273    * Return the radius of the v0
274    */
275   Double_t GetXYRadius() const{return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}
276         
277   /*
278    * Get the opening angle between the two tracks
279    */
280   Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}
281
282   /*
283    * Get the Cos Pointing angle between the two tracks
284    */
285   Double_t GetCosPointingAngle(){return fCurrentV0->GetV0CosineOfPointingAngle();}
286
287   /*
288    * Get the DCA between the two tracks
289    */
290   Double_t GetDcaDaughters(){return fCurrentV0->GetDcaV0Daughters();}
291
292   /*
293    * Get the Normalized DCA between the two tracks
294    */
295   Double_t GetNormDcaDistDaughters(){return fCurrentV0->GetDcaV0Daughters()/fCurrentV0->GetDistSigma();}
296
297   /*
298    * Get the Likelihood for a Conversion
299    */
300   Double_t GetLikelihoodAP(){return fCurrentV0->GetLikelihoodAP(0,0);}
301       
302   /*
303    * Gets the Energy of the negative track.
304    */
305   Double_t GetNegativeTrackEnergy() const{return fCurrentNegativeKFParticle->E();}
306         
307   /*
308    * Gets the Energy of the positive track.
309    */
310   Double_t GetPositiveTrackEnergy() const{return fCurrentPositiveKFParticle->E();}
311         
312   /*
313    * Gets the Energy of the mother candidate.
314    */
315   Double_t GetMotherCandidateEnergy() const{return fCurrentMotherKFCandidate->E();}
316         
317   /*
318    * Gets the Pt of the negative track.
319    */
320   Double_t GetNegativeTrackPt() const{return fNegativeTrackLorentzVector->Pt();}
321         
322   /*
323    * Gets the Pt of the positive track.
324    */
325   Double_t GetPositiveTrackPt() const{return fPositiveTrackLorentzVector->Pt();}
326         
327
328   /*
329    * Gets the Pt of the mother candidate.
330    */
331   Double_t GetMotherCandidatePt() const{return fMotherCandidateLorentzVector->Pt();}
332
333
334   /*
335    * Gets the P of the mother candidate.
336    */
337   Double_t GetMotherCandidateP() const{return fMotherCandidateLorentzVector->P();}
338         
339
340   /*
341    * Gets the Eta of the negative track.
342    */
343   Double_t GetNegativeTrackEta() const{return fNegativeTrackLorentzVector->Eta();}
344   /*
345    * Gets the Eta of the positive track.
346    */
347   Double_t GetPositiveTrackEta() const{return fPositiveTrackLorentzVector->Eta();}
348   /*
349    * Gets the Eta of the mother candidate.
350    */
351   Double_t GetMotherCandidateEta() const{return fMotherCandidateLorentzVector->Eta();}
352         
353   /*
354    * Gets the NDF of the mother candidate.
355    */
356   Double_t GetMotherCandidateNDF() const{return fCurrentMotherKFCandidate->GetNDF();}
357         
358   /*
359    * Gets the Chi2 of the mother candidate.
360    */
361   Double_t GetMotherCandidateChi2() const{return fCurrentMotherKFCandidate->GetChi2();}
362         
363   /*
364    * Gets the Mass of the mother candidate.
365    */
366   Double_t GetMotherCandidateMass() const{return fMotherCandidateKFMass;}
367         
368   /*
369    * Gets the Width of the mother candidate.
370    */
371   Double_t GetMotherCandidateWidth() const{return fMotherCandidateKFWidth;}
372         
373   /*
374    * Gets the Phi of the negative track.
375    */
376   Double_t GetNegativeTrackPhi() const;
377         
378   /*
379    * Gets the Phi of the positive track.
380    */
381   Double_t GetPositiveTrackPhi() const;
382         
383   /*
384    * Gets the Phi of the mother candidate.
385    */
386   Double_t GetMotherCandidatePhi() const;
387         
388   /*
389    * Gets the Rapidity of the mother candidate.
390    */
391   Double_t GetMotherCandidateRapidity() const;
392         
393
394   /*
395    * Gets the P of the negative track.
396    */
397   Double_t GetNegativeTrackP() const{return fNegativeTrackLorentzVector->P();}
398         
399   /*
400    * Gets the P of the positive track.
401    */
402   Double_t GetPositiveTrackP() const{return fPositiveTrackLorentzVector->P();}
403
404   /*
405    * Gets the dE/dx in the TPC of the negative track.
406    */
407   Double_t GetNegativeTrackTPCdEdx() const{return fCurrentNegativeESDTrack->GetTPCsignal();}
408         
409   /*
410    * Gets the dE/dx in the TPC of the positive track.
411    */
412   Double_t GetPositiveTrackTPCdEdx() const{return fCurrentPositiveESDTrack->GetTPCsignal();}
413
414   /*
415    * Gets the Number of the TPC clusters of the negative track.
416    */
417   Int_t GetNegativeTracknTPCClusters() const{return fCurrentNegativeESDTrack->GetNcls(1);}
418
419   /*
420    * Gets the Number of the TPC clusters of the positive track.
421    */
422   Int_t GetPositiveTracknTPCClusters() const{return fCurrentPositiveESDTrack->GetNcls(1);}
423         
424   /*
425    * Gets the Number of the ITS clusters of the negative track.
426    */
427   Int_t GetNegativeTracknITSClusters() const{return fCurrentNegativeESDTrack->GetNcls(0);}
428
429   /*
430    * Gets the Number of the ITS clusters of the positive track.
431    */
432   Int_t GetPositiveTracknITSClusters() const{return fCurrentPositiveESDTrack->GetNcls(0);}
433         
434   /*
435    * Update data which need to be updated every event.
436    */
437   void UpdateEventByEventData();
438
439   /*
440    * Gets the MaxRCut value.
441    */
442   Double_t GetMaxVertexZ() const{return fMaxVertexZ;}
443         
444   /*
445    * Gets the MaxRCut value.
446    */
447   Double_t GetMaxRCut() const{return fMaxR;}
448         
449   /*
450    * Gets the Eta cut value.
451    */
452   Double_t GetEtaCut() const{return fEtaCut;}
453         
454   /*
455    * Gets the Pt cut value.
456    */
457   Double_t GetPtCut() const{return fPtCut;}
458         
459         
460   /*
461    * Gets the MaxZCut value.
462    */
463   Double_t GetMaxZCut() const{return fMaxZ;}
464         
465         
466   /*
467    * Gets the MinClsTPC value.
468    */
469   Double_t GetMinClsTPCCut() const{return fMinClsTPC;}
470         
471
472   /*
473    * Gets the line cut values.
474    */
475   Double_t GetLineCutZRSlope() const{return fLineCutZRSlope;}
476   Double_t GetLineCutZValue() const{return fLineCutZValue;}
477         
478   /*
479    * Gets the Chi2 cut value for the conversions.
480    */
481   Double_t GetChi2CutConversion() const{return fChi2CutConversion;}
482         
483   /*
484    * Gets the Chi2 cut value for the mesons.
485    */
486   Double_t GetChi2CutMeson() const{return fChi2CutMeson;}
487         
488   /*
489    * Gets the alpha cut value for the mesons.
490    */
491   Double_t GetAlphaCutMeson() const{return fAlphaCutMeson;}
492
493
494   Double_t GetPositiveTrackLength() const{return fCurrentPositiveESDTrack->GetIntegratedLength();}
495   Double_t GetNegativeTrackLength() const{return fCurrentNegativeESDTrack->GetIntegratedLength();}
496         
497   Double_t GetPositiveNTPCClusters() const{return fCurrentPositiveESDTrack->GetTPCNcls();}
498   Double_t GetNegativeNTPCClusters() const{return fCurrentNegativeESDTrack->GetTPCNcls();}
499         
500   /*
501    * Sets the MaxVertexZ value.
502    */
503   void SetMaxVertexZ(Double_t maxVertexZ){fMaxVertexZ=maxVertexZ;}
504
505   /*
506    * Sets the MaxRCut value.
507    */
508   void SetMaxRCut(Double_t maxR){fMaxR=maxR;}
509         
510   /*
511    * Sets the EtaCut value.
512    */
513   void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}
514         
515   /*
516    * Sets the PtCut value.
517    */
518   void SetPtCut(Double_t ptCut){fPtCut=ptCut;}
519         
520   /*
521    * Sets the PtCut value.
522    */
523   void SetSinglePtCut(Double_t singleptCut){fSinglePtCut=singleptCut;}
524         
525     
526   /*
527    * Sets the MaxZCut value.
528    */
529   void SetMaxZCut(Double_t maxZ){fMaxZ=maxZ;}
530         
531  /*
532    * Sets the MinClsTPC value.
533    */
534   void SetMinClsTPCCut(Double_t minClsTPC){fMinClsTPC=minClsTPC;}
535         
536   /*
537    * Sets the LineCut values.
538    */
539   void SetLineCutZRSlope(Double_t LineCutZRSlope){fLineCutZRSlope=LineCutZRSlope;}
540   void SetLineCutZValue(Double_t LineCutZValue){fLineCutZValue=LineCutZValue;}
541         
542   /*
543    * Sets the Chi2Cut value for conversions.
544    */
545   void SetChi2CutConversion(Double_t chi2){fChi2CutConversion=chi2;}
546         
547   /*
548    * Sets the Chi2Cut for the mesons.
549    */
550   void SetChi2CutMeson(Double_t chi2){fChi2CutMeson=chi2;}
551         
552   /*
553    * Sets the AlphaCut for the mesons.
554    */
555   void SetAlphaCutMeson(Double_t alpha){fAlphaCutMeson=alpha;}
556         
557
558   /*
559    * Sets the XVertexCut value.
560    */
561   void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}
562         
563   /*
564    * Sets the YVertexCut value.
565    */
566   void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}
567         
568   /*
569    * Sets the ZVertexCut value.
570    */
571   void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}
572         
573   /*
574    * Sets the PIDProbabilityCut value for track particles.
575    */
576   void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}
577         
578   /*
579    * Sets the PIDProbability cut value for the negative track.
580    */
581   void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}
582         
583   /*
584    * Sets the PIDProbability cut value for the positive track.
585    */
586   void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}
587
588   /*
589    * Sets the PIDnSigmaAboveElectron cut value for the tracks.
590    */
591   void SetPIDnSigmaAboveElectronLine(Double_t nSigmaAbove){fPIDnSigmaAboveElectronLine=nSigmaAbove;}
592         
593   /*
594    * Sets the PIDnSigmaBelowElectron cut value for the tracks.
595    */
596   void SetPIDnSigmaBelowElectronLine(Double_t nSigmaBelow){fPIDnSigmaBelowElectronLine=nSigmaBelow;}
597         
598   /*
599    * Sets the PIDnSigmaAbovePion cut value for the tracks.
600    */
601   void SetPIDnSigmaAbovePionLine(Double_t nSigmaAbovePion){fPIDnSigmaAbovePionLine=nSigmaAbovePion;}
602
603   /*
604    * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
605    */
606   void SetPIDMinPnSigmaAbovePionLine(Double_t MinPnSigmaAbovePion){fPIDMinPnSigmaAbovePionLine=MinPnSigmaAbovePion;}
607
608  /*
609    * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
610    */
611   void SetPIDMaxPnSigmaAbovePionLine(Double_t MaxPnSigmaAbovePion){fPIDMaxPnSigmaAbovePionLine=MaxPnSigmaAbovePion;}
612
613   /*
614    * Sets the SigmaMassCut value.
615    */
616   void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}
617         
618   /*
619    * Sets the flag to enable/disable the usage of MC information. 
620    */
621   void SetDoMCTruth(Bool_t doMC){fDoMC = doMC;}
622
623   /*
624    * Sets the flag to enable/disable the usage of MC information. 
625    */
626   Bool_t GetDoMCTruth() const {return fDoMC;}
627         
628   /*
629    * Sets the flag to enable/disable the cut dedx N sigma 
630    */
631
632   void SetDodEdxSigmaCut( Bool_t dodEdxSigmaCut){fDodEdxSigmaCut=dodEdxSigmaCut;}
633
634   /*
635    * Sets the flag to enable/disable the cut dedx N sigma for Kaon Rejection at low p 
636    */
637   void SetDoKaonRejectionLowP( Bool_t doKaonRejectionLowP){fDoKaonRejectionLowP=doKaonRejectionLowP;}
638   /*
639    * Sets the flag to enable/disable the cut dedx N sigma for Proton Rejection at low p 
640    */
641   void SetDoProtonRejectionLowP( Bool_t doProtonRejectionLowP){fDoProtonRejectionLowP=doProtonRejectionLowP;}
642
643   /*
644    * Sets the flag to enable/disable the cut dedx N sigma for Pion Rejection at low p 
645    */
646   void SetDoPionRejectionLowP( Bool_t doPionRejectionLowP){fDoPionRejectionLowP=doPionRejectionLowP;}
647
648   /*
649    * Sets the PIDMinPnSigmaAroundKaon cut value for the tracks.
650    */
651   void SetPIDnSigmaAtLowPAroundKaonLine(Double_t nSigmaAtLowPAroundKaon){fPIDnSigmaAtLowPAroundKaonLine =nSigmaAtLowPAroundKaon;}
652
653   /*
654    * Sets the PIDMinPnSigmaAroundProton cut value for the tracks.
655    */
656   void SetPIDnSigmaAtLowPAroundProtonLine(Double_t nSigmaAtLowPAroundProton){fPIDnSigmaAtLowPAroundProtonLine =nSigmaAtLowPAroundProton;}
657
658   /*
659    * Sets the PIDMinPnSigmaAroundPion cut value for the tracks.
660    */
661   void SetPIDnSigmaAtLowPAroundPionLine(Double_t nSigmaAtLowPAroundPion){fPIDnSigmaAtLowPAroundPionLine =nSigmaAtLowPAroundPion;}
662
663   /*
664    * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
665    */
666   void SetPIDMinPKaonRejectionLowP(Double_t PIDMinPKaonRejectionLowP ){fPIDMinPKaonRejectionLowP=PIDMinPKaonRejectionLowP;}
667
668   /*
669    * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
670    */
671   void SetPIDMinPProtonRejectionLowP(Double_t PIDMinPProtonRejectionLowP ){fPIDMinPProtonRejectionLowP=PIDMinPProtonRejectionLowP;}
672   /*
673    * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
674    */
675   void SetPIDMinPPionRejectionLowP(Double_t PIDMinPPionRejectionLowP ){fPIDMinPPionRejectionLowP=PIDMinPPionRejectionLowP;}
676
677   /*
678    *Set if we want to use Gamma Selection based on Qt from Armenteros
679    */
680   void SetDoQtGammaSelection(Bool_t doQtGammaSelection){fDoQtGammaSelection=doQtGammaSelection;}
681   /*
682    * Sets the MaxQtCut value.
683    */
684   void SetQtMax(Double_t qtMax){fQtMax=qtMax;}
685
686   /*
687    * Updates the V0 information of the current V0.
688    */
689   Bool_t UpdateV0Information();
690         
691   /*
692    * Resets the V0 index.
693    */
694   void ResetV0IndexNumber(){fCurrentV0IndexNumber=0;}
695         
696   /*
697    * Sets the histograms.
698    */
699   void SetHistograms(AliGammaConversionHistograms * const histograms){fHistograms=histograms;}
700         
701   /*
702    * Check for primary vertex.
703    */
704   Bool_t CheckForPrimaryVertex();
705         
706   /*
707    * Check for primary vertex Z.
708    */
709   Bool_t CheckForPrimaryVertexZ();
710
711   /*
712    * Gets a vector of good v0s.
713    */
714   TClonesArray* GetCurrentEventGoodV0s() const{return fCurrentEventGoodV0s;}
715         
716   /*
717    * Gets the vector of previous events v0s (for bacground analysis)
718    */
719   AliGammaConversionKFVector* GetBGGoodV0s(Int_t event);
720   //  vector<AliKFParticle> GetPreviousEventGoodV0s() const{return fPreviousEventGoodV0s;}
721
722   void SetUseOwnXYZCalculation(Bool_t flag){fUseOwnXYZCalculation=flag;}
723
724   Bool_t GetHelixCenter(AliESDtrack* track, Double_t b,Int_t charge, Double_t center[2]);
725         
726   Bool_t GetConvPosXY(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b, Double_t convpos[2]);
727         
728   Double_t GetConvPosZ(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b);
729
730   Bool_t GetArmenterosQtAlfa(AliKFParticle * posKFparticle,AliKFParticle * negKFparticle,AliKFParticle * gamKFparticle,Double_t armenterosQtAlfa[2]);
731
732   void SetDoCF(Bool_t flag){fDoCF = flag;}
733
734   Bool_t CheckV0FinderStatus(Int_t index);
735
736   void SetOnFlyFlag(Bool_t flag){fUseOnFlyV0Finder = flag;}
737   
738   Int_t GetNBGEvents(){return fBGEventHandler->GetNBGEvents();}
739
740   void SetCalculateBackground(Bool_t flag){fCalculateBackground=flag;}
741
742   AliGammaConversionBGHandler* GetBGHandler() const {return fBGEventHandler;}
743
744   Double_t GetVertexZ(){return fESDEvent->GetPrimaryVertex()->GetZ();}
745
746   Int_t GetMultiplicity(){return CountESDTracks();}
747
748   void SetESDtrackCuts(AliESDtrackCuts * const trackCuts){fEsdTrackCuts = trackCuts;}
749
750   void SetNEventsForBG(Int_t nev){nEventsForBGCalculation = nev;}
751
752   Int_t CountESDTracks();
753
754   Int_t GetCurrentV0IndexNumber() const {return fCurrentV0IndexNumber;}
755
756   Bool_t CheckIfPi0IsMother(Int_t label);
757   Bool_t CheckIfEtaIsMother(Int_t label);
758
759   static void InitESDpid(Int_t type=0);
760   static void SetESDpid(AliESDpid * const pid) {fgESDpid=pid;}
761   static AliESDpid* GetESDpid() {return fgESDpid;}
762  
763
764
765  private:
766   AliStack * fMCStack;           // pointer to MonteCarlo particle stack 
767   //  AliMCEventHandler* fMCTruth;   // for CF    pointer to the MC object
768   AliMCEvent *fMCEvent;                 //  for CF      pointer to MC event
769   TChain * fChain;               // pointer to the TChain
770         
771   //  AliESDInputHandler* fESDHandler;      //! pointer to esd object
772   AliESDEvent *fESDEvent;               //! pointer to esd object
773         
774         
775   // for CF
776   AliCFManager *fCFManager; // pointer to the cf manager
777   //  AliCFContainer *container;
778         
779   // for dEdx cut based on nSigma to a particle line
780   //AliESDpid * fESDpid; // esd pid
781         
782   AliGammaConversionHistograms *fHistograms; // pointer to histogram handling class
783         
784   Int_t fCurrentV0IndexNumber;  // the current v0 index number
785   AliESDv0 * fCurrentV0;                //! pointer to the current v0
786   AliKFParticle * fCurrentNegativeKFParticle;  //! pointer to the negative KF particle
787   AliKFParticle * fCurrentPositiveKFParticle;  //! pointer to the positive KF particle
788   AliKFParticle * fCurrentMotherKFCandidate;   //! pointer to the positive KF particle
789         
790   AliESDtrack * fCurrentNegativeESDTrack;      //! pointer to the negative ESD track
791   AliESDtrack * fCurrentPositiveESDTrack;      //! pointer to the positive ESD track
792         
793   TLorentzVector * fNegativeTrackLorentzVector; //! pointer to the negative Track Lorentz Vector
794   TLorentzVector * fPositiveTrackLorentzVector; //! pointer to the positive Track Lorentz Vector
795   TLorentzVector * fMotherCandidateLorentzVector;   //! pointer to the mother candidate Track Lorentz Vector
796         
797   Double_t fCurrentXValue;   // current x value
798   Double_t fCurrentYValue;   // current y value
799   Double_t fCurrentZValue;   // current z value
800         
801   Int_t fPositiveTrackPID;   // positive track pid
802   Int_t fNegativeTrackPID;   // negative track pid
803         
804   TParticle *fNegativeMCParticle;      //!
805   TParticle *fPositiveMCParticle;      //!
806   TParticle *fMotherMCParticle;        //!
807         
808   Double_t fMotherCandidateKFMass;   // mass of mother candidate KF particle
809   Double_t fMotherCandidateKFWidth;  // width of mother candidate KF particle
810         
811   Bool_t fUseKFParticle;   // flag 
812   Bool_t fUseESDTrack;     // flag 
813   Bool_t fDoMC;            // flag 
814
815   //Event Cuts
816   Double_t fMaxVertexZ;
817   //cuts
818   Double_t fMaxR; //r cut
819   Double_t fEtaCut; //eta cut
820   Double_t fPtCut; // pt cut
821   Double_t fSinglePtCut; // pt cut for electron/positron
822   Double_t fMaxZ; //z cut
823   Double_t fMinClsTPC;
824   Double_t fLineCutZRSlope; //linecut
825   Double_t fLineCutZValue; //linecut
826   Double_t fChi2CutConversion; //chi2cut
827   Double_t fChi2CutMeson;  //chi2cut
828   Double_t fAlphaCutMeson;  //alphacut
829   Double_t fPIDProbabilityCutNegativeParticle; //pid cut
830   Double_t fPIDProbabilityCutPositiveParticle; //pid cut
831   Bool_t   fDodEdxSigmaCut; // flag to use the dEdxCut based on sigmas
832   Double_t fPIDnSigmaAboveElectronLine; // sigma cut
833   Double_t fPIDnSigmaBelowElectronLine; // sigma cut
834   Double_t fPIDnSigmaAbovePionLine;     // sigma cut
835   Double_t fPIDMinPnSigmaAbovePionLine; // sigma cut
836   Double_t fPIDMaxPnSigmaAbovePionLine; // sigma cut
837   Double_t fDoKaonRejectionLowP;   // Kaon rejection at low p
838   Double_t fDoProtonRejectionLowP; // Proton rejection at low p
839   Double_t fDoPionRejectionLowP;   // Pion rejection at low p
840   Double_t fPIDnSigmaAtLowPAroundKaonLine; // sigma cut
841   Double_t fPIDnSigmaAtLowPAroundProtonLine; // sigma cut
842   Double_t fPIDnSigmaAtLowPAroundPionLine; // sigma cut
843   Double_t fPIDMinPKaonRejectionLowP; // Momentum limit to apply kaon rejection
844   Double_t fPIDMinPProtonRejectionLowP; // Momentum limit to apply proton rejection
845   Double_t fPIDMinPPionRejectionLowP; // Momentum limit to apply proton rejection
846   Bool_t   fDoQtGammaSelection; // Select gammas using qtMax
847   Double_t fQtMax; // Maximum Qt from Armenteros to select Gammas
848   Double_t fXVertexCut; //vertex cut
849   Double_t fYVertexCut; //vertex cut
850   Double_t fZVertexCut; // vertexcut
851         
852   Double_t fNSigmaMass; //nsigma cut
853         
854   Bool_t fUseImprovedVertex; //flag
855
856   Bool_t fUseOwnXYZCalculation; //flag that determines if we use our own calculation of xyz (markus)
857
858   Bool_t fDoCF; //flag
859
860   Bool_t fUseOnFlyV0Finder; //flag
861
862   Bool_t fUpdateV0AlreadyCalled; //flag
863         
864   TClonesArray* fCurrentEventGoodV0s; //vector of good v0s
865   //  vector<AliKFParticle> fPreviousEventGoodV0s; // vector of good v0s from prevous events
866
867   Bool_t fCalculateBackground; //flag
868   AliGammaConversionBGHandler *fBGEventHandler; // background handler
869   Bool_t fBGEventInitialized; //flag
870         
871   AliESDtrackCuts *fEsdTrackCuts; // track cuts
872   Int_t fNumberOfESDTracks; //track counter
873
874   static AliESDpid* fgESDpid;                 // ESD pid object
875
876   Int_t nEventsForBGCalculation;
877
878   ClassDef(AliV0Reader,12)
879 };
880
881 inline void AliV0Reader::InitESDpid(Int_t type)
882 {
883   //
884   // initialize PID parameters
885   // type=0 is simulation
886   // type=1 is data
887
888   if (!fgESDpid) fgESDpid=new AliESDpid;
889   Double_t alephParameters[5];
890   // simulation
891   alephParameters[0] = 2.15898e+00/50.;
892   alephParameters[1] = 1.75295e+01;
893   alephParameters[2] = 3.40030e-09;
894   alephParameters[3] = 1.96178e+00;
895   alephParameters[4] = 3.91720e+00;
896   fgESDpid->GetTOFResponse().SetTimeResolution(80.);
897
898   // data
899   if (type==1){
900     alephParameters[0] = 0.0283086;
901     alephParameters[1] = 2.63394e+01;
902     alephParameters[2] = 5.04114e-11;
903     alephParameters[3] = 2.12543e+00;
904     alephParameters[4] = 4.88663e+00;
905     fgESDpid->GetTOFResponse().SetTimeResolution(130.);
906     fgESDpid->GetTPCResponse().SetMip(47.9);
907   }
908
909   fgESDpid->GetTPCResponse().SetBetheBlochParameters(
910     alephParameters[0],alephParameters[1],alephParameters[2],
911     alephParameters[3],alephParameters[4]);
912
913   fgESDpid->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
914
915 }
916
917 #endif
918
919
920
921
922