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