]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/patchSSDBadChannelsMap.C
Macro for SDD dE/dx calibration (from Leonardo and Riccardo)
[u/mrichter/AliRoot.git] / ITS / patchSSDBadChannelsMap.C
CommitLineData
a65a7e70 1//
2// Macro adds ether ladder, module or channel to the SSD bad channels map.
3// Input parameter: file name with original SSD bad channels map (optional),
4// equipment identifier (DDL),
5// Slot number, ADC number (optional) and channel number (optional)
6// Author: Oleksandr.Borysov@cern.ch
7//
8
9
10#if !defined(__CINT__) || defined(__MAKECINT__)
11#include <time.h>
12#include <Riostream.h>
13#include "TStyle.h"
14#include "TFile.h"
15#include "TCanvas.h"
16#include "TH2F.h"
17#include "TString.h"
18#include "AliLog.h"
19#include "AliRawReaderRoot.h"
20#include "AliITSRawStreamSSDv1.h"
21#include "AliITSBadChannelsSSDv2.h"
22#include "AliITSModuleDaSSD.h"
23#include "AliITSHandleDaSSD.h"
24#endif
25
26
27#define EQUIPMENTTODDLMASK 0xFF
28#define NUMBEROFSSDMODULESPERSLOT 12
29#define MAXSSDDDLID 15
30#define MINSSDMODULEID 500
31
32Bool_t readBCMap(const Char_t *filename, AliITSBadChannelsSSDv2*& bcl);
33Int_t addModuleBCMap(AliITSBadChannelsSSDv2 *bcl, const Int_t EqId, const Int_t SlotId, const Int_t adc);
34Int_t addLadderBCMap(AliITSBadChannelsSSDv2 *bcl, const Int_t EqId, const Int_t SlotId);
35Int_t addChannelBCMap(AliITSBadChannelsSSDv2 *bcl, const Int_t EqId, const Int_t SlotId, const Int_t adc, const Int_t strp);
36void Usage (void);
37void drawBadChannelMapDAQDB(const char* filename);
38
39//class gives an access to the protected array with SSD_DDL_Map in case the geometry is not initialized
40class AliITSRawStreamSSDtmp : public AliITSRawStreamSSDv1 {
41public:
42 AliITSRawStreamSSDtmp(AliRawReader* rawReader): AliITSRawStreamSSDv1(rawReader) {
43 Setv11HybridDDLMapping();
44 AliInfo("Using SSD DDL Map initialized by AliITSRawStreamSSDv1::Setv11HybridDDLMapping()");
45 }
46 virtual ~AliITSRawStreamSSDtmp() {};
47 Int_t GetModId(Int_t iDDL, Int_t iModule) {return fgkDDLModuleMap[iDDL][iModule];}
48 ClassDef(AliITSRawStreamSSDtmp, 0)
49};
50
51//_________________________________________________________//
52Int_t patchSSDBadChannelsMap(const Char_t *fname = 0,
53 const Int_t EqId = -1,
54 const Int_t SlotId = -1,
55 const Int_t adc = -1,
56 const Int_t strp = -1) {
57 //Macro to patch the bad channel list
58 TString bcfname, pathstr;
59 AliITSBadChannelsSSDv2 *bc = 0;
60 if (EqId < 0) { Usage(); cerr << "Equipment number (that is DDL) must be specified! Exit.\n"; return 0; }
61 else if ((EqId & EQUIPMENTTODDLMASK) > MAXSSDDDLID) {
62 cerr << "Icorrect equipment number (that is DDL)! Exit.\n";
63 return -1;
64 }
65 if (SlotId < 1 || SlotId > 9) {
66 cerr << "Slot number must be specified (in the range 1 - 9)! Exit.\n";
67 Usage();
68 return -2;
69 }
70 if (!fname || fname[0]==0) {
71 cout << "File name with bad channels map is not specified, an empty map is used!\n";
72 bc = new AliITSBadChannelsSSDv2();
73 if (!bc) { cerr << "Error creating the AliITSBadChannelsSSDv2 object! Exit.\n"; return -1; }
74 pathstr = "";
75 }
76 else {
77 pathstr.Form(fname);
78 if (!readBCMap(fname, bc)) {
79 cerr << "Error reading file " << fname << " with Static Bad Channels Map! Exit.\n";
80 return -3;
81 }
82 }
83 if (adc < 0) addLadderBCMap(bc, EqId, SlotId);
84 else if ((adc < 6) || (adc > 7 && (adc - 2) < NUMBEROFSSDMODULESPERSLOT) ) {
85 if (strp < 0) addModuleBCMap(bc, EqId, SlotId, adc);
86 else if (strp < AliITSModuleDaSSD::GetStripsPerModuleConst()) addChannelBCMap(bc, EqId, SlotId, adc, strp);
87 else {
88 cerr << "Incorrect number for Strip. Exit\n";
89 Usage();
90 if (bc) delete bc;
91 return -5;
92 }
93 }
94 else {
95 cerr << "Incorrect number for ADC. Exit\n";
96 if (bc) delete bc;
97 return -4;
98 }
99 bcfname = pathstr(0, pathstr.Last('/')+1);
100 bcfname.Append(Form("ssdbcmap_%i.root", time(NULL)));
101 TFile *bcfile = new TFile (bcfname.Data(),"RECREATE");
102 if (bcfile->IsZombie()) {
103 cerr << "Error open file " << bcfname.Data() << " for writing new bad channels map!\n";
104 if (bc) delete bc;
105 return -1;
106 }
107 bcfile->WriteTObject(bc);
108 bcfile->Close();
109 delete bcfile;
110 cout << "New SSD bad channels map was saved in file " << bcfname.Data() << endl;
111 if (bc) delete bc;
112 return 0;
113}
114
115//_________________________________________________________//
116Bool_t readBCMap(const Char_t *filename, AliITSBadChannelsSSDv2*& bcl) {
117 // Reads Static Bad Channels Map from the file
118 TFile *bcfile;
119 if (!filename) {
120 cout << "No file name is specified for Static Bad Channels Map!\n";
121 return kFALSE;
122 }
123 cout << "Reading SSD Bad Channels Map from the file " << filename << endl;
124 bcfile = new TFile(filename, "READ");
125 if (bcfile->IsZombie()) {
126 cerr << "Error reading file " << filename << " with Static Bad Channels Map!\n";
127 return kFALSE;
128 }
129 bcfile->GetObject("AliITSBadChannelsSSDv2;1", bcl);
130 if (!bcl) {
131 cout << "Error bcl == NULL!\n";
132 bcfile->Close();
133 delete bcfile;
134 return kFALSE;
135 }
136 bcfile->Close();
137 delete bcfile;
138 return kTRUE;
139}
140
141//_________________________________________________________//
142Int_t addModuleBCMap(AliITSBadChannelsSSDv2 *bcl,
143 const Int_t EqId, const Int_t SlotId, const Int_t adc) {
144 // Add module to bad channels map.
145 const Char_t isbad = 3;
146 Int_t ddl, mn, modid;
147 AliRawReaderRoot *rwr = 0;
148 AliITSRawStreamSSDtmp *rst = 0;
149 rwr = new AliRawReaderRoot();
150 rst = new AliITSRawStreamSSDtmp(rwr);
151 ddl = EqId & EQUIPMENTTODDLMASK;
152 mn = (SlotId - 1) * NUMBEROFSSDMODULESPERSLOT + (adc<8 ? adc : adc-2);
153 modid = rst->GetModId(ddl, mn);
154 modid -= MINSSDMODULEID;
155 if (modid < 0) return 0;
156 for (Int_t strind = 0; strind < AliITSModuleDaSSD::GetPNStripsPerModule(); strind++) {
157 bcl->AddBadChannelP(modid, strind, isbad);
158 bcl->AddBadChannelN(modid, strind, isbad);
159 }
160 delete rst;
161 delete rwr;
162 cout << "Module added: ModId = " << modid + MINSSDMODULEID << "; Ddl/Ad/Adc: " << ddl
163 << " / " << SlotId << " / " << adc << endl;
164 return 0;
165}
166
167//_________________________________________________________//
168Int_t addLadderBCMap(AliITSBadChannelsSSDv2 *bcl,
169 const Int_t EqId, const Int_t SlotId) {
170 // Add ladder to bad channels map.
171 for (Int_t adc = 0; adc < NUMBEROFSSDMODULESPERSLOT; adc++)
172 addModuleBCMap(bcl, EqId, SlotId, (adc<6 ? adc : adc+2));
173 return 0;
174}
175
176//_________________________________________________________//
177Int_t addChannelBCMap(AliITSBadChannelsSSDv2 *bcl,
178 const Int_t EqId, const Int_t SlotId,
179 const Int_t adc, const Int_t strp) {
180 // Add strip to bad channels map.
181 const Char_t isbad = 3;
182 Int_t ddl, mn, modid;
183 AliRawReaderRoot *rwr = 0;
184 AliITSRawStreamSSDtmp *rst = 0;
185 rwr = new AliRawReaderRoot();
186 rst = new AliITSRawStreamSSDtmp(rwr);
187 ddl = EqId & EQUIPMENTTODDLMASK;
188 mn = (SlotId - 1) * NUMBEROFSSDMODULESPERSLOT + (adc<8 ? adc : adc-2);
189 modid = rst->GetModId(ddl, mn);
190 modid -= MINSSDMODULEID;
191 if (modid < 0) { cout << "There is no module with given Equipment, Slot, adc.\n" ; return 0; }
192 if (strp < AliITSModuleDaSSD::GetPNStripsPerModule() ) {
193 bcl->AddBadChannelP(modid, strp, isbad);
194 } else {
195 bcl->AddBadChannelN(modid, (AliITSChannelDaSSD::GetMaxStripIdConst() - strp), isbad);
196 }
197 delete rst;
198 delete rwr;
199 cout << "Channel added (ModId/Ddl/Ad/Adc/Strip): " << modid + MINSSDMODULEID << " / " << ddl << " / " << SlotId
200 << " / " << adc << " / " << strp << endl;
201 return 0;
202}
203
204//_________________________________________________________//
205void Usage (void) {
206 //Usage function
207 cout << "Usage: PatchSSDBadChannelsMap(bc_fname /* can be \"\" */, EqipmentId, SlotId, adc /*optional*/, strip /*optional*/)\n";
208}
209
210//_______________________________________//
211void drawBadChannelMapDAQDB(const char* filename) {
212 //Draws the 2D plots of the bad channels maps
213 const Int_t fgkSSDMODULES = 1698;
214 static const Int_t fgkDefaultNStripsSSD = 768;
215 gStyle->SetPalette(1,0);
216
217 TH2F *fHistDAQDBPSideBadChannelMapLayer5 = new TH2F("fHistDAQDBPSideBadChannelMapLayer5",
218 "Layer 5;N_{module};N_{ladder}",
219 22,1,23,
220 34,500,534);
221 fHistDAQDBPSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
222 fHistDAQDBPSideBadChannelMapLayer5->SetStats(kFALSE);
223 fHistDAQDBPSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
224 fHistDAQDBPSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
225 fHistDAQDBPSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
226 fHistDAQDBPSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
227 fHistDAQDBPSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
228 fHistDAQDBPSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
229 fHistDAQDBPSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
230
231 TH2F *fHistDAQDBNSideBadChannelMapLayer5 = new TH2F("fHistDAQDBNSideBadChannelMapLayer5",
232 "Layer 5;N_{module};N_{ladder}",
233 22,1,23,
234 34,500,534);
235 fHistDAQDBNSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
236 fHistDAQDBNSideBadChannelMapLayer5->SetStats(kFALSE);
237 fHistDAQDBNSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
238 fHistDAQDBNSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
239 fHistDAQDBNSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
240 fHistDAQDBNSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
241 fHistDAQDBNSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
242 fHistDAQDBNSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
243 fHistDAQDBNSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
244
245 TH2F *fHistDAQDBPSideBadChannelMapLayer6 = new TH2F("fHistDAQDBPSideBadChannelMapLayer6",
246 "Layer 6;N_{module};N_{ladder}",
247 25,1,26,
248 38,600,638);
249 fHistDAQDBPSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
250 fHistDAQDBPSideBadChannelMapLayer6->SetStats(kFALSE);
251 fHistDAQDBPSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
252 fHistDAQDBPSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
253 fHistDAQDBPSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
254 fHistDAQDBPSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
255 fHistDAQDBPSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
256 fHistDAQDBPSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
257 fHistDAQDBPSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
258
259 TH2F *fHistDAQDBNSideBadChannelMapLayer6 = new TH2F("fHistDAQDBNSideBadChannelMapLayer6",
260 "Layer 6;N_{module};N_{ladder}",
261 25,1,26,
262 38,600,638);
263 fHistDAQDBNSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
264 fHistDAQDBNSideBadChannelMapLayer6->SetStats(kFALSE);
265 fHistDAQDBNSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
266 fHistDAQDBNSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
267 fHistDAQDBNSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
268 fHistDAQDBNSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
269 fHistDAQDBNSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
270 fHistDAQDBNSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
271 fHistDAQDBNSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
272
273 //===============================//
274 TFile *f = TFile::Open(filename);
275 if(!f) {
276 Printf("File poiter not valid");
277 return;
278 }
279
280 if(!f->IsOpen()) {
281 Printf("The file was not found");
282 return;
283 }
284 //===============================//
285
286 AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
287 badChannelsSSD = dynamic_cast<AliITSBadChannelsSSDv2 *>(f->Get("AliITSBadChannelsSSDv2"));
288
289 //_____________________________________________________________________________//
290 Int_t nPSideChannelsTotal = 0, nNSideChannelsTotal = 0;
291 Int_t nBadPSideChannelsTotal = 0, nBadNSideChannelsTotal = 0;
292 Int_t nBadPSideChannels = 0, nBadNSideChannels = 0;
293 Int_t layer = 0, ladder = 0, module = 0;
294 Int_t nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
295 Int_t nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
296 //_____________________________________________________________________________//
297
298 for(Int_t i = 0; i < fgkSSDMODULES; i++) {
299 //for(Int_t i = 0; i < 1; i++) {
300 AliITSgeomTGeo::GetModuleId(i+500,layer,ladder,module);
301 nBadPSideChannels = 0, nBadNSideChannels = 0;
302 nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
303 nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
304
305 Int_t badChannel = 0;
306 for(Int_t j = 0; j < fgkDefaultNStripsSSD; j++) {
307 badChannel = (Int_t)(badChannelsSSD->GetBadChannelP(i,j));
308 //cout<<"Module: "<<i+500<< " Strip: "<<j<<" - "<<badChannel<<endl;
309 if(badChannel != 0) {
310 if(layer == 5)
311 nPSideChannelsLayer5 += 1;
312 if(layer == 6)
313 nPSideChannelsLayer6 += 1;
314 nBadPSideChannels += 1;
315 }
316 badChannel = (Int_t)(badChannelsSSD->GetBadChannelN(i,j));
317 //cout<<"Module: "<<i+500<< " Strip: "<<fgkDefaultNStripsSSD+j+1<<" - "<<badChannel<<endl;
318 if(badChannel != 0) {
319 if(layer == 5)
320 nNSideChannelsLayer5 += 1;
321 if(layer == 6)
322 nNSideChannelsLayer6 += 1;
323 nBadNSideChannels += 1;
324 }
325 }
326 if(layer == 5) {
327 if(nPSideChannelsLayer5 > 0)
328 fHistDAQDBPSideBadChannelMapLayer5->Fill(module,499+ladder,
329 100.*nPSideChannelsLayer5/fgkDefaultNStripsSSD);
330 else fHistDAQDBPSideBadChannelMapLayer5->Fill(module,499+ladder,0.0001);
331 if(nNSideChannelsLayer5 > 0)
332 fHistDAQDBNSideBadChannelMapLayer5->Fill(module,499+ladder,
333 100.*nNSideChannelsLayer5/fgkDefaultNStripsSSD);
334 else fHistDAQDBNSideBadChannelMapLayer5->Fill(module,499+ladder,0.0001);
335 }//layer 5
336 if(layer == 6) {
337 if(nPSideChannelsLayer6 > 0)
338 fHistDAQDBPSideBadChannelMapLayer6->Fill(module,599+ladder,
339 100.*nPSideChannelsLayer6/fgkDefaultNStripsSSD);
340 else fHistDAQDBPSideBadChannelMapLayer6->Fill(module,599+ladder,0.0001);
341 if(nNSideChannelsLayer6 > 0)
342 fHistDAQDBNSideBadChannelMapLayer6->Fill(module,599+ladder,
343 100.*nNSideChannelsLayer6/fgkDefaultNStripsSSD);
344 else fHistDAQDBNSideBadChannelMapLayer6->Fill(module,599+ladder,0.0001);
345 }//layer 6
346
347 nBadPSideChannelsTotal += nBadPSideChannels;
348 nBadNSideChannelsTotal += nBadNSideChannels;
349 nPSideChannelsTotal += fgkDefaultNStripsSSD;
350 nNSideChannelsTotal += fgkDefaultNStripsSSD;
351 }
352
353 cout<<"================================="<<endl;
354 cout<<"Bad p-Side channels: "<<100.*nBadPSideChannelsTotal/nPSideChannelsTotal<<endl;
355 cout<<"Bad n-Side channels: "<<100.*nBadNSideChannelsTotal/nNSideChannelsTotal<<endl;
356 cout<<"================================="<<endl;
357
358 TCanvas *cBadChannelDAQDB = new TCanvas("cBadChannelDAQDB",
359 "Bad channel list - DAQ DB",
360 0,0,900,900);
361 cBadChannelDAQDB->SetHighLightColor(10); cBadChannelDAQDB->SetFillColor(10);
362 cBadChannelDAQDB->Divide(2,2);
363
364 cBadChannelDAQDB->cd(1)->SetBottomMargin(.2);
365 cBadChannelDAQDB->cd(1)->SetLeftMargin(.15);
366 cBadChannelDAQDB->cd(1)->SetRightMargin(.2);
367 cBadChannelDAQDB->cd(1)->SetGridx(); cBadChannelDAQDB->cd(1)->SetGridy();
368 cBadChannelDAQDB->cd(1); fHistDAQDBPSideBadChannelMapLayer5->Draw("colz");
369 cBadChannelDAQDB->cd(2)->SetBottomMargin(.2);
370 cBadChannelDAQDB->cd(2)->SetLeftMargin(.15);
371 cBadChannelDAQDB->cd(2)->SetRightMargin(.2);
372 cBadChannelDAQDB->cd(2)->SetGridx(); cBadChannelDAQDB->cd(2)->SetGridy();
373 cBadChannelDAQDB->cd(2); fHistDAQDBPSideBadChannelMapLayer6->Draw("colz");
374 cBadChannelDAQDB->cd(3)->SetBottomMargin(.2);
375 cBadChannelDAQDB->cd(3)->SetLeftMargin(.15);
376 cBadChannelDAQDB->cd(3)->SetRightMargin(.2);
377 cBadChannelDAQDB->cd(3)->SetGridx(); cBadChannelDAQDB->cd(3)->SetGridy();
378 cBadChannelDAQDB->cd(3); fHistDAQDBNSideBadChannelMapLayer5->Draw("colz");
379 cBadChannelDAQDB->cd(4)->SetBottomMargin(.2);
380 cBadChannelDAQDB->cd(4)->SetLeftMargin(.15);
381 cBadChannelDAQDB->cd(4)->SetRightMargin(.2);
382 cBadChannelDAQDB->cd(4)->SetGridx(); cBadChannelDAQDB->cd(4)->SetGridy();
383 cBadChannelDAQDB->cd(4); fHistDAQDBNSideBadChannelMapLayer6->Draw("colz");
384
385 return;
386}