]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PMD/AliPMDUtility.cxx
Fixed the DA so it picks up the pedestals from the right sample
[u/mrichter/AliRoot.git] / PMD / AliPMDUtility.cxx
index 402775df6d62b26447a9e8c18d40348df46f461a..ccb63332d3e2c28e214816ceba14b3d1b82062ee 100644 (file)
 //-----------------------------------------------------//
 
 #include "Riostream.h"
-#include "AliPMDUtility.h"
 #include "TMath.h"
+#include "TText.h"
+#include "TLine.h"
+
 #include <stdio.h>
 #include <math.h>
 
+#include "AliPMDUtility.h"
+
 
 ClassImp(AliPMDUtility)
 
-AliPMDUtility::AliPMDUtility()
+AliPMDUtility::AliPMDUtility():
+  fPx(0.),
+  fPy(0.),
+  fPz(0.),
+  fTheta(0.),
+  fEta(0.),
+  fPhi(0.),
+  fWriteModule(1)
 {
   // Default constructor
-  fPx    = 0.;
-  fPy    = 0.;
-  fPz    = 0.;
-  fTheta = 0.;
-  fEta   = 0.;
-  fPhi   = 0.;
 }
 
-AliPMDUtility::AliPMDUtility(Float_t px, Float_t py, Float_t pz)
+AliPMDUtility::AliPMDUtility(Float_t px, Float_t py, Float_t pz):
+  fPx(px),
+  fPy(py),
+  fPz(pz),
+  fTheta(0.),
+  fEta(0.),
+  fPhi(0.),
+  fWriteModule(1)
 {
   // Constructor
-  fPx = px;
-  fPy = py;
-  fPz = pz;
-  fTheta = 0.;
-  fEta   = 0.;
-  fPhi   = 0.;
 }
-
+AliPMDUtility::AliPMDUtility(const AliPMDUtility &pmdutil):
+  fPx(pmdutil.fPx),
+  fPy(pmdutil.fPy),
+  fPz(pmdutil.fPz),
+  fTheta(pmdutil.fTheta),
+  fEta(pmdutil.fEta),
+  fPhi(pmdutil.fPhi),
+  fWriteModule(pmdutil.fWriteModule)
+{
+  // copy constructor
+}
+AliPMDUtility & AliPMDUtility::operator=(const AliPMDUtility &pmdutil)
+{
+  // assignment operator
+  if(this != &pmdutil)
+    {
+      fPx = pmdutil.fPx;
+      fPy = pmdutil.fPy;
+      fPz = pmdutil.fPz;
+      fTheta = pmdutil.fTheta;
+      fEta = pmdutil.fEta;
+      fPhi = pmdutil.fPhi;
+      fWriteModule = pmdutil.fWriteModule;
+    }
+  return *this;
+}
 AliPMDUtility::~AliPMDUtility()
 {
   // Default destructor
@@ -236,6 +267,143 @@ void AliPMDUtility::RectGeomCellPos(Int_t ism, Float_t xpad, Float_t ypad, Float
     }
 
 }
