]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerV2.cxx
75a7168de127dc6cb465047aa0a8f47bfdd3c6e9
[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             
263         fBestTrack.PropagateTo(3.,0.0028,65.19);
264         fBestTrack.PropagateToVertex();
265         
266         itsTree.Fill();
267         UseClusters(&fBestTrack);
268         delete itsTracks.RemoveAt(i);
269
270      }
271   }
272
273   itsTree.Write();
274
275   itsTracks.Delete();
276
277   UnloadClusters();
278
279   savedir->cd();
280   cerr<<"Number of TPC tracks: "<<nentr<<endl;
281   cerr<<"Number of prolonged tracks: "<<itsTree.GetEntries()<<endl;
282
283   return 0;
284 }
285
286 Int_t AliITStrackerV2::PropagateBack(const TFile *inp, TFile *out) {
287   //--------------------------------------------------------------------
288   //This functions propagates reconstructed ITS tracks back
289   //--------------------------------------------------------------------
290   TFile *in=(TFile*)inp;
291   TDirectory *savedir=gDirectory; 
292
293   LoadClusters();
294
295   if (!in->IsOpen()) {
296      cerr<<"AliITStrackerV2::PropagateBack(): ";
297      cerr<<"file with ITS tracks is not open !\n";
298      return 1;
299   }
300
301   if (!out->IsOpen()) {
302      cerr<<"AliITStrackerV2::PropagateBack(): ";
303      cerr<<"file for back propagated ITS tracks is not open !\n";
304      return 2;
305   }
306
307   in->cd();
308
309   Char_t tname[100];
310   sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
311   TTree *itsTree=(TTree*)in->Get(tname);
312   if (!itsTree) {
313      cerr<<"AliITStrackerV2::PropagateBack() ";
314      cerr<<"can't get a tree with ITS tracks !\n";
315      return 3;
316   }
317   AliITStrackV2 *itrack=new AliITStrackV2; 
318   itsTree->SetBranchAddress("tracks",&itrack);
319
320   out->cd();
321
322   sprintf(tname,"TreeT_ITSb_%d",GetEventNumber());
323   TTree backTree(tname,"Tree with back propagated ITS tracks");
324   AliTPCtrack *otrack=0;
325   backTree.Branch("tracks","AliTPCtrack",&otrack,32000,2);
326
327   Int_t ntrk=0;
328
329   Int_t nentr=(Int_t)itsTree->GetEntries();
330   for (Int_t i=0; i<nentr; i++) {
331
332     itsTree->GetEvent(i);
333     ResetTrackToFollow(*itrack);
334
335     // propagete to vertex [SR, GSI 17.02.2003]
336     fTrackToFollow.PropagateTo(3.,0.0028,65.19);
337     fTrackToFollow.PropagateToVertex();
338
339     // Start Time measurement [SR, GSI 17.02.2003]
340     fTrackToFollow.StartTimeIntegral();
341     fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
342     //
343
344     fTrackToFollow.ResetCovariance(); 
345     fTrackToFollow.ResetClusters();
346
347     Int_t itsLabel=fTrackToFollow.GetLabel(); //save the ITS track label
348
349 #ifdef DEBUG
350 if (TMath::Abs(itsLabel)!=LAB) continue;
351 cout<<itsLabel<<" *****************\n";
352 #endif
353
354     try {
355        Int_t nc=itrack->GetNumberOfClusters();
356 #ifdef DEBUG
357 for (Int_t k=0; k<nc; k++) {
358     Int_t index=itrack->GetClusterIndex(k);
359     AliITSclusterV2 *c=(AliITSclusterV2*)GetCluster(index);
360     cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl;
361 }
362 #endif       
363        Int_t idx=0, l=0; 
364        const  AliITSclusterV2 *c=0; 
365        if (nc--) {
366           idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
367           c=(AliITSclusterV2*)GetCluster(idx);
368        }
369        for (fI=0; fI<kMaxLayer; fI++) {
370          AliITSlayer &layer=fLayers[fI];
371          Double_t r=layer.GetR();
372          if (fI==2 || fI==4) {             
373             Double_t rs=0.5*(fLayers[fI-1].GetR() + r);
374             Double_t d=0.0034, x0=38.6;
375             if (fI==2) {rs=9.; d=0.0097; x0=42.;}
376             if (!fTrackToFollow.PropagateTo(rs,-d,x0)) throw "";
377          }
378
379          // remember old position [SR, GSI 18.02.2003]
380          Double_t oldX, oldY, oldZ;
381          fTrackToFollow.GetGlobalXYZat(fTrackToFollow.GetX(),oldX,oldY,oldZ);
382          //
383
384          Double_t x,y,z;
385          if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) 
386             throw "AliITStrackerV2::PropagateBack: failed to estimate track !";
387          Double_t phi=TMath::ATan2(y,x);
388          Int_t idet=layer.FindDetectorIndex(phi,z);
389          if (idet<0) 
390          throw "AliITStrackerV2::PropagateBack: failed to find a detector !\n";
391          const AliITSdetector &det=layer.GetDetector(idet);
392          r=det.GetR(); phi=det.GetPhi();
393          if (!fTrackToFollow.Propagate(phi,r)) throw "";
394          fTrackToFollow.SetDetectorIndex(idet);
395
396          const AliITSclusterV2 *cl=0;
397          Int_t index=0;
398          Double_t maxchi2=kMaxChi2;
399
400          if (l==fI) {
401            idet=c->GetDetectorIndex();
402            if (idet != fTrackToFollow.GetDetectorIndex()) {
403              const AliITSdetector &det=layer.GetDetector(idet);
404              r=det.GetR(); phi=det.GetPhi();
405              if (!fTrackToFollow.Propagate(phi,r)) throw "";
406              fTrackToFollow.SetDetectorIndex(idet);
407            }
408            Double_t chi2=fTrackToFollow.GetPredictedChi2(c);
409            if (chi2<kMaxChi2) {
410               cl=c; maxchi2=chi2; index=idx;
411            }
412            if (nc--) {
413               idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
414               c=(AliITSclusterV2*)GetCluster(idx);
415            } 
416          }
417
418          if (fTrackToFollow.GetNumberOfClusters()>2) {
419            Double_t dz=3*TMath::Sqrt(fTrackToFollow.GetSigmaZ2()+kSigmaZ2[fI]);
420            Double_t dy=3*TMath::Sqrt(fTrackToFollow.GetSigmaY2()+kSigmaY2[fI]);
421            Double_t zmin=fTrackToFollow.GetZ() - dz;
422            Double_t zmax=fTrackToFollow.GetZ() + dz;
423            Double_t ymin=fTrackToFollow.GetY() + phi*r - dy;
424            Double_t ymax=fTrackToFollow.GetY() + phi*r + dy;
425            layer.SelectClusters(zmin,zmax,ymin,ymax);
426
427            const AliITSclusterV2 *cc=0; Int_t ci;
428            while ((cc=layer.GetNextCluster(ci))!=0) {
429               idet=cc->GetDetectorIndex();
430               if (idet != fTrackToFollow.GetDetectorIndex()) continue;
431               Double_t chi2=fTrackToFollow.GetPredictedChi2(cc);
432               if (chi2<maxchi2) {
433                  cl=cc; index=(fI<<28)+ci; maxchi2=chi2;
434               }
435            }
436          }
437
438          if (cl) {
439            if (!fTrackToFollow.Update(cl,maxchi2,index)) 
440              cerr<<"AliITStrackerV2::PropagateBack: filtering failed !\n";
441          }
442          {
443            Double_t x0;
444            x=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
445            fTrackToFollow.CorrectForMaterial(-x,x0); 
446          }
447                  
448          // track time update [SR, GSI 17.02.2003]
449          Double_t newX, newY, newZ;
450          fTrackToFollow.GetGlobalXYZat(fTrackToFollow.GetX(),newX,newY,newZ);
451          Double_t dL2 = (oldX-newX)*(oldX-newX)+(oldY-newY)*(oldY-newY)+(oldZ-newZ)*(oldZ-newZ);
452          fTrackToFollow.AddTimeStep(TMath::Sqrt(dL2));
453          //
454
455        }
456
457        fTrackToFollow.PropagateTo(50.,-0.001);
458        //Double_t xk=52.,x,y,z; fTrackToFollow.GetGlobalXYZat(xk,x,y,z);
459        //if (TMath::Abs(y)<7.77) fTrackToFollow.PropagateTo(xk,-0.19,24.); 
460        fTrackToFollow.PropagateTo(61,-0.0053,30);
461        fTrackToFollow.PropagateTo(80.,-0.0053,30);
462
463        fTrackToFollow.SetLabel(itsLabel);
464        otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha()); 
465        backTree.Fill(); delete otrack;
466        UseClusters(&fTrackToFollow);
467        cerr << ++ntrk << "                \r";
468     }
469     catch (const Char_t *msg) {
470        cerr<<msg<<endl;
471     }
472   }
473
474   backTree.Write();
475
476   cerr<<"Number of ITS tracks: "<<nentr<<endl;
477   cerr<<"Number of back propagated ITS tracks: "<<ntrk<<endl;
478
479   delete itrack;
480
481   delete itsTree; //Thanks to Mariana Bondila
482
483   UnloadClusters();
484
485   savedir->cd();
486
487   return 0;
488 }
489
490 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
491   //--------------------------------------------------------------------
492   //       Return pointer to a given cluster
493   //--------------------------------------------------------------------
494   Int_t l=(index & 0xf0000000) >> 28;
495   Int_t c=(index & 0x0fffffff) >> 00;
496   return fLayers[l].GetCluster(c);
497 }
498
499
500 void AliITStrackerV2::FollowProlongation() {
501   //--------------------------------------------------------------------
502   //This function finds a track prolongation 
503   //--------------------------------------------------------------------
504   Int_t tryAgain=kLayersToSkip;
505
506   while (fI) {
507     Int_t i=fI-1;
508 #ifdef DEBUG
509 cout<<i<<' ';
510 #endif
511     AliITSlayer &layer=fLayers[i];
512     AliITStrackV2 &track=fTracks[i];
513
514     Double_t r=layer.GetR();
515
516     if (i==3 || i==1) {
517        Double_t rs=0.5*(fLayers[i+1].GetR() + r);
518        Double_t d=0.0034, x0=38.6;
519        if (i==1) {rs=9.; d=0.0097; x0=42;}
520        if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
521          //cerr<<"AliITStrackerV2::FollowProlongation: "
522          //"propagation failed !\n";
523          break;
524        }
525     }
526
527     //find intersection
528     Double_t x,y,z;  
529     if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
530       //cerr<<"AliITStrackerV2::FollowProlongation: "
531       //"failed to estimate track !\n";
532       break;
533     }
534     Double_t phi=TMath::ATan2(y,x);
535     Int_t idet=layer.FindDetectorIndex(phi,z);
536     if (idet<0) {
537       //cerr<<"AliITStrackerV2::FollowProlongation: "
538       //"failed to find a detector !\n";
539       break;
540     }
541
542     //propagate to the intersection
543     const AliITSdetector &det=layer.GetDetector(idet);
544     phi=det.GetPhi();
545     if (!fTrackToFollow.Propagate(phi,det.GetR())) {
546       //cerr<<"AliITStrackerV2::FollowProlongation: propagation failed !\n";
547       break;
548     }
549     fTrackToFollow.SetDetectorIndex(idet);
550
551     //Select possible prolongations and store the current track estimation
552     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
553     Double_t dz=3*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
554     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
555     if (dz > kMaxRoad) {
556       //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n";
557       break;
558     }
559
560     if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) break;
561
562     Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
563     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
564     if (dy > kMaxRoad) {
565       //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n";
566       break;
567     }
568
569     Double_t zmin=track.GetZ() - dz; 
570     Double_t zmax=track.GetZ() + dz;
571     Double_t ymin=track.GetY() + r*phi - dy;
572     Double_t ymax=track.GetY() + r*phi + dy;
573     layer.SelectClusters(zmin,zmax,ymin,ymax); 
574     fI--;
575
576     //take another prolongation
577     if (!TakeNextProlongation()) if (!tryAgain--) break;
578     tryAgain=kLayersToSkip;
579
580   } 
581
582   //deal with the best track
583   Int_t ncl=fTrackToFollow.GetNumberOfClusters();
584   Int_t nclb=fBestTrack.GetNumberOfClusters();
585   if (ncl)
586   if (ncl >= nclb) {
587      Double_t chi2=fTrackToFollow.GetChi2();
588      if (chi2/ncl < kChi2PerCluster) {        
589         if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
590            ResetBestTrack();
591         }
592      }
593   }
594
595 }
596
597 Int_t AliITStrackerV2::TakeNextProlongation() {
598   //--------------------------------------------------------------------
599   // This function takes another track prolongation 
600   //
601   //  dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
602   //--------------------------------------------------------------------
603   AliITSlayer &layer=fLayers[fI];
604   ResetTrackToFollow(fTracks[fI]);
605
606   Double_t dz=4*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
607   Double_t dy=4*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
608
609   const AliITSclusterV2 *c=0; Int_t ci=-1;
610   Double_t chi2=12345.;
611   while ((c=layer.GetNextCluster(ci))!=0) {
612     //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue; 
613     Int_t idet=c->GetDetectorIndex();
614
615     if (fTrackToFollow.GetDetectorIndex()!=idet) {
616        const AliITSdetector &det=layer.GetDetector(idet);
617        ResetTrackToFollow(fTracks[fI]);
618        if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
619          //cerr<<"AliITStrackerV2::TakeNextProlongation: "
620          //"propagation failed !\n";
621          continue;
622        }
623        fTrackToFollow.SetDetectorIndex(idet);
624        if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
625
626 #ifdef DEBUG
627 cout<<fI<<" change detector !\n";
628 #endif
629
630     }
631
632     if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
633     if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
634
635     chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
636   }
637
638 #ifdef DEBUG
639 cout<<fI<<" chi2="<<chi2<<' '
640     <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
641     <<dy<<' '<<dz<<endl;
642 {
643 Double_t phi=TMath::ATan(fTrackToFollow.GetY()/fTrackToFollow.GetX())
644             +fTrackToFollow.GetAlpha(); phi=phi*180./3.141593;
645 cout<<phi<<endl;
646 }
647 #endif
648
649   if (chi2>=kMaxChi2) return 0;
650   if (!c) return 0;
651
652   if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
653      //cerr<<"AliITStrackerV2::TakeNextProlongation: filtering failed !\n";
654      return 0;
655   }
656
657   if (fTrackToFollow.GetNumberOfClusters()>1)
658   if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
659
660   fTrackToFollow.
661     SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
662
663   {
664  Double_t x0;
665  Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
666    fTrackToFollow.CorrectForMaterial(d,x0);
667   }
668
669   if (fConstraint[fPass]) {
670     Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
671     fTrackToFollow.Improve(d,GetY(),GetZ());
672   }
673
674 #ifdef DEBUG
675 cout<<"accepted lab="<<c->GetLabel(0)<<' '
676     <<fTrackToFollow.GetNumberOfClusters()<<' '
677     <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
678     <<fTrackToFollow.Get1Pt()<<endl<<endl;
679 #endif
680
681   return 1;
682 }
683
684
685 AliITStrackerV2::AliITSlayer::AliITSlayer() {
686   //--------------------------------------------------------------------
687   //default AliITSlayer constructor
688   //--------------------------------------------------------------------
689   fN=0;
690   fDetectors=0;
691 }
692
693 AliITStrackerV2::AliITSlayer::
694 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
695   //--------------------------------------------------------------------
696   //main AliITSlayer constructor
697   //--------------------------------------------------------------------
698   fR=r; fPhiOffset=p; fZOffset=z;
699   fNladders=nl; fNdetectors=nd;
700   fDetectors=new AliITSdetector[fNladders*fNdetectors];
701
702   fN=0;
703   fI=0;
704 }
705
706 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
707   //--------------------------------------------------------------------
708   // AliITSlayer destructor
709   //--------------------------------------------------------------------
710   delete[] fDetectors;
711   for (Int_t i=0; i<fN; i++) delete fClusters[i];
712 }
713
714 void AliITStrackerV2::AliITSlayer::ResetClusters() {
715   //--------------------------------------------------------------------
716   // This function removes loaded clusters
717   //--------------------------------------------------------------------
718   for (Int_t i=0; i<fN; i++) delete fClusters[i];
719   fN=0;
720   fI=0;
721 }
722
723 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
724   //--------------------------------------------------------------------
725   //This function adds a cluster to this layer
726   //--------------------------------------------------------------------
727   if (fN==kMaxClusterPerLayer) {
728      cerr<<"AliITStrackerV2::AliITSlayer::InsertCluster(): "
729            "Too many clusters !\n";
730      return 1;
731   }
732
733   if (fN==0) {fClusters[fN++]=c; return 0;}
734   Int_t i=FindClusterIndex(c->GetZ());
735   memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
736   fClusters[i]=c; fN++;
737
738   return 0;
739 }
740
741 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
742   //--------------------------------------------------------------------
743   // This function returns the index of the nearest cluster 
744   //--------------------------------------------------------------------
745   if (fN==0) return 0;
746   if (z <= fClusters[0]->GetZ()) return 0;
747   if (z > fClusters[fN-1]->GetZ()) return fN;
748   Int_t b=0, e=fN-1, m=(b+e)/2;
749   for (; b<e; m=(b+e)/2) {
750     if (z > fClusters[m]->GetZ()) b=m+1;
751     else e=m; 
752   }
753   return m;
754 }
755
756 void AliITStrackerV2::AliITSlayer::
757 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
758   //--------------------------------------------------------------------
759   // This function sets the "window"
760   //--------------------------------------------------------------------
761   fI=FindClusterIndex(zmin); fZmax=zmax;
762   Double_t circle=2*TMath::Pi()*fR;
763   if (ymax>circle) { ymax-=circle; ymin-=circle; }
764   fYmin=ymin; fYmax=ymax;
765 }
766
767 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
768   //--------------------------------------------------------------------
769   // This function returns clusters within the "window" 
770   //--------------------------------------------------------------------
771   const AliITSclusterV2 *cluster=0;
772   for (Int_t i=fI; i<fN; i++) {
773     const AliITSclusterV2 *c=fClusters[i];
774     if (c->GetZ() > fZmax) break;
775     if (c->IsUsed()) continue;
776     const AliITSdetector &det=GetDetector(c->GetDetectorIndex());    
777     Double_t y=fR*det.GetPhi() + c->GetY();
778
779     if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
780     if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
781
782     if (y<fYmin) continue;
783     if (y>fYmax) continue;
784     cluster=c; ci=i;
785     fI=i+1;
786     break; 
787   }
788
789   return cluster;
790 }
791
792 Int_t AliITStrackerV2::AliITSlayer::
793 FindDetectorIndex(Double_t phi, Double_t z) const {
794   //--------------------------------------------------------------------
795   //This function finds the detector crossed by the track
796   //--------------------------------------------------------------------
797   Double_t dphi=phi-fPhiOffset;
798   if      (dphi <  0) dphi += 2*TMath::Pi();
799   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
800   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
801   if (np>=fNladders) np-=fNladders;
802   if (np<0)          np+=fNladders;
803
804   Double_t dz=fZOffset-z;
805   Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
806   if (nz>=fNdetectors) return -1;
807   if (nz<0)            return -1;
808
809 #ifdef DEBUG
810 cout<<np<<' '<<nz<<endl;
811 #endif
812
813   return np*fNdetectors + nz;
814 }
815
816 Double_t 
817 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
818 const {
819   //--------------------------------------------------------------------
820   //This function returns the layer thickness at this point (units X0)
821   //--------------------------------------------------------------------
822   Double_t d=0.0085;
823   x0=21.82;
824
825   if (43<fR&&fR<45) { //SSD2
826      Double_t dd=0.0034;
827      d=dd;
828      if (TMath::Abs(y-0.00)>3.40) d+=dd;
829      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
830      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
831      for (Int_t i=0; i<12; i++) {
832        if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
833           if (TMath::Abs(y-0.00)>3.40) d+=dd;
834           d+=0.0034; 
835           break;
836        }
837        if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
838           if (TMath::Abs(y-0.00)>3.40) d+=dd;
839           d+=0.0034; 
840           break;
841        }         
842        if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
843        if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
844      }
845   } else 
846   if (37<fR&&fR<41) { //SSD1
847      Double_t dd=0.0034;
848      d=dd;
849      if (TMath::Abs(y-0.00)>3.40) d+=dd;
850      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
851      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
852      for (Int_t i=0; i<11; i++) {
853        if (TMath::Abs(z-3.9*i)<0.15) {
854           if (TMath::Abs(y-0.00)>3.40) d+=dd;
855           d+=dd; 
856           break;
857        }
858        if (TMath::Abs(z+3.9*i)<0.15) {
859           if (TMath::Abs(y-0.00)>3.40) d+=dd;
860           d+=dd; 
861           break;
862        }         
863        if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
864        if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}         
865      }
866   } else
867   if (13<fR&&fR<26) { //SDD
868      Double_t dd=0.0033;
869      d=dd;
870      if (TMath::Abs(y-0.00)>3.30) d+=dd;
871
872      if (TMath::Abs(y-1.80)<0.55) {
873         d+=0.016;
874         for (Int_t j=0; j<20; j++) {
875           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
876           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
877         } 
878      }
879      if (TMath::Abs(y+1.80)<0.55) {
880         d+=0.016;
881         for (Int_t j=0; j<20; j++) {
882           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
883           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
884         } 
885      }
886
887      for (Int_t i=0; i<4; i++) {
888        if (TMath::Abs(z-7.3*i)<0.60) {
889           d+=dd;
890           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
891           break;
892        }
893        if (TMath::Abs(z+7.3*i)<0.60) {
894           d+=dd; 
895           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
896           break;
897        }
898      }
899   } else
900   if (6<fR&&fR<8) {   //SPD2
901      Double_t dd=0.0063; x0=21.5;
902      d=dd;
903      if (TMath::Abs(y-3.08)>0.5) d+=dd;
904      //if (TMath::Abs(y-3.08)>0.45) d+=dd;
905      if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
906   } else
907   if (3<fR&&fR<5) {   //SPD1
908      Double_t dd=0.0063; x0=21.5;
909      d=dd;
910      if (TMath::Abs(y+0.21)>0.6) d+=dd;
911      //if (TMath::Abs(y+0.21)>0.45) d+=dd;
912      if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
913   }
914
915 #ifdef DEBUG
916   cout<<"d="<<d<<endl;
917 #endif
918
919   return d;
920 }
921
922 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
923 {
924   //--------------------------------------------------------------------
925   //Returns the thickness between the current layer and the vertex (units X0)
926   //--------------------------------------------------------------------
927   Double_t d=0.0028*3*3; //beam pipe
928   Double_t x0=0;
929
930   Double_t xn=fLayers[fI].GetR();
931   for (Int_t i=0; i<fI; i++) {
932     Double_t xi=fLayers[i].GetR();
933     d+=fLayers[i].GetThickness(y,z,x0)*xi*xi;
934   }
935
936   if (fI>1) {
937     Double_t xi=9.;
938     d+=0.0097*xi*xi;
939   }
940
941   if (fI>3) {
942     Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
943     d+=0.0034*xi*xi;
944   }
945
946   return d/(xn*xn);
947 }
948
949 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
950   //--------------------------------------------------------------------
951   // This function returns number of clusters within the "window" 
952   //--------------------------------------------------------------------
953   Int_t ncl=0;
954   for (Int_t i=fI; i<fN; i++) {
955     const AliITSclusterV2 *c=fClusters[i];
956     if (c->GetZ() > fZmax) break;
957     if (c->IsUsed()) continue;
958     const AliITSdetector &det=GetDetector(c->GetDetectorIndex());    
959     Double_t y=fR*det.GetPhi() + c->GetY();
960
961     if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
962     if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
963
964     if (y<fYmin) continue;
965     if (y>fYmax) continue;
966     ncl++;
967   }
968   return ncl;
969 }
970
971 Bool_t AliITStrackerV2::RefitAt(Double_t x, AliITStrackV2 *t, Int_t *index) {
972   //--------------------------------------------------------------------
973   // This function refits a track at a given position
974   //--------------------------------------------------------------------
975   Int_t from, to, step;
976   if (x > t->GetX()) {
977       from=0; to=kMaxLayer;
978       step=+1;
979   } else {
980       from=kMaxLayer-1; to=-1;
981       step=-1;
982   }
983
984   for (Int_t i=from; i != to; i += step) {
985      AliITSlayer &layer=fLayers[i];
986      Double_t r=layer.GetR();
987  
988      {
989      Double_t hI=i-0.5*step; 
990      if (hI==1.5 || hI==3.5) {             
991         Double_t rs=0.5*(fLayers[i-step].GetR() + r);
992         Double_t d=0.0034, x0=38.6; 
993         if (hI==1.5) {rs=9.; d=0.0097; x0=42;}
994         if (!t->PropagateTo(rs,d,x0)) {
995           return kFALSE;
996         }
997      }
998      }
999
1000      Double_t x,y,z;
1001      if (!t->GetGlobalXYZat(r,x,y,z)) { 
1002        return kFALSE;
1003      }
1004      Double_t phi=TMath::ATan2(y,x);
1005      Int_t idet=layer.FindDetectorIndex(phi,z);
1006      if (idet<0) { 
1007        return kFALSE;
1008      }
1009      const AliITSdetector &det=layer.GetDetector(idet);
1010      phi=det.GetPhi();
1011      if (!t->Propagate(phi,det.GetR())) {
1012        return kFALSE;
1013      }
1014      t->SetDetectorIndex(idet);
1015
1016      const AliITSclusterV2 *cl=0;
1017      Double_t maxchi2=kMaxChi2;
1018
1019      Int_t idx=index[i];
1020      if (idx>0) {
1021         const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx); 
1022         if (idet != c->GetDetectorIndex()) {
1023            idet=c->GetDetectorIndex();
1024            const AliITSdetector &det=layer.GetDetector(idet);
1025            if (!t->Propagate(det.GetPhi(),det.GetR())) {
1026              return kFALSE;
1027            }
1028            t->SetDetectorIndex(idet);
1029         }
1030         Double_t chi2=t->GetPredictedChi2(c);
1031         if (chi2<maxchi2) { cl=c; maxchi2=chi2; }
1032         else return kFALSE;
1033      }
1034
1035      /*
1036      if (cl==0)
1037      if (t->GetNumberOfClusters()>2) {
1038         Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1039         Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1040         Double_t zmin=t->GetZ() - dz;
1041         Double_t zmax=t->GetZ() + dz;
1042         Double_t ymin=t->GetY() + phi*r - dy;
1043         Double_t ymax=t->GetY() + phi*r + dy;
1044         layer.SelectClusters(zmin,zmax,ymin,ymax);
1045
1046         const AliITSclusterV2 *c=0; Int_t ci=-1;
1047         while ((c=layer.GetNextCluster(ci))!=0) {
1048            if (idet != c->GetDetectorIndex()) continue;
1049            Double_t chi2=t->GetPredictedChi2(c);
1050            if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1051         }
1052      }
1053      */
1054
1055      if (cl) {
1056        if (!t->Update(cl,maxchi2,idx)) {
1057           return kFALSE;
1058        }
1059      }
1060
1061      {
1062      Double_t x0;
1063      Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1064      t->CorrectForMaterial(-step*d,x0);
1065      }
1066
1067   }
1068
1069   if (!t->PropagateTo(x,0.,0.)) return kFALSE;
1070   return kTRUE;
1071 }
1072
1073
1074