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