]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/hough/AliL3HoughDisplay.cxx
Seems to be working properly now.
[u/mrichter/AliRoot.git] / HLT / hough / AliL3HoughDisplay.cxx
1 //$Id$
2
3 // Author: Anders Vestbo <mailto:vestbo@fi.uib.no>
4 //*-- Copyright &copy ASV 
5
6 #include "AliL3StandardIncludes.h"
7
8 #include <TCanvas.h>
9 #include <TView.h>
10 #include <TPolyMarker3D.h>
11 #include <TPolyLine3D.h>
12 #include <TNode.h>
13 #include <TGeometry.h>
14 #include <TShape.h>
15 #include <TFile.h>
16
17 #include "AliL3Logging.h"
18 #include "AliL3Transform.h"
19 #include "AliL3HoughTrack.h"
20 #include "AliL3TrackArray.h"
21 #include "AliL3MemHandler.h"
22 #include "AliL3HoughDisplay.h"
23
24 //_____________________________________________________________
25 // Display class for Hough transform code
26
27 ClassImp(AliL3HoughDisplay)
28
29
30 AliL3HoughDisplay::AliL3HoughDisplay()
31 {
32
33   fTracks = 0;
34   fDigitRowData = 0;
35   fNDigitRowData = 0;
36   fShowSlice = -1;
37   fPatch = -1;
38 }
39
40 AliL3HoughDisplay::~AliL3HoughDisplay()
41 {
42   if(fTracks)
43     delete fTracks;
44 }
45
46 void AliL3HoughDisplay::Init(Char_t *trackfile, Char_t *gfile)
47 {
48   TFile *file = TFile::Open(gfile);
49   if(!file->IsOpen())
50     cerr<<"AliL3HoughDisplay::AliL3HoughDisplay : Geometry file " << gfile << " does not exist"<<endl;
51   fGeom = (TGeometry*)file->Get("AliceGeom");
52   file->Close();
53   
54   fTracks = new AliL3TrackArray();
55   AliL3MemHandler *tfile = new AliL3MemHandler();
56   tfile->SetBinaryInput(trackfile);
57   tfile->Binary2TrackArray(fTracks);
58   tfile->CloseBinaryInput();
59   delete tfile;
60 }
61
62 void AliL3HoughDisplay::GenerateHits(AliL3Track *track,Float_t *x,Float_t *y,Float_t *z,Int_t &n)
63 {
64   n=0;
65   Float_t xyz[3];
66   for(Int_t i=AliL3Transform::GetFirstRow(0); i<AliL3Transform::GetLastRow(5); i++)
67     {
68       if(track->GetCrossingPoint(i,xyz))
69         {
70           AliL3Transform::Local2Global(xyz,0);
71           x[n] = xyz[0];
72           y[n] = xyz[1];
73           z[n] = xyz[2];
74           n++;
75         }
76       else
77         break;
78     }
79 }
80
81 TPolyMarker3D *AliL3HoughDisplay::LoadDigits()
82 {
83   
84   AliL3DigitRowData *tempPt = fDigitRowData;
85   if(!tempPt)
86     {
87       cerr<<"AliL3HoughDisplay::LoadDigits : No data"<<endl;
88       return 0;
89     }
90   
91   UInt_t nrows = AliL3Transform::GetNRows(fPatch);
92   Int_t count=0;
93   for(UInt_t i=0; i<nrows; i++)
94     {
95       count += tempPt->fNDigit;
96       AliL3MemHandler::UpdateRowPointer(tempPt);
97     }
98   tempPt = fDigitRowData;
99   TPolyMarker3D *pm = new TPolyMarker3D(count);
100   Float_t xyz[3];
101   Int_t sector,row;
102   count=0;
103   for(UInt_t i=0; i<nrows; i++)
104     {
105       AliL3DigitData *digPt = tempPt->fDigitData;
106       Int_t padrow = (Int_t)tempPt->fRow;
107       for(UInt_t j=0; j<tempPt->fNDigit; j++)
108         {
109           AliL3Transform::Slice2Sector(fShowSlice,padrow,sector,row);
110           AliL3Transform::Raw2Global(xyz,sector,row,(Int_t)digPt->fPad,(Int_t)digPt->fTime);
111           pm->SetPoint(count,xyz[0],xyz[1],xyz[2]);
112           count++;
113         }
114       AliL3MemHandler::UpdateRowPointer(tempPt);
115     }
116
117   cout<<"Displaying "<<count<<" digits"<<endl;
118   return pm;
119 }
120
121 void AliL3HoughDisplay::DisplayEvent()
122 {
123   //Display the found tracks.
124   
125   if(!fTracks)
126     {
127       cerr<<"AliL3HoughDisplay::DisplayTracks() : No tracks"<<endl;
128       return;
129     }
130   
131   TCanvas *c1 = new TCanvas("c1","",700,700);
132   c1->cd();
133   
134   TView *v = new TView(1);
135   v->SetRange(-430,-560,-430,430,560,1710);
136
137   c1->Clear();
138   c1->SetFillColor(1);
139   c1->SetTheta(90.);
140   c1->SetPhi(0.);
141     
142   Int_t ntracks = fTracks->GetNTracks();
143   TPolyLine3D *line = new TPolyLine3D[ntracks];
144   
145   Int_t n;
146   Float_t x[176],y[176],z[176];
147   for(Int_t j=0; j<ntracks; j++)
148     {
149       AliL3Track *track = fTracks->GetCheckedTrack(j); 
150       if(!track) continue;        
151       track->CalculateHelix();
152       GenerateHits(track,x,y,z,n);
153       TPolyMarker3D *pm = new TPolyMarker3D(n);
154       
155       for(Int_t h=0; h<n; h++)
156         pm->SetPoint(h,x[h],y[h],z[h]);
157       
158       pm->SetMarkerColor(2);
159       pm->Draw();
160       TPolyLine3D *current_line = &(line[j]);
161       current_line = new TPolyLine3D(n,x,y,z,"");
162       current_line->SetLineColor(4);
163       current_line->Draw("same");
164       
165     }
166   
167   if(fShowSlice>=0)
168     {
169       TPolyMarker3D *pm = LoadDigits();
170       pm->SetMarkerColor(2);
171       pm->Draw("same");
172     }
173   
174   cout<<"Displaying...."<<endl;
175   fGeom->Draw("same");
176   c1->x3d();
177 }
178   
179