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