3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
6 ////////////////////////////////////////////////
7 //---------------------------------------------
8 // Class used to do analysis on conversion pairs
9 //---------------------------------------------
10 ////////////////////////////////////////////////
12 // --- ROOT system ---
14 #include "AliMCEvent.h" // for CF
16 #include "AliESDEvent.h"
17 #include "AliKFParticle.h"
18 #include "TParticle.h"
19 #include "AliGammaConversionHistograms.h"
21 #include "AliCFManager.h"
22 #include "AliGammaConversionBGHandler.h"
23 #include "AliESDpid.h"
30 //--- AliRoot system ---
33 class AliMCEvent; // for CF
35 class AliMCEventHandler;
36 class AliESDInputHandler;
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;
46 class AliV0Reader : public TObject {
54 kStepReconstructable = 1,
59 kStepdEdx_electronselection = 6,
60 kStepdEdx_pionrejection = 7,
61 kStepNContributors = 8,
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
81 *Initialize the reader
84 void SetInputAndMCEvent(AliVEvent* esd, AliMCEvent* mc) ;
87 virtual void SetInputEvent(AliVEvent* const input) {fESDEvent = dynamic_cast<AliESDEvent*>(input);}
88 virtual void SetMC(AliMCEvent* const mc) {fMCEvent = mc;}
92 void SetCFManager(AliCFManager * const io){fCFManager = io;};
93 AliCFManager *GetCFManager() const {return fCFManager;}
101 AliESDEvent* GetESDEvent() const{return fESDEvent;}
104 *Returns the number of v0s in the event, no cuts applied.
106 Int_t GetNumberOfV0s() const{return fESDEvent->GetNumberOfV0s();}
109 *Returns the number of contributors to the vertex
111 Int_t GetNumberOfContributorsVtx() const{return fESDEvent->GetPrimaryVertex()->GetNContributors();}
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
122 * Returns the v0 at the given index, no checks are done on the v0.
124 AliESDv0* GetV0(Int_t index);
127 * Returns the current v0
129 AliESDv0* GetCurrentV0() const{return fCurrentV0;}
132 * Returns the negative ESD track which belongs to fCurrentV0
134 // AliESDtrack* GetNegativeESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetNindex());}
135 AliESDtrack* GetNegativeESDTrack(){return fCurrentNegativeESDTrack;}
138 * Returns the positive ESD track which belongs to fCurrentV0
140 // AliESDtrack* GetPositiveESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetPindex());}
141 AliESDtrack* GetPositiveESDTrack(){return fCurrentPositiveESDTrack;}
144 * Returns the negative KF particle which belongs to fCurrentV0
146 AliKFParticle* GetNegativeKFParticle() const{return fCurrentNegativeKFParticle;}
149 * Returns the positive KF particle which belongs to fCurrentV0
151 AliKFParticle* GetPositiveKFParticle() const{return fCurrentPositiveKFParticle;}
154 * Returns the KFParticle object of the 2 tracks.
156 AliKFParticle* GetMotherCandidateKFCombination() const{return fCurrentMotherKFCandidate;}
159 * Checks the probablity that the PID of the particle is what we want it to be.
161 Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);
164 * Checks if the PID of the two particles are within our cuts.
166 void GetPIDProbability(Double_t &negPIDProb, Double_t &posPIDProb);
169 * Checks if the PID of the two particles are within our cuts.
171 void GetPIDProbabilityMuonPion(Double_t &negPIDProb, Double_t &posPIDProb);
174 *Get the negative MC TParticle from the stack
176 TParticle * GetNegativeMCParticle() const{return fNegativeMCParticle;}
179 *Get the positive MC TParticle from the stack
181 TParticle * GetPositiveMCParticle() const{return fPositiveMCParticle;}
184 *Get the mother MC TParticle from the stack
186 TParticle * GetMotherMCParticle() const{return fMotherMCParticle;}
189 * Flag to see if the v0 particles share the same mother
191 Bool_t HasSameMCMother();
195 *Get the PID of the MC mother particle
197 Int_t GetMotherMCParticlePDGCode() const{if(fMotherMCParticle != NULL){ cout<<"MCParticle exists"<<endl;} return fMotherMCParticle->GetPdgCode();}
202 AliStack* GetMCStack() const{return fMCStack;}
206 * Setup AliMCEventHandler
208 // AliMCEventHandler* GetMCTruth() const{return fMCTruth;} // for CF
214 AliMCEvent* GetMCEvent() const{return fMCEvent;} // for CF
218 *Get the magnetic field from the ESD event
220 Double_t GetMagneticField() const{return fESDEvent->GetMagneticField();}
223 *Get the primary vertex from the esd event
225 const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}
228 * Set the PID of the negative track
230 void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}
233 * Set the PID of the positive track
235 void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}
238 * Set the flag to use the kfparticle class. Will also disable the use of esd tracks
240 void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}
243 * Set the flag to use the esd track class. Will also disable the use of kf particles
245 void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}
248 * Set the flag to use improved vertex or not
250 void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}
253 * Return the number in the species array belonging to the negative or positive track pid.
255 Int_t GetSpeciesIndex(Int_t chargeOfTrack);
258 * Return the x coordinate of the v0
260 Double_t GetX() const{return fCurrentXValue;}
263 * Return the y coordinate of the v0
265 Double_t GetY() const{return fCurrentYValue;}
268 * Return the Z coordinate of the v0
270 Double_t GetZ() const{return fCurrentZValue;}
273 * Return the radius of the v0
275 Double_t GetXYRadius() const{return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}
278 * Get the opening angle between the two tracks
280 Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}
283 * Get the Cos Pointing angle between the two tracks
285 Double_t GetCosPointingAngle(){return fCurrentV0->GetV0CosineOfPointingAngle();}
288 * Get the DCA between the two tracks
290 Double_t GetDcaDaughters(){return fCurrentV0->GetDcaV0Daughters();}
293 * Get the Normalized DCA between the two tracks
295 Double_t GetNormDcaDistDaughters(){return fCurrentV0->GetDcaV0Daughters()/fCurrentV0->GetDistSigma();}
298 * Get the Likelihood for a Conversion
300 Double_t GetLikelihoodAP(){return fCurrentV0->GetLikelihoodAP(0,0);}
303 * Gets the Energy of the negative track.
305 Double_t GetNegativeTrackEnergy() const{return fCurrentNegativeKFParticle->E();}
308 * Gets the Energy of the positive track.
310 Double_t GetPositiveTrackEnergy() const{return fCurrentPositiveKFParticle->E();}
313 * Gets the Energy of the mother candidate.
315 Double_t GetMotherCandidateEnergy() const{return fCurrentMotherKFCandidate->E();}
318 * Gets the Pt of the negative track.
320 Double_t GetNegativeTrackPt() const{return fNegativeTrackLorentzVector->Pt();}
323 * Gets the Pt of the positive track.
325 Double_t GetPositiveTrackPt() const{return fPositiveTrackLorentzVector->Pt();}
329 * Gets the Pt of the mother candidate.
331 Double_t GetMotherCandidatePt() const{return fMotherCandidateLorentzVector->Pt();}
335 * Gets the P of the mother candidate.
337 Double_t GetMotherCandidateP() const{return fMotherCandidateLorentzVector->P();}
341 * Gets the Eta of the negative track.
343 Double_t GetNegativeTrackEta() const{return fNegativeTrackLorentzVector->Eta();}
345 * Gets the Eta of the positive track.
347 Double_t GetPositiveTrackEta() const{return fPositiveTrackLorentzVector->Eta();}
349 * Gets the Eta of the mother candidate.
351 Double_t GetMotherCandidateEta() const{return fMotherCandidateLorentzVector->Eta();}
354 * Gets the NDF of the mother candidate.
356 Double_t GetMotherCandidateNDF() const{return fCurrentMotherKFCandidate->GetNDF();}
359 * Gets the Chi2 of the mother candidate.
361 Double_t GetMotherCandidateChi2() const{return fCurrentMotherKFCandidate->GetChi2();}
364 * Gets the Mass of the mother candidate.
366 Double_t GetMotherCandidateMass() const{return fMotherCandidateKFMass;}
369 * Gets the Width of the mother candidate.
371 Double_t GetMotherCandidateWidth() const{return fMotherCandidateKFWidth;}
374 * Gets the Phi of the negative track.
376 Double_t GetNegativeTrackPhi() const;
379 * Gets the Phi of the positive track.
381 Double_t GetPositiveTrackPhi() const;
384 * Gets the Phi of the mother candidate.
386 Double_t GetMotherCandidatePhi() const;
389 * Gets the Rapidity of the mother candidate.
391 Double_t GetMotherCandidateRapidity() const;
395 * Gets the P of the negative track.
397 Double_t GetNegativeTrackP() const{return fNegativeTrackLorentzVector->P();}
400 * Gets the P of the positive track.
402 Double_t GetPositiveTrackP() const{return fPositiveTrackLorentzVector->P();}
405 * Gets the dE/dx in the TPC of the negative track.
407 Double_t GetNegativeTrackTPCdEdx() const{return fCurrentNegativeESDTrack->GetTPCsignal();}
410 * Gets the dE/dx in the TPC of the positive track.
412 Double_t GetPositiveTrackTPCdEdx() const{return fCurrentPositiveESDTrack->GetTPCsignal();}
415 * Gets the Number of the TPC clusters of the negative track.
417 Int_t GetNegativeTracknTPCClusters() const{return fCurrentNegativeESDTrack->GetNcls(1);}
420 * Gets the Number of the TPC clusters of the positive track.
422 Int_t GetPositiveTracknTPCClusters() const{return fCurrentPositiveESDTrack->GetNcls(1);}
425 * Gets the Number of the ITS clusters of the negative track.
427 Int_t GetNegativeTracknITSClusters() const{return fCurrentNegativeESDTrack->GetNcls(0);}
430 * Gets the Number of the ITS clusters of the positive track.
432 Int_t GetPositiveTracknITSClusters() const{return fCurrentPositiveESDTrack->GetNcls(0);}
435 * Update data which need to be updated every event.
437 void UpdateEventByEventData();
440 * Gets the MaxRCut value.
442 Double_t GetMaxVertexZ() const{return fMaxVertexZ;}
445 * Gets the MaxRCut value.
447 Double_t GetMaxRCut() const{return fMaxR;}
450 * Gets the Eta cut value.
452 Double_t GetEtaCut() const{return fEtaCut;}
455 * Gets the Pt cut value.
457 Double_t GetPtCut() const{return fPtCut;}
461 * Gets the MaxZCut value.
463 Double_t GetMaxZCut() const{return fMaxZ;}
467 * Gets the MinClsTPC value.
469 Double_t GetMinClsTPCCut() const{return fMinClsTPC;}
473 * Gets the line cut values.
475 Double_t GetLineCutZRSlope() const{return fLineCutZRSlope;}
476 Double_t GetLineCutZValue() const{return fLineCutZValue;}
479 * Gets the Chi2 cut value for the conversions.
481 Double_t GetChi2CutConversion() const{return fChi2CutConversion;}
484 * Gets the Chi2 cut value for the mesons.
486 Double_t GetChi2CutMeson() const{return fChi2CutMeson;}
489 * Gets the alpha cut value for the mesons.
491 Double_t GetAlphaCutMeson() const{return fAlphaCutMeson;}
494 Double_t GetPositiveTrackLength() const{return fCurrentPositiveESDTrack->GetIntegratedLength();}
495 Double_t GetNegativeTrackLength() const{return fCurrentNegativeESDTrack->GetIntegratedLength();}
497 Double_t GetPositiveNTPCClusters() const{return fCurrentPositiveESDTrack->GetTPCNcls();}
498 Double_t GetNegativeNTPCClusters() const{return fCurrentNegativeESDTrack->GetTPCNcls();}
501 * Sets the MaxVertexZ value.
503 void SetMaxVertexZ(Double_t maxVertexZ){fMaxVertexZ=maxVertexZ;}
506 * Sets the MaxRCut value.
508 void SetMaxRCut(Double_t maxR){fMaxR=maxR;}
511 * Sets the EtaCut value.
513 void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}
516 * Sets the PtCut value.
518 void SetPtCut(Double_t ptCut){fPtCut=ptCut;}
521 * Sets the PtCut value.
523 void SetSinglePtCut(Double_t singleptCut){fSinglePtCut=singleptCut;}
527 * Sets the MaxZCut value.
529 void SetMaxZCut(Double_t maxZ){fMaxZ=maxZ;}
532 * Sets the MinClsTPC value.
534 void SetMinClsTPCCut(Double_t minClsTPC){fMinClsTPC=minClsTPC;}
537 * Sets the LineCut values.
539 void SetLineCutZRSlope(Double_t LineCutZRSlope){fLineCutZRSlope=LineCutZRSlope;}
540 void SetLineCutZValue(Double_t LineCutZValue){fLineCutZValue=LineCutZValue;}
543 * Sets the Chi2Cut value for conversions.
545 void SetChi2CutConversion(Double_t chi2){fChi2CutConversion=chi2;}
548 * Sets the Chi2Cut for the mesons.
550 void SetChi2CutMeson(Double_t chi2){fChi2CutMeson=chi2;}
553 * Sets the AlphaCut for the mesons.
555 void SetAlphaCutMeson(Double_t alpha){fAlphaCutMeson=alpha;}
559 * Sets the XVertexCut value.
561 void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}
564 * Sets the YVertexCut value.
566 void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}
569 * Sets the ZVertexCut value.
571 void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}
574 * Sets the PIDProbabilityCut value for track particles.
576 void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}
579 * Sets the PIDProbability cut value for the negative track.
581 void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}
584 * Sets the PIDProbability cut value for the positive track.
586 void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}
589 * Sets the PIDnSigmaAboveElectron cut value for the tracks.
591 void SetPIDnSigmaAboveElectronLine(Double_t nSigmaAbove){fPIDnSigmaAboveElectronLine=nSigmaAbove;}
594 * Sets the PIDnSigmaBelowElectron cut value for the tracks.
596 void SetPIDnSigmaBelowElectronLine(Double_t nSigmaBelow){fPIDnSigmaBelowElectronLine=nSigmaBelow;}
599 * Sets the PIDnSigmaAbovePion cut value for the tracks.
601 void SetPIDnSigmaAbovePionLine(Double_t nSigmaAbovePion){fPIDnSigmaAbovePionLine=nSigmaAbovePion;}
604 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
606 void SetPIDMinPnSigmaAbovePionLine(Double_t MinPnSigmaAbovePion){fPIDMinPnSigmaAbovePionLine=MinPnSigmaAbovePion;}
609 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
611 void SetPIDMaxPnSigmaAbovePionLine(Double_t MaxPnSigmaAbovePion){fPIDMaxPnSigmaAbovePionLine=MaxPnSigmaAbovePion;}
614 * Sets the SigmaMassCut value.
616 void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}
619 * Sets the flag to enable/disable the usage of MC information.
621 void SetDoMCTruth(Bool_t doMC){fDoMC = doMC;}
624 * Sets the flag to enable/disable the usage of MC information.
626 Bool_t GetDoMCTruth() const {return fDoMC;}
629 * Sets the flag to enable/disable the cut dedx N sigma
632 void SetDodEdxSigmaCut( Bool_t dodEdxSigmaCut){fDodEdxSigmaCut=dodEdxSigmaCut;}
635 * Sets the flag to enable/disable the cut dedx N sigma for Kaon Rejection at low p
637 void SetDoKaonRejectionLowP( Bool_t doKaonRejectionLowP){fDoKaonRejectionLowP=doKaonRejectionLowP;}
639 * Sets the flag to enable/disable the cut dedx N sigma for Proton Rejection at low p
641 void SetDoProtonRejectionLowP( Bool_t doProtonRejectionLowP){fDoProtonRejectionLowP=doProtonRejectionLowP;}
644 * Sets the flag to enable/disable the cut dedx N sigma for Pion Rejection at low p
646 void SetDoPionRejectionLowP( Bool_t doPionRejectionLowP){fDoPionRejectionLowP=doPionRejectionLowP;}
649 * Sets the PIDMinPnSigmaAroundKaon cut value for the tracks.
651 void SetPIDnSigmaAtLowPAroundKaonLine(Double_t nSigmaAtLowPAroundKaon){fPIDnSigmaAtLowPAroundKaonLine =nSigmaAtLowPAroundKaon;}
654 * Sets the PIDMinPnSigmaAroundProton cut value for the tracks.
656 void SetPIDnSigmaAtLowPAroundProtonLine(Double_t nSigmaAtLowPAroundProton){fPIDnSigmaAtLowPAroundProtonLine =nSigmaAtLowPAroundProton;}
659 * Sets the PIDMinPnSigmaAroundPion cut value for the tracks.
661 void SetPIDnSigmaAtLowPAroundPionLine(Double_t nSigmaAtLowPAroundPion){fPIDnSigmaAtLowPAroundPionLine =nSigmaAtLowPAroundPion;}
664 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
666 void SetPIDMinPKaonRejectionLowP(Double_t PIDMinPKaonRejectionLowP ){fPIDMinPKaonRejectionLowP=PIDMinPKaonRejectionLowP;}
669 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
671 void SetPIDMinPProtonRejectionLowP(Double_t PIDMinPProtonRejectionLowP ){fPIDMinPProtonRejectionLowP=PIDMinPProtonRejectionLowP;}
673 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
675 void SetPIDMinPPionRejectionLowP(Double_t PIDMinPPionRejectionLowP ){fPIDMinPPionRejectionLowP=PIDMinPPionRejectionLowP;}
678 *Set if we want to use Gamma Selection based on Qt from Armenteros
680 void SetDoQtGammaSelection(Bool_t doQtGammaSelection){fDoQtGammaSelection=doQtGammaSelection;}
682 * Sets the MaxQtCut value.
684 void SetQtMax(Double_t qtMax){fQtMax=qtMax;}
687 * Updates the V0 information of the current V0.
689 Bool_t UpdateV0Information();
692 * Resets the V0 index.
694 void ResetV0IndexNumber(){fCurrentV0IndexNumber=0;}
697 * Sets the histograms.
699 void SetHistograms(AliGammaConversionHistograms * const histograms){fHistograms=histograms;}
702 * Check for primary vertex.
704 Bool_t CheckForPrimaryVertex();
707 * Check for primary vertex Z.
709 Bool_t CheckForPrimaryVertexZ();
712 * Gets a vector of good v0s.
714 TClonesArray* GetCurrentEventGoodV0s() const{return fCurrentEventGoodV0s;}
717 * Gets the vector of previous events v0s (for bacground analysis)
719 AliGammaConversionKFVector* GetBGGoodV0s(Int_t event);
720 // vector<AliKFParticle> GetPreviousEventGoodV0s() const{return fPreviousEventGoodV0s;}
722 void SetUseOwnXYZCalculation(Bool_t flag){fUseOwnXYZCalculation=flag;}
724 Bool_t GetHelixCenter(AliESDtrack* track, Double_t b,Int_t charge, Double_t center[2]);
726 Bool_t GetConvPosXY(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b, Double_t convpos[2]);
728 Double_t GetConvPosZ(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b);
730 Bool_t GetArmenterosQtAlfa(AliKFParticle * posKFparticle,AliKFParticle * negKFparticle,AliKFParticle * gamKFparticle,Double_t armenterosQtAlfa[2]);
732 void SetDoCF(Bool_t flag){fDoCF = flag;}
734 Bool_t CheckV0FinderStatus(Int_t index);
736 void SetOnFlyFlag(Bool_t flag){fUseOnFlyV0Finder = flag;}
738 Int_t GetNBGEvents(){return fBGEventHandler->GetNBGEvents();}
740 void SetCalculateBackground(Bool_t flag){fCalculateBackground=flag;}
742 AliGammaConversionBGHandler* GetBGHandler() const {return fBGEventHandler;}
744 Double_t GetVertexZ(){return fESDEvent->GetPrimaryVertex()->GetZ();}
746 Int_t GetMultiplicity(){return CountESDTracks();}
748 void SetESDtrackCuts(AliESDtrackCuts * const trackCuts){fEsdTrackCuts = trackCuts;}
750 void SetNEventsForBG(Int_t nev){nEventsForBGCalculation = nev;}
752 Int_t CountESDTracks();
754 Int_t GetCurrentV0IndexNumber() const {return fCurrentV0IndexNumber;}
756 Bool_t CheckIfPi0IsMother(Int_t label);
757 Bool_t CheckIfEtaIsMother(Int_t label);
759 static void InitESDpid(Int_t type=0);
760 static void SetESDpid(AliESDpid * const pid) {fgESDpid=pid;}
761 static AliESDpid* GetESDpid() {return fgESDpid;}
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
771 // AliESDInputHandler* fESDHandler; //! pointer to esd object
772 AliESDEvent *fESDEvent; //! pointer to esd object
776 AliCFManager *fCFManager; // pointer to the cf manager
777 // AliCFContainer *container;
779 // for dEdx cut based on nSigma to a particle line
780 //AliESDpid * fESDpid; // esd pid
782 AliGammaConversionHistograms *fHistograms; // pointer to histogram handling class
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
790 AliESDtrack * fCurrentNegativeESDTrack; //! pointer to the negative ESD track
791 AliESDtrack * fCurrentPositiveESDTrack; //! pointer to the positive ESD track
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
797 Double_t fCurrentXValue; // current x value
798 Double_t fCurrentYValue; // current y value
799 Double_t fCurrentZValue; // current z value
801 Int_t fPositiveTrackPID; // positive track pid
802 Int_t fNegativeTrackPID; // negative track pid
804 TParticle *fNegativeMCParticle; //!
805 TParticle *fPositiveMCParticle; //!
806 TParticle *fMotherMCParticle; //!
808 Double_t fMotherCandidateKFMass; // mass of mother candidate KF particle
809 Double_t fMotherCandidateKFWidth; // width of mother candidate KF particle
811 Bool_t fUseKFParticle; // flag
812 Bool_t fUseESDTrack; // flag
813 Bool_t fDoMC; // flag
816 Double_t fMaxVertexZ;
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
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
852 Double_t fNSigmaMass; //nsigma cut
854 Bool_t fUseImprovedVertex; //flag
856 Bool_t fUseOwnXYZCalculation; //flag that determines if we use our own calculation of xyz (markus)
860 Bool_t fUseOnFlyV0Finder; //flag
862 Bool_t fUpdateV0AlreadyCalled; //flag
864 TClonesArray* fCurrentEventGoodV0s; //vector of good v0s
865 // vector<AliKFParticle> fPreviousEventGoodV0s; // vector of good v0s from prevous events
867 Bool_t fCalculateBackground; //flag
868 AliGammaConversionBGHandler *fBGEventHandler; // background handler
869 Bool_t fBGEventInitialized; //flag
871 AliESDtrackCuts *fEsdTrackCuts; // track cuts
872 Int_t fNumberOfESDTracks; //track counter
874 static AliESDpid* fgESDpid; // ESD pid object
876 Int_t nEventsForBGCalculation;
878 ClassDef(AliV0Reader,12)
881 inline void AliV0Reader::InitESDpid(Int_t type)
884 // initialize PID parameters
885 // type=0 is simulation
888 if (!fgESDpid) fgESDpid=new AliESDpid;
889 Double_t alephParameters[5];
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.);
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);
909 fgESDpid->GetTPCResponse().SetBetheBlochParameters(
910 alephParameters[0],alephParameters[1],alephParameters[2],
911 alephParameters[3],alephParameters[4]);
913 fgESDpid->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);