]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/ShowDriftSpeedSDD.C
b3716f9b4408d12dc0531b4166542b891b1a1b2a
[u/mrichter/AliRoot.git] / ITS / macrosSDD / 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 iRescaledSpeed=0;
115   Int_t iAverSpeed=0;
116   TLatex* tleft=new TLatex(0.2,0.82,"Side 0");
117   tleft->SetNDC();
118   tleft->SetTextColor(1);
119   TLatex* tright=new TLatex(0.2,0.75,"Side 1");
120   tright->SetNDC();
121   tright->SetTextColor(2);
122
123   for(Int_t i=firstmod; i<lastmod; i++){
124     Int_t iMod=i+240;
125     Int_t lay,lad,det;
126     AliITSgeomTGeo::GetModuleId(iMod,lay,lad,det);
127     Int_t i0=2*i;
128     Int_t i1=1+2*i;
129     vdriftarr0=(AliITSDriftSpeedArraySDD*)drspSDD->At(i0);
130     vdriftarr1=(AliITSDriftSpeedArraySDD*)drspSDD->At(i1);
131     AliITSDriftSpeedSDD* vdrift0=0x0;
132     if(vdriftarr0) vdrift0=vdriftarr0->GetDriftSpeedObject(0);
133     AliITSDriftSpeedSDD* vdrift1=0x0;
134     if(vdriftarr1) vdrift1=vdriftarr1->GetDriftSpeedObject(0);
135
136     gvdr0[i]=new TGraph(0);
137     gvdr1[i]=new TGraph(0);
138     gvdr0[i]->SetMarkerStyle(7);
139     gvdr1[i]->SetMarkerStyle(7);
140     gvdr0[i]->SetMarkerColor(kBlack);
141     gvdr0[i]->SetLineColor(kBlack);
142     gvdr0[i]->SetMinimum(5.);
143     gvdr1[i]->SetMinimum(5.);
144     gvdr0[i]->SetMaximum(7.5);
145     gvdr1[i]->SetMaximum(7.5);
146     sprintf(tit,"Mod %d\n",iMod);
147     gvdr0[i]->SetTitle(tit);
148     gvdr1[i]->SetTitle(tit);
149     gvdr1[i]->SetMarkerColor(kRed);
150     gvdr1[i]->SetLineColor(kRed);
151
152     for(Int_t iAn=0; iAn<256; iAn++){
153       Float_t vel0=0;
154       if(vdrift0) vel0=vdrift0->GetDriftSpeedAtAnode(iAn);
155       Float_t vel1=0;
156       if(vdrift1) vel1=vdrift1->GetDriftSpeedAtAnode(iAn);
157       gvdr0[i]->SetPoint(iAn,(Float_t)iAn,vel0);
158       gvdr1[i]->SetPoint(iAn,(Float_t)iAn,vel1);
159     }
160     Int_t statusInj0=vdriftarr0->GetInjectorStatus();
161     Int_t statusInj1=vdriftarr1->GetInjectorStatus();
162     if(statusInj0>1) iGoodInj++;
163     else if(statusInj0==1) iRescaledSpeed++;
164     else iAverSpeed++;
165     if(statusInj1>0) iGoodInj++;
166     else if(statusInj1==1) iRescaledSpeed++;
167     else iAverSpeed++;
168
169     printf(" Mod. %d \tStatusLR=%X %X \t TimeStamp=%d \t v(an 128l)= %f",iMod,statusInj0,statusInj1,vdrift0->GetEventTimestamp(),vdriftarr0->GetDriftSpeed(0,128));
170     printf("        \t v(an 128r)= %f  Degree=%d %d\n",vdriftarr1->GetDriftSpeed(0,128),vdrift0->GetDegreeofPoly(),vdrift1->GetDegreeofPoly());
171
172     Int_t n7=(statusInj0&(0x1F<<25))>>25;
173     Int_t n6=(statusInj0&(0x1F<<20))>>20;
174     Int_t n5=(statusInj0&(0x1F<<15))>>15;
175     Int_t n4=(statusInj0&(0x1F<<5))>>10;
176     Int_t n3=(statusInj0&(0x1F<<5))>>5;
177     Int_t n2=statusInj0&0x1F;
178     Float_t aveStatus0=(7.*n7+6.*n6+5.*n5+4.*n4+3.*n3+2.*n2)/32.;
179     n7=(statusInj1&(0x1F<<25))>>25;
180     n6=(statusInj1&(0x1F<<20))>>20;
181     n5=(statusInj1&(0x1F<<15))>>15;
182     n4=(statusInj1&(0x1F<<5))>>10;
183     n3=(statusInj1&(0x1F<<5))>>5;
184     n2=statusInj1&0x1F;
185     Float_t aveStatus1=(7.*n7+6.*n6+5.*n5+4.*n4+3.*n3+2.*n2)/32.;
186
187     Int_t index=1+(det-1)*2;
188     if(lay==3){ 
189       hlay3->SetBinContent(index,lad,aveStatus0);
190       hlay3->SetBinContent(index+1,lad,aveStatus1);
191     }
192     if(lay==4){ 
193       hlay4->SetBinContent(index,lad,aveStatus0);
194       hlay4->SetBinContent(index+1,lad,aveStatus1);
195     }
196
197
198     if(!kNoDraw){
199       if (i%12==0 ) {
200         c0->cd();
201         c0->Modified();
202         c0->Update();
203         if (i) c0->Print(psnm1.Data());
204         c0->Clear();
205         c0->Divide(3,4);
206         cntpad = 0;
207       }
208       c0->cd(++cntpad);
209       gvdr0[i]->Draw("AP");
210       gvdr0[i]->GetXaxis()->SetTitle("Anode");
211       gvdr0[i]->GetYaxis()->SetTitle("Vdrift (#mum/ns)");      
212       gvdr1[i]->Draw("P same");
213       gvdr1[i]->GetXaxis()->SetTitle("Anode");
214       gvdr1[i]->GetYaxis()->SetTitle("Vdrift (#mum/ns)");
215       tleft->Draw();
216       tright->Draw();
217     }
218
219     Float_t vel0=0;
220     Float_t pd0=0;
221     if(vdrift0){ 
222       vel0=vdrift0->GetDriftSpeedAtAnode(128);
223      pd0=vdrift0->GetDegreeofPoly();
224     }
225     Float_t vel1=0;
226     Float_t pd1=0;
227     if(vdrift1){ 
228       vel1=vdrift1->GetDriftSpeedAtAnode(128);
229       pd1=vdrift1->GetDegreeofPoly();
230     }
231     Float_t Edrift=(1800-45)/291/0.012;  
232     Float_t mob0=vel0*1.E5/Edrift;  
233     Float_t temper0=293.15*TMath::Power((mob0/1350.),-1/2.4); 
234     Float_t mob1=vel1*1.E5/Edrift;  
235     Float_t temper1=293.15*TMath::Power((mob1/1350.),-1/2.4); 
236     tempvsmod0->SetPoint(tempvsmod0->GetN(),(Float_t)iMod,temper0);
237     tempvsmod1->SetPoint(tempvsmod1->GetN(),(Float_t)iMod,temper1);
238     vvsmod0->SetPoint(vvsmod0->GetN(),(Float_t)iMod,vel0);
239     vvsmod1->SetPoint(vvsmod1->GetN(),(Float_t)iMod,vel1);
240     poldegvsmod0->SetPoint(poldegvsmod0->GetN(),(Float_t)iMod,pd0);
241     poldegvsmod1->SetPoint(poldegvsmod1->GetN(),(Float_t)iMod,pd1);
242
243     for(Int_t ipar=0; ipar<=vdrift0->GetDegreeofPoly(); ipar++){
244       fPoly->SetParameter(ipar,vdrift0->GetDriftSpeedParameter(ipar));
245     }
246     if(vdrift0->GetDegreeofPoly()<3){
247       for(Int_t ipar=vdrift0->GetDegreeofPoly()+1; ipar<=3; ipar++) fPoly->SetParameter(ipar,0.);
248     }
249
250     anmaxvsmod0->SetPoint(anmaxvsmod0->GetN(),(Float_t)iMod,fPoly->GetMaximumX(0.,256.));
251     dvcevsmod0->SetPoint(dvcevsmod0->GetN(),(Float_t)iMod,fPoly->Eval(128)-fPoly->Eval(0));
252     dveevsmod0->SetPoint(dveevsmod0->GetN(),(Float_t)iMod,fPoly->Eval(256)-fPoly->Eval(0));
253     
254     for(Int_t ipar=0; ipar<=vdrift1->GetDegreeofPoly(); ipar++){
255       fPoly->SetParameter(ipar,vdrift1->GetDriftSpeedParameter(ipar));
256     }
257     if(vdrift1->GetDegreeofPoly()<3){
258       for(Int_t ipar=vdrift1->GetDegreeofPoly()+1; ipar<=3; ipar++) fPoly->SetParameter(ipar,0.);
259     }
260     anmaxvsmod1->SetPoint(anmaxvsmod1->GetN(),(Float_t)iMod,fPoly->GetMaximumX(0.,256.));
261     dvcevsmod1->SetPoint(dvcevsmod1->GetN(),(Float_t)iMod,fPoly->Eval(128)-fPoly->Eval(0));
262     dveevsmod1->SetPoint(dveevsmod1->GetN(),(Float_t)iMod,fPoly->Eval(256)-fPoly->Eval(0));
263     //    getchar();
264   }
265
266   if(!kNoDraw){
267     c0->cd();
268     c0->Modified();
269     c0->Update();
270     c0->Print(psnm2.Data());
271   }
272   printf("Number of half-modules with drift speed from injectors               = %d\n",iGoodInj);
273   printf("Number of half-modules with drift speed rewscaled from golden module = %d\n",iRescaledSpeed);
274   printf("Number of half-modules with average drift speed                      = %d\n",iAverSpeed);
275
276   gStyle->SetPalette(59);
277
278   TCanvas* clay=new TCanvas("clay","Injector Status",900,600);
279   clay->Divide(2,1);
280   clay->cd(1);
281   hlay3->Draw("colz");
282   TLine** linv3=new TLine*[5];
283   for(Int_t i=0;i<5;i++){
284     linv3[i]=new TLine(i+0.5,-0.5,i+0.5,13.5);
285     linv3[i]->SetLineColor(kGray+1);
286     linv3[i]->Draw();
287   }
288   TLine** linh3=new TLine*[13];
289   for(Int_t i=0;i<13;i++){
290     linh3[i]=new TLine(-0.5,i+0.5,5.5,i+0.5);
291     linh3[i]->SetLineColor(kGray+1);
292     linh3[i]->Draw();
293   }
294   clay->cd(2);
295   hlay4->Draw("colz");
296   TLine** linv4=new TLine*[7];
297   for(Int_t i=0;i<7;i++){
298     linv4[i]=new TLine(i+0.5,-0.5,i+0.5,21.5);
299     linv4[i]->SetLineColor(kGray+1);
300     linv4[i]->Draw();
301   }
302   TLine** linh4=new TLine*[21];
303   for(Int_t i=0;i<21;i++){
304     linh4[i]=new TLine(-0.5,i+0.5,7.5,i+0.5);
305     linh4[i]->SetLineColor(kGray+1);
306     linh4[i]->Draw();
307   }
308   clay->Modified();
309
310
311   TCanvas* c2;
312   c2=new TCanvas("c2","Vdrift vs. mod",1000,700);
313   vvsmod0->SetMarkerStyle(20);
314   vvsmod0->Draw("AP");
315   vvsmod0->GetXaxis()->SetTitle("Module Number");
316   vvsmod0->GetYaxis()->SetTitle("Vdrift (#mum/ns)");
317   vvsmod1->SetMarkerStyle(21);
318   vvsmod1->SetMarkerColor(2);
319   vvsmod1->Draw("SAMEP");
320   tleft->Draw();
321   tright->Draw();
322
323   TCanvas* c2t;
324   c2t=new TCanvas("c2t","Temper vs. mod",1000,700);
325   tempvsmod0->SetMarkerStyle(20);
326   tempvsmod0->Draw("AP");
327   tempvsmod0->GetXaxis()->SetTitle("Module Number");
328   tempvsmod0->GetYaxis()->SetTitle("Temperature (K)");
329   tempvsmod1->SetMarkerStyle(21);
330   tempvsmod1->SetMarkerColor(2);
331   tempvsmod1->Draw("SAMEP");
332   tleft->Draw();
333   tright->Draw();
334
335   TCanvas* c3;
336   c3=new TCanvas("c3","Params vs. mod",900,900);
337   c3->Divide(2,2);
338   
339   c3->cd(1);
340   gPad->SetLeftMargin(0.14);
341   poldegvsmod0->SetMarkerStyle(20);
342   poldegvsmod0->Draw("AP");
343   poldegvsmod0->GetXaxis()->SetTitle("Module Number");
344   poldegvsmod0->GetYaxis()->SetTitle("Degree of Polynomial fit");
345   poldegvsmod0->GetYaxis()->SetTitleOffset(1.4);
346   poldegvsmod1->SetMarkerStyle(21);
347   poldegvsmod1->SetMarkerColor(2);
348   poldegvsmod1->Draw("SAMEP");
349   tleft->Draw();
350   tright->Draw();
351   c3->cd(2);
352   gPad->SetLeftMargin(0.14);
353   anmaxvsmod0->SetMarkerStyle(20);
354   anmaxvsmod0->Draw("AP");
355   anmaxvsmod0->GetXaxis()->SetTitle("Module Number");
356   anmaxvsmod0->GetYaxis()->SetTitle("Anode with max. drift speed");
357   anmaxvsmod0->GetYaxis()->SetTitleOffset(1.4);
358   anmaxvsmod1->SetMarkerStyle(21);
359   anmaxvsmod1->SetMarkerColor(2);
360   anmaxvsmod1->Draw("SAMEP");
361   tleft->Draw();
362   tright->Draw();
363   c3->cd(3);
364   gPad->SetLeftMargin(0.14);
365   dvcevsmod0->SetMarkerStyle(20);
366   dvcevsmod0->Draw("AP");
367   dvcevsmod0->GetXaxis()->SetTitle("Module Number");
368   dvcevsmod0->GetYaxis()->SetTitle("vdrift(anode128)-vdrift(anode0)");
369   dvcevsmod0->GetYaxis()->SetTitleOffset(1.4);
370   dvcevsmod1->SetMarkerStyle(21);
371   dvcevsmod1->SetMarkerColor(2);
372   dvcevsmod1->Draw("SAMEP");
373   tleft->Draw();
374   tright->Draw();
375   c3->cd(4);
376   gPad->SetLeftMargin(0.14);
377   dveevsmod0->SetMarkerStyle(20);
378   dveevsmod0->Draw("AP");
379   dveevsmod0->GetYaxis()->SetTitleOffset(1.4);
380   dveevsmod0->GetXaxis()->SetTitle("Module Number");
381   dveevsmod0->GetYaxis()->SetTitle("vdrift(anode256)-vdrift(anode0)");
382   dveevsmod1->SetMarkerStyle(21);
383   dveevsmod1->SetMarkerColor(2);
384   dveevsmod1->Draw("SAMEP");
385   tleft->Draw();
386   tright->Draw();
387
388   
389 }
390
391
392
393 void ShowDriftSpeedSDD(Int_t nrun, Int_t year=2011, Int_t nv=-1){
394   TGrid::Connect("alien:",0,0,"t");
395   TString cmd=Form("gbbox find \"/alice/data/%d/OCDB/ITS/Calib/DriftSpeedSDD\" \"Run%d*.root\" > run.txt",year,nrun);
396   if(nv>0){
397     cmd.Form("gbbox find \"/alice/data/%d/OCDB/ITS/Calib/DriftSpeedSDD\" \"Run%d_999999999_v%d_s0.root\" > run.txt",year,nrun,nv);  
398   }
399   gSystem->Exec(cmd.Data());
400   
401   Char_t filnam[200],filnamalien[200];
402   FILE* runtxt=fopen("run.txt","r");
403   fscanf(runtxt,"%s\n",filnam);    
404   if(!strstr(filnam,"/alice/data/")){
405     printf("Bad run number\n");
406     gSystem->Exec("rm run.txt");
407     return;
408   }  
409   sprintf(filnamalien,"alien://%s",filnam);
410   
411   printf("Open file: %s\n",filnamalien);
412   ShowDriftSpeedSDD(filnamalien,0,260,nrun);
413   fclose(runtxt);
414   gSystem->Exec("rm run.txt");
415   
416 }