]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ITS/AliITStrackerV2.cxx
Uncertainty of the tracking step included (Yu.Belikov)
[u/mrichter/AliRoot.git] / ITS / AliITStrackerV2.cxx
... / ...
CommitLineData
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// dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
21//-------------------------------------------------------------------------
22#include <TFile.h>
23#include <TTree.h>
24#include <TRandom.h>
25#include <Riostream.h>
26
27#include "AliITSgeom.h"
28#include "AliITSRecPoint.h"
29#include "AliTPCtrack.h"
30#include "AliITSclusterV2.h"
31#include "AliITStrackerV2.h"
32
33//#define DEBUG
34
35#ifdef DEBUG
36Int_t LAB=70201;
37#endif
38
39AliITStrackerV2::AliITSlayer AliITStrackerV2::fLayers[kMaxLayer]; // ITS layers
40
41AliITStrackerV2::
42AliITStrackerV2(const AliITSgeom *geom, Int_t eventn) throw (const Char_t *) :
43AliTracker() {
44 //--------------------------------------------------------------------
45 //This is the AliITStracker constructor
46 //It also reads clusters from a root file
47 //--------------------------------------------------------------------
48 fEventN=eventn; //MI change add event number - used to generate identifier
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 //MI change
94 char cname[100];
95 sprintf(cname,"TreeC_ITS_%d",eventn);
96 TTree *cTree=(TTree*)gDirectory->Get(cname);
97
98 if (!cTree) throw
99 ("AliITStrackerV2::AliITStrackerV2 can't get cTree !\n");
100
101 TBranch *branch=cTree->GetBranch("Clusters");
102 if (!branch) throw
103 ("AliITStrackerV2::AliITStrackerV2 can't get Clusters branch !\n");
104
105 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
106 branch->SetAddress(&clusters);
107
108 Int_t nentr=(Int_t)cTree->GetEntries();
109 for (i=0; i<nentr; i++) {
110 if (!cTree->GetEvent(i)) continue;
111 Int_t lay,lad,det; g->GetModuleId(i,lay,lad,det);
112 Int_t ncl=clusters->GetEntriesFast();
113 while (ncl--) {
114 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
115
116#ifdef DEBUG
117if (c->GetLabel(2)!=LAB)
118if (c->GetLabel(1)!=LAB)
119if (c->GetLabel(0)!=LAB) continue;
120cout<<lay-1<<' '<<lad-1<<' '<<det-1<<' '<<c->GetY()<<' '<<c->GetZ()<<endl;
121#endif
122
123 fLayers[lay-1].InsertCluster(new AliITSclusterV2(*c));
124 }
125 clusters->Delete();
126 }
127 delete cTree; //Thanks to Mariana Bondila
128 }
129 catch (const Char_t *msg) {
130 cerr<<msg<<endl;
131 throw;
132 }
133
134 fI=kMaxLayer;
135
136 fPass=0;
137 fConstraint[0]=1; fConstraint[1]=0;
138}
139
140//#ifdef DEBUG
141static Int_t lbl;
142//#endif
143
144Int_t AliITStrackerV2::Clusters2Tracks(const TFile *inp, TFile *out) {
145 //--------------------------------------------------------------------
146 //This functions reconstructs ITS tracks
147 //--------------------------------------------------------------------
148 TFile *in=(TFile*)inp;
149 TDirectory *savedir=gDirectory;
150
151 if (!in->IsOpen()) {
152 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
153 cerr<<"file with TPC tracks is not open !\n";
154 return 1;
155 }
156
157 if (!out->IsOpen()) {
158 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
159 cerr<<"file for ITS tracks is not open !\n";
160 return 2;
161 }
162
163 in->cd();
164
165 char tname[100];
166 Int_t nentr=0; TObjArray itsTracks(15000);
167
168 {/* Read TPC tracks */
169 sprintf(tname,"TreeT_TPC_%d",fEventN);
170 TTree *tpcTree=(TTree*)in->Get(tname);
171 if (!tpcTree) {
172 cerr<<"AliITStrackerV2::Clusters2Tracks(): "
173 "can't get a tree with TPC tracks !\n";
174 return 3;
175 }
176 AliTPCtrack *itrack=new AliTPCtrack;
177 tpcTree->SetBranchAddress("tracks",&itrack);
178 nentr=(Int_t)tpcTree->GetEntries();
179 for (Int_t i=0; i<nentr; i++) {
180 tpcTree->GetEvent(i);
181 AliITStrackV2 *t=0;
182 try {
183 t=new AliITStrackV2(*itrack);
184 } catch (const Char_t *msg) {
185 cerr<<msg<<endl;
186 delete t;
187 continue;
188 }
189 if (TMath::Abs(t->GetD())>4) continue;
190
191 t->PropagateTo(80.,0.0053);
192 if (TMath::Abs(t->GetY())>13.) t->CorrectForMaterial(0.03);
193 if (TMath::Abs(t->GetZ())<0.2) t->CorrectForMaterial(0.40);
194 t->PropagateTo(61.,0.0052);
195 Double_t xk=52.,x,y,z; t->GetGlobalXYZat(xk,x,y,z);
196 if (TMath::Abs(y)<7.77) t->PropagateTo(xk,0.19,24.);
197 t->PropagateTo(50.,0.001);
198
199 itsTracks.AddLast(t);
200 }
201 delete tpcTree; //Thanks to Mariana Bondila
202 delete itrack;
203 }
204 itsTracks.Sort();
205
206 out->cd();
207
208 sprintf(tname,"TreeT_ITS_%d",fEventN);
209 TTree itsTree(tname,"Tree with ITS tracks");
210 AliITStrackV2 *otrack=&fBestTrack;
211 itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
212
213 for (fPass=0; fPass<2; fPass++) {
214 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
215 for (Int_t i=0; i<nentr; i++) {
216 if (i%10==0) cerr<<nentr-i<<" \r";
217 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
218 if (t==0) continue; //this track has been already tracked
219 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
220
221lbl=tpcLabel;
222#ifdef DEBUG
223lbl=tpcLabel;
224if (TMath::Abs(tpcLabel)!=LAB) continue;
225cout<<tpcLabel<<" *****************\n";
226#endif
227
228 ResetTrackToFollow(*t);
229 ResetBestTrack();
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) continue;
240
241 if (fConstraint[fPass]) {
242 Int_t index[kMaxLayer];
243 Int_t k;
244 for (k=0; k<kMaxLayer; k++) index[k]=-1;
245 Int_t nc=fBestTrack.GetNumberOfClusters();
246 for (k=0; k<nc; k++) {
247 Int_t idx=fBestTrack.GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
248 index[nl]=idx;
249 }
250 fBestTrack.~AliITStrackV2(); new(&fBestTrack) AliITStrackV2(*t);
251 if (!RefitAt(3.7, &fBestTrack, index)) continue;
252 }
253
254 fBestTrack.SetLabel(tpcLabel);
255 fBestTrack.CookdEdx();
256 CookLabel(&fBestTrack,0.); //For comparison only
257 itsTree.Fill();
258 UseClusters(&fBestTrack);
259 delete itsTracks.RemoveAt(i);
260
261 }
262 }
263
264 itsTree.Write();
265
266 itsTracks.Delete();
267 savedir->cd();
268 cerr<<"Number of TPC tracks: "<<nentr<<endl;
269 cerr<<"Number of prolonged tracks: "<<itsTree.GetEntries()<<endl;
270
271 return 0;
272}
273
274Int_t AliITStrackerV2::PropagateBack(const TFile *inp, TFile *out) {
275 //--------------------------------------------------------------------
276 //This functions propagates reconstructed ITS tracks back
277 //--------------------------------------------------------------------
278 TFile *in=(TFile*)inp;
279 TDirectory *savedir=gDirectory;
280
281 if (!in->IsOpen()) {
282 cerr<<"AliITStrackerV2::PropagateBack(): ";
283 cerr<<"file with ITS tracks is not open !\n";
284 return 1;
285 }
286
287 if (!out->IsOpen()) {
288 cerr<<"AliITStrackerV2::PropagateBack(): ";
289 cerr<<"file for back propagated ITS tracks is not open !\n";
290 return 2;
291 }
292
293 in->cd();
294 TTree *itsTree=(TTree*)in->Get("TreeT_ITS_0");
295 if (!itsTree) {
296 cerr<<"AliITStrackerV2::PropagateBack() ";
297 cerr<<"can't get a tree with ITS tracks !\n";
298 return 3;
299 }
300 AliITStrackV2 *itrack=new AliITStrackV2;
301 itsTree->SetBranchAddress("tracks",&itrack);
302
303 out->cd();
304 TTree backTree("TreeT_ITSb_0","Tree with back propagated ITS tracks");
305 AliTPCtrack *otrack=0;
306 backTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
307
308 Int_t ntrk=0;
309
310 Int_t nentr=(Int_t)itsTree->GetEntries();
311 for (Int_t i=0; i<nentr; i++) {
312 itsTree->GetEvent(i);
313 ResetTrackToFollow(*itrack);
314 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
315 Int_t itsLabel=fTrackToFollow.GetLabel(); //save the ITS track label
316
317#ifdef DEBUG
318if (TMath::Abs(itsLabel)!=LAB) continue;
319cout<<itsLabel<<" *****************\n";
320#endif
321
322 try {
323 Int_t nc=itrack->GetNumberOfClusters();
324#ifdef DEBUG
325for (Int_t k=0; k<nc; k++) {
326 Int_t index=itrack->GetClusterIndex(k);
327 AliITSclusterV2 *c=(AliITSclusterV2*)GetCluster(index);
328 cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl;
329}
330#endif
331 Int_t idx=0, l=0;
332 const AliITSclusterV2 *c=0;
333 if (nc--) {
334 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
335 c=(AliITSclusterV2*)GetCluster(idx);
336 }
337 for (fI=0; fI<kMaxLayer; fI++) {
338 AliITSlayer &layer=fLayers[fI];
339 Double_t r=layer.GetR();
340 if (fI==2 || fI==4) {
341 Double_t rs=0.5*(fLayers[fI-1].GetR() + r);
342 Double_t d=0.011; if (fI==4) d=0.0053;
343 if (!fTrackToFollow.PropagateTo(rs,-d)) throw "";
344 }
345
346 Double_t x,y,z;
347 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z))
348 throw "AliITStrackerV2::PropagateBack: failed to estimate track !";
349 Double_t phi=TMath::ATan2(y,x);
350 Int_t idet=layer.FindDetectorIndex(phi,z);
351 if (idet<0)
352 throw "AliITStrackerV2::PropagateBack: failed to find a detector !\n";
353 const AliITSdetector &det=layer.GetDetector(idet);
354 r=det.GetR(); phi=det.GetPhi();
355 if (!fTrackToFollow.Propagate(phi,r)) throw "";
356 fTrackToFollow.SetDetectorIndex(idet);
357
358 const AliITSclusterV2 *cl=0;
359 Int_t index=0;
360 Double_t maxchi2=kMaxChi2;
361
362 if (l==fI) {
363 idet=c->GetDetectorIndex();
364 if (idet != fTrackToFollow.GetDetectorIndex()) {
365 const AliITSdetector &det=layer.GetDetector(idet);
366 r=det.GetR(); phi=det.GetPhi();
367 if (!fTrackToFollow.Propagate(phi,r)) throw "";
368 fTrackToFollow.SetDetectorIndex(idet);
369 }
370 Double_t chi2=fTrackToFollow.GetPredictedChi2(c);
371 if (chi2<kMaxChi2) {
372 cl=c; maxchi2=chi2; index=idx;
373 }
374 if (nc--) {
375 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
376 c=(AliITSclusterV2*)GetCluster(idx);
377 }
378 }
379
380 if (fTrackToFollow.GetNumberOfClusters()>2) {
381 Double_t dz=3*TMath::Sqrt(fTrackToFollow.GetSigmaZ2()+kSigmaZ2[fI]);
382 Double_t dy=3*TMath::Sqrt(fTrackToFollow.GetSigmaY2()+kSigmaY2[fI]);
383 Double_t zmin=fTrackToFollow.GetZ() - dz;
384 Double_t zmax=fTrackToFollow.GetZ() + dz;
385 Double_t ymin=fTrackToFollow.GetY() + phi*r - dy;
386 Double_t ymax=fTrackToFollow.GetY() + phi*r + dy;
387 layer.SelectClusters(zmin,zmax,ymin,ymax);
388
389 const AliITSclusterV2 *cc=0; Int_t ci;
390 while ((cc=layer.GetNextCluster(ci))!=0) {
391 idet=cc->GetDetectorIndex();
392 if (idet != fTrackToFollow.GetDetectorIndex()) continue;
393 Double_t chi2=fTrackToFollow.GetPredictedChi2(cc);
394 if (chi2<maxchi2) {
395 cl=cc; index=(fI<<28)+ci; maxchi2=chi2;
396 }
397 }
398 }
399
400 if (cl) {
401 if (!fTrackToFollow.Update(cl,maxchi2,index))
402 cerr<<"AliITStrackerV2::PropagateBack: filtering failed !\n";
403 }
404
405 x=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ());
406 fTrackToFollow.CorrectForMaterial(-x);
407
408 }
409
410 Double_t xk=52.,x,y,z; fTrackToFollow.GetGlobalXYZat(xk,x,y,z);
411 if (TMath::Abs(y)<7.77)
412 fTrackToFollow.PropagateTo(xk,-0.19,24.);
413 fTrackToFollow.PropagateTo(61,-0.0110);
414 fTrackToFollow.PropagateTo(80.,-0.0053);
415
416 fTrackToFollow.SetLabel(itsLabel);
417 otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha());
418 backTree.Fill(); delete otrack;
419 UseClusters(&fTrackToFollow);
420 cerr << ++ntrk << " \r";
421 }
422 catch (const Char_t *msg) {
423 cerr<<msg<<endl;
424 }
425 }
426
427 backTree.Write();
428 savedir->cd();
429 cerr<<"Number of ITS tracks: "<<nentr<<endl;
430 cerr<<"Number of back propagated ITS tracks: "<<ntrk<<endl;
431
432 delete itrack;
433
434 delete itsTree; //Thanks to Mariana Bondila
435
436 return 0;
437}
438
439
440AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
441 //--------------------------------------------------------------------
442 // Return pointer to a given cluster
443 //--------------------------------------------------------------------
444 Int_t l=(index & 0xf0000000) >> 28;
445 Int_t c=(index & 0x0fffffff) >> 00;
446 return fLayers[l].GetCluster(c);
447}
448
449
450void AliITStrackerV2::FollowProlongation() {
451 //--------------------------------------------------------------------
452 //This function finds a track prolongation
453 //--------------------------------------------------------------------
454 Int_t tryAgain=kLayersToSkip;
455
456 while (fI) {
457 Int_t i=fI-1;
458#ifdef DEBUG
459cout<<i<<' ';
460#endif
461 AliITSlayer &layer=fLayers[i];
462 AliITStrackV2 &track=fTracks[i];
463
464 Double_t r=layer.GetR();
465 if (i==3 || i==1) {
466 Double_t rs=0.5*(fLayers[i+1].GetR() + r);
467 Double_t d=0.011; if (i==3) d=0.0053;
468 if (!fTrackToFollow.PropagateTo(rs,d)) {
469 //cerr<<"AliITStrackerV2::FollowProlongation: "
470 //"propagation failed !\n";
471 break;
472 }
473 }
474
475 //find intersection
476 Double_t x,y,z;
477 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
478 //cerr<<"AliITStrackerV2::FollowProlongation: "
479 //"failed to estimate track !\n";
480 break;
481 }
482 Double_t phi=TMath::ATan2(y,x);
483 Int_t idet=layer.FindDetectorIndex(phi,z);
484 if (idet<0) {
485 //cerr<<"AliITStrackerV2::FollowProlongation: "
486 //"failed to find a detector !\n";
487 break;
488 }
489
490 //propagate to the intersection
491 const AliITSdetector &det=layer.GetDetector(idet);
492 phi=det.GetPhi();
493 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
494 //cerr<<"AliITStrackerV2::FollowProlongation: propagation failed !\n";
495 break;
496 }
497 fTrackToFollow.SetDetectorIndex(idet);
498
499 //Select possible prolongations and store the current track estimation
500 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
501 Double_t dz=3*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
502 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
503 if (dz > kMaxRoad) {
504 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n";
505 break;
506 }
507
508 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) break;
509
510 Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
511 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
512 if (dy > kMaxRoad) {
513 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n";
514 break;
515 }
516
517 Double_t zmin=track.GetZ() - dz;
518 Double_t zmax=track.GetZ() + dz;
519 Double_t ymin=track.GetY() + r*phi - dy;
520 Double_t ymax=track.GetY() + r*phi + dy;
521 layer.SelectClusters(zmin,zmax,ymin,ymax);
522 fI--;
523
524 //take another prolongation
525 if (!TakeNextProlongation()) if (!tryAgain--) break;
526 tryAgain=kLayersToSkip;
527
528 }
529
530 //deal with the best track
531 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
532 Int_t nclb=fBestTrack.GetNumberOfClusters();
533 if (ncl)
534 if (ncl >= nclb) {
535 Double_t chi2=fTrackToFollow.GetChi2();
536 if (chi2/ncl < kChi2PerCluster) {
537 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
538 ResetBestTrack();
539 }
540 }
541 }
542
543}
544
545Int_t AliITStrackerV2::TakeNextProlongation() {
546 //--------------------------------------------------------------------
547 // This function takes another track prolongation
548 //
549 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
550 //--------------------------------------------------------------------
551 AliITSlayer &layer=fLayers[fI];
552 ResetTrackToFollow(fTracks[fI]);
553
554 Double_t dz=4*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
555 Double_t dy=4*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
556
557 const AliITSclusterV2 *c=0; Int_t ci=-1;
558 Double_t chi2=12345.;
559 while ((c=layer.GetNextCluster(ci))!=0) {
560 //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue;
561 Int_t idet=c->GetDetectorIndex();
562
563 if (fTrackToFollow.GetDetectorIndex()!=idet) {
564 const AliITSdetector &det=layer.GetDetector(idet);
565 ResetTrackToFollow(fTracks[fI]);
566 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
567 //cerr<<"AliITStrackerV2::TakeNextProlongation: "
568 //"propagation failed !\n";
569 continue;
570 }
571 fTrackToFollow.SetDetectorIndex(idet);
572 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
573
574#ifdef DEBUG
575cout<<fI<<" change detector !\n";
576#endif
577
578 }
579
580 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
581 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
582
583 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
584 }
585
586#ifdef DEBUG
587cout<<fI<<" chi2="<<chi2<<' '
588 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
589 <<dy<<' '<<dz<<endl;
590#endif
591
592 if (chi2>=kMaxChi2) return 0;
593 if (!c) return 0;
594
595 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
596 //cerr<<"AliITStrackerV2::TakeNextProlongation: filtering failed !\n";
597 return 0;
598 }
599
600 if (fTrackToFollow.GetNumberOfClusters()>1)
601 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
602
603 fTrackToFollow.
604 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
605
606 {
607 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ());
608 fTrackToFollow.CorrectForMaterial(d);
609 }
610
611 if (fConstraint[fPass]) {
612 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
613 fTrackToFollow.Improve(d,GetY(),GetZ());
614 }
615
616#ifdef DEBUG
617cout<<"accepted lab="<<c->GetLabel(0)<<' '
618 <<fTrackToFollow.GetNumberOfClusters()<<' '
619 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
620 <<fTrackToFollow.Get1Pt()<<endl<<endl;
621#endif
622
623 return 1;
624}
625
626
627AliITStrackerV2::AliITSlayer::AliITSlayer() {
628 //--------------------------------------------------------------------
629 //default AliITSlayer constructor
630 //--------------------------------------------------------------------
631 fN=0;
632 fDetectors=0;
633}
634
635AliITStrackerV2::AliITSlayer::
636AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
637 //--------------------------------------------------------------------
638 //main AliITSlayer constructor
639 //--------------------------------------------------------------------
640 fR=r; fPhiOffset=p; fZOffset=z;
641 fNladders=nl; fNdetectors=nd;
642 fDetectors=new AliITSdetector[fNladders*fNdetectors];
643
644 fN=0;
645 fI=0;
646}
647
648AliITStrackerV2::AliITSlayer::~AliITSlayer() {
649 //--------------------------------------------------------------------
650 // AliITSlayer destructor
651 //--------------------------------------------------------------------
652 delete[] fDetectors;
653 for (Int_t i=0; i<fN; i++) delete fClusters[i];
654}
655
656Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
657 //--------------------------------------------------------------------
658 //This function adds a cluster to this layer
659 //--------------------------------------------------------------------
660 if (fN==kMaxClusterPerLayer) {
661 cerr<<"AliITStrackerV2::AliITSlayer::InsertCluster(): "
662 "Too many clusters !\n";
663 return 1;
664 }
665
666 if (fN==0) {fClusters[fN++]=c; return 0;}
667 Int_t i=FindClusterIndex(c->GetZ());
668 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
669 fClusters[i]=c; fN++;
670
671 return 0;
672}
673
674Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
675 //--------------------------------------------------------------------
676 // This function returns the index of the nearest cluster
677 //--------------------------------------------------------------------
678 if (fN==0) return 0;
679 if (z <= fClusters[0]->GetZ()) return 0;
680 if (z > fClusters[fN-1]->GetZ()) return fN;
681 Int_t b=0, e=fN-1, m=(b+e)/2;
682 for (; b<e; m=(b+e)/2) {
683 if (z > fClusters[m]->GetZ()) b=m+1;
684 else e=m;
685 }
686 return m;
687}
688
689void AliITStrackerV2::AliITSlayer::
690SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
691 //--------------------------------------------------------------------
692 // This function sets the "window"
693 //--------------------------------------------------------------------
694 fI=FindClusterIndex(zmin); fZmax=zmax;
695 Double_t circle=2*TMath::Pi()*fR;
696 if (ymax>circle) { ymax-=circle; ymin-=circle; }
697 fYmin=ymin; fYmax=ymax;
698}
699
700const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
701 //--------------------------------------------------------------------
702 // This function returns clusters within the "window"
703 //--------------------------------------------------------------------
704 const AliITSclusterV2 *cluster=0;
705 for (Int_t i=fI; i<fN; i++) {
706 const AliITSclusterV2 *c=fClusters[i];
707 if (c->GetZ() > fZmax) break;
708 if (c->IsUsed()) continue;
709 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
710 Double_t y=fR*det.GetPhi() + c->GetY();
711
712 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
713 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
714
715 if (y<fYmin) continue;
716 if (y>fYmax) continue;
717 cluster=c; ci=i;
718 fI=i+1;
719 break;
720 }
721
722 return cluster;
723}
724
725Int_t AliITStrackerV2::AliITSlayer::
726FindDetectorIndex(Double_t phi, Double_t z) const {
727 //--------------------------------------------------------------------
728 //This function finds the detector crossed by the track
729 //--------------------------------------------------------------------
730 Double_t dphi=phi-fPhiOffset;
731 if (dphi < 0) dphi += 2*TMath::Pi();
732 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
733 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
734 if (np>=fNladders) np-=fNladders;
735 if (np<0) np+=fNladders;
736
737 Double_t dz=fZOffset-z;
738 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
739 if (nz>=fNdetectors) return -1;
740 if (nz<0) return -1;
741
742#ifdef DEBUG
743cout<<np<<' '<<nz<<endl;
744#endif
745
746 return np*fNdetectors + nz;
747}
748
749Double_t
750AliITStrackerV2::AliITSlayer::GetThickness(Double_t y, Double_t z) const {
751 //--------------------------------------------------------------------
752 //This function returns the layer thickness at this point (units X0)
753 //--------------------------------------------------------------------
754 Double_t d=0.0085;
755
756 if (43<fR&&fR<45) { //SSD2
757 d=0.0036;
758 if (TMath::Abs(y-0.00)>3.40) d+=0.0036;
759 if (TMath::Abs(y-2.50)<0.10) d+=(0.02-0.0036);
760 if (TMath::Abs(y+1.90)<0.10) d+=(0.02-0.0036);
761 for (Int_t i=0; i<12; i++) {
762 if (TMath::Abs(z-3.6*(i+0.5))<0.20) {d+=0.0036; break;}
763 if (TMath::Abs(z+3.6*(i+0.5))<0.20) {d+=0.0036; break;}
764 if (TMath::Abs(z-3.6*(i+0.929))<0.50) {d+=(0.02-0.0036); break;}
765 if (TMath::Abs(z+3.6*(i+0.104))<0.50) {d+=(0.02-0.0036); break;}
766 }
767 } else
768 if (37<fR&&fR<41) { //SSD1
769 d=0.0036;
770 if (TMath::Abs(y-0.00)>3.40) d+=0.0036;
771 if (TMath::Abs(y-2.20)<0.10) d+=(0.02-0.0036);
772 if (TMath::Abs(y+2.20)<0.10) d+=(0.02-0.0036);
773 for (Int_t i=0; i<11; i++) {
774 if (TMath::Abs(z-3.6*i)<0.20) {d+=0.0036; break;}
775 if (TMath::Abs(z+3.6*i)<0.20) {d+=0.0036; break;}
776 if (TMath::Abs(z-3.6*(i+0.54))<0.50) {d+=(0.02-0.0036); break;}
777 if (TMath::Abs(z+3.6*(i+0.58))<0.50) {d+=(0.02-0.0036); break;}
778 }
779 } else
780 if (13<fR&&fR<26) { //SDD
781 d=0.0034;
782 if (TMath::Abs(y-0.00)>3.30) d+=0.0034;
783 if (TMath::Abs(y-2.10)<0.20) d+=0.0034*3;
784 if (TMath::Abs(y+2.10)<0.20) d+=0.0034*3;
785 for (Int_t i=0; i<4; i++) {
786 if (TMath::Abs(z-7.3*i)<0.60) {d+=0.0034; break;}
787 if (TMath::Abs(z+7.3*i)<0.60) {d+=0.0034; break;}
788 }
789 } else
790 if (6<fR&&fR<8) { //SPD2
791 d=0.0093;
792 if (TMath::Abs(y-3.08)>0.45) d+=0.0064;
793 if (TMath::Abs(y-3.03)<0.10) d+=0.0192;
794 } else
795 if (3<fR&&fR<5) { //SPD1
796 d=0.0093;
797 if (TMath::Abs(y+0.21)>0.55) d+=0.0064;
798 if (TMath::Abs(y+0.10)<0.10) d+=0.0192;
799 }
800
801 d+=0.002;
802
803 return d;
804}
805
806Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
807{
808 //--------------------------------------------------------------------
809 //Returns the thickness between the current layer and the vertex (units X0)
810 //--------------------------------------------------------------------
811 Double_t d=0.1/65.19*1.848;
812
813 Double_t xn=fLayers[fI].GetR();
814 for (Int_t i=0; i<fI; i++) {
815 Double_t xi=fLayers[i].GetR();
816 d+=fLayers[i].GetThickness(y,z)*xi*xi;
817 }
818
819 if (fI>1) {
820 Double_t xi=0.5*(fLayers[1].GetR()+fLayers[2].GetR());
821 d+=0.011*xi*xi;
822 }
823
824 if (fI>3) {
825 Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
826 d+=0.0053*xi*xi;
827 }
828 return d/(xn*xn);
829}
830
831
832
833Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
834 //--------------------------------------------------------------------
835 // This function returns number of clusters within the "window"
836 //--------------------------------------------------------------------
837 Int_t ncl=0;
838 for (Int_t i=fI; i<fN; i++) {
839 const AliITSclusterV2 *c=fClusters[i];
840 if (c->GetZ() > fZmax) break;
841 if (c->IsUsed()) continue;
842 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
843 Double_t y=fR*det.GetPhi() + c->GetY();
844
845 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
846 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
847
848 if (y<fYmin) continue;
849 if (y>fYmax) continue;
850 ncl++;
851 }
852 return ncl;
853}
854
855Bool_t AliITStrackerV2::RefitAt(Double_t x, AliITStrackV2 *t, Int_t *index) {
856 //--------------------------------------------------------------------
857 // This function refits a track at a given position
858 //--------------------------------------------------------------------
859 Int_t from, to, step;
860 if (x > t->GetX()) {
861 from=0; to=kMaxLayer;
862 step=+1;
863 } else {
864 from=kMaxLayer-1; to=-1;
865 step=-1;
866 }
867
868 for (Int_t i=from; i != to; i += step) {
869 AliITSlayer &layer=fLayers[i];
870 Double_t r=layer.GetR();
871
872 {
873 Double_t hI=i-0.5*step;
874 if (hI==1.5 || hI==3.5) {
875 Double_t rs=0.5*(fLayers[i-step].GetR() + r);
876 Double_t ds=0.011; if (hI==3.5) ds=0.0053;
877 if (!t->PropagateTo(rs,ds)) {
878 return kFALSE;
879 }
880 }
881 }
882
883 Double_t x,y,z;
884 if (!t->GetGlobalXYZat(r,x,y,z)) {
885 return kFALSE;
886 }
887 Double_t phi=TMath::ATan2(y,x);
888 Int_t idet=layer.FindDetectorIndex(phi,z);
889 if (idet<0) {
890 return kFALSE;
891 }
892 const AliITSdetector &det=layer.GetDetector(idet);
893 phi=det.GetPhi();
894 if (!t->Propagate(phi,det.GetR())) {
895 return kFALSE;
896 }
897 t->SetDetectorIndex(idet);
898
899 const AliITSclusterV2 *cl=0;
900 Double_t maxchi2=kMaxChi2;
901
902 Int_t idx=index[i];
903 if (idx>0) {
904 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
905 if (idet != c->GetDetectorIndex()) {
906 idet=c->GetDetectorIndex();
907 const AliITSdetector &det=layer.GetDetector(idet);
908 if (!t->Propagate(det.GetPhi(),det.GetR())) {
909 return kFALSE;
910 }
911 t->SetDetectorIndex(idet);
912 }
913 Double_t chi2=t->GetPredictedChi2(c);
914 if (chi2<maxchi2) { cl=c; maxchi2=chi2; }
915 else return kFALSE;
916 }
917
918 /*
919 if (cl==0)
920 if (t->GetNumberOfClusters()>2) {
921 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
922 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
923 Double_t zmin=t->GetZ() - dz;
924 Double_t zmax=t->GetZ() + dz;
925 Double_t ymin=t->GetY() + phi*r - dy;
926 Double_t ymax=t->GetY() + phi*r + dy;
927 layer.SelectClusters(zmin,zmax,ymin,ymax);
928
929 const AliITSclusterV2 *c=0; Int_t ci=-1;
930 while ((c=layer.GetNextCluster(ci))!=0) {
931 if (idet != c->GetDetectorIndex()) continue;
932 Double_t chi2=t->GetPredictedChi2(c);
933 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
934 }
935 }
936 */
937
938 if (cl) {
939 if (!t->Update(cl,maxchi2,idx)) {
940 return kFALSE;
941 }
942 }
943
944 {
945 Double_t d=layer.GetThickness(t->GetY(),t->GetZ());
946 t->CorrectForMaterial(-step*d);
947 }
948
949 }
950
951 if (!t->PropagateTo(x,0.,0.)) return kFALSE;
952 return kTRUE;
953}
954
955
956