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