]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add haff chamber status to simulation part (Julian)
authorcblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 20 Sep 2012 10:02:46 +0000 (10:02 +0000)
committercblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 20 Sep 2012 10:02:46 +0000 (10:02 +0000)
TRD/AliTRDarrayADC.cxx
TRD/AliTRDdigitizer.cxx
TRD/Macros/AliTRDanadigits.C [new file with mode: 0644]

index c8c48ea7ad2a3dbb279c3864ddb2fa53a83142ca..0a8b3e4a9d474124e7151e58c66e7a77fe19bc9a 100644 (file)
@@ -167,7 +167,7 @@ Short_t AliTRDarrayADC::GetDataBits(Int_t row, Int_t col, Int_t time) const
   // Get the ADC value for a given position: row, col, time
   // Taking bit masking into account
   //
-  // Adapted from code of the class AliTRDdataArrayDigits 
+  // Adapted from code of the class AliTRDclusterizer 
   //
 
   Short_t tempval = GetData(row,col,time);
@@ -193,7 +193,7 @@ UChar_t AliTRDarrayADC::GetPadStatus(Int_t row, Int_t col, Int_t time) const
   //               Bridged Right Masking    8
   //               Not Connected Masking Digits
   //
-  // Adapted from code of the class AliTRDdataArrayDigits
+  // Adapted from code of the class AliTRDclusterizer
   //
 
   UChar_t padstatus = 0;
