]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PMD/AliPMDUtility.cxx
Coding rule violation corrected.
[u/mrichter/AliRoot.git] / PMD / AliPMDUtility.cxx
index 7f1ba5aa14f3fe78212b379f0b81bab84f4cdc0c..b16c3ce223f6358c726585d51b9e260c53aa62d5 100644 (file)
@@ -31,7 +31,6 @@
 
 #include "AliPMDUtility.h"
 
-
 ClassImp(AliPMDUtility)
 
 AliPMDUtility::AliPMDUtility():
@@ -44,6 +43,14 @@ AliPMDUtility::AliPMDUtility():
   fWriteModule(1)
 {
   // Default constructor
+  for (Int_t i = 0; i < 4; i++)
+    {
+      for (Int_t j = 0; j < 3; j++)
+       {
+         fSecTr[i][j] = 0.;
+       }
+    }
+
 }
 
 AliPMDUtility::AliPMDUtility(Float_t px, Float_t py, Float_t pz):
@@ -56,8 +63,17 @@ AliPMDUtility::AliPMDUtility(Float_t px, Float_t py, Float_t pz):
   fWriteModule(1)
 {
   // Constructor
+  for (Int_t i = 0; i < 4; i++)
+    {
+      for (Int_t j = 0; j < 3; j++)
+       {
+         fSecTr[i][j] = 0.;
+       }
+    }
+
 }
 AliPMDUtility::AliPMDUtility(const AliPMDUtility &pmdutil):
+  TObject(pmdutil),
   fPx(pmdutil.fPx),
   fPy(pmdutil.fPy),
   fPz(pmdutil.fPz),
@@ -67,6 +83,14 @@ AliPMDUtility::AliPMDUtility(const AliPMDUtility &pmdutil):
   fWriteModule(pmdutil.fWriteModule)
 {
   // copy constructor
+    for (Int_t i = 0; i < 4; i++)
+    {
+      for (Int_t j = 0; j < 3; j++)
+       {
+         fSecTr[i][j] = pmdutil.fSecTr[i][j];
+       }
+    }
+
 }
 AliPMDUtility & AliPMDUtility::operator=(const AliPMDUtility &pmdutil)
 {
@@ -80,6 +104,14 @@ AliPMDUtility & AliPMDUtility::operator=(const AliPMDUtility &pmdutil)
       fEta = pmdutil.fEta;
       fPhi = pmdutil.fPhi;
       fWriteModule = pmdutil.fWriteModule;
+      for (Int_t i = 0; i < 4; i++)
+       {
+         for (Int_t j = 0; j < 3; j++)
+           {
+             fSecTr[i][j] = pmdutil.fSecTr[i][j];
+           }
+       }
+
     }
   return *this;
 }
@@ -174,9 +206,30 @@ void AliPMDUtility::RectGeomCellPos(Int_t ism, Int_t xpad, Int_t ypad, Float_t &
       ypos = ycorner[ism] + (Float_t) xpad*kCellRadius*2.0 + shift;
       xpos = xcorner[ism] + (Float_t) ypad*kSqroot3*kCellRadius;
     }
+  // Apply the alignment here to the x, y values
+  if(ism < 6)
+    {
+      xpos += fSecTr[0][0];
+      ypos += fSecTr[0][1];
+    }
+  else if(ism >= 6 && ism < 12)
+    {
+      xpos += fSecTr[1][0];
+      ypos += fSecTr[1][1];
+    }
+  else if(ism >=12 && ism < 18)
+    {
+      xpos += fSecTr[2][0];
+      ypos += fSecTr[2][1];
+    }
+  else if(ism >= 18 && ism < 24)
+    {
+      xpos += fSecTr[3][0];
+      ypos += fSecTr[3][1];
+    }
 
 }
