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