]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/GammaConv/AliV0Reader.h
Cleaning up files, added documentation. Added ntuple for v0s.
[u/mrichter/AliRoot.git] / PWG4 / GammaConv / 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   /*\r
95    * Returns the KFParticle object of the 2 tracks.\r
96    */\r
97   AliKFParticle* GetMotherCandidateKFCombination() const{return fCurrentMotherKFCandidate;}\r
98 \r
99   /*\r
100    * Checks the probablity that the PID of the particle is what we want it to be.\r
101    */\r
102   Bool_t CheckPIDProbability(Double_t negProbCut, Double_t posProbCut);\r
103   \r
104   /*\r
105    * Checks if the PID of the two particles are within our cuts.\r
106    */\r
107   void GetPIDProbability(Double_t &negPIDProb, Double_t &posPIDProb);\r
108 \r
109   /*\r
110    *Get the negative MC TParticle from the stack \r
111    */\r
112   TParticle * GetNegativeMCParticle() const{return fNegativeMCParticle;}\r
113 \r
114   /*\r
115    *Get the positive MC TParticle from the stack \r
116    */\r
117   TParticle * GetPositiveMCParticle() const{return fPositiveMCParticle;}\r
118 \r
119   /*\r
120    *Get the mother MC TParticle from the stack \r
121    */\r
122   TParticle * GetMotherMCParticle() const{return fMotherMCParticle;}\r
123 \r
124   /*\r
125    * Flag to see if the v0 particles share the same mother\r
126    */\r
127   Bool_t HasSameMCMother();\r
128 \r
129 \r
130   /*\r
131    *Get the PID of the MC mother particle\r
132    */\r
133   Int_t GetMotherMCParticlePDGCode() const{if(fMotherMCParticle != NULL){ cout<<"MCParticle exists"<<endl;} return fMotherMCParticle->GetPdgCode();}\r
134 \r
135   /*\r
136    *Get the MC stack \r
137    */\r
138   AliStack* GetMCStack() const{return fMCStack;}\r
139 \r
140   /*\r
141    *Get the magnetic field from the ESD event \r
142    */\r
143   Double_t GetMagneticField() const{return fESDEvent->GetMagneticField();}\r
144 \r
145   /*\r
146    *Get the primary vertex from the esd event\r
147    */\r
148   const AliESDVertex *GetPrimaryVertex() const {return fESDEvent->GetPrimaryVertex();}\r
149 \r
150   /*\r
151    * Set the PID of the negative track\r
152    */\r
153   void SetNegativeTrackPID(Int_t negTrackPID){fNegativeTrackPID=negTrackPID;}\r
154 \r
155   /*\r
156    * Set the PID of the positive track\r
157    */\r
158   void SetPositiveTrackPID(Int_t posTrackPID){fPositiveTrackPID=posTrackPID;}\r
159 \r
160   /*\r
161    * Set the flag to use the kfparticle class. Will also disable the use of esd tracks\r
162    */\r
163   void UseKFParticle(){fUseKFParticle = kTRUE; fUseESDTrack = kFALSE;}\r
164 \r
165   /*\r
166    *  Set the flag to use the esd track class. Will also disable the use of kf particles\r
167    */\r
168   void UseESDTrack(){fUseESDTrack = kTRUE; fUseKFParticle = kFALSE;}\r
169 \r
170   /*\r
171    *  Set the flag to use improved vertex or not\r
172    */\r
173   void SetUseImprovedVertex(Bool_t useImprovedVertex){fUseImprovedVertex=useImprovedVertex;}\r
174 \r
175   /*\r
176    * Return the number in the species array belonging to the negative or positive track pid.\r
177    */\r
178   Int_t GetSpeciesIndex(Int_t chargeOfTrack);\r
179 \r
180   /*\r
181    * Return the x coordinate of the v0\r
182    */\r
183   Double_t GetX() const{return fCurrentXValue;}\r
184 \r
185   /*\r
186    * Return the y coordinate of the v0\r
187    */\r
188   Double_t GetY() const{return fCurrentYValue;}\r
189 \r
190   /*\r
191    * Return the Z coordinate of the v0\r
192    */\r
193   Double_t GetZ() const{return fCurrentZValue;}\r
194 \r
195   /*\r
196    * Return the radius of the v0\r
197    */\r
198   Double_t GetXYRadius() const{return sqrt((Double_t)(fCurrentXValue*fCurrentXValue + fCurrentYValue*fCurrentYValue));}\r
199 \r
200   /*\r
201    * Get the opening angle between the two tracks\r
202    */\r
203   Double_t GetOpeningAngle(){return fNegativeTrackLorentzVector->Angle(fPositiveTrackLorentzVector->Vect());}\r
204 \r
205   /*\r
206    * Gets the Energy of the negative track.\r
207    */\r
208   Double_t GetNegativeTrackEnergy() const{return fCurrentNegativeKFParticle->E();}\r
209 \r
210   /*\r
211    * Gets the Energy of the positive track.\r
212    */\r
213   Double_t GetPositiveTrackEnergy() const{return fCurrentPositiveKFParticle->E();}\r
214 \r
215   /*\r
216    * Gets the Energy of the mother candidate.\r
217    */\r
218   Double_t GetMotherCandidateEnergy() const{return fCurrentMotherKFCandidate->E();}\r
219 \r
220   /*\r
221    * Gets the Pt of the negative track.\r
222    */\r
223   Double_t GetNegativeTrackPt() const{return fNegativeTrackLorentzVector->Pt();}\r
224 \r
225   /*\r
226    * Gets the Pt of the positive track.\r
227    */\r
228   Double_t GetPositiveTrackPt() const{return fPositiveTrackLorentzVector->Pt();}\r
229 \r
230   /*\r
231    * Gets the Pt of the mother candidate.\r
232    */\r
233   Double_t GetMotherCandidatePt() const{return fMotherCandidateLorentzVector->Pt();}\r
234 \r
235   /*\r
236    * Gets the Eta of the negative track.\r
237    */\r
238   Double_t GetNegativeTrackEta() const{return fNegativeTrackLorentzVector->Eta();}\r
239   /*\r
240    * Gets the Eta of the positive track.\r
241    */\r
242   Double_t GetPositiveTrackEta() const{return fPositiveTrackLorentzVector->Eta();}\r
243   /*\r
244    * Gets the Eta of the mother candidate.\r
245    */\r
246   Double_t GetMotherCandidateEta() const{return fMotherCandidateLorentzVector->Eta();}\r
247 \r
248   /*\r
249    * Gets the NDF of the mother candidate.\r
250    */\r
251   Double_t GetMotherCandidateNDF() const{return fCurrentMotherKFCandidate->GetNDF();}\r
252 \r
253   /*\r
254    * Gets the Chi2 of the mother candidate.\r
255    */\r
256   Double_t GetMotherCandidateChi2() const{return fCurrentMotherKFCandidate->GetChi2();}\r
257 \r
258   /*\r
259    * Gets the Mass of the mother candidate.\r
260    */\r
261   Double_t GetMotherCandidateMass() const{return fMotherCandidateKFMass;}\r
262 \r
263   /*\r
264    * Gets the Width of the mother candidate.\r
265    */\r
266   Double_t GetMotherCandidateWidth() const{return fMotherCandidateKFWidth;}\r
267 \r
268   /*\r
269    * Gets the Phi of the negative track.\r
270    */\r
271   Double_t GetNegativeTrackPhi() const;\r
272 \r
273   /*\r
274    * Gets the Phi of the positive track.\r
275    */\r
276   Double_t GetPositiveTrackPhi() const;\r
277 \r
278   /*\r
279    * Gets the Phi of the mother candidate.\r
280    */\r
281   Double_t GetMotherCandidatePhi() const;\r
282 \r
283   /*\r
284    * Update data which need to be updated every event.\r
285    */\r
286   void UpdateEventByEventData();\r
287   \r
288   /*\r
289    * Gets the MaxRCut value.\r
290    */\r
291   Double_t GetMaxRCut() const{return fMaxR;}\r
292 \r
293   /*\r
294    * Gets the Eta cut value.\r
295    */\r
296   Double_t GetEtaCut() const{return fEtaCut;}\r
297 \r
298   /*\r
299    * Gets the Pt cut value.\r
300    */\r
301   Double_t GetPtCut() const{return fPtCut;}\r
302 \r
303   /*\r
304    * Gets the Chi2 cut value for the conversions.\r
305    */\r
306   Double_t GetChi2CutConversion() const{return fChi2CutConversion;}\r
307 \r
308   /*\r
309    * Gets the Chi2 cut value for the mesons.\r
310    */\r
311   Double_t GetChi2CutMeson() const{return fChi2CutMeson;}\r
312 \r
313   /*\r
314    * Sets the MaxRCut value.\r
315    */\r
316   void SetMaxRCut(Double_t maxR){fMaxR=maxR;}\r
317 \r
318   /*\r
319    * Sets the EtaCut value.\r
320    */\r
321   void SetEtaCut(Double_t etaCut){fEtaCut=etaCut;}\r
322 \r
323   /*\r
324    * Sets the PtCut value.\r
325    */\r
326   void SetPtCut(Double_t ptCut){fPtCut=ptCut;}\r
327 \r
328   /*\r
329    * Sets the Chi2Cut value for conversions.\r
330    */\r
331   void SetChi2CutConversion(Double_t chi2){fChi2CutConversion=chi2;}\r
332 \r
333   /*\r
334    * Sets the Chi2Cut for the mesons.\r
335    */\r
336   void SetChi2CutMeson(Double_t chi2){fChi2CutMeson=chi2;}\r
337   \r
338   /*\r
339    * Sets the XVertexCut value.\r
340    */\r
341   void SetXVertexCut(Double_t xVtx){fCurrentXValue=xVtx;}\r
342 \r
343   /*\r
344    * Sets the YVertexCut value.\r
345    */\r
346   void SetYVertexCut(Double_t yVtx){fCurrentYValue=yVtx;}\r
347 \r
348   /*\r
349    * Sets the ZVertexCut value.\r
350    */\r
351   void SetZVertexCut(Double_t zVtx){fCurrentZValue=zVtx;}\r
352 \r
353   /*\r
354    * Sets the PIDProbabilityCut value for track particles.\r
355    */\r
356   void SetPIDProbability(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb; fPIDProbabilityCutNegativeParticle=pidProb;}\r
357 \r
358   /*\r
359    * Sets the PIDProbability cut value for the negative track.\r
360    */\r
361   void SetPIDProbabilityNegativeParticle(Double_t pidProb){fPIDProbabilityCutNegativeParticle=pidProb;}\r
362 \r
363   /*\r
364    * Sets the PIDProbability cut value for the positive track.\r
365    */\r
366   void SetPIDProbabilityPositiveParticle(Double_t pidProb){fPIDProbabilityCutPositiveParticle=pidProb;}\r
367 \r
368   /*\r
369    * Sets the SigmaMassCut value.\r
370    */\r
371   void SetSigmaMass(Double_t sigmaMass){fNSigmaMass=sigmaMass;}\r
372 \r
373   /*\r
374    * Sets the flag to enable/disable the usage of MC information. \r
375    */\r
376   void SetDoMCTruth(Bool_t doMC){fDoMC = doMC;}\r
377 \r
378   /*\r
379    * Updates the V0 information of the current V0.\r
380    */\r
381   void UpdateV0Information();\r
382 \r
383   /*\r
384    * Resets the V0 index.\r
385    */\r
386   void ResetV0IndexNumber(){fCurrentV0IndexNumber=0;}\r
387 \r
388   /*\r
389    * Sets the histograms.\r
390    */\r
391   void SetHistograms(AliGammaConversionHistograms *histograms){fHistograms=histograms;}\r
392 \r
393   /*\r
394    * Check for primary vertex.\r
395    */\r
396   Bool_t CheckForPrimaryVertex();\r
397 \r
398   /*\r
399    * Gets a vector of good v0s.\r
400    */\r
401   vector<AliKFParticle> GetCurrentEventGoodV0s() const{return fCurrentEventGoodV0s;}\r
402 \r
403   /*\r
404    * Gets the vector of previous events v0s (for bacground analysis)\r
405    */\r
406   vector<AliKFParticle> GetPreviousEventGoodV0s() const{return fPreviousEventGoodV0s;}\r
407 \r
408  private:\r
409   AliStack * fMCStack;           // pointer to MonteCarlo particle stack \r
410   AliMCEventHandler* fMCTruth;   // pointer to the MC event handler\r
411   TChain * fChain;               // pointer to the TChain\r
412   \r
413   AliESDInputHandler* fESDHandler;      //! pointer to esd object\r
414   AliESDEvent *fESDEvent;               //! pointer to esd object\r
415 \r
416   AliGammaConversionHistograms *fHistograms;\r
417   \r
418   Int_t fCurrentV0IndexNumber;\r
419   AliESDv0 * fCurrentV0;                //! pointer to the current v0\r
420   AliKFParticle * fCurrentNegativeKFParticle;  //! pointer to the negative KF particle\r
421   AliKFParticle * fCurrentPositiveKFParticle;  //! pointer to the positive KF particle\r
422   AliKFParticle * fCurrentMotherKFCandidate;   //! pointer to the positive KF particle\r
423 \r
424   AliESDtrack * fCurrentNegativeESDTrack;      //! pointer to the negative ESD track\r
425   AliESDtrack * fCurrentPositiveESDTrack;      //! pointer to the positive ESD track\r
426  \r
427   TLorentzVector * fNegativeTrackLorentzVector; //! pointer to the negative Track Lorentz Vector\r
428   TLorentzVector * fPositiveTrackLorentzVector; //! pointer to the positive Track Lorentz Vector\r
429   TLorentzVector * fMotherCandidateLorentzVector;   //! pointer to the mother candidate Track Lorentz Vector\r
430 \r
431   Double_t fCurrentXValue;\r
432   Double_t fCurrentYValue;\r
433   Double_t fCurrentZValue;\r
434 \r
435   Int_t fPositiveTrackPID;\r
436   Int_t fNegativeTrackPID;\r
437 \r
438   TParticle *fNegativeMCParticle;      //!\r
439   TParticle *fPositiveMCParticle;      //!\r
440   TParticle *fMotherMCParticle;        //!\r
441 \r
442   Double_t fMotherCandidateKFMass;\r
443   Double_t fMotherCandidateKFWidth;\r
444 \r
445   Bool_t fUseKFParticle;\r
446   Bool_t fUseESDTrack;\r
447   Bool_t fDoMC;\r
448 \r
449   //cuts\r
450   Double_t fMaxR;\r
451   Double_t fEtaCut;\r
452   Double_t fPtCut;\r
453   Double_t fChi2CutConversion;\r
454   Double_t fChi2CutMeson;\r
455   Double_t fPIDProbabilityCutNegativeParticle;\r
456   Double_t fPIDProbabilityCutPositiveParticle;\r
457   Double_t fXVertexCut;\r
458   Double_t fYVertexCut;\r
459   Double_t fZVertexCut;\r
460   \r
461   Double_t fNSigmaMass;\r
462   \r
463   Bool_t fUseImprovedVertex;\r
464   \r
465   vector<AliKFParticle> fCurrentEventGoodV0s;\r
466   vector<AliKFParticle> fPreviousEventGoodV0s;\r
467 \r
468   ClassDef(AliV0Reader,0)\r
469 };\r
470 \r
471 \r
472 #endif\r
473 \r
474 \r
475 \r