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