]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/Macros/AliTRDcompareDigits.C
Adding Domenico Colella as responsible for SPD part in TRI pp
[u/mrichter/AliRoot.git] / TRD / Macros / AliTRDcompareDigits.C
1 // AliTRDcompareDigits.C
2 //
3 // Author
4 //   Ken Oyama
5 //   Mateusz Ploskon
6 //
7 // Example and diagnostics macro to compare digits
8 // before and after the raw data simulation.
9 // Please simulate raw data in root format with keeping
10 // all intermediate data unerased to run the macro.
11
12 class Histograms {
13 public:
14   Histograms() {
15     adc_all = new TH1F("adc_all", "adc_all", 1025, -1.5, 1023.5 );
16   }
17
18   ~Histograms() {
19     delete adc_all;
20   }
21
22   setattr( Color_t c, Color_t fc, Float_t lw ) {
23     if( c  != -1 ) adc_all->SetLineColor( c );
24     if( fc != -1 ) adc_all->SetFillColor( fc );
25     if( lw != -1 ) adc_all->SetLineWidth( lw );
26   }
27
28   TH1F *adc_all; 
29 };
30
31 AliTRDgeometry *gGeo;
32
33 void AliTRDcompareDigits( Int_t id_start = 0 , Int_t id_end = 539 )
34 {
35   AliCDBManager *cdb = AliCDBManager::Instance();
36   cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
37   cdb->SetRun(0);
38   //AliCDBStorage* localStorage = cdb->GetStorage("local://$ALICE_ROOT/OCDB");
39   //  cout << "[I] Is storage set : " << cdb->IsDefaultStorageSet() << endl;
40   cout << endl;
41
42   // Read raw data after zero suppression
43   AliLog::SetClassDebugLevel("AliTRDRawStream", 10);
44   AliLog::SetFileOutput("decodeZS.log");
45   AliRawReaderRoot reader("raw.root", 0);
46   reader.Select("TRD");
47   AliTRDrawData data;
48   AliTRDdigitsManager *dmz = data.Raw2Digits(&reader);
49
50   gGeo = new AliTRDgeometry(); 
51   Ho   = new Histograms(); // Original
52   Hz   = new Histograms(); // Zero suppressed
53   Ho->setattr( 46 , 46, -1);
54   Hz->setattr( -1 , -1, 2 );
55
56   // Read original data from digits tree
57   AliRunLoader *runLoader  = AliRunLoader::Open("galice.root" ,AliConfig::GetDefaultEventFolderName(),"READ");
58   AliLoader *loader = runLoader->GetLoader("TRDLoader");
59   loader->LoadDigits();
60   AliTRDdigitsManager *dmo = new AliTRDdigitsManager();  // Original digits
61   dmo->ReadDigits(loader->TreeD());
62
63   for( int id = id_start ; id <= id_end ; id++ ) {
64     AliTRDdataArrayI *dao = dmo->GetDigits( id );
65     AliTRDdataArrayI *daz = dmz->GetDigits( id );
66     dao->Expand();
67     daz->Expand();
68     analyze_det( id, dao, Ho );
69     analyze_det( id, daz, Hz );
70     compare_signal_position( id, dao, daz );
71   }
72
73   TCanvas *c = new TCanvas("decodeZS", "decodeZS", 800, 500);
74   c->Draw();
75   c->SetLogy();
76   Ho->adc_all->Draw();
77   Hz->adc_all->Draw("SAME");
78 }
79
80 //
81 // Analyze detector digits, printout data and fill histogram
82 //
83 void analyze_det( Int_t det, AliTRDdataArrayI *d, Histograms *h )
84 {
85   Int_t        plan = gGeo->GetPlane( det );   // Plane
86   Int_t        cham = gGeo->GetChamber( det ); // Chamber
87   Int_t        sect = gGeo->GetSector( det );  // Sector (=iDDL)
88   Int_t        nRow = gGeo->GetRowMax( plan, cham, sect );
89   Int_t        nCol = gGeo->GetColMax( plan );
90   Int_t       nTBin = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
91   Int_t       ctype = 0;                       // Chamber type (0:C0, 1:C1)
92   Int_t        adc;
93
94   printf("Det=%03d Sector=%02d Stack=%d Layer=%d\n", det,sect, cham, plan );
95   // printf("Size=%d\n", d->GetSize());
96
97   for (Int_t irow = 0; irow < nRow; irow++ ) {
98     for (Int_t icol = 0; icol < nCol; icol++ ) {
99       // printf(" Row=%02d Col=%03d :", irow, icol);
100       Int_t signal = 0;
101       for (Int_t it = 0; it < nTBin; it++ ) {
102         adc = d->GetData(irow, icol, it);
103         // printf(" % 3d", adc );
104         h->adc_all->Fill( adc );
105         if( adc > 5 ) signal++;
106       }
107       // if( signal != 0 ) printf(" ... signal!");
108       // printf("\n");
109     }
110   } 
111 }
112
113 //
114 // Compare two digits array
115 //
116 void compare_signal_position( Int_t det, AliTRDdataArrayI *dao, AliTRDdataArrayI *daz )
117 {
118   Int_t        plan = gGeo->GetPlane( det );   // Plane
119   Int_t        cham = gGeo->GetChamber( det ); // Chamber
120   Int_t        sect = gGeo->GetSector( det );  // Sector (=iDDL)
121   Int_t        nRow = gGeo->GetRowMax( plan, cham, sect );
122   Int_t        nCol = gGeo->GetColMax( plan );
123   Int_t       nTBin = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
124   Int_t       wrong = 0;
125   Int_t       adc_o;
126   Int_t       adc_z;
127
128   for (Int_t irow = 0; irow < nRow; irow++ ) {
129     for (Int_t icol = 0; icol < nCol; icol++ ) {
130       for (Int_t it = 0; it < nTBin; it++ ) {
131         adc_o = dao->GetData(irow, icol, it);
132         adc_z = daz->GetData(irow, icol, it);
133         if( adc_o > 5 && adc_o != adc_z ) {
134           printf("              can not find adc=%d at Row=%02d Col=%03d TimeBin=%02d in ZS digits.\n", adc_o, irow, icol, it ); 
135           wrong++;
136         }
137       }
138     }
139   } 
140   if( wrong != 0 ) printf("               %d wrong ADC hit(s) found !!!\n", wrong);
141 }