]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/Macros/ReadCDD.C
correcting DDL bit for AD and HLT
[u/mrichter/AliRoot.git] / TRD / Macros / ReadCDD.C
CommitLineData
04bf3ff3 1#if !defined( __CINT__) || defined(__MAKECINT__)
2
3#include <iostream>
4#include <TRandom.h>
5#include <TSystem.h>
6#include <TDatime.h>
7#include <TCanvas.h>
8#include <TH1F.h>
9#include <TH2F.h>
10#include <TStyle.h>
11
12#include <AliCDBManager.h>
13#include <AliTRDcalibDB.h>
14#include <AliCDBStorage.h>
15#include <AliCDBEntry.h>
16#include <AliCDBMetaData.h>
17#include <AliGeometry.h>
18#include <AliPID.h>
19
20#include "../TRD/AliTRDgeometry.h"
21#include "../TRD/Cal/AliTRDCalROC.h"
22#include "../TRD/Cal/AliTRDCalPad.h"
23#include "../TRD/Cal/AliTRDCalDet.h"
24
25#endif
26
27Int_t GetPlane(Int_t d);
28Int_t GetSector(Int_t d);
29Int_t GetChamber(Int_t d);
30void ReadCDD(Bool_t residual = kFALSE, Int_t ndet = 319);
31
32void ReadCDD(Bool_t residual, Int_t ndet){
33
34
35 gStyle->SetPalette(1);
36 gStyle->SetOptStat(1111);
37 gStyle->SetPadBorderMode(0);
38 gStyle->SetCanvasColor(10);
39 gStyle->SetPadLeftMargin(0.13);
40 gStyle->SetPadRightMargin(0.05);
41
42
43 //definition number of detector et histos
44 Int_t Ndet = AliTRDgeometry::kNdet;
45 TH1F *Vdriftdet = new TH1F("Vdriftdet","Vdrift det",Ndet,0,Ndet);
46 TH1F *Gaindet = new TH1F("Gaindet","Gain det",Ndet,0,Ndet);
47 TH1F *Vdriftdets = new TH1F("Vdriftdets","Vdrift",20,1.0,2.0);
48 TH1F *Gaindets = new TH1F("Gaindets","Gain",15,0.0,2.0);
49
50
51 TH1F *T0pad = new TH1F("T0pad","T0pad",900,-1.0,1.0);
52 TH1F *Vdriftpad = new TH1F("Vdriftpad","Vdriftpad",100,0.9,2.1);
53 TH1F *Gainpad = new TH1F("Gainpad","Gainpad",150,0.2,1.8);
54
55 //////////////////
56 //Datebase used
57 //////////////////
58 AliCDBManager *man = AliCDBManager::Instance();
59
60
162637e4 61 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
04bf3ff3 62 man->SetSpecificStorage("TRD/Calib/LocalGainFactor","local://.");
63 man->SetSpecificStorage("TRD/Calib/LocalT0","local://.");
64 man->SetSpecificStorage("TRD/Calib/LocalVdrift","local://.");
65 man->SetSpecificStorage("TRD/Calib/ChamberT0","local://.");
66 if(! residual){
67 man->SetSpecificStorage("TRD/Calib/ChamberVdrift","local://.");
68 man->SetSpecificStorage("TRD/Calib/ChamberGainFactor","local://.");
69 }
70 man->SetRun(0);
71
72
73 AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
74 AliTRDgeometry *geo = new AliTRDgeometry();
75
76
77 ////////////////////////////////////////////
78 // drift velocity and gain for detectors
79 ///////////////////////////////////////////
80
81 //drift velocity par detector
82 for(Int_t k = 0; k < Ndet; k++){
83 Vdriftdet->Fill(k+0.5,cal->GetVdriftAverage(k));
84 Vdriftdets->Fill(cal->GetVdriftAverage(k));
85 //cout << "Vdrift average: " << cal->GetVdriftAverage(k) << endl;
86 }
87
88 TCanvas *c1 = new TCanvas("c1","Calibchamberdrift",50,50,600,800);
89 c1->Divide(1,2);
90 c1->cd(1);
91 Vdriftdet->SetStats(0);
92 Vdriftdet->SetXTitle("Det Number");
93 Vdriftdet->SetYTitle("Mean Drift Velocity [cm/#mu s]");
94 Vdriftdet->Draw();
95
96 //gain par detector
97 for(Int_t k = 0; k < Ndet; k++){
98 Gaindet->Fill(k+0.5,cal->GetGainFactorAverage(k));
99 Gaindets->Fill(cal->GetGainFactorAverage(k));
100 //cout << "Vdrift average: " << cal->GetVdriftAverage(k) << endl;
101 }
102
103 //TCanvas *c1g = new TCanvas("c1g","Calibchambergain",50,50,600,800);
104 c1->cd(2);
105 Gaindet->SetStats(0);
106 Gaindet->SetXTitle("Det Number");
107 Gaindet->SetYTitle("Mean Gain");
108 Gaindet->Draw();
109
110
111 TCanvas *c1u = new TCanvas("c1u","Calibchamberdriftu",50,50,600,800);
112 c1u->Divide(1,2);
113 c1u->cd(1);
114 Vdriftdets->SetStats(0);
115 Vdriftdets->SetXTitle("Mean Drift Velocity [cm/#mu s]");
116 Vdriftdets->SetYTitle("Number of detectors");
117 Vdriftdets->Draw();
118
119 c1u->cd(2);
120 Gaindets->SetStats(0);
121 Gaindets->SetXTitle("Mean Gain");
122 Gaindets->SetYTitle("Number of detectors");
123 Gaindets->Draw();
124
125 /////////////////////////////////
126 // Maps of one detector choosen
127 /////////////////////////////////
128
129 Int_t rowMaxd = geo->GetRowMax(GetPlane(ndet),GetChamber(ndet),GetSector(ndet));
130 Int_t colMaxd = geo->GetColMax(GetPlane(ndet));
131
132 TH2F *hZYvd = new TH2F("hZYvd" ,"Y vs Z", rowMaxd,0,rowMaxd,colMaxd,0,colMaxd);
133 TH2F *hZYg = new TH2F("hZYg" ,"Y vs Z", rowMaxd,0,rowMaxd,colMaxd,0,colMaxd);
134 TH2F *hZYT0 = new TH2F("hZYT0" ,"Y vs Z", rowMaxd,0,rowMaxd,colMaxd,0,colMaxd);
135
136
137 for (Int_t col = 0; col < colMaxd; col++) {
138 for (Int_t row = 0; row < rowMaxd; row++) {
139 hZYvd->Fill(row,col,cal->GetVdrift(ndet,col,row));
140 hZYg->Fill(row,col,cal->GetGainFactor(ndet,col,row));
141 hZYT0->Fill(row,col,cal->GetT0(ndet,col,row));
142 }//boucle row
143 }//boucle col
144
145 ////////////////////////////////
146 // coefficients per pad
147 ///////////////////////////////
148
149 for(Int_t isector = 0; isector < 18; isector++){
150 for(Int_t iplane =0; iplane < 6; iplane++){
151 for(Int_t ichamb =0; ichamb < 5; ichamb++){
152
153 Int_t det = AliTRDgeometry::GetDetector(iplane,ichamb,isector);
154 Int_t rowMax = geo->GetRowMax(iplane,ichamb,isector);
155 Int_t colMax = geo->GetColMax(iplane);
156
157 for (Int_t col = 0; col < colMax; col++) {
158 for (Int_t row = 0; row < rowMax; row++) {
159 Gainpad->Fill(cal->GetGainFactor(det,col,row));
160 T0pad->Fill(cal->GetT0(det,col,row));
161 Vdriftpad->Fill(cal->GetVdrift(det,col,row));
162 }//boucle row
163 }//boucle col
164
165 }//loop chamber
166 }//loop plane
167 }
168
169 /////////////////////////////
170 // Plot the choosen detector
171 /////////////////////////////
172
173 TCanvas *c2vd = new TCanvas("c2vd","choosen detector",50,50,600,800);
174 c2vd->Divide(3,1);
175 c2vd->cd(1);
176 hZYvd->SetXTitle("row Number");
177 hZYvd->SetYTitle("col Number");
178 hZYvd->SetZTitle("Drift velocity");
179 hZYvd->SetStats(0);
180 hZYvd->Draw("LEGO");
181
182 c2vd->cd(2);
183 hZYT0->SetXTitle("row Number");
184 hZYT0->SetYTitle("col Number");
185 hZYT0->SetZTitle("T0 [time bin]");
186 hZYT0->SetStats(0);
187 hZYT0->Draw("SURF1");
188
189 c2vd->cd(3);
190 hZYg->SetXTitle("row Number");
191 hZYg->SetYTitle("col Number");
192 hZYg->SetZTitle("gain factor");
193 hZYg->SetStats(0);
194 hZYg->Draw("SURF1");
195
196 /////////////////////////////
197 // Plot the pad info
198 /////////////////////////////
199
200 TCanvas *c5vd = new TCanvas("c5vd","PadInfo",50,50,600,800);
201 c5vd->Divide(3,1);
202 c5vd->cd(1);
203 Vdriftpad->SetXTitle("v_{d} [cm/#mus]");
204 Vdriftpad->SetYTitle("number of pads");
205 Vdriftpad->SetStats(0);
206 Vdriftpad->Draw();
207
208 c5vd->cd(2);
209 T0pad->SetXTitle("T0 [time bin]");
210 T0pad->SetYTitle("number of pads");
211 T0pad->SetStats(0);
212 T0pad->Draw();
213
214
215 c5vd->cd(3);
216 Gainpad->SetXTitle("gain factor");
217 Gainpad->SetYTitle("number of pads");
218 Gainpad->SetStats(0);
219 Gainpad->Draw();
220
221
222}
223//____________________________________________________________________________________________
224Int_t GetPlane(Int_t d)
225{
226 //
227 // Reconstruct the plane number from the detector number
228 //
229
230 return ((Int_t) (d % 6));
231
232}
233
234//_____________________________________________________________________________
235Int_t GetChamber(Int_t d)
236{
237 //
238 // Reconstruct the chamber number from the detector number
239 //
240 Int_t fgkNplan = 6;
241
242 return ((Int_t) (d % 30) / fgkNplan);
243
244}
245
246//_____________________________________________________________________________
247Int_t GetSector(Int_t d)
248{
249 //
250 // Reconstruct the sector number from the detector number
251 //
252 Int_t fg = 30;
253
254 return ((Int_t) (d / fg));
255
256}