]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemto/AliFemtoEventReaderESDChain.cxx
17109b474bf32ffc76837b4d0d9540ff358fd7a3
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemto / AliFemtoEventReaderESDChain.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 ///                                                                          ///
3 /// AliFemtoEventReaderESDChain - the reader class for the Alice ESD         ///
4 /// tailored for the Task framework                                 ///
5 /// Reads in AliESDfriend to create shared hit/quality information           ///
6 /// Authors: Adam Kisiel kisiel@mps.ohio-state.edu                           ///
7 ///                                                                          ///
8 ////////////////////////////////////////////////////////////////////////////////
9 #include "AliFemtoEventReaderESDChain.h"
10
11 #include "TFile.h"
12 #include "TTree.h"
13 #include "AliESDEvent.h"
14 #include "AliESDtrack.h"
15 #include "AliESDVertex.h"
16 #include "AliMultiplicity.h"
17 #include "AliCentrality.h"
18 #include "AliEventplane.h"
19 #include "AliESDVZERO.h"
20 #include "AliFmPhysicalHelixD.h"
21 #include "AliFmThreeVectorF.h"
22 #include "SystemOfUnits.h"
23 #include "AliFemtoEvent.h"
24 #include "AliFemtoModelHiddenInfo.h"
25 #include "AliPID.h"
26
27 ClassImp(AliFemtoEventReaderESDChain)
28
29 #if !(ST_NO_NAMESPACES)
30   using namespace units;
31 #endif
32
33 using namespace std;
34
35 //____________________________
36 AliFemtoEventReaderESDChain::AliFemtoEventReaderESDChain():
37   fFileName(" "),
38   fConstrained(true),
39   fReadInner(false),
40   fUseTPCOnly(false),
41   fNumberofEvent(0),
42   fCurEvent(0),
43   fCurFile(0),
44   fEvent(0x0),
45   fUsePhysicsSel(kFALSE),
46   fSelect(0x0),
47   fTrackType(kGlobal),
48   fEstEventMult(kITSTPC), 
49   fEventTrig(AliVEvent::kMB), //trigger
50   fESDpid(0),
51   fIsPidOwner(0),
52   fReadV0(0)
53 {
54   //constructor with 0 parameters , look at default settings 
55   //   fClusterPerPadrow = (list<Int_t> **) malloc(sizeof(list<Int_t> *) * AliESDfriendTrack::kMaxTPCcluster);
56   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
57   //     fClusterPerPadrow[tPad] = new list<Int_t>();
58   //   }
59   //   fSharedList = (list<Int_t> **) malloc(sizeof(list<Int_t> *) * AliESDfriendTrack::kMaxTPCcluster);
60   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
61   //     fSharedList[tPad] = new list<Int_t>();
62   //   }
63 }
64
65 //__________________
66 AliFemtoEventReaderESDChain::AliFemtoEventReaderESDChain(const AliFemtoEventReaderESDChain& aReader):
67   AliFemtoEventReader(aReader),
68   fFileName(" "),
69   fConstrained(true),
70   fReadInner(false),
71   fUseTPCOnly(false),
72   fNumberofEvent(0),
73   fCurEvent(0),
74   fCurFile(0),
75   fEvent(0x0),
76   fUsePhysicsSel(kFALSE),
77   fSelect(0x0),
78   fTrackType(kGlobal),
79   fEstEventMult(kITSTPC),
80   fEventTrig(AliVEvent::kMB), //trigger
81   fESDpid(0),
82   fIsPidOwner(0),
83   fReadV0(0)
84 {
85   // Copy constructor
86   fConstrained = aReader.fConstrained;
87   fReadInner = aReader.fReadInner;
88   fUseTPCOnly = aReader.fUseTPCOnly;
89   fNumberofEvent = aReader.fNumberofEvent;
90   fCurEvent = aReader.fCurEvent;
91   fCurFile = aReader.fCurFile;
92   //  fEvent = new AliESD(*aReader.fEvent);
93   fEvent = new AliESDEvent();
94   fUsePhysicsSel = aReader.fUsePhysicsSel;
95   if (aReader.fUsePhysicsSel)
96     fSelect = new AliPhysicsSelection();
97   fTrackType = aReader.fTrackType;
98   fEstEventMult = aReader.fEstEventMult;
99   fEventTrig = aReader.fEventTrig; //trigger
100
101   //   fEventFriend = aReader.fEventFriend;
102   //   fClusterPerPadrow = (list<Int_t> **) malloc(sizeof(list<Int_t> *) * AliESDfriendTrack::kMaxTPCcluster);
103   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
104   //     fClusterPerPadrow[tPad] = new list<Int_t>();
105   //     list<Int_t>::iterator iter;
106   //     for (iter=aReader.fClusterPerPadrow[tPad]->begin(); iter!=aReader.fClusterPerPadrow[tPad]->end(); iter++) {
107   //       fClusterPerPadrow[tPad]->push_back(*iter);
108   //     }
109   //   }
110   //   fSharedList = (list<Int_t> **) malloc(sizeof(list<Int_t> *) * AliESDfriendTrack::kMaxTPCcluster);
111   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
112   //     fSharedList[tPad] = new list<Int_t>();
113   //     list<Int_t>::iterator iter;
114   //     for (iter=aReader.fSharedList[tPad]->begin(); iter!=aReader.fSharedList[tPad]->end(); iter++) {
115   //       fSharedList[tPad]->push_back(*iter);
116   //     }
117   //   }
118 }
119 //__________________
120 AliFemtoEventReaderESDChain::~AliFemtoEventReaderESDChain()
121 {
122   //Destructor
123   delete fEvent;
124
125   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
126   //     fClusterPerPadrow[tPad]->clear();
127   //     delete fClusterPerPadrow[tPad];
128   //   }
129   //   delete [] fClusterPerPadrow;
130   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
131   //     fSharedList[tPad]->clear();
132   //     delete fSharedList[tPad];
133   //   }
134   //   delete [] fSharedList;
135   if (fSelect) delete fSelect;
136 }
137
138 //__________________
139 AliFemtoEventReaderESDChain& AliFemtoEventReaderESDChain::operator=(const AliFemtoEventReaderESDChain& aReader)
140 {
141   // Assignment operator
142   if (this == &aReader)
143     return *this;
144
145   fConstrained = aReader.fConstrained;
146   fReadInner = aReader.fReadInner;
147   fUseTPCOnly = aReader.fUseTPCOnly;
148   fNumberofEvent = aReader.fNumberofEvent;
149   fCurEvent = aReader.fCurEvent;
150   fCurFile = aReader.fCurFile;
151   if (fEvent) delete fEvent;
152   fEvent = new AliESDEvent();
153   fTrackType = aReader.fTrackType;
154   fEstEventMult = aReader.fEstEventMult;
155
156   fUsePhysicsSel = aReader.fUsePhysicsSel;
157   fReadV0 = aReader.fReadV0;
158   if (aReader.fUsePhysicsSel)
159     fSelect = new AliPhysicsSelection();
160   //  fEventFriend = aReader.fEventFriend;
161   
162   //   if (fClusterPerPadrow) {
163   //     for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
164   //       fClusterPerPadrow[tPad]->clear();
165   //       delete fClusterPerPadrow[tPad];
166   //     }
167   //     delete [] fClusterPerPadrow;
168   //   }
169   
170   //   if (fSharedList) {
171   //     for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
172   //       fSharedList[tPad]->clear();
173   //       delete fSharedList[tPad];
174   //     }
175   //     delete [] fSharedList;
176   //   }
177
178   //   fClusterPerPadrow = (list<Int_t> **) malloc(sizeof(list<Int_t> *) * AliESDfriendTrack::kMaxTPCcluster);
179   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
180   //     fClusterPerPadrow[tPad] = new list<Int_t>();
181   //     list<Int_t>::iterator iter;
182   //     for (iter=aReader.fClusterPerPadrow[tPad]->begin(); iter!=aReader.fClusterPerPadrow[tPad]->end(); iter++) {
183   //       fClusterPerPadrow[tPad]->push_back(*iter);
184   //     }
185   //   }
186   //   fSharedList = (list<Int_t> **) malloc(sizeof(list<Int_t> *) * AliESDfriendTrack::kMaxTPCcluster);
187   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
188   //     fSharedList[tPad] = new list<Int_t>();
189   //     list<Int_t>::iterator iter;
190   //     for (iter=aReader.fSharedList[tPad]->begin(); iter!=aReader.fSharedList[tPad]->end(); iter++) {
191   //       fSharedList[tPad]->push_back(*iter);
192   //     }
193   //   }
194   
195   return *this;
196 }
197 //__________________
198 // Simple report
199 AliFemtoString AliFemtoEventReaderESDChain::Report()
200 {
201   AliFemtoString temp = "\n This is the AliFemtoEventReaderESDChain\n";
202   return temp;
203 }
204
205 //__________________
206 void AliFemtoEventReaderESDChain::SetConstrained(const bool constrained)
207 {
208   // Select whether to read constrained or not constrained momentum
209   fConstrained=constrained;
210 }
211 //__________________
212 bool AliFemtoEventReaderESDChain::GetConstrained() const
213 {
214   // Check whether we read constrained or not constrained momentum
215   return fConstrained;
216 }
217 //__________________
218 void AliFemtoEventReaderESDChain::SetReadTPCInner(const bool readinner)
219 {
220   fReadInner=readinner;
221 }
222
223 bool AliFemtoEventReaderESDChain::GetReadTPCInner() const
224 {
225   return fReadInner;
226 }
227
228 //__________________
229 void AliFemtoEventReaderESDChain::SetUseTPCOnly(const bool usetpconly)
230 {
231   fUseTPCOnly=usetpconly;
232 }
233
234 bool AliFemtoEventReaderESDChain::GetUseTPCOnly() const
235 {
236   return fUseTPCOnly;
237 }
238
239 void AliFemtoEventReaderESDChain::SetUsePhysicsSelection(const bool usephysics)
240 {
241   fUsePhysicsSel = usephysics;
242   if (!fSelect) fSelect = new AliPhysicsSelection();
243 }
244
245 void AliFemtoEventReaderESDChain::SetUseMultiplicity(EstEventMult aType)
246 {
247   fEstEventMult = aType;
248 }
249
250 AliFemtoEvent* AliFemtoEventReaderESDChain::ReturnHbtEvent()
251 {
252   // Get the event, read all the relevant information
253   // and fill the AliFemtoEvent class
254   // Returns a valid AliFemtoEvent
255   AliFemtoEvent *hbtEvent = 0;
256   string tFriendFileName;
257
258   // Get the friend information
259   if (Debug()>1) cout<<"starting to read event "<<fCurEvent<<endl;
260   //  fEvent->SetESDfriend(fEventFriend);
261   if(fEvent->GetAliESDOld())fEvent->CopyFromOldESD();
262   
263   hbtEvent = new AliFemtoEvent;
264
265   if (fUsePhysicsSel) {
266     hbtEvent->SetIsCollisionCandidate(fSelect->IsCollisionCandidate(fEvent));
267     if (!(fSelect->IsCollisionCandidate(fEvent)))
268       printf("Event not a collision candidate\n");
269   }
270   else
271     hbtEvent->SetIsCollisionCandidate(kTRUE);
272
273   //setting basic things
274   //  hbtEvent->SetEventNumber(fEvent->GetEventNumber());
275   hbtEvent->SetRunNumber(fEvent->GetRunNumber());
276   //hbtEvent->SetNumberOfTracks(fEvent->GetNumberOfTracks());
277   hbtEvent->SetMagneticField(fEvent->GetMagneticField()*kilogauss);//to check if here is ok
278   hbtEvent->SetZDCN1Energy(fEvent->GetZDCN1Energy());
279   hbtEvent->SetZDCP1Energy(fEvent->GetZDCP1Energy());
280   hbtEvent->SetZDCN2Energy(fEvent->GetZDCN2Energy());
281   hbtEvent->SetZDCP2Energy(fEvent->GetZDCP2Energy());
282   hbtEvent->SetZDCEMEnergy(fEvent->GetZDCEMEnergy());
283   hbtEvent->SetZDCParticipants(fEvent->GetZDCParticipants());
284   hbtEvent->SetTriggerMask(fEvent->GetTriggerMask());
285   //  hbtEvent->SetTriggerCluster(fEvent->GetTriggerCluster());
286
287   if ((fEvent->IsTriggerClassFired("CINT1WU-B-NOPF-ALL")) ||
288       (fEvent->IsTriggerClassFired("CINT1B-ABCE-NOPF-ALL")) ||
289       (fEvent->IsTriggerClassFired("CINT1-B-NOPF-ALLNOTRD")) ||
290       (fEvent->IsTriggerClassFired("CINT1-B-NOPF-FASTNOTRD")))
291     hbtEvent->SetTriggerCluster(1);
292   else if ((fEvent->IsTriggerClassFired("CSH1WU-B-NOPF-ALL")) ||
293            (fEvent->IsTriggerClassFired("CSH1-B-NOPF-ALLNOTRD")))
294     hbtEvent->SetTriggerCluster(2);
295   else 
296     hbtEvent->SetTriggerCluster(0);
297         
298   //Vertex
299   double fV1[3];
300   double fVCov[6];
301   if (fUseTPCOnly) {
302     fEvent->GetPrimaryVertexTPC()->GetXYZ(fV1);
303     fEvent->GetPrimaryVertexTPC()->GetCovMatrix(fVCov);
304     if (!fEvent->GetPrimaryVertexTPC()->GetStatus())
305       fVCov[4] = -1001.0;
306   }
307   else {
308     fEvent->GetPrimaryVertex()->GetXYZ(fV1);
309     fEvent->GetPrimaryVertex()->GetCovMatrix(fVCov);
310     if (!fEvent->GetPrimaryVertex()->GetStatus())
311       fVCov[4] = -1001.0;
312   }
313
314   AliFmThreeVectorF vertex(fV1[0],fV1[1],fV1[2]);
315   hbtEvent->SetPrimVertPos(vertex);
316   hbtEvent->SetPrimVertCov(fVCov);
317   
318   Int_t spdetaonecount = 0;
319   
320   for (int iter=0; iter<fEvent->GetMultiplicity()->GetNumberOfTracklets(); iter++) 
321     if (fabs(fEvent->GetMultiplicity()->GetEta(iter)) < 1.0)
322       spdetaonecount++;
323
324   //  hbtEvent->SetSPDMult(fEvent->GetMultiplicity()->GetNumberOfTracklets());
325   hbtEvent->SetSPDMult(spdetaonecount);
326
327   //starting to reading tracks
328   int nofTracks=0;  //number of reconstructed tracks in event
329   nofTracks=fEvent->GetNumberOfTracks();
330   int realnofTracks=0;//number of track which we use ina analysis
331
332   //   // Clear the shared cluster list
333   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
334   //     fClusterPerPadrow[tPad]->clear();
335   //   }
336   //   for (int tPad=0; tPad<AliESDfriendTrack::kMaxTPCcluster; tPad++) {
337   //     fSharedList[tPad]->clear();
338   //   }
339
340
341   //   for (int i=0;i<nofTracks;i++) {
342   //     const AliESDtrack *esdtrack=fEvent->GetTrack(i);//getting next track
343
344   //     list<Int_t>::iterator tClustIter;
345
346   //     Int_t tTrackIndices[AliESDfriendTrack::kMaxTPCcluster];
347   //     Int_t tNClusters = esdtrack->GetTPCclusters(tTrackIndices);
348   //     for (int tNcl=0; tNcl<AliESDfriendTrack::kMaxTPCcluster; tNcl++) {
349   //       if (tTrackIndices[tNcl] >= 0) {
350   //    tClustIter = find(fClusterPerPadrow[tNcl]->begin(), fClusterPerPadrow[tNcl]->end(), tTrackIndices[tNcl]);
351   //    if (tClustIter == fClusterPerPadrow[tNcl]->end()) {
352   //      fClusterPerPadrow[tNcl]->push_back(tTrackIndices[tNcl]);
353   //    }
354   //    else {
355   //      fSharedList[tNcl]->push_back(tTrackIndices[tNcl]);
356   //    }
357   //       }
358   //     }
359       
360   //   }
361
362   int tNormMult = 0;
363   int tNormMultPos = 0;
364   int tNormMultNeg = 0;
365
366   Float_t tTotalPt = 0.0;
367
368   Float_t b[2];
369   Float_t bCov[3];
370
371   Int_t tTracklet=0, tITSTPC=0, tITSPure=0;
372   
373   fEvent->EstimateMultiplicity(tTracklet, tITSTPC, tITSPure, 1.2);
374   
375   hbtEvent->SetMultiplicityEstimateITSTPC(tITSTPC);
376   hbtEvent->SetMultiplicityEstimateTracklets(tTracklet);
377   //  hbtEvent->SetMultiplicityEstimateITSPure(tITSPure);
378   hbtEvent->SetMultiplicityEstimateITSPure(fEvent->GetMultiplicity()->GetNumberOfITSClusters(1));
379   
380   for (int i=0;i<nofTracks;i++)
381     {
382       bool  tGoodMomentum=true; //flaga to chcek if we can read momentum of this track
383
384       const AliESDtrack *esdtrack=fEvent->GetTrack(i);//getting next track
385       //      const AliESDfriendTrack *tESDfriendTrack = esdtrack->GetFriendTrack();
386
387       if ((esdtrack->GetStatus() & AliESDtrack::kTPCrefit) &&
388           (esdtrack->GetStatus() & AliESDtrack::kITSrefit)) {
389         if (esdtrack->GetTPCNcls() > 70) 
390           if (esdtrack->GetTPCchi2()/esdtrack->GetTPCNcls() < 4.0) {
391             if (TMath::Abs(esdtrack->Eta()) < 1.2) {
392               esdtrack->GetImpactParameters(b,bCov);
393               if ((b[0]<0.2) && (b[1] < 0.25)) {
394                 tNormMult++;
395                 tTotalPt += esdtrack->Pt();
396               }
397             }
398           }
399       }
400       else if (esdtrack->GetStatus() & AliESDtrack::kTPCrefit) {
401         if (esdtrack->GetTPCNcls() > 100) 
402           if (esdtrack->GetTPCchi2()/esdtrack->GetTPCNcls() < 4.0) {
403             if (TMath::Abs(esdtrack->Eta()) < 1.2) {
404               esdtrack->GetImpactParameters(b,bCov);
405               if ((b[0]<2.4) && (b[1] < 3.2)) {
406                 tNormMult++;
407                 tTotalPt += esdtrack->Pt();
408               }
409             }
410           }
411       }
412       
413       hbtEvent->SetZDCEMEnergy(tTotalPt);
414       //       if (esdtrack->GetStatus() & AliESDtrack::kTPCrefit)
415       //        if (esdtrack->GetTPCNcls() > 80) 
416       //          if (esdtrack->GetTPCchi2()/esdtrack->GetTPCNcls() < 6.0) 
417       //            if (esdtrack->GetConstrainedParam())
418       //              if (fabs(esdtrack->GetConstrainedParam()->Eta()) < 0.5)
419       //                if (esdtrack->GetConstrainedParam()->Pt() < 1.0) {
420       //                  if (esdtrack->GetSign() > 0)
421       //                    tNormMultPos++;
422       //                  else if (esdtrack->GetSign() < 0)
423       //                    tNormMultNeg--;
424       //                }
425
426       // If reading ITS-only tracks, reject all with TPC
427       if (fTrackType == kITSOnly) {
428         if (esdtrack->GetStatus() & AliESDtrack::kTPCrefit) continue;
429         if (!(esdtrack->GetStatus() & AliESDtrack::kITSrefit)) continue;
430         if (esdtrack->GetStatus() & AliESDtrack::kTPCin) continue;
431         UChar_t iclm = esdtrack->GetITSClusterMap();
432         Int_t incls = 0;
433         for (int iter=0; iter<6; iter++) if (iclm&(1<<iter)) incls++;
434         if (incls<=3) {
435           if (Debug()>1) cout << "Rejecting track with " << incls << " clusters" << endl;
436           continue;
437         }
438       }
439
440       AliFemtoTrack* trackCopy = new AliFemtoTrack();   
441       trackCopy->SetCharge((short)esdtrack->GetSign());
442
443       //in aliroot we have AliPID 
444       //0-electron 1-muon 2-pion 3-kaon 4-proton 5-photon 6-pi0 7-neutron 8-kaon0 9-eleCon   
445       //we use only 5 first
446       double esdpid[5];
447       //       esdtrack->GetESDpid(esdpid);
448       esdtrack->GetTPCpid(esdpid);
449       trackCopy->SetPidProbElectron(esdpid[0]);
450       trackCopy->SetPidProbMuon(esdpid[1]);
451       trackCopy->SetPidProbPion(esdpid[2]);
452       trackCopy->SetPidProbKaon(esdpid[3]);
453       trackCopy->SetPidProbProton(esdpid[4]);
454       
455       esdpid[0] = -100000.0;
456       esdpid[1] = -100000.0;
457       esdpid[2] = -100000.0;
458       esdpid[3] = -100000.0;
459       esdpid[4] = -100000.0;
460       
461       double tTOF = 0.0;
462
463       if (esdtrack->GetStatus()&AliESDtrack::kTOFpid) {
464         tTOF = esdtrack->GetTOFsignal();
465         esdtrack->GetIntegratedTimes(esdpid);
466       }
467
468       trackCopy->SetTofExpectedTimes(tTOF-esdpid[2], tTOF-esdpid[3], tTOF-esdpid[4]);
469
470       //////  TPC ////////////////////////////////////////////
471
472       float nsigmaTPCK=-1000.;                                                        
473       float nsigmaTPCPi=-1000.;                                                        
474       float nsigmaTPCP=-1000.;                                                        
475           
476   
477       if ((fESDpid) && (esdtrack->IsOn(AliESDtrack::kTPCpid))){
478         nsigmaTPCK = fESDpid->NumberOfSigmasTPC(esdtrack,AliPID::kKaon);
479         nsigmaTPCPi = fESDpid->NumberOfSigmasTPC(esdtrack,AliPID::kPion);
480         nsigmaTPCP = fESDpid->NumberOfSigmasTPC(esdtrack,AliPID::kProton);
481
482       }
483       trackCopy->SetNSigmaTPCPi(nsigmaTPCPi);
484       trackCopy->SetNSigmaTPCK(nsigmaTPCK);
485       trackCopy->SetNSigmaTPCP(nsigmaTPCP);
486
487       ///// TOF ///////////////////////////////////////////////
488
489         float vp=-1000.;
490         float nsigmaTOFPi=-1000.;
491         float nsigmaTOFK=-1000.;
492         float nsigmaTOFP=-1000.;
493
494         if ((esdtrack->GetStatus()&AliESDtrack::kTOFpid) &&
495             (esdtrack->GetStatus()&AliESDtrack::kTOFout) &&
496             (esdtrack->GetStatus()&AliESDtrack::kTIME) &&
497             !(esdtrack->GetStatus()&AliESDtrack::kTOFmismatch))
498           {
499
500             //if ((esdtrack->GetStatus()&AliESDtrack::kTOFpid) &&
501             //(esdtrack->GetStatus()&AliESDtrack::kTOFout) &&
502             //(esdtrack->GetStatus()&AliESDtrack::kTIME)){
503             // collect info from ESDpid class
504           
505             if ((fESDpid) && (esdtrack->IsOn(AliESDtrack::kTOFpid))) {
506               
507               
508               double tZero = fESDpid->GetTOFResponse().GetStartTime(esdtrack->P());
509               
510               nsigmaTOFPi = fESDpid->NumberOfSigmasTOF(esdtrack,AliPID::kPion,tZero);
511               nsigmaTOFK = fESDpid->NumberOfSigmasTOF(esdtrack,AliPID::kKaon,tZero);
512               nsigmaTOFP = fESDpid->NumberOfSigmasTOF(esdtrack,AliPID::kProton,tZero);
513               
514               Double_t len=esdtrack->GetIntegratedLength();
515               Double_t tof=esdtrack->GetTOFsignal();
516               if(tof > 0.) vp=len/tof/0.03;
517             }
518           }
519         
520         trackCopy->SetVTOF(vp);
521         trackCopy->SetNSigmaTOFPi(nsigmaTOFPi);
522         trackCopy->SetNSigmaTOFK(nsigmaTOFK);
523         trackCopy->SetNSigmaTOFP(nsigmaTOFP);
524         
525         double pxyz[3];
526         double rxyz[3];
527         double impact[2];
528         double covimpact[3];
529       
530         if (fUseTPCOnly) {
531           if (!esdtrack->GetTPCInnerParam()) {
532             delete trackCopy;
533             continue;
534           }
535
536
537           AliExternalTrackParam *param = new AliExternalTrackParam(*esdtrack->GetTPCInnerParam());
538           param->GetXYZ(rxyz);
539           param->PropagateToDCA(fEvent->GetPrimaryVertexTPC(), (fEvent->GetMagneticField()), 10000, impact, covimpact);
540           param->GetPxPyPz(pxyz);//reading noconstarined momentum
541
542           if (fReadInner == true) {
543             AliFemtoModelHiddenInfo *tInfo = new AliFemtoModelHiddenInfo();
544             tInfo->SetPDGPid(211);
545             tInfo->SetTrueMomentum(pxyz[0], pxyz[1], pxyz[2]);
546             tInfo->SetMass(0.13957);
547             //    tInfo->SetEmissionPoint(rxyz[0], rxyz[1], rxyz[2], 0.0);
548             //    tInfo->SetEmissionPoint(fV1[0], fV1[1], fV1[2], 0.0);
549             tInfo->SetEmissionPoint(rxyz[0]-fV1[0], rxyz[1]-fV1[1], rxyz[2]-fV1[2], 0.0);
550             trackCopy->SetHiddenInfo(tInfo);
551           }
552         
553           AliFemtoThreeVector v(pxyz[0],pxyz[1],pxyz[2]);
554           if (v.Mag() < 0.0001) {
555             //  cout << "Found 0 momentum ???? " <<endl;
556             delete trackCopy;
557             continue;
558           }
559           trackCopy->SetP(v);//setting momentum
560           trackCopy->SetPt(sqrt(pxyz[0]*pxyz[0]+pxyz[1]*pxyz[1]));
561
562           const AliFmThreeVectorD kP(pxyz[0],pxyz[1],pxyz[2]);
563           const AliFmThreeVectorD kOrigin(fV1[0],fV1[1],fV1[2]);
564           //setting helix I do not if it is ok
565           AliFmPhysicalHelixD helix(kP,kOrigin,(double)(fEvent->GetMagneticField())*kilogauss,(double)(trackCopy->Charge())); 
566           trackCopy->SetHelix(helix);
567
568           //some stuff which could be useful 
569           trackCopy->SetImpactD(impact[0]);
570           trackCopy->SetImpactZ(impact[1]);
571           trackCopy->SetCdd(covimpact[0]);
572           trackCopy->SetCdz(covimpact[1]);
573           trackCopy->SetCzz(covimpact[2]);
574           trackCopy->SetSigmaToVertex(GetSigmaToVertex(impact, covimpact));     
575
576           delete param;
577         }
578         else {
579           if (fReadInner == true) {
580           
581             if (esdtrack->GetTPCInnerParam()) {
582               AliExternalTrackParam *param = new AliExternalTrackParam(*esdtrack->GetTPCInnerParam());
583               param->GetXYZ(rxyz);
584               //            param->PropagateToDCA(fEvent->GetPrimaryVertex(), (fEvent->GetMagneticField()), 10000);
585               param->GetPxPyPz(pxyz);//reading noconstarined momentum
586               delete param;
587             
588               AliFemtoModelHiddenInfo *tInfo = new AliFemtoModelHiddenInfo();
589               tInfo->SetPDGPid(211);
590               tInfo->SetTrueMomentum(pxyz[0], pxyz[1], pxyz[2]);
591               tInfo->SetMass(0.13957);
592               //            tInfo->SetEmissionPoint(rxyz[0], rxyz[1], rxyz[2], 0.0);
593               //tInfo->SetEmissionPoint(fV1[0], fV1[1], fV1[2], 0.0);
594               tInfo->SetEmissionPoint(rxyz[0]-fV1[0], rxyz[1]-fV1[1], rxyz[2]-fV1[2], 0.0);
595               trackCopy->SetHiddenInfo(tInfo);
596             }
597           }
598
599           if (fTrackType == kGlobal) {
600             if (fConstrained==true)                 
601               tGoodMomentum=esdtrack->GetConstrainedPxPyPz(pxyz); //reading constrained momentum
602             else
603               tGoodMomentum=esdtrack->GetPxPyPz(pxyz);//reading noconstarined momentum
604           }
605           else if (fTrackType == kTPCOnly) {
606             if (esdtrack->GetTPCInnerParam())
607               esdtrack->GetTPCInnerParam()->GetPxPyPz(pxyz);
608             else {
609               delete trackCopy;
610               continue;
611             }
612           }
613           else if (fTrackType == kITSOnly) {
614             if (fConstrained==true)                 
615               tGoodMomentum=esdtrack->GetConstrainedPxPyPz(pxyz); //reading constrained momentum
616             else
617               tGoodMomentum=esdtrack->GetPxPyPz(pxyz);//reading noconstarined momentum
618           }
619
620
621           AliFemtoThreeVector v(pxyz[0],pxyz[1],pxyz[2]);
622           if (v.Mag() < 0.0001) {
623             //  cout << "Found 0 momentum ???? " <<endl;
624             delete trackCopy;
625             continue;
626           }
627           trackCopy->SetP(v);//setting momentum
628           trackCopy->SetPt(sqrt(pxyz[0]*pxyz[0]+pxyz[1]*pxyz[1]));
629           const AliFmThreeVectorD kP(pxyz[0],pxyz[1],pxyz[2]);
630           const AliFmThreeVectorD kOrigin(fV1[0],fV1[1],fV1[2]);
631           //setting helix I do not if it is ok
632           AliFmPhysicalHelixD helix(kP,kOrigin,(double)(fEvent->GetMagneticField())*kilogauss,(double)(trackCopy->Charge())); 
633           trackCopy->SetHelix(helix);
634
635           //some stuff which could be useful 
636           float imp[2];
637           float cim[3];
638           // if (fTrackType == kTPCOnly) {
639           //   esdtrack->GetTPCInnerParam()->GetImpactParameters(imp,cim);
640           // }
641           // else {
642           esdtrack->GetImpactParameters(imp,cim);
643           // }
644
645           impact[0] = imp[0];
646           impact[1] = imp[1];
647           covimpact[0] = cim[0];
648           covimpact[1] = cim[1];
649           covimpact[2] = cim[2];
650
651           trackCopy->SetImpactD(impact[0]);
652           trackCopy->SetImpactZ(impact[1]);
653           trackCopy->SetCdd(covimpact[0]);
654           trackCopy->SetCdz(covimpact[1]);
655           trackCopy->SetCzz(covimpact[2]);
656           trackCopy->SetSigmaToVertex(GetSigmaToVertex(impact,covimpact));
657         }
658
659         trackCopy->SetTrackId(esdtrack->GetID());
660         trackCopy->SetFlags(esdtrack->GetStatus());
661         trackCopy->SetLabel(esdtrack->GetLabel());
662                 
663         trackCopy->SetITSchi2(esdtrack->GetITSchi2());    
664         if (esdtrack->GetITSFakeFlag())
665           trackCopy->SetITSncls(-esdtrack->GetNcls(0));     
666         else
667           trackCopy->SetITSncls(esdtrack->GetNcls(0));     
668         trackCopy->SetTPCchi2(esdtrack->GetTPCchi2());       
669         trackCopy->SetTPCncls(esdtrack->GetTPCNcls());       
670         trackCopy->SetTPCnclsF(esdtrack->GetTPCNclsF());      
671         trackCopy->SetTPCsignal(esdtrack->GetTPCsignal());
672         trackCopy->SetTPCsignalN((short)esdtrack->GetTPCsignalN()); //due to bug in aliesdtrack class   
673         trackCopy->SetTPCsignalS(esdtrack->GetTPCsignalSigma()); 
674
675         trackCopy->SetTPCClusterMap(esdtrack->GetTPCClusterMap());
676         trackCopy->SetTPCSharedMap(esdtrack->GetTPCSharedMap());
677
678         double xtpc[3];
679         esdtrack->GetInnerXYZ(xtpc);
680         xtpc[2] -= fV1[2];
681         trackCopy->SetNominalTPCEntrancePoint(xtpc);
682
683         esdtrack->GetOuterXYZ(xtpc);
684         xtpc[2] -= fV1[2];
685         trackCopy->SetNominalTPCExitPoint(xtpc);
686
687         int indexes[3];
688         for (int ik=0; ik<3; ik++) {
689           indexes[ik] = esdtrack->GetKinkIndex(ik);
690         }
691         trackCopy->SetKinkIndexes(indexes);
692         //decision if we want this track
693         //if we using diffrent labels we want that this label was use for first time 
694         //if we use hidden info we want to have match between sim data and ESD
695         if (tGoodMomentum==true)
696           {
697             hbtEvent->TrackCollection()->push_back(trackCopy);//adding track to analysis
698             realnofTracks++;//real number of tracks
699             //    delete trackCopy;
700           }
701         else
702           {
703             delete  trackCopy;
704           }
705                 
706     }
707
708   hbtEvent->SetNumberOfTracks(realnofTracks);//setting number of track which we read in event   
709
710   AliCentrality *cent = fEvent->GetCentrality();
711   if (cent) {
712     hbtEvent->SetCentralityV0(cent->GetCentralityPercentile("V0M"));
713     //    hbtEvent->SetCentralityFMD(cent->GetCentralityPercentile("FMD"));
714     hbtEvent->SetCentralitySPD1(cent->GetCentralityPercentile("CL1"));
715     //    hbtEvent->SetCentralityTrk(cent->GetCentralityPercentile("TRK"));
716
717     if (Debug()>1) printf("  FemtoReader Got Event with %f %f %f %f\n", cent->GetCentralityPercentile("V0M"), 0.0, cent->GetCentralityPercentile("CL1"), 0.0);
718   }
719
720   if (fEstEventMult == kGlobalCount) 
721     hbtEvent->SetNormalizedMult(tNormMult);
722   else if (fEstEventMult == kTracklet)
723     hbtEvent->SetNormalizedMult(tTracklet);
724   else if (fEstEventMult == kITSTPC)
725     hbtEvent->SetNormalizedMult(tITSTPC);
726   else if (fEstEventMult == kITSPure)
727     hbtEvent->SetNormalizedMult(tITSPure);
728   else if (fEstEventMult == kSPDLayer1)
729     hbtEvent->SetNormalizedMult(fEvent->GetMultiplicity()->GetNumberOfITSClusters(1));
730   else if (fEstEventMult == kV0Centrality) {
731     // centrality between 0 (central) and 1 (very peripheral)
732
733     if (cent) {
734       if (cent->GetCentralityPercentile("V0M") < 0.00001)
735         hbtEvent->SetNormalizedMult(-1);
736       else
737         hbtEvent->SetNormalizedMult(lrint(10.0*cent->GetCentralityPercentile("V0M")));
738       if (Debug()>1) printf ("Set Centrality %i %f %li\n", hbtEvent->UncorrectedNumberOfPrimaries(), 
739                              10.0*cent->GetCentralityPercentile("V0M"), lrint(10.0*cent->GetCentralityPercentile("V0M")));
740     }
741   }
742
743   if (tNormMultPos > tNormMultNeg)
744     hbtEvent->SetZDCParticipants(tNormMultPos);
745   else
746     hbtEvent->SetZDCParticipants(tNormMultNeg);
747   
748   AliEventplane* ep = fEvent->GetEventplane();
749   if (ep) {
750     hbtEvent->SetEP(ep);
751     hbtEvent->SetReactionPlaneAngle(ep->GetEventplane("Q"));
752   }
753
754   //V0
755   if(fReadV0)
756     {
757       for (Int_t i = 0; i < fEvent->GetNumberOfV0s(); i++) {
758         AliESDv0* esdv0 = fEvent->GetV0(i);
759         if (!esdv0) continue;
760         //if(esdv0->GetNDaughters()>2) continue;
761         //if(esdv0->GetNProngs()>2) continue;
762         if(esdv0->Charge()!=0) continue;
763         AliESDtrack *trackPos = fEvent->GetTrack(esdv0->GetPindex());
764         if(!trackPos) continue;
765         AliESDtrack *trackNeg = fEvent->GetTrack(esdv0->GetNindex());
766         if(!trackNeg) continue;
767         if(trackPos->Charge()==trackNeg->Charge()) continue;
768         AliFemtoV0* trackCopyV0 = new AliFemtoV0();
769         CopyESDtoFemtoV0(esdv0, trackCopyV0, fEvent);
770         hbtEvent->V0Collection()->push_back(trackCopyV0);
771         //cout<<"Pushback v0 to v0collection"<<endl;
772       }
773     }
774
775
776   fCurEvent++;  
777   if (Debug()>1) cout<<"end of reading nt "<<nofTracks<<" real number "<<realnofTracks<<endl;
778
779   return hbtEvent; 
780 }
781 //___________________
782 void AliFemtoEventReaderESDChain::SetESDSource(AliESDEvent *aESD)
783 {
784   // The chain loads the ESD for us
785   // You must provide the address where it can be found
786   fEvent = aESD;
787 }
788 //___________________
789 // void AliFemtoEventReaderESDChain::SetESDfriendSource(AliESDfriend *aFriend)
790 // {
791 //   // We need the ESD tree to obtain
792 //   // information about the friend file location
793 //   fEventFriend = aFriend;
794 // }
795
796 //____________________________________________________________________
797 Float_t AliFemtoEventReaderESDChain::GetSigmaToVertex(double *impact, double *covar)
798 {
799   // Calculates the number of sigma to the vertex.
800
801   Float_t b[2];
802   Float_t bRes[2];
803   Float_t bCov[3];
804
805   b[0] = impact[0];
806   b[1] = impact[1];
807   bCov[0] = covar[0];
808   bCov[1] = covar[1];
809   bCov[2] = covar[2];
810
811   bRes[0] = TMath::Sqrt(bCov[0]);
812   bRes[1] = TMath::Sqrt(bCov[2]);
813
814   // -----------------------------------
815   // How to get to a n-sigma cut?
816   //
817   // The accumulated statistics from 0 to d is
818   //
819   // ->  Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
820   // ->  1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
821   //
822   // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
823   // Can this be expressed in a different way?
824
825   if (bRes[0] == 0 || bRes[1] ==0)
826     return -1;
827
828   Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
829
830   // stupid rounding problem screws up everything:
831   // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
832   if (TMath::Exp(-d * d / 2) < 1e-10)
833     return 1000;
834
835   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
836   return d;
837 }
838
839
840
841
842 void AliFemtoEventReaderESDChain::CopyESDtoFemtoV0(AliESDv0 *tESDv0, AliFemtoV0 *tFemtoV0, AliESDEvent *tESDevent)
843 {
844   AliESDtrack *pTrack = tESDevent->GetTrack(tESDv0->GetPindex());
845   AliESDtrack *nTrack = tESDevent->GetTrack(tESDv0->GetNindex());
846
847   //tFemtoV0->SetdecayLengthV0(tESDv0->DecayLengthV0()); //wrocic do tego
848   //tFemtoV0->SetdecayVertexV0X(tESDv0->DecayVertexV0X());
849   //tFemtoV0->SetdecayVertexV0Y(tESDv0->DecayVertexV0Y());
850   //tFemtoV0->SetdecayVertexV0Z(tESDv0->DecayVertexV0Z()); //nie ma w AliESDv0
851   //AliFemtoThreeVector decayvertex(tESDv0->DecayVertexV0X(),tESDv0->DecayVertexV0Y(),tESDv0->DecayVertexV0Z());
852   //tFemtoV0->SetdecayVertexV0(decayvertex);
853   tFemtoV0->SetdcaV0Daughters(tESDv0->GetDcaV0Daughters());
854   double fPrimaryVtxPosition[3];
855   tESDevent->GetPrimaryVertex()->GetXYZ(fPrimaryVtxPosition);
856   tFemtoV0->SetdcaV0ToPrimVertex(tESDv0->GetD(fPrimaryVtxPosition[0],fPrimaryVtxPosition[1],fPrimaryVtxPosition[2]));
857   tFemtoV0->SetdcaPosToPrimVertex(TMath::Abs(pTrack->GetD(fPrimaryVtxPosition[0],fPrimaryVtxPosition[1],tESDevent->GetMagneticField())));
858   tFemtoV0->SetdcaNegToPrimVertex(TMath::Abs(pTrack->GetD(fPrimaryVtxPosition[0],fPrimaryVtxPosition[1],tESDevent->GetMagneticField())));
859   double MomPos[3];
860   pTrack->PxPyPz(MomPos);
861   tFemtoV0->SetmomPosX(MomPos[0]);
862   tFemtoV0->SetmomPosY(MomPos[1]);
863   tFemtoV0->SetmomPosZ(MomPos[2]);
864   AliFemtoThreeVector mompos(MomPos[0],MomPos[1],MomPos[2]);
865   tFemtoV0->SetmomPos(mompos);
866   double MomNeg[3];
867   nTrack->PxPyPz(MomNeg);
868   tFemtoV0->SetmomNegX(MomNeg[0]);
869   tFemtoV0->SetmomNegY(MomNeg[1]);
870   tFemtoV0->SetmomNegZ(MomNeg[2]);
871   AliFemtoThreeVector momneg(MomNeg[0],MomNeg[1],MomNeg[2]);
872   tFemtoV0->SetmomNeg(momneg);
873
874   //jest cos takiego w AliFemtoV0.h czego nie ma w AliAODv0.h
875   //void SettpcHitsPos(const int& i);      
876   //void SettpcHitsNeg(const int& i);      
877
878   //void SetTrackTopologyMapPos(unsigned int word, const unsigned long& m);
879   //void SetTrackTopologyMapNeg(unsigned int word, const unsigned long& m);
880
881   tFemtoV0->SetmomV0X(tESDv0->Px());
882   tFemtoV0->SetmomV0Y(tESDv0->Py());
883   tFemtoV0->SetmomV0Z(tESDv0->Pz());
884   AliFemtoThreeVector momv0(tESDv0->Px(),tESDv0->Py(),tESDv0->Pz());
885   tFemtoV0->SetmomV0(momv0);
886   tFemtoV0->SetalphaV0(tESDv0->AlphaV0());
887   tFemtoV0->SetptArmV0(tESDv0->PtArmV0());
888   //tFemtoV0->SeteLambda(tESDv0->ELambda());
889   //tFemtoV0->SeteK0Short(tESDv0->EK0Short());
890   //tFemtoV0->SetePosProton(tESDv0->EPosProton());
891   //tFemtoV0->SeteNegProton(tESDv0->ENegProton());
892   tFemtoV0->SetmassLambda(tESDv0->GetEffMass(4,2));
893   tFemtoV0->SetmassAntiLambda(tESDv0->GetEffMass(2,4));
894   tFemtoV0->SetmassK0Short(tESDv0->GetEffMass(2,2));
895   //tFemtoV0->SetrapLambda(tESDv0->RapLambda());
896   //tFemtoV0->SetrapK0Short(tESDv0->RapK0Short());
897   
898   //void SetcTauLambda( float x);   
899   //void SetcTauK0Short( float x); 
900   
901   tFemtoV0->SetptV0(tESDv0->Pt());
902   tFemtoV0->SetptotV0(tESDv0->P());
903   tFemtoV0->SetptPos(pTrack->Pt());
904   tFemtoV0->SetptotPos(pTrack->P());
905   tFemtoV0->SetptNeg(nTrack->Pt());
906   tFemtoV0->SetptotNeg(nTrack->P());
907   
908   tFemtoV0->SetidNeg(nTrack->GetID());
909   //cout<<"tESDv0->GetNegID(): "<<tESDv0->GetNegID()<<endl;
910   //cout<<"tFemtoV0->IdNeg(): "<<tFemtoV0->IdNeg()<<endl;
911   tFemtoV0->SetidPos(pTrack->GetID());
912
913   tFemtoV0->SetEtaV0(tESDv0->Eta());
914   tFemtoV0->SetPhiV0(tESDv0->Phi());
915   tFemtoV0->SetCosPointingAngle(tESDv0->GetV0CosineOfPointingAngle(fPrimaryVtxPosition[0],fPrimaryVtxPosition[1], fPrimaryVtxPosition[2]));
916   //tFemtoV0->SetYV0(tESDv0->Y());
917
918   //void SetdedxNeg(float x);
919   //void SeterrdedxNeg(float x);//Gael 04Fev2002
920   //void SetlendedxNeg(float x);//Gael 04Fev2002
921   //void SetdedxPos(float x);
922   //void SeterrdedxPos(float x);//Gael 04Fev2002
923   //void SetlendedxPos(float x);//Gael 04Fev2002
924
925   tFemtoV0->SetEtaPos(pTrack->Eta());
926   tFemtoV0->SetEtaNeg(nTrack->Eta());
927
928   //AliAODTrack *trackpos = (AliAODTrack*)tESDv0->GetDaughter(0);
929   //AliAODTrack *trackneg = (AliAODTrack*)tESDv0->GetDaughter(1);
930
931   if(pTrack && nTrack)
932     {
933       tFemtoV0->SetEtaPos(pTrack->Eta()); //tESDv0->PseudoRapPos()
934       tFemtoV0->SetEtaNeg(nTrack->Eta()); //tESDv0->PseudoRapNeg()
935       tFemtoV0->SetTPCNclsPos(pTrack->GetTPCNcls());
936       tFemtoV0->SetTPCNclsNeg(nTrack->GetTPCNcls());
937       tFemtoV0->SetTPCclustersPos(pTrack->GetTPCClusterMap());
938       tFemtoV0->SetTPCclustersNeg(nTrack->GetTPCClusterMap());
939       tFemtoV0->SetTPCsharingPos(pTrack->GetTPCSharedMap());
940       tFemtoV0->SetTPCsharingNeg(nTrack->GetTPCSharedMap());
941       tFemtoV0->SetNdofPos(pTrack->GetTPCchi2()/pTrack->GetTPCNcls());
942       tFemtoV0->SetNdofNeg(nTrack->GetTPCchi2()/nTrack->GetTPCNcls());
943       tFemtoV0->SetStatusPos(pTrack->GetStatus());
944       tFemtoV0->SetStatusNeg(nTrack->GetStatus());
945
946       tFemtoV0->SetPosNSigmaTPCK(fESDpid->NumberOfSigmasTPC(pTrack,AliPID::kKaon));
947       tFemtoV0->SetNegNSigmaTPCK(fESDpid->NumberOfSigmasTPC(nTrack,AliPID::kKaon));
948       tFemtoV0->SetPosNSigmaTPCP(fESDpid->NumberOfSigmasTPC(pTrack,AliPID::kProton));
949       tFemtoV0->SetNegNSigmaTPCP(fESDpid->NumberOfSigmasTPC(nTrack,AliPID::kProton));
950       tFemtoV0->SetPosNSigmaTPCPi(fESDpid->NumberOfSigmasTPC(pTrack,AliPID::kPion));
951       tFemtoV0->SetNegNSigmaTPCPi(fESDpid->NumberOfSigmasTPC(nTrack,AliPID::kPion));
952
953
954       if((tFemtoV0->StatusPos()&AliESDtrack::kTOFpid)==0 || (tFemtoV0->StatusPos()&AliESDtrack::kTIME)==0 || (tFemtoV0->StatusPos()&AliESDtrack::kTOFout)==0 || (tFemtoV0->StatusPos()&AliESDtrack::kTOFmismatch)>0)
955         {
956           if((tFemtoV0->StatusNeg()&AliESDtrack::kTOFpid)==0 || (tFemtoV0->StatusNeg()&AliESDtrack::kTIME)==0 || (tFemtoV0->StatusNeg()&AliESDtrack::kTOFout)==0 || (tFemtoV0->StatusNeg()&AliESDtrack::kTOFmismatch)>0)
957             {
958               tFemtoV0->SetPosNSigmaTOFK(-9999);
959               tFemtoV0->SetNegNSigmaTOFK(-9999);
960               tFemtoV0->SetPosNSigmaTOFP(-9999);
961               tFemtoV0->SetNegNSigmaTOFP(-9999);
962               tFemtoV0->SetPosNSigmaTOFPi(-9999);
963               tFemtoV0->SetNegNSigmaTOFPi(-9999);
964             }
965         }
966       else
967         {
968           tFemtoV0->SetPosNSigmaTOFK(fESDpid->NumberOfSigmasTOF(pTrack,AliPID::kKaon));
969           tFemtoV0->SetNegNSigmaTOFK(fESDpid->NumberOfSigmasTOF(nTrack,AliPID::kKaon));
970           tFemtoV0->SetPosNSigmaTOFP(fESDpid->NumberOfSigmasTOF(pTrack,AliPID::kProton));
971           tFemtoV0->SetNegNSigmaTOFP(fESDpid->NumberOfSigmasTOF(nTrack,AliPID::kProton));
972           tFemtoV0->SetPosNSigmaTOFPi(fESDpid->NumberOfSigmasTOF(pTrack,AliPID::kPion));
973           tFemtoV0->SetNegNSigmaTOFPi(fESDpid->NumberOfSigmasTOF(nTrack,AliPID::kPion));
974         }
975     }
976   else
977     {
978       tFemtoV0->SetStatusPos(999);
979       tFemtoV0->SetStatusNeg(999);
980     }
981   tFemtoV0->SetOnFlyStatusV0(tESDv0->GetOnFlyStatus());
982 }
983
984
985
986
987 void AliFemtoEventReaderESDChain::SetReadTrackType(ReadTrackType aType)
988 {
989   fTrackType = aType;
990 }
991
992 //trigger
993 void AliFemtoEventReaderESDChain::SetEventTrigger(UInt_t eventtrig)
994 {
995   fEventTrig = eventtrig;
996 }
997
998 //V0 reading
999 void AliFemtoEventReaderESDChain::SetReadV0(bool a)
1000 {
1001   fReadV0 = a;
1002 }
1003
1004
1005
1006
1007