]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCtracker.cxx
- Replace the TClonesArray of AliMUONTrackParam by a TObjArray in
[u/mrichter/AliRoot.git] / TPC / AliTPCtracker.cxx
CommitLineData
73042f01 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
0fa7dfa7 16/* $Id$ */
73042f01 17
18//-------------------------------------------------------
19// Implementation of the TPC tracker
20//
21// Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
22//-------------------------------------------------------
73042f01 23#include <TObjArray.h>
c630aafd 24#include <TError.h>
73042f01 25#include <TFile.h>
be5b9287 26#include <TTree.h>
c630aafd 27
af885e0f 28#include "AliESDEvent.h"
aad72f45 29#include "AliESDtrack.h"
be5b9287 30
31#include "AliTPCtracker.h"
32#include "AliTPCcluster.h"
b9de75e1 33#include "AliTPCParam.h"
c630aafd 34#include "AliClusters.h"
35
36ClassImp(AliTPCtracker)
83d73500 37
6c94f330 38//_____________________________________________________________________________
39AliTPCtracker::AliTPCtracker():
40 AliTracker(),
41 fkNIS(0),
42 fInnerSec(0),
43 fkNOS(0),
44 fOuterSec(0),
45 fN(0),
46 fSectors(0),
47 fParam(0),
48 fSeeds(0)
49{
50 //
51 // The default TPC tracker constructor
52 //
53}
54
b9de75e1 55//_____________________________________________________________________________
c630aafd 56AliTPCtracker::AliTPCtracker(const AliTPCParam *par):
6c94f330 57 AliTracker(),
58 fkNIS(par->GetNInnerSector()/2),
59 fInnerSec(new AliTPCSector[fkNIS]),
60 fkNOS(par->GetNOuterSector()/2),
61 fOuterSec(new AliTPCSector[fkNOS]),
62 fN(0),
63 fSectors(0),
64 fParam((AliTPCParam*) par),
65 fSeeds(0)
b9de75e1 66{
67 //---------------------------------------------------------------------
68 // The main TPC tracker constructor
69 //---------------------------------------------------------------------
b9de75e1 70
71 Int_t i;
72 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
73 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
74
b9de75e1 75}
76
77//_____________________________________________________________________________
78AliTPCtracker::~AliTPCtracker() {
79 //------------------------------------------------------------------
80 // TPC tracker destructor
81 //------------------------------------------------------------------
82 delete[] fInnerSec;
83 delete[] fOuterSec;
1c8d9aad 84 if (fSeeds) {
85 fSeeds->Delete();
86 delete fSeeds;
87 }
c630aafd 88}
89
73042f01 90//_____________________________________________________________________________
bdbd0f7a 91Double_t f1(Double_t x1,Double_t y1,
73042f01 92 Double_t x2,Double_t y2,
93 Double_t x3,Double_t y3)
94{
95 //-----------------------------------------------------------------
96 // Initial approximation of the track curvature
97 //-----------------------------------------------------------------
98 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
99 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
100 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
101 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
102 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
103
104 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
105
106 return -xr*yr/sqrt(xr*xr+yr*yr);
107}
108
109
110//_____________________________________________________________________________
bdbd0f7a 111Double_t f2(Double_t x1,Double_t y1,
73042f01 112 Double_t x2,Double_t y2,
113 Double_t x3,Double_t y3)
114{
115 //-----------------------------------------------------------------
116 // Initial approximation of the track curvature times center of curvature
117 //-----------------------------------------------------------------
118 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
119 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
120 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
121 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
122 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
123
124 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
125
126 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
127}
128
129//_____________________________________________________________________________
bdbd0f7a 130Double_t f3(Double_t x1,Double_t y1,
73042f01 131 Double_t x2,Double_t y2,
132 Double_t z1,Double_t z2)
133{
134 //-----------------------------------------------------------------
135 // Initial approximation of the tangent of the track dip angle
136 //-----------------------------------------------------------------
137 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
138}
139
83d73500 140//_____________________________________________________________________________
c630aafd 141Int_t AliTPCtracker::LoadClusters(TTree *cTree) {
b9de75e1 142 //-----------------------------------------------------------------
c630aafd 143 // This function loads TPC clusters.
b9de75e1 144 //-----------------------------------------------------------------
c630aafd 145 TBranch *branch=cTree->GetBranch("Segment");
146 if (!branch) {
147 Error("LoadClusters","Can't get the branch !");
148 return 1;
b9de75e1 149 }
150
c630aafd 151 AliClusters carray, *addr=&carray;
e7de6656 152 addr = new AliClusters("AliTPCcluster");
c630aafd 153 branch->SetAddress(&addr);
b9de75e1 154
c630aafd 155 Int_t nentr=(Int_t)cTree->GetEntries();
83d73500 156
c630aafd 157 for (Int_t i=0; i<nentr; i++) {
158 cTree->GetEvent(i);
bcc04e2a 159
c630aafd 160 Int_t ncl=carray.GetArray()->GetEntriesFast();
161
162 Int_t nir=fInnerSec->GetNRows(), nor=fOuterSec->GetNRows();
163 Int_t id=carray.GetID();
164 if ((id<0) || (id>2*(fkNIS*nir + fkNOS*nor))) {
df02fd67 165 Fatal("LoadClusters","Wrong index !");
bcc04e2a 166 }
c630aafd 167 Int_t outindex = 2*fkNIS*nir;
168 if (id<outindex) {
169 Int_t sec = id/nir;
170 Int_t row = id - sec*nir;
171 sec %= fkNIS;
172 AliTPCRow &padrow=fInnerSec[sec][row];
173 while (ncl--) {
174 AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
175 padrow.InsertCluster(c,sec,row);
176 }
177 } else {
178 id -= outindex;
179 Int_t sec = id/nor;
180 Int_t row = id - sec*nor;
181 sec %= fkNOS;
182 AliTPCRow &padrow=fOuterSec[sec][row];
183 while (ncl--) {
184 AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
185 padrow.InsertCluster(c,sec+fkNIS,row);
186 }
187 }
188 carray.GetArray()->Clear();
bcc04e2a 189 }
88cb7938 190
c630aafd 191 return 0;
b9de75e1 192}
193
194//_____________________________________________________________________________
c630aafd 195void AliTPCtracker::UnloadClusters() {
b9de75e1 196 //-----------------------------------------------------------------
c630aafd 197 // This function unloads TPC clusters.
b9de75e1 198 //-----------------------------------------------------------------
c630aafd 199 Int_t i;
200 for (i=0; i<fkNIS; i++) {
201 Int_t nr=fInnerSec->GetNRows();
202 for (Int_t n=0; n<nr; n++) fInnerSec[i][n].ResetClusters();
203 }
204 for (i=0; i<fkNOS; i++) {
205 Int_t nr=fOuterSec->GetNRows();
206 for (Int_t n=0; n<nr; n++) fOuterSec[i][n].ResetClusters();
b9de75e1 207 }
b9de75e1 208}
209
210//_____________________________________________________________________________
211Int_t AliTPCtracker::FollowProlongation(AliTPCseed& t, Int_t rf) {
73042f01 212 //-----------------------------------------------------------------
213 // This function tries to find a track prolongation.
214 //-----------------------------------------------------------------
b9de75e1 215 Double_t xt=t.GetX();
216 const Int_t kSKIP=(t.GetNumberOfClusters()<10) ? kRowsToSkip :
217 Int_t(0.5*fSectors->GetNRows());
73042f01 218 Int_t tryAgain=kSKIP;
73042f01 219
b9de75e1 220 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
221 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
222 if (alpha < 0. ) alpha += 2.*TMath::Pi();
223 Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
224
f03e3423 225 Int_t nrows=fSectors->GetRowNumber(xt)-1;
226 for (Int_t nr=nrows; nr>=rf; nr--) {
b9de75e1 227 Double_t x=fSectors->GetX(nr), ymax=fSectors->GetMaxY(nr);
73042f01 228 if (!t.PropagateTo(x)) return 0;
229
230 AliTPCcluster *cl=0;
231 UInt_t index=0;
b9de75e1 232 Double_t maxchi2=kMaxCHI2;
233 const AliTPCRow &krow=fSectors[s][nr];
6c23ffed 234 Double_t pt=t.GetSignedPt();
6d75e4b6 235 Double_t sy2=AliTPCcluster::SigmaY2(t.GetX(),t.GetTgl(),pt);
236 Double_t sz2=AliTPCcluster::SigmaZ2(t.GetX(),t.GetTgl());
73042f01 237 Double_t road=4.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
238
b9de75e1 239 if (road>kMaxROAD) {
73042f01 240 if (t.GetNumberOfClusters()>4)
c630aafd 241 Warning("FindProlongation","Too broad road !");
73042f01 242 return 0;
243 }
244
245 if (krow) {
246 for (Int_t i=krow.Find(y-road); i<krow; i++) {
c630aafd 247 AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
248 if (c->GetY() > y+road) break;
249 if (c->IsUsed()) continue;
250 if ((c->GetZ()-z)*(c->GetZ()-z) > 16.*(t.GetSigmaZ2()+sz2)) continue;
251 Double_t chi2=t.GetPredictedChi2(c);
252 if (chi2 > maxchi2) continue;
253 maxchi2=chi2;
254 cl=c;
73042f01 255 index=krow.GetIndex(i);
256 }
257 }
258 if (cl) {
b9de75e1 259 Float_t l=fSectors->GetPadPitchWidth();
024a7fe9 260 Float_t corr=1.; if (nr>63) corr=0.67; // new (third) pad response !
261 t.SetSampledEdx(cl->GetQ()/l*corr,t.GetNumberOfClusters());
be9c5115 262 if (!t.Update(cl,maxchi2,index)) {
263 if (!tryAgain--) return 0;
264 } else tryAgain=kSKIP;
73042f01 265 } else {
266 if (tryAgain==0) break;
267 if (y > ymax) {
b9de75e1 268 s = (s+1) % fN;
269 if (!t.Rotate(fSectors->GetAlpha())) return 0;
73042f01 270 } else if (y <-ymax) {
b9de75e1 271 s = (s-1+fN) % fN;
272 if (!t.Rotate(-fSectors->GetAlpha())) return 0;
73042f01 273 }
274 tryAgain--;
275 }
276 }
277
278 return 1;
b9de75e1 279}
c630aafd 280//_____________________________________________________________________________
281
282Int_t AliTPCtracker::FollowRefitInward(AliTPCseed *seed, AliTPCtrack *track) {
283 //
284 // This function propagates seed inward TPC using old clusters
285 // from the track.
286 //
287 // Sylwester Radomski, GSI
288 // 26.02.2003
289 //
290
291 // loop over rows
292
293 Int_t nRows = fSectors->GetNRows();
0fa7dfa7 294 for (Int_t iRow = nRows-1; iRow >= 0; iRow--) {
c630aafd 295
296 Double_t x = fSectors->GetX(iRow);
297 if (!seed->PropagateTo(x)) return 0;
298
299 // try to find an assigned cluster in this row
300
301 AliTPCcluster* cluster = NULL;
302 Int_t idx = -1;
303 Int_t sec = -1;
304 for (Int_t iCluster = 0; iCluster < track->GetNumberOfClusters(); iCluster++) {
305 idx = track->GetClusterIndex(iCluster);
306 sec = (idx&0xff000000)>>24;
307 Int_t row = (idx&0x00ff0000)>>16;
308 if (((fSectors == fInnerSec) && (sec >= fkNIS)) ||
309 ((fSectors == fOuterSec) && (sec < fkNIS))) continue;
310 if (row == iRow) {
311 cluster = (AliTPCcluster*) GetCluster(idx);
312 break;
313 }
314 }
315
316 // update the track seed with the found cluster
317
318 if (cluster) {
319 Double_t dAlpha = fParam->GetAngle(sec) - seed->GetAlpha();
320 if (TMath::Abs(dAlpha) > 0.0001) {
321 if (!seed->Rotate(dAlpha)) return 0;
322 if (!seed->PropagateTo(x)) return 0;
323 }
324
325 seed->Update(cluster, seed->GetPredictedChi2(cluster), idx);
326 }
327 }
328
329 return 1;
330}
b9de75e1 331
332//_____________________________________________________________________________
333Int_t AliTPCtracker::FollowBackProlongation
334(AliTPCseed& seed, const AliTPCtrack &track) {
335 //-----------------------------------------------------------------
336 // This function propagates tracks back through the TPC
337 //-----------------------------------------------------------------
338 Double_t alpha=seed.GetAlpha() - fSectors->GetAlphaShift();
339 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
340 if (alpha < 0. ) alpha += 2.*TMath::Pi();
341 Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
342
343 Int_t idx=-1, sec=-1, row=-1;
6c94f330 344 Int_t nc=track.GetNumberOfClusters();
c630aafd 345
b9de75e1 346 if (nc--) {
347 idx=track.GetClusterIndex(nc);
348 sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
349 }
c630aafd 350 if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
351 else { if (sec < fkNIS) row=-1; }
b9de75e1 352
353 Int_t nr=fSectors->GetNRows();
354 for (Int_t i=0; i<nr; i++) {
355 Double_t x=fSectors->GetX(i), ymax=fSectors->GetMaxY(i);
6c94f330 356 Double_t y;
357 if (!seed.GetYAt(x,GetBz(),y)) return 0;
c630aafd 358
b9de75e1 359 if (y > ymax) {
360 s = (s+1) % fN;
361 if (!seed.Rotate(fSectors->GetAlpha())) return 0;
362 } else if (y <-ymax) {
363 s = (s-1+fN) % fN;
364 if (!seed.Rotate(-fSectors->GetAlpha())) return 0;
365 }
366
c630aafd 367 if (!seed.PropagateTo(x)) return 0;
368
b9de75e1 369 AliTPCcluster *cl=0;
370 Int_t index=0;
371 Double_t maxchi2=kMaxCHI2;
6c23ffed 372 Double_t pt=seed.GetSignedPt();
6d75e4b6 373 Double_t sy2=AliTPCcluster::SigmaY2(seed.GetX(),seed.GetTgl(),pt);
374 Double_t sz2=AliTPCcluster::SigmaZ2(seed.GetX(),seed.GetTgl());
b9de75e1 375 Double_t road=4.*sqrt(seed.GetSigmaY2() + sy2), z=seed.GetZ();
376 if (road>kMaxROAD) {
c630aafd 377 Warning("FollowBackProlongation","Too broad road !");
b9de75e1 378 return 0;
379 }
380
b9de75e1 381 Int_t accepted=seed.GetNumberOfClusters();
382 if (row==i) {
383 //try to accept already found cluster
384 AliTPCcluster *c=(AliTPCcluster*)GetCluster(idx);
385 Double_t chi2;
386 if ((chi2=seed.GetPredictedChi2(c))<maxchi2 || accepted<27) {
387 index=idx; cl=c; maxchi2=chi2;
388 } //else cerr<<"AliTPCtracker::FollowBackProlongation: oulier !\n";
389
390 if (nc--) {
391 idx=track.GetClusterIndex(nc);
392 sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
393 }
c630aafd 394 if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
395 else { if (sec < fkNIS) row=-1; }
73042f01 396
b9de75e1 397 }
398 if (!cl) {
399 //try to fill the gap
400 const AliTPCRow &krow=fSectors[s][i];
401 if (accepted>27)
402 if (krow) {
51678892 403 for (Int_t icl=krow.Find(y-road); icl<krow; icl++) {
404 AliTPCcluster *c=(AliTPCcluster*)(krow[icl]);
c630aafd 405 if (c->GetY() > y+road) break;
406 if (c->IsUsed()) continue;
407 if ((c->GetZ()-z)*(c->GetZ()-z)>16.*(seed.GetSigmaZ2()+sz2)) continue;
408 Double_t chi2=seed.GetPredictedChi2(c);
409 if (chi2 > maxchi2) continue;
410 maxchi2=chi2;
411 cl=c;
51678892 412 index=krow.GetIndex(icl);
b9de75e1 413 }
414 }
415 }
416
417 if (cl) {
418 Float_t l=fSectors->GetPadPitchWidth();
024a7fe9 419 Float_t corr=1.; if (i>63) corr=0.67; // new (third) pad response !
420 seed.SetSampledEdx(cl->GetQ()/l*corr,seed.GetNumberOfClusters());
b9de75e1 421 seed.Update(cl,maxchi2,index);
422 }
423
424 }
425
b9de75e1 426 return 1;
73042f01 427}
428
429//_____________________________________________________________________________
b9de75e1 430void AliTPCtracker::MakeSeeds(Int_t i1, Int_t i2) {
73042f01 431 //-----------------------------------------------------------------
432 // This function creates track seeds.
433 //-----------------------------------------------------------------
434 Double_t x[5], c[15];
435
c630aafd 436 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
73042f01 437 Double_t cs=cos(alpha), sn=sin(alpha);
438
c630aafd 439 Double_t x1 =fSectors->GetX(i1);
440 Double_t xx2=fSectors->GetX(i2);
73042f01 441
c630aafd 442 for (Int_t ns=0; ns<fN; ns++) {
443 Int_t nl=fSectors[(ns-1+fN)%fN][i2];
444 Int_t nm=fSectors[ns][i2];
445 Int_t nu=fSectors[(ns+1)%fN][i2];
446 const AliTPCRow& kr1=fSectors[ns][i1];
73042f01 447 for (Int_t is=0; is < kr1; is++) {
448 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
449 for (Int_t js=0; js < nl+nm+nu; js++) {
c630aafd 450 const AliTPCcluster *kcl;
73042f01 451 Double_t x2, y2, z2;
7f6ddf58 452 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
73042f01 453
c630aafd 454 if (js<nl) {
455 const AliTPCRow& kr2=fSectors[(ns-1+fN)%fN][i2];
456 kcl=kr2[js];
73042f01 457 y2=kcl->GetY(); z2=kcl->GetZ();
458 x2= xx2*cs+y2*sn;
459 y2=-xx2*sn+y2*cs;
c630aafd 460 } else
461 if (js<nl+nm) {
462 const AliTPCRow& kr2=fSectors[ns][i2];
463 kcl=kr2[js-nl];
73042f01 464 x2=xx2; y2=kcl->GetY(); z2=kcl->GetZ();
c630aafd 465 } else {
466 const AliTPCRow& kr2=fSectors[(ns+1)%fN][i2];
467 kcl=kr2[js-nl-nm];
73042f01 468 y2=kcl->GetY(); z2=kcl->GetZ();
469 x2=xx2*cs-y2*sn;
470 y2=xx2*sn+y2*cs;
c630aafd 471 }
73042f01 472
7f6ddf58 473 Double_t zz=z1 - (z1-z3)/(x1-x3)*(x1-x2);
73042f01 474 if (TMath::Abs(zz-z2)>5.) continue;
475
476 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
c630aafd 477 if (d==0.) {
478 Warning("MakeSeeds","Straight seed !");
479 continue;
480 }
481 x[0]=y1;
482 x[1]=z1;
483 x[4]=f1(x1,y1,x2,y2,x3,y3);
484 if (TMath::Abs(x[4]) >= 0.0066) continue;
485 x[2]=f2(x1,y1,x2,y2,x3,y3);
486 //if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
487 x[3]=f3(x1,y1,x2,y2,z1,z2);
488 if (TMath::Abs(x[3]) > 1.2) continue;
489 Double_t a=asin(x[2]);
490 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
491 if (TMath::Abs(zv-z3)>10.) continue;
73042f01 492
493 Double_t sy1=kr1[is]->GetSigmaY2(), sz1=kr1[is]->GetSigmaZ2();
494 Double_t sy2=kcl->GetSigmaY2(), sz2=kcl->GetSigmaZ2();
024a7fe9 495 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
496 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
73042f01 497
c630aafd 498 Double_t f40=(f1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
499 Double_t f42=(f1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
500 Double_t f43=(f1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
501 Double_t f20=(f2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
502 Double_t f22=(f2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
503 Double_t f23=(f2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
504 Double_t f30=(f3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
505 Double_t f31=(f3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
506 Double_t f32=(f3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
507 Double_t f34=(f3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
73042f01 508
509 c[0]=sy1;
510 c[1]=0.; c[2]=sz1;
b9de75e1 511 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
512 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
513 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
514 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
515 c[13]=f30*sy1*f40+f32*sy2*f42;
516 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
73042f01 517
6c94f330 518 Int_t index=kr1.GetIndex(is);
519 AliTPCseed *track=new AliTPCseed(x1, ns*alpha+shift, x, c, index);
c630aafd 520 Float_t l=fSectors->GetPadPitchWidth();
73042f01 521 track->SetSampledEdx(kr1[is]->GetQ()/l,0);
522
b9de75e1 523 Int_t rc=FollowProlongation(*track, i2);
be9c5115 524 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track;
c630aafd 525 else fSeeds->AddLast(track);
73042f01 526 }
527 }
528 }
529}
530
531//_____________________________________________________________________________
b9de75e1 532Int_t AliTPCtracker::ReadSeeds(const TFile *inp) {
73042f01 533 //-----------------------------------------------------------------
b9de75e1 534 // This function reades track seeds.
73042f01 535 //-----------------------------------------------------------------
536 TDirectory *savedir=gDirectory;
537
b9de75e1 538 TFile *in=(TFile*)inp;
539 if (!in->IsOpen()) {
c630aafd 540 Error("ReadSeeds","Input file has not been open !");
be5b9287 541 return 1;
73042f01 542 }
543
b9de75e1 544 in->cd();
545 TTree *seedTree=(TTree*)in->Get("Seeds");
546 if (!seedTree) {
c630aafd 547 Error("ReadSeeds","Can't get a tree with track seeds !");
b9de75e1 548 return 2;
549 }
550 AliTPCtrack *seed=new AliTPCtrack;
551 seedTree->SetBranchAddress("tracks",&seed);
552
553 if (fSeeds==0) fSeeds=new TObjArray(15000);
73042f01 554
b9de75e1 555 Int_t n=(Int_t)seedTree->GetEntries();
556 for (Int_t i=0; i<n; i++) {
557 seedTree->GetEvent(i);
a91e73bc 558 seed->ResetClusters();
bbc6cd2c 559 fSeeds->AddLast(new AliTPCseed(*seed));
b9de75e1 560 }
561
562 delete seed;
f38c8ae5 563
564 delete seedTree; //Thanks to Mariana Bondila
565
b9de75e1 566 savedir->cd();
73042f01 567
b9de75e1 568 return 0;
569}
73042f01 570
cab3f0a3 571//_____________________________________________________________________________
af885e0f 572Int_t AliTPCtracker::Clusters2Tracks(AliESDEvent *event) {
b9de75e1 573 //-----------------------------------------------------------------
574 // This is a track finder.
c630aafd 575 // The clusters must be already loaded !
b9de75e1 576 //-----------------------------------------------------------------
b9de75e1 577
73042f01 578 //find track seeds
b9de75e1 579 Int_t nup=fOuterSec->GetNRows(), nlow=fInnerSec->GetNRows();
73042f01 580 Int_t nrows=nlow+nup;
b9de75e1 581 if (fSeeds==0) {
582 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
c630aafd 583 fSectors=fOuterSec; fN=fkNOS;
584 fSeeds=new TObjArray(15000);
b9de75e1 585 MakeSeeds(nup-1, nup-1-gap);
586 MakeSeeds(nup-1-shift, nup-1-shift-gap);
c630aafd 587 }
b9de75e1 588 fSeeds->Sort();
73042f01 589
b9de75e1 590 Int_t nseed=fSeeds->GetEntriesFast();
c630aafd 591 for (Int_t i=0; i<nseed; i++) {
592 //tracking in the outer sectors
593 fSectors=fOuterSec; fN=fkNOS;
594
b9de75e1 595 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
c630aafd 596 if (!FollowProlongation(t)) {
597 delete fSeeds->RemoveAt(i);
73042f01 598 continue;
599 }
c630aafd 600
601 //tracking in the inner sectors
602 fSectors=fInnerSec; fN=fkNIS;
603
604 Double_t alpha=t.GetAlpha() - fInnerSec->GetAlphaShift();
605 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
606 if (alpha < 0. ) alpha += 2.*TMath::Pi();
607 Int_t ns=Int_t(alpha/fInnerSec->GetAlpha())%fkNIS;
608
609 alpha=ns*fInnerSec->GetAlpha()+fInnerSec->GetAlphaShift()-t.GetAlpha();
610
611 if (t.Rotate(alpha)) {
612 if (FollowProlongation(t)) {
613 if (t.GetNumberOfClusters() >= Int_t(0.4*nrows)) {
614 t.CookdEdx();
615 CookLabel(pt,0.1); //For comparison only
616 pt->PropagateTo(fParam->GetInnerRadiusLow());
617 AliESDtrack iotrack;
618 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
619
620 event->AddTrack(&iotrack);
621
622 UseClusters(&t);
623 }
624 }
625 }
88cb7938 626 delete fSeeds->RemoveAt(i);
c630aafd 627 }
628
629 Info("Clusters2Tracks","Number of found tracks : %d",
630 event->GetNumberOfTracks());
631
632 fSeeds->Clear(); delete fSeeds; fSeeds=0;
633
634 return 0;
635}
636
c630aafd 637//_____________________________________________________________________________
af885e0f 638Int_t AliTPCtracker::RefitInward(AliESDEvent* event) {
c630aafd 639 //
640 // The function propagates tracks throught TPC inward
641 // using already associated clusters.
642 // The clusters must be already loaded !
643 //
644
645 Int_t nTracks = event->GetNumberOfTracks();
646 Int_t nRefited = 0;
647
648 for (Int_t i = 0; i < nTracks; i++) {
649 AliESDtrack* track = event->GetTrack(i);
650 ULong_t status = track->GetStatus();
651
4ec88b9c 652 if ( (status & AliESDtrack::kTPCrefit) != 0 ) continue;
c630aafd 653 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
4ec88b9c 654
655 if ( (status & AliESDtrack::kTRDout ) != 0 )
656 if ( (status & AliESDtrack::kTRDrefit ) == 0 ) continue;
657
c630aafd 658 AliTPCtrack* tpcTrack = new AliTPCtrack(*track);
a91e73bc 659 AliTPCseed* seed=new AliTPCseed(*tpcTrack); seed->ResetClusters();
08dba16b 660
6c94f330 661 if ( (status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance(10.);
4ec88b9c 662
c630aafd 663 fSectors = fOuterSec;
664
665 Int_t res = FollowRefitInward(seed, tpcTrack);
666 UseClusters(seed);
667 Int_t nc = seed->GetNumberOfClusters();
668
669 fSectors = fInnerSec;
670
671 res = FollowRefitInward(seed, tpcTrack);
672 UseClusters(seed, nc);
673
674 if (res) {
675 seed->PropagateTo(fParam->GetInnerRadiusLow());
676 seed->SetLabel(tpcTrack->GetLabel());
677 seed->SetdEdx(tpcTrack->GetdEdx());
678 track->UpdateTrackParams(seed, AliESDtrack::kTPCrefit);
679 nRefited++;
680 }
681
682 delete seed;
683 delete tpcTrack;
684 }
685
686 Info("RefitInward","Number of refitted tracks : %d",nRefited);
687
688 return 0;
689}
690
af885e0f 691Int_t AliTPCtracker::PropagateBack(AliESDEvent *event) {
c630aafd 692 //-----------------------------------------------------------------
693 // This function propagates tracks back through the TPC.
694 // The clusters must be already loaded !
695 //-----------------------------------------------------------------
696 Int_t nentr=event->GetNumberOfTracks();
697 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
698
699 Int_t ntrk=0;
700 for (Int_t i=0; i<nentr; i++) {
701 AliESDtrack *esd=event->GetTrack(i);
702 ULong_t status=esd->GetStatus();
703
704 if ( (status & AliESDtrack::kTPCin ) == 0 ) continue;
705 if ( (status & AliESDtrack::kTPCout) != 0 ) continue;
14718359 706 if ( (status & AliESDtrack::kITSin) != 0 )
707 if ( (status & AliESDtrack::kITSout) == 0 ) continue;
c630aafd 708
709 const AliTPCtrack t(*esd);
a91e73bc 710 AliTPCseed s(t); s.ResetClusters();
c630aafd 711
6c94f330 712 if ( (status & AliESDtrack::kITSout) == 0 ) s.ResetCovariance(10.);
c630aafd 713
714 //inner sectors
715 fSectors=fInnerSec; fN=fkNIS;
716
717 Double_t alpha=s.GetAlpha() - fSectors->GetAlphaShift();
718 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
719 if (alpha < 0. ) alpha += 2.*TMath::Pi();
720 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
721 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
722 alpha-=s.GetAlpha();
723
724 if (!s.Rotate(alpha)) continue;
725 if (!FollowBackProlongation(s,t)) continue;
726
727 UseClusters(&s);
728
729 //outer sectors
730 fSectors=fOuterSec; fN=fkNOS;
731
6c94f330 732 Int_t nc=s.GetNumberOfClusters();
c630aafd 733
734 alpha=s.GetAlpha() - fSectors->GetAlphaShift();
735 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
736 if (alpha < 0. ) alpha += 2.*TMath::Pi();
737 ns=Int_t(alpha/fSectors->GetAlpha())%fN;
738
739 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
740 alpha-=s.GetAlpha();
741
742 if (!s.Rotate(alpha)) continue;
743 if (!FollowBackProlongation(s,t)) continue;
744 {
745 Int_t nrows=fOuterSec->GetNRows()+fInnerSec->GetNRows();
746 if (s.GetNumberOfClusters() < Int_t(0.4*nrows)) continue;
747 }
748 s.PropagateTo(fParam->GetOuterRadiusUp());
749 s.CookdEdx();
750 CookLabel(&s,0.1); //For comparison only
751 UseClusters(&s,nc);
752 esd->UpdateTrackParams(&s,AliESDtrack::kTPCout);
753 ntrk++;
754 }
755 Info("PropagateBack","Number of back propagated tracks: %d",ntrk);
08dba16b 756
757 return 0;
758}
b9de75e1 759
b9de75e1 760//_________________________________________________________________________
761AliCluster *AliTPCtracker::GetCluster(Int_t index) const {
762 //--------------------------------------------------------------------
763 // Return pointer to a given cluster
764 //--------------------------------------------------------------------
765 Int_t sec=(index&0xff000000)>>24;
766 Int_t row=(index&0x00ff0000)>>16;
767 Int_t ncl=(index&0x0000ffff)>>00;
768
c630aafd 769 const AliTPCcluster *cl=0;
770
771 if (sec<fkNIS) {
772 cl=fInnerSec[sec][row].GetUnsortedCluster(ncl);
773 } else {
774 sec-=fkNIS;
775 cl=fOuterSec[sec][row].GetUnsortedCluster(ncl);
776 }
777
778 return (AliCluster*)cl;
b9de75e1 779}
780
781//__________________________________________________________________________
782void AliTPCtracker::CookLabel(AliKalmanTrack *t, Float_t wrong) const {
783 //--------------------------------------------------------------------
784 //This function "cooks" a track label. If label<0, this track is fake.
785 //--------------------------------------------------------------------
786 Int_t noc=t->GetNumberOfClusters();
787 Int_t *lb=new Int_t[noc];
788 Int_t *mx=new Int_t[noc];
789 AliCluster **clusters=new AliCluster*[noc];
790
791 Int_t i;
792 for (i=0; i<noc; i++) {
793 lb[i]=mx[i]=0;
794 Int_t index=t->GetClusterIndex(i);
795 clusters[i]=GetCluster(index);
796 }
797
798 Int_t lab=123456789;
799 for (i=0; i<noc; i++) {
800 AliCluster *c=clusters[i];
801 lab=TMath::Abs(c->GetLabel(0));
802 Int_t j;
803 for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
804 lb[j]=lab;
805 (mx[j])++;
806 }
807
808 Int_t max=0;
809 for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
810
811 for (i=0; i<noc; i++) {
812 AliCluster *c=clusters[i];
813 if (TMath::Abs(c->GetLabel(1)) == lab ||
814 TMath::Abs(c->GetLabel(2)) == lab ) max++;
815 }
816
817 if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
818
819 else {
820 Int_t tail=Int_t(0.10*noc);
821 max=0;
822 for (i=1; i<=tail; i++) {
823 AliCluster *c=clusters[noc-i];
824 if (lab == TMath::Abs(c->GetLabel(0)) ||
825 lab == TMath::Abs(c->GetLabel(1)) ||
826 lab == TMath::Abs(c->GetLabel(2))) max++;
827 }
828 if (max < Int_t(0.5*tail)) lab=-lab;
829 }
830
831 t->SetLabel(lab);
832
833 delete[] lb;
834 delete[] mx;
835 delete[] clusters;
836}
837
838//_________________________________________________________________________
839void AliTPCtracker::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) {
840 //-----------------------------------------------------------------------
841 // Setup inner sector
842 //-----------------------------------------------------------------------
843 if (f==0) {
844 fAlpha=par->GetInnerAngle();
845 fAlphaShift=par->GetInnerAngleShift();
846 fPadPitchWidth=par->GetInnerPadPitchWidth();
024a7fe9 847 f1PadPitchLength=par->GetInnerPadPitchLength();
848 f2PadPitchLength=f1PadPitchLength;
b9de75e1 849
850 fN=par->GetNRowLow();
851 fRow=new AliTPCRow[fN];
852 for (Int_t i=0; i<fN; i++) fRow[i].SetX(par->GetPadRowRadiiLow(i));
853 } else {
854 fAlpha=par->GetOuterAngle();
855 fAlphaShift=par->GetOuterAngleShift();
856 fPadPitchWidth=par->GetOuterPadPitchWidth();
f03e3423 857 f1PadPitchLength=par->GetOuter1PadPitchLength();
858 f2PadPitchLength=par->GetOuter2PadPitchLength();
b9de75e1 859
860 fN=par->GetNRowUp();
861 fRow=new AliTPCRow[fN];
f03e3423 862 for (Int_t i=0; i<fN; i++){
863 fRow[i].SetX(par->GetPadRowRadiiUp(i));
864 }
b9de75e1 865 }
866}
867
73042f01 868//_________________________________________________________________________
c630aafd 869void AliTPCtracker::
870AliTPCRow::InsertCluster(const AliTPCcluster* c, Int_t sec, Int_t row) {
73042f01 871 //-----------------------------------------------------------------------
872 // Insert a cluster into this pad row in accordence with its y-coordinate
873 //-----------------------------------------------------------------------
b9de75e1 874 if (fN==kMaxClusterPerRow) {
c630aafd 875 ::Error("InsertCluster","Too many clusters !");
876 return;
877 }
878
879 Int_t index=(((sec<<8)+row)<<16)+fN;
880
881 if (fN==0) {
882 fSize=kMaxClusterPerRow/8;
883 fClusterArray=new AliTPCcluster[fSize];
884 fIndex[0]=index;
885 fClusterArray[0]=*c; fClusters[fN++]=fClusterArray;
886 return;
73042f01 887 }
c630aafd 888
889 if (fN==fSize) {
890 Int_t size=fSize*2;
891 AliTPCcluster *buff=new AliTPCcluster[size];
892 memcpy(buff,fClusterArray,fSize*sizeof(AliTPCcluster));
893 for (Int_t i=0; i<fN; i++)
894 fClusters[i]=buff+(fClusters[i]-fClusterArray);
895 delete[] fClusterArray;
896 fClusterArray=buff;
897 fSize=size;
898 }
899
73042f01 900 Int_t i=Find(c->GetY());
901 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCcluster*));
902 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
c630aafd 903 fIndex[i]=index;
904 fClusters[i]=fClusterArray+fN; fClusterArray[fN++]=*c;
73042f01 905}
906
907//___________________________________________________________________
908Int_t AliTPCtracker::AliTPCRow::Find(Double_t y) const {
909 //-----------------------------------------------------------------------
910 // Return the index of the nearest cluster
911 //-----------------------------------------------------------------------
11c7b548 912 if(fN<=0) return 0;
73042f01 913 if (y <= fClusters[0]->GetY()) return 0;
914 if (y > fClusters[fN-1]->GetY()) return fN;
915 Int_t b=0, e=fN-1, m=(b+e)/2;
916 for (; b<e; m=(b+e)/2) {
917 if (y > fClusters[m]->GetY()) b=m+1;
918 else e=m;
919 }
920 return m;
921}
922
923//_____________________________________________________________________________
924void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) {
925 //-----------------------------------------------------------------
926 // This funtion calculates dE/dX within the "low" and "up" cuts.
927 //-----------------------------------------------------------------
928 Int_t i;
be9c5115 929 Int_t nc=GetNumberOfClusters();
c630aafd 930
931 Int_t swap;//stupid sorting
932 do {
933 swap=0;
934 for (i=0; i<nc-1; i++) {
935 if (fdEdxSample[i]<=fdEdxSample[i+1]) continue;
936 Float_t tmp=fdEdxSample[i];
937 fdEdxSample[i]=fdEdxSample[i+1]; fdEdxSample[i+1]=tmp;
938 swap++;
939 }
940 } while (swap);
73042f01 941
be9c5115 942 Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
73042f01 943 Float_t dedx=0;
c630aafd 944 for (i=nl; i<=nu; i++) dedx += fdEdxSample[i];
73042f01 945 dedx /= (nu-nl+1);
946 SetdEdx(dedx);
7f6ddf58 947
14718359 948 return;
73042f01 949}
950
73042f01 951
952