@@ -234,7 +234,7 @@ void AliTRDarrayADC::SetPadStatus(Int_t row, Int_t col, Int_t time, UChar_t stat
   //               Bridged Left masking:    Bit 11(0), Bit 12(1)
   //               Bridged Right masking:   Bit 11(1), Bit 12(1)
   // 
-  // Adapted from code of the class AliTRDdataArrayDigits
+  // Adapted from code of the class AliTRDclusterizer
   //
 
   Short_t signal = GetData(row,col,time);
@@ -280,7 +280,7 @@ Bool_t AliTRDarrayADC::IsPadCorrupted(Int_t row, Int_t col, Int_t time)
   // 
   // Checks if the pad has any masking as corrupted (Bit 10 in signal set)
   // 
-  // Adapted from code of the class AliTRDdataArrayDigits
+  // Adapted from code of the class AliTRDclusterizer
   //
 
   Short_t signal = GetData(row,col,time);
index b1f4dc3e35b6484157c6547fb1b76e969b137b96..931e34b6d605ce761d92b6941c525d7ba34153f6 100644 (file)
@@ -665,7 +665,7 @@ Bool_t AliTRDdigitizer::MakeDigits()
         (nhit[det] > 0)) {
 
       signals = new AliTRDarraySignal();
-         
+
       // Convert the hits of the current detector to detector signals
       if (!ConvertHits(det,hits[det],nhit[det],signals)) {
        AliError(Form("Conversion of hits failed for detector=%d",det));
@@ -1032,7 +1032,7 @@ Bool_t AliTRDdigitizer::ConvertHits(Int_t det
       colE       = padPlane->GetPadColNumber(locC+offsetTilt);
       if (colE < 0) continue;         
       colOffset  = padPlane->GetPadColOffset(colE,locC+offsetTilt);
-         
+
       // Also re-retrieve drift velocity because col and row may have changed
       driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
       Float_t t0    = calT0DetValue     + calT0ROC->GetValue(colE,rowE);
@@ -1190,7 +1190,7 @@ Bool_t AliTRDdigitizer::ConvertSignals(Int_t det, AliTRDarraySignal *signals)
   }
 
   // Compress the arrays
-  CompressOutputArrays(det);   
+  CompressOutputArrays(det);
 
   return kTRUE;
 
@@ -1291,6 +1291,13 @@ Bool_t AliTRDdigitizer::Signal2ADC(Int_t det, AliTRDarraySignal *signals)
   for (row  = 0; row  <  nRowMax; row++ ) {
     for (col  = 0; col  <  nColMax; col++ ) {
 
+      // halfchamber masking
+      Int_t iMcm            = (Int_t)(col/18); // current group of 18 col pads
+      Int_t halfchamberside = (iMcm>3 ? 1 : 0); // 0=Aside, 1=Bside
+      // Halfchambers that are switched off, masked by calibration
+      if (calibration->IsHalfChamberNoData(det, halfchamberside)) 
+       continue;
+
       // Check whether pad is masked
       // Bridged pads are not considered yet!!!
       if (calibration->IsPadMasked(det,col,row) || 
@@ -1383,12 +1390,20 @@ Bool_t AliTRDdigitizer::Signal2SDigits(Int_t det, AliTRDarraySignal *signals)
   // Create the sdigits for this chamber
   for (row  = 0; row  <  nRowMax; row++ ) {
     for (col  = 0; col  <  nColMax; col++ ) {
-      for (time = 0; time < nTimeTotal; time++) {         
+
+      // halfchamber masking
+      Int_t iMcm            = (Int_t)(col/18); // current group of 18 col pads
+      Int_t halfchamberside = (iMcm>3 ? 1 : 0); // 0=Aside, 1=Bside
+      // Halfchambers that are switched off, masked by calibration
+      if (calibration->IsHalfChamberNoData(det, halfchamberside))
+       continue;
+
+      for (time = 0; time < nTimeTotal; time++) {
         digits->SetData(row,col,time,signals->GetData(row,col,time));
       } // for: time
     } // for: col
   } // for: row
-  
+
   return kTRUE;
 
 }
@@ -1962,7 +1977,7 @@ void AliTRDdigitizer::RunDigitalProcessing(Int_t det)
   for (Int_t side = 0; side <= 1; side++) {
     for(Int_t rob = side; rob < digits->GetNrow() / 2; rob += 2) {
       for(Int_t mcm = 0; mcm < 16; mcm++) {
-       fMcmSim->Init(det, rob, mcm); 
+       fMcmSim->Init(det, rob, mcm);
        fMcmSim->SetDataByPad(digits, fDigitsManager);
        fMcmSim->Filter();
        if (feeParam->GetTracklet()) {
diff --git a/TRD/Macros/AliTRDanadigits.C b/TRD/Macros/AliTRDanadigits.C
new file mode 100644 (file)
index 0000000..db527ef
--- /dev/null
@@ -0,0 +1,149 @@
+#if !defined(__CINT__) || defined(__MAKECINT__)
+
+
+#include <iostream>
+
+// root
+#include <TStyle.h>
+#include <TH2.h>
+#include <TCanvas.h>
+#include <TStyle.h>
+
+
+
+// aliroot
+#include "AliRun.h"
+#include "AliTRDgeometry.h"
+#include "AliRunLoader.h"
+#include "AliLoader.h"
+#include "AliTRDdigitsManager.h"
+#include "AliTRDarrayADC.h"
+#include "AliTRDfeeParam.h"
+#include "AliTRDdigitsParam.h"
+
+
+
+#endif
+
+
+
+void AliTRDanadigits()
+{
+  gStyle->SetNdivisions(005,"X");
+
+  Int_t iev1 = 1;
+  Int_t iev2 = 20;
+  //
+  // Analyzes the digits
+  //
+  const Int_t kNdet = AliTRDgeometry::Ndet();
+
+  //open run loader and load gAlice, kinematics and header
+  AliRunLoader* rl = AliRunLoader::Open("galice.root");
+  if (!rl)   return;
+  rl->LoadgAlice();
+  gAlice = rl->GetAliRun();
+  if (!gAlice) return;
+
+  // Import the Trees for the event nEvent in the file
+  rl->LoadKinematics();
+  rl->GetEvent(0);
+  rl->GetHeader();
+
+  AliLoader* loader = rl->GetLoader("TRDLoader");
+  if (!loader) {
+    cout << "<AliTRDanalyzeHits> No TRDLoader found" << endl;
+    return;
+  }
+
+  //TH2D
+  Int_t detA = 240;//AliTRDgeometry::GetDetector(0,4,1); //layer,stac,sm
+  Int_t detB = 250;//AliTRDgeometry::GetDetector(5,4,1);
+  printf("detA %d \t detB %d \n",detA,detB);
+  TH2D *chA = new TH2D("chA", "chA", kNdet, 0, kNdet, 200, -50, 50);
+  TH2D *chB = new TH2D("chB", "chB", kNdet, 0, kNdet, 200, -50, 50);
+
+  AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
+
+  for (Int_t ievent = iev1; ievent < iev2; ievent++) {
+    printf("Process event %d\n",ievent);
+    rl->GetEvent(ievent);
+    if (!loader->TreeD()){
+      printf("loader Loading Digits ... \n");
+      loader->LoadDigits();
+    }
+
+    // Read the digits from the file
+    if (!(digitsManager->ReadDigits(loader->TreeD()))) {
+      cout << "<AliTRDanalyzeDigits> Cannot read the digits" << endl;
+      return;
+    }
+    else {
+      printf("digitsManager Read Digits Done\n");
+    }
+
+    // loop over selected chambers
+    for(Int_t det = detA; det < detB; det++) {
+
+      printf("try in %d :\t",det);
+      Int_t        plan = AliTRDgeometry::GetLayer( det );   // Plane
+      Int_t        cham = AliTRDgeometry::GetStack( det ); // Chamber
+      Int_t        sect = AliTRDgeometry::GetSector( det );  // Sector (=iDDL)
+      Int_t        nRow = AliTRDgeometry::GetRowMax( plan, cham, sect );
+      Int_t        nCol = AliTRDgeometry::GetColMax( plan );
+      printf(" sm %d \t lay %d \t stac %d \n",sect,plan,cham);
+
+
+      AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(det); //mod
+      digits->Expand();
+      if (digits->HasData()) { //array
+
+       printf("data\n");
+
+       Int_t nbtimebin  = digitsManager->GetDigitsParam()->GetNTimeBins(det);
+
+       AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();  
+       // Signal for regular pads
+       for (Int_t irow = 0; irow < nRow; irow++ ) {
+         for (Int_t icol = 0; icol < nCol; icol++ ) {
+
+           Int_t rob = feeParam->GetROBfromPad(irow,icol);
+
+           for(Int_t k = 0; k < nbtimebin; k++){
+             //              printf(" row %d \t col %d \t timebin %d \t --> signal %d \n",irow,icol,k,signal);
+             Short_t signal = 999;
+             signal = (Short_t) digits->GetData(irow,icol,k);
+
+             if(signal==999) continue;
+             if(rob%2==0) { chB->Fill(det,signal); }
+             else {chA->Fill(det,signal); }
+           } // loop: ntimebin
+         } // loop: nCol
+       } //loop: nRow
+      } // if: hasdata
+    } // loop: det
+    loader->UnloadDigits();  
+  }// event
+
+  // Display the detector matrix
+  TCanvas *c = new TCanvas("c","",50,50,900,600);
+  c->Divide(2,1);
+
+  chA->SetAxisRange(detA,detB-1,"X");
+  chA->SetTitle("HalfChamberSide A");
+  chA->SetXTitle("Detector number");
+  chA->SetYTitle("charge deposit [a.u]");
+  chB->SetAxisRange(detA,detB-1,"X");
+  chB->SetTitle("HalfChamberSide B");
+  chB->SetXTitle("Detector number");
+  chB->SetYTitle("charge deposit [a.u]");
+
+  c->cd(1);
+  chA->Draw("COLZ");
+
+  c->cd(2);
+  chB->Draw("COLZ");
+
+  //c->SaveAs("ChargeDeposit.pdf");
+}
+