]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/ShowCalibrationSDD.C
Added loop for extraction of clusters really attached to its track.
[u/mrichter/AliRoot.git] / ITS / macrosSDD / ShowCalibrationSDD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TH1F.h>
4 #include <TH2I.h>
5 #include <TGraph.h>
6 #include <TExec.h>
7 #include <TStyle.h>
8 #include <TString.h>
9 #include <TSystem.h>
10 #include <TGrid.h>
11 #include <TLine.h>
12 #include <TCanvas.h>
13 #include <TObjArray.h>
14 #include "AliCDBEntry.h"
15 #include "AliITSCalibrationSDD.h"
16 #include "AliITSgeomTGeo.h"
17 #endif
18
19 /* $Id$ */
20
21 // Macro to plot the calibration parameters from the OCDB file 
22 // created from PEDESTAL and PULSER runs (OCDB/ITS/Calib/CalibSDD)
23 // Two methods ShowCalibrationSDD:
24 //  - the first takes the name of the file to be displayed
25 //  - the second builds the alien path+name from run number and file version
26 //
27 // Origin: F. Prino (prino@to.infn.it)
28
29 void MakePalette(){
30   Int_t palette[3]={kGray,2,3};  
31   gStyle->SetPalette(3,palette);
32 }
33
34 void ShowCalibrationSDD(Char_t *filnam="$ALICE_ROOT/OCDB/ITS/Calib/CalibSDD/Run0_9999999_v0_s0.root", Int_t iMod=0){
35
36
37   TFile *f=TFile::Open(filnam);
38   AliCDBEntry *ent=(AliCDBEntry*)f->Get("AliCDBEntry");
39   TH2I* hlay3=new TH2I("hlay3","Layer 3",12,-0.5,5.5,14,-0.5,13.5);
40   hlay3->GetXaxis()->SetTitle("Detector");
41   hlay3->GetYaxis()->SetTitle("Ladder");
42   hlay3->GetXaxis()->SetTickLength(0);
43   hlay3->GetYaxis()->SetTickLength(0);
44   hlay3->SetStats(0);
45   hlay3->SetMinimum(-1);
46   TH2I* hlay4=new TH2I("hlay4","Layer 4",16,-0.5,7.5,22,-0.5,21.5);
47   hlay4->GetXaxis()->SetTitle("Detector");
48   hlay4->GetYaxis()->SetTitle("Ladder");
49   hlay4->GetXaxis()->SetTickLength(0);
50   hlay4->GetYaxis()->SetTickLength(0);
51   hlay4->GetYaxis()->SetTitle("Ladder");
52   hlay4->SetStats(0);
53   hlay4->SetMinimum(-1);
54   TH2I* hdeadlay3=new TH2I("hdlay3","Layer 3",6,-0.5,5.5,14,-0.5,13.5);
55   hdeadlay3->GetXaxis()->SetTitle("Detector");
56   hdeadlay3->GetYaxis()->SetTitle("Ladder");
57   hdeadlay3->GetXaxis()->SetTickLength(0);
58   hdeadlay3->GetYaxis()->SetTickLength(0);
59   hdeadlay3->SetStats(0);
60   hdeadlay3->SetMinimum(-1.);
61   TH2I* hdeadlay4=new TH2I("hdlay4","Layer 4",8,-0.5,7.5,22,-0.5,21.5);
62   hdeadlay4->GetXaxis()->SetTitle("Detector");
63   hdeadlay4->GetYaxis()->SetTitle("Ladder");
64   hdeadlay4->GetXaxis()->SetTickLength(0);
65   hdeadlay4->GetYaxis()->SetTickLength(0);
66   hdeadlay4->GetYaxis()->SetTitle("Ladder");
67   hdeadlay4->SetStats(0);
68   hdeadlay4->SetMinimum(-1.);
69
70   TObjArray *calSDD = (TObjArray *)ent->GetObject();
71   printf("Entries in array=%d\n",calSDD->GetEntriesFast());
72   TH1F* hmodstatus=new TH1F("hmodstatus","",260,0.5,260.5);
73   TH1F* hnbadch=new TH1F("hnbadch","",260,0.5,260.5);
74   TH1F* hbase=new TH1F("hbase","",60,0.5,120.5);
75   TH2F* hbasemod=new TH2F("hbasemod","",260,239.5,499.5,50,0.,100.);
76   TH1F* hnoise=new TH1F("hnoise","",100,0.,7.);
77   TH2F* hnoisemod=new TH2F("hnoisemod","",260,239.5,499.5,50,0.,10.);
78   TH1F* hgain=new TH1F("hgain","",100,0.,4.);
79   TH2F* hgainmod=new TH2F("hgainmod","",260,239.5,499.5,50,0.,4.);
80   TH1F* hchstatus=new TH1F("hchstatus","",2,-0.5,1.5);
81
82   AliITSDDLModuleMapSDD* dmap=new AliITSDDLModuleMapSDD();
83   dmap->SetJun09Map();
84   TH2I *hddlcarlos=new TH2I("hddlcarlos","",24,-0.5,11.5,24,-0.5,23.5);
85   hddlcarlos->GetXaxis()->SetTitle("Module");
86   hddlcarlos->GetYaxis()->SetTitle("DDL");
87   hddlcarlos->GetXaxis()->SetTickLength(0);
88   hddlcarlos->GetYaxis()->SetTickLength(0);
89   hddlcarlos->SetStats(0);
90   hddlcarlos->SetMinimum(-1.);
91
92   AliITSCalibrationSDD *cal;
93   Int_t badModCounter3=0;
94   Int_t badModCounter4=0;
95   Int_t badHybridCounter3=0;
96   Int_t badHybridCounter4=0;
97   Int_t badAnodeCounter3=0;
98   Int_t badAnodeCounter4=0;
99   Int_t badAnodeCounterGoodMod3=0;
100   Int_t badAnodeCounterGoodMod4=0;
101   Int_t badAnodeCounterGoodHybrid3=0;
102   Int_t badAnodeCounterGoodHybrid4=0;
103   Int_t badAnodeCounterGoodModAndChip3=0;
104   Int_t badAnodeCounterGoodModAndChip4=0;
105   Int_t badChipCounter3=0;
106   Int_t badChipCounter4=0;
107   for(Int_t i=0; i<260; i++){
108     cal=(AliITSCalibrationSDD*)calSDD->At(i);
109     if(cal==0) continue;
110     printf("Module %d (%d)   status = ",i,i+240);
111     Int_t lay,lad,det;
112     AliITSgeomTGeo::GetModuleId(i+240,lay,lad,det);
113     Int_t index=1+(det-1)*2;
114     Int_t ddl,carlos;
115     dmap->FindInDDLMap(i+240,ddl,carlos);
116     Int_t index2=1+carlos*2;
117     ddl+=1;
118     if(cal->IsBad()){ 
119       printf("BAD\t");
120       hddlcarlos->SetBinContent(index2,ddl,0);
121       hddlcarlos->SetBinContent(index2+1,ddl,0);
122       if(lay==3){ 
123         badModCounter3++;
124         badHybridCounter3+=2;
125         hlay3->SetBinContent(index,lad,0);
126         hlay3->SetBinContent(index+1,lad,0);
127       }else if(lay==4){ 
128         badModCounter4++;
129         badHybridCounter4+=2;
130         hlay4->SetBinContent(index,lad,0);
131         hlay4->SetBinContent(index+1,lad,0);
132       }
133       hmodstatus->SetBinContent(i+1,0);
134     }else{ 
135       printf("OK\t");
136       hmodstatus->SetBinContent(i+1,1);
137       if(lay==3){ 
138         badAnodeCounterGoodMod3+=cal->GetDeadChannels();
139         if(cal->IsChipBad(0) && cal->IsChipBad(1) && cal->IsChipBad(2) && cal->IsChipBad(3)){
140           hlay3->SetBinContent(index,lad,0);
141           hddlcarlos->SetBinContent(index2,ddl,0);
142           badHybridCounter3++;
143         }else{
144           hlay3->SetBinContent(index,lad,1);
145           hddlcarlos->SetBinContent(index2,ddl,1);
146           for(Int_t iAn=0; iAn<256; iAn++){
147             if(cal->IsBadChannel(iAn)) badAnodeCounterGoodHybrid3++;
148           }
149         }
150         if(cal->IsChipBad(4) && cal->IsChipBad(5) && cal->IsChipBad(6) && cal->IsChipBad(7)){
151           hlay3->SetBinContent(index+1,lad,0);
152           hddlcarlos->SetBinContent(index2+1,ddl,0);
153           badHybridCounter3++;
154         }else{
155           hlay3->SetBinContent(index+1,lad,1);
156           hddlcarlos->SetBinContent(index2+1,ddl,1);
157           for(Int_t iAn=256; iAn<512; iAn++){
158             if(cal->IsBadChannel(iAn)) badAnodeCounterGoodHybrid3++;
159           }
160         }
161       }else{ 
162         badAnodeCounterGoodMod4+=cal->GetDeadChannels();
163         if(cal->IsChipBad(0) && cal->IsChipBad(1) && cal->IsChipBad(2) && cal->IsChipBad(3)){
164           hlay4->SetBinContent(index,lad,0);
165           hddlcarlos->SetBinContent(index2,ddl,0);
166           badHybridCounter4++;
167         }else{
168           hlay4->SetBinContent(index,lad,1);
169           hddlcarlos->SetBinContent(index2,ddl,1);
170           for(Int_t iAn=0; iAn<256; iAn++){
171             if(cal->IsBadChannel(iAn)) badAnodeCounterGoodHybrid4++;
172           }
173         }
174         if(cal->IsChipBad(4) && cal->IsChipBad(5) && cal->IsChipBad(6) && cal->IsChipBad(7)){
175           hlay4->SetBinContent(index+1,lad,0);
176           hddlcarlos->SetBinContent(index2+1,ddl,0);
177           badHybridCounter4++;
178         }else{
179           hlay4->SetBinContent(index+1,lad,1);
180           hddlcarlos->SetBinContent(index2+1,ddl,1);
181           for(Int_t iAn=256; iAn<512; iAn++){
182             if(cal->IsBadChannel(iAn)) badAnodeCounterGoodHybrid4++;
183           }
184         }
185       }
186     }
187     printf("   Chip Status (0=OK, 1=BAD): ");  
188     for(Int_t ic=0; ic<8;ic++){ 
189       printf("%d ",cal->IsChipBad(ic));
190       if(cal->IsChipBad(ic) && !cal->IsBad()){ 
191         if(i<84) badChipCounter3++;
192         else badChipCounter4++;
193       }
194     }
195     printf(" # bad anodes = %d  ",cal->GetDeadChannels());
196     if(cal->IsAMAt20MHz()) printf("      20 MHz sampling");
197     else printf("      40 MHz sampling");
198     printf(" Threshold L %d %d H %d %d\n",cal->GetZSLowThreshold(0),cal->GetZSLowThreshold(1),cal->GetZSHighThreshold(0),cal->GetZSHighThreshold(1));
199     if(i<84) badAnodeCounter3+=cal->GetDeadChannels();
200     else badAnodeCounter4+=cal->GetDeadChannels();
201     hnbadch->SetBinContent(i+1,cal->GetDeadChannels());
202     if(lay==3) hdeadlay3->SetBinContent(det,lad,cal->GetDeadChannels());
203     if(lay==4) hdeadlay4->SetBinContent(det,lad,cal->GetDeadChannels());
204     for(Int_t iAn=0; iAn<512; iAn++){
205       Int_t ic=cal->GetChip(iAn);
206       if(!cal->IsChipBad(ic) && !cal->IsBad() && cal->IsBadChannel(iAn)){ 
207         if(i<84) badAnodeCounterGoodModAndChip3++;
208         else badAnodeCounterGoodModAndChip4++;
209       }
210       Float_t base=cal->GetBaseline(iAn);
211       Float_t noise=cal->GetNoiseAfterElectronics(iAn);
212       Float_t gain=cal->GetChannelGain(iAn);
213       if(cal->IsBadChannel(iAn)) hchstatus->Fill(0);
214       if(!cal->IsBadChannel(iAn) && !cal->IsChipBad(ic) && !cal->IsBad() ){
215         hbase->Fill(base);
216         hbasemod->Fill(i+240,base);
217         hchstatus->Fill(1);
218         hnoise->Fill(noise);
219         hnoisemod->Fill(i+240,noise);
220         hgain->Fill(gain);
221         hgainmod->Fill(i+240,gain);
222       }
223     }
224   }
225   Int_t totbad3=badModCounter3*512+badChipCounter3*64+badAnodeCounterGoodModAndChip3;
226   Int_t tot3=6*14*512;
227   Float_t fracbad3=(Float_t)totbad3/(Float_t)tot3;
228   Float_t fracgm3=(Float_t)(84.-badModCounter3)/84.;
229   Float_t fracgm4=(Float_t)(176.-badModCounter4)/176.;
230   Float_t fracgh3=(Float_t)(84.*2.-badHybridCounter3)/84./2.;
231   Float_t fracgh4=(Float_t)(176.*2-badHybridCounter4)/176./2.;
232   Float_t fraccgm3=0.;
233   if(badModCounter3!=84){
234     fraccgm3=1.-(Float_t)(badAnodeCounterGoodModAndChip3+badChipCounter3*64)/(512.*(Float_t)(84.-badModCounter3));
235   }
236   Float_t fraccgm4=0.;
237   if(badModCounter4!=176){
238     fraccgm4=1.-(Float_t)(badAnodeCounterGoodModAndChip4+badChipCounter4*64)/(512.*(Float_t)(176.-badModCounter4));
239   }
240   Float_t fraccgh3=0.;
241   if(badHybridCounter3!=(84*2)){
242     fraccgh3=1.-(Float_t)badAnodeCounterGoodHybrid3/(256.*(84.*2.-badHybridCounter3));
243   }
244   Float_t fraccgh4=0.;
245   if(badHybridCounter4!=(176*2)){
246     fraccgh4=1.-(Float_t)badAnodeCounterGoodHybrid4/(256.*(176.*2.-badHybridCounter4));
247   }
248   Int_t totbad4=badModCounter4*512+badChipCounter4*64+badAnodeCounterGoodModAndChip4;
249   Int_t tot4=8*22*512;
250   Float_t fracbad4=(Float_t)totbad4/(Float_t)tot4;
251   Float_t fractot=(Float_t)(totbad3+totbad4)/(Float_t)(tot3+tot4);
252   printf("----------------------Summary----------------------\n");
253   printf("---- Layer 3 ----\n");
254   printf("# of bad modules                      = %d\n",badModCounter3);
255   printf("# of bad chips in good modules        = %d\n",badChipCounter3);
256   printf("# of bad hybrids                      = %d\n",badHybridCounter3);
257   printf("# of bad anodes in good hybrids       = %d\n",badAnodeCounterGoodHybrid3);
258   printf("Fraction of Good modules=%f\n",fracgm3);
259   printf("Fraction of Good hybrids=%f\n",fracgh3);
260   printf("Fraction of good anodes in good modules = %f\n",fraccgm3);
261   printf("Fraction of good anodes in good hybrids = %f\n",fraccgh3);
262   printf("Fraction of bads (anodes+chips+mod)     = %f\n",fracbad3);
263   printf("---- Layer 4 ----\n");
264   printf("# of bad modules                      = %d\n",badModCounter4);
265   printf("# of bad chips in good modules        = %d\n",badChipCounter4);
266   printf("# of bad hybrids                      = %d\n",badHybridCounter4);
267   printf("# of bad anodes in good hybrids       = %d\n",badAnodeCounterGoodHybrid4);
268   printf("Fraction of Good modules=%f\n",fracgm4);
269   printf("Fraction of Good hybrids=%f\n",fracgh4);
270   printf("Fraction of good anodes in good modules = %f\n",fraccgm4);
271   printf("Fraction of good anodes in good hybrids = %f\n",fraccgh4);
272   printf("Fraction of bads (anodes+chips+mod)     = %f\n",fracbad4);
273   printf("---- Total   ----\n");
274   printf("# of bad modules                      = %d\n",badModCounter3+badModCounter4);
275   printf("# of bad chips in good modules        = %d\n",badChipCounter3+badChipCounter4);
276   printf("# of bad anodes in good modules+chips = %d\n",badAnodeCounterGoodModAndChip3+badAnodeCounterGoodModAndChip4);
277   printf("Fraction of bads (anodes+chips+mod)   = %f\n",fractot);
278   printf("---------------------------------------------------\n");
279   
280
281   TLine* lin=new TLine(0,0,0,23);  
282   TExec *ex1 = new TExec("ex1","MakePalette();");
283   TExec *ex2 = new TExec("ex2","gStyle->SetPalette(1);");
284
285   TCanvas* clay=new TCanvas("clay","Layer status",900,600);
286   clay->Divide(2,1);
287   clay->cd(1);
288   hlay3->Draw("col");
289   ex1->Draw();
290   hlay3->DrawCopy("col same");
291   for(Int_t i=0;i<6;i++){
292     lin->SetY1(-0.5);
293     lin->SetY2(13.5);
294     lin->SetX1(i+0.5);
295     lin->SetX2(i+0.5);
296     lin->DrawClone();
297   }
298   for(Int_t i=0;i<14;i++){
299     lin->SetX1(-0.5);
300     lin->SetX2(5.5);
301     lin->SetY1(i+0.5);
302     lin->SetY2(i+0.5);
303     lin->DrawClone();
304   }
305   clay->cd(2);
306   hlay4->DrawCopy("col");
307   for(Int_t i=0;i<8;i++){
308     lin->SetY1(-0.5);
309     lin->SetY2(21.5);
310     lin->SetX1(i+0.5);
311     lin->SetX2(i+0.5);
312     lin->DrawClone();
313   }
314   for(Int_t i=0;i<22;i++){
315     lin->SetX1(-0.5);
316     lin->SetX2(7.5);
317     lin->SetY1(i+0.5);
318     lin->SetY2(i+0.5);
319     lin->DrawClone();
320   }
321
322
323   TCanvas* cddl=new TCanvas("cddl","DDL status",800,800);
324   hddlcarlos->Draw("col");
325   ex1->Draw();
326   hddlcarlos->DrawCopy("col same");
327   for(Int_t i=0;i<12;i++){
328     lin->SetY1(-0.5);
329     lin->SetY2(23.5);
330     lin->SetX1(i+0.5);
331     lin->SetX2(i+0.5);
332     lin->DrawClone();
333   }
334   for(Int_t i=0;i<24;i++){
335     lin->SetX1(-0.5);
336     lin->SetX2(11.5);
337     lin->SetY1(i+0.5);
338     lin->SetY2(i+0.5);
339     lin->DrawClone();
340   }
341  
342
343   TCanvas *c0b=new TCanvas("c0b","Bad Channels",900,600);
344   c0b->Divide(2,1);
345   c0b->cd(1);
346   hdeadlay3->DrawCopy("colz");
347   ex2->Draw();
348   hdeadlay3->DrawCopy("colz same");
349   for(Int_t i=0;i<6;i++){
350     lin->SetY1(-0.5);
351     lin->SetY2(13.5);
352     lin->SetX1(i+0.5);
353     lin->SetX2(i+0.5);
354     lin->DrawClone();
355   }
356   for(Int_t i=0;i<14;i++){
357     lin->SetX1(-0.5);
358     lin->SetX2(5.5);
359     lin->SetY1(i+0.5);
360     lin->SetY2(i+0.5);
361     lin->DrawClone();
362   }
363   c0b->cd(2);
364   hdeadlay4->DrawCopy("colz");
365   ex2->Draw();
366   hdeadlay4->DrawCopy("colz same");
367   for(Int_t i=0;i<8;i++){
368     lin->SetY1(-0.5);
369     lin->SetY2(21.5);
370     lin->SetX1(i+0.5);
371     lin->SetX2(i+0.5);
372     lin->DrawClone();
373   }
374   for(Int_t i=0;i<22;i++){
375     lin->SetX1(-0.5);
376     lin->SetX2(7.5);
377     lin->SetY1(i+0.5);
378     lin->SetY2(i+0.5);
379     lin->DrawClone();
380   }
381
382
383   
384   TCanvas *c1=new TCanvas("c1","Anode calibration",800,800);
385   c1->Divide(2,2);
386   c1->cd(1);
387   hbase->Draw();
388   hbase->GetXaxis()->SetTitle("Baseline after equalization");
389   hbase->GetXaxis()->CenterTitle();  
390   c1->cd(2);
391   hnoise->Draw(); 
392   hnoise->GetXaxis()->SetTitle("Noise");
393   hnoise->GetXaxis()->CenterTitle();  
394   c1->cd(3);
395   hgain->Draw();
396   hgain->GetXaxis()->SetTitle("Gain");
397   hgain->GetXaxis()->CenterTitle();  
398   c1->cd(4);
399   hchstatus->Draw();
400   hchstatus->GetXaxis()->SetTitle("Anode status (0=bad, 1=OK)");
401   hchstatus->GetXaxis()->CenterTitle();
402
403   TCanvas *c1m=new TCanvas("c1m","Calib. vs. mod",1000,800);
404   c1m->Divide(2,2);
405   c1m->cd(1);
406   gPad->SetRightMargin(0.14);
407   hbasemod->SetStats(0);
408   hbasemod->Draw("colz"); 
409   hbasemod->GetXaxis()->SetTitle("Module Number");
410   hbasemod->GetYaxis()->SetTitle("Baseline");
411   c1m->cd(2);
412   gPad->SetRightMargin(0.14);
413   hnoisemod->SetStats(0);
414   hnoisemod->Draw("colz"); 
415   hnoisemod->GetXaxis()->SetTitle("Module Number");
416   hnoisemod->GetYaxis()->SetTitle("Noise");
417   c1m->cd(3);
418   gPad->SetRightMargin(0.14);
419   hgainmod->SetStats(0);
420   hgainmod->Draw("colz");
421   hgainmod->GetXaxis()->SetTitle("Module Number");
422   hgainmod->GetYaxis()->SetTitle("Gain");
423   c1m->cd(4);
424   hnbadch->Scale(1/512.);
425   hnbadch->SetMarkerStyle(20);
426   hnbadch->SetMarkerSize(0.8);
427   hnbadch->SetStats(0);
428   hnbadch->Draw("P");
429   hnbadch->GetXaxis()->SetTitle("Module number");   
430   hnbadch->GetYaxis()->SetTitle("Fraction of bad anodes");
431
432
433
434
435
436
437
438
439   // Plot quantities for specified module
440
441   cal=(AliITSCalibrationSDD*)calSDD->At(iMod);
442   if(cal==0) return;
443   printf("-----------------------------------\n");
444   printf("Module %d    status = ",iMod);
445   if(cal->IsBad()) printf("BAD\n");
446   else printf("OK\n");
447   printf("   Chip Status (0=OK, 1=BAD): ");  
448   for(Int_t ic=0; ic<8;ic++) printf("%d ",cal->IsChipBad(ic));
449   printf("\n");
450   printf("   Number of bad anodes =%d\n",cal->GetDeadChannels());
451   printf("-----------------------------------\n");
452   Int_t ipt=0;
453   TGraph *gbad=new TGraph(0);
454   gbad->SetTitle("Bad Channels");
455   TGraph *gbase=new TGraph(0);
456   gbase->SetTitle("Baselines");
457   TGraph *gnoi=new TGraph(0);
458   gnoi->SetTitle("Noise");
459   TGraph *ggain=new TGraph(0);
460   ggain->SetTitle("Gain");
461   for(Int_t iAn=0; iAn<512; iAn++){
462     Float_t bad=1;
463     if(cal->IsBadChannel(iAn)) bad=0;
464     Float_t base=cal->GetBaseline(iAn);
465     Float_t noise=cal->GetNoiseAfterElectronics(iAn);
466     Float_t gain=cal->GetChannelGain(iAn);
467     gbad->SetPoint(ipt,(Float_t)iAn,bad);
468     gbase->SetPoint(ipt,(Float_t)iAn,base);
469     ggain->SetPoint(ipt,(Float_t)iAn,gain);
470     gnoi->SetPoint(ipt,(Float_t)iAn,noise);
471     ipt++;
472   }
473   Char_t ctit[100];
474   sprintf(ctit,"Module %d",iMod);
475
476   TCanvas *c2=new TCanvas("c2",ctit,1200,800);
477   c2->Divide(2,2);
478   
479   c2->cd(1);
480   gbase->SetMarkerStyle(7);
481   gbase->Draw("AP");
482   gbase->GetXaxis()->SetTitle("Anode Number");
483   gbase->GetYaxis()->SetTitle("Baseline after equalization");  
484   c2->cd(2);
485   gnoi->SetMarkerStyle(7);
486   gnoi->Draw("AP");
487   gnoi->GetXaxis()->SetTitle("Anode Number");
488   gnoi->GetYaxis()->SetTitle("Noise");  
489   c2->cd(3);
490   ggain->SetMarkerStyle(7);
491   ggain->Draw("AP");
492   ggain->GetXaxis()->SetTitle("Anode Number");
493   ggain->GetYaxis()->SetTitle("Gain");
494   c2->cd(4);
495   gbad->SetMarkerStyle(7);
496   gbad->Draw("AP");
497   gbad->SetMinimum(-0.1);
498   gbad->GetXaxis()->SetTitle("Anode Number");
499   gbad->GetYaxis()->SetTitle("Anode Status (1=OK, 0=bad)");
500 }
501
502 void ShowCalibrationSDD(Int_t nrun, Int_t year=2012, Int_t nmod=0){
503   TGrid::Connect("alien:",0,0,"t");
504   TString cmd=Form("gbbox find \"/alice/data/%d/OCDB/ITS/Calib/CalibSDD\" \"Run%d*.root\" > run.txt",year,nrun);
505   gSystem->Exec(cmd.Data());
506   Char_t filnam[200],filnamalien[200];
507   FILE* runtxt=fopen("run.txt","r");
508   Bool_t found=kFALSE;
509   while(!feof(runtxt)){
510     fscanf(runtxt,"%s\n",filnam);    
511     if(strstr(filnam,"/alice/data/")){
512       found=kTRUE;
513       break;
514     }  
515   }
516   if(!found){
517     printf("Bad run number\n");
518     //    gSystem->Exec("rm run.txt");
519     return;
520   }
521   sprintf(filnamalien,"alien://%s",filnam);
522   
523   printf("Open file: %s\n",filnamalien);
524   ShowCalibrationSDD(filnamalien,nmod);
525   fclose(runtxt);
526   gSystem->Exec("rm run.txt");
527 }