-
+// ---------------------------------------------------------- 
 void AliPMDUtility::RectGeomCellPos(Int_t ism, Float_t xpad, Float_t ypad, Float_t &xpos, Float_t &ypos)
 {
   // If the xpad and ypad inputs are float, then 0.5 is added to it
@@ -266,6 +319,151 @@ void AliPMDUtility::RectGeomCellPos(Int_t ism, Float_t xpad, Float_t ypad, Float
       xpos = xcorner[ism] + ypad*kSqroot3*kCellRadius;
     }
 
+  // Apply the alignment here to the x, y values
+  if(ism < 6)
+    {
+      xpos += fSecTr[0][0];
+      ypos += fSecTr[0][1];
+    }
+  else if(ism >= 6 && ism < 12)
+    {
+      xpos += fSecTr[1][0];
+      ypos += fSecTr[1][1];
+    }
+  else if(ism >=12 && ism < 18)
+    {
+      xpos += fSecTr[2][0];
+      ypos += fSecTr[2][1];
+    }
+  else if(ism >= 18 && ism < 24)
+    {
+      xpos += fSecTr[3][0];
+      ypos += fSecTr[3][1];
+    }
+
+}
+
+// -------------------------------------------------------- //
+
+void AliPMDUtility::RectGeomCellPos(Int_t ism, Float_t xpad,
+                                   Float_t ypad, Float_t &xpos,
+                                   Float_t &ypos, Float_t & zpos)
+{
+  // If the xpad and ypad inputs are float, then 0.5 is added to it
+  // to find the layer which is shifted.
+  // This routine finds the cell eta,phi for the new PMD rectangular 
+  // geometry in ALICE
+  // Authors : Bedanga Mohanty and Dipak Mishra - 29.4.2003
+  // modified by B. K. Nnadi for change of coordinate sys
+  //
+  // SMA  ---> Supermodule Type A           ( SM - 0)
+  // SMAR ---> Supermodule Type A ROTATED   ( SM - 1)
+  // SMB  ---> Supermodule Type B           ( SM - 2)
+  // SMBR ---> Supermodule Type B ROTATED   ( SM - 3)
+  //
+  // ism   : Serial Module number from 0 to 23 for each plane
+
+  // Corner positions (x,y) of the 24 unit moudles in ALICE PMD
+
+  double xcorner[24] =
+    {
+      74.8833,  53.0045, 31.1255,    //Type-A
+      74.8833,  53.0045, 31.1255,    //Type-A
+      -74.8833, -53.0044, -31.1255,  //Type-AR
+      -74.8833, -53.0044, -31.1255,  //Type-AR
+      8.9165, -33.7471,            //Type-B
+      8.9165, -33.7471,            //Type-B
+      8.9165, -33.7471,            //Type-B
+      -8.9165, 33.7471,            //Type-BR
+      -8.9165, 33.7471,            //Type-BR
+      -8.9165, 33.7471,            //Type-BR
+    };
+
+  
+
+  double ycorner[24] =
+    {
+      86.225,  86.225,  86.225,      //Type-A
+      37.075,  37.075,  37.075,      //Type-A
+      -86.225, -86.225, -86.225,     //Type-AR
+      -37.075, -37.075, -37.075,     //Type-AR
+      86.225,  86.225,               //Type-B
+      61.075,  61.075,               //Type-B
+      35.925,  35.925,               //Type-B
+      -86.225, -86.225,              //Type-BR
+      -61.075, -61.075,              //Type-BR
+      -35.925, -35.925               //Type-BR
+    };
+
+
+  const Float_t kSqroot3    = 1.73205;  // sqrt(3.);
+  const Float_t kCellRadius = 0.25;
+  
+  //
+  //Every even row of cells is shifted and placed
+  //in geant so this condition
+  //
+  Float_t cellRadius = 0.25;
+  Float_t shift = 0.0;
+  Int_t iirow = (Int_t) (xpad+0.5);
+  if(iirow%2 == 0)
+    {
+      shift = -cellRadius/2.0;
+    }
+  else
+    {
+      shift = 0.0;
+    }
+
+  if(ism < 6)
+    {
+      ypos = ycorner[ism] - xpad*kCellRadius*2.0 + shift;
+      xpos = xcorner[ism] - ypad*kSqroot3*kCellRadius;
+    }
+  else if(ism >=6 && ism < 12)
+    {
+      ypos = ycorner[ism] + xpad*kCellRadius*2.0 + shift;
+      xpos = xcorner[ism] + ypad*kSqroot3*kCellRadius;
+    }
+  else if(ism >= 12 && ism < 18)
+    {
+      ypos = ycorner[ism] - xpad*kCellRadius*2.0 + shift;
+      xpos = xcorner[ism] - ypad*kSqroot3*kCellRadius;
+    }
+  else if(ism >= 18 && ism < 24)
+    {
+      ypos = ycorner[ism] + xpad*kCellRadius*2.0 + shift;
+      xpos = xcorner[ism] + ypad*kSqroot3*kCellRadius;
+    }
+
+  // Apply the alignment here to the x, y, and z values
+  if(ism < 6)
+    {
+      xpos += fSecTr[0][0];
+      ypos += fSecTr[0][1];
+      zpos += fSecTr[0][2];
+    }
+  else if(ism >= 6 && ism < 12)
+    {
+      xpos += fSecTr[1][0];
+      ypos += fSecTr[1][1];
+      zpos += fSecTr[1][2];
+    }
+  else if(ism >=12 && ism < 18)
+    {
+      xpos += fSecTr[2][0];
+      ypos += fSecTr[2][1];
+      zpos += fSecTr[2][2];
+    }
+  else if(ism >= 18 && ism < 24)
+    {
+      xpos += fSecTr[3][0];
+      ypos += fSecTr[3][1];
+      zpos += fSecTr[3][2];
+    }
+
+
+
 }
 // -------------------------------------------------------- //
 
