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