]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/util/AliHLTGlobalBarrelTrack.cxx
f81f3f112afe646563538d2eac7a8096acfd5e6d
[u/mrichter/AliRoot.git] / HLT / BASE / util / AliHLTGlobalBarrelTrack.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTGlobalBarrelTrack.cxx
20     @author Matthias Richter
21     @date   2009-06-24
22     @brief  An AliKalmanTrack implementation for global HLT barrel tracks.
23 */
24
25 // see header file for class documentation
26 // or
27 // refer to README to build package
28 // or
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30
31 #include <cassert>
32 #include <memory>
33 #include <iostream>
34 #include "AliHLTGlobalBarrelTrack.h"
35 #include "AliHLTSpacePointContainer.h"
36 #include "AliHLTTrackGeometry.h"
37 #include "AliHLTMisc.h"
38 #include "TClonesArray.h"
39 #include "TFile.h"
40 #include "TArrayC.h"
41 #include "TMath.h"
42 #include "TMarker.h"
43 #include "TArc.h"
44
45 /** ROOT macro for the implementation of ROOT specific class methods */
46 ClassImp(AliHLTGlobalBarrelTrack)
47
48 AliHLTGlobalBarrelTrack::AliHLTGlobalBarrelTrack()
49 : AliKalmanTrack()
50   , fPoints()
51   , fLastX(0.0)
52   , fLastY(0.0)
53   , fTrackID(-1)
54   , fHelixRadius(0.0)
55   , fHelixCenterX(0.0)
56   , fHelixCenterY(0.0)
57   , fSpacePoints(NULL)
58   , fTrackPoints(NULL)
59 {
60   // see header file for class documentation
61   // or
62   // refer to README to build package
63   // or
64   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
65 }
66
67 AliHLTGlobalBarrelTrack::AliHLTGlobalBarrelTrack(const AliHLTGlobalBarrelTrack& t)
68   : AliKalmanTrack(t)
69   , fPoints()
70   , fLastX(t.GetLastPointX())
71   , fLastY(t.GetLastPointY())
72   , fTrackID(t.TrackID())
73   , fHelixRadius(t.fHelixRadius)
74   , fHelixCenterX(t.fHelixCenterX)
75   , fHelixCenterY(t.fHelixCenterY)
76   , fSpacePoints(NULL)
77   , fTrackPoints(NULL)
78 {
79   // see header file for class documentation
80   fPoints.assign(t.fPoints.begin(), t.fPoints.end());
81 }
82
83 AliHLTGlobalBarrelTrack::AliHLTGlobalBarrelTrack(const AliHLTExternalTrackParam& p )
84   : AliKalmanTrack()
85   , fPoints()
86   , fLastX(p.fLastX)
87   , fLastY(p.fLastY)
88   , fTrackID(p.fTrackID)
89   , fHelixRadius(0.0)
90   , fHelixCenterX(0.0)
91   , fHelixCenterY(0.0)
92   , fSpacePoints(NULL)
93   , fTrackPoints(NULL)
94 {
95   // see header file for class documentation
96
97   // the 5 track parameters are named in the AliHLTExternalTrackParam
98   // while AliExternalTrackParam just uses an array[5]
99   // the members have the same order, fY is the first one
100   Set(p.fX, p.fAlpha, &p.fY, p.fC);
101   SetPoints(p.fPointIDs, p.fNPoints);
102   SetNumberOfClusters(p.fNPoints);
103   //SetIntegratedLength(GetPathLengthTo( GetLastPointX(), b);
104 }
105
106 AliHLTGlobalBarrelTrack::AliHLTGlobalBarrelTrack(const AliExternalTrackParam& p )
107   : AliKalmanTrack()
108   , fPoints()
109   , fLastX(0)
110   , fLastY(0)
111   , fTrackID(0)
112   , fHelixRadius(0.0)
113   , fHelixCenterX(0.0)
114   , fHelixCenterY(0.0)
115   , fSpacePoints(NULL)
116   , fTrackPoints(NULL)
117 {
118   // see header file for class documentation
119   *(dynamic_cast<AliExternalTrackParam*>(this))=p;
120 }
121
122 AliHLTGlobalBarrelTrack::~AliHLTGlobalBarrelTrack()
123 {
124   // see header file for class documentation
125 }
126
127
128 Double_t AliHLTGlobalBarrelTrack::GetPathLengthTo( Double_t x, Double_t b ) const
129 {
130   // calculate the trajectory length for dx step
131     
132   Double_t dx = x - GetX();
133   Double_t ey = GetSnp();
134   if( TMath::Abs( ey )>=kAlmost1 ) return 0;
135
136   Double_t ex = TMath::Sqrt(1-ey*ey);
137   Double_t k  = GetC(b);
138
139   Double_t ey1 = k * dx + ey;
140
141   // check for intersection with X=x
142
143   if ( TMath::Abs( ey1 ) >= kAlmost1  ) return 0;
144
145   Double_t ex1 = TMath::Sqrt(1-ey1*ey1);
146
147   Double_t ss = ey + ey1;
148   Double_t cc = ex + ex1;
149
150   if ( TMath::Abs( cc ) < 1.e-4  ) return 0;
151
152   Double_t tg = ss / cc; 
153   Double_t dl = dx * TMath::Sqrt( 1 + tg * tg );
154   Double_t dSin = dl * k / 2;
155   if ( dSin > 1 ) dSin = 1;
156   if ( dSin < -1 ) dSin = -1;
157   Double_t dS = ( TMath::Abs( k ) > 1.e-4 )  ? ( 2 * TMath::ASin( dSin ) / k ) : dl;
158
159   return dS*TMath::Sqrt(1 + GetTgl()*GetTgl() );
160 }
161
162
163
164
165 int AliHLTGlobalBarrelTrack::ConvertTrackDataArray(const AliHLTTracksData* pTracks, unsigned sizeInByte, vector<AliHLTGlobalBarrelTrack> &tgtArray )
166 {
167   // see header file for class documentation
168   int iResult=0;
169   if (!pTracks || sizeInByte<sizeof(AliHLTTracksData) || pTracks->fCount==0) return 0;
170
171   const AliHLTUInt8_t* pEnd=reinterpret_cast<const AliHLTUInt8_t*>(pTracks);
172   pEnd+=sizeInByte;
173
174   tgtArray.resize(pTracks->fCount + tgtArray.size());
175   const AliHLTUInt8_t* pCurrent=reinterpret_cast<const AliHLTUInt8_t*>(pTracks->fTracklets);
176   for (unsigned i=0; i<pTracks->fCount; i++) {
177     if (pCurrent+sizeof(AliHLTExternalTrackParam)>pEnd) {
178       iResult=-EINVAL; break;
179     }
180     const AliHLTExternalTrackParam* track=reinterpret_cast<const AliHLTExternalTrackParam*>(pCurrent);
181     if (pCurrent+sizeof(AliHLTExternalTrackParam)+track->fNPoints*sizeof(UInt_t)>pEnd) {
182       iResult=-EINVAL; break;
183     }
184     tgtArray[i]=*track;
185     pCurrent+=sizeof(AliHLTExternalTrackParam)+track->fNPoints*sizeof(UInt_t);
186   }
187   if (iResult<0) tgtArray.clear();
188   else iResult=tgtArray.size();
189   return iResult;
190 }
191
192 UInt_t AliHLTGlobalBarrelTrack::GetNumberOfPoints() const
193 {
194   // see header file for class documentation
195   return fPoints.size();
196 }
197
198 const UInt_t* AliHLTGlobalBarrelTrack::GetPoints() const
199 {
200   // see header file for class documentation
201   if (fPoints.size()==0) return NULL;
202   return &fPoints[0];
203 }
204
205 int AliHLTGlobalBarrelTrack::SetPoints(const UInt_t* pArray, UInt_t arraySize)
206 {
207   // see header file for class documentation
208   if (!pArray || arraySize==0) return 0;
209   fPoints.resize(arraySize);
210   for (unsigned i=0; i<arraySize; i++) fPoints[i]=pArray[i];
211   return fPoints.size();
212 }
213
214 int AliHLTGlobalBarrelTrack::CalculateHelixParams()
215 {
216   // calculate radius and center of the helix
217   float bfield=AliHLTMisc::Instance().GetBz();
218   if (TMath::Abs(bfield)<kAlmost0) {
219     // no magnetic field -> straight lines
220     fHelixRadius=kVeryBig;
221     fHelixCenterX=kVeryBig;
222     fHelixCenterY=kVeryBig;
223   } else {
224     fHelixRadius = GetSignedPt()/(-kB2C*bfield);
225     Double_t trackPhi = Phi()-GetAlpha();
226
227     //cout << "Helix: phi=" << trackPhi << " x=" << GetX() << " y=" << GetY() << endl;
228     fHelixCenterX = GetX() + fHelixRadius *  sin(trackPhi);
229     fHelixCenterY = GetY() - fHelixRadius *  cos(trackPhi);
230     //cout << "Helix: center" << " x=" << fHelixCenterX << " y=" << fHelixCenterY << endl;
231   }
232   return 0;
233 }
234
235 int AliHLTGlobalBarrelTrack::CalculateCrossingPoint(float xPlane, float alphaPlane, float& y, float& z)
236 {
237   // calculate crossing point of helix with a plane in yz
238   // in the local coordinates of the plane
239   int iResult=0;
240   if (TMath::Abs(fHelixRadius)<kAlmost0 &&
241       (iResult=CalculateHelixParams())<0) {
242     return iResult;
243   }
244
245   if (TMath::Abs(fHelixRadius)>=kVeryBig) {
246     // no magnetic field -> straight lines
247   } else {
248     // rotate helix center to local coordinates of the plane reference frame
249     float cosa=TMath::Cos(alphaPlane-GetAlpha());
250     float sina=TMath::Sin(alphaPlane-GetAlpha());
251     float cx= fHelixCenterX * cosa + fHelixCenterY * sina;
252     float cy=-fHelixCenterX * sina + fHelixCenterY * cosa;
253
254     // crossing point of helix with plane
255     Double_t aa = (cx - xPlane)*(cx - xPlane);
256     Double_t r2 = fHelixRadius*fHelixRadius;
257     if(aa > r2) // no crossing
258       return 0;
259
260     Double_t aa2 = sqrt(r2 - aa);
261     Double_t y1 = cy + aa2;
262     Double_t y2 = cy - aa2;
263     y = y1;
264     if(TMath::Abs(y2) < TMath::Abs(y1)) y = y2;
265   
266     // calculate the arc length between (x,y) and (x0,y0) with radius (cx,cy)
267     // reference point is (x0,y0) rotated by the diffence of the plane angle and
268     // the track reference frame angle alpha
269     // 1) angle of (x,y)
270     Double_t angle1 = atan2((y - cy),(xPlane - cx));
271     if(angle1 < 0) angle1 += TMath::TwoPi();
272
273     // 2) angle of (x0,y0)
274     float x0= GetX() * cosa + GetY() * sina;
275     float y0=-GetX() * sina + GetY() * cosa;
276     Double_t angle2 = atan2((y0 - cy),(x0 - cx));
277     if(angle2 < 0) angle2 += TMath::TwoPi();
278
279     // 3) angle between (x,y) and (x0,y0)
280     Double_t diffangle = angle1 - angle2;
281     diffangle = fmod(diffangle,TMath::TwoPi());
282
283     // 4) arc length
284     Double_t arclength = TMath::Abs(diffangle*fHelixRadius);
285
286     // 5) direction depending on whether going outwards or inwards
287     int direction=GetX()>xPlane?-1:1;
288     z = GetZ() + direction*arclength*GetTgl();
289
290     //cout << "x=" << xPlane << " y=" << y << " cx=" << cx << " cy=" << cy << " a1=" << angle1 << " a2=" << angle2 << " diffa=" << diffangle << " s=" << arclength << " z=" << z << endl;
291   }
292   return 1;
293 }
294
295 void AliHLTGlobalBarrelTrack::Print(Option_t* option) const
296 {
297   // see header file for class documentation
298   cout << "********* Track Id: " << fTrackID << " *******************" << endl;
299   AliExternalTrackParam::Print(option);
300 //   cout << "  Alpha "     << GetAlpha();
301 //   cout << "  X "         << GetX();
302 //   cout << "  Y "         << GetY();
303 //   cout << "  Z "         << GetZ() << endl;
304 //   cout << "  Snp "       << GetSnp();
305 //   cout << "  Tgl "       << GetTgl();
306 //   cout << "  Signed1Pt " << GetSigned1Pt() << endl;
307 }
308
309 void AliHLTGlobalBarrelTrack::Draw(Option_t *option)
310 {
311   /// Inherited from TObject, draw the track
312   float scale=250;
313   float center[2]={0.5,0.5};
314
315   if (TMath::Abs(fHelixRadius)<kAlmost0 &&
316       (CalculateHelixParams())<0) {
317     return;
318   }
319
320   TString strOption(option);
321   if (strOption.IsNull()) strOption="spacepoints trackarc";
322   std::auto_ptr<TObjArray> tokens(strOption.Tokenize(" "));
323   if (!tokens.get()) return;
324   for (int i=0; i<tokens->GetEntriesFast(); i++) {
325     if (!tokens->At(i)) continue;
326     const char* key="";
327     TString arg=tokens->At(i)->GetName();
328
329     key="scale=";
330     if (arg.BeginsWith(key)) {
331       arg.ReplaceAll(key, "");
332       scale=arg.Atof();
333       continue;
334     }
335     key="centerx=";
336     if (arg.BeginsWith(key)) {
337       arg.ReplaceAll(key, "");
338       center[0]=arg.Atof();
339       continue;
340     }
341     key="centery=";
342     if (arg.BeginsWith(key)) {
343       arg.ReplaceAll(key, "");
344       center[1]=arg.Atof();
345       continue;
346     }
347     key="spacepoints";
348     if (arg.CompareTo(key)==0) {
349       if (fSpacePoints) DrawProjXYSpacePoints(option, fSpacePoints, scale, center);
350       continue;
351     }
352     key="trackarc";
353     if (arg.CompareTo(key)==0) {
354       DrawProjXYTrack(option, scale, center);
355       continue;
356     }
357   }
358 }
359
360 int AliHLTGlobalBarrelTrack::DrawProjXYSpacePoints(Option_t */*option*/, const AliHLTSpacePointContainer* spacePoints, const float scale, float center[2])
361 {
362   /// draw space points
363   int markerColor=3;
364
365   if (!spacePoints) return -EINVAL;
366
367   const UInt_t* pointids=GetPoints();
368   for (unsigned i=0; i<GetNumberOfPoints() && pointids; i++) {
369     float clusterphi=spacePoints->GetPhi(pointids[i]);
370     float cosphi=TMath::Cos(clusterphi);
371     float sinphi=TMath::Sin(clusterphi);
372     float clusterx=spacePoints->GetX(pointids[i]);
373     float clustery=spacePoints->GetY(pointids[i]);
374     // rotate
375     float pointx= clusterx*sinphi + clustery*cosphi;
376     float pointy=-clusterx*cosphi + clustery*sinphi;
377
378     // FIXME: cleanup of marker objects
379     TMarker* m=new TMarker(pointx/(2*scale)+center[0], pointy/(2*scale)+center[1], 3);
380     m->SetMarkerColor(markerColor);
381     m->Draw("same");
382   }
383   return 0;
384 }
385
386 // FIXME: make this a general geometry definition
387 // through an abstract class interface
388 const Double_t gkTPCX[159] = {
389   85.195,
390   85.945,
391   86.695,
392   87.445,
393   88.195,
394   88.945,
395   89.695,
396   90.445,
397   91.195,
398   91.945,
399   92.695,
400   93.445,
401   94.195,
402   94.945,
403   95.695,
404   96.445,
405   97.195,
406   97.945,
407   98.695,
408   99.445,
409   100.195,
410   100.945,
411   101.695,
412   102.445,
413   103.195,
414   103.945,
415   104.695,
416   105.445,
417   106.195,
418   106.945,
419   107.695,
420   108.445,
421   109.195,
422   109.945,
423   110.695,
424   111.445,
425   112.195,
426   112.945,
427   113.695,
428   114.445,
429   115.195,
430   115.945,
431   116.695,
432   117.445,
433   118.195,
434   118.945,
435   119.695,
436   120.445,
437   121.195,
438   121.945,
439   122.695,
440   123.445,
441   124.195,
442   124.945,
443   125.695,
444   126.445,
445   127.195,
446   127.945,
447   128.695,
448   129.445,
449   130.195,
450   130.945,
451   131.695,
452   135.180,
453   136.180,
454   137.180,
455   138.180,
456   139.180,
457   140.180,
458   141.180,
459   142.180,
460   143.180,
461   144.180,
462   145.180,
463   146.180,
464   147.180,
465   148.180,
466   149.180,
467   150.180,
468   151.180,
469   152.180,
470   153.180,
471   154.180,
472   155.180,
473   156.180,
474   157.180,
475   158.180,
476   159.180,
477   160.180,
478   161.180,
479   162.180,
480   163.180,
481   164.180,
482   165.180,
483   166.180,
484   167.180,
485   168.180,
486   169.180,
487   170.180,
488   171.180,
489   172.180,
490   173.180,
491   174.180,
492   175.180,
493   176.180,
494   177.180,
495   178.180,
496   179.180,
497   180.180,
498   181.180,
499   182.180,
500   183.180,
501   184.180,
502   185.180,
503   186.180,
504   187.180,
505   188.180,
506   189.180,
507   190.180,
508   191.180,
509   192.180,
510   193.180,
511   194.180,
512   195.180,
513   196.180,
514   197.180,
515   198.180,
516   199.430,
517   200.930,
518   202.430,
519   203.930,
520   205.430,
521   206.930,
522   208.430,
523   209.930,
524   211.430,
525   212.930,
526   214.430,
527   215.930,
528   217.430,
529   218.930,
530   220.430,
531   221.930,
532   223.430,
533   224.930,
534   226.430,
535   227.930,
536   229.430,
537   230.930,
538   232.430,
539   233.930,
540   235.430,
541   236.930,
542   238.430,
543   239.930,
544   241.430,
545   242.930,
546   244.430,
547   245.930
548 };
549
550 int AliHLTGlobalBarrelTrack::DrawProjXYTrack(Option_t *option, const float scale, float center[2])
551 {
552   /// draw track
553   bool bDrawArc=false; // draw TArc
554   bool bNoTrackPoints=false; // don't draw track points
555   TString strOption(option);
556   if (strOption.IsNull()) strOption="spacepoints trackarc";
557   std::auto_ptr<TObjArray> tokens(strOption.Tokenize(" "));
558   if (!tokens.get()) return 0;
559   for (int i=0; i<tokens->GetEntriesFast(); i++) {
560     if (!tokens->At(i)) continue;
561     const char* key="";
562     TString arg=tokens->At(i)->GetName();
563
564     key="drawarc";
565     if (arg.BeginsWith(key)) {
566       bDrawArc=true;
567       continue;
568     }
569     key="notrackpoints";
570     if (arg.BeginsWith(key)) {
571       bNoTrackPoints=true;
572       continue;
573     }
574   }
575
576   float cosa=TMath::Cos(GetAlpha());
577   float sina=TMath::Sin(GetAlpha());
578
579   // first point
580   float firstpoint[2];
581   firstpoint[0]= GetX()*sina + GetY()*cosa;
582   firstpoint[1]=-GetX()*cosa + GetY()*sina;
583   {
584     //cout << " first point alpha=" << GetAlpha() << " x: " << firstpoint[0] << " y: " << firstpoint[1] << endl;
585     // FIXME: cleanup of marker objects
586     TMarker* m=new TMarker(firstpoint[0]/(2*scale)+center[0], firstpoint[1]/(2*scale)+center[1], 29);
587     m->SetMarkerSize(2);
588     m->SetMarkerColor(2);
589     m->Draw("same");
590   }
591
592   // draw points in step width and remember the last point
593   float lastpoint[2]={0.0, 0.0};
594   int firstpadrow=0;
595   for (; firstpadrow<159 && gkTPCX[firstpadrow]<GetX(); firstpadrow++);
596   DrawProjXYTrackPoints(option, scale, center, firstpadrow, -1, firstpoint);
597   DrawProjXYTrackPoints(option, scale, center, firstpadrow, 1, lastpoint);
598
599   if (bDrawArc) {
600     if (TMath::Abs(fHelixRadius)>=kVeryBig) {
601       // no magnetic field -> straight lines
602     } else {
603       // rotate helix center to local coordinates of the plane reference frame
604       float cx= fHelixCenterX * sina + fHelixCenterY * cosa;
605       float cy=-fHelixCenterX * cosa + fHelixCenterY * sina;
606
607       float diffx=cx-firstpoint[0];
608       float diffy=cy-firstpoint[1];
609       float phimin=0.0;
610       float phimax=2*TMath::Pi();
611       if (TMath::Abs(diffx)<kAlmost0) {
612         phimin=TMath::Pi()/2;
613       } else {
614         phimin=TMath::ATan(diffy/diffx);
615       }
616       if (diffx>0) phimin+=TMath::Pi();
617
618       diffx=cx-lastpoint[0];
619       diffy=cy-lastpoint[1];
620       if (TMath::Abs(diffx)<kAlmost0) {
621         phimax=TMath::Pi()/2;
622       } else {
623         phimax=TMath::ATan(diffy/diffx);
624       }
625       //cout << "diffx=" << diffx << " diffy=" << diffy << " phimin=" << phimin << " phimax=" << phimax << endl;
626       if (diffx>0) phimax+=TMath::Pi();
627       if (phimax<0 && phimin>=0 && 
628           phimax+TMath::TwoPi()-phimin<TMath::Pi()) phimax+=TMath::TwoPi();
629       //if (phimax<0 && TMath::Abs(phimax-phimin)<TMath::Pi()) phimax+=TMath::TwoPi();
630       if (phimax-phimin>TMath::Pi()) phimax-=TMath::TwoPi();
631
632       if (0/*phimin>phimax*/) {
633         float tmp=phimin;
634         phimin=phimax;
635         phimax=tmp;
636       }
637       phimin*=360.0/(2*TMath::Pi());
638       phimax*=360.0/(2*TMath::Pi());
639       //cout << " cx=" << cx << " cy=" << cy << " r=" << fHelixRadius << " phimin=" << phimin << " phimax=" << phimax << endl;
640       // FIXME: cleanup of graphics objects
641       TArc* tarc=new TArc(cx/(2*scale)+center[0], cy/(2*scale)+center[1], TMath::Abs(fHelixRadius)/(2*scale), phimin, phimax);
642       tarc->SetNoEdges();
643       tarc->SetFillStyle(0);
644       tarc->Draw("same");
645     }
646   }
647   return 0;
648 }
649
650 int AliHLTGlobalBarrelTrack::DrawProjXYTrackPoints(Option_t */*option*/, const float scale, const float center[2], int firstpadrow, int step, float point[2])
651 {
652   // draw points in step width and return the last point
653   float offsetAlpha=0.0;
654   float cosa=TMath::Cos(GetAlpha());
655   float sina=TMath::Sin(GetAlpha());
656   int markerColor=1;
657   for (int padrow=firstpadrow; padrow>=0 && padrow<159; padrow+=step) {
658     float x=gkTPCX[padrow];
659     float y=0.0;
660     float z=0.0;
661
662     int maxshift=9;
663     int shift=0;
664     int result=0;
665     do {
666       if ((result=CalculateCrossingPoint(x, GetAlpha()-offsetAlpha, y, z))<1) break;
667       float pointAlpha=TMath::ATan(y/x);
668       if (TMath::Abs(pointAlpha)>TMath::Pi()/18) {
669         offsetAlpha+=(pointAlpha>0?-1:1)*TMath::Pi()/9;
670         result=0;
671         markerColor++;
672         cosa=TMath::Cos(GetAlpha()-offsetAlpha);
673         sina=TMath::Sin(GetAlpha()-offsetAlpha);
674       }
675     } while (result==0 && shift++<maxshift);
676     if (result<1) continue;
677     point[0]= x*sina + y*cosa;
678     point[1]=-x*cosa + y*sina;
679
680     //cout << x << " : x=" << x << " y=" << y << endl;
681     // FIXME: cleanup of TMarker objects?
682     TMarker* m=new TMarker(point[0]/(2*scale)+center[0], point[1]/(2*scale)+center[1], z>=0?2:5);
683     m->SetMarkerColor(markerColor);
684     m->Draw("same");
685   }
686   return 0;
687 }
688
689 void AliHLTGlobalBarrelTrack::SetTrackGeometry(AliHLTTrackGeometry* points)
690 {
691   /// set the instance to the track points container
692   fTrackPoints=points; 
693   if (!fTrackPoints) return; 
694   fTrackPoints->SetTrackId(GetID());
695 }
696
697 int AliHLTGlobalBarrelTrack::AssociateSpacePoints(AliHLTTrackGeometry* trackpoints, AliHLTSpacePointContainer& spacepoints) const
698 {
699   /// associate the track space points to the calculated track points
700   AliHLTTrackGeometry* instance=trackpoints;
701   if (!instance) instance=fTrackPoints;
702   if (!instance) return 0;
703
704   UInt_t nofIds=GetNumberOfPoints();
705   const UInt_t* ids=GetPoints();
706   int result=instance->AssociateSpacePoints(ids, nofIds, spacepoints);
707   return result;
708 }
709
710 int AliHLTGlobalBarrelTrack::ReadTracks(const char* filename, TClonesArray& tgt, AliHLTComponentDataType /*dt*/, unsigned /*specification*/)
711 {
712   // open block from file and add to collection
713   if (!filename) return -EINVAL;
714   
715   TString input=filename;
716   input+="?filetype=raw";
717   std::auto_ptr<TFile> pFile(new TFile(input));
718   if (!pFile.get()) return -ENOMEM;
719   if (pFile->IsZombie()) return -ENOENT;
720
721   int iResult=0;
722   pFile->Seek(0);
723   std::auto_ptr<TArrayC> buffer(new TArrayC);
724   if (!buffer.get()) return -ENOMEM;
725
726   buffer->Set(pFile->GetSize());
727   if (pFile->ReadBuffer(buffer->GetArray(), buffer->GetSize())==0) {
728     const AliHLTTracksData* pTracks=reinterpret_cast<const AliHLTTracksData*>(buffer->GetArray());
729     vector<AliHLTGlobalBarrelTrack> tracks;
730     iResult=ConvertTrackDataArray(pTracks, buffer->GetSize(), tracks);
731     if (iResult>=0) {
732       int offset=tgt.GetEntriesFast();
733       tgt.ExpandCreate(offset+tracks.size());
734       for (unsigned i=0; i<tracks.size(); i++) {
735         new (tgt[offset+i]) AliHLTGlobalBarrelTrack(tracks[i]);
736       }
737       iResult=tracks.size();
738     } else {
739       //HLTError("failed to convert tracks from file %s size %d byte(s) ", filename, pFile->GetSize());
740     }
741   } else {
742     //HLTError("failed reading %d byte(s) from file %s", pFile->GetSize(), filename);
743     iResult=-ENODATA;
744   }
745
746   return iResult;
747 }
748
749 int AliHLTGlobalBarrelTrack::ReadTrackList(const char* listfile, TClonesArray& tgt, AliHLTComponentDataType dt, unsigned specification)
750 {
751   // open blank separated list of files and read tracks
752   ifstream list(listfile);
753   if (!list.good()) return -ENOENT;
754
755   int count=0;
756   TString file;
757   while (file.ReadLine(list)) {
758     //HLTInfo("adding tracks from file %s", file.Data());
759     int iResult=ReadTracks(file.Data(), tgt, dt, specification);
760     if (iResult<0) {
761       //HLTInfo("failed to add data from file %s: error %d", file.Data(), iResult);
762       return iResult;
763     }
764     count+=iResult;
765   }
766
767   return count;
768 }