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