+// -------------------------------------------------------- //
+
+void AliPMDUtility::GenerateBoundaryPoints(Int_t ism, Float_t &x1ism, 
+                                          Float_t &y1ism, Float_t &x2ism,
+                                          Float_t &y2ism)
+{
+  // Generate bounding-box.
+
+
+    Float_t xism = 0, yism = 0;
+    Float_t dxism = 0., dyism = 0.;
+
+    const Float_t kRad     = 0.25;
+    const Float_t kSqRoot3 = 1.732050808;
+    const Float_t kDia     = 0.50;
+
+
+  const Double_t kXcorner[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
+    };
+
+
+  const Double_t kYcorner[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
+    };
+
+
+  if (ism > 23) ism -= 24;
+
+
+  if (ism < 6)
+    {
+      xism  = kXcorner[ism] + kRad;
+      yism  = kYcorner[ism] + kRad;
+      dxism = -kRad*kSqRoot3*48.;
+      dyism = -kDia*96. - kRad;
+  }
+  if (ism >= 6 && ism < 12)
+    {
+      xism  = kXcorner[ism] - kRad;
+      yism  = kYcorner[ism] - kRad;
+      dxism = kRad*kSqRoot3*48.;
+      dyism = kDia*96. + kRad;
+  }
+  if (ism >= 12 && ism < 18)
+    {
+      xism  = kXcorner[ism] + kRad;
+      yism  = kYcorner[ism] + kRad;
+      dxism = -kRad*kSqRoot3*96.;
+      dyism = -kDia*48. - kRad;
+  }
+  if (ism >= 18 && ism < 24)
+    {
+      xism  = kXcorner[ism] - kRad;
+      yism  = kYcorner[ism] - kRad;
+      dxism = kRad*kSqRoot3*96.;
+      dyism = kDia*48. + kRad;
+  }
+
+  x1ism = xism;
+  x2ism = xism + dxism;
+  y1ism = yism;
+  y2ism = yism + dyism;
+
+}
+// ------------------------------------------------------------------- //
+
+void AliPMDUtility::DrawPMDModule()
+{
+
+    Float_t x1ism, x2ism, y1ism, y2ism;
+    Float_t deltaX, deltaY;
+    
+    //TH2F *h2 = new TH2F("h2","Y vs. X",200,-100.,100.,200,-100.,100.);
+    //h2->Draw();
+
+    TLine t;
+    t.SetLineColor(2);
+
+    TText tt;
+    tt.SetTextColor(4);
+
+    Char_t smnumber[10];
+
+    for(Int_t ism=0; ism < 24; ism++)
+    {
+       GenerateBoundaryPoints(ism, x1ism, y1ism, x2ism, y2ism);
+       deltaX = (x2ism - x1ism)/2.;
+       deltaY = (y2ism - y1ism)/2.;
+       if (fWriteModule == 1)
+       {
+           sprintf(smnumber,"%d",ism);
+           tt.DrawText(x1ism+deltaX,y1ism+deltaY,smnumber);
+       }
+       t.DrawLine(x1ism, y1ism, x1ism, y2ism);
+       t.DrawLine(x1ism, y1ism, x2ism, y1ism);
+       t.DrawLine(x2ism, y1ism, x2ism, y2ism);
+       t.DrawLine(x1ism, y2ism, x2ism, y2ism);
+    }
+
+}
+
+// ------------------------------------------------------------------- //
+
+
+void AliPMDUtility::ApplyVertexCorrection(Float_t vertex[], Float_t xpos,
+                                         Float_t ypos, Float_t zpos)
+{
+  // Not implemented
+  fPx = xpos - vertex[0];
+  fPy = ypos - vertex[1];
+  fPz = zpos - vertex[2];
+}
+void AliPMDUtility::ApplyAlignment()
+{
+  // Not implemented
+}
 
 void AliPMDUtility::SetPxPyPz(Float_t px, Float_t py, Float_t pz)
 {
@@ -250,6 +418,10 @@ void AliPMDUtility::SetXYZ(Float_t xpos, Float_t ypos, Float_t zpos)
   fPy = ypos;
   fPz = zpos;
 }
+void AliPMDUtility::SetWriteModule(Int_t wrmod)
+{
+    fWriteModule = wrmod;
+}
 void AliPMDUtility::CalculateEta()
 {
   Float_t rpxpy, theta, eta;
@@ -295,7 +467,7 @@ void AliPMDUtility::CalculateEtaPhi()
   theta = TMath::ATan2(rpxpy,fPz);
   eta   = -TMath::Log(TMath::Tan(0.5*theta));
   
-  if(fPx==0)
+  if(fPx == 0)
     {
       if(fPy>0) phi = 90.;
       if(fPy<0) phi = 270.;
@@ -317,6 +489,22 @@ void AliPMDUtility::CalculateEtaPhi()
   fEta   = eta;
   fPhi   = phi;
 }
+void AliPMDUtility::CalculateXY(Float_t eta, Float_t phi, Float_t zpos)
+{
+  // Not implemented
+
+  //  eta   = -TMath::Log(TMath::Tan(0.5*theta));
+
+  Float_t xpos = 0., ypos = 0.;
+
+  //  Float_t theta = 2.0*TMath::ATan(TMath::Log(-eta));
+
+  fEta = eta;
+  fPhi = phi;
+  fPx  = xpos;
+  fPy  = ypos;
+  fPz  = zpos;
+}
 Float_t AliPMDUtility::GetTheta() const
 {
   return fTheta;
@@ -329,4 +517,17 @@ Float_t AliPMDUtility::GetPhi() const
 {
   return fPhi;
 }
+Float_t AliPMDUtility::GetX() const
+{
+  return fPx;
+}
+Float_t AliPMDUtility::GetY() const
+{
+  return fPy;
+}
+Float_t AliPMDUtility::GetZ() const
+{
+  return fPz;
+}
+