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"
29 //--- AliRoot system ---
32 class AliMCEvent; // for CF
34 class AliMCEventHandler;
35 class AliESDInputHandler;
40 class AliCFManager; // for CF
41 class AliCFContainer; // for CF
42 class AliTPCpidESD; // for dEdx cut based on nSigma to particle lines
45 class AliV0Reader : public TObject {
53 kStepReconstructable = 1,
58 kStepNContributors = 6,
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
77 *Initialize the reader
83 void SetCFManager(AliCFManager * const io){fCFManager = io;};
84 AliCFManager *GetCFManager() const {return fCFManager;}
92 AliESDEvent* GetESDEvent() const{return fESDEvent;}
95 *Returns the number of v0s in the event, no cuts applied.
97 Int_t GetNumberOfV0s() const{return fESDEvent->GetNumberOfV0s();}
100 *Returns the number of contributors to the vertex
102 Int_t GetNumberOfContributorsVtx() const{return fESDEvent->GetPrimaryVertex()->GetNContributors();}
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
113 * Returns the v0 at the given index, no checks are done on the v0.
115 AliESDv0* GetV0(Int_t index);
118 * Returns the current v0
120 AliESDv0* GetCurrentV0() const{return fCurrentV0;}
123 * Returns the negative ESD track which belongs to fCurrentV0
125 AliESDtrack* GetNegativeESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetNindex());}
128 * Returns the positive ESD track which belongs to fCurrentV0
130 AliESDtrack* GetPositiveESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetPindex());}
133 * Returns the negative KF particle which belongs to fCurrentV0
135 AliKFParticle* GetNegativeKFParticle() const{return fCurrentNegativeKFParticle;}
138 * Returns the positive KF particle which belongs to fCurrentV0
140 AliKFParticle* GetPositiveKFParticle() const{return fCurrentPositiveKFParticle;}
143 * Returns the KFParticle object of the 2 tracks.
145 AliKFParticle* GetMotherCandidateKFCombination() const{return fCurrentMotherKFCandidate;}
148 * Checks the probablity that the PID of the particle is what we want it to be.
150 Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);
153 * Checks if the PID of the two particles are within our cuts.
155 void GetPIDProbability(Double_t &negPIDProb, Double_t &posPIDProb);
158 *Get the negative MC TParticle from the stack
160 TParticle * GetNegativeMCParticle() const{return fNegativeMCParticle;}
163 *Get the positive MC TParticle from the stack
165 TParticle * GetPositiveMCParticle() const{return fPositiveMCParticle;}
168 *Get the mother MC TParticle from the stack
170 TParticle * GetMotherMCParticle() const{return fMotherMCParticle;}
173 * Flag to see if the v0 particles share the same mother
175 Bool_t HasSameMCMother();
179 *Get the PID of the MC mother particle
181 Int_t GetMotherMCParticlePDGCode() const{if(fMotherMCParticle != NULL){ cout<<"MCParticle exists"<<endl;} return fMotherMCParticle->GetPdgCode();}
186 AliStack* GetMCStack() const{return fMCStack;}
190 * Setup AliMCEventHandler
192 AliMCEventHandler* GetMCTruth() const{return fMCTruth;} // for CF
198 AliMCEvent* GetMCEvent() const{return fMCEvent;} // for CF
202 *Get the magnetic field from the ESD event
204 Double_t GetMagneticField() const{return fESDEvent->GetMagneticField();}
207 *Get the primary vertex from the esd event
209 const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}
212 * Set the PID of the negative track
214 void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}
217 * Set the PID of the positive track
219 void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}
222 * Set the flag to use the kfparticle class. Will also disable the use of esd tracks
224 void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}
227 * Set the flag to use the esd track class. Will also disable the use of kf particles
229 void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}
232 * Set the flag to use improved vertex or not
234 void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}
237 * Return the number in the species array belonging to the negative or positive track pid.
239 Int_t GetSpeciesIndex(Int_t chargeOfTrack);
242 * Return the x coordinate of the v0
244 Double_t GetX() const{return fCurrentXValue;}
247 * Return the y coordinate of the v0
249 Double_t GetY() const{return fCurrentYValue;}
252 * Return the Z coordinate of the v0
254 Double_t GetZ() const{return fCurrentZValue;}
257 * Return the radius of the v0
259 Double_t GetXYRadius() const{return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}
262 * Get the opening angle between the two tracks
264 Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}
267 * Get the Cos Pointing angle between the two tracks
269 Double_t GetCosPointingAngle(){return fCurrentV0->GetV0CosineOfPointingAngle();}
272 * Get the DCA between the two tracks
274 Double_t GetDcaDaughters(){return fCurrentV0->GetDcaV0Daughters();}
277 * Get the Normalized DCA between the two tracks
279 Double_t GetNormDcaDistDaughters(){return fCurrentV0->GetDcaV0Daughters()/fCurrentV0->GetDistSigma();}
282 * Get the Likelihood for a Conversion
284 Double_t GetLikelihoodAP(){return fCurrentV0->GetLikelihoodAP(0,0);}
287 * Gets the Energy of the negative track.
289 Double_t GetNegativeTrackEnergy() const{return fCurrentNegativeKFParticle->E();}
292 * Gets the Energy of the positive track.
294 Double_t GetPositiveTrackEnergy() const{return fCurrentPositiveKFParticle->E();}
297 * Gets the Energy of the mother candidate.
299 Double_t GetMotherCandidateEnergy() const{return fCurrentMotherKFCandidate->E();}
302 * Gets the Pt of the negative track.
304 Double_t GetNegativeTrackPt() const{return fNegativeTrackLorentzVector->Pt();}
307 * Gets the Pt of the positive track.
309 Double_t GetPositiveTrackPt() const{return fPositiveTrackLorentzVector->Pt();}
313 * Gets the Pt of the mother candidate.
315 Double_t GetMotherCandidatePt() const{return fMotherCandidateLorentzVector->Pt();}
319 * Gets the P of the mother candidate.
321 Double_t GetMotherCandidateP() const{return fMotherCandidateLorentzVector->P();}
325 * Gets the Eta of the negative track.
327 Double_t GetNegativeTrackEta() const{return fNegativeTrackLorentzVector->Eta();}
329 * Gets the Eta of the positive track.
331 Double_t GetPositiveTrackEta() const{return fPositiveTrackLorentzVector->Eta();}
333 * Gets the Eta of the mother candidate.
335 Double_t GetMotherCandidateEta() const{return fMotherCandidateLorentzVector->Eta();}
338 * Gets the NDF of the mother candidate.
340 Double_t GetMotherCandidateNDF() const{return fCurrentMotherKFCandidate->GetNDF();}
343 * Gets the Chi2 of the mother candidate.
345 Double_t GetMotherCandidateChi2() const{return fCurrentMotherKFCandidate->GetChi2();}
348 * Gets the Mass of the mother candidate.
350 Double_t GetMotherCandidateMass() const{return fMotherCandidateKFMass;}
353 * Gets the Width of the mother candidate.
355 Double_t GetMotherCandidateWidth() const{return fMotherCandidateKFWidth;}
358 * Gets the Phi of the negative track.
360 Double_t GetNegativeTrackPhi() const;
363 * Gets the Phi of the positive track.
365 Double_t GetPositiveTrackPhi() const;
368 * Gets the Phi of the mother candidate.
370 Double_t GetMotherCandidatePhi() const;
373 * Gets the Rapidity of the mother candidate.
375 Double_t GetMotherCandidateRapidity() const;
379 * Gets the P of the negative track.
381 Double_t GetNegativeTrackP() const{return fNegativeTrackLorentzVector->P();}
384 * Gets the P of the positive track.
386 Double_t GetPositiveTrackP() const{return fPositiveTrackLorentzVector->P();}
389 * Gets the dE/dx in the TPC of the negative track.
391 Double_t GetNegativeTrackTPCdEdx() const{return fCurrentNegativeESDTrack->GetTPCsignal();}
394 * Gets the dE/dx in the TPC of the positive track.
396 Double_t GetPositiveTrackTPCdEdx() const{return fCurrentPositiveESDTrack->GetTPCsignal();}
399 * Update data which need to be updated every event.
401 void UpdateEventByEventData();
404 * Gets the MaxRCut value.
406 Double_t GetMaxRCut() const{return fMaxR;}
409 * Gets the Eta cut value.
411 Double_t GetEtaCut() const{return fEtaCut;}
414 * Gets the Pt cut value.
416 Double_t GetPtCut() const{return fPtCut;}
420 * Gets the MaxZCut value.
422 Double_t GetMaxZCut() const{return fMaxZ;}
426 * Gets the line cut values.
428 Double_t GetLineCutZRSlope() const{return fLineCutZRSlope;}
429 Double_t GetLineCutZValue() const{return fLineCutZValue;}
432 * Gets the Chi2 cut value for the conversions.
434 Double_t GetChi2CutConversion() const{return fChi2CutConversion;}
437 * Gets the Chi2 cut value for the mesons.
439 Double_t GetChi2CutMeson() const{return fChi2CutMeson;}
441 Double_t GetPositiveTrackLength() const{return fCurrentPositiveESDTrack->GetIntegratedLength();}
442 Double_t GetNegativeTrackLength() const{return fCurrentNegativeESDTrack->GetIntegratedLength();}
444 Double_t GetPositiveNTPCClusters() const{return fCurrentPositiveESDTrack->GetTPCNcls();}
445 Double_t GetNegativeNTPCClusters() const{return fCurrentNegativeESDTrack->GetTPCNcls();}
448 * Sets the MaxRCut value.
450 void SetMaxRCut(Double_t maxR){fMaxR=maxR;}
453 * Sets the EtaCut value.
455 void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}
458 * Sets the PtCut value.
460 void SetPtCut(Double_t ptCut){fPtCut=ptCut;}
464 * Sets the MaxZCut value.
466 void SetMaxZCut(Double_t maxZ){fMaxZ=maxZ;}
470 * Sets the LineCut values.
472 void SetLineCutZRSlope(Double_t LineCutZRSlope){fLineCutZRSlope=LineCutZRSlope;}
473 void SetLineCutZValue(Double_t LineCutZValue){fLineCutZValue=LineCutZValue;}
476 * Sets the Chi2Cut value for conversions.
478 void SetChi2CutConversion(Double_t chi2){fChi2CutConversion=chi2;}
481 * Sets the Chi2Cut for the mesons.
483 void SetChi2CutMeson(Double_t chi2){fChi2CutMeson=chi2;}
486 * Sets the XVertexCut value.
488 void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}
491 * Sets the YVertexCut value.
493 void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}
496 * Sets the ZVertexCut value.
498 void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}
501 * Sets the PIDProbabilityCut value for track particles.
503 void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}
506 * Sets the PIDProbability cut value for the negative track.
508 void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}
511 * Sets the PIDProbability cut value for the positive track.
513 void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}
516 * Sets the PIDnSigmaAboveElectron cut value for the tracks.
518 void SetPIDnSigmaAboveElectronLine(Double_t nSigmaAbove){fPIDnSigmaAboveElectronLine=nSigmaAbove;}
521 * Sets the PIDnSigmaBelowElectron cut value for the tracks.
523 void SetPIDnSigmaBelowElectronLine(Double_t nSigmaBelow){fPIDnSigmaBelowElectronLine=nSigmaBelow;}
526 * Sets the PIDnSigmaAbovePion cut value for the tracks.
528 void SetPIDnSigmaAbovePionLine(Double_t nSigmaAbovePion){fPIDnSigmaAbovePionLine=nSigmaAbovePion;}
531 * Sets the PIDMinPnSigmaAbovePion cut value for the tracks.
533 void SetPIDMinPnSigmaAbovePionLine(Double_t MinPnSigmaAbovePion){fPIDMinPnSigmaAbovePionLine=MinPnSigmaAbovePion;}
536 * Sets the SigmaMassCut value.
538 void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}
541 * Sets the flag to enable/disable the usage of MC information.
543 void SetDoMCTruth(Bool_t doMC){fDoMC = doMC;}
546 * Sets the flag to enable/disable the usage of MC information.
548 Bool_t GetDoMCTruth(){return fDoMC;}
551 * Sets the flag to enable/disable the cut dedx N sigma
554 void SetDodEdxSigmaCut( Bool_t dodEdxSigmaCut){fDodEdxSigmaCut=dodEdxSigmaCut;}
558 * Updates the V0 information of the current V0.
560 Bool_t UpdateV0Information();
563 * Resets the V0 index.
565 void ResetV0IndexNumber(){fCurrentV0IndexNumber=0;}
568 * Sets the histograms.
570 void SetHistograms(AliGammaConversionHistograms * const histograms){fHistograms=histograms;}
573 * Check for primary vertex.
575 Bool_t CheckForPrimaryVertex();
578 * Gets a vector of good v0s.
580 TClonesArray* GetCurrentEventGoodV0s() const{return fCurrentEventGoodV0s;}
583 * Gets the vector of previous events v0s (for bacground analysis)
585 AliGammaConversionKFVector* GetBGGoodV0s(Int_t event);
586 // vector<AliKFParticle> GetPreviousEventGoodV0s() const{return fPreviousEventGoodV0s;}
588 void SetUseOwnXYZCalculation(Bool_t flag){fUseOwnXYZCalculation=flag;}
590 Bool_t GetHelixCenter(AliESDtrack* track, Double_t b,Int_t charge, Double_t center[2]);
592 Bool_t GetConvPosXY(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b, Double_t convpos[2]);
594 Double_t GetConvPosZ(AliESDtrack* ptrack,AliESDtrack* ntrack, Double_t b);
596 void SetDoCF(Bool_t flag){fDoCF = flag;}
598 Bool_t CheckV0FinderStatus(Int_t index);
600 void SetOnFlyFlag(Bool_t flag){fUseOnFlyV0Finder = flag;}
602 Int_t GetNBGEvents(){return fBGEventHandler->GetNBGEvents();}
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
610 AliESDInputHandler* fESDHandler; //! pointer to esd object
611 AliESDEvent *fESDEvent; //! pointer to esd object
615 AliCFManager *fCFManager; // pointer to the cf manager
616 // AliCFContainer *container;
618 // for dEdx cut based on nSigma to a particle line
619 AliTPCpidESD * fTPCpid;
621 AliGammaConversionHistograms *fHistograms; //! pointer to histogram handling class
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
629 AliESDtrack * fCurrentNegativeESDTrack; //! pointer to the negative ESD track
630 AliESDtrack * fCurrentPositiveESDTrack; //! pointer to the positive ESD track
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
636 Double_t fCurrentXValue; // current x value
637 Double_t fCurrentYValue; // current y value
638 Double_t fCurrentZValue; // current z value
640 Int_t fPositiveTrackPID; // positive track pid
641 Int_t fNegativeTrackPID; // negative track pid
643 TParticle *fNegativeMCParticle; //!
644 TParticle *fPositiveMCParticle; //!
645 TParticle *fMotherMCParticle; //!
647 Double_t fMotherCandidateKFMass; // mass of mother candidate KF particle
648 Double_t fMotherCandidateKFWidth; // width of mother candidate KF particle
650 Bool_t fUseKFParticle; // flag
651 Bool_t fUseESDTrack; // flag
652 Bool_t fDoMC; // flag
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
674 Double_t fNSigmaMass; //nsigma cut
676 Bool_t fUseImprovedVertex; //flag
678 Bool_t fUseOwnXYZCalculation; //flag that determines if we use our own calculation of xyz (markus)
682 Bool_t fUseOnFlyV0Finder;
684 Bool_t fUpdateV0AlreadyCalled; //flag
686 TClonesArray* fCurrentEventGoodV0s; //vector of good v0s
687 // vector<AliKFParticle> fPreviousEventGoodV0s; // vector of good v0s from prevous events
689 AliGammaConversionBGHandler *fBGEventHandler;
690 Bool_t fBGEventInitialized;
692 ClassDef(AliV0Reader,8)