]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITStrackerV2.cxx
Faster version of the tracker V2
[u/mrichter/AliRoot.git] / ITS / AliITStrackerV2.cxx
CommitLineData
006b5f7f 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
16//-------------------------------------------------------------------------
17// Implementation of the ITS tracker class
18856a77 18// It reads AliITSclusterV2 clusters and creates AliITStrackV2 tracks
19// and fills with them the ESD
006b5f7f 20// Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
a9a2d814 21// dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
006b5f7f 22//-------------------------------------------------------------------------
1f9a65c4 23
9c369bfe 24#include <new>
25
006b5f7f 26#include <TFile.h>
27#include <TTree.h>
28#include <TRandom.h>
006b5f7f 29
30#include "AliITSgeom.h"
31#include "AliITSRecPoint.h"
c630aafd 32#include "AliESD.h"
006b5f7f 33#include "AliITSclusterV2.h"
34#include "AliITStrackerV2.h"
35
8676d691 36ClassImp(AliITStrackerV2)
006b5f7f 37
18856a77 38AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[kMaxLayer]; // ITS layers
006b5f7f 39
c630aafd 40AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
006b5f7f 41 //--------------------------------------------------------------------
61ab8ea8 42 //This is the AliITStrackerV2 constructor
006b5f7f 43 //--------------------------------------------------------------------
006b5f7f 44 AliITSgeom *g=(AliITSgeom*)geom;
45
c219fd63 46 Float_t x,y,z; Int_t i;
c630aafd 47 for (i=1; i<kMaxLayer+1; i++) {
006b5f7f 48 Int_t nlad=g->GetNladders(i);
49 Int_t ndet=g->GetNdetectors(i);
50
51 g->GetTrans(i,1,1,x,y,z);
52 Double_t r=TMath::Sqrt(x*x + y*y);
53 Double_t poff=TMath::ATan2(y,x);
54 Double_t zoff=z;
55
56 g->GetTrans(i,1,2,x,y,z);
57 r += TMath::Sqrt(x*x + y*y);
58 g->GetTrans(i,2,1,x,y,z);
59 r += TMath::Sqrt(x*x + y*y);
60 g->GetTrans(i,2,2,x,y,z);
61 r += TMath::Sqrt(x*x + y*y);
62 r*=0.25;
63
18856a77 64 new (fgLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
006b5f7f 65
66 for (Int_t j=1; j<nlad+1; j++) {
67 for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
68 Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift);
69 Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
70
f363d377 71 Double_t phi=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
72 phi+=TMath::Pi()/2;
73 if (i==1) phi+=TMath::Pi();
c219fd63 74
75 if (phi<0) phi+=TMath::TwoPi();
76 else if (phi>=TMath::TwoPi()) phi-=TMath::TwoPi();
77
f363d377 78 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
79 Double_t r=x*cp+y*sp;
006b5f7f 80
18856a77 81 AliITSdetector &det=fgLayers[i-1].GetDetector((j-1)*ndet + k-1);
006b5f7f 82 new(&det) AliITSdetector(r,phi);
83 }
84 }
c630aafd 85
006b5f7f 86 }
87
61ab8ea8 88 fI=kMaxLayer;
743a19f2 89
61ab8ea8 90 fPass=0;
c219fd63 91 fConstraint[0]=1; fConstraint[1]=0;
8676d691 92
93 Double_t xyz[]={kXV,kYV,kZV}, ers[]={kSigmaXV,kSigmaYV,kSigmaZV};
94 SetVertex(xyz,ers);
c0dd5278 95
96 for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=kLayersNotToSkip[i];
97 fLastLayerToTrackTo=kLastLayerToTrackTo;
98
99}
100
101void AliITStrackerV2::SetLayersNotToSkip(Int_t *l) {
102 //--------------------------------------------------------------------
103 //This function set masks of the layers which must be not skipped
104 //--------------------------------------------------------------------
105 for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=l[i];
61ab8ea8 106}
006b5f7f 107
c630aafd 108Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
61ab8ea8 109 //--------------------------------------------------------------------
110 //This function loads ITS clusters
111 //--------------------------------------------------------------------
61ab8ea8 112 TBranch *branch=cTree->GetBranch("Clusters");
113 if (!branch) {
c630aafd 114 Error("LoadClusters"," can't get the branch !\n");
8676d691 115 return 1;
61ab8ea8 116 }
006b5f7f 117
61ab8ea8 118 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
119 branch->SetAddress(&clusters);
006b5f7f 120
61ab8ea8 121 Int_t j=0;
122 for (Int_t i=0; i<kMaxLayer; i++) {
18856a77 123 Int_t ndet=fgLayers[i].GetNdetectors();
124 Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
c219fd63 125
126 Double_t r=fgLayers[i].GetR();
127 Double_t circ=TMath::TwoPi()*r;
128
61ab8ea8 129 for (; j<jmax; j++) {
130 if (!cTree->GetEvent(j)) continue;
131 Int_t ncl=clusters->GetEntriesFast();
132 while (ncl--) {
133 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
c219fd63 134
135 Int_t idx=c->GetDetectorIndex();
136 Double_t y=r*fgLayers[i].GetDetector(idx).GetPhi()+c->GetY();
137 if (y>circ) y-=circ; else if (y<0) y+=circ;
138 c->SetPhiR(y);
139
18856a77 140 fgLayers[i].InsertCluster(new AliITSclusterV2(*c));
61ab8ea8 141 }
142 clusters->Delete();
143 }
18856a77 144 fgLayers[i].ResetRoad(); //road defined by the cluster density
006b5f7f 145 }
c630aafd 146
8676d691 147 return 0;
61ab8ea8 148}
006b5f7f 149
61ab8ea8 150void AliITStrackerV2::UnloadClusters() {
151 //--------------------------------------------------------------------
152 //This function unloads ITS clusters
153 //--------------------------------------------------------------------
18856a77 154 for (Int_t i=0; i<kMaxLayer; i++) fgLayers[i].ResetClusters();
006b5f7f 155}
156
880f41b9 157static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
158 //--------------------------------------------------------------------
159 // Correction for the material between the TPC and the ITS
160 // (should it belong to the TPC code ?)
161 //--------------------------------------------------------------------
162 Double_t riw=80., diw=0.0053, x0iw=30; // TPC inner wall ?
163 Double_t rcd=61., dcd=0.0053, x0cd=30; // TPC "central drum" ?
164 Double_t yr=12.8, dr=0.03; // rods ?
165 Double_t zm=0.2, dm=0.40; // membrane
166 //Double_t rr=52., dr=0.19, x0r=24., yyr=7.77; //rails
167 Double_t rs=50., ds=0.001; // something belonging to the ITS (screen ?)
168
169 if (t->GetX() > riw) {
170 if (!t->PropagateTo(riw,diw,x0iw)) return 1;
171 if (TMath::Abs(t->GetY())>yr) t->CorrectForMaterial(dr);
172 if (TMath::Abs(t->GetZ())<zm) t->CorrectForMaterial(dm);
173 if (!t->PropagateTo(rcd,dcd,x0cd)) return 1;
174 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
175 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,dr,x0r);
176 if (!t->PropagateTo(rs,ds)) return 1;
177 } else if (t->GetX() < rs) {
178 if (!t->PropagateTo(rs,-ds)) return 1;
179 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
180 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,-dr,x0r);
181 if (!t->PropagateTo(rcd,-dcd,x0cd)) return 1;
ba1a3e1a 182 if (!t->PropagateTo(riw+0.001,-diw,x0iw)) return 1;
880f41b9 183 } else {
184 ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
185 return 1;
186 }
187
188 return 0;
189}
190
c630aafd 191Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
006b5f7f 192 //--------------------------------------------------------------------
c630aafd 193 // This functions reconstructs ITS tracks
194 // The clusters must be already loaded !
006b5f7f 195 //--------------------------------------------------------------------
c630aafd 196 TObjArray itsTracks(15000);
006b5f7f 197
c630aafd 198 {/* Read ESD tracks */
199 Int_t nentr=event->GetNumberOfTracks();
200 Info("Clusters2Tracks", "Number of ESD tracks: %d\n", nentr);
201 while (nentr--) {
202 AliESDtrack *esd=event->GetTrack(nentr);
203
ba1a3e1a 204 if ((esd->GetStatus()&AliESDtrack::kTPCin)==0) continue;
205 if (esd->GetStatus()&AliESDtrack::kTPCout) continue;
206 if (esd->GetStatus()&AliESDtrack::kITSin) continue;
c630aafd 207
208 AliITStrackV2 *t=0;
209 try {
210 t=new AliITStrackV2(*esd);
211 } catch (const Char_t *msg) {
212 Warning("Clusters2Tracks",msg);
213 delete t;
214 continue;
215 }
15dd636f 216 if (TMath::Abs(t->GetD())>4) {
ebd8369c 217 delete t;
218 continue;
219 }
006b5f7f 220
c630aafd 221 if (CorrectForDeadZoneMaterial(t)!=0) {
222 Warning("Clusters2Tracks",
223 "failed to correct for the material in the dead zone !\n");
224 delete t;
225 continue;
226 }
227 itsTracks.AddLast(t);
228 }
229 } /* End Read ESD tracks */
230
231 itsTracks.Sort();
232 Int_t nentr=itsTracks.GetEntriesFast();
15dd636f 233
c630aafd 234 Int_t ntrk=0;
235 for (fPass=0; fPass<2; fPass++) {
236 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
237 for (Int_t i=0; i<nentr; i++) {
238 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
15dd636f 239 if (t==0) continue; //this track has been already tracked
c630aafd 240 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
241
242 ResetTrackToFollow(*t);
243 ResetBestTrack();
244
245 for (FollowProlongation(); fI<kMaxLayer; fI++) {
15dd636f 246 while (TakeNextProlongation()) FollowProlongation();
c630aafd 247 }
248
de5d3595 249 if (fBestTrack.GetNumberOfClusters() == 0) continue;
15dd636f 250
c630aafd 251 if (fConstraint[fPass]) {
15dd636f 252 ResetTrackToFollow(*t);
253 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
254 ResetBestTrack();
255 }
256
c630aafd 257 fBestTrack.SetLabel(tpcLabel);
258 fBestTrack.CookdEdx();
259 CookLabel(&fBestTrack,0.); //For comparison only
260 fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
15dd636f 261 UseClusters(&fBestTrack);
c630aafd 262 delete itsTracks.RemoveAt(i);
263 ntrk++;
264 }
265 }
babd135a 266
c630aafd 267 itsTracks.Delete();
2257f27e 268
c630aafd 269 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
270
271 return 0;
272}
273
c630aafd 274Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
275 //--------------------------------------------------------------------
276 // This functions propagates reconstructed ITS tracks back
277 // The clusters must be loaded !
278 //--------------------------------------------------------------------
279 Int_t nentr=event->GetNumberOfTracks();
280 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
281
282 Int_t ntrk=0;
283 for (Int_t i=0; i<nentr; i++) {
284 AliESDtrack *esd=event->GetTrack(i);
285
ba1a3e1a 286 if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
287 if (esd->GetStatus()&AliESDtrack::kITSout) continue;
c630aafd 288
289 AliITStrackV2 *t=0;
290 try {
291 t=new AliITStrackV2(*esd);
292 } catch (const Char_t *msg) {
293 Warning("PropagateBack",msg);
294 delete t;
295 continue;
296 }
297
298 ResetTrackToFollow(*t);
299
300 // propagete to vertex [SR, GSI 17.02.2003]
a9459f9c 301 // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
302 if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
303 if (fTrackToFollow.PropagateToVertex()) {
304 fTrackToFollow.StartTimeIntegral();
305 }
306 fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
307 }
c630aafd 308
309 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
310 if (RefitAt(49.,&fTrackToFollow,t)) {
311 if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
312 Warning("PropagateBack",
313 "failed to correct for the material in the dead zone !\n");
314 delete t;
315 continue;
316 }
317 fTrackToFollow.SetLabel(t->GetLabel());
ba1a3e1a 318 //fTrackToFollow.CookdEdx();
c630aafd 319 CookLabel(&fTrackToFollow,0.); //For comparison only
320 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
15dd636f 321 UseClusters(&fTrackToFollow);
c630aafd 322 ntrk++;
323 }
324 delete t;
325 }
326
327 Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
61ab8ea8 328
83d73500 329 return 0;
330}
331
c630aafd 332Int_t AliITStrackerV2::RefitInward(AliESD *event) {
83d73500 333 //--------------------------------------------------------------------
334 // This functions refits ITS tracks using the
335 // "inward propagated" TPC tracks
c630aafd 336 // The clusters must be loaded !
83d73500 337 //--------------------------------------------------------------------
c630aafd 338 Int_t nentr=event->GetNumberOfTracks();
339 Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
340
341 Int_t ntrk=0;
342 for (Int_t i=0; i<nentr; i++) {
343 AliESDtrack *esd=event->GetTrack(i);
344
ba1a3e1a 345 if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
346 if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
347 if (esd->GetStatus()&AliESDtrack::kTPCout)
15dd636f 348 if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
c630aafd 349
350 AliITStrackV2 *t=0;
351 try {
352 t=new AliITStrackV2(*esd);
353 } catch (const Char_t *msg) {
354 Warning("RefitInward",msg);
355 delete t;
356 continue;
357 }
358
359 if (CorrectForDeadZoneMaterial(t)!=0) {
360 Warning("RefitInward",
361 "failed to correct for the material in the dead zone !\n");
362 delete t;
363 continue;
364 }
365
366 ResetTrackToFollow(*t);
367 fTrackToFollow.ResetClusters();
368
369 //Refitting...
370 if (RefitAt(3.7, &fTrackToFollow, t)) {
371 fTrackToFollow.SetLabel(t->GetLabel());
372 fTrackToFollow.CookdEdx();
15dd636f 373 CookLabel(&fTrackToFollow,0.); //For comparison only
18856a77 374
67c3dcbe 375 if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe
376 Double_t a=fTrackToFollow.GetAlpha();
377 Double_t cs=TMath::Cos(a),sn=TMath::Sin(a);
378 Double_t xv= GetX()*cs + GetY()*sn;
379 Double_t yv=-GetX()*sn + GetY()*cs;
380
381 Double_t c=fTrackToFollow.GetC(), snp=fTrackToFollow.GetSnp();
382 Double_t x=fTrackToFollow.GetX(), y=fTrackToFollow.GetY();
383 Double_t tgfv=-(c*(x-xv)-snp)/(c*(y-yv) + TMath::Sqrt(1.-snp*snp));
384 Double_t fv=TMath::ATan(tgfv);
385
386 cs=TMath::Cos(fv); sn=TMath::Sin(fv);
387 x = xv*cs + yv*sn;
388 yv=-xv*sn + yv*cs; xv=x;
389
390 if (fTrackToFollow.Propagate(fv+a,xv)) {
df02fd67 391 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
399fb957 392 Float_t d=fTrackToFollow.GetD(GetX(),GetY());
393 Float_t z=fTrackToFollow.GetZ()-GetZ();
394 fTrackToFollow.GetESDtrack()->SetImpactParameters(d,z);
15dd636f 395 UseClusters(&fTrackToFollow);
67c3dcbe 396 {
397 AliITSclusterV2 c; c.SetY(yv); c.SetZ(GetZ());
398 c.SetSigmaY2(GetSigmaY()*GetSigmaY());
399 c.SetSigmaZ2(GetSigmaZ()*GetSigmaZ());
400 Double_t chi2=fTrackToFollow.GetPredictedChi2(&c);
401 if (chi2<kMaxChi2)
402 if (fTrackToFollow.Update(&c,-chi2,0))
403 fTrackToFollow.SetConstrainedESDtrack(chi2);
404 }
df02fd67 405 ntrk++;
406 }
67c3dcbe 407 }
c630aafd 408 }
409 delete t;
410 }
411
412 Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
83d73500 413
c630aafd 414 return 0;
415}
83d73500 416
006b5f7f 417AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
418 //--------------------------------------------------------------------
419 // Return pointer to a given cluster
420 //--------------------------------------------------------------------
421 Int_t l=(index & 0xf0000000) >> 28;
422 Int_t c=(index & 0x0fffffff) >> 00;
18856a77 423 return fgLayers[l].GetCluster(c);
006b5f7f 424}
425
426
427void AliITStrackerV2::FollowProlongation() {
428 //--------------------------------------------------------------------
429 //This function finds a track prolongation
430 //--------------------------------------------------------------------
c0dd5278 431 while (fI>fLastLayerToTrackTo) {
a9a2d814 432 Int_t i=fI-1;
8676d691 433
18856a77 434 AliITSlayer &layer=fgLayers[i];
a9a2d814 435 AliITStrackV2 &track=fTracks[i];
006b5f7f 436
14825d5a 437 Double_t r=layer.GetR();
61ab8ea8 438
a9a2d814 439 if (i==3 || i==1) {
18856a77 440 Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
61ab8ea8 441 Double_t d=0.0034, x0=38.6;
442 if (i==1) {rs=9.; d=0.0097; x0=42;}
443 if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
8676d691 444 //Warning("FollowProlongation","propagation failed !\n");
c0dd5278 445 return;
a9a2d814 446 }
006b5f7f 447 }
448
449 //find intersection
450 Double_t x,y,z;
14825d5a 451 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
8676d691 452 //Warning("FollowProlongation","failed to estimate track !\n");
c0dd5278 453 return;
006b5f7f 454 }
455 Double_t phi=TMath::ATan2(y,x);
f363d377 456
006b5f7f 457 Int_t idet=layer.FindDetectorIndex(phi,z);
458 if (idet<0) {
8676d691 459 //Warning("FollowProlongation","failed to find a detector !\n");
c0dd5278 460 return;
006b5f7f 461 }
462
463 //propagate to the intersection
464 const AliITSdetector &det=layer.GetDetector(idet);
14825d5a 465 phi=det.GetPhi();
a9a2d814 466 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
8676d691 467 //Warning("FollowProlongation","propagation failed !\n");
c0dd5278 468 return;
006b5f7f 469 }
470 fTrackToFollow.SetDetectorIndex(idet);
471
472 //Select possible prolongations and store the current track estimation
473 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
8676d691 474 Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
475 Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
476 Double_t road=layer.GetRoad();
477 if (dz*dy>road*road) {
478 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
479 dz=road*scz; dy=road*scy;
480 }
b87bd7cd 481
482 //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
14825d5a 483 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
7f6ddf58 484 if (dz > kMaxRoad) {
8676d691 485 //Warning("FollowProlongation","too broad road in Z !\n");
c0dd5278 486 return;
006b5f7f 487 }
a9a2d814 488
15dd636f 489 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
a9a2d814 490
b87bd7cd 491 //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
14825d5a 492 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
7f6ddf58 493 if (dy > kMaxRoad) {
8676d691 494 //Warning("FollowProlongation","too broad road in Y !\n");
c0dd5278 495 return;
006b5f7f 496 }
a9a2d814 497
c219fd63 498 fI--;
499
7f6ddf58 500 Double_t zmin=track.GetZ() - dz;
006b5f7f 501 Double_t zmax=track.GetZ() + dz;
14825d5a 502 Double_t ymin=track.GetY() + r*phi - dy;
503 Double_t ymax=track.GetY() + r*phi + dy;
c219fd63 504 if (layer.SelectClusters(zmin,zmax,ymin,ymax)==0)
505 if (fLayersNotToSkip[fI]) return;
006b5f7f 506
15dd636f 507 if (!TakeNextProlongation())
508 if (fLayersNotToSkip[fI]) return;
509
006b5f7f 510 }
511
512 //deal with the best track
513 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
514 Int_t nclb=fBestTrack.GetNumberOfClusters();
15dd636f 515 if (ncl)
516 if (ncl >= nclb) {
517 Double_t chi2=fTrackToFollow.GetChi2();
518 if (chi2/ncl < kChi2PerCluster) {
519 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
520 ResetBestTrack();
006b5f7f 521 }
15dd636f 522 }
006b5f7f 523 }
15dd636f 524
006b5f7f 525}
526
006b5f7f 527Int_t AliITStrackerV2::TakeNextProlongation() {
528 //--------------------------------------------------------------------
a9a2d814 529 // This function takes another track prolongation
530 //
531 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
006b5f7f 532 //--------------------------------------------------------------------
18856a77 533 AliITSlayer &layer=fgLayers[fI];
4ab260d6 534 ResetTrackToFollow(fTracks[fI]);
006b5f7f 535
b87bd7cd 536 Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
537 Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
8676d691 538 Double_t road=layer.GetRoad();
539 if (dz*dy>road*road) {
540 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
541 dz=road*scz; dy=road*scy;
542 }
006b5f7f 543
544 const AliITSclusterV2 *c=0; Int_t ci=-1;
c219fd63 545 const AliITSclusterV2 *cc=0; Int_t cci=-1;
546 Double_t chi2=kMaxChi2;
006b5f7f 547 while ((c=layer.GetNextCluster(ci))!=0) {
006b5f7f 548 Int_t idet=c->GetDetectorIndex();
15dd636f 549
4ab260d6 550 if (fTrackToFollow.GetDetectorIndex()!=idet) {
006b5f7f 551 const AliITSdetector &det=layer.GetDetector(idet);
4ab260d6 552 ResetTrackToFollow(fTracks[fI]);
553 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
8676d691 554 //Warning("TakeNextProlongation","propagation failed !\n");
006b5f7f 555 continue;
556 }
4ab260d6 557 fTrackToFollow.SetDetectorIndex(idet);
558 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
006b5f7f 559 }
560
4ab260d6 561 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
562 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
006b5f7f 563
c219fd63 564 Double_t ch2=fTrackToFollow.GetPredictedChi2(c);
565 if (ch2 > chi2) continue;
566 chi2=ch2;
567 cc=c; cci=ci;
568 break;
006b5f7f 569 }
570
c219fd63 571 if (!cc) return 0;
006b5f7f 572
c219fd63 573 if (!fTrackToFollow.Update(cc,chi2,(fI<<28)+cci)) {
8676d691 574 //Warning("TakeNextProlongation","filtering failed !\n");
006b5f7f 575 return 0;
576 }
15dd636f 577
4ab260d6 578 if (fTrackToFollow.GetNumberOfClusters()>1)
579 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
580
a9a2d814 581 fTrackToFollow.
c219fd63 582 SetSampledEdx(cc->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
a9a2d814 583
584 {
880f41b9 585 Double_t x0;
61ab8ea8 586 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
880f41b9 587 fTrackToFollow.CorrectForMaterial(d,x0);
a9a2d814 588 }
4ab260d6 589
a9a2d814 590 if (fConstraint[fPass]) {
4ab260d6 591 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
8676d691 592 Double_t xyz[]={GetX(),GetY(),GetZ()};
593 Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
594 fTrackToFollow.Improve(d,xyz,ers);
a9a2d814 595 }
006b5f7f 596
006b5f7f 597 return 1;
598}
599
600
006b5f7f 601AliITStrackerV2::AliITSlayer::AliITSlayer() {
602 //--------------------------------------------------------------------
603 //default AliITSlayer constructor
604 //--------------------------------------------------------------------
c219fd63 605 fR=0.; fPhiOffset=0.; fZOffset=0.;
606 fNladders=0; fNdetectors=0;
006b5f7f 607 fDetectors=0;
c219fd63 608
609 for (Int_t i=0; i<kNsector; i++) fN[i]=0;
610 fNsel=0;
611
612 fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
006b5f7f 613}
614
615AliITStrackerV2::AliITSlayer::
616AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
617 //--------------------------------------------------------------------
618 //main AliITSlayer constructor
619 //--------------------------------------------------------------------
620 fR=r; fPhiOffset=p; fZOffset=z;
621 fNladders=nl; fNdetectors=nd;
622 fDetectors=new AliITSdetector[fNladders*fNdetectors];
623
c219fd63 624 for (Int_t i=0; i<kNsector; i++) fN[i]=0;
625 fNsel=0;
626
627 for (Int_t i=0; i<kMaxClusterPerLayer; i++) fClusters[i]=0;
b87bd7cd 628
629 fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
006b5f7f 630}
631
632AliITStrackerV2::AliITSlayer::~AliITSlayer() {
633 //--------------------------------------------------------------------
634 // AliITSlayer destructor
635 //--------------------------------------------------------------------
636 delete[] fDetectors;
c219fd63 637 ResetClusters();
006b5f7f 638}
639
61ab8ea8 640void AliITStrackerV2::AliITSlayer::ResetClusters() {
641 //--------------------------------------------------------------------
642 // This function removes loaded clusters
643 //--------------------------------------------------------------------
c219fd63 644 for (Int_t s=0; s<kNsector; s++) {
645 Int_t &n=fN[s];
646 while (n) {
647 n--;
648 delete fClusters[s*kMaxClusterPerSector+n];
649 }
650 }
61ab8ea8 651}
652
b87bd7cd 653void AliITStrackerV2::AliITSlayer::ResetRoad() {
654 //--------------------------------------------------------------------
655 // This function calculates the road defined by the cluster density
656 //--------------------------------------------------------------------
657 Int_t n=0;
c219fd63 658 for (Int_t s=0; s<kNsector; s++) {
659 Int_t i=fN[s];
660 while (i--)
661 if (TMath::Abs(fClusters[s*kMaxClusterPerSector+i]->GetZ())<fR) n++;
b87bd7cd 662 }
663 if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
664}
665
006b5f7f 666Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
667 //--------------------------------------------------------------------
c219fd63 668 // This function inserts a cluster to this layer in increasing
669 // order of the cluster's fZ
006b5f7f 670 //--------------------------------------------------------------------
c219fd63 671 Float_t circ=TMath::TwoPi()*fR;
672 Int_t sec=Int_t(kNsector*c->GetPhiR()/circ);
673 Int_t &n=fN[sec];
674 if (n>=kMaxClusterPerSector) {
675 ::Error("InsertCluster","Too many clusters !\n");
676 return 1;
006b5f7f 677 }
c219fd63 678 if (n==0) fClusters[sec*kMaxClusterPerSector]=c;
679 else {
680 Int_t i=FindClusterIndex(c->GetZ(),sec);
681 Int_t k=n-i+sec*kMaxClusterPerSector;
682 memmove(fClusters+i+1 ,fClusters+i,k*sizeof(AliITSclusterV2*));
683 fClusters[i]=c;
684 }
685 n++;
006b5f7f 686 return 0;
687}
688
c219fd63 689Int_t
690AliITStrackerV2::AliITSlayer::FindClusterIndex(Float_t z,Int_t s) const {
006b5f7f 691 //--------------------------------------------------------------------
c219fd63 692 // For the sector "s", this function returns the index of the first
693 // with its fZ >= "z".
006b5f7f 694 //--------------------------------------------------------------------
c219fd63 695 Int_t nc=fN[s];
696 if (nc==0) return kMaxClusterPerSector*s;
697
698 Int_t b=kMaxClusterPerSector*s;
699 if (z <= fClusters[b]->GetZ()) return b;
700
701 Int_t e=b+nc-1;
702 if (z > fClusters[e]->GetZ()) return e+1;
703
704 Int_t m=(b+e)/2;
006b5f7f 705 for (; b<e; m=(b+e)/2) {
706 if (z > fClusters[m]->GetZ()) b=m+1;
707 else e=m;
708 }
709 return m;
710}
711
c219fd63 712Int_t AliITStrackerV2::AliITSlayer::
713SelectClusters(Float_t zmin,Float_t zmax,Float_t ymin, Float_t ymax) {
006b5f7f 714 //--------------------------------------------------------------------
c219fd63 715 // This function selects clusters within the "window"
006b5f7f 716 //--------------------------------------------------------------------
c219fd63 717 Float_t circ=fR*TMath::TwoPi();
718
719 if (ymin>circ) ymin-=circ; else if (ymin<0) ymin+=circ;
720 if (ymax>circ) ymax-=circ; else if (ymax<0) ymax+=circ;
721
722 Int_t i1=Int_t(kNsector*ymin/circ);
723 if (fN[i1]!=0) {
724 Float_t ym = (ymax<ymin) ? ymax+circ : ymax;
725 Int_t i=FindClusterIndex(zmin,i1), imax=i1*kMaxClusterPerSector+fN[i1];
726 for (; i<imax; i++) {
727 AliITSclusterV2 *c=fClusters[i];
728 if (c->IsUsed()) continue;
729 if (c->GetZ()>zmax) break;
730 if (c->GetPhiR()<=ymin) continue;
731 if (c->GetPhiR()>ym) continue;
732 fIndex[fNsel++]=i;
733 }
734 }
735
736 Int_t i2=Int_t(kNsector*ymax/circ);
737 if (i2==i1) return fNsel;
738
739 if (fN[i2]!=0) {
740 Float_t ym = (ymin>ymax) ? ymin-circ : ymin;
741 Int_t i=FindClusterIndex(zmin,i2), imax=i2*kMaxClusterPerSector+fN[i2];
742 for (; i<imax; i++) {
743 AliITSclusterV2 *c=fClusters[i];
744 if (c->IsUsed()) continue;
745 if (c->GetZ()>zmax) break;
746 if (c->GetPhiR()<=ym) continue;
747 if (c->GetPhiR()>ymax) continue;
748 fIndex[fNsel++]=i;
749 }
750 }
751
752 return fNsel;
006b5f7f 753}
754
755const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
756 //--------------------------------------------------------------------
757 // This function returns clusters within the "window"
758 //--------------------------------------------------------------------
c219fd63 759 AliITSclusterV2 *c=0;
760 ci=-1;
761 if (fNsel) {
762 fNsel--;
763 ci=fIndex[fNsel];
764 c=fClusters[ci];
006b5f7f 765 }
c219fd63 766 return c;
767}
7f6ddf58 768
c219fd63 769Int_t AliITStrackerV2::AliITSlayer::GetNumberOfClusters() const {
770 Int_t n=0;
771 for (Int_t s=0; s<kNsector; s++) n+=fN[s];
772 return n;
006b5f7f 773}
774
c219fd63 775Int_t
776AliITStrackerV2::AliITSlayer::FindDetectorIndex(Double_t phi,Double_t z)const {
006b5f7f 777 //--------------------------------------------------------------------
778 //This function finds the detector crossed by the track
779 //--------------------------------------------------------------------
f363d377 780 Double_t dphi=-(phi-fPhiOffset);
006b5f7f 781 if (dphi < 0) dphi += 2*TMath::Pi();
782 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
783 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
3e3e8427 784 if (np>=fNladders) np-=fNladders;
785 if (np<0) np+=fNladders;
006b5f7f 786
787 Double_t dz=fZOffset-z;
788 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
3e3e8427 789 if (nz>=fNdetectors) return -1;
790 if (nz<0) return -1;
006b5f7f 791
006b5f7f 792 return np*fNdetectors + nz;
793}
794
795Double_t
61ab8ea8 796AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
797const {
a9a2d814 798 //--------------------------------------------------------------------
799 //This function returns the layer thickness at this point (units X0)
800 //--------------------------------------------------------------------
801 Double_t d=0.0085;
61ab8ea8 802 x0=21.82;
a9a2d814 803
804 if (43<fR&&fR<45) { //SSD2
61ab8ea8 805 Double_t dd=0.0034;
806 d=dd;
807 if (TMath::Abs(y-0.00)>3.40) d+=dd;
808 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
809 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
a9a2d814 810 for (Int_t i=0; i<12; i++) {
61ab8ea8 811 if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
812 if (TMath::Abs(y-0.00)>3.40) d+=dd;
813 d+=0.0034;
814 break;
815 }
816 if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
817 if (TMath::Abs(y-0.00)>3.40) d+=dd;
818 d+=0.0034;
819 break;
820 }
821 if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
822 if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
a9a2d814 823 }
824 } else
825 if (37<fR&&fR<41) { //SSD1
61ab8ea8 826 Double_t dd=0.0034;
827 d=dd;
828 if (TMath::Abs(y-0.00)>3.40) d+=dd;
829 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
830 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
a9a2d814 831 for (Int_t i=0; i<11; i++) {
61ab8ea8 832 if (TMath::Abs(z-3.9*i)<0.15) {
833 if (TMath::Abs(y-0.00)>3.40) d+=dd;
834 d+=dd;
835 break;
836 }
837 if (TMath::Abs(z+3.9*i)<0.15) {
838 if (TMath::Abs(y-0.00)>3.40) d+=dd;
839 d+=dd;
840 break;
841 }
842 if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
843 if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
a9a2d814 844 }
845 } else
846 if (13<fR&&fR<26) { //SDD
61ab8ea8 847 Double_t dd=0.0033;
848 d=dd;
849 if (TMath::Abs(y-0.00)>3.30) d+=dd;
850
851 if (TMath::Abs(y-1.80)<0.55) {
852 d+=0.016;
853 for (Int_t j=0; j<20; j++) {
854 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
855 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
856 }
857 }
858 if (TMath::Abs(y+1.80)<0.55) {
859 d+=0.016;
860 for (Int_t j=0; j<20; j++) {
861 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
862 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
863 }
864 }
865
866 for (Int_t i=0; i<4; i++) {
867 if (TMath::Abs(z-7.3*i)<0.60) {
868 d+=dd;
869 if (TMath::Abs(y-0.00)>3.30) d+=dd;
870 break;
871 }
872 if (TMath::Abs(z+7.3*i)<0.60) {
873 d+=dd;
874 if (TMath::Abs(y-0.00)>3.30) d+=dd;
875 break;
876 }
a9a2d814 877 }
878 } else
879 if (6<fR&&fR<8) { //SPD2
61ab8ea8 880 Double_t dd=0.0063; x0=21.5;
881 d=dd;
882 if (TMath::Abs(y-3.08)>0.5) d+=dd;
883 //if (TMath::Abs(y-3.08)>0.45) d+=dd;
884 if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
a9a2d814 885 } else
886 if (3<fR&&fR<5) { //SPD1
61ab8ea8 887 Double_t dd=0.0063; x0=21.5;
888 d=dd;
889 if (TMath::Abs(y+0.21)>0.6) d+=dd;
890 //if (TMath::Abs(y+0.21)>0.45) d+=dd;
891 if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
a9a2d814 892 }
893
a9a2d814 894 return d;
895}
006b5f7f 896
a9a2d814 897Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
006b5f7f 898{
899 //--------------------------------------------------------------------
a9a2d814 900 //Returns the thickness between the current layer and the vertex (units X0)
006b5f7f 901 //--------------------------------------------------------------------
61ab8ea8 902 Double_t d=0.0028*3*3; //beam pipe
903 Double_t x0=0;
006b5f7f 904
18856a77 905 Double_t xn=fgLayers[fI].GetR();
006b5f7f 906 for (Int_t i=0; i<fI; i++) {
18856a77 907 Double_t xi=fgLayers[i].GetR();
908 d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
006b5f7f 909 }
910
911 if (fI>1) {
61ab8ea8 912 Double_t xi=9.;
913 d+=0.0097*xi*xi;
006b5f7f 914 }
915
916 if (fI>3) {
18856a77 917 Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
61ab8ea8 918 d+=0.0034*xi*xi;
006b5f7f 919 }
61ab8ea8 920
006b5f7f 921 return d/(xn*xn);
922}
923
880f41b9 924Bool_t
c630aafd 925AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,const AliITStrackV2 *c) {
4ab260d6 926 //--------------------------------------------------------------------
c630aafd 927 // This function refits the track "t" at the position "x" using
928 // the clusters from "c"
4ab260d6 929 //--------------------------------------------------------------------
880f41b9 930 Int_t index[kMaxLayer];
931 Int_t k;
932 for (k=0; k<kMaxLayer; k++) index[k]=-1;
c630aafd 933 Int_t nc=c->GetNumberOfClusters();
880f41b9 934 for (k=0; k<nc; k++) {
c630aafd 935 Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
880f41b9 936 index[nl]=idx;
937 }
880f41b9 938
4ab260d6 939 Int_t from, to, step;
ee34184f 940 if (xx > t->GetX()) {
4ab260d6 941 from=0; to=kMaxLayer;
942 step=+1;
943 } else {
944 from=kMaxLayer-1; to=-1;
945 step=-1;
946 }
947
948 for (Int_t i=from; i != to; i += step) {
18856a77 949 AliITSlayer &layer=fgLayers[i];
4ab260d6 950 Double_t r=layer.GetR();
951
952 {
953 Double_t hI=i-0.5*step;
880f41b9 954 if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {
18856a77 955 Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
61ab8ea8 956 Double_t d=0.0034, x0=38.6;
880f41b9 957 if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
958 if (!t->PropagateTo(rs,-step*d,x0)) {
4ab260d6 959 return kFALSE;
960 }
961 }
962 }
963
880f41b9 964 // remember old position [SR, GSI 18.02.2003]
965 Double_t oldX=0., oldY=0., oldZ=0.;
966 if (t->IsStartedTimeIntegral() && step==1) {
967 t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
968 }
969 //
970
4ab260d6 971 Double_t x,y,z;
972 if (!t->GetGlobalXYZat(r,x,y,z)) {
973 return kFALSE;
974 }
975 Double_t phi=TMath::ATan2(y,x);
976 Int_t idet=layer.FindDetectorIndex(phi,z);
977 if (idet<0) {
978 return kFALSE;
979 }
980 const AliITSdetector &det=layer.GetDetector(idet);
981 phi=det.GetPhi();
982 if (!t->Propagate(phi,det.GetR())) {
983 return kFALSE;
984 }
985 t->SetDetectorIndex(idet);
986
987 const AliITSclusterV2 *cl=0;
988 Double_t maxchi2=kMaxChi2;
989
990 Int_t idx=index[i];
991 if (idx>0) {
992 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
993 if (idet != c->GetDetectorIndex()) {
994 idet=c->GetDetectorIndex();
995 const AliITSdetector &det=layer.GetDetector(idet);
996 if (!t->Propagate(det.GetPhi(),det.GetR())) {
997 return kFALSE;
998 }
999 t->SetDetectorIndex(idet);
1000 }
1001 Double_t chi2=t->GetPredictedChi2(c);
c630aafd 1002 if (chi2<maxchi2) {
1003 cl=c;
1004 maxchi2=chi2;
1005 } else {
1006 return kFALSE;
1007 }
4ab260d6 1008 }
4ab260d6 1009 /*
1010 if (cl==0)
1011 if (t->GetNumberOfClusters()>2) {
1012 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1013 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1014 Double_t zmin=t->GetZ() - dz;
1015 Double_t zmax=t->GetZ() + dz;
1016 Double_t ymin=t->GetY() + phi*r - dy;
1017 Double_t ymax=t->GetY() + phi*r + dy;
1018 layer.SelectClusters(zmin,zmax,ymin,ymax);
1019
1020 const AliITSclusterV2 *c=0; Int_t ci=-1;
1021 while ((c=layer.GetNextCluster(ci))!=0) {
1022 if (idet != c->GetDetectorIndex()) continue;
1023 Double_t chi2=t->GetPredictedChi2(c);
1024 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1025 }
1026 }
1027 */
4ab260d6 1028 if (cl) {
1029 if (!t->Update(cl,maxchi2,idx)) {
1030 return kFALSE;
1031 }
880f41b9 1032 t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
4ab260d6 1033 }
1034
1035 {
61ab8ea8 1036 Double_t x0;
1037 Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1038 t->CorrectForMaterial(-step*d,x0);
4ab260d6 1039 }
880f41b9 1040
1041 // track time update [SR, GSI 17.02.2003]
1042 if (t->IsStartedTimeIntegral() && step==1) {
1043 Double_t newX, newY, newZ;
1044 t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1045 Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) +
1046 (oldZ-newZ)*(oldZ-newZ);
1047 t->AddTimeStep(TMath::Sqrt(dL2));
1048 }
1049 //
4ab260d6 1050
1051 }
1052
ee34184f 1053 if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
4ab260d6 1054 return kTRUE;
1055}
1056
b87bd7cd 1057void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1058 //--------------------------------------------------------------------
1059 // This function marks clusters assigned to the track
1060 //--------------------------------------------------------------------
1061 AliTracker::UseClusters(t,from);
880f41b9 1062
b87bd7cd 1063 AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
c219fd63 1064 if (c->GetSigmaZ2()>0.1) c->UnUse();
b87bd7cd 1065 c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
c219fd63 1066 if (c->GetSigmaZ2()>0.1) c->UnUse();
880f41b9 1067
b87bd7cd 1068}