]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/PlotCalibSDDVsTime.C
First set of modifications to speed-up AliITStrackerSA::SearchClusters (Marco Van...
[u/mrichter/AliRoot.git] / ITS / PlotCalibSDDVsTime.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TH1F.h>
4 #include <TGraph.h>
5 #include <TGraphErrors.h>
6 #include <TStyle.h>
7 #include <TGrid.h>
8 #include <TCanvas.h>
9 #include <TSystem.h>
10 #include <TLatex.h>
11 #include <TLegend.h>
12 #include <TLegendEntry.h>
13 #include <TObjArray.h>
14 #include "AliCDBEntry.h"
15 #include "AliITSCalibrationSDD.h"
16 #endif
17
18 /*  $Id$    */
19
20 // Macro to plot the calibration parameters from the OCDB files 
21 // created from PEDESTAL and PULSER runs vs. Time
22 // Origin: F. Prino (prino@to.infn.it)
23
24 void PlotCalibSDDVsTime(Int_t year=2010, Int_t firstRun=77677, 
25                         Int_t lastRun=999999999){
26
27   gStyle->SetOptTitle(0);
28   gStyle->SetOptStat(0);
29   gStyle->SetPadLeftMargin(0.14);
30   gStyle->SetTitleOffset(1.4,"Y");  
31
32
33   TGrid::Connect("alien:",0,0,"t");
34   gSystem->Exec(Form("gbbox find \"/alice/data/%d/OCDB/ITS/Calib/CalibSDD\" \"Run*.root\" > runCalibAlien.txt",year));
35   FILE* listruns=fopen("runCalibAlien.txt","r");
36
37   TH1F* hbase=new TH1F("hbase","",60,0.5,120.5);
38   TH1F* hnoise=new TH1F("hnoise","",100,0.,7.);
39   TH1F* hgain=new TH1F("hgain","",100,0.,4.);
40   TH1F* hchstatus=new TH1F("hchstatus","",2,-0.5,1.5);
41   TH1F* hchstatus3=new TH1F("hchstatus3","",2,-0.5,1.5);
42   TH1F* hchstatus4=new TH1F("hchstatus4","",2,-0.5,1.5);
43   TGraphErrors* gbasevstim=new TGraphErrors(0);
44   TGraphErrors* gnoisevstim=new TGraphErrors(0);
45   TGraphErrors* ggainvstim=new TGraphErrors(0);
46   TGraphErrors* gstatvstim=new TGraphErrors(0);
47   TGraphErrors* gfracvstim=new TGraphErrors(0);
48   TGraphErrors* gfrac3vstim=new TGraphErrors(0);
49   TGraphErrors* gfrac4vstim=new TGraphErrors(0);
50   gbasevstim->SetName("gbasevstim");
51   gnoisevstim->SetName("gnoisevstim");
52   ggainvstim->SetName("ggainvstim");
53   gstatvstim->SetName("gstatvstim");
54   gfracvstim->SetName("gfracvstim");
55   gfrac3vstim->SetName("gfrac3vstim");
56   gfrac4vstim->SetName("gfrac4vstim");
57   gbasevstim->SetTitle("Baseline vs. run");
58   gnoisevstim->SetTitle("Noise vs. run");
59   ggainvstim->SetTitle("Gain vs. run");
60   gstatvstim->SetTitle("Good Anodes vs. run");
61   gfracvstim->SetTitle("Fraction of Good Anodes vs. run");
62   gfrac3vstim->SetTitle("Fraction of Good Anodes vs. run");
63   gfrac4vstim->SetTitle("Fraction of Good Anodes vs. run");
64
65
66   Char_t filnam[200],filnamalien[200];
67   Int_t iPoint=0;
68   Int_t nrun,nrun2,nv,ns;
69
70   while(!feof(listruns)){
71     hbase->Reset();
72     hnoise->Reset();
73     hgain->Reset();
74     hchstatus->Reset();
75     hchstatus3->Reset();
76     hchstatus4->Reset();
77     fscanf(listruns,"%s\n",filnam);    
78     Char_t directory[100];
79     sprintf(directory,"/alice/data/%d",year);
80     if(!strstr(filnam,directory)) continue;
81     sscanf(filnam,"/alice/data/%d/OCDB/ITS/Calib/CalibSDD/Run%d_%d_v%d_s%d.root",&year,&nrun,&nrun2,&nv,&ns);
82     if(year==2009 && (nrun<85639 && nrun2> 85639)) continue; // protection for files with swapped ladders 4-5 of layer 3 
83     if(year==2009 && (nrun>100000 && nv< 184)) continue; // protection for files with swapped ladder 0-1 of layer 4
84     if(year==2010 && (nrun>=114603 && nv< 98)) continue; // protection for files without treatment of masked hybrids 
85     if(nrun<firstRun) continue;
86     if(nrun>lastRun) continue;
87     sprintf(filnamalien,"alien://%s",filnam);
88     printf("Open file: %s\n",filnam);
89     TFile *f= TFile::Open(filnamalien);  
90     AliCDBEntry *ent=(AliCDBEntry*)f->Get("AliCDBEntry");
91     TObjArray *calSDD = (TObjArray *)ent->GetObject();
92     printf("Run %d Entries in array=%d \n",nrun,calSDD->GetEntriesFast());
93
94
95     AliITSCalibrationSDD *cal;
96     for(Int_t iMod=0; iMod<260;iMod++){
97       cal=(AliITSCalibrationSDD*)calSDD->At(iMod);
98       if(cal==0) continue;
99       for(Int_t iAn=0; iAn<512; iAn++){
100         Int_t ic=cal->GetChip(iAn);
101         Float_t base=cal->GetBaseline(iAn);
102         Float_t noise=cal->GetNoiseAfterElectronics(iAn);
103         Float_t gain=cal->GetChannelGain(iAn);
104         if(cal->IsBadChannel(iAn)){
105           hchstatus->Fill(0);
106           if(iMod<84) hchstatus3->Fill(0);
107           else hchstatus4->Fill(0);
108         }
109         if(!cal->IsBadChannel(iAn) && !cal->IsChipBad(ic) && !cal->IsBad() ){
110           hbase->Fill(base);
111           hchstatus->Fill(1);
112           if(iMod<84) hchstatus3->Fill(1);
113           else hchstatus4->Fill(1);
114           hnoise->Fill(noise);
115           hgain->Fill(gain);
116         }
117       } 
118     }
119     printf("Run %d <Base> = %f <Noise> =%f Entries = %d\n",nrun,hbase->GetMean(),hnoise->GetMean(),(Int_t)hbase->GetEntries());
120     if((Int_t)hbase->GetEntries()==0) continue;
121     gbasevstim->SetPoint(iPoint,(Double_t)nrun,hbase->GetMean());
122     gbasevstim->SetPointError(iPoint,0.,hbase->GetRMS());
123     gnoisevstim->SetPoint(iPoint,(Double_t)nrun,hnoise->GetMean());
124     gnoisevstim->SetPointError(iPoint,0.,hnoise->GetRMS());
125     ggainvstim->SetPoint(iPoint,(Double_t)nrun,hgain->GetMean());
126     ggainvstim->SetPointError(iPoint,0.,hgain->GetRMS());
127     gstatvstim->SetPoint(iPoint,(Double_t)nrun,hchstatus->GetBinContent(2));
128     gfracvstim->SetPoint(iPoint,(Double_t)nrun,hchstatus->GetBinContent(2)/260./512.);
129     gfrac3vstim->SetPoint(iPoint,(Double_t)nrun,hchstatus3->GetBinContent(2)/84./512.);
130     gfrac4vstim->SetPoint(iPoint,(Double_t)nrun,hchstatus4->GetBinContent(2)/176./512.);
131     iPoint++;
132     f->Close();
133   }
134
135   TFile *ofil=new TFile(Form("Calib%dVsTime.root",year),"recreate");
136   gbasevstim->Write();
137   gnoisevstim->Write();
138   ggainvstim->Write();
139   gstatvstim->Write();
140   ofil->Close();
141
142   TCanvas* cbase=new TCanvas("cbase","Baselines");
143   gbasevstim->SetFillColor(kOrange-2);
144   gbasevstim->SetMarkerStyle(20);
145   gbasevstim->Draw("AP3");
146   gbasevstim->Draw("PLXSAME");
147   gbasevstim->SetMinimum(0.);
148   gbasevstim->SetMaximum(70.);  
149   gbasevstim->GetXaxis()->SetTitle("Run number");
150   gbasevstim->GetYaxis()->SetTitle("<Baseline> (ADC counts)");
151   cbase->SaveAs(Form("BaseRun%d.gif",year));
152
153   TCanvas* cnoise=new TCanvas("cnoise","Noise");
154   gnoisevstim->SetFillColor(kOrange-2);
155   gnoisevstim->SetMarkerStyle(20);
156   gnoisevstim->Draw("AP3");
157   gnoisevstim->Draw("PLXSAME");
158   gnoisevstim->SetMinimum(0.);
159   gnoisevstim->SetMaximum(4.);
160   gnoisevstim->GetXaxis()->SetTitle("Run number");
161   gnoisevstim->GetYaxis()->SetTitle("<Noise> (ADC counts)");
162   cnoise->SaveAs(Form("NoiseRun%d.gif",year));
163
164   TCanvas* cgain=new TCanvas("cgain","Gain");
165   ggainvstim->SetFillColor(kOrange-2);
166   ggainvstim->SetMarkerStyle(20);
167   ggainvstim->Draw("AP3");
168   ggainvstim->Draw("PLXSAME");
169   ggainvstim->SetMinimum(0.);
170   ggainvstim->SetMaximum(4.);
171   ggainvstim->GetXaxis()->SetTitle("Run number");
172   ggainvstim->GetYaxis()->SetTitle("<Gain> (ADC/DAC)");
173   cgain->SaveAs(Form("GainRun%d.gif",year));
174
175   TCanvas* cstatus=new TCanvas("cstatus","Good channels");
176   gstatvstim->SetFillColor(kOrange-2);
177   gstatvstim->SetMarkerStyle(20);
178   gstatvstim->Draw("AP3");
179   gstatvstim->Draw("PLXSAME");
180   gstatvstim->SetMinimum(100000.);
181   gstatvstim->SetMaximum(133000.);
182   gstatvstim->GetXaxis()->SetTitle("Run number");
183   gstatvstim->GetYaxis()->SetTitle("Number of good anodes in acquisition");
184   cstatus->SaveAs(Form("GoodAnodesRun%d.gif",year));
185
186   TCanvas* cfrac=new TCanvas("cfrac","Fraction of Good");
187   gfracvstim->SetMarkerStyle(20);
188   gfrac3vstim->SetMarkerStyle(22);
189   gfrac3vstim->SetMarkerColor(2);
190   gfrac3vstim->SetLineColor(2);
191   gfrac4vstim->SetMarkerStyle(23);
192   gfrac4vstim->SetMarkerColor(4);
193   gfrac4vstim->SetLineColor(4);
194   gfracvstim->Draw("APL");
195   gfrac3vstim->Draw("PLSAME");
196   gfrac4vstim->Draw("PLSAME");
197   gfracvstim->SetMinimum(0.7);
198   gfracvstim->SetMaximum(1.05);
199   gfracvstim->GetXaxis()->SetTitle("Run number");
200   gfracvstim->GetYaxis()->SetTitle("Fraction of good anodes in acquisition");
201   TLegend* leg=new TLegend(0.2,0.15,0.45,0.35);
202   leg->SetFillColor(0);
203   TLegendEntry* entr=leg->AddEntry(gfrac3vstim,"Layer 3","P");
204   entr->SetTextColor(2);
205   entr=leg->AddEntry(gfrac4vstim,"Layer 4","P");
206   entr->SetTextColor(4);
207   entr=leg->AddEntry(gfracvstim,"All","P");
208   entr->SetTextColor(1);
209   leg->Draw();
210   cfrac->SaveAs(Form("FractionGoodRun%d.gif",year));
211
212 }