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