]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/exa/AliTPCComparison.C
Introduction of a fast version of the AliITSVertexerZ for HLT
[u/mrichter/AliRoot.git] / HLT / exa / AliTPCComparison.C
1 /****************************************************************************
2  *           Very important, delicate and rather obscure macro.             *
3  *                                                                          *
4  *               Creates list of "trackable" tracks,                        *
5  *             calculates efficiency, resolutions etc.                      *
6  *     There is a possibility to run this macro over several events.         *
7  *                                                                          *
8  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
9  * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de        *
10  ****************************************************************************/
11
12 #if !defined(__CINT__) || defined(__MAKECINT__)
13   #include <TMath.h>
14   #include <TError.h>
15   #include <Riostream.h>
16   #include <TH1F.h>
17   #include <TH2F.h>
18   #include <TTree.h>
19   #include <TParticle.h>
20   #include <TPad.h>
21   #include <TCanvas.h>
22   #include <TLine.h>
23   #include <TText.h>
24   #include <TBenchmark.h>
25   #include <TStyle.h>
26   #include <TKey.h>
27   #include <TROOT.h>
28
29   #include "AliStack.h"
30   #include "AliHeader.h"
31   #include "AliTrackReference.h"
32   #include "AliRunLoader.h"
33   #include "AliMagF.h"
34   #include "AliRun.h"
35   #include "AliESD.h"
36
37   #include "AliSimDigits.h"
38   #include "AliTPC.h"
39   #include "AliTPCParamSR.h"
40   #include "AliTPCClustersArray.h"
41   #include "AliTPCClustersRow.h"
42   #include "AliTPCcluster.h"
43   #include "AliTPCLoader.h"
44   #include "AliConfig.h"
45 #endif
46
47 struct GoodTrackTPC {
48   Int_t lab;
49   Int_t code;
50   Float_t px,py,pz;
51   Float_t x,y,z;
52 };
53
54 Int_t good_tracks_tpc(GoodTrackTPC *gt, const Int_t max, 
55 const char* evfoldname=AliConfig::GetDefaultEventFolderName());
56
57 extern AliRun *gAlice;
58 extern TBenchmark *gBenchmark;
59 extern TROOT *gROOT;
60
61 Int_t AliTPCComparison(Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
62    gBenchmark->Start("AliTPCComparison");
63
64    ::Info("AliTPCComparison.C","Doing comparison...");
65
66    if (gAlice) { 
67        delete gAlice->GetRunLoader();
68        delete gAlice;//if everything was OK here it is already NULL
69        gAlice = 0x0;
70     }
71
72    Char_t fname[100];
73    sprintf(fname,"%s/galice.root",dir);
74
75    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
76    if (!rl) {
77        ::Fatal("AliTPCComparison.C","Can't start session !");
78    }
79
80    /* Generate a list of "good" tracks */
81    const Int_t MAX=20000;
82    GoodTrackTPC gt[MAX];
83    Int_t ngood=0;
84    sprintf(fname,"%s/good_tracks_tpc",dir);
85    ifstream in(fname);
86    if (in) {
87       ::Info("AliTPCComparison.C","Reading good tracks...");
88       while (in>>gt[ngood].lab>>gt[ngood].code>>
89                  gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
90                  gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
91           ngood++;
92           if (ngood==MAX) {
93             ::Warning("AliTPCComparison.C","Too many good tracks !");
94             break;
95           }     
96       }
97       if (!in.eof()) 
98         ::Fatal("AliTPCComparison.C","Read error (good_tracks_tpc) !");
99    } else {
100       ::Info
101       ("AliTPCComparison","Marking good tracks (this will take a while)...");
102       ngood=good_tracks_tpc(gt,MAX,"COMPARISON");
103       ofstream out(fname);
104       if (out) {
105          for (Int_t ngd=0; ngd<ngood; ngd++)            
106             out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
107                  gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
108                  gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
109       } else 
110          ::Fatal("AliTPCComparison.C","Can not open file (good_tracks_tpc) !");
111       out.close();
112    }
113
114    return 0;
115    AliESD *event=0;
116    { /*Load tracks*/
117    sprintf(fname,"%s/AliESDs.root",dir);
118    TFile *ef=TFile::Open(fname);
119    if ((!ef)||(!ef->IsOpen())) {
120       sprintf(fname,"%s/AliESDtpc.root",dir);
121       ef=TFile::Open(fname);
122       if ((!ef)||(!ef->IsOpen()))
123          ::Fatal("AliTPCComparison.C","Can't open AliESDtpc.root !");
124    }
125    TKey *key=0;
126    TIter next(ef->GetListOfKeys());
127    if ((key=(TKey*)next())!=0) event=(AliESD*)key->ReadObj();
128    ef->Close();
129    }
130    Int_t nentr=event->GetNumberOfTracks();
131
132    TH1F *hp=(TH1F*)gROOT->FindObject("hp");
133    if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.); 
134    hp->SetFillColor(4);
135
136    TH1F *hl=(TH1F*)gROOT->FindObject("hl");
137    if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
138    hl->SetFillColor(4);
139
140    TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
141    if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
142    hpt->SetFillColor(2);
143  
144    TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
145    if (!hmpt) 
146       hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60); 
147    hmpt->SetFillColor(6);
148
149
150
151    TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
152    if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
153     
154    TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
155    if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
156
157    TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
158    if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
159
160    TH1F *hg=(TH1F*)gROOT->FindObject("hg");
161    if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
162    hg->SetLineColor(4); hg->SetLineWidth(2);
163
164    TH1F *hf=(TH1F*)gROOT->FindObject("hf");
165    if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
166    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
167
168    TH1F *he=(TH1F*)gROOT->FindObject("he");
169    if (!he) 
170       he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
171
172    TH2F *hep=(TH2F*)gROOT->FindObject("hep");
173    if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
174    hep->SetMarkerStyle(8);
175    hep->SetMarkerSize(0.4);
176
177    //MI change
178    Int_t mingood=ngood;
179    Int_t track_notfound[MAX], itrack_notfound=0;
180    Int_t track_fake[MAX], itrack_fake=0;
181    Int_t track_multifound[MAX], track_multifound_n[MAX], itrack_multifound=0;
182
183    Int_t i;
184    while (ngood--) {
185       Int_t lab=gt[ngood].lab,tlab=-1;
186       Float_t ptg=
187       TMath::Sqrt(gt[ngood].px*gt[ngood].px + gt[ngood].py*gt[ngood].py);
188
189       if (ptg<1e-33) continue; // for those not crossing 0 pad row
190       if (ptg<ptcutl) continue;
191       if (ptg>ptcuth) continue;
192       hgood->Fill(ptg);
193
194       AliESDtrack *track=0;
195       for (i=0; i<nentr; i++) {
196           track=event->GetTrack(i);
197           tlab=track->GetTPCLabel();
198           if (lab==TMath::Abs(tlab)) break;
199       }
200       if (i==nentr) {
201         track_notfound[itrack_notfound++]=lab;
202         continue;
203       }
204       
205       //MI change  - addition
206       Int_t micount=0;
207       Int_t mi;
208       AliESDtrack * mitrack;
209       for (mi=0; mi<nentr; mi++) {
210         mitrack=event->GetTrack(mi);          
211         if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
212       }
213       if (micount>1) {
214         track_multifound[itrack_multifound]=lab;
215         track_multifound_n[itrack_multifound]=micount;
216         itrack_multifound++;
217       }
218       if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
219       if (lab==tlab) hfound->Fill(ptg);
220       else { 
221         track_fake[itrack_fake++]=lab;
222         hfake->Fill(ptg); 
223       }
224       
225       Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
226       Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
227       if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
228       if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
229       Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
230       Float_t lam=TMath::ATan2(pxpypz[2],pt); 
231       Float_t pt_1=1/pt;
232
233       if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
234          hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
235       } else {
236          Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
237          hp->Fill((phi - phig)*1000.);
238
239          Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
240          hl->Fill((lam - lamg)*1000.);
241
242          hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
243       }
244
245       Float_t mom=pt/TMath::Cos(lam);
246       Float_t dedx=track->GetTPCsignal();
247       hep->Fill(mom,dedx,1.);
248       if (TMath::Abs(gt[ngood].code)==211)
249          if (mom>0.4 && mom<0.5) {
250             he->Fill(dedx,1.);
251          }
252
253    }
254
255    cout<<"\nList of Not found tracks :\n";
256    for ( i = 0; i< itrack_notfound; i++){
257      cout<<track_notfound[i]<<"\t";
258      if ((i%5)==4) cout<<"\n";
259    }
260    cout<<"\nList of fake  tracks :\n";
261    for ( i = 0; i< itrack_fake; i++){
262      cout<<track_fake[i]<<"\t";
263      if ((i%5)==4) cout<<"\n";
264    }
265     cout<<"\nList of multiple found tracks :\n";
266    for ( i=0; i<itrack_multifound; i++) {
267        cout<<"id.   "<<track_multifound[i]
268            <<"     found - "<<track_multifound_n[i]<<"times\n";
269    }
270
271    Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
272    if (ng!=0) cerr<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
273    cout<<"Total number of found tracks ="<<nentr<<endl;
274    cout<<"Total number of \"good\" tracks ="
275        <<mingood<<"   (selected for comparison: "<<ng<<')'<<endl<<endl;
276
277
278    gStyle->SetOptStat(111110);
279    gStyle->SetOptFit(1);
280
281    TCanvas *c1=new TCanvas("c1","",0,0,700,850);
282
283    Int_t minc=33; 
284
285    TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
286    p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
287    hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); 
288    if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
289
290    TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
291    p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
292    hl->SetXTitle("(mrad)");
293    if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
294
295    TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
296    p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
297    hpt->SetXTitle("(%)");
298    if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
299
300    TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
301    p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
302    hmpt->SetXTitle("(%)");
303    if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();
304
305    TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); 
306    p5->SetFillColor(41); p5->SetFrameFillColor(10);
307    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
308    hg->Divide(hfound,hgood,1,1.,"b");
309    hf->Divide(hfake,hgood,1,1.,"b");
310    hg->SetMaximum(1.4);
311    hg->SetYTitle("Tracking efficiency");
312    hg->SetXTitle("Pt (GeV/c)");
313    hg->Draw();
314
315    TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
316    line1->Draw("same");
317    TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
318    line2->Draw("same");
319
320    hf->SetFillColor(1);
321    hf->SetFillStyle(3013);
322    hf->SetLineColor(2);
323    hf->SetLineWidth(2);
324    hf->Draw("histsame");
325    TText *text = new TText(0.461176,0.248448,"Fake tracks");
326    text->SetTextSize(0.05);
327    text->Draw();
328    text = new TText(0.453919,1.11408,"Good tracks");
329    text->SetTextSize(0.05);
330    text->Draw();
331
332
333
334    TCanvas *c2=new TCanvas("c2","",320,32,530,590);
335
336    TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
337    p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); 
338    he->SetFillColor(2); he->SetFillStyle(3005);  
339    he->SetXTitle("Arbitrary Units"); 
340    if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
341
342    TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); 
343    p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
344    hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); 
345    hep->Draw(); c1->cd();
346
347    gBenchmark->Stop("AliTPCComparison");
348    gBenchmark->Show("AliTPCComparison");
349
350    TFile fc("AliTPCComparison.root","RECREATE");
351    c1->Write();
352    c2->Write();
353    fc.Close();
354
355    delete rl;
356    return 0;
357 }
358
359
360 Int_t 
361 good_tracks_tpc(GoodTrackTPC *gt, const Int_t max, const char* evfoldname) {
362    Int_t nt=0;
363
364    AliRunLoader *rl = AliRunLoader::GetRunLoader(evfoldname);
365    if (rl == 0x0) {
366       ::Fatal("AliTPCComparison.C::good_tracks_tpc",
367               "Can not find Run Loader in Folder Named %s",
368               evfoldname);
369    }
370    AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
371    if (tpcl == 0x0) {
372        ::Fatal("AliTPCHits2Digits.C","Can not find TPCLoader !");
373    }
374    
375    rl->LoadgAlice();
376    
377    AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
378    Int_t ver = TPC->IsVersion(); 
379    cerr<<"TPC version "<<ver<<" has been found !\n";
380
381    rl->CdGAFile();
382    AliTPCParamSR *digp=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
383    if (!digp) { 
384      ::Fatal("AliTPCHits2Digits.C","TPC parameters have not been found !"); 
385    }
386    TPC->SetParam(digp);
387
388    rl->LoadHeader();
389
390    Int_t np = rl->GetHeader()->GetNtrack();
391       
392
393    Int_t nrow_up=digp->GetNRowUp();
394    Int_t nrows=digp->GetNRowLow()+nrow_up;
395    Int_t zero=digp->GetZeroSup();
396    Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
397    Int_t good_number=Int_t(0.4*nrows);
398
399    Int_t *good=new Int_t[np];
400    Int_t i;
401    for (i=0; i<np; i++) good[i]=0;
402
403
404    switch (ver) {
405    case 1:
406      {
407       tpcl->LoadRecPoints();      
408       AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
409       ca.Setup(digp);
410       ca.SetClusterType("AliTPCcluster");
411       ca.ConnectTree(tpcl->TreeR());
412       Int_t nrows=Int_t(ca.GetTree()->GetEntries());
413       for (Int_t n=0; n<nrows; n++) {
414           AliSegmentID *s=ca.LoadEntry(n);
415           Int_t sec,row;
416           digp->AdjustSectorRow(s->GetID(),sec,row);
417           AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
418           Int_t ncl=clrow.GetArray()->GetEntriesFast();
419           while (ncl--) {
420               AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
421               Int_t lab=c->GetLabel(0);
422               if (lab<0) continue; //noise cluster
423               lab=TMath::Abs(lab);
424
425               if (sec>=digp->GetNInnerSector())
426                  if (row==nrow_up-1) good[lab]|=0x4000;
427               if (sec>=digp->GetNInnerSector())
428                  if (row==nrow_up-1-gap) good[lab]|=0x1000;
429
430               if (sec>=digp->GetNInnerSector())
431                  if (row==nrow_up-1-shift) good[lab]|=0x2000;
432               if (sec>=digp->GetNInnerSector())
433                  if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
434
435               good[lab]++;
436           }
437           ca.ClearRow(sec,row);
438       }
439       delete pca;
440       tpcl->UnloadRecPoints();
441      }
442      break;
443    case 2:
444      {
445       tpcl->LoadDigits();
446       TTree *TD=tpcl->TreeD();
447       
448       AliSimDigits da, *digits=&da;
449       TD->GetBranch("Segment")->SetAddress(&digits);
450
451       Int_t *count = new Int_t[np];
452       Int_t i;
453       for (i=0; i<np; i++) count[i]=0;
454
455       Int_t sectors_by_rows=(Int_t)TD->GetEntries();
456       for (i=0; i<sectors_by_rows; i++) {
457           if (!TD->GetEvent(i)) continue;
458           Int_t sec,row;
459           digp->AdjustSectorRow(digits->GetID(),sec,row);
460           //PH          cerr<<sec<<' '<<row<<"                                     \r";
461           digits->First();
462           do { //Many thanks to J.Chudoba who noticed this
463               Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
464               Short_t dig = digits->GetDigit(it,ip);
465               Int_t idx0=digits->GetTrackID(it,ip,0); 
466               Int_t idx1=digits->GetTrackID(it,ip,1);
467               Int_t idx2=digits->GetTrackID(it,ip,2);
468               if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
469               if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
470               if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
471           } while (digits->Next());
472           for (Int_t j=0; j<np; j++) {
473               if (count[j]>1) {
474                  if (sec>=digp->GetNInnerSector())
475                    if (row==nrow_up-1    ) good[j]|=0x4000;
476                  if (sec>=digp->GetNInnerSector())
477                    if (row==nrow_up-1-gap) good[j]|=0x1000;
478
479                  if (sec>=digp->GetNInnerSector())
480                    if (row==nrow_up-1-shift) good[j]|=0x2000;
481                  if (sec>=digp->GetNInnerSector())
482                    if (row==nrow_up-1-gap-shift) good[j]|=0x800;
483                  good[j]++;
484               }
485               count[j]=0;
486           }
487       }
488       delete[] count;
489       tpcl->UnloadDigits();
490      }
491       break;
492    default:
493       ::Fatal("AliTPCComparison.C","Invalid TPC version !");
494    }
495
496    rl->LoadKinematics();
497    AliStack* stack = rl->Stack();
498    /** select tracks which are "good" enough **/
499    for (i=0; i<np; i++) {
500       if ((good[i]&0x5000) != 0x5000)
501       if ((good[i]&0x2800) != 0x2800) continue;
502       if ((good[i]&0x7FF ) < good_number) continue;
503
504       TParticle *p = (TParticle*)stack->Particle(i);
505       if (p == 0x0)
506        {
507          cerr<<"Can not get particle "<<i<<endl;
508          continue;
509        }
510       if (p->Pt()<0.100) continue;
511       if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
512
513       Double_t vx=p->Vx(),vy=p->Vy();
514       if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
515
516       /*
517       Int_t j=p->GetFirstMother();
518       if (j>=0) {
519         TParticle *pp = (TParticle*)stack->Particle(j);
520         if (pp == 0x0)
521          {
522            cerr<<"Can not get particle "<<j<<endl;
523            continue;
524          }
525         //if (pp->GetFirstMother()>=0) continue;//only one decay is allowed
526         //  for cascade hyperons only
527         Int_t jj=pp->GetFirstMother();
528         if (jj>=0) {
529           TParticle *ppp = (TParticle*)stack->Particle(jj);
530           if (ppp->GetFirstMother()>=0) continue;//two decays are allowed
531         }
532
533       }
534       */
535
536       gt[nt].lab=i;
537       gt[nt].code=p->GetPdgCode();
538       gt[nt].px=0.; gt[nt].py=0.; gt[nt].pz=0.;
539       gt[nt].x=0.; gt[nt].y=0.; gt[nt].z=0.;
540       nt++;
541       if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
542       //      cerr<<np-i<<"                \r";
543    }
544    rl->UnloadKinematics();
545
546    /** check if there is also information at the entrance of the TPC **/
547    
548    rl->LoadTrackRefs();
549    TTree *TR=rl->TreeTR();
550    TBranch *branch=TR->GetBranch("TPC");
551    if (branch==0) {
552       ::Fatal("AliTPCComparison.C::good_tracks_tpc","No track references !");
553    }
554    TClonesArray *references=new TClonesArray("AliTrackReference",1000);
555    branch->SetAddress(&references);
556    np=(Int_t)TR->GetEntries();
557    for (i=0; i<np; i++) {
558       references->Clear();
559       TR->GetEvent(i);
560       Int_t nref=references->GetEntriesFast();
561       if (nref==0) continue;
562       AliTrackReference *ref=(AliTrackReference*)references->UncheckedAt(0);
563
564       Int_t j, lab=ref->Label();
565       for (j=0; j<nt; j++) {if (gt[j].lab==lab) break;}
566       if (j==nt) continue;         
567
568       // (px,py,pz) - in global coordinate system, (x,y,z) - in local !
569       gt[j].px=ref->Px(); gt[j].py=ref->Py(); gt[j].pz=ref->Pz();
570       gt[j].x = ref->LocalX();
571       gt[j].y = ref->LocalY();
572       gt[j].z = ref->Z();
573    }
574
575    delete references;
576    delete[] good;
577    
578    tpcl->UnloadHits();
579    rl->UnloadTrackRefs();
580    rl->UnloadgAlice();
581
582    gBenchmark->Stop("AliTPCComparison");
583    gBenchmark->Show("AliTPCComparison");   
584    return nt;
585 }