]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/AliV0Reader.h
Updated compilation options
[u/mrichter/AliRoot.git] / PWG4 / 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 // --- ROOT system ---
8 #include "TObject.h" 
9 #include "TChain.h"
10 #include "AliESDv0.h"
11 #include "AliESDVertex.h"
12 #include "AliESDInputHandler.h"
13 #include "AliMCEventHandler.h"
14 #include "AliESDEvent.h"
15 #include "AliESDtrack.h"
16 #include "AliKFParticle.h"
17 #include "TParticle.h"
18 #include "AliStack.h"
19 #include "AliGammaConversionHistograms.h"
20 #include <vector>
21
22 class TClonesArray ; 
23 class TFormula ;
24 class TParticle ; 
25 class Riostream ;
26 class TChain;
27 //--- AliRoot system ---
28
29 class AliStack ;
30 class AliESDEvent ; 
31 class AliMCEventHandler;
32 class AliLog ;
33
34 class AliV0Reader : public TObject {
35
36  public: 
37
38   AliV0Reader();                                        //constructor
39   AliV0Reader(const AliV0Reader & g);                   //copy constructor
40   AliV0Reader & operator = (const AliV0Reader & g);     //assignment operator
41   virtual ~AliV0Reader() {;}                            //virtual destructor
42   /*
43    *Initialize the reader
44    */
45   void Initialize();
46
47   /*
48    *Returns the number of v0s in the event, no cuts applied.
49    */
50   Int_t GetNumberOfV0s(){return fESDEvent->GetNumberOfV0s();}
51
52   /*
53    * Check if there are any more good v0s left in the v0 stack
54    * if so, fCurrent v0 is set to this v0 and can be retrieved
55    * by GetCurrentV0 function.
56    * returns kFALSE if there is no more good v0s in the v0 stack
57    */
58   Bool_t NextV0();
59   
60   /*
61    * Returns the v0 at the given index, no checks are done on the v0. 
62    */
63   AliESDv0* GetV0(Int_t index);
64
65   /*
66    * Returns the current v0
67    */
68   AliESDv0* GetCurrentV0(){return fCurrentV0;}
69
70   /*
71    * Returns the negative ESD track which belongs to fCurrentV0
72    */
73   AliESDtrack* GetNegativeESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetNindex());}
74
75   /*
76    * Returns the positive ESD track which belongs to fCurrentV0
77    */
78   AliESDtrack* GetPositiveESDTrack(){return fESDEvent->GetTrack(fCurrentV0->GetPindex());}
79
80   /*
81    * Returns the negative KF particle which belongs to fCurrentV0
82    */
83   AliKFParticle* GetNegativeKFParticle();
84
85   /*
86    * Returns the positive KF particle which belongs to fCurrentV0
87    */
88   AliKFParticle* GetPositiveKFParticle();
89   /*
90    * Returns the KFParticle object of the 2 tracks.
91    */
92   AliKFParticle* GetMotherCandidateKFCombination();
93   /*
94    * Checks the probablity that the PID of the particle is what we want it to be.
95    */
96   Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);
97
98   /*
99    *Get the negative MC TParticle from the stack 
100    */
101   TParticle * GetNegativeMCParticle(){return fNegativeMCParticle;}//fMCStack->Particle(TMath::Abs(fESDEvent->GetTrack(fCurrentV0->GetNindex())->GetLabel()));}
102
103   /*
104    *Get the positive MC TParticle from the stack 
105    */
106   TParticle * GetPositiveMCParticle(){return fPositiveMCParticle;}//{return fMCStack->Particle(TMath::Abs(fESDEvent->GetTrack(fCurrentV0->GetPindex())->GetLabel()));}
107
108   /*
109    *Get the mother MC TParticle from the stack 
110    */
111   TParticle * GetMotherMCParticle(){return fMotherMCParticle;}
112
113   Bool_t HasSameMCMother();
114
115   /*
116    *Get the MC stack 
117    */
118   AliStack* GetMCStack(){return fMCStack;}
119
120   /*
121    *Get the magnetic field from the ESD event 
122    */
123   Double_t GetMagneticField(){return fESDEvent->GetMagneticField();}
124
125   /*
126    *Get the primary vertex from the esd event
127    */
128   const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}
129
130   /*
131    * Set the PID of the negative track
132    */
133   void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}
134
135   /*
136    * Set the PID of the positive track
137    */
138   void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}
139
140   /*
141    * Set the flag to use the kfparticle class. Will also disable the use of esd tracks
142    */
143   void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}
144
145   /*
146    *  Set the flag to use the esd track class. Will also disable the use of kf particles
147    */
148   void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}
149
150   /*
151    *  Set the flag to use improved vertex or not
152    */
153   void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}
154
155   /*
156    * Return the number in the species array belonging to the negative or positive track pid.
157    */
158   Int_t GetSpeciesIndex(Int_t chargeOfTrack);
159
160   /*
161    * Return the x coordinate of the v0
162    */
163   Double_t GetX(){return fCurrentXValue;}
164
165   /*
166    * Return the y coordinate of the v0
167    */
168   Double_t GetY(){return fCurrentYValue;}
169
170   /*
171    * Return the Z coordinate of the v0
172    */
173   Double_t GetZ(){return fCurrentZValue;}
174
175   /*
176    * Return the radius of the v0
177    */
178   Double_t GetXYRadius(){return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}
179
180   /*
181    * Get the opening angle between the two tracks
182    */
183   Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}
184
185   Double_t GetNegativeTrackEnergy(){return fCurrentNegativeKFParticle->E();}
186   Double_t GetPositiveTrackEnergy(){return fCurrentPositiveKFParticle->E();}
187   Double_t GetMotherCandidateEnergy(){return fCurrentMotherKFCandidate->E();}
188
189   Double_t GetNegativeTrackPt(){return fNegativeTrackLorentzVector->Pt();}
190   Double_t GetPositiveTrackPt(){return fPositiveTrackLorentzVector->Pt();}
191   Double_t GetMotherCandidatePt(){return fMotherCandidateLorentzVector->Pt();}
192
193   Double_t GetNegativeTrackEta(){return fNegativeTrackLorentzVector->Eta();}
194   Double_t GetPositiveTrackEta(){return fPositiveTrackLorentzVector->Eta();}
195   Double_t GetMotherCandidateEta(){return fMotherCandidateLorentzVector->Eta();}
196
197   Double_t GetMotherCandidateNDF(){return fCurrentMotherKFCandidate->GetNDF();}
198   Double_t GetMotherCandidateChi2(){return fCurrentMotherKFCandidate->GetChi2();}
199   Double_t GetMotherCandidateMass(){return fMotherCandidateKFMass;}
200   Double_t GetMotherCandidateWidth(){return fMotherCandidateKFWidth;}
201
202   Double_t GetNegativeTrackPhi();
203   Double_t GetPositiveTrackPhi();
204   Double_t GetMotherCandidatePhi();
205
206   void UpdateEventByEventData();
207   
208   Double_t GetMaxRCut(){return fMaxR;}
209   Double_t GetEtaCut(){return fEtaCut;}
210   Double_t GetPtCut(){return fPtCut;}
211   Double_t GetChi2Cut(){return fChi2Cut;}
212
213   void SetMaxRCut(Double_t maxR){fMaxR=maxR;}
214   void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}
215   void SetPtCut(Double_t ptCut){fPtCut=ptCut;}
216   void SetChi2Cut(Double_t chi2){fChi2Cut=chi2;}
217   
218   void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}
219   void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}
220   void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}
221   void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}
222   void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}
223   void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}
224   void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}
225   void UpdateV0Information();
226
227   void SetHistograms(AliGammaConversionHistograms *histograms){fHistograms=histograms;}
228
229   std::vector<AliKFParticle> fCurrentEventGoodV0s;
230   std::vector<AliKFParticle> fPreviousEventGoodV0s;
231
232  private:
233   AliStack * fMCStack;           // pointer to MonteCarlo particle stack 
234   AliMCEventHandler* fMCTruth;   // pointer to the MC event handler
235   TChain * fChain;               // pointer to the TChain
236   
237   AliESDInputHandler* fESDHandler;      //! pointer to esd object
238   AliESDEvent *fESDEvent;               //! pointer to esd object
239
240   AliGammaConversionHistograms *fHistograms;
241   
242   Int_t fCurrentV0IndexNumber;
243   AliESDv0 * fCurrentV0;                //! pointer to the current v0
244   AliKFParticle * fCurrentNegativeKFParticle;  //! pointer to the negative KF particle
245   AliKFParticle * fCurrentPositiveKFParticle;  //! pointer to the positive KF particle
246   AliKFParticle * fCurrentMotherKFCandidate;   //! pointer to the positive KF particle
247
248   AliESDtrack * fCurrentNegativeESDTrack;      //! pointer to the negative ESD track
249   AliESDtrack * fCurrentPositiveESDTrack;      //! pointer to the positive ESD track
250  
251   TLorentzVector * fNegativeTrackLorentzVector; //! pointer to the negative Track Lorentz Vector
252   TLorentzVector * fPositiveTrackLorentzVector; //! pointer to the positive Track Lorentz Vector
253   TLorentzVector * fMotherCandidateLorentzVector;   //! pointer to the mother candidate Track Lorentz Vector
254
255   Double_t fCurrentXValue;
256   Double_t fCurrentYValue;
257   Double_t fCurrentZValue;
258
259   Int_t fPositiveTrackPID;
260   Int_t fNegativeTrackPID;
261
262   TParticle *fNegativeMCParticle;      //!
263   TParticle *fPositiveMCParticle;      //!
264   TParticle *fMotherMCParticle;        //!
265
266   Double_t fMotherCandidateKFMass;
267   Double_t fMotherCandidateKFWidth;
268
269   Bool_t fUseKFParticle;
270   Bool_t fUseESDTrack;
271   Bool_t fDoMC;
272
273   //cuts
274   Double_t fMaxR;
275   Double_t fEtaCut;
276   Double_t fPtCut;
277   Double_t fChi2Cut;
278   Double_t fPIDProbabilityCutNegativeParticle;
279   Double_t fPIDProbabilityCutPositiveParticle;
280   Double_t fXVertexCut;
281   Double_t fYVertexCut;
282   Double_t fZVertexCut;
283   
284   Double_t fNSigmaMass;
285   
286   Bool_t fUseImprovedVertex;
287   
288   ClassDef(AliV0Reader,0)
289 };
290
291
292 #endif
293
294
295