]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/GammaConv/AliV0Reader.h
Added new background scheeme, did some cleanup. added new bethe block parameters...
[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
24 class TClonesArray; 
25 class TFormula;
26 class Riostream;
27 class TChain;
28
29 //--- AliRoot system ---
30
31 class AliStack;
32 class AliMCEvent;       // for CF
33 class AliESDEvent; 
34 class AliMCEventHandler;
35 class AliESDInputHandler;
36 class AliESDVertex;
37 class AliLog;
38 class TChain;
39 class TChain;
40 class AliCFManager;   // for CF
41 class AliCFContainer;  // for CF
42 class AliTPCpidESD; // for dEdx cut based on nSigma to particle lines 
43
44
45 class AliV0Reader : public TObject {
46         
47  public: 
48         
49         
50   // for CF
51   enum{
52     kStepGenerated = 0,
53     kStepReconstructable = 1, 
54     kStepLikeSign = 2,
55     kStepTPCRefit = 3,
56     kStepKinks = 4,
57     kStepGetOnFly = 5,
58     kStepNContributors = 6,
59     kStepTPCPID = 7,
60     kStepR = 8,
61     kStepLine = 9,
62     kStepZ = 10,
63     kStepNDF = 11,
64     kStepChi2 = 12,
65     kStepEta = 13,
66     kStepPt = 14
67   };
68         
69         
70         
71   AliV0Reader();                                        //constructor
72   AliV0Reader(const AliV0Reader & g);                   //copy constructor
73   AliV0Reader & operator = (const AliV0Reader & g);     //assignment operator
74   //  virtual ~AliV0Reader() {;}                            //virtual destructor
75   virtual ~AliV0Reader();                            //virtual destructor
76   /*
77    *Initialize the reader
78    */
79   void Initialize();
80         
81         
82   // for CF
83   void SetCFManager(AliCFManager * const io){fCFManager = io;};
84   AliCFManager *GetCFManager() const {return fCFManager;}
85         
86         
87         
88         
89   /*
90    * Returns  AliESDEvent
91    */                   
92   AliESDEvent* GetESDEvent() const{return fESDEvent;}   
93         
94   /*
95    *Returns the number of v0s in the event, no cuts applied.
96    */
97   Int_t GetNumberOfV0s() const{return fESDEvent->GetNumberOfV0s();}
98         
99   /*
100    *Returns the number of contributors to the vertex
101    */
102   Int_t GetNumberOfContributorsVtx() const{return fESDEvent->GetPrimaryVertex()->GetNContributors();}
103   
104   /*
105    * Check if there are any more good v0s left in the v0 stack
106    * if so, fCurrent v0 is set to this v0 and can be retrieved
107    * by GetCurrentV0 function.
108    * returns kFALSE if there is no more good v0s in the v0 stack
109    */
110   Bool_t NextV0();
111         
112   /*
113    * Returns the v0 at the given index, no checks are done on the v0. 
114    */
115   AliESDv0* GetV0(Int_t index);
116         
117   /*
118    * Returns the current v0
119    */
120   AliESDv0* GetCurrentV0() const{return fCurrentV0;}
121         
122   /*
123    * Returns the negative ESD track which belongs to fCurrentV0
124    */
125   AliESDtrack* GetNegativeESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetNindex());}
126         
127   /*
128    * Returns the positive ESD track which belongs to fCurrentV0
129    */
130   AliESDtrack* GetPositiveESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetPindex());}
131         
132   /*
133    * Returns the negative KF particle which belongs to fCurrentV0
134    */
135   AliKFParticle* GetNegativeKFParticle() const{return fCurrentNegativeKFParticle;}
136         
137   /*
138    * Returns the positive KF particle which belongs to fCurrentV0
139    */
140   AliKFParticle* GetPositiveKFParticle() const{return fCurrentPositiveKFParticle;}
141         
142   /*
143    * Returns the KFParticle object of the 2 tracks.
144    */
145   AliKFParticle* GetMotherCandidateKFCombination() const{return fCurrentMotherKFCandidate;}
146         
147   /*
148    * Checks the probablity that the PID of the particle is what we want it to be.
149    */
150   Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);
151         
152   /*
153    * Checks if the PID of the two particles are within our cuts.
154    */
155   void GetPIDProbability(Double_t &negPIDProb, Double_t &posPIDProb);
156         
157   /*
158    *Get the negative MC TParticle from the stack 
159    */
160   TParticle * GetNegativeMCParticle() const{return fNegativeMCParticle;}
161         
162   /*
163    *Get the positive MC TParticle from the stack 
164    */
165   TParticle * GetPositiveMCParticle() const{return fPositiveMCParticle;}
166         
167   /*
168    *Get the mother MC TParticle from the stack 
169    */
170   TParticle * GetMotherMCParticle() const{return fMotherMCParticle;}
171         
172   /*
173    * Flag to see if the v0 particles share the same mother
174    */
175   Bool_t HasSameMCMother();
176         
177         
178   /*
179    *Get the PID of the MC mother particle
180    */
181   Int_t GetMotherMCParticlePDGCode() const{if(fMotherMCParticle != NULL){ cout<<"MCParticle exists"<<endl;} return fMotherMCParticle->GetPdgCode();}
182         
183   /*
184    *Get the MC stack 
185    */
186   AliStack* GetMCStack() const{return fMCStack;}
187         
188         
189   /*
190    * Setup  AliMCEventHandler
191    */                   
192   AliMCEventHandler* GetMCTruth() const{return fMCTruth;}       // for CF
193         
194         
195   /*
196    *Get the MC stack 
197    */
198   AliMCEvent* GetMCEvent() const{return fMCEvent;}   // for CF
199         
200         
201   /*
202    *Get the magnetic field from the ESD event 
203    */
204   Double_t GetMagneticField() const{return fESDEvent->GetMagneticField();}
205         
206   /*
207    *Get the primary vertex from the esd event
208    */
209   const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}
210         
211   /*
212    * Set the PID of the negative track
213    */
214   void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}
215         
216   /*
217    * Set the PID of the positive track
218    */
219   void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}
220         
221   /*
222    * Set the flag to use the kfparticle class. Will also disable the use of esd tracks
223    */
224   void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}
225         
226   /*
227    *  Set the flag to use the esd track class. Will also disable the use of kf particles
228    */
229   void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}
230         
231   /*
232    *  Set the flag to use improved vertex or not
233    */
234   void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}
235         
236   /*
237    * Return the number in the species array belonging to the negative or positive track pid.
238    */
239   Int_t GetSpeciesIndex(Int_t chargeOfTrack);
240         
241   /*
242    * Return the x coordinate of the v0
243    */
244   Double_t GetX() const{return fCurrentXValue;}
245         
246   /*
247    * Return the y coordinate of the v0
248    */
249   Double_t GetY() const{return fCurrentYValue;}
250         
251   /*
252    * Return the Z coordinate of the v0
253    */
254   Double_t GetZ() const{return fCurrentZValue;}
255         
256   /*
257    * Return the radius of the v0
258    */
259   Double_t GetXYRadius() const{return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}
260         
261   /*
262    * Get the opening angle between the two tracks
263    */
264   Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}
265
266   /*
267    * Get the Cos Pointing angle between the two tracks
268    */
269   Double_t GetCosPointingAngle(){return fCurrentV0->GetV0CosineOfPointingAngle();}
270
271   /*
272    * Get the DCA between the two tracks
273    */
274   Double_t GetDcaDaughters(){return fCurrentV0->GetDcaV0Daughters();}
275
276   /*
277    * Get the Normalized DCA between the two tracks
278    */
279   Double_t GetNormDcaDistDaughters(){return fCurrentV0->GetDcaV0Daughters()/fCurrentV0->GetDistSigma();}
280
281   /*
282    * Get the Likelihood for a Conversion
283    */
284   Double_t GetLikelihoodAP(){return fCurrentV0->GetLikelihoodAP(0,0);}
285       
286   /*
287    * Gets the Energy of the negative track.
288    */
289   Double_t GetNegativeTrackEnergy() const{return fCurrentNegativeKFParticle->E();}
290         
291   /*
292    * Gets the Energy of the positive track.
293    */
294   Double_t GetPositiveTrackEnergy() const{return fCurrentPositiveKFParticle->E();}
295         
296   /*
297    * Gets the Energy of the mother candidate.
298    */
299   Double_t GetMotherCandidateEnergy() const{return fCurrentMotherKFCandidate->E();}
300         
301   /*
302    * Gets the Pt of the negative track.
303    */
304   Double_t GetNegativeTrackPt() const{return fNegativeTrackLorentzVector->Pt();}
305         
306   /*
307    * Gets the Pt of the positive track.
308    */
309   Double_t GetPositiveTrackPt() const{return fPositiveTrackLorentzVector->Pt();}
310         
311
312   /*
313    * Gets the Pt of the mother candidate.
314    */
315   Double_t GetMotherCandidatePt() const{return fMotherCandidateLorentzVector->Pt();}
316
317
318   /*
319    * Gets the P of the mother candidate.
320    */
321   Double_t GetMotherCandidateP() const{return fMotherCandidateLorentzVector->P();}
322         
323
324   /*
325    * Gets the Eta of the negative track.
326    */
327   Double_t GetNegativeTrackEta() const{return fNegativeTrackLorentzVector->Eta();}
328   /*
329    * Gets the Eta of the positive track.
330    */
331   Double_t GetPositiveTrackEta() const{return fPositiveTrackLorentzVector->Eta();}
332   /*
333    * Gets the Eta of the mother candidate.
334    */
335   Double_t GetMotherCandidateEta() const{return fMotherCandidateLorentzVector->Eta();}
336         
337   /*
338    * Gets the NDF of the mother candidate.
339    */
340   Double_t GetMotherCandidateNDF() const{return fCurrentMotherKFCandidate->GetNDF();}
341         
342   /*
343    * Gets the Chi2 of the mother candidate.
344    */
345   Double_t GetMotherCandidateChi2() const{return fCurrentMotherKFCandidate->GetChi2();}
346         
347   /*
348    * Gets the Mass of the mother candidate.
349    */
350   Double_t GetMotherCandidateMass() const{return fMotherCandidateKFMass;}
351         
352   /*
353    * Gets the Width of the mother candidate.
354    */
355   Double_t GetMotherCandidateWidth() const{return fMotherCandidateKFWidth;}
356         
357   /*
358    * Gets the Phi of the negative track.
359    */
360   Double_t GetNegativeTrackPhi() const;
361         
362   /*
363    * Gets the Phi of the positive track.
364    */
365   Double_t GetPositiveTrackPhi() const;
366         
367   /*
368    * Gets the Phi of the mother candidate.
369    */
370   Double_t GetMotherCandidatePhi() const;
371         
372   /*
373    * Gets the Rapidity of the mother candidate.
374    */
375   Double_t GetMotherCandidateRapidity() const;
376         
377
378   /*
379    * Gets the P of the negative track.
380    */
381   Double_t GetNegativeTrackP() const{return fNegativeTrackLorentzVector->P();}
382         
383   /*
384    * Gets the P of the positive track.
385    */
386   Double_t GetPositiveTrackP() const{return fPositiveTrackLorentzVector->P();}
387
388   /*
389    * Gets the dE/dx in the TPC of the negative track.
390    */
391   Double_t GetNegativeTrackTPCdEdx() const{return fCurrentNegativeESDTrack->GetTPCsignal();}
392         
393   /*
394    * Gets the dE/dx in the TPC of the positive track.
395    */
396   Double_t GetPositiveTrackTPCdEdx() const{return fCurrentPositiveESDTrack->GetTPCsignal();}
397         
398   /*
399    * Update data which need to be updated every event.
400    */
401   void UpdateEventByEventData();
402         
403   /*
404    * Gets the MaxRCut value.
405    */
406   Double_t GetMaxRCut() const{return fMaxR;}
407         
408   /*
409    * Gets the Eta cut value.
410    */
411   Double_t GetEtaCut() const{return fEtaCut;}
412         
413   /*
414    * Gets the Pt cut value.
415    */
416   Double_t GetPtCut() const{return fPtCut;}
417         
418         
419   /*
420    * Gets the MaxZCut value.
421    */
422   Double_t GetMaxZCut() const{return fMaxZ;}
423         
424         
425   /*
426    * Gets the line cut values.
427    */
428   Double_t GetLineCutZRSlope() const{return fLineCutZRSlope;}
429   Double_t GetLineCutZValue() const{return fLineCutZValue;}
430         
431   /*
432    * Gets the Chi2 cut value for the conversions.
433    */
434   Double_t GetChi2CutConversion() const{return fChi2CutConversion;}
435         
436   /*
437    * Gets the Chi2 cut value for the mesons.
438    */
439   Double_t GetChi2CutMeson() const{return fChi2CutMeson;}
440         
441   Double_t GetPositiveTrackLength() const{return fCurrentPositiveESDTrack->GetIntegratedLength();}
442   Double_t GetNegativeTrackLength() const{return fCurrentNegativeESDTrack->GetIntegratedLength();}
443         
444   Double_t GetPositiveNTPCClusters() const{return fCurrentPositiveESDTrack->GetTPCNcls();}
445   Double_t GetNegativeNTPCClusters() const{return fCurrentNegativeESDTrack->GetTPCNcls();}
446         
447   /*
448    * Sets the MaxRCut value.
449    */
450   void SetMaxRCut(Double_t maxR){fMaxR=maxR;}
451         
452   /*
453    * Sets the EtaCut value.
454    */
455   void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}
456         
457   /*
458    * Sets the PtCut value.
459    */
460   void SetPtCut(Double_t ptCut){fPtCut=ptCut;}
461         
462     
463   /*
464    * Sets the MaxZCut value.
465    */
466   void SetMaxZCut(Double_t maxZ){fMaxZ=maxZ;}
467         
468         
469   /*
470    * Sets the LineCut values.
471    */
472   void SetLineCutZRSlope(Double_t LineCutZRSlope){fLineCutZRSlope=LineCutZRSlope;}
473   void SetLineCutZValue(Double_t LineCutZValue){fLineCutZValue=LineCutZValue;}
474         
475   /*
476    * Sets the Chi2Cut value for conversions.
477    */
478   void SetChi2CutConversion(Double_t chi2){fChi2CutConversion=chi2;}
479         
480   /*
481    * Sets the Chi2Cut for the mesons.
482    */
483   void SetChi2CutMeson(Double_t chi2){fChi2CutMeson=chi2;}
484         
485   /*
486    * Sets the XVertexCut value.
487    */
488   void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}
489         
490   /*
491    * Sets the YVertexCut value.
492    */
493   void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}
494         
495   /*
496    * Sets the ZVertexCut value.
497    */
498   void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}
499         
500   /*
501    * Sets the PIDProbabilityCut value for track particles.
502    */
503   void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}
504         
505   /*
506    * Sets the PIDProbability cut value for the negative track.
507    */
508   void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}
509         
510   /*
511    * Sets the PIDProbability cut value for the positive track.
512    */
513   void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}
514
515   /*
516    * Sets the PIDnSigmaAboveElectron cut value for the tracks.
517    */
518   void SetPIDnSigmaAboveElectronLine(Double_t nSigmaAbove){fPIDnSigmaAboveElectronLine=nSigmaAbove;}
519         
520   /*
521    * Sets the PIDnSigmaBelowElectron cut value for the tracks.
522    */
523   void SetPIDnSigmaBelowElectronLine(Double_t nSigmaBelow){fPIDnSigmaBelowElectronLine=nSigmaBelow;}
524         
525   /*
526    * Sets the PIDnSigmaAbovePion cut value for the tracks.
527    */
528   void SetPIDnSigmaAbovePionLine(Double_t nSigmaAbovePion){fPIDnSigmaAbovePionLine=nSigmaAbovePion;}
529
530   /*
531    * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
532    */
533   void SetPIDMinPnSigmaAbovePionLine(Double_t MinPnSigmaAbovePion){fPIDMinPnSigmaAbovePionLine=MinPnSigmaAbovePion;}
534
535   /*
536    * Sets the SigmaMassCut value.
537    */
538   void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}
539         
540   /*
541    * Sets the flag to enable/disable the usage of MC information. 
542    */
543   void SetDoMCTruth(Bool_t doMC){fDoMC = doMC;}
544
545   /*
546    * Sets the flag to enable/disable the usage of MC information. 
547    */
548   Bool_t GetDoMCTruth(){return fDoMC;}
549         
550   /*
551    * Sets the flag to enable/disable the cut dedx N sigma 
552    */
553
554   void SetDodEdxSigmaCut( Bool_t dodEdxSigmaCut){fDodEdxSigmaCut=dodEdxSigmaCut;}
555
556
557   /*
558    * Updates the V0 information of the current V0.
559    */
560   Bool_t UpdateV0Information();
561         
562   /*
563    * Resets the V0 index.
564    */
565   void ResetV0IndexNumber(){fCurrentV0IndexNumber=0;}
566         
567   /*
568    * Sets the histograms.
569    */
570   void SetHistograms(AliGammaConversionHistograms * const histograms){fHistograms=histograms;}
571         
572   /*
573    * Check for primary vertex.
574    */
575   Bool_t CheckForPrimaryVertex();
576         
577   /*
578    * Gets a vector of good v0s.
579    */
580   TClonesArray* GetCurrentEventGoodV0s() const{return fCurrentEventGoodV0s;}
581         
582   /*
583    * Gets the vector of previous events v0s (for bacground analysis)
584    */
585   AliGammaConversionKFVector* GetBGGoodV0s(Int_t event);
586   //  vector<AliKFParticle> GetPreviousEventGoodV0s() const{return fPreviousEventGoodV0s;}
587
588   void SetUseOwnXYZCalculation(Bool_t flag){fUseOwnXYZCalculation=flag;}
589
590   Bool_t GetHelixCenter(AliESDtrack* track, Double_t b,Int_t charge, Double_t center[2]);
591         
592   Bool_t GetConvPosXY(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b, Double_t convpos[2]);
593         
594   Double_t GetConvPosZ(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b);
595
596   void SetDoCF(Bool_t flag){fDoCF = flag;}
597
598   Bool_t CheckV0FinderStatus(Int_t index);
599
600   void SetOnFlyFlag(Bool_t flag){fUseOnFlyV0Finder = flag;}
601
602  Int_t GetNBGEvents(){return fBGEventHandler->GetNBGEvents();}
603
604  private:
605   AliStack * fMCStack;           // pointer to MonteCarlo particle stack 
606   AliMCEventHandler* fMCTruth;   // for CF    pointer to the MC object
607   AliMCEvent *fMCEvent;                 //  for CF      pointer to MC event
608   TChain * fChain;               // pointer to the TChain
609         
610   AliESDInputHandler* fESDHandler;      //! pointer to esd object
611   AliESDEvent *fESDEvent;               //! pointer to esd object
612         
613         
614   // for CF
615   AliCFManager *fCFManager; // pointer to the cf manager
616   //  AliCFContainer *container;
617         
618   // for dEdx cut based on nSigma to a particle line
619   AliTPCpidESD * fTPCpid; 
620         
621   AliGammaConversionHistograms *fHistograms; //! pointer to histogram handling class
622         
623   Int_t fCurrentV0IndexNumber;  // the current v0 index number
624   AliESDv0 * fCurrentV0;                //! pointer to the current v0
625   AliKFParticle * fCurrentNegativeKFParticle;  //! pointer to the negative KF particle
626   AliKFParticle * fCurrentPositiveKFParticle;  //! pointer to the positive KF particle
627   AliKFParticle * fCurrentMotherKFCandidate;   //! pointer to the positive KF particle
628         
629   AliESDtrack * fCurrentNegativeESDTrack;      //! pointer to the negative ESD track
630   AliESDtrack * fCurrentPositiveESDTrack;      //! pointer to the positive ESD track
631         
632   TLorentzVector * fNegativeTrackLorentzVector; //! pointer to the negative Track Lorentz Vector
633   TLorentzVector * fPositiveTrackLorentzVector; //! pointer to the positive Track Lorentz Vector
634   TLorentzVector * fMotherCandidateLorentzVector;   //! pointer to the mother candidate Track Lorentz Vector
635         
636   Double_t fCurrentXValue;   // current x value
637   Double_t fCurrentYValue;   // current y value
638   Double_t fCurrentZValue;   // current z value
639         
640   Int_t fPositiveTrackPID;   // positive track pid
641   Int_t fNegativeTrackPID;   // negative track pid
642         
643   TParticle *fNegativeMCParticle;      //!
644   TParticle *fPositiveMCParticle;      //!
645   TParticle *fMotherMCParticle;        //!
646         
647   Double_t fMotherCandidateKFMass;   // mass of mother candidate KF particle
648   Double_t fMotherCandidateKFWidth;  // width of mother candidate KF particle
649         
650   Bool_t fUseKFParticle;   // flag 
651   Bool_t fUseESDTrack;     // flag 
652   Bool_t fDoMC;            // flag 
653         
654   //cuts
655   Double_t fMaxR; //r cut
656   Double_t fEtaCut; //eta cut
657   Double_t fPtCut; // pt cut
658   Double_t fMaxZ; //z cut  
659   Double_t fLineCutZRSlope; //linecut
660   Double_t fLineCutZValue; //linecut
661   Double_t fChi2CutConversion; //chi2cut
662   Double_t fChi2CutMeson;  //chi2cut
663   Double_t fPIDProbabilityCutNegativeParticle; //pid cut
664   Double_t fPIDProbabilityCutPositiveParticle; //pid cut
665   Bool_t   fDodEdxSigmaCut; // flag to use the dEdxCut based on sigmas
666   Double_t fPIDnSigmaAboveElectronLine;
667   Double_t fPIDnSigmaBelowElectronLine;
668   Double_t fPIDnSigmaAbovePionLine;
669   Double_t fPIDMinPnSigmaAbovePionLine;
670   Double_t fXVertexCut; //vertex cut
671   Double_t fYVertexCut; //vertex cut
672   Double_t fZVertexCut; // vertexcut
673         
674   Double_t fNSigmaMass; //nsigma cut
675         
676   Bool_t fUseImprovedVertex; //flag
677
678   Bool_t fUseOwnXYZCalculation; //flag that determines if we use our own calculation of xyz (markus)
679
680   Bool_t fDoCF; //flag
681
682   Bool_t fUseOnFlyV0Finder;
683
684   Bool_t fUpdateV0AlreadyCalled; //flag
685         
686   TClonesArray* fCurrentEventGoodV0s; //vector of good v0s
687   //  vector<AliKFParticle> fPreviousEventGoodV0s; // vector of good v0s from prevous events
688
689   AliGammaConversionBGHandler *fBGEventHandler;
690   Bool_t fBGEventInitialized;
691         
692   ClassDef(AliV0Reader,8)
693 };
694 #endif
695
696
697