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