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