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