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