]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ITS/AliITStrackerV2.cxx
- fixing warnings/coverity
[u/mrichter/AliRoot.git] / ITS / AliITStrackerV2.cxx
index 9f40d657d23b0cc996b33039bd0ef5cdc0e6019d..758c7641c155ba2733a867ff64d1e29b789d0777 100644 (file)
@@ -12,7 +12,7 @@
  * about the suitability of this software for any purpose. It is          *
  * provided "as is" without express or implied warranty.                  *
  **************************************************************************/
-
+/* $Id$ */
 //-------------------------------------------------------------------------
 //               Implementation of the ITS tracker class
 //    It reads AliITSRecPoint clusters and creates AliITStrackV2 tracks
 
 #include <new>
 
+#include <TError.h>
 #include <TFile.h>
 #include <TTree.h>
 #include <TRandom.h>
+#include <TGeoMatrix.h>
 
-#include "AliITSgeom.h"
+#include "AliITSgeomTGeo.h"
+#include "AliAlignObj.h"
 #include "AliITSRecPoint.h"
-#include "AliESD.h"
+#include "AliESDEvent.h"
+#include "AliESDtrack.h"
 #include "AliITSRecPoint.h"
+#include "AliITSReconstructor.h"
 #include "AliITStrackerV2.h"
 
 ClassImp(AliITStrackerV2)
 
-AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[kMaxLayer]; // ITS layers
+AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[AliITSgeomTGeo::kNLayers]; //ITS layers
+
+AliITStrackerV2::AliITStrackerV2(): 
+  AliTracker(), 
+  fI(AliITSgeomTGeo::GetNLayers()),
+  fBestTrack(),
+  fTrackToFollow(),
+  fPass(0),
+  fLastLayerToTrackTo(AliITSRecoParam::GetLastLayerToTrackTo())
+{
+  //--------------------------------------------------------------------
+  //This is the AliITStrackerV2 default constructor
+  //--------------------------------------------------------------------
+
+  for (Int_t i=1; i<AliITSgeomTGeo::GetNLayers()+1; i++) new(fgLayers+i-1) AliITSlayer();
+
+  fConstraint[0]=1; fConstraint[1]=0;
+
+  Double_t xyz[]={AliITSReconstructor::GetRecoParam()->GetXVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetYVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetZVdef()}; 
+  Double_t ers[]={AliITSReconstructor::GetRecoParam()->GetSigmaXVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetSigmaYVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetSigmaZVdef()}; 
+  SetVertex(xyz,ers);
+
+  for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=AliITSRecoParam::GetLayersNotToSkip(i);
+
+}
+
+AliITStrackerV2::AliITStrackerV2(const AliITStrackerV2 &t): 
+  AliTracker(t), 
+  fI(t.fI),
+  fBestTrack(t.fBestTrack),
+  fTrackToFollow(t.fTrackToFollow),
+  fPass(t.fPass),
+  fLastLayerToTrackTo(t.fLastLayerToTrackTo)
+{
+  //--------------------------------------------------------------------
+  //This is the AliITStrackerV2 copy constructor
+  //--------------------------------------------------------------------
 
-AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
+  //for (Int_t i=1; i<AliITSgeomTGeo::GetNLayers()+1; i++) new(fgLayers+i-1) AliITSlayer();
+
+  fConstraint[0]=t.fConstraint[0]; fConstraint[1]=t.fConstraint[1];
+
+  Double_t xyz[]={AliITSReconstructor::GetRecoParam()->GetXVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetYVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetZVdef()}; 
+  Double_t ers[]={AliITSReconstructor::GetRecoParam()->GetSigmaXVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetSigmaYVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetSigmaZVdef()}; 
+  xyz[0]=t.GetX(); xyz[1]=t.GetY(); xyz[2]=t.GetZ(); 
+  ers[0]=t.GetSigmaX(); ers[1]=t.GetSigmaY(); ers[2]=t.GetSigmaZ(); 
+  SetVertex(xyz,ers);
+
+  for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=t.fLayersNotToSkip[i];
+
+}
+
+AliITStrackerV2::AliITStrackerV2(const Char_t *geom) : 
+  AliTracker(), 
+  fI(AliITSgeomTGeo::GetNLayers()),
+  fBestTrack(),
+  fTrackToFollow(),
+  fPass(0),
+  fLastLayerToTrackTo(AliITSRecoParam::GetLastLayerToTrackTo())
+{
   //--------------------------------------------------------------------
   //This is the AliITStrackerV2 constructor
   //--------------------------------------------------------------------
-  AliITSgeom *g=(AliITSgeom*)geom;
+  if (geom) {
+    AliWarning("\"geom\" is actually a dummy argument !");
+  }
 
-  Float_t x,y,z;  Int_t i;
-  for (i=1; i<kMaxLayer+1; i++) {
-    Int_t nlad=g->GetNladders(i);
-    Int_t ndet=g->GetNdetectors(i);
+  for (Int_t i=1; i<AliITSgeomTGeo::GetNLayers()+1; i++) {
+    Int_t nlad=AliITSgeomTGeo::GetNLadders(i);
+    Int_t ndet=AliITSgeomTGeo::GetNDetectors(i);
 
-    g->GetTrans(i,1,1,x,y,z); 
-    Double_t r=TMath::Sqrt(x*x + y*y);
+    Double_t xyz[3], &x=xyz[0], &y=xyz[1], &z=xyz[2];
+    AliITSgeomTGeo::GetOrigTranslation(i,1,1,xyz); 
     Double_t poff=TMath::ATan2(y,x);
     Double_t zoff=z;
+    Double_t r=TMath::Sqrt(x*x + y*y);
 
-    g->GetTrans(i,1,2,x,y,z);
+    AliITSgeomTGeo::GetOrigTranslation(i,1,2,xyz);
     r += TMath::Sqrt(x*x + y*y);
-    g->GetTrans(i,2,1,x,y,z);
+    AliITSgeomTGeo::GetOrigTranslation(i,2,1,xyz);
     r += TMath::Sqrt(x*x + y*y);
-    g->GetTrans(i,2,2,x,y,z);
+    AliITSgeomTGeo::GetOrigTranslation(i,2,2,xyz);
     r += TMath::Sqrt(x*x + y*y);
     r*=0.25;
 
@@ -65,19 +137,18 @@ AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
 
     for (Int_t j=1; j<nlad+1; j++) {
       for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
-        Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift); 
-        Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
-
-        Double_t phi=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
-        phi+=TMath::Pi()/2;
-        if (i==1) phi+=TMath::Pi();
+        TGeoHMatrix m; AliITSgeomTGeo::GetOrigMatrix(i,j,k,m);
+        const TGeoHMatrix *tm=AliITSgeomTGeo::GetTracking2LocalMatrix(i,j,k);
+        m.Multiply(tm);
+        Double_t txyz[3]={0.}; 
+       xyz[0]=0.; xyz[1]=0.; xyz[2]=0.;
+        m.LocalToMaster(txyz,xyz);
+        r=TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]);
+        Double_t phi=TMath::ATan2(xyz[1],xyz[0]);
 
         if (phi<0) phi+=TMath::TwoPi();
         else if (phi>=TMath::TwoPi()) phi-=TMath::TwoPi();
 
-        Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
-        Double_t r=x*cp+y*sp;
-
         AliITSdetector &det=fgLayers[i-1].GetDetector((j-1)*ndet + k-1); 
         new(&det) AliITSdetector(r,phi); 
       } 
@@ -85,16 +156,17 @@ AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
 
   }
 
-  fI=kMaxLayer;
-
-  fPass=0;
   fConstraint[0]=1; fConstraint[1]=0;
 
-  Double_t xyz[]={kXV,kYV,kZV}, ers[]={kSigmaXV,kSigmaYV,kSigmaZV}; 
+  Double_t xyz[]={AliITSReconstructor::GetRecoParam()->GetXVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetYVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetZVdef()}; 
+  Double_t ers[]={AliITSReconstructor::GetRecoParam()->GetSigmaXVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetSigmaYVdef(),
+                 AliITSReconstructor::GetRecoParam()->GetSigmaZVdef()}; 
   SetVertex(xyz,ers);
 
