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