]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/hough/AliHLTHoughTest.cxx
compilation warnings corrected
[u/mrichter/AliRoot.git] / HLT / hough / AliHLTHoughTest.cxx
CommitLineData
3e87ef69 1// @(#) $Id$
7a1f1432 2
3// Author: Anders Vestbo <mailto:vestbo@fi.uib.no>
3e87ef69 4//*-- Copyright &copy ALICE HLT Group
7a1f1432 5
4aa41877 6#include "AliHLTStandardIncludes.h"
7#include "AliHLTHoughTest.h"
8#include "AliHLTModelTrack.h"
9#include "AliHLTTransform.h"
10#include "AliHLTHistogram.h"
11#include "AliHLTTrackArray.h"
12#include "AliHLTHoughTrack.h"
3e87ef69 13
14#include <TRandom.h>
15#include <TMath.h>
16#include <TH2.h>
17#include <TH3.h>
18#include <TAxis.h>
7a1f1432 19
5929c18d 20#if __GNUC__ >= 3
7a1f1432 21using namespace std;
22#endif
23
24//_____________________________________________________________
4aa41877 25// AliHLTHoughTest
7a1f1432 26
27
4aa41877 28ClassImp(AliHLTHoughTest)
7a1f1432 29
4aa41877 30AliHLTHoughTest::AliHLTHoughTest()
7a1f1432 31{
c62b480b 32 //ctor
7a1f1432 33 fData=0;
34}
35
36
4aa41877 37AliHLTHoughTest::~AliHLTHoughTest()
7a1f1432 38{
c62b480b 39 //dtor
7a1f1432 40 if(fData)
41 delete [] fData;
42}
43
4aa41877 44Bool_t AliHLTHoughTest::GenerateTrackData(Double_t pt,Double_t psi,Double_t tgl,Int_t sign,Int_t patch,Int_t minhits)
7a1f1432 45{
c62b480b 46 //Generate digits according to given track parameters?
7a1f1432 47 fCurrentPatch=patch;
48 if(fData)
49 delete fData;
4aa41877 50 fData = new AliHLTSimData[AliHLTTransform::GetNRows(patch)];
51 memset(fData,0,AliHLTTransform::GetNRows(patch)*sizeof(AliHLTSimData));
7a1f1432 52
4aa41877 53 AliHLTModelTrack *track = new AliHLTModelTrack();
7a1f1432 54 track->Init(0,patch);
55 track->SetPt(pt);
56 track->SetPsi(psi);
c3853fab 57 track->SetTgl(tgl);
d7ab65da 58 track->SetCharge(sign);
59 track->SetFirstPoint(0,0,0);
7a1f1432 60 track->CalculateHelix();
7a1f1432 61
62 Int_t temp[200];
4aa41877 63 // Int_t temp2[AliHLTTransform::GetNTimeBins()];
64 Int_t * temp2 = new Int_t[AliHLTTransform::GetNTimeBins()];
7a1f1432 65 Int_t entries=100;
66 Int_t clustercharge=100;
d7ab65da 67 Int_t hitcounter=0;
4aa41877 68 for(Int_t i=AliHLTTransform::GetFirstRow(patch); i<=AliHLTTransform::GetLastRow(patch); i++)
7a1f1432 69 {
7a1f1432 70
3e87ef69 71 Float_t xyz[3];
7a1f1432 72 if(!track->GetCrossingPoint(i,xyz))
73 continue;
3e87ef69 74
4aa41877 75 Int_t rowindex = i - AliHLTTransform::GetFirstRow(patch);
7a1f1432 76 Int_t sector,row;
4aa41877 77 AliHLTTransform::Slice2Sector(0,i,sector,row);
78 AliHLTTransform::Local2Raw(xyz,sector,row);
7a1f1432 79
4aa41877 80 if(xyz[1] < 0 || xyz[1] >= AliHLTTransform::GetNPads(i) || xyz[2] < 0 || xyz[2] >= AliHLTTransform::GetNTimeBins())
7a1f1432 81 continue;
d7ab65da 82 hitcounter++;
83 track->SetPadHit(i,xyz[1]);
c3853fab 84 track->SetTimeHit(i,xyz[2]);
3e87ef69 85
86 Double_t beta = track->GetCrossingAngle(i);
87 track->SetCrossingAngleLUT(i,beta);
88 track->CalculateClusterWidths(i);
c3853fab 89
7a1f1432 90 memset(temp,0,200*sizeof(Int_t));
4aa41877 91 memset(temp2,0,AliHLTTransform::GetNTimeBins()*sizeof(Int_t));
7a1f1432 92 Double_t xysigma = sqrt(track->GetParSigmaY2(i));
c3853fab 93 Double_t zsigma = sqrt(track->GetParSigmaZ2(i));
3e87ef69 94
7a1f1432 95 Int_t minpad=200,j;
c3853fab 96 Int_t mintime = 1000;
7a1f1432 97 for(j=0; j<entries; j++)
98 {
99 Int_t pad = TMath::Nint(gRandom->Gaus(xyz[1],xysigma));
c3853fab 100 Int_t time = TMath::Nint(gRandom->Gaus(xyz[2],zsigma));
4aa41877 101 if(pad < 0 || pad >= AliHLTTransform::GetNPads(i) || time < 0 || time >= AliHLTTransform::GetNTimeBins())
7a1f1432 102 continue;
103 temp[pad]++;
c3853fab 104 temp2[time]++;
7a1f1432 105 if(pad < minpad)
106 minpad=pad;
c3853fab 107 if(time < mintime)
108 mintime=time;
7a1f1432 109 }
c3853fab 110 Int_t npads=0;
7a1f1432 111 for(j=0; j<200; j++)
112 {
113 if(temp[j]==0) continue;
c3853fab 114
7a1f1432 115 Int_t index = j - minpad;
c3853fab 116
7a1f1432 117 if(index < 0 || index >= 10)
118 {
4aa41877 119 cerr<<"AliHLTHoughTest::GenerateTrackData : Wrong index "<<index<<endl;
7a1f1432 120 exit(5);
121 }
c3853fab 122 npads++;
c62b480b 123 Int_t seqcharge = clustercharge*temp[j]/entries;
c3853fab 124 Int_t ntimes=0;
4aa41877 125 for(Int_t k=0; k<AliHLTTransform::GetNTimeBins(); k++)
c3853fab 126 {
127 if(temp2[k]==0) continue;
128 Int_t tindex = k - mintime;
129 if(tindex < 0 || tindex >= 10)
130 {
4aa41877 131 cerr<<"AliHLTHoughTest::GenerateTrackData : Wrong timeindex "<<tindex<<" "<<k<<" "<<mintime<<endl;
c3853fab 132 exit(5);
133 }
c62b480b 134 Int_t charge = seqcharge*temp2[k]/entries;
c3853fab 135 if(charge < 3 )
136 continue;
137 if(charge > 1023)
138 charge=1023;
139 //cout<<"row "<<i<<" pad "<<j<<" time "<<k<<" charge "<<charge<<endl;
140 ntimes++;
c62b480b 141 fData[rowindex].fPads[index][tindex]=charge;
c3853fab 142 }
7a1f1432 143 }
c62b480b 144 fData[rowindex].fMinpad=minpad;
145 fData[rowindex].fMintime=mintime;
146 fData[rowindex].fNPads=npads;
7a1f1432 147 }
148 delete track;
62bb4b3d 149 delete [] temp2;
d7ab65da 150 if(hitcounter < minhits)
151 return kFALSE;
152 return kTRUE;
7a1f1432 153}
154
4aa41877 155void AliHLTHoughTest::Transform2Circle(AliHLTHistogram *hist)
7a1f1432 156{
c62b480b 157 //Hough trasnform
7a1f1432 158 if(!fData)
159 {
4aa41877 160 cerr<<"AliHLTHoughTest::Transform : No data"<<endl;
7a1f1432 161 return;
162 }
c62b480b 163 Float_t r,phi,phi0,kappa,xyz[3];
c3853fab 164 Int_t pad,time,charge,sector,row;
4aa41877 165 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
7a1f1432 166 {
4aa41877 167 Int_t rowindex = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
168 AliHLTTransform::Slice2Sector(0,i,sector,row);
c62b480b 169 for(Int_t j=0; j<fData[rowindex].fNPads; j++)
7a1f1432 170 {
c62b480b 171 pad = j + fData[rowindex].fMinpad;
c3853fab 172 for(Int_t k=0; k<10; k++)
7a1f1432 173 {
c62b480b 174 time = k + fData[rowindex].fMintime;
175 charge = fData[rowindex].fPads[j][k];
c3853fab 176 if(charge == 0) continue;
4aa41877 177 AliHLTTransform::Raw2Local(xyz,sector,row,pad,time);
c3853fab 178
c62b480b 179 r = sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]);
c3853fab 180
4aa41877 181 phi = AliHLTTransform::GetPhi(xyz);
c3853fab 182
183 for(Int_t k=hist->GetFirstYbin(); k<=hist->GetLastYbin(); k++)
184 {
185 phi0 = hist->GetBinCenterY(k);
c62b480b 186 kappa = 2*sin(phi-phi0)/r;
c3853fab 187 hist->Fill(kappa,phi0,charge);
188 }
7a1f1432 189 }
190 }
191 }
192}
193
4aa41877 194void AliHLTHoughTest::Transform2CircleC(AliHLTHistogram *hist)
7a1f1432 195{
c62b480b 196 //Hough transform
7a1f1432 197 if(!fData)
198 {
4aa41877 199 cerr<<"AliHLTHoughTest::TransformC : No data"<<endl;
7a1f1432 200 return;
201 }
c3853fab 202 Int_t pad1,pad2,time1,time2,charge1,charge2,sector,row;
c62b480b 203 Float_t r1,r2,phi1,phi2,phi0,kappa,hit[3],hit2[3];
4aa41877 204 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
7a1f1432 205 {
4aa41877 206 Int_t rowindex1 = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 207 for(Int_t d1=0; d1<fData[rowindex1].fNPads; d1++)
7a1f1432 208 {
c62b480b 209 pad1 = d1 + fData[rowindex1].fMinpad;
c3853fab 210 for(Int_t j=0; j<10; j++)
7a1f1432 211 {
c62b480b 212 time1 = j + fData[rowindex1].fMintime;
213 charge1 = fData[rowindex1].fPads[d1][j];
c3853fab 214 if(charge1==0) continue;
4aa41877 215 AliHLTTransform::Slice2Sector(0,i,sector,row);
216 AliHLTTransform::Raw2Local(hit,sector,row,pad1,time1);
c3853fab 217 r1 = sqrt(hit[0]*hit[0]+hit[1]*hit[1]);
218 phi1 = atan2(hit[1],hit[0]);
219
4aa41877 220 for(Int_t j=i+1; j<=AliHLTTransform::GetLastRow(fCurrentPatch); j++)
7a1f1432 221 {
4aa41877 222 Int_t rowindex2 = j - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 223 for(Int_t d2=0; d2<fData[rowindex2].fNPads; d2++)
c3853fab 224 {
c62b480b 225 pad2 = d2 + fData[rowindex2].fMinpad;
c3853fab 226 for(Int_t k=0; k<10; k++)
227 {
c62b480b 228 time2 = k + fData[rowindex2].fMintime;
229 charge2 = fData[rowindex2].fPads[d2][k];
c3853fab 230 if(charge2==0) continue;
4aa41877 231 AliHLTTransform::Slice2Sector(0,j,sector,row);
232 AliHLTTransform::Raw2Local(hit2,sector,row,pad2,time2);
c3853fab 233 r2 = sqrt(hit2[0]*hit2[0]+hit2[1]*hit2[1]);
234 phi2 = atan2(hit2[1],hit2[0]);
c62b480b 235 phi0 = atan( (r2*sin(phi1) - r1*sin(phi2)) / (r2*cos(phi1) - r1*cos(phi2)) );
c3853fab 236
c62b480b 237 kappa = 2*sin(phi1-phi0) / r1;
238 hist->Fill(kappa,phi0,charge1+charge2);
c3853fab 239 }
240 }
7a1f1432 241 }
3e87ef69 242 return;
7a1f1432 243 }
244 }
245 }
246}
247
4aa41877 248void AliHLTHoughTest::Transform2CircleF(AliHLTHistogram *hist)
3e87ef69 249{
250 //Fix one point in the middle of the tpc
251
252 if(!fData)
253 {
4aa41877 254 cerr<<"AliHLTHoughTest::TransformF : No data"<<endl;
3e87ef69 255 return;
256 }
257 Int_t pad1,pad2,time1,time2,charge1,charge2,sector,row;
c62b480b 258 Float_t r1,r2,phi1,phi2,phi0,kappa,hit[3],hit2[3];
3e87ef69 259
4aa41877 260 Int_t rowindex1 = 80 - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 261 for(Int_t d1=1; d1<fData[rowindex1].fNPads; d1++)
3e87ef69 262 {
c62b480b 263 pad1 = d1 + fData[rowindex1].fMinpad;
3e87ef69 264
4aa41877 265 AliHLTTransform::Slice2Sector(0,80,sector,row);
266 AliHLTTransform::Raw2Local(hit,sector,row,pad1,0);
3e87ef69 267 r1 = sqrt(hit[0]*hit[0]+hit[1]*hit[1]);
268 phi1 = atan2(hit[1],hit[0]);
269
270 for(Int_t j=0; j<10; j++)
271 {
c62b480b 272 time1 = j + fData[rowindex1].fMintime;
273 charge1 = fData[rowindex1].fPads[d1][j];
3e87ef69 274 if(charge1==0) continue;
275
4aa41877 276 for(Int_t j=0; j<=AliHLTTransform::GetLastRow(fCurrentPatch); j++)
3e87ef69 277 {
278 if(j==80) continue;
4aa41877 279 Int_t rowindex2 = j - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 280 for(Int_t d2=0; d2<fData[rowindex2].fNPads; d2++)
3e87ef69 281 {
c62b480b 282 pad2 = d2 + fData[rowindex2].fMinpad;
3e87ef69 283 for(Int_t k=0; k<10; k++)
284 {
c62b480b 285 time2 = k + fData[rowindex2].fMintime;
286 charge2 = fData[rowindex2].fPads[d2][k];
3e87ef69 287 if(charge2==0) continue;
4aa41877 288 AliHLTTransform::Slice2Sector(0,j,sector,row);
289 AliHLTTransform::Raw2Local(hit2,sector,row,pad2,time2);
3e87ef69 290 r2 = sqrt(hit2[0]*hit2[0]+hit2[1]*hit2[1]);
291 phi2 = atan2(hit2[1],hit2[0]);
c62b480b 292 phi0 = atan( (r2*sin(phi1) - r1*sin(phi2)) / (r2*cos(phi1) - r1*cos(phi2)) );
3e87ef69 293
c62b480b 294 kappa = 2*sin(phi1-phi0) / r1;
295 hist->Fill(kappa,phi0,charge1+charge2);
296 //cout<<"Filling "<<kappa<<" psi "<<phi0<<" charge "<<charge1<<endl;
3e87ef69 297 }
298 }
299 }
300 }
3e87ef69 301 }
e33f3609 302 return;
3e87ef69 303}
304
4aa41877 305void AliHLTHoughTest::Transform2Line(AliHLTHistogram *hist,Int_t *rowrange)
d7ab65da 306{
c62b480b 307 //Hough transform
d7ab65da 308 if(!fData)
309 {
4aa41877 310 cerr<<"AliHLTHoughTest::Transform2Line : No data"<<endl;
d7ab65da 311 return;
312 }
313
c3853fab 314 Int_t pad,time,charge,sector,row;
d7ab65da 315 Float_t hit[3],theta,rho;
4aa41877 316 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
d7ab65da 317 {
3e87ef69 318 if(i<rowrange[0])
319 continue;
320 if(i>rowrange[1])
321 break;
4aa41877 322 Int_t rowindex = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 323 for(Int_t d=0; d<fData[rowindex].fNPads; d++)
d7ab65da 324 {
c62b480b 325 pad = d + fData[rowindex].fMinpad;
c3853fab 326 for(Int_t j=0; j<10; j++)
d7ab65da 327 {
c62b480b 328 time = j + fData[rowindex].fMintime;
329 charge = fData[rowindex].fPads[d][j];
c3853fab 330 if(charge==0) continue;
4aa41877 331 AliHLTTransform::Slice2Sector(0,i,sector,row);
332 AliHLTTransform::Raw2Local(hit,sector,row,pad,time);
c3853fab 333
4aa41877 334 hit[0] = hit[0] - AliHLTTransform::Row2X(rowrange[0]);
c3853fab 335
336 for(Int_t xbin=hist->GetFirstXbin(); xbin<hist->GetLastXbin(); xbin++)
337 {
338 theta = hist->GetBinCenterX(xbin);
339 rho = hit[0]*cos(theta) + hit[1]*sin(theta);
340 hist->Fill(theta,rho,charge);
341 }
d7ab65da 342 }
343 }
344 }
345}
346
4aa41877 347void AliHLTHoughTest::Transform2LineC(AliHLTHistogram *hist,Int_t *rowrange)
d7ab65da 348{
c62b480b 349 //Hough transform?
d7ab65da 350 if(!fData)
351 {
4aa41877 352 cerr<<"AliHLTHoughTest::Transform2Line : No data"<<endl;
d7ab65da 353 return;
354 }
355
c3853fab 356 Int_t pad1,pad2,time1,time2,charge1,charge2,sector,row;
d7ab65da 357 Float_t theta,rho,hit[3],hit2[3];
358 for(Int_t i=rowrange[0]; i<=rowrange[1]; i++)
359 {
4aa41877 360 Int_t rowindex1 = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 361 for(Int_t d1=0; d1<fData[rowindex1].fNPads; d1++)
d7ab65da 362 {
c62b480b 363 pad1 = d1 + fData[rowindex1].fMinpad;
c3853fab 364 for(Int_t j=0; j<10; j++)
d7ab65da 365 {
c62b480b 366 time1 = j + fData[rowindex1].fMintime;
367 charge1 = fData[rowindex1].fPads[d1][j];
c3853fab 368 if(charge1==0) continue;
4aa41877 369 AliHLTTransform::Slice2Sector(0,i,sector,row);
370 AliHLTTransform::Raw2Local(hit,sector,row,pad1,time1);
c3853fab 371
4aa41877 372 hit[0] = hit[0] - AliHLTTransform::Row2X(rowrange[0]);
c3853fab 373
3e87ef69 374 for(Int_t i2=i+1; i2<=rowrange[1]; i2++)
d7ab65da 375 {
4aa41877 376 Int_t rowindex2 = i2 - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 377 for(Int_t d2=0; d2<fData[rowindex2].fNPads; d2++)
c3853fab 378 {
c62b480b 379 pad2 = d2 + fData[rowindex2].fMinpad;
c3853fab 380 for(Int_t k=0; k<10; k++)
381 {
c62b480b 382 time2 = k + fData[rowindex2].fMintime;
383 charge2 = fData[rowindex2].fPads[d2][k];
c3853fab 384 if(charge2==0) continue;
4aa41877 385 AliHLTTransform::Slice2Sector(0,i2,sector,row);
386 AliHLTTransform::Raw2Local(hit2,sector,row,pad2,time2);
c3853fab 387
4aa41877 388 hit2[0] = hit2[0] - AliHLTTransform::Row2X(rowrange[0]);
c3853fab 389
390 theta = atan2(hit2[0]-hit[0],hit[1]-hit2[1]);
391 rho = hit[0]*cos(theta)+hit[1]*sin(theta);
392 hist->Fill(theta,rho,charge1+charge2);
393 }
394 }
d7ab65da 395 }
396 }
397 }
398 }
399}
400
4aa41877 401void AliHLTHoughTest::FillImage(TH2 *hist,Int_t displayrow)
7a1f1432 402{
c62b480b 403 //Draw digits data
7a1f1432 404 if(!fData)
405 {
4aa41877 406 cerr<<"AliHLTHoughTest::FillImage : No data to fill"<<endl;
7a1f1432 407 return;
408 }
c3853fab 409
4aa41877 410 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
7a1f1432 411 {
4aa41877 412 Int_t rowindex = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
3e87ef69 413 if(displayrow >=0)
414 if(i != displayrow) continue;
c3853fab 415
c62b480b 416 //cout<<"row "<<i<<" fNPads "<<fData[rowindex].fNPads<<endl;
417 for(Int_t j=0; j<fData[rowindex].fNPads; j++)
7a1f1432 418 {
c62b480b 419 Int_t pad = j + fData[rowindex].fMinpad;
c3853fab 420 for(Int_t k=0; k<10; k++)
421 {
c62b480b 422 Int_t time = k + fData[rowindex].fMintime;
423 Int_t charge = fData[rowindex].fPads[j][k];
c3853fab 424 if(charge==0) continue;
425 //cout<<i<<" "<<pad<<" "<<time<<" "<<charge<<endl;
426 Float_t xyz[3];
427 Int_t sector,row;
4aa41877 428 AliHLTTransform::Slice2Sector(0,i,sector,row);
429 AliHLTTransform::Raw2Local(xyz,sector,row,pad,time);
3e87ef69 430 if(displayrow>=0)
c3853fab 431 hist->Fill(pad,time,charge);
432 else
433 hist->Fill(xyz[0],xyz[1],charge);
434 }
7a1f1432 435 }
3e87ef69 436 if(displayrow>=0)
c3853fab 437 break;
7a1f1432 438 }
439}
c3853fab 440
4aa41877 441void AliHLTHoughTest::Transform2Line3D(TH3 *hist,Int_t *rowrange,Float_t *phirange)
c3853fab 442{
c62b480b 443 //HT in 3D?
c3853fab 444 if(!fData)
445 {
4aa41877 446 cerr<<"AliHLTHoughTest::Transform2Line : No data"<<endl;
c3853fab 447 return;
448 }
449
450 Int_t pad,time,charge,sector,row;
c62b480b 451 Float_t hit[3],theta,rho,r,delta;
4aa41877 452 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
c3853fab 453 {
3e87ef69 454
455 if(i<rowrange[0])
456 continue;
457 if(i>rowrange[1])
458 break;
459
4aa41877 460 Int_t rowindex = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 461 for(Int_t d=0; d<fData[rowindex].fNPads; d++)
c3853fab 462 {
c62b480b 463 pad = d + fData[rowindex].fMinpad;
c3853fab 464 for(Int_t j=0; j<10; j++)
465 {
c62b480b 466 time = j + fData[rowindex].fMintime;
467 charge = fData[rowindex].fPads[d][j];
c3853fab 468 if(charge==0) continue;
4aa41877 469 AliHLTTransform::Slice2Sector(0,i,sector,row);
470 AliHLTTransform::Raw2Local(hit,sector,row,pad,time);
c3853fab 471
4aa41877 472 Float_t phi = AliHLTTransform::GetPhi(hit);
3e87ef69 473 if(phi < phirange[0] || phi > phirange[1])
474 continue;
475
4aa41877 476 hit[0] = hit[0] - AliHLTTransform::Row2X(rowrange[0]);
477 Float_t x = hit[0] + AliHLTTransform::Row2X(rowrange[0]);
c62b480b 478 r = sqrt(x*x + hit[1]*hit[1]);
479 delta = atan(hit[2]/r);
c3853fab 480
481 for(Int_t xbin=hist->GetXaxis()->GetFirst(); xbin<=hist->GetXaxis()->GetLast(); xbin++)
482 {
483 theta = hist->GetXaxis()->GetBinCenter(xbin);
484 rho = hit[0]*cos(theta) + hit[1]*sin(theta);
c3853fab 485 hist->Fill(theta,rho,delta,charge);
486 }
487 }
488 }
489 }
490}
491
4aa41877 492void AliHLTHoughTest::Transform2LineC3D(TH3 *hist,Int_t *rowrange)
c3853fab 493{
c62b480b 494 //HT in 3D?
c3853fab 495 if(!fData)
496 {
4aa41877 497 cerr<<"AliHLTHoughTest::Transform2Line : No data"<<endl;
c3853fab 498 return;
499 }
500
501 Int_t pad1,pad2,time1,time2,charge1,charge2,sector,row;
c62b480b 502 Float_t theta,rho,hit[3],hit2[3],r1,r2,delta,delta1,delta2;
4aa41877 503 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
c3853fab 504 {
3e87ef69 505 if(i<rowrange[0])
506 continue;
507 if(i>rowrange[1])
508 break;
4aa41877 509 Int_t rowindex1 = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 510 for(Int_t d1=0; d1<fData[rowindex1].fNPads; d1++)
c3853fab 511 {
c62b480b 512 pad1 = d1 + fData[rowindex1].fMinpad;
c3853fab 513 for(Int_t j=0; j<10; j++)
514 {
c62b480b 515 time1 = j + fData[rowindex1].fMintime;
516 charge1 = fData[rowindex1].fPads[d1][j];
c3853fab 517 if(charge1==0) continue;
4aa41877 518 AliHLTTransform::Slice2Sector(0,i,sector,row);
519 AliHLTTransform::Raw2Local(hit,sector,row,pad1,time1);
c62b480b 520 r1 = sqrt(hit[0]*hit[0]+hit[1]*hit[1]);
521 delta1 = atan(hit[2]/r1);
4aa41877 522 hit[0] = hit[0] - AliHLTTransform::Row2X(rowrange[0]);
c3853fab 523
3e87ef69 524 for(Int_t i2=i+1; i2<=rowrange[1]; i2++)
c3853fab 525 {
4aa41877 526 Int_t rowindex2 = i2 - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 527 for(Int_t d2=0; d2<fData[rowindex2].fNPads; d2++)
c3853fab 528 {
c62b480b 529 pad2 = d2 + fData[rowindex2].fMinpad;
c3853fab 530 for(Int_t k=0; k<10; k++)
531 {
c62b480b 532 time2 = k + fData[rowindex2].fMintime;
533 charge2 = fData[rowindex2].fPads[d2][k];
c3853fab 534 if(charge2==0) continue;
4aa41877 535 AliHLTTransform::Slice2Sector(0,i2,sector,row);
536 AliHLTTransform::Raw2Local(hit2,sector,row,pad2,time2);
c62b480b 537 r2 = sqrt(hit2[0]*hit2[0]+hit2[1]*hit2[1]);
538 delta2 = atan(hit2[2]/r2);
c3853fab 539 delta = (charge1*delta1 + charge2*delta2)/(charge1+charge2);
4aa41877 540 hit2[0] = hit2[0] - AliHLTTransform::Row2X(rowrange[0]);
c3853fab 541
542 theta = atan2(hit2[0]-hit[0],hit[1]-hit2[1]);
543 rho = hit[0]*cos(theta)+hit[1]*sin(theta);
544 hist->Fill(theta,rho,delta,charge1+charge2);
545 }
546 }
547 }
548 }
549 }
550 }
551}
3e87ef69 552
4aa41877 553void AliHLTHoughTest::TransformLines2Circle(TH3 *hist,AliHLTTrackArray *tracks)
3e87ef69 554{
c62b480b 555 //Another HT?
3e87ef69 556 for(Int_t i=0; i<tracks->GetNTracks(); i++)
557 {
4aa41877 558 AliHLTHoughTrack *tr = (AliHLTHoughTrack*)tracks->GetCheckedTrack(i);
3e87ef69 559 if(!tr) continue;
560 Int_t middlerow = (Int_t)(tr->GetFirstRow()+(tr->GetLastRow()-tr->GetFirstRow())/2);
561 Float_t hit[3];
562 tr->GetLineCrossingPoint(middlerow,hit);
4aa41877 563 hit[0] += AliHLTTransform::Row2X(tr->GetFirstRow());
c62b480b 564 Float_t r = sqrt(hit[0]*hit[0] + hit[1]*hit[1]);
565 hit[2] = r*tr->GetTgl();
3e87ef69 566 Float_t phi = atan2(hit[1],hit[0]);
4aa41877 567 Float_t theta = tr->GetPsiLine() - AliHLTTransform::Pi()/2;
3e87ef69 568 Float_t psi = 2*phi - theta;
c62b480b 569 Float_t kappa = 2/r*sin(phi-psi);
570 Float_t delta = atan(hit[2]/r);
3e87ef69 571 hist->Fill(kappa,psi,delta,tr->GetWeight());
572
573 }
574}
575
4aa41877 576void AliHLTHoughTest::Transform2Center(AliHLTHistogram *hist)
3e87ef69 577{
578 //Choose parameter space to be center of curvature.
579
580 if(!fData)
581 {
4aa41877 582 cerr<<"AliHLTHoughTest::TransformC : No data"<<endl;
3e87ef69 583 return;
584 }
585 Int_t pad1,pad2,time1,time2,charge1,charge2,sector,row;
586 Float_t r1,r2,phi1,phi2,hit[3],hit2[3];
587 //Float_t phi_0,kappa;
4aa41877 588 for(Int_t i=AliHLTTransform::GetFirstRow(fCurrentPatch); i<=AliHLTTransform::GetLastRow(fCurrentPatch); i++)
3e87ef69 589 {
4aa41877 590 Int_t rowindex1 = i - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 591 for(Int_t d1=0; d1<fData[rowindex1].fNPads; d1++)
3e87ef69 592 {
c62b480b 593 pad1 = d1 + fData[rowindex1].fMinpad;
3e87ef69 594 for(Int_t j=0; j<10; j++)
595 {
c62b480b 596 time1 = j + fData[rowindex1].fMintime;
597 charge1 = fData[rowindex1].fPads[d1][j];
3e87ef69 598 if(charge1==0) continue;
4aa41877 599 AliHLTTransform::Slice2Sector(0,i,sector,row);
600 AliHLTTransform::Raw2Local(hit,sector,row,pad1,time1);
3e87ef69 601 r1 = sqrt(hit[0]*hit[0]+hit[1]*hit[1])/2;
602 phi1 = atan2(hit[1],hit[0]);
603
4aa41877 604 for(Int_t j=i+1; j<=AliHLTTransform::GetLastRow(fCurrentPatch); j++)
3e87ef69 605 {
4aa41877 606 Int_t rowindex2 = j - AliHLTTransform::GetFirstRow(fCurrentPatch);
c62b480b 607 for(Int_t d2=0; d2<fData[rowindex2].fNPads; d2++)
3e87ef69 608 {
c62b480b 609 pad2 = d2 + fData[rowindex2].fMinpad;
3e87ef69 610 for(Int_t k=0; k<10; k++)
611 {
c62b480b 612 time2 = k + fData[rowindex2].fMintime;
613 charge2 = fData[rowindex2].fPads[d2][k];
3e87ef69 614 if(charge2==0) continue;
4aa41877 615 AliHLTTransform::Slice2Sector(0,j,sector,row);
616 AliHLTTransform::Raw2Local(hit2,sector,row,pad2,time2);
3e87ef69 617 r2 = sqrt(hit2[0]*hit2[0]+hit2[1]*hit2[1])/2;
618 phi2 = atan2(hit2[1],hit2[0]);
619 Float_t yc = (r2-(r1/cos(phi1))*cos(phi2))/(sin(phi2)-tan(phi1)*cos(phi2));
620 Float_t xc = r1/cos(phi1) - (r2*tan(phi1)-r1*sin(phi1)*cos(phi2))/(sin(phi2)-tan(phi1)*cos(phi2));
621 hist->Fill(xc,yc,charge1+charge2);
622 }
623 }
624 }
625 }
626 }
627 }
628}
629
630
4aa41877 631void AliHLTHoughTest::FindAbsMaxima(TH3 *hist,Int_t zsearch,Float_t &maxx,Float_t &maxy,Float_t &maxz,Int_t &maxvalue) const
3e87ef69 632{
c62b480b 633 //Find peaks in the Hough space
3e87ef69 634 TH1 *h1 = hist->Project3D("z");
635
636 TAxis *z = hist->GetZaxis();
637 Float_t zpeak[50];
638 Int_t n=0,i,zbin[50];
639 for(i=z->GetFirst()+1; i<z->GetLast()-1; i++)
640 {
641 int bin1 = h1->GetBin(i-1);
642 int bin2 = h1->GetBin(i);
643 int bin3 = h1->GetBin(i+1);
644 if(h1->GetBinContent(bin2) > h1->GetBinContent(bin1) && h1->GetBinContent(bin2) > h1->GetBinContent(bin3))
645 {
646 zbin[n]=bin2;
647 zpeak[n++]=h1->GetBinCenter(bin2);
648 }
649 }
650
651 Int_t zrange[2] = {z->GetFirst(),z->GetLast()};
652 z->SetRange(zbin[0]-zsearch,zbin[0]+zsearch);
653
654 TH1 *h2 = hist->Project3D("yx");
655 z->SetRange(zrange[0],zrange[1]);
656
657 TAxis *x = h2->GetXaxis();
658 TAxis *y = h2->GetYaxis();
659 maxvalue=0;
660 for(i=0; i<n; i++)
661 {
662 for(Int_t xbin=x->GetFirst(); xbin<=x->GetLast(); xbin++)
663 {
664 Float_t xvalue = x->GetBinCenter(xbin);
665 for(Int_t ybin=y->GetFirst(); ybin<=y->GetLast(); ybin++)
666 {
667 Float_t yvalue = y->GetBinCenter(ybin);
668 Int_t value = (Int_t)h2->GetBinContent(xbin,ybin);
669 if(value > maxvalue)
670 {
671 maxvalue=value;
c62b480b 672 maxx = xvalue;
673 maxy = yvalue;
674 maxz = zpeak[i];
3e87ef69 675 }
676 }
677 }
678 }
679
680 delete h1;
681 delete h2;
682}
683