-  for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=kLayersNotToSkip[i];
-  fLastLayerToTrackTo=kLastLayerToTrackTo;
+  for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=AliITSRecoParam::GetLayersNotToSkip(i);
 
 }
 
@@ -102,7 +174,7 @@ void AliITStrackerV2::SetLayersNotToSkip(Int_t *l) {
   //--------------------------------------------------------------------
   //This function set masks of the layers which must be not skipped
   //--------------------------------------------------------------------
-  for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=l[i];
+  for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=l[i];
 }
 
 Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
@@ -119,7 +191,7 @@ Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
   branch->SetAddress(&clusters);
 
   Int_t j=0;
-  for (Int_t i=0; i<kMaxLayer; i++) {
+  for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) {
     Int_t ndet=fgLayers[i].GetNdetectors();
     Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
 
@@ -129,11 +201,16 @@ Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
     for (; j<jmax; j++) {           
       if (!cTree->GetEvent(j)) continue;
       Int_t ncl=clusters->GetEntriesFast();
       while (ncl--) {
         AliITSRecPoint *c=(AliITSRecPoint*)clusters->UncheckedAt(ncl);
 
+       if (!c->Misalign()) AliWarning("Can't misalign this cluster !");
+
         Int_t idx=c->GetDetectorIndex();
-        Double_t y=r*fgLayers[i].GetDetector(idx).GetPhi()+c->GetY();
+        AliITSdetector &det=fgLayers[i].GetDetector(idx);
+   
+        Double_t y=r*det.GetPhi()+c->GetY();
         if (y>circ) y-=circ; else if (y<0) y+=circ;
         c->SetPhiR(y);
 
@@ -151,7 +228,7 @@ void AliITStrackerV2::UnloadClusters() {
   //--------------------------------------------------------------------
   //This function unloads ITS clusters
   //--------------------------------------------------------------------
-  for (Int_t i=0; i<kMaxLayer; i++) fgLayers[i].ResetClusters();
+  for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fgLayers[i].ResetClusters();
 }
 
 static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
@@ -188,7 +265,7 @@ static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
   return 0;
 }
 
-Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
+Int_t AliITStrackerV2::Clusters2Tracks(AliESDEvent *event) {
   //--------------------------------------------------------------------
   // This functions reconstructs ITS tracks
   // The clusters must be already loaded !
@@ -205,15 +282,9 @@ Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
       if (esd->GetStatus()&AliESDtrack::kTPCout) continue;
       if (esd->GetStatus()&AliESDtrack::kITSin) continue;
 
-      AliITStrackV2 *t=0;
-      try {
-        t=new AliITStrackV2(*esd);
-      } catch (const Char_t *msg) {
-        Warning("Clusters2Tracks",msg);
-        delete t;
-        continue;
-      }
-      if (TMath::Abs(t->GetD())>4) {
+      AliITStrackV2 *t = new AliITStrackV2(*esd);
+
+      if (TMath::Abs(t->GetD(GetX(),GetY()))>4) {
        delete t;
        continue;
       }
@@ -242,7 +313,7 @@ Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
        ResetTrackToFollow(*t);
        ResetBestTrack();
 
-       for (FollowProlongation(); fI<kMaxLayer; fI++) {
+       for (FollowProlongation(); fI<AliITSgeomTGeo::GetNLayers(); fI++) {
           while (TakeNextProlongation()) FollowProlongation();
        }
 
@@ -271,7 +342,7 @@ Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
   return 0;
 }
 
-Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
+Int_t AliITStrackerV2::PropagateBack(AliESDEvent *event) {
   //--------------------------------------------------------------------
   // This functions propagates reconstructed ITS tracks back
   // The clusters must be loaded !
@@ -286,27 +357,20 @@ Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
      if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
      if (esd->GetStatus()&AliESDtrack::kITSout) continue;
 
-     AliITStrackV2 *t=0;
-     try {
-        t=new AliITStrackV2(*esd);
-     } catch (const Char_t *msg) {
-        Warning("PropagateBack",msg);
-        delete t;
-        continue;
-     }
+     AliITStrackV2 *t = new AliITStrackV2(*esd);
 
      ResetTrackToFollow(*t);
 
      // propagete to vertex [SR, GSI 17.02.2003]
      // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
      if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
-       if (fTrackToFollow.PropagateToVertex()) {
+       if (fTrackToFollow.PropagateToVertex(event->GetVertex())) {
           fTrackToFollow.StartTimeIntegral();
        }
        fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
      }
 
-     fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
+     fTrackToFollow.ResetCovariance(10.); fTrackToFollow.ResetClusters();
      if (RefitAt(49.,&fTrackToFollow,t)) {
         if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
           Warning("PropagateBack",
@@ -329,7 +393,7 @@ Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
   return 0;
 }
 
-Int_t AliITStrackerV2::RefitInward(AliESD *event) {
+Int_t AliITStrackerV2::RefitInward(AliESDEvent *event) {
   //--------------------------------------------------------------------
   // This functions refits ITS tracks using the 
   // "inward propagated" TPC tracks
@@ -347,14 +411,7 @@ Int_t AliITStrackerV2::RefitInward(AliESD *event) {
     if (esd->GetStatus()&AliESDtrack::kTPCout)
     if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
 
-    AliITStrackV2 *t=0;
-    try {
-        t=new AliITStrackV2(*esd);
-    } catch (const Char_t *msg) {
-        Warning("RefitInward",msg);
-        delete t;
-        continue;
-    }
+    AliITStrackV2 *t = new AliITStrackV2(*esd);
 
     if (CorrectForDeadZoneMaterial(t)!=0) {
        Warning("RefitInward",
@@ -373,9 +430,9 @@ Int_t AliITStrackerV2::RefitInward(AliESD *event) {
        CookLabel(&fTrackToFollow,0.); //For comparison only
 
        if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe 
+        fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
          AliESDtrack  *esdTrack =fTrackToFollow.GetESDtrack();
-         esdTrack->UpdateTrackParams(&fTrackToFollow,AliESDtrack::kITSrefit);
-         Float_t r[3]={0.,0.,0.};
+         Double_t r[3]={0.,0.,0.};
          Double_t maxD=3.;
         esdTrack->RelateToVertex(event->GetVertex(),GetBz(r),maxD);
          ntrk++;
@@ -422,12 +479,11 @@ void AliITStrackerV2::FollowProlongation() {
     }
 
     //find intersection
-    Double_t x,y,z;  
-    if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
+    Double_t phi,z;  
+    if (!fTrackToFollow.GetPhiZat(r,phi,z)) {
       //Warning("FollowProlongation","failed to estimate track !\n");
       return;
     }
-    Double_t phi=TMath::ATan2(y,x);
 
     Int_t idet=layer.FindDetectorIndex(phi,z);
     if (idet<0) {
@@ -446,26 +502,26 @@ void AliITStrackerV2::FollowProlongation() {
 
     //Select possible prolongations and store the current track estimation
     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
-    Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
-    Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
+    Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
+    Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
     Double_t road=layer.GetRoad();
     if (dz*dy>road*road) {
        Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
        dz=road*scz; dy=road*scy;
     } 
 
-    //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
+    //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
-    if (dz > kMaxRoad) {
+    if (dz > AliITSReconstructor::GetRecoParam()->GetMaxRoad()) {
       //Warning("FollowProlongation","too broad road in Z !\n");
       return;
     }
 
     if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
 
-    //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
+    //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
-    if (dy > kMaxRoad) {
+    if (dy > AliITSReconstructor::GetRecoParam()->GetMaxRoad()) {
       //Warning("FollowProlongation","too broad road in Y !\n");
       return;
     }
@@ -490,7 +546,7 @@ void AliITStrackerV2::FollowProlongation() {
   if (ncl)
   if (ncl >= nclb) {
      Double_t chi2=fTrackToFollow.GetChi2();
-     if (chi2/ncl < kChi2PerCluster) {        
+     if (chi2/ncl < AliITSReconstructor::GetRecoParam()->GetChi2PerCluster()) {        
         if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
            ResetBestTrack();
         }
@@ -508,8 +564,8 @@ Int_t AliITStrackerV2::TakeNextProlongation() {
   AliITSlayer &layer=fgLayers[fI];
   ResetTrackToFollow(fTracks[fI]);
 
-  Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
-  Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
+  Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(fI));
+  Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(fI));
   Double_t road=layer.GetRoad();
   if (dz*dy>road*road) {
      Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
@@ -518,7 +574,7 @@ Int_t AliITStrackerV2::TakeNextProlongation() {
 
   const AliITSRecPoint *c=0; Int_t ci=-1;
   const AliITSRecPoint *cc=0; Int_t cci=-1;
-  Double_t chi2=kMaxChi2;
+  Double_t chi2=AliITSReconstructor::GetRecoParam()->GetMaxChi2();
   while ((c=layer.GetNextCluster(ci))!=0) {
     Int_t idet=c->GetDetectorIndex();
 
@@ -545,16 +601,20 @@ Int_t AliITStrackerV2::TakeNextProlongation() {
 
   if (!cc) return 0;
 
+  {// Take into account the mis-alignment
+    Double_t x = fTrackToFollow.GetX() + cc->GetX();
+    if (!fTrackToFollow.PropagateTo(x,0.,0.)) return 0;
+  }
   if (!fTrackToFollow.Update(cc,chi2,(fI<<28)+cci)) {
      //Warning("TakeNextProlongation","filtering failed !\n");
      return 0;
   }
 
   if (fTrackToFollow.GetNumberOfClusters()>1)
-  if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
+    if (TMath::Abs(fTrackToFollow.GetD(GetX(),GetY()))>4) return 0;
 
   fTrackToFollow.
-    SetSampledEdx(cc->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
+    SetSampledEdx(cc->GetQ(),fI-2); //b.b.
 
   {
   Double_t x0;
@@ -573,35 +633,48 @@ Int_t AliITStrackerV2::TakeNextProlongation() {
 }
 
 
-AliITStrackerV2::AliITSlayer::AliITSlayer() {
+AliITStrackerV2::AliITSlayer::AliITSlayer():
+  fR(0.),
+  fPhiOffset(0.),
+  fNladders(0),
+  fZOffset(0.),
+  fNdetectors(0),
+  fDetectors(0),
+  fNsel(0),
+  fRoad(2*fR*TMath::Sqrt(3.14/1.)) //assuming that there's only one cluster
+{
   //--------------------------------------------------------------------
   //default AliITSlayer constructor
   //--------------------------------------------------------------------
-  fR=0.; fPhiOffset=0.; fZOffset=0.;
-  fNladders=0; fNdetectors=0;
-  fDetectors=0;
   
   for (Int_t i=0; i<kNsector; i++) fN[i]=0;
-  fNsel=0;
-
-  fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
+  for (Int_t i=0; i<AliITSRecoParam::fgkMaxClusterPerLayer; i++){
+    fClusters[i]=0;
+    fIndex[i]=0;
+  }
 }
 
 AliITStrackerV2::AliITSlayer::
-AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
+AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd): 
+  fR(r), 
+  fPhiOffset(p), 
+  fNladders(nl),
+  fZOffset(z),
+  fNdetectors(nd),
+  fDetectors(new AliITSdetector[nl*nd]),
+  fNsel(0),
+  fRoad(2*r*TMath::Sqrt(3.14/1.)) //assuming that there's only one cluster
+{
   //--------------------------------------------------------------------
   //main AliITSlayer constructor
   //--------------------------------------------------------------------
-  fR=r; fPhiOffset=p; fZOffset=z;
-  fNladders=nl; fNdetectors=nd;
-  fDetectors=new AliITSdetector[fNladders*fNdetectors];
 
   for (Int_t i=0; i<kNsector; i++) fN[i]=0;
-  fNsel=0;
 
-  for (Int_t i=0; i<kMaxClusterPerLayer; i++) fClusters[i]=0;
-
-  fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
+  for (Int_t i=0; i<AliITSRecoParam::fgkMaxClusterPerLayer; i++){
+    fClusters[i]=0;
+    fIndex[i]=0;
+  }
 }
 
 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
@@ -756,7 +829,12 @@ AliITStrackerV2::AliITSlayer::FindDetectorIndex(Double_t phi,Double_t z)const {
   //--------------------------------------------------------------------
   //This function finds the detector crossed by the track
   //--------------------------------------------------------------------
-  Double_t dphi=-(phi-fPhiOffset);
+  Double_t dphi;
+  if (fZOffset<0)            // old geometry
+    dphi = -(phi-fPhiOffset);
+  else                       // new geometry
+    dphi = phi-fPhiOffset;
+
   if      (dphi <  0) dphi += 2*TMath::Pi();
   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
@@ -908,9 +986,9 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
   // If "extra"==kTRUE, 
   //    the clusters from overlapped modules get attached to "t" 
   //--------------------------------------------------------------------
-  Int_t index[kMaxLayer];
+  Int_t index[AliITSgeomTGeo::kNLayers];
   Int_t k;
-  for (k=0; k<kMaxLayer; k++) index[k]=-1;
+  for (k=0; k<AliITSgeomTGeo::GetNLayers(); k++) index[k]=-1;
   Int_t nc=c->GetNumberOfClusters();
   for (k=0; k<nc; k++) { 
     Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
@@ -919,10 +997,10 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
 
   Int_t from, to, step;
   if (xx > t->GetX()) {
-      from=0; to=kMaxLayer;
+      from=0; to=AliITSgeomTGeo::GetNLayers();
       step=+1;
   } else {
-      from=kMaxLayer-1; to=-1;
+      from=AliITSgeomTGeo::GetNLayers()-1; to=-1;
       step=-1;
   }
 
@@ -932,13 +1010,21 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
  
      {
      Double_t hI=i-0.5*step; 
-     if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {             
-        Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
-        Double_t d=0.0034, x0=38.6; 
-        if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
-        if (!t->PropagateTo(rs,-step*d,x0)) {
-          return kFALSE;
-        }
+     if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {  
+       Int_t iLay = i-step;
+       Double_t rs = 0.;
+       if(iLay<0 || iLay>= AliITSgeomTGeo::kNLayers){
+        AliError(Form("Invalid layer %d ",iLay));
+        return kFALSE;
+       }
+       else{
+        rs=0.5*(fgLayers[i-step].GetR() + r);
+       }
+       Double_t d=0.0034, x0=38.6; 
+       if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
+       if (!t->PropagateTo(rs,-step*d,x0)) {
+        return kFALSE;
+       }
      }
      }
 
@@ -949,11 +1035,11 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
      }
      //
 
-     Double_t x,y,z;
-     if (!t->GetGlobalXYZat(r,x,y,z)) { 
+     Double_t phi,z;
+     if (!t->GetPhiZat(r,phi,z)) { 
        return kFALSE;
      }
-     Double_t phi=TMath::ATan2(y,x);
+
      Int_t idet=layer.FindDetectorIndex(phi,z);
      if (idet<0) { 
        return kFALSE;
@@ -966,22 +1052,22 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
      t->SetDetectorIndex(idet);
 
      const AliITSRecPoint *cl=0;
-     Double_t maxchi2=kMaxChi2;
+     Double_t maxchi2=AliITSReconstructor::GetRecoParam()->GetMaxChi2();
 
      Int_t idx=index[i];
-     if (idx>0) {
-        const AliITSRecPoint *c=(AliITSRecPoint *)GetCluster(idx); 
-        if (idet != c->GetDetectorIndex()) {
-           idet=c->GetDetectorIndex();
-           const AliITSdetector &det=layer.GetDetector(idet);
-           if (!t->Propagate(det.GetPhi(),det.GetR())) {
+     if (idx>=0) {
+        const AliITSRecPoint *ccc=(AliITSRecPoint *)GetCluster(idx); 
+        if (idet != ccc->GetDetectorIndex()) {
+           idet=ccc->GetDetectorIndex();
+           const AliITSdetector &det2=layer.GetDetector(idet);
+           if (!t->Propagate(det2.GetPhi(),det2.GetR())) {
              return kFALSE;
            }
            t->SetDetectorIndex(idet);
         }
-        Double_t chi2=t->GetPredictedChi2(c);
+        Double_t chi2=t->GetPredictedChi2(ccc);
         if (chi2<maxchi2) { 
-         cl=c; 
+         cl=ccc
          maxchi2=chi2; 
        } else {
          return kFALSE;
@@ -989,10 +1075,13 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
      }
  
      if (cl) {
+       // Take into account the mis-alignment
+       Double_t x=t->GetX()+cl->GetX();
+       if (!t->PropagateTo(x,0.,0.)) return kFALSE;
        if (!t->Update(cl,maxchi2,idx)) {
           return kFALSE;
        }
-       t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
+       t->SetSampledEdx(cl->GetQ(),i-2);
      }
 
      {
@@ -1003,9 +1092,9 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
                  
      if (extra) { //search for extra clusters
         AliITStrackV2 tmp(*t);
-        Double_t dz=4*TMath::Sqrt(tmp.GetSigmaZ2()+kSigmaZ2[i]);
+        Double_t dz=4*TMath::Sqrt(tmp.GetSigmaZ2()+AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
         if (dz < 0.5*TMath::Abs(tmp.GetTgl())) dz=0.5*TMath::Abs(tmp.GetTgl());
-        Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
+        Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
         if (dy < 0.5*TMath::Abs(tmp.GetSnp())) dy=0.5*TMath::Abs(tmp.GetSnp());
         Double_t zmin=t->GetZ() - dz;
         Double_t zmax=t->GetZ() + dz;
@@ -1013,19 +1102,20 @@ Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
         Double_t ymax=t->GetY() + phi*r + dy;
         layer.SelectClusters(zmin,zmax,ymin,ymax);
 
-        const AliITSRecPoint *c=0; Int_t ci=-1,cci=-1;
-        Double_t maxchi2=1000.*kMaxChi2, tolerance=0.1;
-        while ((c=layer.GetNextCluster(ci))!=0) {
-           if (idet == c->GetDetectorIndex()) continue;
+        const AliITSRecPoint *cx=0; Int_t ci=-1,cci=-1;
+        maxchi2=1000.*AliITSReconstructor::GetRecoParam()->GetMaxChi2();
+       Double_t tolerance=0.1;
+        while ((cx=layer.GetNextCluster(ci))!=0) {
+           if (idet == cx->GetDetectorIndex()) continue;
 
-          const AliITSdetector &det=layer.GetDetector(c->GetDetectorIndex());
+          const AliITSdetector &detx=layer.GetDetector(cx->GetDetectorIndex());
 
-          if (!tmp.Propagate(det.GetPhi(),det.GetR())) continue;
+          if (!tmp.Propagate(detx.GetPhi(),detx.GetR())) continue;
            
-          if (TMath::Abs(tmp.GetZ() - c->GetZ()) > tolerance) continue;
-           if (TMath::Abs(tmp.GetY() - c->GetY()) > tolerance) continue;
+          if (TMath::Abs(tmp.GetZ() - cx->GetZ()) > tolerance) continue;
+           if (TMath::Abs(tmp.GetY() - cx->GetY()) > tolerance) continue;
 
-           Double_t chi2=tmp.GetPredictedChi2(c);
+           Double_t chi2=tmp.GetPredictedChi2(cx);
            if (chi2<maxchi2) { maxchi2=chi2; cci=ci; }
         }
         if (cci>=0) t->SetExtraCluster(i,(i<<28)+cci);