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