]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ShowDriftSpeedSDD.C
Coverity fix
[u/mrichter/AliRoot.git] / ITS / ShowDriftSpeedSDD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TF1.h>
3 #include <TFile.h>
4 #include <TLine.h>
5 #include <TH1F.h>
6 #include <TH2F.h>
7 #include <TGraph.h>
8 #include <TStyle.h>
9 #include <TSystem.h>
10 #include <TString.h>
11 #include <TGrid.h>
12 #include <TCanvas.h>
13 #include <TLatex.h>
14 #include <TObjArray.h>
15 #include "AliCDBEntry.h"
16 #include "AliITSDriftSpeedArraySDD.h"
17 #include "AliITSDriftSpeedSDD.h"
18 #include "AliITSgeomTGeo.h"
19 #endif
20
21 // Macro to plot the calibration parameters from the OCDB file 
22 // created from an INJECTOR run (OCDB/ITS/Calib/DriftSpeedSDD)
23 // Two methods ShowDriftSpeedSDD:
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 Bool_t kNoDraw = kFALSE; // set to kTRUE to eliminate module dependent plots
30
31
32 void ShowDriftSpeedSDD(Char_t filnam[150]="$ALICE_ROOT/ITS/Calib/DriftSpeedSDD/Run0_9999999_v0_s0.root", Int_t firstmod=0, Int_t lastmod=260,Int_t nrun=0){
33   TFile *f= TFile::Open(filnam);
34   AliCDBEntry *ent=(AliCDBEntry*)f->Get("AliCDBEntry");
35   TObjArray *drspSDD = (TObjArray *)ent->GetObject();
36   AliITSDriftSpeedArraySDD *vdriftarr0;
37   AliITSDriftSpeedArraySDD *vdriftarr1;
38
39   TH2F* hlay3=new TH2F("hlay3","Injector Status 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(-0.01);
46   hlay3->SetMaximum(7.);
47   TH2F* hlay4=new TH2F("hlay4","Injector Status Layer 4",16,-0.5,7.5,22,-0.5,21.5);
48   hlay4->GetXaxis()->SetTitle("Detector");
49   hlay4->GetYaxis()->SetTitle("Ladder");
50   hlay4->GetXaxis()->SetTickLength(0);
51   hlay4->GetYaxis()->SetTickLength(0);
52   hlay4->GetYaxis()->SetTitle("Ladder");
53   hlay4->SetStats(0);
54   hlay4->SetMinimum(-0.01);
55   hlay4->SetMaximum(7.);
56
57   TGraph **gvdr0=new TGraph*[260];
58   TGraph **gvdr1=new TGraph*[260];
59   TCanvas *c0=0x0;
60   if(!kNoDraw)c0=new TCanvas("c0","Module Drift Speed",700,1000);
61
62   TGraph *tempvsmod0=new TGraph(0);
63   TGraph *tempvsmod1=new TGraph(0);
64   TGraph *vvsmod0=new TGraph(0);
65   TGraph *vvsmod1=new TGraph(0);
66   TGraph *poldegvsmod0=new TGraph(0); 
67   TGraph *poldegvsmod1=new TGraph(0); 
68   TGraph *anmaxvsmod0=new TGraph(0); 
69   TGraph *anmaxvsmod1=new TGraph(0); 
70   TGraph *dvcevsmod0=new TGraph(0);
71   TGraph *dvcevsmod1=new TGraph(0);
72   TGraph *dveevsmod0=new TGraph(0);
73   TGraph *dveevsmod1=new TGraph(0);
74
75   char tit0[100];
76   sprintf(tit0,"Temperature vs. mod. number");
77   if(nrun!=0)sprintf(tit0,"Temperature vs. mod. number - Run %d",nrun);
78   tempvsmod0->SetTitle(tit0);
79   tempvsmod1->SetTitle(tit0);
80
81   sprintf(tit0,"Drift Speed vs. mod. number");
82   if(nrun!=0)sprintf(tit0,"Drift Speed vs. mod. number - Run %d",nrun);
83   vvsmod0->SetTitle(tit0);
84   vvsmod1->SetTitle(tit0);
85
86   sprintf(tit0,"Degree of poly fit vs. mod. number");
87   if(nrun!=0)sprintf(tit0,"Degree of poly fit vs. mod. number - Run %d",nrun);
88   poldegvsmod0->SetTitle(tit0);
89   poldegvsmod1->SetTitle(tit0);
90
91   sprintf(tit0,"Anode with max. vdrift vs. mod. number");
92   if(nrun!=0)sprintf(tit0,"Anode with max. vdrift vs. mod. number - Run %d",nrun);
93   anmaxvsmod0->SetTitle(tit0);
94   anmaxvsmod1->SetTitle(tit0);
95
96   sprintf(tit0,"Delta Vdrift 128-0 vs. mod. number");
97   if(nrun!=0)sprintf(tit0,"Delta Vdrift 128-0 vs. mod. number - Run %d",nrun);
98   dvcevsmod0->SetTitle(tit0);
99   dvcevsmod1->SetTitle(tit0);
100
101   sprintf(tit0,"Delta Vdrift 256-0 vs. mod. number");
102   if(nrun!=0)sprintf(tit0,"Delta Vdrift 256-0 vs. mod. number - Run %d",nrun);
103   dveevsmod0->SetTitle(tit0);
104   dveevsmod1->SetTitle(tit0);
105
106   TF1* fPoly=new TF1("fPoly","pol3",0.,256.);
107   Char_t tit[100];
108   TString psnm0 = "vdriftSDD.ps[";
109   TString psnm1 = "vdriftSDD.ps";
110   TString psnm2 = "vdriftSDD.ps]";
111   if(!kNoDraw) c0->Print(psnm0.Data());
112   Int_t cntpad = 0;
113   Int_t iGoodInj=0;
114   Int_t iAverSpeed=0;
115   TLatex* tleft=new TLatex(0.2,0.82,"Side 0");
116   tleft->SetNDC();
117   tleft->SetTextColor(1);
118   TLatex* tright=new TLatex(0.2,0.75,"Side 1");
119   tright->SetNDC();
120   tright->SetTextColor(2);
121
122   for(Int_t i=firstmod; i<lastmod; i++){
123     Int_t iMod=i+240;
124     Int_t lay,lad,det;
125     AliITSgeomTGeo::GetModuleId(iMod,lay,lad,det);
126     Int_t i0=2*i;
127     Int_t i1=1+2*i;
128     vdriftarr0=(AliITSDriftSpeedArraySDD*)drspSDD->At(i0);
129     vdriftarr1=(AliITSDriftSpeedArraySDD*)drspSDD->At(i1);
130     AliITSDriftSpeedSDD* vdrift0=0x0;
131     if(vdriftarr0) vdrift0=vdriftarr0->GetDriftSpeedObject(0);
132     AliITSDriftSpeedSDD* vdrift1=0x0;
133     if(vdriftarr1) vdrift1=vdriftarr1->GetDriftSpeedObject(0);
134
135     gvdr0[i]=new TGraph(0);
136     gvdr1[i]=new TGraph(0);
137     gvdr0[i]->SetMarkerStyle(7);
138     gvdr1[i]->SetMarkerStyle(7);
139     gvdr0[i]->SetMarkerColor(kBlack);
140     gvdr0[i]->SetLineColor(kBlack);
141     gvdr0[i]->SetMinimum(5.);
142     gvdr1[i]->SetMinimum(5.);
143     gvdr0[i]->SetMaximum(7.5);
144     gvdr1[i]->SetMaximum(7.5);
145     sprintf(tit,"Mod %d\n",iMod);
146     gvdr0[i]->SetTitle(tit);
147     gvdr1[i]->SetTitle(tit);
148     gvdr1[i]->SetMarkerColor(kRed);
149     gvdr1[i]->SetLineColor(kRed);
150
151     for(Int_t iAn=0; iAn<256; iAn++){
152       Float_t vel0=0;
153       if(vdrift0) vel0=vdrift0->GetDriftSpeedAtAnode(iAn);
154       Float_t vel1=0;
155       if(vdrift1) vel1=vdrift1->GetDriftSpeedAtAnode(iAn);
156       gvdr0[i]->SetPoint(iAn,(Float_t)iAn,vel0);
157       gvdr1[i]->SetPoint(iAn,(Float_t)iAn,vel1);
158     }
159     Int_t statusInj0=vdriftarr0->GetInjectorStatus();
160     Int_t statusInj1=vdriftarr1->GetInjectorStatus();
161     if(statusInj0>0) iGoodInj++;
162     else iAverSpeed++;
163     if(statusInj1>0) iGoodInj++;
164     else iAverSpeed++;
165
166     printf(" Mod. %d \tStatusLR=%X %X \t v(an 128l)= %f",iMod,statusInj0,statusInj1,vdriftarr0->GetDriftSpeed(0,128));
167     printf("        \t v(an 128r)= %f  Degree=%d %d\n",vdriftarr1->GetDriftSpeed(0,128),vdrift0->GetDegreeofPoly(),vdrift1->GetDegreeofPoly());
168
169     Int_t n7=(statusInj0&(0x1F<<25))>>25;
170     Int_t n6=(statusInj0&(0x1F<<20))>>20;
171     Int_t n5=(statusInj0&(0x1F<<15))>>15;
172     Int_t n4=(statusInj0&(0x1F<<5))>>10;
173     Int_t n3=(statusInj0&(0x1F<<5))>>5;
174     Int_t n2=statusInj0&0x1F;
175     Float_t aveStatus0=(7.*n7+6.*n6+5.*n5+4.*n4+3.*n3+2.*n2)/32.;
176     n7=(statusInj1&(0x1F<<25))>>25;
177     n6=(statusInj1&(0x1F<<20))>>20;
178     n5=(statusInj1&(0x1F<<15))>>15;
179     n4=(statusInj1&(0x1F<<5))>>10;
180     n3=(statusInj1&(0x1F<<5))>>5;
181     n2=statusInj1&0x1F;
182     Float_t aveStatus1=(7.*n7+6.*n6+5.*n5+4.*n4+3.*n3+2.*n2)/32.;
183
184     Int_t index=1+(det-1)*2;
185     if(lay==3){ 
186       hlay3->SetBinContent(index,lad,aveStatus0);
187       hlay3->SetBinContent(index+1,lad,aveStatus1);
188     }
189     if(lay==4){ 
190       hlay4->SetBinContent(index,lad,aveStatus0);
191       hlay4->SetBinContent(index+1,lad,aveStatus1);
192     }
193
194
195     if(!kNoDraw){
196       if (i%12==0 ) {
197         c0->cd();
198         c0->Modified();
199         c0->Update();
200         if (i) c0->Print(psnm1.Data());
201         c0->Clear();
202         c0->Divide(3,4);
203         cntpad = 0;
204       }
205       c0->cd(++cntpad);
206       gvdr0[i]->Draw("AP");
207       gvdr0[i]->GetXaxis()->SetTitle("Anode");
208       gvdr0[i]->GetYaxis()->SetTitle("Vdrift (#mum/ns)");      
209       gvdr1[i]->Draw("P same");
210       gvdr1[i]->GetXaxis()->SetTitle("Anode");
211       gvdr1[i]->GetYaxis()->SetTitle("Vdrift (#mum/ns)");
212       tleft->Draw();
213       tright->Draw();
214     }
215
216     Float_t vel0=0;
217     Float_t pd0=0;
218     if(vdrift0){ 
219       vel0=vdrift0->GetDriftSpeedAtAnode(128);
220      pd0=vdrift0->GetDegreeofPoly();
221     }
222     Float_t vel1=0;
223     Float_t pd1=0;
224     if(vdrift1){ 
225       vel1=vdrift1->GetDriftSpeedAtAnode(128);
226       pd1=vdrift1->GetDegreeofPoly();
227     }
228     Float_t Edrift=(1800-45)/291/0.012;  
229     Float_t mob0=vel0*1.E5/Edrift;  
230     Float_t temper0=293.15*TMath::Power((mob0/1350.),-1/2.4); 
231     Float_t mob1=vel1*1.E5/Edrift;  
232     Float_t temper1=293.15*TMath::Power((mob1/1350.),-1/2.4); 
233     tempvsmod0->SetPoint(tempvsmod0->GetN(),(Float_t)iMod,temper0);
234     tempvsmod1->SetPoint(tempvsmod1->GetN(),(Float_t)iMod,temper1);
235     vvsmod0->SetPoint(vvsmod0->GetN(),(Float_t)iMod,vel0);
236     vvsmod1->SetPoint(vvsmod1->GetN(),(Float_t)iMod,vel1);
237     poldegvsmod0->SetPoint(poldegvsmod0->GetN(),(Float_t)iMod,pd0);
238     poldegvsmod1->SetPoint(poldegvsmod1->GetN(),(Float_t)iMod,pd1);
239
240     for(Int_t ipar=0; ipar<=vdrift0->GetDegreeofPoly(); ipar++){
241       fPoly->SetParameter(ipar,vdrift0->GetDriftSpeedParameter(ipar));
242     }
243     if(vdrift0->GetDegreeofPoly()<3){
244       for(Int_t ipar=vdrift0->GetDegreeofPoly()+1; ipar<=3; ipar++) fPoly->SetParameter(ipar,0.);
245     }
246
247     anmaxvsmod0->SetPoint(anmaxvsmod0->GetN(),(Float_t)iMod,fPoly->GetMaximumX(0.,256.));
248     dvcevsmod0->SetPoint(dvcevsmod0->GetN(),(Float_t)iMod,fPoly->Eval(128)-fPoly->Eval(0));
249     dveevsmod0->SetPoint(dveevsmod0->GetN(),(Float_t)iMod,fPoly->Eval(256)-fPoly->Eval(0));
250     
251     for(Int_t ipar=0; ipar<=vdrift1->GetDegreeofPoly(); ipar++){
252       fPoly->SetParameter(ipar,vdrift1->GetDriftSpeedParameter(ipar));
253     }
254     if(vdrift1->GetDegreeofPoly()<3){
255       for(Int_t ipar=vdrift1->GetDegreeofPoly()+1; ipar<=3; ipar++) fPoly->SetParameter(ipar,0.);
256     }
257     anmaxvsmod1->SetPoint(anmaxvsmod1->GetN(),(Float_t)iMod,fPoly->GetMaximumX(0.,256.));
258     dvcevsmod1->SetPoint(dvcevsmod1->GetN(),(Float_t)iMod,fPoly->Eval(128)-fPoly->Eval(0));
259     dveevsmod1->SetPoint(dveevsmod1->GetN(),(Float_t)iMod,fPoly->Eval(256)-fPoly->Eval(0));
260     //    getchar();
261   }
262
263   if(!kNoDraw){
264     c0->cd();
265     c0->Modified();
266     c0->Update();
267     c0->Print(psnm2.Data());
268   }
269   printf("Number of half-modules with drift speed from injectors = %d\n",iGoodInj);
270   printf("Number of half-modules with average drift speed        = %d\n",iAverSpeed);
271
272   gStyle->SetPalette(59);
273
274   TCanvas* clay=new TCanvas("clay","Injector Status",900,600);
275   clay->Divide(2,1);
276   clay->cd(1);
277   hlay3->Draw("colz");
278   TLine** linv3=new TLine*[5];
279   for(Int_t i=0;i<5;i++){
280     linv3[i]=new TLine(i+0.5,-0.5,i+0.5,13.5);
281     linv3[i]->SetLineColor(kGray+1);
282     linv3[i]->Draw();
283   }
284   TLine** linh3=new TLine*[13];
285   for(Int_t i=0;i<13;i++){
286     linh3[i]=new TLine(-0.5,i+0.5,5.5,i+0.5);
287     linh3[i]->SetLineColor(kGray+1);
288     linh3[i]->Draw();
289   }
290   clay->cd(2);
291   hlay4->Draw("colz");
292   TLine** linv4=new TLine*[7];
293   for(Int_t i=0;i<7;i++){
294     linv4[i]=new TLine(i+0.5,-0.5,i+0.5,21.5);
295     linv4[i]->SetLineColor(kGray+1);
296     linv4[i]->Draw();
297   }
298   TLine** linh4=new TLine*[21];
299   for(Int_t i=0;i<21;i++){
300     linh4[i]=new TLine(-0.5,i+0.5,7.5,i+0.5);
301     linh4[i]->SetLineColor(kGray+1);
302     linh4[i]->Draw();
303   }
304   clay->Modified();
305
306
307   TCanvas* c2;
308   c2=new TCanvas("c2","Vdrift vs. mod",1000,700);
309   vvsmod0->SetMarkerStyle(20);
310   vvsmod0->Draw("AP");
311   vvsmod0->GetXaxis()->SetTitle("Module Number");
312   vvsmod0->GetYaxis()->SetTitle("Vdrift (#mum/ns)");
313   vvsmod1->SetMarkerStyle(21);
314   vvsmod1->SetMarkerColor(2);
315   vvsmod1->Draw("SAMEP");
316   tleft->Draw();
317   tright->Draw();
318
319   TCanvas* c2t;
320   c2t=new TCanvas("c2t","Temper vs. mod",1000,700);
321   tempvsmod0->SetMarkerStyle(20);
322   tempvsmod0->Draw("AP");
323   tempvsmod0->GetXaxis()->SetTitle("Module Number");
324   tempvsmod0->GetYaxis()->SetTitle("Temperature (K)");
325   tempvsmod1->SetMarkerStyle(21);
326   tempvsmod1->SetMarkerColor(2);
327   tempvsmod1->Draw("SAMEP");
328   tleft->Draw();
329   tright->Draw();
330
331   TCanvas* c3;
332   c3=new TCanvas("c3","Params vs. mod",900,900);
333   c3->Divide(2,2);
334   
335   c3->cd(1);
336   gPad->SetLeftMargin(0.14);
337   poldegvsmod0->SetMarkerStyle(20);
338   poldegvsmod0->Draw("AP");
339   poldegvsmod0->GetXaxis()->SetTitle("Module Number");
340   poldegvsmod0->GetYaxis()->SetTitle("Degree of Polynomial fit");
341   poldegvsmod0->GetYaxis()->SetTitleOffset(1.4);
342   poldegvsmod1->SetMarkerStyle(21);
343   poldegvsmod1->SetMarkerColor(2);
344   poldegvsmod1->Draw("SAMEP");
345   tleft->Draw();
346   tright->Draw();
347   c3->cd(2);
348   gPad->SetLeftMargin(0.14);
349   anmaxvsmod0->SetMarkerStyle(20);
350   anmaxvsmod0->Draw("AP");
351   anmaxvsmod0->GetXaxis()->SetTitle("Module Number");
352   anmaxvsmod0->GetYaxis()->SetTitle("Anode with max. drift speed");
353   anmaxvsmod0->GetYaxis()->SetTitleOffset(1.4);
354   anmaxvsmod1->SetMarkerStyle(21);
355   anmaxvsmod1->SetMarkerColor(2);
356   anmaxvsmod1->Draw("SAMEP");
357   tleft->Draw();
358   tright->Draw();
359   c3->cd(3);
360   gPad->SetLeftMargin(0.14);
361   dvcevsmod0->SetMarkerStyle(20);
362   dvcevsmod0->Draw("AP");
363   dvcevsmod0->GetXaxis()->SetTitle("Module Number");
364   dvcevsmod0->GetYaxis()->SetTitle("vdrift(anode128)-vdrift(anode0)");
365   dvcevsmod0->GetYaxis()->SetTitleOffset(1.4);
366   dvcevsmod1->SetMarkerStyle(21);
367   dvcevsmod1->SetMarkerColor(2);
368   dvcevsmod1->Draw("SAMEP");
369   tleft->Draw();
370   tright->Draw();
371   c3->cd(4);
372   gPad->SetLeftMargin(0.14);
373   dveevsmod0->SetMarkerStyle(20);
374   dveevsmod0->Draw("AP");
375   dveevsmod0->GetYaxis()->SetTitleOffset(1.4);
376   dveevsmod0->GetXaxis()->SetTitle("Module Number");
377   dveevsmod0->GetYaxis()->SetTitle("vdrift(anode256)-vdrift(anode0)");
378   dveevsmod1->SetMarkerStyle(21);
379   dveevsmod1->SetMarkerColor(2);
380   dveevsmod1->Draw("SAMEP");
381   tleft->Draw();
382   tright->Draw();
383
384   
385 }
386
387
388
389 void ShowDriftSpeedSDD(Int_t nrun, Int_t year=2010, Int_t nv=-1){
390   TGrid::Connect("alien:",0,0,"t");
391   TString cmd=Form("gbbox find \"/alice/data/%d/OCDB/ITS/Calib/DriftSpeedSDD\" \"Run%d*.root\" > run.txt",year,nrun);
392   if(nv>0){
393     cmd.Form("gbbox find \"/alice/data/%d/OCDB/ITS/Calib/DriftSpeedSDD\" \"Run%d_999999999_v%d_s0.root\" > run.txt",year,nrun,nv);  
394   }
395   gSystem->Exec(cmd.Data());
396   
397   Char_t filnam[200],filnamalien[200];
398   FILE* runtxt=fopen("run.txt","r");
399   fscanf(runtxt,"%s\n",filnam);    
400   if(!strstr(filnam,"/alice/data/")){
401     printf("Bad run number\n");
402     gSystem->Exec("rm run.txt");
403     return;
404   }  
405   sprintf(filnamalien,"alien://%s",filnam);
406   
407   printf("Open file: %s\n",filnamalien);
408   ShowDriftSpeedSDD(filnamalien,0,260,nrun);
409   fclose(runtxt);
410   gSystem->Exec("rm run.txt");
411   
412 }