]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTrackPoints.cxx
df9bdd83ea2b7690470e7749caff0995e0c8c8ad
[u/mrichter/AliRoot.git] / ANALYSIS / AliTrackPoints.cxx
1 #include "AliTrackPoints.h"
2 //_________________________________
3 ////////////////////////////////////////////////////////////
4 //                                                        //
5 // class AliTrackPoints                                //
6 //                                                        //
7 // used by Anti-Merging cut                               //
8 // contains set of poits the lay on track trajectory      //
9 // according to reconstructed track parameters -          //
10 // NOT CLUSTERS POSITIONS!!!                              //
11 // Anti-Merging cut is applied only on tracks coming from //
12 // different events (that are use to fill deniminators)   //
13 //                                                        //
14 ////////////////////////////////////////////////////////////
15
16 #include <TClonesArray.h>
17 #include <TFile.h>
18 #include <TMath.h>
19
20 #include "AliESDtrack.h"
21 #include "AliTPCtrack.h"
22 #include "AliTrackReference.h"
23 #include "AliITStrackV2.h"
24
25 ClassImp(AliTrackPoints)
26
27 Int_t AliTrackPoints::fgDebug = 0;
28
29 AliTrackPoints::AliTrackPoints():
30  fN(0),
31  fX(0x0),
32  fY(0x0),
33  fZ(0x0)
34 {
35   //constructor
36 }
37 /***************************************************************/
38
39 AliTrackPoints::AliTrackPoints(AliTrackPoints::ETypes type, AliESDtrack* track, Float_t mf):
40  fN(0),
41  fX(0x0),
42  fY(0x0),
43  fZ(0x0)
44 {
45   //constructor 
46   //tupe -  what kind of track points should be calculated
47   //mf - magnetic field in [kG] = [T]*10.0
48   switch (type)
49    {
50      case kITS:
51        fN = 6;
52        fX = new Float_t[fN];
53        fY = new Float_t[fN];
54        fZ = new Float_t[fN];
55        MakeITSPoints(track);
56        break;
57
58      case kITSInnerFromVertexOuterFromTPC:
59        fN = 6;
60        fX = new Float_t[fN];
61        fY = new Float_t[fN];
62        fZ = new Float_t[fN];
63        MakeITSPointsInnerFromVertexOuterFromTPC(track,mf);
64        break;
65
66      default:
67        Info("AliTrackPoints","Not recognized type");
68    }
69    
70 }
71 /***************************************************************/
72
73 AliTrackPoints::AliTrackPoints(Int_t n, AliESDtrack* track, Float_t mf, Float_t dr, Float_t r0):
74  fN(n),
75  fX(new Float_t[fN]),
76  fY(new Float_t[fN]),
77  fZ(new Float_t[fN])
78 {
79   //constructor
80   //mf - magnetic field in kG - needed to calculate curvature out of Pt
81   //r0 - starting radius
82   //dr - calculate points every dr cm, default every 30cm
83   if (track == 0x0)
84    {
85      Error("AliTrackPoints","ESD track is null");
86      fN = 0;
87      delete [] fX;
88      delete [] fY;
89      delete [] fZ;
90      fX = fY = fZ = 0x0;
91      return;
92    }
93
94   if ( ((track->GetStatus() & AliESDtrack::kTPCrefit) == kFALSE)&& 
95        ((track->GetStatus() & AliESDtrack::kTPCin) == kFALSE)  )
96    {
97      //could happend: its stand alone tracking
98      if (GetDebug() > 3) 
99        Warning("AliTrackPoints","This ESD track does not contain TPC information");
100        
101      fN = 0;
102      delete [] fX;
103      delete [] fY;
104      delete [] fZ;
105      fX = fY = fZ = 0x0;
106      
107      return;
108    }
109   
110   Double_t x;
111   Double_t par[5];
112   track->GetInnerExternalParameters(x,par);     //get properties of the track
113   if (par[4] == 0)
114    {
115      Error("AliTrackPoints","This ESD track seem not to contain TPC information (curv is 0)");
116      return;
117    }
118
119   if (mf == 0.0)
120    {
121      Error("AliTrackPoints","Zero Magnetic field passed as parameter.");
122      return;
123    }
124    
125   Double_t alpha = track->GetInnerAlpha();
126   Double_t cc = 1000./0.299792458/mf;//conversion constant
127   Double_t c=par[4]/cc;
128   
129   MakePoints(dr,r0,x,par,c,alpha);
130   
131 }
132 /***************************************************************/
133
134 AliTrackPoints::AliTrackPoints(Int_t n, AliTPCtrack* track, Float_t dr, Float_t r0):
135  fN(n),
136  fX(new Float_t[fN]),
137  fY(new Float_t[fN]),
138  fZ(new Float_t[fN])
139 {
140   //constructor
141   //r0 starting radius
142   //dr - calculate points every dr cm, default every 30cm
143   if (track == 0x0)
144    {
145      Error("AliTrackPoints","TPC track is null");
146      fN = 0;
147      delete [] fX;
148      delete [] fY;
149      delete [] fZ;
150      fX = fY = fZ = 0x0;
151      return;
152    }
153   track->PropagateTo(r0);
154
155  //* This formation is now fixed in the following way:                         *
156  //*      external param0:   local Y-coordinate of a track (cm)                *
157  //*      external param1:   local Z-coordinate of a track (cm)                *
158  //*      external param2:   local sine of the track momentum azimuth angle    *
159  //*      external param3:   tangent of the track momentum dip angle           *
160  //*      external param4:   1/pt (1/(GeV/c))                                  *
161     
162   Double_t x = 0;
163   Double_t par[5];
164   track->GetExternalParameters(x,par);     //get properties of the track
165   
166   Double_t alpha = track->GetAlpha();
167   Double_t c=track->GetC();
168   MakePoints(dr,r0,x,par,c,alpha);
169 }  
170 /***************************************************************/
171
172 AliTrackPoints::~AliTrackPoints()
173 {
174   //destructor
175   delete [] fX;
176   delete [] fY;
177   delete [] fZ;
178 }
179 /***************************************************************/
180
181 void AliTrackPoints::MakePoints( Float_t dr, Float_t r0, Double_t x, Double_t* par, Double_t c, Double_t alpha)
182 {
183   //Calculates points starting at radius r0
184   //spacing every dr (in radial direction)
185   // according to track parameters
186   // x - position in sector local reference frame. x i parallel to R and sector is symmetric with respect to x
187   // par - track external parameters; array with 5 elements;  look at AliTPCtrack.h or AliESDtrack.h for their meaning
188   // c - track curvature
189   // alpha - sector's rotation angle (phi) == angle needed for local to global transformation
190   
191   Double_t y = par[0];
192   Double_t z0 = par[1];
193       
194   Double_t phi0local = TMath::ATan2(y,x);
195   Double_t phi0global = phi0local + alpha;
196   
197   if (phi0local<0) phi0local+=2*TMath::Pi();
198   if (phi0local>=2.*TMath::Pi()) phi0local-=2*TMath::Pi();
199   
200   if (phi0global<0) phi0global+=2*TMath::Pi();
201   if (phi0global>=2.*TMath::Pi()) phi0global-=2*TMath::Pi();
202   
203   Double_t r = TMath::Hypot(x,y);
204   
205
206   if (GetDebug() > 9) 
207     Info("AliTrackPoints","Radius0 %f, Real Radius %f",r0,r); 
208   
209   if (GetDebug() > 5) 
210     Info("AliTrackPoints","Phi Global at first padraw %f, Phi locat %f",phi0global,phi0local);
211     
212   Double_t eta = x*c - par[2] ;//par[2] = fX*C - eta; eta==fP2 ; C==fP4
213   
214   //this calculattions are assuming new (current) model 
215   Double_t tmp = par[2];
216   tmp = 1. - tmp*tmp;
217   tmp = c*y + TMath::Sqrt(tmp);
218   Double_t dca=(TMath::Hypot(eta,tmp) - 1. )/TMath::Abs(c);
219
220   //Here is old model Cold=Cnew/2.
221   Double_t dcasq = dca*dca;
222   Double_t c2 = c/2.;
223   Double_t cst1 = (1.+c2*dca)*dca;//first constant
224   Double_t cst2 = 1. + 2.*c2*dca;//second constant
225       
226   Double_t factorPhi0 = TMath::ASin((c2*r + cst1/r)/cst2);
227   Double_t factorZ0 = TMath::ASin(c2*TMath::Sqrt((r*r-dcasq)/cst2))*par[3]/c2;
228  
229   for(Int_t i = 0; i<fN; i++)
230    {
231      Double_t rc = r0 + i*dr;
232      Double_t ftmp = (c2*rc + cst1/rc)/cst2;
233      if (ftmp > 1.0) 
234       {
235         if (GetDebug() > 1) 
236           Warning("AliTrackPoints","ASin argument > 1 %f:",ftmp);
237         ftmp=1.0;
238       }
239      else if (ftmp < -1.0) 
240       {
241         if (GetDebug() > 1) 
242           Warning("AliTrackPoints","ASin argument < -1 %f:",ftmp);
243         ftmp=-1.0;
244       }
245       
246      Double_t factorPhi = TMath::ASin( ftmp );//factor phi od rc
247      Double_t phi = phi0global + factorPhi - factorPhi0;
248      
249      ftmp = (rc*rc-dcasq)/cst2;
250      if (ftmp < 0.0) 
251       {
252         if (GetDebug() > 1) 
253           Warning("AliTrackPoints","Sqrt argument < 0: %f",ftmp);
254         ftmp=0.0;
255       }
256      
257      ftmp = c2*TMath::Sqrt(ftmp);
258      if (ftmp > 1.0) 
259       {
260         if (GetDebug() > 1) 
261           Warning("AliTrackPoints","ASin argument > 1: %f",ftmp);
262         ftmp=1.0;
263       }
264      else if (ftmp < -1.0) 
265       {
266         if (GetDebug() > 1) 
267           Warning("AliTrackPoints","ASin argument < -1: %f",ftmp);
268         ftmp=-1.0;
269       }
270      Double_t factorZ = TMath::ASin(ftmp)*par[3]/c2;
271      fZ[i] = z0 + factorZ - factorZ0;
272      fX[i] = rc*TMath::Cos(phi);
273      fY[i] = rc*TMath::Sin(phi);
274      
275      if ( GetDebug() > 2 )
276       {
277         Info("AliTrackPoints","X %f Y %f Z %f R asked %f R obtained %f",
278              fX[i],fY[i],fZ[i],rc,TMath::Hypot(fX[i],fY[i]));
279       }
280    }
281 }
282 /***************************************************************/
283
284 void AliTrackPoints::MakeITSPoints(AliESDtrack* track)
285 {
286 //Calculates points in ITS
287 // z=R*Pz/Pt
288  AliITStrackV2 itstrack(*track,kTRUE);
289  Double_t x,y,z;
290  static const Double_t r[6] = {4.0, 7.0, 14.9, 23.8, 39.1, 43.6};
291  for (Int_t i = 0; i < 6; i++)
292   {
293     itstrack.GetGlobalXYZat(r[i],x,y,z);
294     fX[i] = x;
295     fY[i] = y;
296     fZ[i] = z;
297 //    Info("MakeITSPoints","X %f Y %f Z %f R asked %f R obtained %f",
298 //             fX[i],fY[i],fZ[i],r[i],TMath::Hypot(fX[i],fY[i]));
299   }   
300  
301 }
302
303 /***************************************************************/
304 void AliTrackPoints::MakeITSPointsInnerFromVertexOuterFromTPC(AliESDtrack* track, Float_t mf)
305 {
306 //makes trackpoints for ITS  
307 //for 3 inner layers calculates out of the vector at vertex
308 //for 3 outer ---------------//------------------ at inner TPC  
309
310  static const Double_t r[6] = {4.0, 7.0, 14.9, 23.8, 39.1, 43.6};
311  AliITStrackV2 itstrack(*track,kTRUE);
312  Double_t x,y,z;
313  for (Int_t i = 0; i < 3; i++)
314   {
315     itstrack.GetGlobalXYZat(r[i],x,y,z);
316     fX[i] = x;
317     fY[i] = y;
318     fZ[i] = z;
319     Info("MakeITSPoints","X %f Y %f Z %f R asked %f R obtained %f",
320              fX[i],fY[i],fZ[i],r[i],TMath::Hypot(fX[i],fY[i]));
321   }   
322  
323  for (Int_t i = 3; i < 6; i++)
324   {
325     Float_t ax,ay,az;
326     AliTrackPoints tmptp(1,track,mf,0,r[i]);
327     tmptp.PositionAt(0,ax,ay,az);
328     fX[i] = ax;
329     fY[i] = ay;
330     fZ[i] = az;
331     Info("MakeITSPoints","X %f Y %f Z %f R asked %f R obtained %f",
332              fX[i],fY[i],fZ[i],r[i],TMath::Hypot(fX[i],fY[i]));
333   }
334  
335 }
336     
337
338 /***************************************************************/
339
340 void AliTrackPoints::PositionAt(Int_t n, Float_t &x,Float_t &y,Float_t &z)
341 {
342   //returns position at point n
343   if ((n<0) || (n>fN))
344    {
345      Error("PositionAt","Point %d out of range",n);
346      return;
347    }
348   
349   x = fX[n];
350   y = fY[n];
351   z = fZ[n];
352   if ( GetDebug() > 1 )
353     {
354       Info("AliTrackPoints","n %d; X %f; Y %f; Z %f",n,x,y,z);
355     }
356 }
357 /***************************************************************/
358
359 void AliTrackPoints::Move(Float_t x, Float_t y, Float_t z)
360 {
361 //Moves all points about vector
362  for (Int_t i = 0; i<fN; i++)
363    {
364      fX[i]+=x;
365      fY[i]+=y;
366      fZ[i]+=z;
367    }   
368 }
369 /***************************************************************/
370
371 Double_t AliTrackPoints::AvarageDistance(const AliTrackPoints& tr)
372 {
373   //returns the aritmethic avarage distance between two tracks
374 //  Info("AvarageDistance","Entered");
375   if ( (fN <= 0) || (tr.fN <=0) )
376    {
377      if (GetDebug()) Warning("AvarageDistance","One of tracks is empty");
378      return -1;
379    }
380
381   if (fN != tr.fN)
382    {
383      Warning("AvarageDistance","Number of points is not equal");
384      return -1;
385    }
386    
387   Double_t sum = 0;
388   for (Int_t i = 0; i<fN; i++)
389    {
390      if (GetDebug()>9)
391       {
392 //       Float_t r1sq = fX[i]*fX[i]+fY[i]*fY[i];
393 //       Float_t r2sq = tr.fX[i]*tr.fX[i]+tr.fY[i]*tr.fY[i];
394        Float_t r1sq = TMath::Hypot(fX[i],fY[i]);
395        Float_t r2sq = TMath::Hypot(tr.fX[i],tr.fY[i]);
396        Info("AvarageDistance","radii: %f %f",r1sq,r2sq);
397       } 
398
399       
400      Double_t dx = fX[i]-tr.fX[i];
401      Double_t dy = fY[i]-tr.fY[i];
402      Double_t dz = fZ[i]-tr.fZ[i];
403      sum+=TMath::Sqrt(dx*dx + dy*dy + dz*dz);
404      
405      if (GetDebug()>1)
406       {
407        Info("AvarageDistance","Diff: x ,y z: %f , %f, %f",dx,dy,dz);
408        Info("AvarageDistance","xxyyzz %f %f %f %f %f %f",
409             fX[i],tr.fX[i],fY[i],tr.fY[i],fZ[i],tr.fZ[i]);
410       } 
411    }
412    
413   Double_t retval = sum/((Double_t)fN);
414   if ( GetDebug() )
415     {
416       Info("AvarageDistance","Avarage distance is %f.",retval);
417     }
418   return retval;
419 }
420 /***************************************************************/
421 /***************************************************************/
422 /***************************************************************/
423 /***************************************************************/
424 /***************************************************************/
425 /***************************************************************/
426
427 #include "AliRun.h"
428 #include "AliESD.h"
429 #include "AliRunLoader.h"
430 #include "AliTPCtrack.h"
431 #include "TTree.h"
432 #include "TBranch.h"
433 #include "TH2D.h"
434 #include "TCanvas.h"
435 #include "AliMagF.h"
436
437
438
439
440 void AliTrackPoints::Testesd(Int_t entr,const char* fname )
441 {
442   delete gAlice;
443   gAlice = 0x0;
444   AliRunLoader* rl = AliRunLoader::Open();
445   rl->LoadgAlice();
446   
447   Float_t mf = rl->GetAliRun()->Field()->SolenoidField();
448  
449   
450   TFile* fFile = TFile::Open(fname);
451   
452   if (fFile == 0x0)
453    {
454      printf("testesd: There is no suche a ESD file\n");
455      return;
456    }
457   AliESD* esd = dynamic_cast<AliESD*>(fFile->Get("0"));
458   AliESDtrack *t = esd->GetTrack(entr);
459   if (t == 0x0)
460    {
461      ::Error("testesd","Can not get track %d",entr);
462      return;
463    }
464
465   
466   Int_t N = 170;
467   AliTrackPoints* tp = new AliTrackPoints(N,t,mf,1.);
468   
469   Float_t xmin = -250;
470   Float_t xmax =  250;
471   
472   Float_t ymin = -250;
473   Float_t ymax =  250;
474
475   Float_t zmin = -250;
476   Float_t zmax =  250;
477   
478   TH2D* hxy = new TH2D("hxy","hxy",1000,xmin,xmax,1000,ymin,ymax);
479   TH2D* hxyt = new TH2D("hxyt","hxyt",1000,xmin,xmax,1000,ymin,ymax);
480   TH2D* hxyTR = new TH2D("hxyTR","hxyTR",1000,xmin,xmax,1000,ymin,ymax);
481
482   TH2D* hxz = new TH2D("hxz","hxz",1000,xmin,xmax,1000,zmin,zmax);
483   TH2D* hxzt = new TH2D("hxzt","hxzt",1000,xmin,xmax,1000,zmin,zmax);
484   TH2D* hxzTR = new TH2D("hxzTR","hxzTR",1000,xmin,xmax,1000,zmin,zmax);
485
486   hxyt->SetDirectory(0x0);   
487   hxy->SetDirectory(0x0);   
488   hxyTR->SetDirectory(0x0);   
489   
490   hxzt->SetDirectory(0x0);   
491   hxz->SetDirectory(0x0);   
492   hxzTR->SetDirectory(0x0);   
493
494   Float_t x,y,z;
495    
496   for (Int_t i = 0;i<N;i++)
497    {  
498       Double_t r = 84.1+i;
499       tp->PositionAt(i,x,y,z);
500       hxy->Fill(x,y);
501       hxz->Fill(x,z);
502       printf("Rdemanded %f\n",r);
503       printf("tpx %f tpy %f tpz %f Rt =%f\n", x,y,z,TMath::Hypot(x,y));
504       
505    }
506
507   rl->LoadTrackRefs();
508   TTree* treeTR = rl->TreeTR();
509   TBranch* b = treeTR->GetBranch("TPC");
510   
511   TClonesArray* trackrefs = new TClonesArray("AliTrackReference", 100);
512   AliTrackReference* tref;
513   b->SetAddress(&trackrefs);
514   
515   Int_t tlab = TMath::Abs(t->GetLabel());
516
517   Int_t netr = (Int_t)treeTR->GetEntries();
518   printf("Found %d entries in TR tree\n",netr);
519   
520   for (Int_t e = 0; e < netr; e++)
521    { 
522      treeTR->GetEntry(e);
523      tref = (AliTrackReference*)trackrefs->At(0);
524      if (tref == 0x0) continue;
525      if (tref->GetTrack() != tlab) continue;
526     
527      printf("Found %d entries in TR array\n",trackrefs->GetEntries());
528      
529      for (Int_t i = 0; i < trackrefs->GetEntries(); i++)
530       {
531         tref = (AliTrackReference*)trackrefs->At(i);
532         if (tref->GetTrack() != tlab) continue;
533         x = tref->X();
534         y = tref->Y();
535         z = tref->Z();
536         printf("Track Ref: x %f y %f z %f\n",tref->X(),tref->Y(),tref->Z());
537
538         hxzTR->Fill(x,z);
539         hxyTR->Fill(x,y);
540         for (Int_t j = 1; j < 10; j++)
541          {
542             hxyTR->Fill(x, y+j*0.1);
543             hxyTR->Fill(x, y-j*0.1);
544             hxyTR->Fill(x+j*0.1,y);
545             hxyTR->Fill(x-j*0.1,y);
546
547             hxzTR->Fill(x,z-j*0.1);
548             hxzTR->Fill(x,z+j*0.1);
549             hxzTR->Fill(x-j*0.1,z);
550             hxzTR->Fill(x+j*0.1,z);
551          }
552       }
553       break; 
554     }  
555   hxy->Draw("");
556 //  hxzt->Draw("same");
557   hxyTR->Draw("same");
558   
559   delete rl;
560 }
561
562 /***************************************************************/
563 /***************************************************************/
564 /***************************************************************/
565
566 void AliTrackPoints::Testtpc(Int_t entr)
567 {
568   delete gAlice;
569   gAlice = 0x0;
570   AliRunLoader* rl = AliRunLoader::Open();
571   AliLoader* l = rl->GetLoader("TPCLoader");
572   rl->LoadgAlice();
573   AliKalmanTrack::SetConvConst(100/0.299792458/0.2/rl->GetAliRun()->Field()->Factor());
574   l->LoadTracks();
575   AliTPCtrack* t = new AliTPCtrack();
576   TBranch* b=l->TreeT()->GetBranch("tracks");
577   b->SetAddress(&t);
578   l->TreeT()->GetEntry(entr);  
579   Int_t N = 160;
580   AliTrackPoints* tp = new AliTrackPoints(N,t,1.);
581   
582   Float_t xmin = -250;
583   Float_t xmax =  250;
584   
585   Float_t ymin = -250;
586   Float_t ymax =  250;
587
588   Float_t zmin = -250;
589   Float_t zmax =  250;
590   
591   TH2D* hxy = new TH2D("hxy","hxy",1000,xmin,xmax,1000,ymin,ymax);
592   TH2D* hxyt = new TH2D("hxyt","hxyt",1000,xmin,xmax,1000,ymin,ymax);
593   TH2D* hxyTR = new TH2D("hxyTR","hxyTR",1000,xmin,xmax,1000,ymin,ymax);
594
595   TH2D* hxz = new TH2D("hxz","hxz",1000,xmin,xmax,1000,zmin,zmax);
596   TH2D* hxzt = new TH2D("hxzt","hxzt",1000,xmin,xmax,1000,zmin,zmax);
597   TH2D* hxzTR = new TH2D("hxzTR","hxzTR",1000,xmin,xmax,1000,zmin,zmax);
598
599   hxyt->SetDirectory(0x0);   
600   hxy->SetDirectory(0x0);   
601   hxyTR->SetDirectory(0x0);   
602   
603   hxzt->SetDirectory(0x0);   
604   hxz->SetDirectory(0x0);   
605   hxzTR->SetDirectory(0x0);   
606
607   Float_t x,y,z;
608    
609   for (Int_t i = 0;i<N;i++)
610    {  
611       Double_t r = 84.1+i;
612       tp->PositionAt(i,x,y,z);
613       hxy->Fill(x,y);
614       hxz->Fill(x,z);
615       printf("Rdemanded %f\n",r);
616       printf("tpx %f tpy %f tpz %f Rt =%f\n", x,y,z,TMath::Hypot(x,y));
617       
618       //BUT they are local!!!!
619       t->PropagateTo(r);
620 //      Double_t phi = t->Phi();
621       Double_t rl = TMath::Hypot(t->GetX(),t->GetY());//real radius 
622       
623       Double_t alpha = t->GetAlpha();
624       Double_t salpha = TMath::Sin(alpha);
625       Double_t calpha = TMath::Cos(alpha);
626       x = t->GetX()*calpha - t->GetY()*salpha;
627       y = t->GetX()*salpha + t->GetY()*calpha;
628       z = t->GetZ();
629       
630       printf("tx %f ty %f tz %f Rt = %f R from XY %f\n",x,y,z,TMath::Hypot(x,y),rl);
631       
632       printf("tpz - tz %f\n",z-t->GetZ());
633       printf("\n");
634       hxyt->Fill(x,y);
635       hxzt->Fill(x,z);
636       
637    }
638
639   rl->LoadTrackRefs();
640   TTree* treeTR = rl->TreeTR();
641   b = treeTR->GetBranch("TPC");
642   
643   TClonesArray* trackrefs = new TClonesArray("AliTrackReference", 100);
644   AliTrackReference* tref;
645   b->SetAddress(&trackrefs);
646   
647   Int_t tlab = TMath::Abs(t->GetLabel());
648
649   Int_t netr = (Int_t)treeTR->GetEntries();
650   printf("Found %d entries in TR tree\n",netr);
651   
652   for (Int_t e = 0; e < netr; e++)
653    { 
654      treeTR->GetEntry(e);
655      tref = (AliTrackReference*)trackrefs->At(0);
656      if (tref == 0x0) continue;
657      if (tref->GetTrack() != tlab) continue;
658     
659      printf("Found %d entries in TR array\n",trackrefs->GetEntries());
660      
661      for (Int_t i = 0; i < trackrefs->GetEntries(); i++)
662       {
663         tref = (AliTrackReference*)trackrefs->At(i);
664         if (tref->GetTrack() != tlab) continue;
665         x = tref->X();
666         y = tref->Y();
667         z = tref->Z();
668         printf("Track Ref: x %f y %f z %f\n",tref->X(),tref->Y(),tref->Z());
669
670         hxzTR->Fill(x,z);
671         hxyTR->Fill(x,y);
672         for (Int_t j = 1; j < 10; j++)
673          {
674             hxyTR->Fill(x, y+j*0.1);
675             hxyTR->Fill(x, y-j*0.1);
676             hxyTR->Fill(x+j*0.1,y);
677             hxyTR->Fill(x-j*0.1,y);
678
679             hxzTR->Fill(x,z-j*0.1);
680             hxzTR->Fill(x,z+j*0.1);
681             hxzTR->Fill(x-j*0.1,z);
682             hxzTR->Fill(x+j*0.1,z);
683          }
684       }
685       break; 
686     }  
687   hxz->Draw("");
688 //  hxzt->Draw("same");
689   hxzTR->Draw("same");
690   
691   delete rl;
692 }
693