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