]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/hough/AliL3HoughTest.cxx
Several bugfixes
[u/mrichter/AliRoot.git] / HLT / hough / AliL3HoughTest.cxx
1 //$Id$
2
3 // Author: Anders Vestbo <mailto:vestbo@fi.uib.no>
4 //*-- Copyright &copy ASV 
5
6 #include "AliL3StandardIncludes.h"
7 #include "AliL3HoughTest.h"
8 #include "AliL3ModelTrack.h"
9 #include "AliL3Transform.h"
10 #include "AliL3Histogram.h"
11 #include "TRandom.h"
12 #include "TMath.h"
13 #include "TH2.h"
14
15 #if GCCVERSION == 3
16 using namespace std;
17 #endif
18
19 //_____________________________________________________________
20 // AliL3HoughTest
21
22
23 ClassImp(AliL3HoughTest)
24
25 AliL3HoughTest::AliL3HoughTest()
26 {
27   fData=0;
28 }
29
30
31 AliL3HoughTest::~AliL3HoughTest()
32 {
33   if(fData)
34     delete [] fData;
35 }
36
37 Bool_t AliL3HoughTest::GenerateTrackData(Double_t pt,Double_t psi,Int_t sign,Int_t patch,Int_t minhits)
38 {
39   fCurrentPatch=patch;
40   if(fData)
41     delete fData;
42   fData = new SimData[AliL3Transform::GetNRows(patch)];
43   memset(fData,0,AliL3Transform::GetNRows(patch)*sizeof(SimData));
44   
45   AliL3ModelTrack *track = new AliL3ModelTrack();
46   track->Init(0,patch);
47   track->SetPt(pt);
48   track->SetPsi(psi);
49   track->SetCharge(sign);
50   track->SetFirstPoint(0,0,0);
51   track->CalculateHelix();
52   
53
54   Int_t temp[200];
55   Int_t entries=100;
56   Int_t clustercharge=100;
57   Int_t hitcounter=0;
58   for(Int_t i=AliL3Transform::GetFirstRow(patch); i<=AliL3Transform::GetLastRow(patch); i++)
59     {
60       Float_t xyz[3];
61       
62       if(!track->GetCrossingPoint(i,xyz))
63         continue;
64       
65       Int_t rowindex = i - AliL3Transform::GetFirstRow(patch);
66       Int_t sector,row;
67       AliL3Transform::Slice2Sector(0,i,sector,row);
68       AliL3Transform::Local2Raw(xyz,sector,row);
69
70       if(xyz[1] < 0 || xyz[1] >= AliL3Transform::GetNPads(i))
71         continue;
72       hitcounter++;
73       track->SetPadHit(i,xyz[1]);
74       //cout<<i<<" "<<xyz[1]<<endl;
75       memset(temp,0,200*sizeof(Int_t));
76       Double_t xysigma = sqrt(track->GetParSigmaY2(i));
77       Int_t minpad=200,j;
78       for(j=0; j<entries; j++)
79         {
80           Int_t pad = TMath::Nint(gRandom->Gaus(xyz[1],xysigma));
81           if(pad < 0 || pad >= AliL3Transform::GetNPads(i))
82             continue;
83           temp[pad]++;
84           if(pad < minpad)
85             minpad=pad;
86         }
87       for(j=0; j<200; j++)
88         {
89           if(temp[j]==0) continue;
90           Int_t index = j - minpad;
91           if(index < 0 || index >= 10)
92             {
93               cerr<<"AliL3HoughTest::GenerateTrackData : Wrong index "<<index<<endl;
94               exit(5);
95             }
96           Int_t charge = clustercharge*temp[j]/entries;
97           if(charge < 3 ) 
98             continue;
99           if(charge > 1023)
100             charge=1023;
101           fData[rowindex].pads[index]=charge;
102           fData[rowindex].npads++;
103         }
104       fData[rowindex].minpad=minpad;
105     }
106   delete track;
107   if(hitcounter < minhits)
108     return kFALSE;
109   return kTRUE;
110 }
111
112 void AliL3HoughTest::Transform2Circle(AliL3Histogram *hist)
113 {
114   if(!fData)
115     {
116       cerr<<"AliL3HoughTest::Transform : No data"<<endl;
117       return;
118     }
119   Float_t R,phi,phi0,kappa,xyz[3];
120   Int_t pad,charge,sector,row;
121   for(Int_t i=AliL3Transform::GetFirstRow(fCurrentPatch); i<=AliL3Transform::GetLastRow(fCurrentPatch); i++)
122     {
123       Int_t rowindex = i - AliL3Transform::GetFirstRow(fCurrentPatch);
124       AliL3Transform::Slice2Sector(0,i,sector,row);
125       for(Int_t j=0; j<fData[rowindex].npads; j++)
126         {
127           pad = j + fData[rowindex].minpad;
128           charge = fData[rowindex].pads[j];
129           
130           AliL3Transform::Raw2Local(xyz,sector,row,pad,0);
131
132           R = sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]);
133
134           phi = AliL3Transform::GetPhi(xyz);
135           
136           for(Int_t k=hist->GetFirstYbin(); k<=hist->GetLastYbin(); k++)
137             {
138               phi0 = hist->GetBinCenterY(k);
139               kappa = 2*sin(phi-phi0)/R;
140               hist->Fill(kappa,phi0,charge);
141             }
142         }
143     }
144 }
145
146 void AliL3HoughTest::Transform2CircleC(AliL3Histogram *hist)
147 {
148   if(!fData)
149     {
150       cerr<<"AliL3HoughTest::TransformC : No data"<<endl;
151       return;
152     }
153   Int_t pad1,pad2,charge1,charge2,sector,row;
154   Float_t r1,r2,phi1,phi2,phi_0,kappa,hit[3],hit2[3];
155   for(Int_t i=AliL3Transform::GetFirstRow(fCurrentPatch); i<=AliL3Transform::GetLastRow(fCurrentPatch); i++)
156     {
157       Int_t rowindex1 = i - AliL3Transform::GetFirstRow(fCurrentPatch);
158       for(Int_t d1=0; d1<fData[rowindex1].npads; d1++)
159         {
160           pad1 = d1 + fData[rowindex1].minpad;
161           charge1 = fData[rowindex1].pads[d1];
162           AliL3Transform::Slice2Sector(0,i,sector,row);
163           AliL3Transform::Raw2Local(hit,sector,row,pad1,0);
164           r1 = sqrt(hit[0]*hit[0]+hit[1]*hit[1]);
165           phi1 = atan2(hit[1],hit[0]);
166           
167           for(Int_t j=i+1; j<=AliL3Transform::GetLastRow(fCurrentPatch); j++)
168             {
169               Int_t rowindex2 = j - AliL3Transform::GetFirstRow(fCurrentPatch);
170               for(Int_t d2=0; d2<fData[rowindex2].npads; d2++)
171                 {
172                   pad2 = d2 + fData[rowindex2].minpad;
173                   charge2 = fData[rowindex2].pads[d2];
174                   AliL3Transform::Slice2Sector(0,j,sector,row);
175                   AliL3Transform::Raw2Local(hit2,sector,row,pad2,0);
176                   r2 = sqrt(hit2[0]*hit2[0]+hit2[1]*hit2[1]);
177                   phi2 = atan2(hit2[1],hit2[0]);
178                   phi_0 = atan( (r2*sin(phi1) - r1*sin(phi2)) / (r2*cos(phi1) - r1*cos(phi2)) );
179
180                   kappa = 2*sin(phi1-phi_0) / r1;
181                   hist->Fill(kappa,phi_0,charge1+charge2);
182                 }
183             }
184         }
185     }
186 }
187
188 void AliL3HoughTest::Transform2Line(AliL3Histogram *hist,Int_t *rowrange)
189 {
190   if(!fData)
191     {
192       cerr<<"AliL3HoughTest::Transform2Line : No data"<<endl;
193       return;
194     }
195   
196   Int_t pad,charge,sector,row;
197   Float_t hit[3],theta,rho;
198   for(Int_t i=rowrange[0]; i<=rowrange[1]; i++)
199     {
200       Int_t rowindex = i - AliL3Transform::GetFirstRow(fCurrentPatch);
201       for(Int_t d=0; d<fData[rowindex].npads; d++)
202         {
203           pad = d + fData[rowindex].minpad;
204           charge = fData[rowindex].pads[d];
205           AliL3Transform::Slice2Sector(0,i,sector,row);
206           AliL3Transform::Raw2Local(hit,sector,row,pad,0);
207           
208           hit[0] = hit[0] - AliL3Transform::Row2X(rowrange[0]);
209           
210           for(Int_t xbin=hist->GetFirstXbin(); xbin<hist->GetLastXbin(); xbin++)
211             {
212               theta = hist->GetBinCenterX(xbin);
213               rho = hit[0]*cos(theta) + hit[1]*sin(theta);
214               hist->Fill(theta,rho,charge);
215             }
216         }
217     }
218 }
219
220 void AliL3HoughTest::Transform2LineC(AliL3Histogram *hist,Int_t *rowrange)
221 {
222   if(!fData)
223     {
224       cerr<<"AliL3HoughTest::Transform2Line : No data"<<endl;
225       return;
226     }
227   
228   Int_t pad1,pad2,charge1,charge2,sector,row;
229   Float_t theta,rho,hit[3],hit2[3];
230   for(Int_t i=rowrange[0]; i<=rowrange[1]; i++)
231     {
232       Int_t rowindex1 = i - AliL3Transform::GetFirstRow(fCurrentPatch);
233       for(Int_t d1=0; d1<fData[rowindex1].npads; d1++)
234         {
235           pad1 = d1 + fData[rowindex1].minpad;
236           charge1 = fData[rowindex1].pads[d1];
237           AliL3Transform::Slice2Sector(0,i,sector,row);
238           AliL3Transform::Raw2Local(hit,sector,row,pad1,0);
239           
240           hit[0] = hit[0] - AliL3Transform::Row2X(rowrange[0]);
241           
242           for(Int_t j=i+1; j<=rowrange[1]; j++)
243             {
244               Int_t rowindex2 = j - AliL3Transform::GetFirstRow(fCurrentPatch);
245               for(Int_t d2=0; d2<fData[rowindex2].npads; d2++)
246                 {
247                   pad2 = d2 + fData[rowindex2].minpad;
248                   charge2 = fData[rowindex2].pads[d2];
249                   AliL3Transform::Slice2Sector(0,j,sector,row);
250                   AliL3Transform::Raw2Local(hit2,sector,row,pad2,0);
251                   
252                   hit2[0] = hit2[0] - AliL3Transform::Row2X(rowrange[0]);
253                   
254                   theta = atan2(hit2[0]-hit[0],hit[1]-hit2[1]);
255                   rho = hit[0]*cos(theta)+hit[1]*sin(theta);
256                   hist->Fill(theta,rho,charge1+charge2);
257                 }
258             }
259         }
260     }
261 }
262
263
264 void AliL3HoughTest::FillImage(TH2 *hist)
265 {
266   if(!fData)
267     {
268       cerr<<"AliL3HoughTest::FillImage : No data to fill"<<endl;
269       return;
270     }
271   for(Int_t i=AliL3Transform::GetFirstRow(fCurrentPatch); i<=AliL3Transform::GetLastRow(fCurrentPatch); i++)
272     {
273       Int_t rowindex = i - AliL3Transform::GetFirstRow(fCurrentPatch);
274       for(Int_t j=0; j<fData[rowindex].npads; j++)
275         {
276           Int_t pad = j + fData[rowindex].minpad;
277           Int_t charge = fData[rowindex].pads[j];
278           Float_t xyz[3];
279           Int_t sector,row;
280           AliL3Transform::Slice2Sector(0,i,sector,row);
281           AliL3Transform::Raw2Local(xyz,sector,row,pad,0);
282           hist->Fill(xyz[0],xyz[1],charge);
283         }
284     }
285 }