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();}
112 Int_t GetNumberOfContributorsVtx();
115 * Check if there are any more good v0s left in the v0 stack
116 * if so, fCurrent v0 is set to this v0 and can be retrieved
117 * by GetCurrentV0 function.
118 * returns kFALSE if there is no more good v0s in the v0 stack
123 * Returns the v0 at the given index, no checks are done on the v0.
125 AliESDv0* GetV0(Int_t index);
128 * Returns the current v0
130 AliESDv0* GetCurrentV0() const{return fCurrentV0;}
133 * Returns the negative ESD track which belongs to fCurrentV0
135 // AliESDtrack* GetNegativeESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetNindex());}
136 AliESDtrack* GetNegativeESDTrack(){return fCurrentNegativeESDTrack;}
139 * Returns the positive ESD track which belongs to fCurrentV0
141 // AliESDtrack* GetPositiveESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetPindex());}
142 AliESDtrack* GetPositiveESDTrack(){return fCurrentPositiveESDTrack;}
145 * Returns the negative KF particle which belongs to fCurrentV0
147 AliKFParticle* GetNegativeKFParticle() const{return fCurrentNegativeKFParticle;}
150 * Returns the positive KF particle which belongs to fCurrentV0
152 AliKFParticle* GetPositiveKFParticle() const{return fCurrentPositiveKFParticle;}
155 * Returns the KFParticle object of the 2 tracks.
157 AliKFParticle* GetMotherCandidateKFCombination() const{return fCurrentMotherKFCandidate;}
160 * Checks the probablity that the PID of the particle is what we want it to be.
162 Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);
165 * Checks if the PID of the two particles are within our cuts.
167 void GetPIDProbability(Double_t &negPIDProb, Double_t &posPIDProb);
170 * Checks if the PID of the two particles are within our cuts.
172 void GetPIDProbabilityMuonPion(Double_t &negPIDProb, Double_t &posPIDProb);
175 *Get the negative MC TParticle from the stack
177 TParticle * GetNegativeMCParticle() const{return fNegativeMCParticle;}
180 *Get the positive MC TParticle from the stack
182 TParticle * GetPositiveMCParticle() const{return fPositiveMCParticle;}
185 *Get the mother MC TParticle from the stack
187 TParticle * GetMotherMCParticle() const{return fMotherMCParticle;}
190 * Flag to see if the v0 particles share the same mother
192 Bool_t HasSameMCMother();
196 *Get the PID of the MC mother particle
198 Int_t GetMotherMCParticlePDGCode() const{if(fMotherMCParticle != NULL){ cout<<"MCParticle exists"<<endl;} return fMotherMCParticle->GetPdgCode();}
203 AliStack* GetMCStack() const{return fMCStack;}
207 * Setup AliMCEventHandler
209 // AliMCEventHandler* GetMCTruth() const{return fMCTruth;} // for CF
215 AliMCEvent* GetMCEvent() const{return fMCEvent;} // for CF
219 *Get the magnetic field from the ESD event
221 Double_t GetMagneticField() const{return fESDEvent->GetMagneticField();}
224 *Get the primary vertex from the esd event
226 const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}
229 * Set the PID of the negative track
231 void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}
234 * Set the PID of the positive track
236 void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}
239 * Set the flag to use the kfparticle class. Will also disable the use of esd tracks
241 void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}
244 * Set the flag to use the esd track class. Will also disable the use of kf particles
246 void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}
249 * Set the flag to use improved vertex or not
251 void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}
254 * Return the number in the species array belonging to the negative or positive track pid.
256 Int_t GetSpeciesIndex(Int_t chargeOfTrack);
259 * Return the x coordinate of the v0
261 Double_t GetX() const{return fCurrentXValue;}
264 * Return the y coordinate of the v0
266 Double_t GetY() const{return fCurrentYValue;}
269 * Return the Z coordinate of the v0
271 Double_t GetZ() const{return fCurrentZValue;}
274 * Return the radius of the v0
276 Double_t GetXYRadius() const{return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}
279 * Get the opening angle between the two tracks
281 Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}
284 * Get the Cos Pointing angle between the two tracks
286 Double_t GetCosPointingAngle(){return fCurrentV0->GetV0CosineOfPointingAngle();}
289 * Get the DCA between the two tracks
291 Double_t GetDcaDaughters(){return fCurrentV0->GetDcaV0Daughters();}
294 * Get the Normalized DCA between the two tracks
296 Double_t GetNormDcaDistDaughters(){return fCurrentV0->GetDcaV0Daughters()/fCurrentV0->GetDistSigma();}
299 * Get the Likelihood for a Conversion
301 Double_t GetLikelihoodAP(){return fCurrentV0->GetLikelihoodAP(0,0);}
304 * Gets the Energy of the negative track.
306 Double_t GetNegativeTrackEnergy() const{return fCurrentNegativeKFParticle->E();}
309 * Gets the Energy of the positive track.
311 Double_t GetPositiveTrackEnergy() const{return fCurrentPositiveKFParticle->E();}
314 * Gets the Energy of the mother candidate.
316 Double_t GetMotherCandidateEnergy() const{return fCurrentMotherKFCandidate->E();}
319 * Gets the Pt of the negative track.
321 Double_t GetNegativeTrackPt() const{return fNegativeTrackLorentzVector->Pt();}
324 * Gets the Pt of the positive track.
326 Double_t GetPositiveTrackPt() const{return fPositiveTrackLorentzVector->Pt();}
330 * Gets the Pt of the mother candidate.
332 Double_t GetMotherCandidatePt() const{return fMotherCandidateLorentzVector->Pt();}
336 * Gets the P of the mother candidate.
338 Double_t GetMotherCandidateP() const{return fMotherCandidateLorentzVector->P();}
342 * Gets the Eta of the negative track.
344 Double_t GetNegativeTrackEta() const{return fNegativeTrackLorentzVector->Eta();}
346 * Gets the Eta of the positive track.
348 Double_t GetPositiveTrackEta() const{return fPositiveTrackLorentzVector->Eta();}
350 * Gets the Eta of the mother candidate.
352 Double_t GetMotherCandidateEta() const{return fMotherCandidateLorentzVector->Eta();}
355 * Gets the NDF of the mother candidate.
357 Double_t GetMotherCandidateNDF() const{return fCurrentMotherKFCandidate->GetNDF();}
360 * Gets the Chi2 of the mother candidate.
362 Double_t GetMotherCandidateChi2() const{return fCurrentMotherKFCandidate->GetChi2();}
365 * Gets the Mass of the mother candidate.
367 Double_t GetMotherCandidateMass() const{return fMotherCandidateKFMass;}
370 * Gets the Width of the mother candidate.
372 Double_t GetMotherCandidateWidth() const{return fMotherCandidateKFWidth;}
375 * Gets the Phi of the negative track.
377 Double_t GetNegativeTrackPhi() const;
380 * Gets the Phi of the positive track.
382 Double_t GetPositiveTrackPhi() const;
385 * Gets the Phi of the mother candidate.
387 Double_t GetMotherCandidatePhi() const;
390 * Gets the Rapidity of the mother candidate.
392 Double_t GetMotherCandidateRapidity() const;
396 * Gets the P of the negative track.
398 Double_t GetNegativeTrackP() const{return fNegativeTrackLorentzVector->P();}
401 * Gets the P of the positive track.
403 Double_t GetPositiveTrackP() const{return fPositiveTrackLorentzVector->P();}
406 * Gets the dE/dx in the TPC of the negative track.
408 Double_t GetNegativeTrackTPCdEdx() const{return fCurrentNegativeESDTrack->GetTPCsignal();}
411 * Gets the dE/dx in the TPC of the positive track.
413 Double_t GetPositiveTrackTPCdEdx() const{return fCurrentPositiveESDTrack->GetTPCsignal();}
416 * Gets the Number of the TPC clusters of the negative track.
418 Int_t GetNegativeTracknTPCClusters() const{return fCurrentNegativeESDTrack->GetNcls(1);}
421 * Gets the Number of the TPC clusters of the positive track.
423 Int_t GetPositiveTracknTPCClusters() const{return fCurrentPositiveESDTrack->GetNcls(1);}
427 * Gets the Number of the TPC findable clusters of the negative track.
429 Int_t GetNegativeTracknTPCFClusters() const{return fCurrentNegativeESDTrack->GetTPCNclsF();}
432 * Gets the Number of the TPC findable clusters of the positive track.
434 Int_t GetPositiveTracknTPCFClusters() const{return fCurrentPositiveESDTrack->GetTPCNclsF();}
437 * Gets the Number of the ITS clusters of the negative track.
439 Int_t GetNegativeTracknITSClusters() const{return fCurrentNegativeESDTrack->GetNcls(0);}
442 * Gets the Number of the ITS clusters of the positive track.
444 Int_t GetPositiveTracknITSClusters() const{return fCurrentPositiveESDTrack->GetNcls(0);}
447 * Gets the chi2 of the TPC negative track.
449 Double_t GetNegativeTrackTPCchi2() const{return fCurrentNegativeESDTrack->GetTPCchi2();}
452 * Gets the chi2 of the TPC the positive track.
454 Double_t GetPositiveTrackTPCchi2() const{return fCurrentPositiveESDTrack->GetTPCchi2();}
457 * Update data which need to be updated every event.
459 void UpdateEventByEventData();
462 * Gets the MaxRCut value.
464 Double_t GetMaxVertexZ() const{return fMaxVertexZ;}
467 * Gets the MaxRCut value.
469 Double_t GetMaxRCut() const{return fMaxR;}
472 * Gets the MinRCut value.
474 Double_t GetMinRCut() const{return fMinR;}
477 * Gets the Eta cut value.
479 Double_t GetEtaCut() const{return fEtaCut;}
482 * Gets the Rapidity Meson cut value.
484 Double_t GetRapidityMesonCut() const{return fRapidityMesonCut;}
487 * Gets the Pt cut value.
489 Double_t GetPtCut() const{return fPtCut;}
493 * Gets the MaxZCut value.
495 Double_t GetMaxZCut() const{return fMaxZ;}
499 * Gets the MinClsTPC value.
501 Double_t GetMinClsTPCCut() const{return fMinClsTPC;}
504 * Gets the MinClsTPC value.
506 Double_t GetMinClsTPCCutToF() const{return fMinClsTPCToF;}
511 * Gets the line cut values.
513 Double_t GetLineCutZRSlope() const{return fLineCutZRSlope;}
514 Double_t GetLineCutZValue() const{return fLineCutZValue;}
517 * Gets the Chi2 cut value for the conversions.
519 Double_t GetChi2CutConversion() const{return fChi2CutConversion;}
522 * Gets the Chi2 cut value for the mesons.
524 Double_t GetChi2CutMeson() const{return fChi2CutMeson;}
527 * Gets the alpha cut value for the mesons.
529 Double_t GetAlphaCutMeson() const{return fAlphaCutMeson;}
532 * Gets the Minimum alpha cut value for the mesons.
534 Double_t GetAlphaMinCutMeson() const{return fAlphaMinCutMeson;}
536 Double_t GetPositiveTrackLength() const{return fCurrentPositiveESDTrack->GetIntegratedLength();}
537 Double_t GetNegativeTrackLength() const{return fCurrentNegativeESDTrack->GetIntegratedLength();}
539 Double_t GetPositiveNTPCClusters() const{return fCurrentPositiveESDTrack->GetTPCNcls();}
540 Double_t GetNegativeNTPCClusters() const{return fCurrentNegativeESDTrack->GetTPCNcls();}
543 * Sets the MaxVertexZ value.
545 void SetMaxVertexZ(Double_t maxVertexZ){fMaxVertexZ=maxVertexZ;}
548 * Sets the MaxRCut value.
550 void SetMaxRCut(Double_t maxR){fMaxR=maxR;}
552 * Sets the MinRCut value.
554 void SetMinRCut(Double_t minR){fMinR=minR;}
557 * Sets the EtaCut value.
559 void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}
562 * Sets the Rapidity Meson Cut value.
564 void SetRapidityMesonCut(Double_t RapidityMesonCut){fRapidityMesonCut=RapidityMesonCut;}
567 * Sets the PtCut value.
569 void SetPtCut(Double_t ptCut){fPtCut=ptCut;}
572 * Sets the PtCut value.
574 void SetSinglePtCut(Double_t singleptCut){fSinglePtCut=singleptCut;}
578 * Sets the MaxZCut value.
580 void SetMaxZCut(Double_t maxZ){fMaxZ=maxZ;}
583 * Sets the MinClsTPC value.
585 void SetMinClsTPCCut(Double_t minClsTPC){fMinClsTPC=minClsTPC;}
588 * Sets the MinClsTPC value.
590 void SetMinClsTPCCutToF(Double_t minClsTPCToF){fMinClsTPCToF=minClsTPCToF;}
594 * Sets the LineCut values.
596 void SetLineCutZRSlope(Double_t LineCutZRSlope){fLineCutZRSlope=LineCutZRSlope;}
597 void SetLineCutZValue(Double_t LineCutZValue){fLineCutZValue=LineCutZValue;}
600 * Sets the Chi2Cut value for conversions.
602 void SetChi2CutConversion(Double_t chi2){fChi2CutConversion=chi2;}
605 * Sets the Chi2Cut for the mesons.
607 void SetChi2CutMeson(Double_t chi2){fChi2CutMeson=chi2;}
610 * Sets the AlphaCut for the mesons.
612 void SetAlphaCutMeson(Double_t alpha){fAlphaCutMeson=alpha;}
616 * Sets the AlphaCut for the mesons.
618 void SetAlphaMinCutMeson(Double_t alpha){fAlphaMinCutMeson=alpha;}
622 * Sets the XVertexCut value.
624 void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}
627 * Sets the YVertexCut value.
629 void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}
632 * Sets the ZVertexCut value.
634 void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}
637 * Sets the PIDProbabilityCut value for track particles.
639 void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}
642 * Sets the PIDProbability cut value for the negative track.
644 void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}
647 * Sets the PIDProbability cut value for the positive track.
649 void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}
652 * Sets the PIDnSigmaAboveElectron cut value for the tracks.
654 void SetPIDnSigmaAboveElectronLine(Double_t nSigmaAbove){fPIDnSigmaAboveElectronLine=nSigmaAbove;}
657 * Sets the PIDnSigmaBelowElectron cut value for the tracks.
659 void SetPIDnSigmaBelowElectronLine(Double_t nSigmaBelow){fPIDnSigmaBelowElectronLine=nSigmaBelow;}
662 * Sets the PIDnSigmaAbovePion cut value for the tracks.
664 void SetPIDnSigmaAbovePionLine(Double_t nSigmaAbovePion){fPIDnSigmaAbovePionLine=nSigmaAbovePion;}
667 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
669 void SetPIDMinPnSigmaAbovePionLine(Double_t MinPnSigmaAbovePion){fPIDMinPnSigmaAbovePionLine=MinPnSigmaAbovePion;}
672 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
674 void SetPIDMaxPnSigmaAbovePionLine(Double_t MaxPnSigmaAbovePion){fPIDMaxPnSigmaAbovePionLine=MaxPnSigmaAbovePion;}
677 * Sets the SigmaMassCut value.
679 void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}
682 * Sets the flag to enable/disable the usage of MC information.
684 void SetDoMCTruth(Bool_t doMC){fDoMC = doMC;}
687 * Sets the flag to enable/disable the usage of MC information.
689 Bool_t GetDoMCTruth() const {return fDoMC;}
692 * Sets the flag to enable/disable the cut dedx N sigma
695 void SetDodEdxSigmaCut( Bool_t dodEdxSigmaCut){fDodEdxSigmaCut=dodEdxSigmaCut;}
698 * Sets the flag to enable/disable the cut dedx N sigma for Kaon Rejection at low p
700 void SetDoKaonRejectionLowP( Bool_t doKaonRejectionLowP){fDoKaonRejectionLowP=doKaonRejectionLowP;}
702 * Sets the flag to enable/disable the cut dedx N sigma for Proton Rejection at low p
704 void SetDoProtonRejectionLowP( Bool_t doProtonRejectionLowP){fDoProtonRejectionLowP=doProtonRejectionLowP;}
707 * Sets the flag to enable/disable the cut dedx N sigma for Pion Rejection at low p
709 void SetDoPionRejectionLowP( Bool_t doPionRejectionLowP){fDoPionRejectionLowP=doPionRejectionLowP;}
712 * Sets the PIDMinPnSigmaAroundKaon cut value for the tracks.
714 void SetPIDnSigmaAtLowPAroundKaonLine(Double_t nSigmaAtLowPAroundKaon){fPIDnSigmaAtLowPAroundKaonLine =nSigmaAtLowPAroundKaon;}
717 * Sets the PIDMinPnSigmaAroundProton cut value for the tracks.
719 void SetPIDnSigmaAtLowPAroundProtonLine(Double_t nSigmaAtLowPAroundProton){fPIDnSigmaAtLowPAroundProtonLine =nSigmaAtLowPAroundProton;}
722 * Sets the PIDMinPnSigmaAroundPion cut value for the tracks.
724 void SetPIDnSigmaAtLowPAroundPionLine(Double_t nSigmaAtLowPAroundPion){fPIDnSigmaAtLowPAroundPionLine =nSigmaAtLowPAroundPion;}
727 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
729 void SetPIDMinPKaonRejectionLowP(Double_t PIDMinPKaonRejectionLowP ){fPIDMinPKaonRejectionLowP=PIDMinPKaonRejectionLowP;}
732 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
734 void SetPIDMinPProtonRejectionLowP(Double_t PIDMinPProtonRejectionLowP ){fPIDMinPProtonRejectionLowP=PIDMinPProtonRejectionLowP;}
736 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
738 void SetPIDMinPPionRejectionLowP(Double_t PIDMinPPionRejectionLowP ){fPIDMinPPionRejectionLowP=PIDMinPPionRejectionLowP;}
741 *Set if we want to use Gamma Selection based on Qt from Armenteros
743 void SetDoQtGammaSelection(Bool_t doQtGammaSelection){fDoQtGammaSelection=doQtGammaSelection;}
745 * Sets the MaxQtCut value.
747 void SetQtMax(Double_t qtMax){fQtMax=qtMax;}
750 * Updates the V0 information of the current V0.
752 Bool_t UpdateV0Information();
755 * Resets the V0 index.
757 void ResetV0IndexNumber(){fCurrentV0IndexNumber=0;}
761 * Returns number of good v0s in the event
763 Int_t GetNGoodV0s() const {return fNumberOfGoodV0s;}
766 * Sets the histograms.
768 void SetHistograms(AliGammaConversionHistograms * const histograms){fHistograms=histograms;}
771 * Check for primary vertex.
773 Bool_t CheckForPrimaryVertex();
776 * Check for primary vertex Z.
778 Bool_t CheckForPrimaryVertexZ();
781 * Gets a vector of good v0s.
783 TClonesArray* GetCurrentEventGoodV0s() const{return fCurrentEventGoodV0s;}
786 * Gets the vector of previous events v0s (for bacground analysis)
788 AliGammaConversionKFVector* GetBGGoodV0s(Int_t event);
789 // vector<AliKFParticle> GetPreviousEventGoodV0s() const{return fPreviousEventGoodV0s;}
791 void SetUseOwnXYZCalculation(Bool_t flag){fUseOwnXYZCalculation=flag;}
793 void SetUseConstructGamma(Bool_t flag){fUseConstructGamma=flag;}
795 Bool_t GetHelixCenter(AliESDtrack* track, Double_t b,Int_t charge, Double_t center[2]);
797 Bool_t GetConvPosXY(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b, Double_t convpos[2]);
799 Double_t GetConvPosZ(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b);
801 Bool_t GetArmenterosQtAlfa(AliKFParticle * posKFparticle,AliKFParticle * negKFparticle,AliKFParticle * gamKFparticle,Double_t armenterosQtAlfa[2]);
803 void SetDoCF(Bool_t flag){fDoCF = flag;}
805 Bool_t CheckV0FinderStatus(Int_t index);
807 void SetOnFlyFlag(Bool_t flag){fUseOnFlyV0Finder = flag;}
809 Int_t GetNBGEvents(){return fBGEventHandler->GetNBGEvents();}
811 void SetCalculateBackground(Bool_t flag){fCalculateBackground=flag;}
813 AliGammaConversionBGHandler* GetBGHandler() const {return fBGEventHandler;}
815 Double_t GetVertexZ(){return fESDEvent->GetPrimaryVertex()->GetZ();}
817 Int_t GetMultiplicity(){return CountESDTracks();}
819 void SetESDtrackCuts(AliESDtrackCuts * const trackCuts){fEsdTrackCuts = trackCuts;}
821 void SetNEventsForBG(Int_t nev){nEventsForBGCalculation = nev;}
823 Int_t CountESDTracks();
825 Int_t GetCurrentV0IndexNumber() const {return fCurrentV0IndexNumber;}
827 Bool_t CheckIfPi0IsMother(Int_t label);
828 Bool_t CheckIfEtaIsMother(Int_t label);
830 static void InitESDpid(Int_t type=0);
831 static void SetESDpid(AliESDpid * const pid) {fgESDpid=pid;}
832 static AliESDpid* GetESDpid() {return fgESDpid;}
834 void SetUseChargedTracksMultiplicityForBG(Bool_t flag){fUseChargedTrackMultiplicityForBG = flag;}
836 Int_t GetPindex(Int_t i) {return fV0Pindex.at(i);}
837 Int_t GetNindex(Int_t i) {return fV0Nindex.at(i);}
839 void ResetNGoodV0s(){fNumberOfGoodV0s=0;}
843 AliStack * fMCStack; // pointer to MonteCarlo particle stack
844 // AliMCEventHandler* fMCTruth; // for CF pointer to the MC object
845 AliMCEvent *fMCEvent; // for CF pointer to MC event
846 TChain * fChain; // pointer to the TChain
848 // AliESDInputHandler* fESDHandler; //! pointer to esd object
849 AliESDEvent *fESDEvent; //! pointer to esd object
853 AliCFManager *fCFManager; // pointer to the cf manager
854 // AliCFContainer *container;
856 // for dEdx cut based on nSigma to a particle line
857 //AliESDpid * fESDpid; // esd pid
859 AliGammaConversionHistograms *fHistograms; // pointer to histogram handling class
861 Int_t fCurrentV0IndexNumber; // the current v0 index number
862 AliESDv0 * fCurrentV0; //! pointer to the current v0
863 AliKFParticle * fCurrentNegativeKFParticle; //! pointer to the negative KF particle
864 AliKFParticle * fCurrentPositiveKFParticle; //! pointer to the positive KF particle
865 AliKFParticle * fCurrentMotherKFCandidate; //! pointer to the positive KF particle
867 AliESDtrack * fCurrentNegativeESDTrack; //! pointer to the negative ESD track
868 AliESDtrack * fCurrentPositiveESDTrack; //! pointer to the positive ESD track
870 TLorentzVector * fNegativeTrackLorentzVector; //! pointer to the negative Track Lorentz Vector
871 TLorentzVector * fPositiveTrackLorentzVector; //! pointer to the positive Track Lorentz Vector
872 TLorentzVector * fMotherCandidateLorentzVector; //! pointer to the mother candidate Track Lorentz Vector
874 Double_t fCurrentXValue; // current x value
875 Double_t fCurrentYValue; // current y value
876 Double_t fCurrentZValue; // current z value
878 Int_t fPositiveTrackPID; // positive track pid
879 Int_t fNegativeTrackPID; // negative track pid
881 TParticle *fNegativeMCParticle; //!
882 TParticle *fPositiveMCParticle; //!
883 TParticle *fMotherMCParticle; //!
885 Double_t fMotherCandidateKFMass; // mass of mother candidate KF particle
886 Double_t fMotherCandidateKFWidth; // width of mother candidate KF particle
888 Bool_t fUseKFParticle; // flag
889 Bool_t fUseESDTrack; // flag
890 Bool_t fDoMC; // flag
893 Double_t fMaxVertexZ;
895 Double_t fMaxR; //r cut
896 Double_t fMinR; //r cut
897 Double_t fEtaCut; //eta cut
898 Double_t fRapidityMesonCut; //rapidity for meson cut
899 Double_t fPtCut; // pt cut
900 Double_t fSinglePtCut; // pt cut for electron/positron
901 Double_t fMaxZ; //z cut
903 Double_t fMinClsTPCToF;
904 Double_t fLineCutZRSlope; //linecut
905 Double_t fLineCutZValue; //linecut
906 Double_t fChi2CutConversion; //chi2cut
907 Double_t fChi2CutMeson; //chi2cut
908 Double_t fAlphaCutMeson; //alphacut
909 Double_t fAlphaMinCutMeson; //alphacut
910 Double_t fPIDProbabilityCutNegativeParticle; //pid cut
911 Double_t fPIDProbabilityCutPositiveParticle; //pid cut
912 Bool_t fDodEdxSigmaCut; // flag to use the dEdxCut based on sigmas
913 Double_t fPIDnSigmaAboveElectronLine; // sigma cut
914 Double_t fPIDnSigmaBelowElectronLine; // sigma cut
915 Double_t fPIDnSigmaAbovePionLine; // sigma cut
916 Double_t fPIDMinPnSigmaAbovePionLine; // sigma cut
917 Double_t fPIDMaxPnSigmaAbovePionLine; // sigma cut
918 Double_t fDoKaonRejectionLowP; // Kaon rejection at low p
919 Double_t fDoProtonRejectionLowP; // Proton rejection at low p
920 Double_t fDoPionRejectionLowP; // Pion rejection at low p
921 Double_t fPIDnSigmaAtLowPAroundKaonLine; // sigma cut
922 Double_t fPIDnSigmaAtLowPAroundProtonLine; // sigma cut
923 Double_t fPIDnSigmaAtLowPAroundPionLine; // sigma cut
924 Double_t fPIDMinPKaonRejectionLowP; // Momentum limit to apply kaon rejection
925 Double_t fPIDMinPProtonRejectionLowP; // Momentum limit to apply proton rejection
926 Double_t fPIDMinPPionRejectionLowP; // Momentum limit to apply proton rejection
927 Bool_t fDoQtGammaSelection; // Select gammas using qtMax
928 Double_t fQtMax; // Maximum Qt from Armenteros to select Gammas
929 Double_t fXVertexCut; //vertex cut
930 Double_t fYVertexCut; //vertex cut
931 Double_t fZVertexCut; // vertexcut
933 Double_t fNSigmaMass; //nsigma cut
935 Bool_t fUseImprovedVertex; //flag
937 Bool_t fUseOwnXYZCalculation; //flag that determines if we use our own calculation of xyz (markus)
939 Bool_t fUseConstructGamma; //flag that determines if we use ConstructGamma method from AliKF
943 Bool_t fUseOnFlyV0Finder; //flag
945 Bool_t fUpdateV0AlreadyCalled; //flag
947 TClonesArray* fCurrentEventGoodV0s; //vector of good v0s
949 vector<Int_t> fV0Pindex;
950 vector<Int_t> fV0Nindex;
951 // vector<AliKFParticle> fPreviousEventGoodV0s; // vector of good v0s from prevous events
953 Bool_t fCalculateBackground; //flag
954 AliGammaConversionBGHandler *fBGEventHandler; // background handler
955 Bool_t fBGEventInitialized; //flag
957 AliESDtrackCuts *fEsdTrackCuts; // track cuts
958 Int_t fNumberOfESDTracks; //track counter
960 static AliESDpid* fgESDpid; // ESD pid object
962 Int_t nEventsForBGCalculation;
964 Bool_t fUseChargedTrackMultiplicityForBG;
965 Int_t fNumberOfGoodV0s;
967 ClassDef(AliV0Reader,16)
970 inline void AliV0Reader::InitESDpid(Int_t type)
973 // initialize PID parameters
974 // type=0 is simulation
977 if (!fgESDpid) fgESDpid=new AliESDpid;
978 Double_t alephParameters[5];
980 alephParameters[0] = 2.15898e+00/50.;
981 alephParameters[1] = 1.75295e+01;
982 alephParameters[2] = 3.40030e-09;
983 alephParameters[3] = 1.96178e+00;
984 alephParameters[4] = 3.91720e+00;
985 fgESDpid->GetTOFResponse().SetTimeResolution(80.);
989 alephParameters[0] = 0.0283086/0.97;
990 alephParameters[1] = 2.63394e+01;
991 alephParameters[2] = 5.04114e-11;
992 alephParameters[3] = 2.12543e+00;
993 alephParameters[4] = 4.88663e+00;
994 fgESDpid->GetTOFResponse().SetTimeResolution(130.);
995 fgESDpid->GetTPCResponse().SetMip(50.);
998 fgESDpid->GetTPCResponse().SetBetheBlochParameters(
999 alephParameters[0],alephParameters[1],alephParameters[2],
1000 alephParameters[3],alephParameters[4]);
1002 fgESDpid->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);