]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerV2.cxx
Enable multievent run reconstruction (M.Ivanov)
[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 //-------------------------------------------------------------------------
21 #include <TFile.h>
22 #include <TTree.h>
23 #include <TRandom.h>
24 #include <iostream.h>
25
26 #include "AliITSgeom.h"
27 #include "AliITSRecPoint.h"
28 #include "AliTPCtrack.h"
29 #include "AliITSclusterV2.h"
30 #include "AliITStrackerV2.h"
31
32 //#define DEBUG
33
34 #ifdef DEBUG
35 Int_t LAB=236;
36 #endif
37
38 extern TRandom *gRandom;
39
40 AliITStrackerV2::AliITSlayer AliITStrackerV2::fLayers[kMaxLayer]; // ITS layers
41
42 AliITStrackerV2::
43 AliITStrackerV2(const AliITSgeom *geom, Int_t eventn) throw (const Char_t *) {
44   //--------------------------------------------------------------------
45   //This is the AliITStracker constructor
46   //It also reads clusters from a root file
47   //--------------------------------------------------------------------
48  
49   fEventN = eventn;  //MI change add event number  - used to generate identifier 
50   fYV=fZV=0.;
51
52   AliITSgeom *g=(AliITSgeom*)geom;
53
54   Float_t x,y,z;
55   Int_t i;
56   for (i=1; i<kMaxLayer+1; i++) {
57     Int_t nlad=g->GetNladders(i);
58     Int_t ndet=g->GetNdetectors(i);
59
60     g->GetTrans(i,1,1,x,y,z); 
61     Double_t r=TMath::Sqrt(x*x + y*y);
62     Double_t poff=TMath::ATan2(y,x);
63     Double_t zoff=z;
64
65     g->GetTrans(i,1,2,x,y,z);
66     r += TMath::Sqrt(x*x + y*y);
67     g->GetTrans(i,2,1,x,y,z);
68     r += TMath::Sqrt(x*x + y*y);
69     g->GetTrans(i,2,2,x,y,z);
70     r += TMath::Sqrt(x*x + y*y);
71     r*=0.25;
72
73     new (fLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
74
75     for (Int_t j=1; j<nlad+1; j++) {
76       for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
77         Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift); 
78         Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
79
80         Double_t r =-x*rot[1] + y*rot[0];         if (i==1) r=-r;
81         Double_t phi=TMath::ATan2(rot[1],rot[0]); if (i==1) phi-=3.1415927;
82         phi+=0.5*TMath::Pi(); 
83         AliITSdetector &det=fLayers[i-1].GetDetector((j-1)*ndet + k-1); 
84
85 if (phi<0) phi += 2*TMath::Pi();
86
87         new(&det) AliITSdetector(r,phi); 
88       } 
89     }  
90
91   }
92
93   try {
94      //Read clusters
95     //MI change 
96     char   cname[100]; 
97     sprintf(cname,"TreeC_ITS_%d",eventn);
98      TTree *cTree=(TTree*)gDirectory->Get(cname);
99     //     TTree *cTree=(TTree*)gDirectory->Get("cTree");
100
101      if (!cTree) throw 
102         ("AliITStrackerV2::AliITStrackerV2 can't get cTree !\n");
103
104      TBranch *branch=cTree->GetBranch("Clusters");
105      if (!branch) throw 
106         ("AliITStrackerV2::AliITStrackerV2 can't get Clusters branch !\n");
107
108      TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
109      branch->SetAddress(&clusters);
110
111      Int_t nentr=(Int_t)cTree->GetEntries();
112      for (i=0; i<nentr; i++) {
113        if (!cTree->GetEvent(i)) continue;
114        //Int_t lay,lad,det; g->GetModuleId(i-1,lay,lad,det);
115        Int_t lay,lad,det; g->GetModuleId(i,lay,lad,det);
116        Int_t ncl=clusters->GetEntriesFast();
117        while (ncl--) {
118          AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
119
120 #ifdef DEBUG
121 if (c->GetLabel(2)!=LAB)
122 if (c->GetLabel(1)!=LAB)
123 if (c->GetLabel(0)!=LAB) continue;
124 cout<<lay-1<<' '<<lad-1<<' '<<det-1<<' '<<c->GetY()<<' '<<c->GetZ()<<endl;
125 #endif
126
127          AliITSdetector &det=fLayers[lay-1].GetDetector(c->GetDetectorIndex());
128          Double_t r=det.GetR(); r=TMath::Sqrt(r*r + c->GetY()*c->GetY());
129          if (r > TMath::Abs(c->GetZ())-0.5)
130             fLayers[lay-1].InsertCluster(new AliITSclusterV2(*c));
131        }
132        clusters->Delete();
133      }
134   }
135   catch (const Char_t *msg) {
136     cerr<<msg<<endl;
137     throw;
138   }
139
140   fI=kMaxLayer;
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   if (!in->IsOpen()) {
155      cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
156      cerr<<"file with TPC tracks is not open !\n";
157      return 1;
158   }
159
160   if (!out->IsOpen()) {
161      cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
162      cerr<<"file for ITS tracks is not open !\n";
163      return 2;
164   }
165
166   in->cd(); 
167   //
168   char   tname[100];
169   sprintf(tname,"TreeT_TPC_%d",fEventN);
170   TTree *tpcTree=(TTree*)in->Get(tname);
171   //TTree *tpcTree=(TTree*)in->Get("TPCf");
172
173   if (!tpcTree) {
174      cerr<<"AliITStrackerV2::Clusters2Tracks() ";
175      cerr<<"can't get a tree with TPC tracks !\n";
176      return 3;
177   }
178
179   AliTPCtrack *itrack=new AliTPCtrack; 
180   tpcTree->SetBranchAddress("tracks",&itrack);
181
182   out->cd();
183   TTree itsTree("ITSf","Tree with ITS tracks");
184   AliITStrackV2 *otrack=&fBestTrack;
185   itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
186
187   Int_t ntrk=0;
188
189   Int_t nentr=(Int_t)tpcTree->GetEntries();
190   for (Int_t i=0; i<nentr; i++) {
191
192     if (!tpcTree->GetEvent(i)) continue;
193     Int_t tpcLabel=itrack->GetLabel(); //save the TPC track label
194
195 #ifdef DEBUG
196 lbl=tpcLabel;
197 if (TMath::Abs(tpcLabel)!=LAB) continue;
198 cout<<tpcLabel<<" *****************\n";
199 #endif
200
201     try {
202       ResetTrackToFollow(AliITStrackV2(*itrack));
203     } catch (const Char_t *msg) {
204       cerr<<msg<<endl;
205       continue;
206     }
207     ResetBestTrack();
208
209     fYV=gRandom->Gaus(0.,kSigmaYV); fZV=gRandom->Gaus(0.,kSigmaZV);
210
211     Double_t r2=fTrackToFollow.GetX()*fTrackToFollow.GetX();
212     Double_t x0=0.082/21.82*2.33*(45*45+40*40)/r2+2.000/41.28*0.03*75*75/r2;
213     fTrackToFollow.Improve(x0,fYV,fZV);
214
215     Double_t xk=80.;
216     fTrackToFollow.PropagateTo(xk,0.,0.); //Ne if it's still there
217
218     xk-=0.005;
219     fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
220     xk-=0.02;
221     fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);   //Kevlar
222        xk-=2.0;
223        fTrackToFollow.PropagateTo(xk, 2.0/41.28*0.029, 2.0*0.029);//Nomex
224     xk-=0.02;
225     fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);   //Kevlar
226     xk-=0.005;
227     fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
228
229     xk=61.;
230     fTrackToFollow.PropagateTo(xk,0.,0.); //C02
231
232     xk -=0.005;
233     fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7);    //Al    
234     xk -=0.005;
235     fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71);  //Tedlar
236     xk -=0.02;
237     fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);    //Kevlar
238        xk -=0.5;
239        fTrackToFollow.PropagateTo(xk, 0.5/41.28*.029, 0.5*0.029);  //Nomex    
240     xk -=0.02;
241     fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);    //Kevlar
242     xk -=0.005;
243     fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71);  //Tedlar
244     xk -=0.005;
245     fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7);    //Al    
246
247
248     for (FollowProlongation(); fI<kMaxLayer; fI++) {
249        while (TakeNextProlongation()) FollowProlongation();
250     }
251
252 #ifdef DEBUG
253 cout<<fBestTrack.GetNumberOfClusters()<<" number of clusters\n\n";
254 #endif
255
256     if (fBestTrack.GetNumberOfClusters() >= kMaxLayer-kLayersToSkip) {
257        cerr << ++ntrk << "                \r";
258        fBestTrack.SetLabel(tpcLabel);
259        UseClusters(&fBestTrack);
260        itsTree.Fill();
261     }
262
263   }
264   sprintf(tname,"TreeT_ITS_%d",fEventN);
265   itsTree.Write(tname);
266   savedir->cd();
267   cerr<<"Number of TPC tracks: "<<nentr<<endl;
268   cerr<<"Number of prolonged tracks: "<<ntrk<<endl;
269
270   delete itrack;
271
272   return 0;
273 }
274
275 Int_t AliITStrackerV2::PropagateBack(const TFile *inp, TFile *out) {
276   //--------------------------------------------------------------------
277   //This functions propagates reconstructed ITS tracks back
278   //--------------------------------------------------------------------
279   TFile *in=(TFile*)inp;
280   TDirectory *savedir=gDirectory; 
281
282   if (!in->IsOpen()) {
283      cerr<<"AliITStrackerV2::PropagateBack(): ";
284      cerr<<"file with ITS tracks is not open !\n";
285      return 1;
286   }
287
288   if (!out->IsOpen()) {
289      cerr<<"AliITStrackerV2::PropagateBack(): ";
290      cerr<<"file for back propagated ITS tracks is not open !\n";
291      return 2;
292   }
293
294   in->cd();
295   TTree *itsTree=(TTree*)in->Get("ITSf");
296   if (!itsTree) {
297      cerr<<"AliITStrackerV2::PropagateBack() ";
298      cerr<<"can't get a tree with ITS tracks !\n";
299      return 3;
300   }
301   AliITStrackV2 *itrack=new AliITStrackV2; 
302   itsTree->SetBranchAddress("tracks",&itrack);
303
304   out->cd();
305   TTree backTree("ITSb","Tree with back propagated ITS tracks");
306   AliTPCtrack *otrack=0;
307   backTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
308
309   Int_t ntrk=0;
310
311   Int_t nentr=(Int_t)itsTree->GetEntries();
312   for (Int_t i=0; i<nentr; i++) {
313     itsTree->GetEvent(i);
314     ResetTrackToFollow(*itrack);
315     fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
316     Int_t itsLabel=fTrackToFollow.GetLabel(); //save the ITS track label
317
318 #ifdef DEBUG
319 if (TMath::Abs(itsLabel)!=LAB) continue;
320 cout<<itsLabel<<" *****************\n";
321 #endif
322
323     try {
324        Int_t nc=itrack->GetNumberOfClusters();
325 #ifdef DEBUG
326 for (Int_t k=0; k<nc; k++) {
327     Int_t index=itrack->GetClusterIndex(k);
328     AliITSclusterV2 *c=(AliITSclusterV2*)GetCluster(index);
329     cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl;
330 }
331 #endif       
332        Int_t idx=0, l=0; 
333        const  AliITSclusterV2 *c=0; 
334        if (nc--) {
335           idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
336           c=(AliITSclusterV2*)GetCluster(idx);
337        }
338        for (fI=0; fI<kMaxLayer; fI++) {
339          AliITSlayer &layer=fLayers[fI];
340          Double_t r=layer.GetR();
341          if (fI==2 || fI==4) {             
342             Double_t rs=0.5*(fLayers[fI-1].GetR() + r);
343             Double_t ds=0.034; if (fI==4) ds=0.039;
344             Double_t dx0r=ds/21.82*2.33, dr=ds*2.33;
345             if (!fTrackToFollow.PropagateTo(rs,1*dx0r,dr)) throw "";
346          }
347
348          Double_t x,y,z;
349          if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) 
350             throw "AliITStrackerV2::PropagateBack: failed to estimate track !";
351          Double_t phi=TMath::ATan2(y,x);
352          Double_t d=layer.GetThickness(phi,z);
353          Int_t idet=layer.FindDetectorIndex(phi,z);
354          if (idet<0) 
355          throw "AliITStrackerV2::PropagateBack: failed to find a detector !\n";
356          const AliITSdetector &det=layer.GetDetector(idet);
357          r=det.GetR(); phi=det.GetPhi();
358          if (!fTrackToFollow.Propagate(phi,r,d/21.82*2.33,d*2.33)) throw "";
359
360          const AliITSclusterV2 *cl=0;
361          Int_t index=0;
362          Double_t maxchi2=kMaxChi2;
363
364          if (l==fI) {
365            if (idet != c->GetDetectorIndex()) {
366              idet=c->GetDetectorIndex();
367              const AliITSdetector &det=layer.GetDetector(idet);
368              r=det.GetR(); phi=det.GetPhi();
369              if (!fTrackToFollow.Propagate(phi,r,0.,0.)) throw "";
370            }
371            Double_t chi2=fTrackToFollow.GetPredictedChi2(c);
372            if (chi2<kMaxChi2) {
373               cl=c; maxchi2=chi2; index=idx;
374            }
375            if (nc--) {
376               idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
377               c=(AliITSclusterV2*)GetCluster(idx);
378            } 
379          }
380
381          if (fTrackToFollow.GetNumberOfClusters()>2) {
382            Double_t dz=3*TMath::Sqrt(fTrackToFollow.GetSigmaZ2()+kSigmaZ2[fI]);
383            Double_t dy=3*TMath::Sqrt(fTrackToFollow.GetSigmaY2()+kSigmaY2[fI]);
384            Double_t zmin=fTrackToFollow.GetZ() - dz;
385            Double_t zmax=fTrackToFollow.GetZ() + dz;
386            Double_t ymin=fTrackToFollow.GetY() + phi*r - dy;
387            Double_t ymax=fTrackToFollow.GetY() + phi*r + dy;
388            layer.SelectClusters(zmin,zmax,ymin,ymax);
389
390            const AliITSclusterV2 *cc=0; Int_t ci;
391            while ((cc=layer.GetNextCluster(ci))!=0) {
392               if (idet != cc->GetDetectorIndex()) continue;
393               Double_t chi2=fTrackToFollow.GetPredictedChi2(cc);
394               if (chi2<maxchi2) {
395                  cl=cc; index=(fI<<28)+ci; maxchi2=chi2;
396               }
397            }
398          }
399
400          if (cl) {
401             if (!fTrackToFollow.Update(cl,maxchi2,index)) 
402               cerr<<"AliITStrackerV2::PropagateBack: filtering failed !\n";
403               continue;
404          }
405        }
406
407        Double_t xk=61.;
408        fTrackToFollow.PropagateTo(xk,0.,0.); //Air
409
410        xk +=0.005;
411        fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7);    //Al    
412        xk +=0.005;
413        fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71);  //Tedlar
414        xk +=0.02;
415        fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);    //Kevlar
416           xk +=0.5;
417           fTrackToFollow.PropagateTo(xk, 0.5/41.28*.029, 0.5*0.029);  //Nomex 
418        xk +=0.02;
419        fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);    //Kevlar
420        xk +=0.005;
421        fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71);  //Tedlar
422        xk +=0.005;
423        fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7);    //Al    
424
425        xk=80.;
426        fTrackToFollow.PropagateTo(xk,0.,0.); //CO2
427
428        xk+=0.005;
429        fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
430        xk+=0.02;
431        fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);   //Kevlar
432           xk+=2.0;
433           fTrackToFollow.PropagateTo(xk, 2.0/41.28*0.029, 2.0*0.029);//Nomex
434        xk+=0.02;
435        fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45);   //Kevlar
436        xk+=0.005;
437        fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar
438
439        fTrackToFollow.SetLabel(itsLabel);
440        otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha()); 
441        backTree.Fill(); delete otrack;
442        UseClusters(&fTrackToFollow);
443        cerr << ++ntrk << "                \r";
444     }
445     catch (const Char_t *msg) {
446        cerr<<msg<<endl;
447     }
448   }
449
450   backTree.Write();
451   savedir->cd();
452   cerr<<"Number of ITS tracks: "<<nentr<<endl;
453   cerr<<"Number of back propagated ITS tracks: "<<ntrk<<endl;
454
455   delete itrack;
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     fI--;
478 #ifdef DEBUG
479 cout<<fI<<' ';
480 #endif
481     AliITSlayer &layer=fLayers[fI];
482     AliITStrackV2 &track=fTracks[fI];
483
484     Double_t r=layer.GetR();
485     if (fI==3 || fI==1) {
486       Double_t rs=0.5*(fLayers[fI+1].GetR() + r);
487       Double_t ds=0.034; if (fI==3) ds=0.039;
488       Double_t dx0r=ds/21.82*2.33, dr=ds*2.33;
489       fTrackToFollow.Propagate(fTrackToFollow.GetAlpha(),rs,1*dx0r,dr);
490     }
491
492     //find intersection
493     Double_t x,y,z;  
494     if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
495      cerr<<"AliITStrackerV2::FollowProlongation: failed to estimate track !\n";
496      break;
497     }
498     Double_t phi=TMath::ATan2(y,x);
499     Double_t d=layer.GetThickness(phi,z);
500     Int_t idet=layer.FindDetectorIndex(phi,z);
501     if (idet<0) {
502     cerr<<"AliITStrackerV2::FollowProlongation: failed to find a detector !\n";
503       break;
504     }
505
506     //propagate to the intersection
507     const AliITSdetector &det=layer.GetDetector(idet);
508     phi=det.GetPhi();
509     if (!fTrackToFollow.Propagate(phi,det.GetR(),1*d/21.82*2.33,d*2.33)) {
510       cerr<<"AliITStrackerV2::FollowProlongation: propagation failed !\n";
511       break;
512     }
513     fTrackToFollow.SetDetectorIndex(idet);
514
515     //Select possible prolongations and store the current track estimation
516     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
517     Double_t dz=3*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[fI]);
518     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
519     if (dz > kMaxRoad/4) {
520       //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n";
521       break;
522     }
523     Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[fI]);
524     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
525     if (dy > kMaxRoad/4) {
526       //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n";
527       break;
528     }
529
530     Double_t zmin=track.GetZ() - dz;
531     Double_t zmax=track.GetZ() + dz;
532     Double_t ymin=track.GetY() + r*phi - dy;
533     Double_t ymax=track.GetY() + r*phi + dy;
534     layer.SelectClusters(zmin,zmax,ymin,ymax);
535
536     //take another prolongation
537     if (!TakeNextProlongation()) if (!tryAgain--) break;
538     tryAgain=kLayersToSkip;
539
540   } 
541
542   //deal with the best track
543   Int_t ncl=fTrackToFollow.GetNumberOfClusters();
544   Int_t nclb=fBestTrack.GetNumberOfClusters();
545   if (ncl)
546   if (ncl >= nclb) {
547      Double_t chi2=fTrackToFollow.GetChi2();
548      if (chi2/ncl < kChi2PerCluster) {        
549         if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
550            ResetBestTrack();
551         }
552      }
553   }
554
555   if (fI) fI++;
556 }
557
558
559 Int_t AliITStrackerV2::TakeNextProlongation() {
560   //--------------------------------------------------------------------
561   //This function takes another track prolongation 
562   //--------------------------------------------------------------------
563   //Double_t m[20];
564   Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
565
566   AliITSlayer &layer=fLayers[fI];
567   AliITStrackV2 &t=fTracks[fI];
568
569   Double_t dz=4*TMath::Sqrt(t.GetSigmaZ2() + kSigmaZ2[fI]);
570   Double_t dy=4*TMath::Sqrt(t.GetSigmaY2() + kSigmaY2[fI]);
571
572   const AliITSclusterV2 *c=0; Int_t ci=-1;
573   Double_t chi2=12345.;
574   while ((c=layer.GetNextCluster(ci))!=0) {
575
576 #ifdef DEBUG
577 //if (fI==0)
578 //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue; 
579 #endif
580
581     Int_t idet=c->GetDetectorIndex();
582
583     if (t.GetDetectorIndex()!=idet) {
584        const AliITSdetector &det=layer.GetDetector(idet);
585        if (!t.Propagate(det.GetPhi(),det.GetR(),0.,0.)) {
586          cerr<<"AliITStrackerV2::TakeNextProlongation: propagation failed !\n";
587          continue;
588        }
589        t.SetDetectorIndex(idet);
590
591 #ifdef DEBUG
592 cout<<fI<<" change detector !\n";
593 #endif
594
595     }
596
597     if (TMath::Abs(t.GetZ() - c->GetZ()) > dz) continue;
598     if (TMath::Abs(t.GetY() - c->GetY()) > dy) continue;
599
600     //m[0]=fYV; m[1]=fZV;
601     //chi2=t.GetPredictedChi2(c,m,d/21.82*2.33);
602     chi2=t.GetPredictedChi2(c);
603
604     if (chi2<kMaxChi2) break;
605   }
606
607 #ifdef DEBUG
608 cout<<fI<<" chi2="<<chi2<<' '<<t.GetY()<<' '<<t.GetZ()<<' '<<dy<<' '<<dz<<endl;
609 #endif
610
611   if (chi2>=kMaxChi2) return 0;
612   if (!c) return 0;
613
614   ResetTrackToFollow(t);
615
616   //if (!fTrackToFollow.Update(m,chi2,(fI<<28)+ci)) {
617   if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
618      cerr<<"AliITStrackerV2::TakeNextProlongation: filtering failed !\n";
619      return 0;
620   }
621   fTrackToFollow.Improve(d/21.82*2.33,fYV,fZV);
622
623 #ifdef DEBUG
624 cout<<"accepted lab="<<c->GetLabel(0)<<' '
625     <<fTrackToFollow.GetNumberOfClusters()<<' '
626     <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<endl<<endl;
627 #endif
628
629   return 1;
630 }
631
632
633
634
635 AliITStrackerV2::AliITSlayer::AliITSlayer() {
636   //--------------------------------------------------------------------
637   //default AliITSlayer constructor
638   //--------------------------------------------------------------------
639   fN=0;
640   fDetectors=0;
641 }
642
643 AliITStrackerV2::AliITSlayer::
644 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
645   //--------------------------------------------------------------------
646   //main AliITSlayer constructor
647   //--------------------------------------------------------------------
648   fR=r; fPhiOffset=p; fZOffset=z;
649   fNladders=nl; fNdetectors=nd;
650   fDetectors=new AliITSdetector[fNladders*fNdetectors];
651
652   fN=0;
653   fI=0;
654 }
655
656 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
657   //--------------------------------------------------------------------
658   // AliITSlayer destructor
659   //--------------------------------------------------------------------
660   delete[] fDetectors;
661   for (Int_t i=0; i<fN; i++) delete fClusters[i];
662 }
663
664 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
665   //--------------------------------------------------------------------
666   //This function adds a cluster to this layer
667   //--------------------------------------------------------------------
668   if (fN==kMaxClusterPerLayer) {
669  cerr<<"AliITStrackerV2::AliITSlayer::InsertCluster(): Too many clusters !\n"; 
670    return 1;
671   }
672
673   if (fN==0) {fClusters[fN++]=c; return 0;}
674   Int_t i=FindClusterIndex(c->GetZ());
675   memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
676   fClusters[i]=c; fN++;
677
678   return 0;
679 }
680
681 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
682   //--------------------------------------------------------------------
683   // This function returns the index of the nearest cluster 
684   //--------------------------------------------------------------------
685   if (fN==0) return 0;
686   if (z <= fClusters[0]->GetZ()) return 0;
687   if (z > fClusters[fN-1]->GetZ()) return fN;
688   Int_t b=0, e=fN-1, m=(b+e)/2;
689   for (; b<e; m=(b+e)/2) {
690     if (z > fClusters[m]->GetZ()) b=m+1;
691     else e=m; 
692   }
693   return m;
694 }
695
696 void AliITStrackerV2::AliITSlayer::
697 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
698   //--------------------------------------------------------------------
699   // This function sets the "window"
700   //--------------------------------------------------------------------
701   fI=FindClusterIndex(zmin); fZmax=zmax;
702   Double_t circle=2*TMath::Pi()*fR;
703   if (ymax>circle) { ymax-=circle; ymin-=circle; }
704   fYmin=ymin; fYmax=ymax;
705 }
706
707 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
708   //--------------------------------------------------------------------
709   // This function returns clusters within the "window" 
710   //--------------------------------------------------------------------
711   const AliITSclusterV2 *cluster=0;
712   for (Int_t i=fI; i<fN; i++) {
713     const AliITSclusterV2 *c=fClusters[i];
714     if (c->GetZ() > fZmax) break;
715     if (c->IsUsed()) continue;
716     const AliITSdetector &det=GetDetector(c->GetDetectorIndex());    
717     Double_t y=fR*det.GetPhi() + c->GetY();
718
719     if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
720     if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
721
722     if (y<fYmin) continue;
723     if (y>fYmax) continue;
724     cluster=c; ci=i;
725     fI=i+1;
726     break; 
727   }
728   return cluster;
729 }
730
731 Int_t AliITStrackerV2::AliITSlayer::
732 FindDetectorIndex(Double_t phi, Double_t z) const {
733   //--------------------------------------------------------------------
734   //This function finds the detector crossed by the track
735   //--------------------------------------------------------------------
736   Double_t dphi=phi-fPhiOffset;
737   if      (dphi <  0) dphi += 2*TMath::Pi();
738   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
739   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
740
741   Double_t dz=fZOffset-z;
742   Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
743
744   if (np>=fNladders) np-=fNladders;
745   if (np<0)          np+=fNladders;
746
747 #ifdef DEBUG
748 cout<<np<<' '<<nz<<endl;
749 #endif
750
751   return np*fNdetectors + nz;
752 }
753
754 Double_t 
755 AliITStrackerV2::AliITSlayer::GetThickness(Double_t phi, Double_t z) const {
756   //--------------------------------------------------------------------
757   //This function returns the thickness of this layer
758   //--------------------------------------------------------------------
759   //-pi<phi<+pi
760   if (3 <fR&&fR<8 ) return 1.1*0.096;
761   if (13<fR&&fR<26) return 1.1*0.088;
762   if (37<fR&&fR<41) return 1.1*0.085;
763   return 1.1*0.081;
764 }
765
766
767 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t phi,Double_t z) const
768 {
769   //--------------------------------------------------------------------
770   //Returns the thickness between the current layer and the vertex
771   //--------------------------------------------------------------------
772   //-pi<phi<+pi
773   Double_t d=0.1;
774
775   Double_t xn=fLayers[fI].GetR();
776   for (Int_t i=0; i<fI; i++) {
777     Double_t xi=fLayers[i].GetR();
778     d+=fLayers[i].GetThickness(phi,z)*xi*xi;
779   }
780
781   if (fI>1) {
782     Double_t xi=0.5*(fLayers[1].GetR()+fLayers[2].GetR());
783     d+=0.034*xi*xi;
784   }
785
786   if (fI>3) {
787     Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
788     d+=0.039*xi*xi;
789   }
790   return d/(xn*xn);
791 }
792
793
794
795 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
796   //--------------------------------------------------------------------
797   // This function returns number of clusters within the "window" 
798   //--------------------------------------------------------------------
799   Int_t ncl=0;
800   for (Int_t i=fI; i<fN; i++) {
801     const AliITSclusterV2 *c=fClusters[i];
802     if (c->GetZ() > fZmax) break;
803     //if (c->IsUsed()) continue;
804     const AliITSdetector &det=GetDetector(c->GetDetectorIndex());    
805     Double_t y=fR*det.GetPhi() + c->GetY();
806
807     if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
808     if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
809
810     if (y<fYmin) continue;
811     if (y>fYmax) continue;
812     ncl++;
813   }
814   return ncl;
815 }
816
817
818