@@ -357,8 +555,8 @@ void AliPMDUtility::GenerateBoundaryPoints(Int_t ism, Float_t &x1ism,
 void AliPMDUtility::DrawPMDModule(Int_t idet)
 {
 
-    Float_t x1ism, x2ism, y1ism, y2ism;
-    Float_t deltaX, deltaY;
+    Float_t x1ism = 0., x2ism = 0., y1ism = 0., y2ism = 0.;
+    Float_t deltaX = 0., deltaY = 0.;
     
     //TH2F *h2 = new TH2F("h2","Y vs. X",200,-100.,100.,200,-100.,100.);
     //h2->Draw();
@@ -380,11 +578,11 @@ void AliPMDUtility::DrawPMDModule(Int_t idet)
        {
          if(idet == 0)
            {
-             sprintf(smnumber,"%d",ism);
+             snprintf(smnumber,10,"%d",ism);
            }
          else if (idet == 1)
            {
-             sprintf(smnumber,"%d",24+ism);
+             snprintf(smnumber,10,"%d",24+ism);
            }
            tt.DrawText(x1ism+deltaX,y1ism+deltaY,smnumber);
        }
@@ -407,9 +605,17 @@ void AliPMDUtility::ApplyVertexCorrection(Float_t vertex[], Float_t xpos,
   fPy = ypos - vertex[1];
   fPz = zpos - vertex[2];
 }
-void AliPMDUtility::ApplyAlignment()
+void AliPMDUtility::ApplyAlignment(Double_t sectr[][3])
 {
-  // Not implemented
+  // Get the alignment stuff here
+
+  for (Int_t isector=0; isector<4; isector++)
+    {
+      for(Int_t ixyz=0; ixyz < 3; ixyz++)
+       {
+         fSecTr[isector][ixyz] = (Float_t) sectr[isector][ixyz];
+       }
+    }
 }
 
 void AliPMDUtility::SetPxPyPz(Float_t px, Float_t py, Float_t pz)
@@ -431,17 +637,15 @@ void AliPMDUtility::SetWriteModule(Int_t wrmod)
 }
 void AliPMDUtility::CalculateEta()
 {
-  Float_t rpxpy, theta, eta;
-
-  rpxpy  = TMath::Sqrt(fPx*fPx + fPy*fPy);
-  theta  = TMath::ATan2(rpxpy,fPz);
-  eta    = -TMath::Log(TMath::Tan(0.5*theta));
+  Float_t rpxpy  = TMath::Sqrt(fPx*fPx + fPy*fPy);
+  Float_t theta  = TMath::ATan2(rpxpy,fPz);
+  Float_t eta    = -TMath::Log(TMath::Tan(0.5*theta));
   fTheta = theta;
   fEta   = eta;
 }
 void AliPMDUtility::CalculatePhi()
 {
-  Float_t pybypx, phi = 0., phi1;
+  Float_t pybypx = 0., phi = 0., phi1 = 0.;
 
   if(fPx==0)
     {
@@ -467,12 +671,11 @@ void AliPMDUtility::CalculatePhi()
 }
 void AliPMDUtility::CalculateEtaPhi()
 {
-  Float_t rpxpy, theta, eta;
-  Float_t pybypx, phi = 0., phi1;
+  Float_t pybypx = 0., phi = 0., phi1 = 0.;
 
-  rpxpy = TMath::Sqrt(fPx*fPx + fPy*fPy);
-  theta = TMath::ATan2(rpxpy,fPz);
-  eta   = -TMath::Log(TMath::Tan(0.5*theta));
+  Float_t rpxpy = TMath::Sqrt(fPx*fPx + fPy*fPy);
+  Float_t theta = TMath::ATan2(rpxpy,fPz);
+  Float_t eta   = -TMath::Log(TMath::Tan(0.5*theta));
   
   if(fPx == 0)
     {
@@ -536,5 +739,4 @@ Float_t AliPMDUtility::GetZ() const
 {
   return fPz;
 }
-
-
+//--------------------------------------------------------------------//