]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - MFT/AliMFTPlane.cxx
Error changed to info
[u/mrichter/AliRoot.git] / MFT / AliMFTPlane.cxx
index 0e9cc07a7bd92e4179628cb45ef126b4d283d7f3..419744775a9cd4f1030505c2fb386a482ab52a06 100644 (file)
@@ -63,7 +63,8 @@ AliMFTPlane::AliMFTPlane():
   fEquivalentSiliconBeforeBack(0),
   fActiveElements(0),
   fReadoutElements(0),
-  fSupportElements(0)
+  fSupportElements(0),
+  fHasPixelRectangularPatternAlongY(kFALSE)
 {
 
   // default constructor
@@ -91,7 +92,8 @@ AliMFTPlane::AliMFTPlane(const Char_t *name, const Char_t *title):
   fEquivalentSiliconBeforeBack(0),
   fActiveElements(new TClonesArray("THnSparseC")),
   fReadoutElements(new TClonesArray("THnSparseC")),
-  fSupportElements(new TClonesArray("THnSparseC"))
+  fSupportElements(new TClonesArray("THnSparseC")),
+  fHasPixelRectangularPatternAlongY(kFALSE)
 {
 
   // constructor
@@ -119,7 +121,8 @@ AliMFTPlane::AliMFTPlane(const AliMFTPlane& plane):
   fEquivalentSiliconBeforeBack(plane.fEquivalentSiliconBeforeBack),
   fActiveElements(new TClonesArray("THnSparseC")),
   fReadoutElements(new TClonesArray("THnSparseC")),
-  fSupportElements(new TClonesArray("THnSparseC"))
+  fSupportElements(new TClonesArray("THnSparseC")),
+  fHasPixelRectangularPatternAlongY(plane.fHasPixelRectangularPatternAlongY)
 {
 
   // copy constructor
@@ -142,24 +145,25 @@ AliMFTPlane& AliMFTPlane::operator=(const AliMFTPlane& plane) {
     // base class assignement
     TNamed::operator=(plane);
     
-    fPlaneNumber                  = plane.fPlaneNumber;
-    fZCenter                      = plane.fZCenter; 
-    fRMinSupport                  = plane.fRMinSupport; 
-    fRMax                         = plane.fRMax;
-    fRMaxSupport                  = plane.fRMaxSupport;
-    fPixelSizeX                   = plane.fPixelSizeX;
-    fPixelSizeY                   = plane.fPixelSizeY; 
-    fThicknessActive              = plane.fThicknessActive; 
-    fThicknessSupport             = plane.fThicknessSupport; 
-    fThicknessReadout             = plane.fThicknessReadout;
-    fZCenterActiveFront           = plane.fZCenterActiveFront;
-    fZCenterActiveBack            = plane.fZCenterActiveBack;
-    fEquivalentSilicon            = plane.fEquivalentSilicon;
-    fEquivalentSiliconBeforeFront = plane.fEquivalentSiliconBeforeFront;
-    fEquivalentSiliconBeforeBack  = plane.fEquivalentSiliconBeforeBack;
-    *fActiveElements              = *plane.fActiveElements;
-    *fReadoutElements             = *plane.fReadoutElements;
-    *fSupportElements             = *plane.fSupportElements;
+    fPlaneNumber                      = plane.fPlaneNumber;
+    fZCenter                          = plane.fZCenter; 
+    fRMinSupport                      = plane.fRMinSupport; 
+    fRMax                             = plane.fRMax;
+    fRMaxSupport                      = plane.fRMaxSupport;
+    fPixelSizeX                       = plane.fPixelSizeX;
+    fPixelSizeY                       = plane.fPixelSizeY; 
+    fThicknessActive                  = plane.fThicknessActive; 
+    fThicknessSupport                 = plane.fThicknessSupport; 
+    fThicknessReadout                 = plane.fThicknessReadout;
+    fZCenterActiveFront               = plane.fZCenterActiveFront;
+    fZCenterActiveBack                = plane.fZCenterActiveBack;
+    fEquivalentSilicon                = plane.fEquivalentSilicon;
+    fEquivalentSiliconBeforeFront     = plane.fEquivalentSiliconBeforeFront;
+    fEquivalentSiliconBeforeBack      = plane.fEquivalentSiliconBeforeBack;
+    *fActiveElements                  = *plane.fActiveElements;
+    *fReadoutElements                 = *plane.fReadoutElements;
+    *fSupportElements                 = *plane.fSupportElements;
+    fHasPixelRectangularPatternAlongY = plane.fHasPixelRectangularPatternAlongY;
   }
 
   return *this;
@@ -176,7 +180,8 @@ Bool_t AliMFTPlane::Init(Int_t    planeNumber,
                         Double_t pixelSizeY, 
                         Double_t thicknessActive, 
                         Double_t thicknessSupport, 
-                        Double_t thicknessReadout) {
+                        Double_t thicknessReadout,
+                        Bool_t   hasPixelRectangularPatternAlongY) {
 
   AliDebug(1, Form("Initializing Plane Structure for Plane %s", GetName()));
 
@@ -190,6 +195,8 @@ Bool_t AliMFTPlane::Init(Int_t    planeNumber,
   fThicknessSupport = thicknessSupport;
   fThicknessReadout = thicknessReadout;
 
+  fHasPixelRectangularPatternAlongY = hasPixelRectangularPatternAlongY;
+
   fZCenterActiveFront = fZCenter - 0.5*fThicknessSupport - 0.5*fThicknessActive;
   fZCenterActiveBack  = fZCenter + 0.5*fThicknessSupport + 0.5*fThicknessActive;
 
@@ -403,8 +410,6 @@ void AliMFTPlane::DrawPlane(Option_t *opt) {
     TCanvas *cnv = new TCanvas("cnv", GetName(), 900, 900);
     cnv->Draw();
 
-    //    printf("Created Canvas\n");
-
     TH2D *h = new TH2D("tmp", GetName(), 
                       1, 1.1*GetSupportElement(0)->GetAxis(0)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(0)->GetXmax(), 
                       1, 1.1*GetSupportElement(0)->GetAxis(1)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(1)->GetXmax());
@@ -420,10 +425,7 @@ void AliMFTPlane::DrawPlane(Option_t *opt) {
     supportExt -> Draw("same");
     supportInt -> Draw("same");
 
-    //    printf("Created Ellipses\n");
-
     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
-      //      printf("Active element %d\n", iEl);
       if (!IsFront(GetActiveElement(iEl))) continue;
       TPave *pave = new TPave(GetActiveElement(iEl)->GetAxis(0)->GetXmin(), 
                              GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
@@ -434,7 +436,6 @@ void AliMFTPlane::DrawPlane(Option_t *opt) {
     }
 
     for (Int_t iEl=0; iEl<GetNReadoutElements(); iEl++) {
-      //      printf("Readout element %d\n", iEl);
       if (!IsFront(GetReadoutElement(iEl))) continue;
       TPave *pave = new TPave(GetReadoutElement(iEl)->GetAxis(0)->GetXmin(), 
                              GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
@@ -619,3 +620,28 @@ void AliMFTPlane::DrawPlane(Option_t *opt) {
 
 //====================================================================================================================================================
 
+Int_t AliMFTPlane::GetNumberOfChips(Option_t *opt) {
+
+  Int_t nChips = 0;
+
+  if (!strcmp(opt, "front")) {
+    for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
+      if (!IsFront(GetActiveElement(iEl))) continue;
+      Double_t length = GetActiveElement(iEl)->GetAxis(0)->GetXmax() - GetActiveElement(iEl)->GetAxis(0)->GetXmin();
+      nChips += Int_t (length/AliMFTConstants::fWidthChip) + 1;
+    }
+  }
+
+  else if (!strcmp(opt, "back")) {
+    for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
+      if (IsFront(GetActiveElement(iEl))) continue;
+      Double_t length = GetActiveElement(iEl)->GetAxis(0)->GetXmax() - GetActiveElement(iEl)->GetAxis(0)->GetXmin();
+      nChips += Int_t (length/AliMFTConstants::fWidthChip) + 1;
+    }
+  }
+
+  return nChips;
+
+}
+
+//====================================================================================================================================================