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