]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemto/AliFemtoEventReaderAOD.cxx
Reaction Plane analysis updates
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemto / AliFemtoEventReaderAOD.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoEventReaderAOD - the reader class for the Alice AOD                //
4 // Reads in AOD information and converts it into internal AliFemtoEvent       //
5 // Authors: Marek Chojnacki mchojnacki@knf.pw.edu.pl                          //
6 //          Adam Kisiel kisiel@mps.ohio-state.edu                             //
7 //                                                                            //
8 ////////////////////////////////////////////////////////////////////////////////
9
10 #include "AliFemtoEventReaderAOD.h"
11
12 #include "TFile.h"
13 #include "TTree.h"
14 #include "AliAODEvent.h"
15 #include "AliAODTrack.h"
16 #include "AliAODVertex.h"
17 #include "AliAODMCHeader.h"
18 #include "AliESDtrack.h"
19
20 #include "AliFmPhysicalHelixD.h"
21 #include "AliFmThreeVectorF.h"
22
23 #include "SystemOfUnits.h"
24
25 #include "AliFemtoEvent.h"
26 #include "AliFemtoModelHiddenInfo.h"
27 #include "AliFemtoModelGlobalHiddenInfo.h"
28 #include "AliPID.h"
29
30 #include "AliAODpidUtil.h"
31
32 ClassImp(AliFemtoEventReaderAOD)
33
34 #if !(ST_NO_NAMESPACES)
35   using namespace units;
36 #endif
37
38 using namespace std;
39
40 double fV1[3];
41
42 //____________________________
43 //constructor with 0 parameters , look at default settings 
44 AliFemtoEventReaderAOD::AliFemtoEventReaderAOD():
45   fNumberofEvent(0),
46   fCurEvent(0),
47   fEvent(0x0),
48   fAllTrue(160),
49   fAllFalse(160),
50   fFilterBit(0),
51   //  fPWG2AODTracks(0x0),
52   fReadMC(0),
53   fReadV0(0),
54   fUsePreCent(0),
55   fEstEventMult(kCentrality),
56   fAODpidUtil(0),
57   fAODheader(0),
58   fInputFile(" "),
59   fFileName(" "),
60   fTree(0x0),
61   fAodFile(0x0),
62     fMagFieldSign(1),
63     fisEPVZ(kTRUE)
64 {
65   // default constructor
66   fAllTrue.ResetAllBits(kTRUE);
67   fAllFalse.ResetAllBits(kFALSE);
68   fCentRange[0] = 0;
69   fCentRange[1] = 1000;
70 }
71
72 AliFemtoEventReaderAOD::AliFemtoEventReaderAOD(const AliFemtoEventReaderAOD &aReader) :
73   AliFemtoEventReader(),
74   fNumberofEvent(0),
75   fCurEvent(0),
76   fEvent(0x0),
77   fAllTrue(160),
78   fAllFalse(160),
79   fFilterBit(0),
80   //  fPWG2AODTracks(0x0),
81   fReadMC(0),
82   fReadV0(0),
83   fUsePreCent(0),
84   fEstEventMult(kCentrality),
85   fAODpidUtil(0),
86   fAODheader(0),
87   fInputFile(" "),
88   fFileName(" "),
89   fTree(0x0),
90   fAodFile(0x0),
91     fMagFieldSign(1),
92     fisEPVZ(kTRUE)
93 {
94   // copy constructor
95   fReadMC = aReader.fReadMC;
96   fReadV0 = aReader.fReadV0;
97   fInputFile = aReader.fInputFile;
98   fFileName  = aReader.fFileName;
99   fNumberofEvent = aReader.fNumberofEvent;
100   fCurEvent = aReader.fCurEvent;
101   fEvent = new AliAODEvent();
102   fAodFile = new TFile(aReader.fAodFile->GetName());
103   fAllTrue.ResetAllBits(kTRUE);
104   fAllFalse.ResetAllBits(kFALSE);
105   fFilterBit = aReader.fFilterBit;
106   //  fPWG2AODTracks = aReader.fPWG2AODTracks;
107   fAODpidUtil = aReader.fAODpidUtil;
108   fAODheader = aReader.fAODheader;
109   fCentRange[0] = aReader.fCentRange[0];
110   fCentRange[1] = aReader.fCentRange[1];
111   fEstEventMult = aReader.fEstEventMult;
112   fUsePreCent = aReader.fUsePreCent;
113 }
114 //__________________
115 //Destructor
116 AliFemtoEventReaderAOD::~AliFemtoEventReaderAOD()
117 {
118   // destructor
119   delete fTree;
120   delete fEvent;
121   delete fAodFile;
122 //   if (fPWG2AODTracks) {
123 //     fPWG2AODTracks->Delete();
124 //     delete fPWG2AODTracks;
125 //   }
126 }
127
128 //__________________
129 AliFemtoEventReaderAOD& AliFemtoEventReaderAOD::operator=(const AliFemtoEventReaderAOD& aReader)
130 {
131   // assignment operator
132   if (this == &aReader)
133    return *this;
134
135   fInputFile = aReader.fInputFile;
136   fFileName  = aReader.fFileName;
137   fNumberofEvent = aReader.fNumberofEvent;
138   fCurEvent = aReader.fCurEvent;
139   if (fTree) delete fTree;
140   if (fEvent) delete fEvent;
141   fEvent = new AliAODEvent();
142   if (fAodFile) delete fAodFile;
143   fAodFile = new TFile(aReader.fAodFile->GetName());
144   fAllTrue.ResetAllBits(kTRUE);
145   fAllFalse.ResetAllBits(kFALSE);
146   fFilterBit = aReader.fFilterBit;
147   //  fPWG2AODTracks = aReader.fPWG2AODTracks;
148   fAODpidUtil = aReader.fAODpidUtil;
149   fAODheader = aReader.fAODheader;
150   fCentRange[0] = aReader.fCentRange[0];
151   fCentRange[1] = aReader.fCentRange[1];
152   fUsePreCent = aReader.fUsePreCent;
153   fEstEventMult = aReader.fEstEventMult;
154
155   return *this;
156 }
157 //__________________
158 AliFemtoString AliFemtoEventReaderAOD::Report()
159 {
160   // create reader report
161   AliFemtoString temp = "\n This is the AliFemtoEventReaderAOD\n";
162   return temp;
163 }
164
165 //__________________
166 void AliFemtoEventReaderAOD::SetInputFile(const char* inputFile)
167 {
168   //setting the name of file where names of AOD file are written 
169   //it takes only this files which have good trees
170   char buffer[256];
171   fInputFile=string(inputFile);
172   ifstream infile(inputFile);
173
174   fTree = new TChain("aodTree");
175
176   if(infile.good()==true)
177     { 
178       //checking if all give files have good tree inside
179       while (infile.eof()==false)
180         {
181           infile.getline(buffer,256);
182           TFile *aodFile=TFile::Open(buffer,"READ");
183           if (aodFile!=0x0)
184             {   
185               TTree* tree = (TTree*) aodFile->Get("aodTree");
186               if (tree!=0x0)
187                 {
188                   //              cout<<"putting file  "<<string(buffer)<<" into analysis"<<endl;
189                   fTree->AddFile(buffer);
190                   delete tree;
191                 }
192               aodFile->Close(); 
193             }
194           delete aodFile;
195         }
196     }
197 }
198
199 AliFemtoEvent* AliFemtoEventReaderAOD::ReturnHbtEvent()
200 {
201   // read in a next hbt event from the chain
202   // convert it to AliFemtoEvent and return
203   // for further analysis
204   AliFemtoEvent *hbtEvent = 0;
205     // cout<<"reader"<<endl;
206   if (fCurEvent==fNumberofEvent)//open next file  
207     {
208       if(fNumberofEvent==0)     
209         {
210           fEvent=new AliAODEvent();
211           fEvent->ReadFromTree(fTree);
212
213           // Check for the existence of the additional information
214 //        fPWG2AODTracks = (TClonesArray *) fEvent->GetList()->FindObject("pwg2aodtracks");
215
216 //        if (fPWG2AODTracks) {
217 //          cout << "Found additional PWG2 specific information in the AOD!" << endl;
218 //          cout << "Reading only tracks with the additional information" << endl;
219 //        }
220
221           fNumberofEvent=fTree->GetEntries();
222           //      cout<<"Number of Entries in file "<<fNumberofEvent<<endl;
223           fCurEvent=0;
224         }
225       else //no more data to read
226         {
227             // cout<<"no more files "<<hbtEvent<<endl;
228           fReaderStatus=1;
229           return hbtEvent; 
230         }
231     }           
232
233     // cout<<"starting to read event "<<fCurEvent<<endl;
234   fTree->GetEvent(fCurEvent);//getting next event
235   //  cout << "Read event " << fEvent << " from file " << fTree << endl;
236         
237   hbtEvent = new AliFemtoEvent;
238
239   CopyAODtoFemtoEvent(hbtEvent);
240   fCurEvent++;
241
242
243   return hbtEvent; 
244 }
245
246 void AliFemtoEventReaderAOD::CopyAODtoFemtoEvent(AliFemtoEvent *tEvent)
247 {
248
249   // A function that reads in the AOD event
250   // and transfers the neccessary information into
251   // the internal AliFemtoEvent
252
253   // setting global event characteristics
254   tEvent->SetRunNumber(fEvent->GetRunNumber());
255   tEvent->SetMagneticField(fEvent->GetMagneticField()*kilogauss);//to check if here is ok
256   tEvent->SetZDCN1Energy(fEvent->GetZDCN1Energy());
257   tEvent->SetZDCP1Energy(fEvent->GetZDCP1Energy());
258   tEvent->SetZDCN2Energy(fEvent->GetZDCN2Energy());
259   tEvent->SetZDCP2Energy(fEvent->GetZDCP2Energy());
260   tEvent->SetZDCEMEnergy(fEvent->GetZDCEMEnergy(0));
261   tEvent->SetZDCParticipants(0);
262   tEvent->SetTriggerMask(fEvent->GetTriggerMask());
263   tEvent->SetTriggerCluster(fEvent->GetTriggerCluster());
264   
265   // Attempt to access MC header
266   AliAODMCHeader *mcH;
267   TClonesArray *mcP=0;
268   if (fReadMC) {
269     mcH = (AliAODMCHeader *) fEvent->FindListObject(AliAODMCHeader::StdBranchName());
270     if (!mcH) {
271       cout << "AOD MC information requested, but no header found!" << endl;
272     }
273
274     mcP = (TClonesArray *) fEvent->FindListObject(AliAODMCParticle::StdBranchName());
275     if (!mcP) {
276       cout << "AOD MC information requested, but no particle array found!" << endl;
277     }
278   }
279
280   tEvent->SetReactionPlaneAngle(fEvent->GetHeader()->GetQTheta(0)/2.0);
281
282   Int_t *motherids=0;
283   if (mcP) {
284     motherids = new Int_t[((AliAODMCParticle *) mcP->At(mcP->GetEntries()-1))->GetLabel()];
285     for (int ip=0; ip<mcP->GetEntries(); ip++) motherids[ip] = 0;
286
287     // Read in mother ids
288     AliAODMCParticle *motherpart;
289     for (int ip=0; ip<mcP->GetEntries(); ip++) {
290       motherpart = (AliAODMCParticle *) mcP->At(ip);
291       if (motherpart->GetDaughter(0) > 0)
292         motherids[motherpart->GetDaughter(0)] = ip;
293       if (motherpart->GetDaughter(1) > 0)
294         motherids[motherpart->GetDaughter(1)] = ip;
295     }
296   }
297
298   // Primary Vertex position
299   //  double fV1[3];
300   fEvent->GetPrimaryVertex()->GetPosition(fV1);
301
302   AliFmThreeVectorF vertex(fV1[0],fV1[1],fV1[2]);
303   tEvent->SetPrimVertPos(vertex);
304         
305   //starting to reading tracks
306   int nofTracks=0;  //number of reconstructed tracks in event
307
308   // Check to see whether the additional info exists
309   //   if (fPWG2AODTracks)
310   //     nofTracks=fPWG2AODTracks->GetEntries();
311   //   else
312   nofTracks=fEvent->GetNumberOfTracks();
313   
314   AliEventplane *ep = fEvent->GetEventplane();
315   if (ep) {
316     tEvent->SetEP(ep);
317         if (fisEPVZ)
318             tEvent->SetReactionPlaneAngle(ep->GetEventplane("V0",fEvent,2));
319         else
320     tEvent->SetReactionPlaneAngle(ep->GetEventplane("Q"));
321   }
322
323   AliCentrality *cent = fEvent->GetCentrality();
324
325   if (!fEstEventMult && cent && fUsePreCent) {
326     if ((cent->GetCentralityPercentile("V0M")*10 < fCentRange[0]) ||
327         (cent->GetCentralityPercentile("V0M")*10 > fCentRange[1]))
328       {
329             // cout << "Centrality " << cent->GetCentralityPercentile("V0M") << " outside of preselection range " << fCentRange[0] << " - " << fCentRange[1] << endl;
330         
331         return;
332       }
333   }
334
335   int realnofTracks=0;   // number of track which we use in a analysis
336   int tracksPrim=0;     
337
338   int labels[20000];    
339   for (int il=0; il<20000; il++) labels[il] = -1;
340
341   // looking for global tracks and saving their numbers to copy from them PID information to TPC-only tracks in the main loop over tracks
342   for (int i=0;i<nofTracks;i++) {
343     const AliAODTrack *aodtrack=fEvent->GetTrack(i);
344     if (!aodtrack->TestFilterBit(fFilterBit)) {
345       if(aodtrack->GetID() < 0) continue;
346       labels[aodtrack->GetID()] = i;
347     }
348   }
349
350   int tNormMult = 0;
351   for (int i=0;i<nofTracks;i++)
352     {
353       AliFemtoTrack* trackCopy = new AliFemtoTrack();   
354
355 //       if (fPWG2AODTracks) {
356 //      // Read tracks from the additional pwg2 specific AOD part
357 //      // if they exist
358 //      // Note that in that case all the AOD tracks without the 
359 //      // additional information will be ignored !
360 //      AliPWG2AODTrack *pwg2aodtrack = (AliPWG2AODTrack *) fPWG2AODTracks->At(i);
361
362 //      // Getting the AOD track through the ref of the additional info
363 //      AliAODTrack *aodtrack = pwg2aodtrack->GetRefAODTrack(); 
364 //      if (!aodtrack->TestFilterBit(fFilterBit)) {
365 //        delete trackCopy;
366 //        continue;
367 //      }
368
369  
370 //      if (aodtrack->IsOn(AliESDtrack::kTPCrefit))
371 //        if (aodtrack->Chi2perNDF() < 6.0) 
372 //          if (aodtrack->Eta() < 0.9)
373 //            tNormMult++;
374
375
376 //      CopyAODtoFemtoTrack(aodtrack, trackCopy, pwg2aodtrack);
377         
378 //      if (mcP) {
379 //        // Fill the hidden information with the simulated data
380 //        //      Int_t pLabel = aodtrack->GetLabel();
381 //        AliAODMCParticle *tPart = GetParticleWithLabel(mcP, (TMath::Abs(aodtrack->GetLabel())));
382
383 //        // Check the mother information
384           
385 //        // Using the new way of storing the freeze-out information
386 //        // Final state particle is stored twice on the stack
387 //        // one copy (mother) is stored with original freeze-out information
388 //        //   and is not tracked
389 //        // the other one (daughter) is stored with primary vertex position
390 //        //   and is tracked
391           
392 //        // Freeze-out coordinates
393 //        double fpx=0.0, fpy=0.0, fpz=0.0, fpt=0.0;
394 //        fpx = tPart->Xv() - fV1[0];
395 //        fpy = tPart->Yv() - fV1[1];
396 //        fpz = tPart->Zv() - fV1[2];
397 //        fpt = tPart->T();
398
399 //        AliFemtoModelGlobalHiddenInfo *tInfo = new AliFemtoModelGlobalHiddenInfo();
400 //        tInfo->SetGlobalEmissionPoint(fpx, fpy, fpz);
401
402 //        fpx *= 1e13;
403 //        fpy *= 1e13;
404 //        fpz *= 1e13;
405 //        fpt *= 1e13;
406           
407 //        //      cout << "Looking for mother ids " << endl;
408 //        if (motherids[TMath::Abs(aodtrack->GetLabel())]>0) {
409 //          //  cout << "Got mother id" << endl;
410 //          AliAODMCParticle *mother = GetParticleWithLabel(mcP, motherids[TMath::Abs(aodtrack->GetLabel())]);
411 //          // Check if this is the same particle stored twice on the stack
412 //          if ((mother->GetPdgCode() == tPart->GetPdgCode() || (mother->Px() == tPart->Px()))) {
413 //            // It is the same particle
414 //            // Read in the original freeze-out information
415 //            // and convert it from to [fm]
416               
417 //            // EPOS style 
418 //            //          fpx = mother->Xv()*1e13*0.197327;
419 //            //          fpy = mother->Yv()*1e13*0.197327;
420 //            //          fpz = mother->Zv()*1e13*0.197327;
421 //            //          fpt = mother->T() *1e13*0.197327*0.5;
422               
423               
424 //            // Therminator style 
425 //            fpx = mother->Xv()*1e13;
426 //            fpy = mother->Yv()*1e13;
427 //            fpz = mother->Zv()*1e13;
428 //            fpt = mother->T() *1e13*3e10;
429               
430 //          }
431 //        }
432           
433 //        //       if (fRotateToEventPlane) {
434 //        //    double tPhi = TMath::ATan2(fpy, fpx);
435 //        //    double tRad = TMath::Hypot(fpx, fpy);
436         
437 //        //    fpx = tRad*TMath::Cos(tPhi - tReactionPlane);
438 //        //    fpy = tRad*TMath::Sin(tPhi - tReactionPlane);
439 //        //       }
440
441 //        tInfo->SetPDGPid(tPart->GetPdgCode());
442
443 //        //      if (fRotateToEventPlane) {
444 //        //        double tPhi = TMath::ATan2(tPart->Py(), tPart->Px());
445 //        //        double tRad = TMath::Hypot(tPart->Px(), tPart->Py());
446             
447 //        //        tInfo->SetTrueMomentum(tRad*TMath::Cos(tPhi - tReactionPlane),
448 //        //                               tRad*TMath::Sin(tPhi - tReactionPlane),
449 //        //                               tPart->Pz());
450 //        //      }
451 //        //       else
452 //        tInfo->SetTrueMomentum(tPart->Px(), tPart->Py(), tPart->Pz());
453 //        Double_t mass2 = (tPart->E() *tPart->E() -
454 //                          tPart->Px()*tPart->Px() -
455 //                          tPart->Py()*tPart->Py() -
456 //                          tPart->Pz()*tPart->Pz());
457 //        if (mass2>0.0)
458 //          tInfo->SetMass(TMath::Sqrt(mass2));
459 //        else 
460 //          tInfo->SetMass(0.0);
461           
462 //        tInfo->SetEmissionPoint(fpx, fpy, fpz, fpt);
463 //        trackCopy->SetHiddenInfo(tInfo);
464
465 //      }
466
467 //      double pxyz[3];
468 //      aodtrack->PxPyPz(pxyz);//reading noconstarined momentum
469 //      const AliFmThreeVectorD ktP(pxyz[0],pxyz[1],pxyz[2]);
470 //      // Check the sanity of the tracks - reject zero momentum tracks
471 //      if (ktP.Mag() == 0) {
472 //        delete trackCopy;
473 //        continue;
474 //      }
475 //       }
476 //       else {
477         // No additional information exists
478         // Read in the normal AliAODTracks 
479
480         //      const AliAODTrack *aodtrack=fEvent->GetTrack(i); // getting the AODtrack directly
481         AliAODTrack *aodtrack=fEvent->GetTrack(i); // getting the AODtrack directly
482         
483         if (aodtrack->IsPrimaryCandidate()) tracksPrim++;
484         
485         if (!aodtrack->TestFilterBit(fFilterBit)) {
486           delete trackCopy;
487           continue;
488         }
489         
490         //counting particles to set multiplicity
491         double impact[2];
492         double covimpact[3];
493         if (aodtrack->PropagateToDCA(fEvent->GetPrimaryVertex(),fEvent->GetMagneticField(),10000,impact,covimpact)) {
494           if(impact[0]<0.2 && TMath::Abs(impact[1]+fV1[2])<2.0)
495             //if (aodtrack->IsPrimaryCandidate()) //? instead of kinks?
496               if (aodtrack->Chi2perNDF() < 4.0) 
497                 if (aodtrack->Pt() > 0.15 && aodtrack->Pt() < 20) 
498                   if (aodtrack->GetTPCNcls() > 70)
499                     if (aodtrack->Eta() < 0.8)
500                       tNormMult++;
501
502         } 
503
504         CopyAODtoFemtoTrack(aodtrack, trackCopy);
505
506         // copying PID information from the correspondent track
507         //      const AliAODTrack *aodtrackpid = fEvent->GetTrack(labels[-1-fEvent->GetTrack(i)->GetID()]);
508
509
510         AliAODTrack *aodtrackpid;
511         if(fFilterBit ==  (1 << (7))) //for TPC Only tracks we have to copy PID information from corresponding global tracks
512           aodtrackpid = fEvent->GetTrack(labels[-1-fEvent->GetTrack(i)->GetID()]);
513         else
514           aodtrackpid = fEvent->GetTrack(i);
515         CopyPIDtoFemtoTrack(aodtrackpid, trackCopy);
516         
517         if (mcP) {
518           // Fill the hidden information with the simulated data
519           //      Int_t pLabel = aodtrack->GetLabel();
520           AliAODMCParticle *tPart = GetParticleWithLabel(mcP, (TMath::Abs(aodtrack->GetLabel())));
521           
522           AliFemtoModelGlobalHiddenInfo *tInfo = new AliFemtoModelGlobalHiddenInfo();
523           double fpx=0.0, fpy=0.0, fpz=0.0, fpt=0.0;
524           if (!tPart) {
525             fpx = fV1[0];
526             fpy = fV1[1];
527             fpz = fV1[2];
528             tInfo->SetGlobalEmissionPoint(fpx, fpy, fpz);
529             tInfo->SetPDGPid(0);
530             tInfo->SetTrueMomentum(0.0, 0.0, 0.0);
531             tInfo->SetEmissionPoint(0.0, 0.0, 0.0, 0.0);
532             tInfo->SetMass(0);
533           }
534           else {
535             // Check the mother information
536           
537             // Using the new way of storing the freeze-out information
538             // Final state particle is stored twice on the stack
539             // one copy (mother) is stored with original freeze-out information
540             //   and is not tracked
541             // the other one (daughter) is stored with primary vertex position
542             //   and is tracked
543             
544             // Freeze-out coordinates
545             fpx = tPart->Xv() - fV1[0];
546             fpy = tPart->Yv() - fV1[1];
547             fpz = tPart->Zv() - fV1[2];
548             //    fpt = tPart->T();
549             
550             tInfo->SetGlobalEmissionPoint(fpx, fpy, fpz);
551             
552             fpx *= 1e13;
553             fpy *= 1e13;
554             fpz *= 1e13;
555             //    fpt *= 1e13;
556             
557             //      cout << "Looking for mother ids " << endl;
558             if (motherids[TMath::Abs(aodtrack->GetLabel())]>0) {
559               //        cout << "Got mother id" << endl;
560               AliAODMCParticle *mother = GetParticleWithLabel(mcP, motherids[TMath::Abs(aodtrack->GetLabel())]);
561               // Check if this is the same particle stored twice on the stack
562               if (mother) {
563                 if ((mother->GetPdgCode() == tPart->GetPdgCode() || (mother->Px() == tPart->Px()))) {
564                   // It is the same particle
565                   // Read in the original freeze-out information
566                   // and convert it from to [fm]
567                   
568                   // EPOS style 
569                   //      fpx = mother->Xv()*1e13*0.197327;
570                   //      fpy = mother->Yv()*1e13*0.197327;
571                   //      fpz = mother->Zv()*1e13*0.197327;
572                   //      fpt = mother->T() *1e13*0.197327*0.5;
573                   
574                   
575                   // Therminator style 
576                   fpx = mother->Xv()*1e13;
577                   fpy = mother->Yv()*1e13;
578                   fpz = mother->Zv()*1e13;
579                   //          fpt = mother->T() *1e13*3e10;
580                   
581                 }
582               }
583             }
584             
585             //       if (fRotateToEventPlane) {
586             //  double tPhi = TMath::ATan2(fpy, fpx);
587             //  double tRad = TMath::Hypot(fpx, fpy);
588             
589             //  fpx = tRad*TMath::Cos(tPhi - tReactionPlane);
590             //  fpy = tRad*TMath::Sin(tPhi - tReactionPlane);
591             //       }
592             
593             tInfo->SetPDGPid(tPart->GetPdgCode());
594             
595             //    if (fRotateToEventPlane) {
596             //      double tPhi = TMath::ATan2(tPart->Py(), tPart->Px());
597             //      double tRad = TMath::Hypot(tPart->Px(), tPart->Py());
598             
599             //      tInfo->SetTrueMomentum(tRad*TMath::Cos(tPhi - tReactionPlane),
600             //                             tRad*TMath::Sin(tPhi - tReactionPlane),
601             //                             tPart->Pz());
602             //    }
603             //       else
604             tInfo->SetTrueMomentum(tPart->Px(), tPart->Py(), tPart->Pz());
605             Double_t mass2 = (tPart->E() *tPart->E() -
606                               tPart->Px()*tPart->Px() -
607                               tPart->Py()*tPart->Py() -
608                               tPart->Pz()*tPart->Pz());
609             if (mass2>0.0)
610               tInfo->SetMass(TMath::Sqrt(mass2));
611             else 
612               tInfo->SetMass(0.0);
613             
614             tInfo->SetEmissionPoint(fpx, fpy, fpz, fpt);
615           }
616           trackCopy->SetHiddenInfo(tInfo);
617         }
618
619         double pxyz[3];
620
621         //AliExternalTrackParam *param = new AliExternalTrackParam(*aodtrack->GetInnerParam());
622         trackCopy->SetInnerMomentum(aodtrack->GetTPCmomentum());
623
624         aodtrack->PxPyPz(pxyz);//reading noconstarined momentum
625         const AliFmThreeVectorD ktP(pxyz[0],pxyz[1],pxyz[2]);
626         // Check the sanity of the tracks - reject zero momentum tracks
627         if (ktP.Mag() == 0) {
628           delete trackCopy;
629           continue;
630         }
631         //    }
632   
633   
634         tEvent->TrackCollection()->push_back(trackCopy);//adding track to analysis
635         realnofTracks++;//real number of tracks         
636     }
637   
638   tEvent->SetNumberOfTracks(realnofTracks);//setting number of track which we read in event     
639   tEvent->SetNormalizedMult(tracksPrim);
640
641
642   if (fEstEventMult==kCentrality) {
643     //  AliCentrality *cent = fEvent->GetCentrality();
644     if (cent) tEvent->SetNormalizedMult(lrint(10*cent->GetCentralityPercentile("V0M")));
645     //  if (cent) tEvent->SetNormalizedMult((int) cent->GetCentralityPercentile("V0M"));
646     
647     if (cent) {
648       tEvent->SetCentralityV0(cent->GetCentralityPercentile("V0M"));
649       //    tEvent->SetCentralityFMD(cent->GetCentralityPercentile("FMD"));
650       tEvent->SetCentralitySPD1(cent->GetCentralityPercentile("CL1"));
651       //    tEvent->SetCentralityTrk(cent->GetCentralityPercentile("TRK"));
652     }
653   }
654   else if(fEstEventMult==kGlobalCount){
655     tEvent->SetNormalizedMult(tNormMult); //particles counted in the loop, trying to reproduce GetReferenceMultiplicity. If better (default) method appears it should be changed
656   }
657   else if(fEstEventMult==kReference)
658     {
659       tEvent->SetNormalizedMult(fAODheader->GetRefMultiplicity());
660     }
661   else if(fEstEventMult==kTPCOnlyRef)
662     {
663       tEvent->SetNormalizedMult(fAODheader->GetTPConlyRefMultiplicity());
664     }
665
666   if (mcP) delete [] motherids;
667
668     // cout<<"end of reading nt "<<nofTracks<<" real number "<<realnofTracks<<endl;
669
670   if(fReadV0)
671     {
672       int count_pass = 0;
673       for (Int_t i = 0; i < fEvent->GetNumberOfV0s(); i++) {
674         AliAODv0* aodv0 = fEvent->GetV0(i);
675         if (!aodv0) continue;
676         if(aodv0->GetNDaughters()>2) continue;
677         if(aodv0->GetNProngs()>2) continue;
678         if(aodv0->GetCharge()!=0) continue;
679         if(aodv0->ChargeProng(0)==aodv0->ChargeProng(1)) continue;
680         if(aodv0->CosPointingAngle(fV1)<0.998) continue;
681         AliFemtoV0* trackCopyV0 = new AliFemtoV0();
682         count_pass++;
683         CopyAODtoFemtoV0(aodv0, trackCopyV0);
684         tEvent->V0Collection()->push_back(trackCopyV0);
685         //cout<<"Pushback v0 to v0collection"<<endl;
686       }
687     }
688
689 }
690
691 void AliFemtoEventReaderAOD::CopyAODtoFemtoTrack(AliAODTrack *tAodTrack, 
692                                                  AliFemtoTrack *tFemtoTrack 
693                                                  //                                              AliPWG2AODTrack *tPWG2AODTrack
694                                                  )
695 {
696   // Copy the track information from the AOD into the internal AliFemtoTrack
697   // If it exists, use the additional information from the PWG2 AOD
698
699   // Primary Vertex position
700   
701   fEvent->GetPrimaryVertex()->GetPosition(fV1);
702   //  fEvent->GetPrimaryVertex()->GetXYZ(fV1);
703
704   tFemtoTrack->SetCharge(tAodTrack->Charge());
705   
706   double pxyz[3];
707   tAodTrack->PxPyPz(pxyz);//reading noconstrained momentum
708   AliFemtoThreeVector v(pxyz[0],pxyz[1],pxyz[2]);
709   tFemtoTrack->SetP(v);//setting momentum
710   tFemtoTrack->SetPt(sqrt(pxyz[0]*pxyz[0]+pxyz[1]*pxyz[1]));
711   const AliFmThreeVectorD kOrigin(fV1[0],fV1[1],fV1[2]);
712   //setting track helix 
713   const AliFmThreeVectorD ktP(pxyz[0],pxyz[1],pxyz[2]);
714   AliFmPhysicalHelixD helix(ktP,kOrigin,(double)(fEvent->GetMagneticField())*kilogauss,(double)(tFemtoTrack->Charge())); 
715   tFemtoTrack->SetHelix(helix);
716                 
717   // Flags
718   tFemtoTrack->SetTrackId(tAodTrack->GetID());
719   tFemtoTrack->SetFlags(tAodTrack->GetFlags());
720   tFemtoTrack->SetLabel(tAodTrack->GetLabel());
721                 
722   // Track quality information 
723   float covmat[6];
724   tAodTrack->GetCovMatrix(covmat);  
725
726         // ! DCA information is done in CopyPIDtoFemtoTrack()
727   
728         // double impact[2];
729         // double covimpact[3];
730
731         // if (!tAodTrack->PropagateToDCA(fEvent->GetPrimaryVertex(),fEvent->GetMagneticField(),10000,impact,covimpact)) {
732         //   //cout << "sth went wrong with dca propagation" << endl;
733         //   tFemtoTrack->SetImpactD(-1000.0);
734         //   tFemtoTrack->SetImpactZ(-1000.0);
735
736         // }
737         // else {
738         //   tFemtoTrack->SetImpactD(impact[0]);
739         //   tFemtoTrack->SetImpactZ(impact[1]+fV1[2]);
740         // }
741
742   //   if (TMath::Abs(tAodTrack->Xv()) > 0.00000000001)
743   //     tFemtoTrack->SetImpactD(TMath::Hypot(tAodTrack->Xv(), tAodTrack->Yv())*(tAodTrack->Xv()/TMath::Abs(tAodTrack->Xv())));
744   //   else
745   //     tFemtoTrack->SetImpactD(0.0);
746   //   tFemtoTrack->SetImpactD(tAodTrack->DCA());
747     
748   //   tFemtoTrack->SetImpactZ(tAodTrack->ZAtDCA());
749
750
751   //   tFemtoTrack->SetImpactD(TMath::Hypot(tAodTrack->Xv() - fV1[0], tAodTrack->Yv() - fV1[1]));
752   //   tFemtoTrack->SetImpactZ(tAodTrack->Zv() - fV1[2]);
753
754
755   //   cout 
756     //    << "dca" << TMath::Hypot(tAodTrack->Xv() - fV1[0], tAodTrack->Yv() - fV1[1]) 
757     //    << "xv - fv10 = "<< tAodTrack->Xv() - fV1[0] 
758     //    << tAodTrack->Yv() - fV1[1] 
759 //     << "xv = " << tAodTrack->Xv() << endl 
760 //     << "fv1[0] = " << fV1[0]  << endl 
761 //     << "yv = " << tAodTrack->Yv()  << endl 
762 //     << "fv1[1] = " << fV1[1]  << endl 
763 //     << "zv = " << tAodTrack->Zv()  << endl 
764 //     << "fv1[2] = " << fV1[2]  << endl 
765 //     << "impact[0] = " << impact[0]  << endl 
766 //     << "impact[1] = " << impact[1]  << endl 
767 //     << endl << endl ;
768
769   tFemtoTrack->SetCdd(covmat[0]);
770   tFemtoTrack->SetCdz(covmat[1]);
771   tFemtoTrack->SetCzz(covmat[2]);
772   tFemtoTrack->SetITSchi2(tAodTrack->Chi2perNDF());
773   tFemtoTrack->SetITSncls(tAodTrack->GetITSNcls());
774   tFemtoTrack->SetTPCchi2(tAodTrack->Chi2perNDF());
775   tFemtoTrack->SetTPCncls(tAodTrack->GetTPCNcls());
776   tFemtoTrack->SetTPCnclsF(tAodTrack->GetTPCNcls());
777   tFemtoTrack->SetTPCsignalN(1); 
778   tFemtoTrack->SetTPCsignalS(1); 
779   tFemtoTrack->SetTPCsignal(tAodTrack->GetTPCsignal());
780
781 //   if (tPWG2AODTrack) {
782 //     // Copy the PWG2 specific information if it exists
783 //     tFemtoTrack->SetTPCClusterMap(tPWG2AODTrack->GetTPCClusterMap());
784 //     tFemtoTrack->SetTPCSharedMap(tPWG2AODTrack->GetTPCSharedMap());
785     
786 //     double xtpc[3] = {0,0,0};
787 //     tPWG2AODTrack->GetTPCNominalEntrancePoint(xtpc);
788 //     tFemtoTrack->SetNominalTPCEntrancePoint(xtpc);
789 //     tPWG2AODTrack->GetTPCNominalExitPoint(xtpc);
790 //     tFemtoTrack->SetNominalTPCExitPoint(xtpc);
791 //   }
792 //   else {
793     // If not use dummy values
794   tFemtoTrack->SetTPCClusterMap(tAodTrack->GetTPCClusterMap());
795   tFemtoTrack->SetTPCSharedMap(tAodTrack->GetTPCSharedMap());
796   
797
798   float globalPositionsAtRadii[9][3];
799   float bfield = 5*fMagFieldSign;
800   GetGlobalPositionAtGlobalRadiiThroughTPC(tAodTrack,bfield,globalPositionsAtRadii);
801   double tpcEntrance[3]={globalPositionsAtRadii[0][0],globalPositionsAtRadii[0][1],globalPositionsAtRadii[0][2]};
802   double **tpcPositions;
803   tpcPositions = new double*[9];
804   for(int i=0;i<9;i++)
805     tpcPositions[i] = new double[3];
806   double tpcExit[3]={globalPositionsAtRadii[8][0],globalPositionsAtRadii[8][1],globalPositionsAtRadii[8][2]};
807   for(int i=0;i<9;i++)
808     {
809       tpcPositions[i][0] = globalPositionsAtRadii[i][0];
810       tpcPositions[i][1] = globalPositionsAtRadii[i][1];
811       tpcPositions[i][2] = globalPositionsAtRadii[i][2];
812     }
813   tFemtoTrack->SetNominalTPCEntrancePoint(tpcEntrance);
814   tFemtoTrack->SetNominalTPCPoints(tpcPositions);
815   tFemtoTrack->SetNominalTPCExitPoint(tpcExit);
816
817   //   }
818   
819   //   //  cout << "Track has " << TMath::Hypot(tAodTrack->Xv(), tAodTrack->Yv()) << "  " << tAodTrack->Zv() << "  " << tAodTrack->GetTPCNcls() << endl;
820   
821   
822   int indexes[3];
823   for (int ik=0; ik<3; ik++) {
824     indexes[ik] = 0;
825   }
826   tFemtoTrack->SetKinkIndexes(indexes);
827 }
828
829 void AliFemtoEventReaderAOD::CopyAODtoFemtoV0(AliAODv0 *tAODv0, AliFemtoV0 *tFemtoV0)
830 {
831   tFemtoV0->SetdecayLengthV0(tAODv0->DecayLength(fV1));
832   tFemtoV0->SetdecayVertexV0X(tAODv0->DecayVertexV0X());
833   tFemtoV0->SetdecayVertexV0Y(tAODv0->DecayVertexV0Y());
834   tFemtoV0->SetdecayVertexV0Z(tAODv0->DecayVertexV0Z());
835   AliFemtoThreeVector decayvertex(tAODv0->DecayVertexV0X(),tAODv0->DecayVertexV0Y(),tAODv0->DecayVertexV0Z());
836   tFemtoV0->SetdecayVertexV0(decayvertex);
837   tFemtoV0->SetdcaV0Daughters(tAODv0->DcaV0Daughters());
838   tFemtoV0->SetdcaV0ToPrimVertex(tAODv0->DcaV0ToPrimVertex());
839   tFemtoV0->SetdcaPosToPrimVertex(tAODv0->DcaPosToPrimVertex());
840   tFemtoV0->SetdcaNegToPrimVertex(tAODv0->DcaNegToPrimVertex());
841   tFemtoV0->SetmomPosX(tAODv0->MomPosX());
842   tFemtoV0->SetmomPosY(tAODv0->MomPosY());
843   tFemtoV0->SetmomPosZ(tAODv0->MomPosZ());
844   AliFemtoThreeVector mompos(tAODv0->MomPosX(),tAODv0->MomPosY(),tAODv0->MomPosZ());
845   tFemtoV0->SetmomPos(mompos);
846   tFemtoV0->SetmomNegX(tAODv0->MomNegX());
847   tFemtoV0->SetmomNegY(tAODv0->MomNegY());
848   tFemtoV0->SetmomNegZ(tAODv0->MomNegZ());
849   AliFemtoThreeVector momneg(tAODv0->MomNegX(),tAODv0->MomNegY(),tAODv0->MomNegZ());
850   tFemtoV0->SetmomNeg(momneg);
851
852   //jest cos takiego w AliFemtoV0.h czego nie ma w AliAODv0.h
853   //void SettpcHitsPos(const int& i);      
854   //void SettpcHitsNeg(const int& i);      
855
856   //void SetTrackTopologyMapPos(unsigned int word, const unsigned long& m);
857   //void SetTrackTopologyMapNeg(unsigned int word, const unsigned long& m);
858
859   tFemtoV0->SetmomV0X(tAODv0->MomV0X());
860   tFemtoV0->SetmomV0Y(tAODv0->MomV0Y());
861   tFemtoV0->SetmomV0Z(tAODv0->MomV0Z());
862   AliFemtoThreeVector momv0(tAODv0->MomV0X(),tAODv0->MomV0Y(),tAODv0->MomV0Z());
863   tFemtoV0->SetmomV0(momv0);
864   tFemtoV0->SetalphaV0(tAODv0->AlphaV0());
865   tFemtoV0->SetptArmV0(tAODv0->PtArmV0());
866   tFemtoV0->SeteLambda(tAODv0->ELambda());
867   tFemtoV0->SeteK0Short(tAODv0->EK0Short());
868   tFemtoV0->SetePosProton(tAODv0->EPosProton());
869   tFemtoV0->SeteNegProton(tAODv0->ENegProton());
870   tFemtoV0->SetmassLambda(tAODv0->MassLambda());
871   tFemtoV0->SetmassAntiLambda(tAODv0->MassAntiLambda());
872   tFemtoV0->SetmassK0Short(tAODv0->MassK0Short());
873   tFemtoV0->SetrapLambda(tAODv0->RapLambda());
874   tFemtoV0->SetrapK0Short(tAODv0->RapK0Short());
875   
876   //void SetcTauLambda( float x);   
877   //void SetcTauK0Short( float x); 
878   
879   //tFemtoV0->SetptV0(::sqrt(tAODv0->Pt2V0())); //!
880   tFemtoV0->SetptV0(tAODv0->Pt());
881   tFemtoV0->SetptotV0(::sqrt(tAODv0->Ptot2V0()));
882   //tFemtoV0->SetptPos(::sqrt(tAODv0->MomPosX()*tAODv0->MomPosX()+tAODv0->MomPosY()*tAODv0->MomPosY()));
883   //tFemtoV0->SetptotPos(::sqrt(tAODv0->Ptot2Pos()));
884   //tFemtoV0->SetptNeg(::sqrt(tAODv0->MomNegX()*tAODv0->MomNegX()+tAODv0->MomNegY()*tAODv0->MomNegY()));
885   //tFemtoV0->SetptotNeg(::sqrt(tAODv0->Ptot2Neg()));
886   
887   tFemtoV0->SetidNeg(tAODv0->GetNegID());
888   //cout<<"tAODv0->GetNegID(): "<<tAODv0->GetNegID()<<endl;
889   //cout<<"tFemtoV0->IdNeg(): "<<tFemtoV0->IdNeg()<<endl;
890   tFemtoV0->SetidPos(tAODv0->GetPosID());
891
892   tFemtoV0->SetEtaV0(tAODv0->Eta());
893   tFemtoV0->SetPhiV0(tAODv0->Phi());
894   tFemtoV0->SetCosPointingAngle(tAODv0->CosPointingAngle(fV1));
895   //tFemtoV0->SetYV0(tAODv0->Y());
896
897   //void SetdedxNeg(float x);
898   //void SeterrdedxNeg(float x);//Gael 04Fev2002
899   //void SetlendedxNeg(float x);//Gael 04Fev2002
900   //void SetdedxPos(float x);
901   //void SeterrdedxPos(float x);//Gael 04Fev2002
902   //void SetlendedxPos(float x);//Gael 04Fev2002
903
904   //tFemtoV0->SetEtaPos(tAODv0->PseudoRapPos());
905   //tFemtoV0->SetEtaNeg(tAODv0->PseudoRapNeg());
906
907   AliAODTrack *trackpos = (AliAODTrack*)tAODv0->GetDaughter(0);
908   AliAODTrack *trackneg = (AliAODTrack*)tAODv0->GetDaughter(1);
909
910   if(trackpos && trackneg)
911     {
912       tFemtoV0->SetEtaPos(trackpos->Eta());
913       tFemtoV0->SetEtaNeg(trackneg->Eta());
914       tFemtoV0->SetptotPos(tAODv0->PProng(0));
915       tFemtoV0->SetptotNeg(tAODv0->PProng(1));
916       tFemtoV0->SetptPos(trackpos->Pt());
917       tFemtoV0->SetptNeg(trackneg->Pt());
918
919       //tFemtoV0->SetEtaPos(trackpos->Eta()); //tAODv0->PseudoRapPos()
920       //tFemtoV0->SetEtaNeg(trackneg->Eta()); //tAODv0->PseudoRapNeg()
921       tFemtoV0->SetTPCNclsPos(trackpos->GetTPCNcls());
922       tFemtoV0->SetTPCNclsNeg(trackneg->GetTPCNcls());
923       tFemtoV0->SetTPCclustersPos(trackpos->GetTPCClusterMap());
924       tFemtoV0->SetTPCclustersNeg(trackneg->GetTPCClusterMap());
925       tFemtoV0->SetTPCsharingPos(trackpos->GetTPCSharedMap());
926       tFemtoV0->SetTPCsharingNeg(trackneg->GetTPCSharedMap());
927       tFemtoV0->SetNdofPos(trackpos->Chi2perNDF());
928       tFemtoV0->SetNdofNeg(trackneg->Chi2perNDF());
929       tFemtoV0->SetStatusPos(trackpos->GetStatus());
930       tFemtoV0->SetStatusNeg(trackneg->GetStatus());
931
932       tFemtoV0->SetPosNSigmaTPCK(fAODpidUtil->NumberOfSigmasTPC(trackpos,AliPID::kKaon));
933       tFemtoV0->SetNegNSigmaTPCK(fAODpidUtil->NumberOfSigmasTPC(trackneg,AliPID::kKaon));
934       tFemtoV0->SetPosNSigmaTPCP(fAODpidUtil->NumberOfSigmasTPC(trackpos,AliPID::kProton));
935       tFemtoV0->SetNegNSigmaTPCP(fAODpidUtil->NumberOfSigmasTPC(trackneg,AliPID::kProton));
936       tFemtoV0->SetPosNSigmaTPCPi(fAODpidUtil->NumberOfSigmasTPC(trackpos,AliPID::kPion));
937       tFemtoV0->SetNegNSigmaTPCPi(fAODpidUtil->NumberOfSigmasTPC(trackneg,AliPID::kPion));
938
939
940       float bfield = 5*fMagFieldSign;
941       float globalPositionsAtRadiiPos[9][3];
942       GetGlobalPositionAtGlobalRadiiThroughTPC(trackpos,bfield,globalPositionsAtRadiiPos);
943       double tpcEntrancePos[3]={globalPositionsAtRadiiPos[0][0],globalPositionsAtRadiiPos[0][1],globalPositionsAtRadiiPos[0][2]};
944       double tpcExitPos[3]={globalPositionsAtRadiiPos[8][0],globalPositionsAtRadiiPos[8][1],globalPositionsAtRadiiPos[8][2]};
945
946       float globalPositionsAtRadiiNeg[9][3];
947       GetGlobalPositionAtGlobalRadiiThroughTPC(trackneg,bfield,globalPositionsAtRadiiNeg);
948       double tpcEntranceNeg[3]={globalPositionsAtRadiiNeg[0][0],globalPositionsAtRadiiNeg[0][1],globalPositionsAtRadiiNeg[0][2]};
949       double tpcExitNeg[3]={globalPositionsAtRadiiNeg[8][0],globalPositionsAtRadiiNeg[8][1],globalPositionsAtRadiiNeg[8][2]};
950
951       AliFemtoThreeVector tmpVec;
952       tmpVec.SetX(tpcEntrancePos[0]); tmpVec.SetX(tpcEntrancePos[1]); tmpVec.SetX(tpcEntrancePos[2]);
953       tFemtoV0->SetNominalTpcEntrancePointPos(tmpVec);
954
955       tmpVec.SetX(tpcExitPos[0]); tmpVec.SetX(tpcExitPos[1]); tmpVec.SetX(tpcExitPos[2]);
956       tFemtoV0->SetNominalTpcExitPointPos(tmpVec);
957
958       tmpVec.SetX(tpcEntranceNeg[0]); tmpVec.SetX(tpcEntranceNeg[1]); tmpVec.SetX(tpcEntranceNeg[2]);
959       tFemtoV0->SetNominalTpcEntrancePointNeg(tmpVec);
960
961       tmpVec.SetX(tpcExitNeg[0]); tmpVec.SetX(tpcExitNeg[1]); tmpVec.SetX(tpcExitNeg[2]);
962       tFemtoV0->SetNominalTpcExitPointNeg(tmpVec);
963
964       AliFemtoThreeVector vecTpcPos[9];
965       AliFemtoThreeVector vecTpcNeg[9];
966       for(int i=0;i<9;i++)
967         {
968           vecTpcPos[i].SetX(globalPositionsAtRadiiPos[i][0]); vecTpcPos[i].SetY(globalPositionsAtRadiiPos[i][1]); vecTpcPos[i].SetZ(globalPositionsAtRadiiPos[i][2]);
969           vecTpcNeg[i].SetX(globalPositionsAtRadiiNeg[i][0]); vecTpcNeg[i].SetY(globalPositionsAtRadiiNeg[i][1]); vecTpcNeg[i].SetZ(globalPositionsAtRadiiNeg[i][2]);
970         }
971       tFemtoV0->SetNominalTpcPointPos(vecTpcPos);
972       tFemtoV0->SetNominalTpcPointNeg(vecTpcNeg);
973
974       tFemtoV0->SetTPCMomentumPos(trackpos->GetTPCmomentum());
975       tFemtoV0->SetTPCMomentumNeg(trackneg->GetTPCmomentum());
976
977       tFemtoV0->SetdedxPos(trackpos->GetTPCsignal());
978       tFemtoV0->SetdedxNeg(trackneg->GetTPCsignal());
979
980       if((tFemtoV0->StatusPos()&AliESDtrack::kTOFpid)==0 || (tFemtoV0->StatusPos()&AliESDtrack::kTIME)==0 || (tFemtoV0->StatusPos()&AliESDtrack::kTOFout)==0 || (tFemtoV0->StatusPos()&AliESDtrack::kTOFmismatch)>0)
981         {
982           if((tFemtoV0->StatusNeg()&AliESDtrack::kTOFpid)==0 || (tFemtoV0->StatusNeg()&AliESDtrack::kTIME)==0 || (tFemtoV0->StatusNeg()&AliESDtrack::kTOFout)==0 || (tFemtoV0->StatusNeg()&AliESDtrack::kTOFmismatch)>0)
983             {
984               tFemtoV0->SetPosNSigmaTOFK(-1000);
985               tFemtoV0->SetNegNSigmaTOFK(-1000);
986               tFemtoV0->SetPosNSigmaTOFP(-1000);
987               tFemtoV0->SetNegNSigmaTOFP(-1000);
988               tFemtoV0->SetPosNSigmaTOFPi(-1000);
989               tFemtoV0->SetNegNSigmaTOFPi(-1000);
990
991               tFemtoV0->SetTOFProtonTimePos(-1000);
992               tFemtoV0->SetTOFPionTimePos(-1000);
993               tFemtoV0->SetTOFKaonTimePos(-1000);
994               tFemtoV0->SetTOFProtonTimeNeg(-1000);
995               tFemtoV0->SetTOFPionTimeNeg(-1000);
996               tFemtoV0->SetTOFKaonTimeNeg(-1000);
997             }
998         }
999       else
1000         {
1001           tFemtoV0->SetPosNSigmaTOFK(fAODpidUtil->NumberOfSigmasTOF(trackpos,AliPID::kKaon));
1002           tFemtoV0->SetNegNSigmaTOFK(fAODpidUtil->NumberOfSigmasTOF(trackneg,AliPID::kKaon));
1003           tFemtoV0->SetPosNSigmaTOFP(fAODpidUtil->NumberOfSigmasTOF(trackpos,AliPID::kProton));
1004           tFemtoV0->SetNegNSigmaTOFP(fAODpidUtil->NumberOfSigmasTOF(trackneg,AliPID::kProton));
1005           tFemtoV0->SetPosNSigmaTOFPi(fAODpidUtil->NumberOfSigmasTOF(trackpos,AliPID::kPion));
1006           tFemtoV0->SetNegNSigmaTOFPi(fAODpidUtil->NumberOfSigmasTOF(trackneg,AliPID::kPion));
1007
1008           double TOFSignalPos = trackpos->GetTOFsignal();
1009           double TOFSignalNeg = trackneg->GetTOFsignal();
1010           double pidPos[5];
1011           double pidNeg[5];
1012           trackpos->GetIntegratedTimes(pidPos);
1013           trackneg->GetIntegratedTimes(pidNeg);
1014
1015           tFemtoV0->SetTOFPionTimePos(TOFSignalPos-pidPos[2]);
1016           tFemtoV0->SetTOFKaonTimePos(TOFSignalPos-pidPos[3]);
1017           tFemtoV0->SetTOFProtonTimePos(TOFSignalPos-pidPos[4]);
1018           tFemtoV0->SetTOFPionTimeNeg(TOFSignalNeg-pidNeg[2]);
1019           tFemtoV0->SetTOFKaonTimeNeg(TOFSignalNeg-pidNeg[3]);
1020           tFemtoV0->SetTOFProtonTimeNeg(TOFSignalNeg-pidNeg[4]);
1021         }
1022     }
1023   else
1024     {
1025       tFemtoV0->SetStatusPos(999);
1026       tFemtoV0->SetStatusNeg(999);
1027     }
1028   tFemtoV0->SetOnFlyStatusV0(tAODv0->GetOnFlyStatus());
1029 }
1030
1031 void AliFemtoEventReaderAOD::SetFilterBit(UInt_t ibit)
1032 {
1033   fFilterBit = (1 << (ibit));
1034 }
1035
1036 void AliFemtoEventReaderAOD::SetReadMC(unsigned char a)
1037 {
1038   fReadMC = a;
1039 }
1040
1041
1042 void AliFemtoEventReaderAOD::SetReadV0(unsigned char a)
1043 {
1044   fReadV0 = a;
1045 }
1046
1047 void AliFemtoEventReaderAOD::SetUseMultiplicity(EstEventMult aType)
1048 {
1049   fEstEventMult = aType;
1050 }
1051
1052 AliAODMCParticle* AliFemtoEventReaderAOD::GetParticleWithLabel(TClonesArray *mcP, Int_t aLabel)
1053 {
1054   if (aLabel < 0) return 0;
1055   AliAODMCParticle *aodP;
1056   Int_t posstack = 0;
1057   if (aLabel > mcP->GetEntries())
1058     posstack = mcP->GetEntries();
1059   else
1060     posstack = aLabel;
1061
1062   aodP = (AliAODMCParticle *) mcP->At(posstack);
1063   if (aodP->GetLabel() > posstack) {
1064     do {
1065       aodP = (AliAODMCParticle *) mcP->At(posstack);
1066       if (aodP->GetLabel() == aLabel) return aodP;
1067       posstack--;
1068     }
1069     while (posstack > 0);
1070   }
1071   else {
1072     do {
1073       aodP = (AliAODMCParticle *) mcP->At(posstack);
1074       if (aodP->GetLabel() == aLabel) return aodP;
1075       posstack++;
1076     }
1077     while (posstack < mcP->GetEntries());
1078   }
1079   
1080   return 0;
1081 }
1082
1083 void AliFemtoEventReaderAOD::CopyPIDtoFemtoTrack(AliAODTrack *tAodTrack, 
1084                                                  AliFemtoTrack *tFemtoTrack)
1085 {
1086
1087         // copying DCA information (taking it from global tracks gives better resolution than from TPC-only)
1088
1089         double impact[2];
1090         double covimpact[3];
1091
1092         if (!tAodTrack->PropagateToDCA(fEvent->GetPrimaryVertex(),fEvent->GetMagneticField(),10000,impact,covimpact)) {
1093                 //cout << "sth went wrong with dca propagation" << endl;
1094                 tFemtoTrack->SetImpactD(-1000.0);
1095                 tFemtoTrack->SetImpactZ(-1000.0);
1096
1097         }
1098         else {
1099                 tFemtoTrack->SetImpactD(impact[0]);
1100                 tFemtoTrack->SetImpactZ(impact[1]);
1101         }
1102
1103   double aodpid[10];
1104   tAodTrack->GetPID(aodpid);
1105   tFemtoTrack->SetPidProbElectron(aodpid[0]);
1106   tFemtoTrack->SetPidProbMuon(aodpid[1]);
1107   tFemtoTrack->SetPidProbPion(aodpid[2]);
1108   tFemtoTrack->SetPidProbKaon(aodpid[3]);
1109   tFemtoTrack->SetPidProbProton(aodpid[4]);
1110
1111   aodpid[0] = -100000.0;
1112   aodpid[1] = -100000.0;
1113   aodpid[2] = -100000.0;
1114   aodpid[3] = -100000.0;
1115   aodpid[4] = -100000.0;
1116                 
1117   double tTOF = 0.0;
1118
1119   if (tAodTrack->GetStatus() & AliESDtrack::kTOFpid) {  //AliESDtrack::kTOFpid=0x8000
1120     tTOF = tAodTrack->GetTOFsignal();
1121     tAodTrack->GetIntegratedTimes(aodpid);
1122
1123     tTOF -= fAODpidUtil->GetTOFResponse().GetStartTime(tAodTrack->P());
1124   }
1125
1126   tFemtoTrack->SetTofExpectedTimes(tTOF-aodpid[2], tTOF-aodpid[3], tTOF-aodpid[4]);
1127  
1128   //////  TPC ////////////////////////////////////////////
1129
1130   float nsigmaTPCK=-1000.;                                                  
1131   float nsigmaTPCPi=-1000.;                                                 
1132   float nsigmaTPCP=-1000.;                                                  
1133           
1134   //   cout<<"in reader fESDpid"<<fESDpid<<endl;
1135
1136   if (tAodTrack->IsOn(AliESDtrack::kTPCpid)){ //AliESDtrack::kTPCpid=0x0080
1137     nsigmaTPCK = fAODpidUtil->NumberOfSigmasTPC(tAodTrack,AliPID::kKaon);
1138     nsigmaTPCPi = fAODpidUtil->NumberOfSigmasTPC(tAodTrack,AliPID::kPion);
1139     nsigmaTPCP = fAODpidUtil->NumberOfSigmasTPC(tAodTrack,AliPID::kProton);
1140   }
1141
1142   tFemtoTrack->SetNSigmaTPCPi(nsigmaTPCPi);
1143   tFemtoTrack->SetNSigmaTPCK(nsigmaTPCK);
1144   tFemtoTrack->SetNSigmaTPCP(nsigmaTPCP);
1145
1146   tFemtoTrack->SetTPCchi2(tAodTrack->Chi2perNDF());       
1147   tFemtoTrack->SetTPCncls(tAodTrack->GetTPCNcls());       
1148   tFemtoTrack->SetTPCnclsF(tAodTrack->GetTPCNcls());      
1149   
1150   tFemtoTrack->SetTPCsignalN(1); 
1151   tFemtoTrack->SetTPCsignalS(1); 
1152   tFemtoTrack->SetTPCsignal(tAodTrack->GetTPCsignal());
1153  
1154   ///////TOF//////////////////////
1155
1156     float vp=-1000.;
1157     float nsigmaTOFPi=-1000.;
1158     float nsigmaTOFK=-1000.;
1159     float nsigmaTOFP=-1000.;
1160
1161     if ((tAodTrack->GetStatus() & AliESDtrack::kTOFpid) && //AliESDtrack::kTOFpid=0x8000
1162         (tAodTrack->GetStatus() & AliESDtrack::kTOFout) && //AliESDtrack::kTOFout=0x2000
1163         (tAodTrack->GetStatus() & AliESDtrack::kTIME) && //AliESDtrack::kTIME=0x80000000
1164         !(tAodTrack->GetStatus() & AliESDtrack::kTOFmismatch)) //AliESDtrack::kTOFmismatch=0x100000
1165       {
1166         if(tAodTrack->IsOn(AliESDtrack::kTOFpid)) //AliESDtrack::kTOFpid=0x8000
1167           {
1168
1169             nsigmaTOFPi = fAODpidUtil->NumberOfSigmasTOF(tAodTrack,AliPID::kPion);
1170             nsigmaTOFK = fAODpidUtil->NumberOfSigmasTOF(tAodTrack,AliPID::kKaon);
1171             nsigmaTOFP = fAODpidUtil->NumberOfSigmasTOF(tAodTrack,AliPID::kProton);
1172
1173             Double_t len=200;// esdtrack->GetIntegratedLength(); !!!!!
1174             Double_t tof=tAodTrack->GetTOFsignal();
1175             if(tof > 0.) vp=len/tof/0.03;
1176           }
1177       }
1178     tFemtoTrack->SetVTOF(vp);
1179     tFemtoTrack->SetNSigmaTOFPi(nsigmaTOFPi);
1180     tFemtoTrack->SetNSigmaTOFK(nsigmaTOFK);
1181     tFemtoTrack->SetNSigmaTOFP(nsigmaTOFP);
1182
1183     
1184     //////////////////////////////////////
1185
1186 }
1187
1188 void AliFemtoEventReaderAOD::SetCentralityPreSelection(double min, double max)
1189 {
1190   fCentRange[0] = min; fCentRange[1] = max;
1191   fUsePreCent = 1; 
1192   fEstEventMult = kCentrality;
1193 }
1194
1195 void AliFemtoEventReaderAOD::SetNoCentrality(bool anocent)
1196 {
1197   if(anocent==false) {fEstEventMult=kCentrality;}
1198   else {fEstEventMult=kReference; fUsePreCent = 0; }
1199 }
1200
1201 void AliFemtoEventReaderAOD::SetAODpidUtil(AliAODpidUtil *aAODpidUtil)
1202 {
1203   fAODpidUtil = aAODpidUtil;
1204   //  printf("fAODpidUtil: %x\n",fAODpidUtil);
1205 }
1206
1207 void AliFemtoEventReaderAOD::SetAODheader(AliAODHeader *aAODheader)
1208 {
1209   fAODheader = aAODheader;
1210 }
1211
1212
1213 void AliFemtoEventReaderAOD::SetMagneticFieldSign(int s)
1214 {
1215   if(s>0)
1216     fMagFieldSign = 1;
1217   else if(s<0)
1218     fMagFieldSign = -1;
1219   else
1220     fMagFieldSign = 0;
1221 }
1222
1223 void AliFemtoEventReaderAOD::SetEPVZERO(Bool_t iepvz)
1224 {
1225     fisEPVZ = iepvz;
1226 }
1227
1228 void AliFemtoEventReaderAOD::GetGlobalPositionAtGlobalRadiiThroughTPC(AliAODTrack *track, Float_t bfield, Float_t globalPositionsAtRadii[9][3])
1229 {
1230   // Gets the global position of the track at nine different radii in the TPC
1231   // track is the track you want to propagate
1232   // bfield is the magnetic field of your event
1233   // globalPositionsAtRadii is the array of global positions in the radii and xyz
1234
1235   // Initialize the array to something indicating there was no propagation
1236   for(Int_t i=0;i<9;i++){
1237     for(Int_t j=0;j<3;j++){
1238       globalPositionsAtRadii[i][j]=-9999.;
1239     }
1240   }
1241
1242   // Make a copy of the track to not change parameters of the track
1243   AliExternalTrackParam etp; etp.CopyFromVTrack(track);
1244   //printf("\nAfter CopyFromVTrack\n");
1245   //etp.Print();
1246
1247   // The global position of the the track
1248   Double_t xyz[3]={-9999.,-9999.,-9999.};
1249
1250   // Counter for which radius we want
1251   Int_t iR=0;
1252   // The radii at which we get the global positions
1253   // IROC (OROC) from 84.1 cm to 132.1 cm (134.6 cm to 246.6 cm)
1254   Float_t Rwanted[9]={85.,105.,125.,145.,165.,185.,205.,225.,245.};
1255   // The global radius we are at
1256   Float_t globalRadius=0;
1257
1258   // Propagation is done in local x of the track
1259   for (Float_t x = etp.GetX();x<247.;x+=1.){ // GetX returns local coordinates
1260     // Starts at the tracks fX and goes outwards. x = 245 is the outer radial limit
1261     // of the TPC when the track is straight, i.e. has inifinite pt and doesn't get bent.
1262     // If the track's momentum is smaller than infinite, it will develop a y-component, which
1263     // adds to the global radius
1264
1265     // Stop if the propagation was not succesful. This can happen for low pt tracks
1266     // that don't reach outer radii
1267     if(!etp.PropagateTo(x,bfield))break;
1268     etp.GetXYZ(xyz); // GetXYZ returns global coordinates
1269     globalRadius = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]); //Idea to speed up: compare squared radii
1270
1271     // Roughly reached the radius we want
1272     if(globalRadius > Rwanted[iR]){
1273
1274       // Bigger loop has bad precision, we're nearly one centimeter too far, go back in small steps.
1275       while (globalRadius>Rwanted[iR]){
1276         x-=.1;
1277         //      printf("propagating to x %5.2f\n",x);
1278         if(!etp.PropagateTo(x,bfield))break;
1279         etp.GetXYZ(xyz); // GetXYZ returns global coordinates
1280         globalRadius = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]); //Idea to speed up: compare squared radii
1281       }
1282       //printf("At Radius:%05.2f (local x %5.2f). Setting position to x %4.1f y %4.1f z %4.1f\n",globalRadius,x,xyz[0],xyz[1],xyz[2]);
1283       globalPositionsAtRadii[iR][0]=xyz[0];
1284       globalPositionsAtRadii[iR][1]=xyz[1];
1285       globalPositionsAtRadii[iR][2]=xyz[2];
1286       // Indicate we want the next radius
1287       iR+=1;
1288     }
1289     if(iR>=8){
1290       // TPC edge reached
1291       return;
1292     }
1293   }
1294 }
1295
1296