]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerV2.cxx
Changes done in order to remove compilation warnings
[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 AliITSclusterV2 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 "AliITSclusterV2.h"
34 #include "AliITStrackerV2.h"
35
36 ClassImp(AliITStrackerV2)
37
38 AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[kMaxLayer]; // ITS layers
39
40 AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
41   //--------------------------------------------------------------------
42   //This is the AliITStrackerV2 constructor
43   //--------------------------------------------------------------------
44   AliITSgeom *g=(AliITSgeom*)geom;
45
46   Float_t x,y,z;  Int_t i;
47   for (i=1; i<kMaxLayer+1; i++) {
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
64     new (fgLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
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
71         Double_t phi=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
72         phi+=TMath::Pi()/2;
73         if (i==1) phi+=TMath::Pi();
74
75         if (phi<0) phi+=TMath::TwoPi();
76         else if (phi>=TMath::TwoPi()) phi-=TMath::TwoPi();
77
78         Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
79         Double_t r=x*cp+y*sp;
80
81         AliITSdetector &det=fgLayers[i-1].GetDetector((j-1)*ndet + k-1); 
82         new(&det) AliITSdetector(r,phi); 
83       } 
84     }  
85
86   }
87
88   fI=kMaxLayer;
89
90   fPass=0;
91   fConstraint[0]=1; fConstraint[1]=0;
92
93   Double_t xyz[]={kXV,kYV,kZV}, ers[]={kSigmaXV,kSigmaYV,kSigmaZV}; 
94   SetVertex(xyz,ers);
95
96   for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=kLayersNotToSkip[i];
97   fLastLayerToTrackTo=kLastLayerToTrackTo;
98
99 }
100
101 void 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];
106 }
107
108 Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
109   //--------------------------------------------------------------------
110   //This function loads ITS clusters
111   //--------------------------------------------------------------------
112   TBranch *branch=cTree->GetBranch("Clusters");
113   if (!branch) { 
114     Error("LoadClusters"," can't get the branch !\n");
115     return 1;
116   }
117
118   TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
119   branch->SetAddress(&clusters);
120
121   Int_t j=0;
122   for (Int_t i=0; i<kMaxLayer; i++) {
123     Int_t ndet=fgLayers[i].GetNdetectors();
124     Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
125
126     Double_t r=fgLayers[i].GetR();
127     Double_t circ=TMath::TwoPi()*r;
128
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);
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
140         fgLayers[i].InsertCluster(new AliITSclusterV2(*c));
141       }
142       clusters->Delete();
143     }
144     fgLayers[i].ResetRoad(); //road defined by the cluster density
145   }
146
147   return 0;
148 }
149
150 void AliITStrackerV2::UnloadClusters() {
151   //--------------------------------------------------------------------
152   //This function unloads ITS clusters
153   //--------------------------------------------------------------------
154   for (Int_t i=0; i<kMaxLayer; i++) fgLayers[i].ResetClusters();
155 }
156
157 static 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;
182      if (!t->PropagateTo(riw+0.001,-diw,x0iw)) return 1;
183   } else {
184   ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
185     return 1;
186   }
187   
188   return 0;
189 }
190
191 Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
192   //--------------------------------------------------------------------
193   // This functions reconstructs ITS tracks
194   // The clusters must be already loaded !
195   //--------------------------------------------------------------------
196   TObjArray itsTracks(15000);
197
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
204       if ((esd->GetStatus()&AliESDtrack::kTPCin)==0) continue;
205       if (esd->GetStatus()&AliESDtrack::kTPCout) continue;
206       if (esd->GetStatus()&AliESDtrack::kITSin) continue;
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       }
216       if (TMath::Abs(t->GetD())>4) {
217         delete t;
218         continue;
219       }
220
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();
233
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);
239        if (t==0) continue;           //this track has been already tracked
240        Int_t tpcLabel=t->GetLabel(); //save the TPC track label
241
242        ResetTrackToFollow(*t);
243        ResetBestTrack();
244
245        for (FollowProlongation(); fI<kMaxLayer; fI++) {
246           while (TakeNextProlongation()) FollowProlongation();
247        }
248
249        if (fBestTrack.GetNumberOfClusters() == 0) continue;
250
251        if (fConstraint[fPass]) {
252           ResetTrackToFollow(*t);
253           if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
254           ResetBestTrack();
255        }
256
257        fBestTrack.SetLabel(tpcLabel);
258        fBestTrack.CookdEdx();
259        CookLabel(&fBestTrack,0.); //For comparison only
260        fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
261        UseClusters(&fBestTrack);
262        delete itsTracks.RemoveAt(i);
263        ntrk++;
264      }
265   }
266
267   itsTracks.Delete();
268
269   Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
270
271   return 0;
272 }
273
274 Int_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
286      if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
287      if (esd->GetStatus()&AliESDtrack::kITSout) continue;
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]
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      }
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());
318         //fTrackToFollow.CookdEdx();
319         CookLabel(&fTrackToFollow,0.); //For comparison only
320         fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
321         UseClusters(&fTrackToFollow);
322         ntrk++;
323      }
324      delete t;
325   }
326
327   Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
328
329   return 0;
330 }
331
332 Int_t AliITStrackerV2::RefitInward(AliESD *event) {
333   //--------------------------------------------------------------------
334   // This functions refits ITS tracks using the 
335   // "inward propagated" TPC tracks
336   // The clusters must be loaded !
337   //--------------------------------------------------------------------
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
345     if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
346     if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
347     if (esd->GetStatus()&AliESDtrack::kTPCout)
348     if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
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();
373        CookLabel(&fTrackToFollow,0.); //For comparison only
374
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)) {
391             fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
392             Float_t d=fTrackToFollow.GetD(GetX(),GetY());
393             Float_t z=fTrackToFollow.GetZ()-GetZ();
394             fTrackToFollow.GetESDtrack()->SetImpactParameters(d,z);
395             UseClusters(&fTrackToFollow);
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             }
405             ntrk++;
406          }
407        }
408     }
409     delete t;
410   }
411
412   Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
413
414   return 0;
415 }
416
417 AliCluster *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;
423   return fgLayers[l].GetCluster(c);
424 }
425
426
427 void AliITStrackerV2::FollowProlongation() {
428   //--------------------------------------------------------------------
429   //This function finds a track prolongation 
430   //--------------------------------------------------------------------
431   while (fI>fLastLayerToTrackTo) {
432     Int_t i=fI-1;
433
434     AliITSlayer &layer=fgLayers[i];
435     AliITStrackV2 &track=fTracks[i];
436
437     Double_t r=layer.GetR();
438
439     if (i==3 || i==1) {
440        Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
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)) {
444          //Warning("FollowProlongation","propagation failed !\n");
445          return;
446        }
447     }
448
449     //find intersection
450     Double_t x,y,z;  
451     if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
452       //Warning("FollowProlongation","failed to estimate track !\n");
453       return;
454     }
455     Double_t phi=TMath::ATan2(y,x);
456
457     Int_t idet=layer.FindDetectorIndex(phi,z);
458     if (idet<0) {
459       //Warning("FollowProlongation","failed to find a detector !\n");
460       return;
461     }
462
463     //propagate to the intersection
464     const AliITSdetector &det=layer.GetDetector(idet);
465     phi=det.GetPhi();
466     if (!fTrackToFollow.Propagate(phi,det.GetR())) {
467       //Warning("FollowProlongation","propagation failed !\n");
468       return;
469     }
470     fTrackToFollow.SetDetectorIndex(idet);
471
472     //Select possible prolongations and store the current track estimation
473     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
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     } 
481
482     //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
483     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
484     if (dz > kMaxRoad) {
485       //Warning("FollowProlongation","too broad road in Z !\n");
486       return;
487     }
488
489     if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
490
491     //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
492     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
493     if (dy > kMaxRoad) {
494       //Warning("FollowProlongation","too broad road in Y !\n");
495       return;
496     }
497
498     fI--;
499
500     Double_t zmin=track.GetZ() - dz; 
501     Double_t zmax=track.GetZ() + dz;
502     Double_t ymin=track.GetY() + r*phi - dy;
503     Double_t ymax=track.GetY() + r*phi + dy;
504     if (layer.SelectClusters(zmin,zmax,ymin,ymax)==0) 
505        if (fLayersNotToSkip[fI]) return;  
506
507     if (!TakeNextProlongation()) 
508        if (fLayersNotToSkip[fI]) return;
509
510   } 
511
512   //deal with the best track
513   Int_t ncl=fTrackToFollow.GetNumberOfClusters();
514   Int_t nclb=fBestTrack.GetNumberOfClusters();
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();
521         }
522      }
523   }
524
525 }
526
527 Int_t AliITStrackerV2::TakeNextProlongation() {
528   //--------------------------------------------------------------------
529   // This function takes another track prolongation 
530   //
531   //  dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
532   //--------------------------------------------------------------------
533   AliITSlayer &layer=fgLayers[fI];
534   ResetTrackToFollow(fTracks[fI]);
535
536   Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
537   Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
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   } 
543
544   const AliITSclusterV2 *c=0; Int_t ci=-1;
545   const AliITSclusterV2 *cc=0; Int_t cci=-1;
546   Double_t chi2=kMaxChi2;
547   while ((c=layer.GetNextCluster(ci))!=0) {
548     Int_t idet=c->GetDetectorIndex();
549
550     if (fTrackToFollow.GetDetectorIndex()!=idet) {
551        const AliITSdetector &det=layer.GetDetector(idet);
552        ResetTrackToFollow(fTracks[fI]);
553        if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
554          //Warning("TakeNextProlongation","propagation failed !\n");
555          continue;
556        }
557        fTrackToFollow.SetDetectorIndex(idet);
558        if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
559     }
560
561     if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
562     if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
563
564     Double_t ch2=fTrackToFollow.GetPredictedChi2(c); 
565     if (ch2 > chi2) continue;
566     chi2=ch2;
567     cc=c; cci=ci;
568     break;
569   }
570
571   if (!cc) return 0;
572
573   if (!fTrackToFollow.Update(cc,chi2,(fI<<28)+cci)) {
574      //Warning("TakeNextProlongation","filtering failed !\n");
575      return 0;
576   }
577
578   if (fTrackToFollow.GetNumberOfClusters()>1)
579   if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
580
581   fTrackToFollow.
582     SetSampledEdx(cc->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
583
584   {
585   Double_t x0;
586  Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
587   fTrackToFollow.CorrectForMaterial(d,x0);
588   }
589
590   if (fConstraint[fPass]) {
591     Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
592     Double_t xyz[]={GetX(),GetY(),GetZ()};
593     Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
594     fTrackToFollow.Improve(d,xyz,ers);
595   }
596
597   return 1;
598 }
599
600
601 AliITStrackerV2::AliITSlayer::AliITSlayer() {
602   //--------------------------------------------------------------------
603   //default AliITSlayer constructor
604   //--------------------------------------------------------------------
605   fR=0.; fPhiOffset=0.; fZOffset=0.;
606   fNladders=0; fNdetectors=0;
607   fDetectors=0;
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
613 }
614
615 AliITStrackerV2::AliITSlayer::
616 AliITSlayer(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
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;
628
629   fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
630 }
631
632 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
633   //--------------------------------------------------------------------
634   // AliITSlayer destructor
635   //--------------------------------------------------------------------
636   delete[] fDetectors;
637   ResetClusters();
638 }
639
640 void AliITStrackerV2::AliITSlayer::ResetClusters() {
641   //--------------------------------------------------------------------
642   // This function removes loaded clusters
643   //--------------------------------------------------------------------
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    }
651 }
652
653 void AliITStrackerV2::AliITSlayer::ResetRoad() {
654   //--------------------------------------------------------------------
655   // This function calculates the road defined by the cluster density
656   //--------------------------------------------------------------------
657   Int_t n=0;
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++;
662   }
663   if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
664 }
665
666 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
667   //--------------------------------------------------------------------
668   // This function inserts a cluster to this layer in increasing
669   // order of the cluster's fZ
670   //--------------------------------------------------------------------
671   Float_t circ=TMath::TwoPi()*fR;
672   Int_t sec=Int_t(kNsector*c->GetPhiR()/circ);
673   if (sec>=kNsector) {
674      ::Error("InsertCluster","Wrong sector !\n");
675      return 1;
676   }
677   Int_t &n=fN[sec];
678   if (n>=kMaxClusterPerSector) {
679      ::Error("InsertCluster","Too many clusters !\n");
680      return 1;
681   }
682   if (n==0) fClusters[sec*kMaxClusterPerSector]=c;
683   else {
684      Int_t i=FindClusterIndex(c->GetZ(),sec);
685      Int_t k=n-i+sec*kMaxClusterPerSector;
686      memmove(fClusters+i+1 ,fClusters+i,k*sizeof(AliITSclusterV2*));
687      fClusters[i]=c;
688   }
689   n++;
690   return 0;
691 }
692
693 Int_t 
694 AliITStrackerV2::AliITSlayer::FindClusterIndex(Float_t z,Int_t s) const {
695   //--------------------------------------------------------------------
696   // For the sector "s", this function returns the index of the first 
697   // with its fZ >= "z". 
698   //--------------------------------------------------------------------
699   Int_t nc=fN[s];
700   if (nc==0) return kMaxClusterPerSector*s;
701
702   Int_t b=kMaxClusterPerSector*s;
703   if (z <= fClusters[b]->GetZ()) return b;
704
705   Int_t e=b+nc-1;
706   if (z > fClusters[e]->GetZ()) return e+1;
707
708   Int_t m=(b+e)/2;
709   for (; b<e; m=(b+e)/2) {
710     if (z > fClusters[m]->GetZ()) b=m+1;
711     else e=m; 
712   }
713   return m;
714 }
715
716 Int_t AliITStrackerV2::AliITSlayer::
717 SelectClusters(Float_t zmin,Float_t zmax,Float_t ymin, Float_t ymax) {
718   //--------------------------------------------------------------------
719   // This function selects clusters within the "window"
720   //--------------------------------------------------------------------
721     Float_t circ=fR*TMath::TwoPi();
722
723     if (ymin>circ) ymin-=circ; else if (ymin<0) ymin+=circ;
724     if (ymax>circ) ymax-=circ; else if (ymax<0) ymax+=circ;
725
726     Int_t i1=Int_t(kNsector*ymin/circ); if (i1==kNsector) i1--;
727     if (fN[i1]!=0) {
728        Float_t ym = (ymax<ymin) ? ymax+circ : ymax;
729        Int_t i=FindClusterIndex(zmin,i1), imax=i1*kMaxClusterPerSector+fN[i1];
730        for (; i<imax; i++) {
731            AliITSclusterV2 *c=fClusters[i];
732            if (c->IsUsed()) continue;
733            if (c->GetZ()>zmax) break;
734            if (c->GetPhiR()<=ymin) continue;
735            if (c->GetPhiR()>ym) continue;
736            fIndex[fNsel++]=i;
737        }
738     }
739
740     Int_t i2=Int_t(kNsector*ymax/circ); if (i2==kNsector) i2--;
741     if (i2==i1) return fNsel;
742
743     if (fN[i2]!=0) {
744        Float_t ym = (ymin>ymax) ? ymin-circ : ymin;
745        Int_t i=FindClusterIndex(zmin,i2), imax=i2*kMaxClusterPerSector+fN[i2];
746        for (; i<imax; i++) {
747            AliITSclusterV2 *c=fClusters[i];
748            if (c->IsUsed()) continue;
749            if (c->GetZ()>zmax) break;
750            if (c->GetPhiR()<=ym) continue;
751            if (c->GetPhiR()>ymax) continue;
752            fIndex[fNsel++]=i;
753        }
754     }
755
756     return fNsel;
757 }
758
759 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
760   //--------------------------------------------------------------------
761   // This function returns clusters within the "window" 
762   //--------------------------------------------------------------------
763   AliITSclusterV2 *c=0;
764   ci=-1;
765   if (fNsel) {
766      fNsel--;
767      ci=fIndex[fNsel]; 
768      c=fClusters[ci];
769   }
770   return c; 
771 }
772
773 Int_t AliITStrackerV2::AliITSlayer::GetNumberOfClusters() const {
774   Int_t n=0;
775   for (Int_t s=0; s<kNsector; s++) n+=fN[s];
776   return n; 
777 }
778
779 Int_t 
780 AliITStrackerV2::AliITSlayer::FindDetectorIndex(Double_t phi,Double_t z)const {
781   //--------------------------------------------------------------------
782   //This function finds the detector crossed by the track
783   //--------------------------------------------------------------------
784   Double_t dphi=-(phi-fPhiOffset);
785   if      (dphi <  0) dphi += 2*TMath::Pi();
786   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
787   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
788   if (np>=fNladders) np-=fNladders;
789   if (np<0)          np+=fNladders;
790
791   Double_t dz=fZOffset-z;
792   Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
793   if (nz>=fNdetectors) return -1;
794   if (nz<0)            return -1;
795
796   return np*fNdetectors + nz;
797 }
798
799 Double_t 
800 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
801 const {
802   //--------------------------------------------------------------------
803   //This function returns the layer thickness at this point (units X0)
804   //--------------------------------------------------------------------
805   Double_t d=0.0085;
806   x0=21.82;
807
808   if (43<fR&&fR<45) { //SSD2
809      Double_t dd=0.0034;
810      d=dd;
811      if (TMath::Abs(y-0.00)>3.40) d+=dd;
812      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
813      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
814      for (Int_t i=0; i<12; i++) {
815        if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
816           if (TMath::Abs(y-0.00)>3.40) d+=dd;
817           d+=0.0034; 
818           break;
819        }
820        if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
821           if (TMath::Abs(y-0.00)>3.40) d+=dd;
822           d+=0.0034; 
823           break;
824        }         
825        if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
826        if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
827      }
828   } else 
829   if (37<fR&&fR<41) { //SSD1
830      Double_t dd=0.0034;
831      d=dd;
832      if (TMath::Abs(y-0.00)>3.40) d+=dd;
833      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
834      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
835      for (Int_t i=0; i<11; i++) {
836        if (TMath::Abs(z-3.9*i)<0.15) {
837           if (TMath::Abs(y-0.00)>3.40) d+=dd;
838           d+=dd; 
839           break;
840        }
841        if (TMath::Abs(z+3.9*i)<0.15) {
842           if (TMath::Abs(y-0.00)>3.40) d+=dd;
843           d+=dd; 
844           break;
845        }         
846        if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
847        if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}         
848      }
849   } else
850   if (13<fR&&fR<26) { //SDD
851      Double_t dd=0.0033;
852      d=dd;
853      if (TMath::Abs(y-0.00)>3.30) d+=dd;
854
855      if (TMath::Abs(y-1.80)<0.55) {
856         d+=0.016;
857         for (Int_t j=0; j<20; j++) {
858           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
859           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
860         } 
861      }
862      if (TMath::Abs(y+1.80)<0.55) {
863         d+=0.016;
864         for (Int_t j=0; j<20; j++) {
865           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
866           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
867         } 
868      }
869
870      for (Int_t i=0; i<4; i++) {
871        if (TMath::Abs(z-7.3*i)<0.60) {
872           d+=dd;
873           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
874           break;
875        }
876        if (TMath::Abs(z+7.3*i)<0.60) {
877           d+=dd; 
878           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
879           break;
880        }
881      }
882   } else
883   if (6<fR&&fR<8) {   //SPD2
884      Double_t dd=0.0063; x0=21.5;
885      d=dd;
886      if (TMath::Abs(y-3.08)>0.5) d+=dd;
887      //if (TMath::Abs(y-3.08)>0.45) d+=dd;
888      if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
889   } else
890   if (3<fR&&fR<5) {   //SPD1
891      Double_t dd=0.0063; x0=21.5;
892      d=dd;
893      if (TMath::Abs(y+0.21)>0.6) d+=dd;
894      //if (TMath::Abs(y+0.21)>0.45) d+=dd;
895      if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
896   }
897
898   return d;
899 }
900
901 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
902 {
903   //--------------------------------------------------------------------
904   //Returns the thickness between the current layer and the vertex (units X0)
905   //--------------------------------------------------------------------
906   Double_t d=0.0028*3*3; //beam pipe
907   Double_t x0=0;
908
909   Double_t xn=fgLayers[fI].GetR();
910   for (Int_t i=0; i<fI; i++) {
911     Double_t xi=fgLayers[i].GetR();
912     d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
913   }
914
915   if (fI>1) {
916     Double_t xi=9.;
917     d+=0.0097*xi*xi;
918   }
919
920   if (fI>3) {
921     Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
922     d+=0.0034*xi*xi;
923   }
924
925   return d/(xn*xn);
926 }
927
928 Bool_t 
929 AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,const AliITStrackV2 *c) {
930   //--------------------------------------------------------------------
931   // This function refits the track "t" at the position "x" using
932   // the clusters from "c"
933   //--------------------------------------------------------------------
934   Int_t index[kMaxLayer];
935   Int_t k;
936   for (k=0; k<kMaxLayer; k++) index[k]=-1;
937   Int_t nc=c->GetNumberOfClusters();
938   for (k=0; k<nc; k++) { 
939     Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
940     index[nl]=idx; 
941   }
942
943   Int_t from, to, step;
944   if (xx > t->GetX()) {
945       from=0; to=kMaxLayer;
946       step=+1;
947   } else {
948       from=kMaxLayer-1; to=-1;
949       step=-1;
950   }
951
952   for (Int_t i=from; i != to; i += step) {
953      AliITSlayer &layer=fgLayers[i];
954      Double_t r=layer.GetR();
955  
956      {
957      Double_t hI=i-0.5*step; 
958      if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {             
959         Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
960         Double_t d=0.0034, x0=38.6; 
961         if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
962         if (!t->PropagateTo(rs,-step*d,x0)) {
963           return kFALSE;
964         }
965      }
966      }
967
968      // remember old position [SR, GSI 18.02.2003]
969      Double_t oldX=0., oldY=0., oldZ=0.;
970      if (t->IsStartedTimeIntegral() && step==1) {
971         t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
972      }
973      //
974
975      Double_t x,y,z;
976      if (!t->GetGlobalXYZat(r,x,y,z)) { 
977        return kFALSE;
978      }
979      Double_t phi=TMath::ATan2(y,x);
980      Int_t idet=layer.FindDetectorIndex(phi,z);
981      if (idet<0) { 
982        return kFALSE;
983      }
984      const AliITSdetector &det=layer.GetDetector(idet);
985      phi=det.GetPhi();
986      if (!t->Propagate(phi,det.GetR())) {
987        return kFALSE;
988      }
989      t->SetDetectorIndex(idet);
990
991      const AliITSclusterV2 *cl=0;
992      Double_t maxchi2=kMaxChi2;
993
994      Int_t idx=index[i];
995      if (idx>0) {
996         const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx); 
997         if (idet != c->GetDetectorIndex()) {
998            idet=c->GetDetectorIndex();
999            const AliITSdetector &det=layer.GetDetector(idet);
1000            if (!t->Propagate(det.GetPhi(),det.GetR())) {
1001              return kFALSE;
1002            }
1003            t->SetDetectorIndex(idet);
1004         }
1005         Double_t chi2=t->GetPredictedChi2(c);
1006         if (chi2<maxchi2) { 
1007           cl=c; 
1008           maxchi2=chi2; 
1009         } else {
1010           return kFALSE;
1011         }
1012      }
1013      /*
1014      if (cl==0)
1015      if (t->GetNumberOfClusters()>2) {
1016         Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1017         Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1018         Double_t zmin=t->GetZ() - dz;
1019         Double_t zmax=t->GetZ() + dz;
1020         Double_t ymin=t->GetY() + phi*r - dy;
1021         Double_t ymax=t->GetY() + phi*r + dy;
1022         layer.SelectClusters(zmin,zmax,ymin,ymax);
1023
1024         const AliITSclusterV2 *c=0; Int_t ci=-1;
1025         while ((c=layer.GetNextCluster(ci))!=0) {
1026            if (idet != c->GetDetectorIndex()) continue;
1027            Double_t chi2=t->GetPredictedChi2(c);
1028            if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1029         }
1030      }
1031      */
1032      if (cl) {
1033        if (!t->Update(cl,maxchi2,idx)) {
1034           return kFALSE;
1035        }
1036        t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
1037      }
1038
1039      {
1040      Double_t x0;
1041      Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1042      t->CorrectForMaterial(-step*d,x0);
1043      }
1044                  
1045      // track time update [SR, GSI 17.02.2003]
1046      if (t->IsStartedTimeIntegral() && step==1) {
1047         Double_t newX, newY, newZ;
1048         t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1049         Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) + 
1050                        (oldZ-newZ)*(oldZ-newZ);
1051         t->AddTimeStep(TMath::Sqrt(dL2));
1052      }
1053      //
1054
1055   }
1056
1057   if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1058   return kTRUE;
1059 }
1060
1061 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1062   //--------------------------------------------------------------------
1063   // This function marks clusters assigned to the track
1064   //--------------------------------------------------------------------
1065   AliTracker::UseClusters(t,from);
1066
1067   AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
1068   if (c->GetSigmaZ2()>0.1) c->UnUse();
1069   c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
1070   if (c->GetSigmaZ2()>0.1) c->UnUse();
1071
1072 }