]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITStrackerV2.cxx
First commit of ITS tracking version 2 (Yu.Belikov)
[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
18//
19// Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
20//-------------------------------------------------------------------------
21#include <TFile.h>
22#include <TTree.h>
23#include <TRandom.h>
24#include <iostream.h>
25
26#include "AliITSgeom.h"
27#include "AliITSRecPoint.h"
28#include "../TPC/AliTPCtrack.h"
29#include "AliITSclusterV2.h"
30#include "AliITStrackerV2.h"
31
32//#define DEBUG
33
34#ifdef DEBUG
35Int_t LAB=236;
36#endif
37
38extern TRandom *gRandom;
39
40AliITStrackerV2::AliITSlayer AliITStrackerV2::fLayers[kMaxLayer]; // ITS layers
41
42AliITStrackerV2::
43AliITStrackerV2(const AliITSgeom *geom) throw (const Char_t *) {
44 //--------------------------------------------------------------------
45 //This is the AliITStracker constructor
46 //It also reads clusters from a root file
47 //--------------------------------------------------------------------
48 fYV=fZV=0.;
49
50 AliITSgeom *g=(AliITSgeom*)geom;
51
52 Float_t x,y,z;
53 Int_t i;
54 for (i=1; i<kMaxLayer+1; i++) {
55 Int_t nlad=g->GetNladders(i);
56 Int_t ndet=g->GetNdetectors(i);
57
58 g->GetTrans(i,1,1,x,y,z);
59 Double_t r=TMath::Sqrt(x*x + y*y);
60 Double_t poff=TMath::ATan2(y,x);
61 Double_t zoff=z;
62
63 g->GetTrans(i,1,2,x,y,z);
64 r += TMath::Sqrt(x*x + y*y);
65 g->GetTrans(i,2,1,x,y,z);
66 r += TMath::Sqrt(x*x + y*y);
67 g->GetTrans(i,2,2,x,y,z);
68 r += TMath::Sqrt(x*x + y*y);
69 r*=0.25;
70
71 new (fLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
72
73 for (Int_t j=1; j<nlad+1; j++) {
74 for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
75 Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift);
76 Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
77
78 Double_t r =-x*rot[1] + y*rot[0]; if (i==1) r=-r;
79 Double_t phi=TMath::ATan2(rot[1],rot[0]); if (i==1) phi-=3.1415927;
80 phi+=0.5*TMath::Pi();
81 AliITSdetector &det=fLayers[i-1].GetDetector((j-1)*ndet + k-1);
82
83if (phi<0) phi += 2*TMath::Pi();
84
85 new(&det) AliITSdetector(r,phi);
86 }
87 }
88
89 }
90
91 try {
92 //Read clusters
93 TTree *cTree=(TTree*)gDirectory->Get("cTree");
94 if (!cTree) throw
95 ("AliITStrackerV2::AliITStrackerV2 can't get cTree !\n");
96
97 TBranch *branch=cTree->GetBranch("Clusters");
98 if (!branch) throw
99 ("AliITStrackerV2::AliITStrackerV2 can't get Clusters branch !\n");
100
101 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
102 branch->SetAddress(&clusters);
103
104 Int_t nentr=(Int_t)cTree->GetEntries();
105 for (i=0; i<nentr; i++) {
106 if (!cTree->GetEvent(i)) continue;
107 Int_t lay,lad,det; g->GetModuleId(i-1,lay,lad,det);
108 Int_t ncl=clusters->GetEntriesFast();
109 while (ncl--) {
110 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
111
112#ifdef DEBUG
113if (c->GetLabel(2)!=LAB)
114if (c->GetLabel(1)!=LAB)
115if (c->GetLabel(0)!=LAB) continue;
116cout<<lay-1<<' '<<lad-1<<' '<<det-1<<' '<<c->GetY()<<' '<<c->GetZ()<<endl;
117#endif
118
119 fLayers[lay-1].InsertCluster(new AliITSclusterV2(*c));
120 }
121 clusters->Delete();
122 }
123 }
124 catch (const Char_t *msg) {
125 cerr<<msg<<endl;
126 throw;
127 }
128
129 fI=kMaxLayer;
130}
131
132static Int_t lbl;
133
134Int_t AliITStrackerV2::Clusters2Tracks(const TFile *inp, TFile *out) {
135 //--------------------------------------------------------------------
136 //This functions reconstructs ITS tracks
137 //--------------------------------------------------------------------
138 TFile *in=(TFile*)inp;
139 TDirectory *savedir=gDirectory;
140
141 if (!in->IsOpen()) {
142 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
143 cerr<<"file with TPC tracks is not open !\n";
144 return 1;
145 }
146
147 if (!out->IsOpen()) {
148 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
149 cerr<<"file for ITS tracks is not open !\n";
150 return 2;
151 }
152
153 in->cd();
154 TTree *tpcTree=(TTree*)in->Get("TreeT");
155 if (!tpcTree) {
156 cerr<<"AliITStrackerV2::Clusters2Tracks() ";
157 cerr<<"can't get a tree with TPC tracks !\n";
158 return 3;
159 }
160
161 AliTPCtrack *itrack=new AliTPCtrack;
162 tpcTree->SetBranchAddress("tracks",&itrack);
163
164 out->cd();
165 TTree itsTree("TreeT","Tree with ITS tracks");
166 AliITStrackV2 *otrack=&fBestTrack;
167 itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
168
169 Int_t ntrk=0;
170
171 Int_t nentr=(Int_t)tpcTree->GetEntries();
172 for (Int_t i=0; i<nentr; i++) {
173
174 if (!tpcTree->GetEvent(i)) continue;
175 Int_t tpcLabel=itrack->GetLabel(); //save the TPC track label
176lbl=tpcLabel;
177
178#ifdef DEBUG
179if (TMath::Abs(tpcLabel)!=LAB) continue;
180cout<<tpcLabel<<" *****************\n";
181#endif
182
183 try {
184 ResetTrackToFollow(AliITStrackV2(*itrack));
185 } catch (const Char_t *msg) {
186 cerr<<msg<<endl;
187 continue;
188 }
189 ResetBestTrack();
190
191 fYV=gRandom->Gaus(0.,kSigmaYV); fZV=gRandom->Gaus(0.,kSigmaZV);
192
193 Double_t r2=fTrackToFollow.GetX()*fTrackToFollow.GetX();
194 Double_t x0=0.082/21.82*2.33*(45*45+40*40)/r2+2.000/41.28*0.03*75*75/r2;
195 fTrackToFollow.Improve(x0,fYV,fZV);
196
197 //Double_t xk=77.2;
198 Double_t xk=88.;
199 fTrackToFollow.PropagateTo(xk,0.,0.); //Ne if it's still there
200
201 xk-=0.005;
202 fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
203 xk-=0.02;
204 fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar
205 xk-=2.0;
206 fTrackToFollow.PropagateTo(xk, 2.0/41.28*0.029, 2.0*0.029);//Nomex
207 xk-=0.02;
208 fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar
209 xk-=0.005;
210 fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
211
212 xk-=14.5;
213 fTrackToFollow.PropagateTo(xk,0.,0.); //C02
214
215 xk -=0.005;
216 fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7); //Al
217 xk -=0.005;
218 fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
219 xk -=0.02;
220 fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar
221 xk -=0.5;
222 fTrackToFollow.PropagateTo(xk, 0.5/41.28*.029, 0.5*0.029); //Nomex
223 xk -=0.02;
224 fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar
225 xk -=0.005;
226 fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
227 xk -=0.005;
228 fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7); //Al
229
230
231 for (FollowProlongation(); fI<kMaxLayer; fI++) {
232 while (TakeNextProlongation()) FollowProlongation();
233 }
234
235#ifdef DEBUG
236cout<<fBestTrack.GetNumberOfClusters()<<" number of clusters\n\n";
237#endif
238
239 if (fBestTrack.GetNumberOfClusters() >= kMaxLayer-kLayersToSkip) {
240 cerr << ++ntrk << " \r";
241 fBestTrack.SetLabel(tpcLabel);
242 UseClusters(&fBestTrack);
243 itsTree.Fill();
244 }
245
246 }
247
248 itsTree.Write();
249 savedir->cd();
250 cerr<<"Number of TPC tracks: "<<nentr<<endl;
251 cerr<<"Number of prolonged tracks: "<<ntrk<<endl;
252
253 delete itrack;
254
255 return 0;
256}
257
258
259AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
260 //--------------------------------------------------------------------
261 // Return pointer to a given cluster
262 //--------------------------------------------------------------------
263 Int_t l=(index & 0xf0000000) >> 28;
264 Int_t c=(index & 0x0fffffff) >> 00;
265 return fLayers[l].GetCluster(c);
266}
267
268
269void AliITStrackerV2::FollowProlongation() {
270 //--------------------------------------------------------------------
271 //This function finds a track prolongation
272 //--------------------------------------------------------------------
273 Int_t tryAgain=kLayersToSkip;
274
275 while (fI) {
276 fI--;
277
278#ifdef DEBUG
279cout<<fI<<' ';
280#endif
281
282 AliITSlayer &layer=fLayers[fI];
283 AliITStrackV2 &track=fTracks[fI];
284
285 if (fI==3 || fI==1) {
286 Double_t rs=0.5*(fLayers[fI+1].GetR() + layer.GetR());
287 Double_t ds=0.034; if (fI==3) ds=0.039;
288 Double_t dx0r=ds/21.82*2.33, dr=ds*2.33;
289 fTrackToFollow.Propagate(fTrackToFollow.GetAlpha(),rs,1*dx0r,dr);
290 }
291
292 //find intersection
293 Double_t x,y,z;
294 if (!fTrackToFollow.GetGlobalXYZat(layer.GetR(),x,y,z)) {
295 cerr<<"AliITStrackerV2::FollowProlongation: failed to estimate track !\n";
296 break;
297 }
298 Double_t phi=TMath::ATan2(y,x);
299 Double_t d=layer.GetThickness(phi,z);
300 Int_t idet=layer.FindDetectorIndex(phi,z);
301 if (idet<0) {
302 cerr<<"AliITStrackerV2::FollowProlongation: failed to find a detector !\n";
303 break;
304 }
305
306 //propagate to the intersection
307 const AliITSdetector &det=layer.GetDetector(idet);
308 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR(),1*d/21.82*2.33,d*2.33))
309 {
310 cerr<<"AliITStrackerV2::FollowProlongation: propagation failed !\n";
311 break;
312 }
313 fTrackToFollow.SetDetectorIndex(idet);
314
315 //Select possible prolongations and store the current track estimation
316 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
317 Double_t dz=3*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[fI]);
318 if (dz<0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
319 if (dz > kMaxRoad/4) {
320 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n";
321 break;
322 }
323 Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[fI]);
324 if (dy<0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
325 if (dy > kMaxRoad/4) {
326 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n";
327 break;
328 }
329
330 Double_t zmin=track.GetZ() - dz;
331 Double_t zmax=track.GetZ() + dz;
332 Double_t ymin=track.GetY() + layer.GetR()*det.GetPhi() - dy;
333 Double_t ymax=track.GetY() + layer.GetR()*det.GetPhi() + dy;
334 if (ymax>layer.GetR()*2*TMath::Pi()) {
335 ymax-=layer.GetR()*2*TMath::Pi();
336 ymin-=layer.GetR()*2*TMath::Pi();
337 }
338 layer.SelectClusters(zmin,zmax,ymin,ymax);
339
340 //if (1/TMath::Abs(track.Get1Pt())<0.2)
341 //cout<<fI<<' '<<1/TMath::Abs(track.Get1Pt())<<' '
342 // <<dy<<' '<<dz<<' '<<layer.InRoad()<<endl;
343
344 //take another prolongation
345 if (!TakeNextProlongation()) if (!tryAgain--) break;
346 tryAgain=kLayersToSkip;
347
348 }
349
350 //deal with the best track
351 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
352 Int_t nclb=fBestTrack.GetNumberOfClusters();
353 if (ncl)
354 if (ncl >= nclb) {
355 Double_t chi2=fTrackToFollow.GetChi2();
356 if (chi2/ncl < kChi2PerCluster) {
357 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
358 ResetBestTrack();
359 }
360 }
361 }
362
363 if (fI) fI++;
364}
365
366
367Int_t AliITStrackerV2::TakeNextProlongation() {
368 //--------------------------------------------------------------------
369 //This function takes another track prolongation
370 //--------------------------------------------------------------------
371 //Double_t m[20];
372 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
373
374 AliITSlayer &layer=fLayers[fI];
375 AliITStrackV2 &t=fTracks[fI];
376
377 Double_t dz=4*TMath::Sqrt(t.GetSigmaZ2() + kSigmaZ2[fI]);
378 Double_t dy=4*TMath::Sqrt(t.GetSigmaY2() + kSigmaY2[fI]);
379
380 const AliITSclusterV2 *c=0; Int_t ci=-1;
381 Double_t chi2=12345.;
382 while ((c=layer.GetNextCluster(ci))!=0) {
383 //if (fI==0)
384 //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue; //88888888888888888888
385 Int_t idet=c->GetDetectorIndex();
386
387 if (t.GetDetectorIndex()!=idet) {
388 const AliITSdetector &det=layer.GetDetector(idet);
389 if (!t.Propagate(det.GetPhi(),det.GetR(),0.,0.)) {
390 cerr<<"AliITStrackerV2::TakeNextProlongation: propagation failed !\n";
391 continue;
392 }
393 t.SetDetectorIndex(idet);
394
395#ifdef DEBUG
396cout<<fI<<" change detector !\n";
397#endif
398
399 }
400
401 if (TMath::Abs(t.GetZ() - c->GetZ()) > dz) continue;
402 if (TMath::Abs(t.GetY() - c->GetY()) > dy) continue;
403
404 //m[0]=fYV; m[1]=fZV;
405 //chi2=t.GetPredictedChi2(c,m,d/21.82*2.33);
406 chi2=t.GetPredictedChi2(c);
407
408 if (chi2<kMaxChi2) break;
409 }
410
411#ifdef DEBUG
412cout<<fI<<" chi2="<<chi2<<' '<<t.GetY()<<' '<<t.GetZ()<<' '<<dy<<' '<<dz<<endl;
413#endif
414
415 if (chi2>=kMaxChi2) return 0;
416 if (!c) return 0;
417
418 ResetTrackToFollow(t);
419
420 //if (!fTrackToFollow.Update(m,chi2,(fI<<28)+ci)) {
421 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
422 cerr<<"AliITStrackerV2::TakeNextProlongation: filtering failed !\n";
423 return 0;
424 }
425 fTrackToFollow.Improve(d/21.82*2.33,fYV,fZV);
426
427#ifdef DEBUG
428cout<<"accepted lab="<<c->GetLabel(0)<<' '
429 <<fTrackToFollow.GetNumberOfClusters()<<' '
430 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<endl<<endl;
431#endif
432
433 return 1;
434}
435
436
437
438
439AliITStrackerV2::AliITSlayer::AliITSlayer() {
440 //--------------------------------------------------------------------
441 //default AliITSlayer constructor
442 //--------------------------------------------------------------------
443 fN=0;
444 fDetectors=0;
445}
446
447AliITStrackerV2::AliITSlayer::
448AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
449 //--------------------------------------------------------------------
450 //main AliITSlayer constructor
451 //--------------------------------------------------------------------
452 fR=r; fPhiOffset=p; fZOffset=z;
453 fNladders=nl; fNdetectors=nd;
454 fDetectors=new AliITSdetector[fNladders*fNdetectors];
455
456 fN=0;
457 fI=0;
458}
459
460AliITStrackerV2::AliITSlayer::~AliITSlayer() {
461 //--------------------------------------------------------------------
462 // AliITSlayer destructor
463 //--------------------------------------------------------------------
464 delete[] fDetectors;
465 for (Int_t i=0; i<fN; i++) delete fClusters[i];
466}
467
468Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
469 //--------------------------------------------------------------------
470 //This function adds a cluster to this layer
471 //--------------------------------------------------------------------
472 if (fN==kMaxClusterPerLayer) {
473 cerr<<"AliITStrackerV2::AliITSlayer::InsertCluster(): Too many clusters !\n";
474 return 1;
475 }
476
477 if (fN==0) {fClusters[fN++]=c; return 0;}
478 Int_t i=FindClusterIndex(c->GetZ());
479 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
480 fClusters[i]=c; fN++;
481
482 return 0;
483}
484
485Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
486 //--------------------------------------------------------------------
487 // This function returns the index of the nearest cluster
488 //--------------------------------------------------------------------
489 if (fN==0) return 0;
490 if (z <= fClusters[0]->GetZ()) return 0;
491 if (z > fClusters[fN-1]->GetZ()) return fN;
492 Int_t b=0, e=fN-1, m=(b+e)/2;
493 for (; b<e; m=(b+e)/2) {
494 if (z > fClusters[m]->GetZ()) b=m+1;
495 else e=m;
496 }
497 return m;
498}
499
500void AliITStrackerV2::AliITSlayer::
501SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
502 //--------------------------------------------------------------------
503 // This function sets the "window"
504 //--------------------------------------------------------------------
505 fI=FindClusterIndex(zmin); fZmax=zmax;
506 fYmin=ymin; fYmax=ymax;
507}
508
509const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
510 //--------------------------------------------------------------------
511 // This function returns clusters within the "window"
512 //--------------------------------------------------------------------
513 const AliITSclusterV2 *cluster=0;
514 for (Int_t i=fI; i<fN; i++) {
515 const AliITSclusterV2 *c=fClusters[i];
516 if (c->GetZ() > fZmax) break;
517 if (c->IsUsed()) continue;
518 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
519 Double_t y=fR*det.GetPhi() + c->GetY();
520
521 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
522 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
523
524 if (y<fYmin) continue;
525 if (y>fYmax) continue;
526 cluster=c; ci=i;
527 fI=i+1;
528 break;
529 }
530 return cluster;
531}
532
533Int_t AliITStrackerV2::AliITSlayer::
534FindDetectorIndex(Double_t phi, Double_t z) const {
535 //--------------------------------------------------------------------
536 //This function finds the detector crossed by the track
537 //--------------------------------------------------------------------
538 Double_t dphi=phi-fPhiOffset;
539 if (dphi < 0) dphi += 2*TMath::Pi();
540 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
541 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
542
543 Double_t dz=fZOffset-z;
544 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
545
546 if (np>=fNladders) np-=fNladders;
547 if (np<0) np+=fNladders;
548
549#ifdef DEBUG
550cout<<np<<' '<<nz<<endl;
551#endif
552
553 return np*fNdetectors + nz;
554}
555
556Double_t
557AliITStrackerV2::AliITSlayer::GetThickness(Double_t phi, Double_t z) const {
558 //--------------------------------------------------------------------
559 //This function returns the thickness of this layer
560 //--------------------------------------------------------------------
561 //-pi<phi<+pi
562 if (3 <fR&&fR<8 ) return 1.1*0.096;
563 if (13<fR&&fR<26) return 1.1*0.088;
564 if (37<fR&&fR<41) return 1.1*0.085;
565 return 1.1*0.081;
566}
567
568
569Double_t AliITStrackerV2::GetEffectiveThickness(Double_t phi,Double_t z) const
570{
571 //--------------------------------------------------------------------
572 //Returns the thickness between the current layer and the vertex
573 //--------------------------------------------------------------------
574 //-pi<phi<+pi
575 Double_t d=0.1;
576
577 Double_t xn=fLayers[fI].GetR();
578 for (Int_t i=0; i<fI; i++) {
579 Double_t xi=fLayers[i].GetR();
580 d+=fLayers[i].GetThickness(phi,z)*xi*xi;
581 }
582
583 if (fI>1) {
584 Double_t xi=0.5*(fLayers[1].GetR()+fLayers[2].GetR());
585 d+=0.034*xi*xi;
586 }
587
588 if (fI>3) {
589 Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
590 d+=0.039*xi*xi;
591 }
592 return d/(xn*xn);
593}
594
595
596
597Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
598 //--------------------------------------------------------------------
599 // This function returns number of clusters within the "window"
600 //--------------------------------------------------------------------
601 Int_t ncl=0;
602 for (Int_t i=fI; i<fN; i++) {
603 const AliITSclusterV2 *c=fClusters[i];
604 if (c->GetZ() > fZmax) break;
605 //if (c->IsUsed()) continue;
606 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
607 Double_t y=fR*det.GetPhi() + c->GetY();
608
609 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
610 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
611
612 if (y<fYmin) continue;
613 if (y>fYmax) continue;
614 ncl++;
615 }
616 return ncl;
617}
618
619
620