]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/Macros/AliTRDanadigits.C
- coverity fixes for unused but set variables in various places
[u/mrichter/AliRoot.git] / TRD / Macros / AliTRDanadigits.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2
3
4 #include <iostream>
5
6 // root
7 #include <TStyle.h>
8 #include <TH2.h>
9 #include <TCanvas.h>
10 #include <TStyle.h>
11
12
13
14 // aliroot
15 #include "AliRun.h"
16 #include "AliTRDgeometry.h"
17 #include "AliRunLoader.h"
18 #include "AliLoader.h"
19 #include "AliTRDdigitsManager.h"
20 #include "AliTRDarrayADC.h"
21 #include "AliTRDfeeParam.h"
22 #include "AliTRDdigitsParam.h"
23
24
25
26 #endif
27
28
29
30 void AliTRDanadigits()
31 {
32   gStyle->SetNdivisions(005,"X");
33
34   Int_t iev1 = 1;
35   Int_t iev2 = 20;
36   //
37   // Analyzes the digits
38   //
39   const Int_t kNdet = AliTRDgeometry::Ndet();
40
41   //open run loader and load gAlice, kinematics and header
42   AliRunLoader* rl = AliRunLoader::Open("galice.root");
43   if (!rl)   return;
44   rl->LoadgAlice();
45   gAlice = rl->GetAliRun();
46   if (!gAlice) return;
47
48   // Import the Trees for the event nEvent in the file
49   rl->LoadKinematics();
50   rl->GetEvent(0);
51   rl->GetHeader();
52
53   AliLoader* loader = rl->GetLoader("TRDLoader");
54   if (!loader) {
55     cout << "<AliTRDanalyzeHits> No TRDLoader found" << endl;
56     return;
57   }
58
59   //TH2D
60   Int_t detA = 240;//AliTRDgeometry::GetDetector(0,4,1); //layer,stac,sm
61   Int_t detB = 250;//AliTRDgeometry::GetDetector(5,4,1);
62   printf("detA %d \t detB %d \n",detA,detB);
63   TH2D *chA = new TH2D("chA", "chA", kNdet, 0, kNdet, 200, -50, 50);
64   TH2D *chB = new TH2D("chB", "chB", kNdet, 0, kNdet, 200, -50, 50);
65
66   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
67
68   for (Int_t ievent = iev1; ievent < iev2; ievent++) {
69     printf("Process event %d\n",ievent);
70     rl->GetEvent(ievent);
71     if (!loader->TreeD()){
72       printf("loader Loading Digits ... \n");
73       loader->LoadDigits();
74     }
75
76     // Read the digits from the file
77     if (!(digitsManager->ReadDigits(loader->TreeD()))) {
78       cout << "<AliTRDanalyzeDigits> Cannot read the digits" << endl;
79       return;
80     }
81     else {
82       printf("digitsManager Read Digits Done\n");
83     }
84
85     // loop over selected chambers
86     for(Int_t det = detA; det < detB; det++) {
87
88       printf("try in %d :\t",det);
89       Int_t        plan = AliTRDgeometry::GetLayer( det );   // Plane
90       Int_t        cham = AliTRDgeometry::GetStack( det ); // Chamber
91       Int_t        sect = AliTRDgeometry::GetSector( det );  // Sector (=iDDL)
92       Int_t        nRow = AliTRDgeometry::GetRowMax( plan, cham, sect );
93       Int_t        nCol = AliTRDgeometry::GetColMax( plan );
94       printf(" sm %d \t lay %d \t stac %d \n",sect,plan,cham);
95
96
97       AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(det); //mod
98       digits->Expand();
99       if (digits->HasData()) { //array
100
101         printf("data\n");
102
103         Int_t nbtimebin  = digitsManager->GetDigitsParam()->GetNTimeBins(det);
104
105         AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();  
106         // Signal for regular pads
107         for (Int_t irow = 0; irow < nRow; irow++ ) {
108           for (Int_t icol = 0; icol < nCol; icol++ ) {
109
110             Int_t rob = feeParam->GetROBfromPad(irow,icol);
111
112             for(Int_t k = 0; k < nbtimebin; k++){
113               //              printf(" row %d \t col %d \t timebin %d \t --> signal %d \n",irow,icol,k,signal);
114               Short_t signal = 999;
115               signal = (Short_t) digits->GetData(irow,icol,k);
116
117               if(signal==999) continue;
118               if(rob%2==0) { chB->Fill(det,signal); }
119               else {chA->Fill(det,signal); }
120             } // loop: ntimebin
121           } // loop: nCol
122         } //loop: nRow
123       } // if: hasdata
124     } // loop: det
125     loader->UnloadDigits();  
126   }// event
127
128   // Display the detector matrix
129   TCanvas *c = new TCanvas("c","",50,50,900,600);
130   c->Divide(2,1);
131
132   chA->SetAxisRange(detA,detB-1,"X");
133   chA->SetTitle("HalfChamberSide A");
134   chA->SetXTitle("Detector number");
135   chA->SetYTitle("charge deposit [a.u]");
136   chB->SetAxisRange(detA,detB-1,"X");
137   chB->SetTitle("HalfChamberSide B");
138   chB->SetXTitle("Detector number");
139   chB->SetYTitle("charge deposit [a.u]");
140
141   c->cd(1);
142   chA->Draw("COLZ");
143
144   c->cd(2);
145   chB->Draw("COLZ");
146
147   //c->SaveAs("ChargeDeposit.pdf");
148 }
149