]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/PlotCalibSDDVsTime.C
Bayesian PID: new parametrization and code update (E. Biolcati - F.Prino)
[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 <TObjArray.h>
12 #include "AliCDBEntry.h"
13 #include "AliITSCalibrationSDD.h"
14 #endif
15
16 /*  $Id$    */
17
18 // Macro to plot the calibration parameters from the OCDB files 
19 // created from PEDESTAL and PULSER runs vs. Time
20 // Origin: F. Prino (prino@to.infn.it)
21
22 void PlotCalibSDDVsTime(Int_t firstRun=77677, Int_t lastRun=999999999){
23
24   gStyle->SetOptTitle(0);
25   gStyle->SetOptStat(0);
26   gStyle->SetPadLeftMargin(0.14);
27   gStyle->SetTitleOffset(1.4,"Y");  
28
29
30   TGrid::Connect("alien:",0,0,"t");
31   gSystem->Exec("gbbox find \"/alice/data/2009/OCDB/ITS/Calib/CalibSDD\" \"Run*.root\" > runCalibAlien.txt");
32   FILE* listruns=fopen("runCalibAlien.txt","r");
33
34   TH1F* hbase=new TH1F("hbase","",60,0.5,120.5);
35   TH1F* hnoise=new TH1F("hnoise","",100,0.,7.);
36   TH1F* hgain=new TH1F("hgain","",100,0.,4.);
37   TH1F* hchstatus=new TH1F("hchstatus","",2,-0.5,1.5);
38   TGraphErrors* gbasevstim=new TGraphErrors(0);
39   TGraphErrors* gnoisevstim=new TGraphErrors(0);
40   TGraphErrors* ggainvstim=new TGraphErrors(0);
41   TGraphErrors* gstatvstim=new TGraphErrors(0);
42   gbasevstim->SetName("gbasevstim");
43   gnoisevstim->SetName("gnoisevstim");
44   ggainvstim->SetName("ggainvstim");
45   gstatvstim->SetName("gstatvstim");
46   gbasevstim->SetTitle("Baseline vs. run");
47   gnoisevstim->SetTitle("Noise vs. run");
48   ggainvstim->SetTitle("Gain vs. run");
49   gstatvstim->SetTitle("Good Anodes vs. run");
50
51
52   Char_t filnam[200],filnamalien[200];
53   Int_t iPoint=0;
54   Int_t nrun,nrun2,nv,ns;
55
56   while(!feof(listruns)){
57     hbase->Reset();
58     hnoise->Reset();
59     hgain->Reset();
60     hchstatus->Reset();
61     fscanf(listruns,"%s\n",filnam);    
62     if(!strstr(filnam,"/alice/data/2009")) continue;
63     sscanf(filnam,"/alice/data/2009/OCDB/ITS/Calib/CalibSDD/Run%d_%d_v%d_s%d.root",&nrun,&nrun2,&nv,&ns);
64     if(nrun<85639 && nrun2> 85639) continue; // protection for files with swapped ladders 4-5 of layer 3 
65     if(nrun>100000 && nv< 184) continue; // protection for files with swapped ladder 0-1 of layer 4
66     if(nrun<firstRun) continue;
67     if(nrun>lastRun) continue;
68     sprintf(filnamalien,"alien://%s",filnam);
69     printf("Open file: %s\n",filnam);
70     TFile *f= TFile::Open(filnamalien);  
71     AliCDBEntry *ent=(AliCDBEntry*)f->Get("AliCDBEntry");
72     TObjArray *calSDD = (TObjArray *)ent->GetObject();
73     printf("Run %d Entries in array=%d Entries in Histo=%d\n",nrun,calSDD->GetEntriesFast(),int(hbase->GetEntries()));
74
75
76     AliITSCalibrationSDD *cal;
77     for(Int_t iMod=0; iMod<260;iMod++){
78       cal=(AliITSCalibrationSDD*)calSDD->At(iMod);
79       if(cal==0) continue;
80       for(Int_t iAn=0; iAn<512; iAn++){
81         Int_t ic=cal->GetChip(iAn);
82         Float_t base=cal->GetBaseline(iAn);
83         Float_t noise=cal->GetNoiseAfterElectronics(iAn);
84         Float_t gain=cal->GetChannelGain(iAn);
85         if(cal->IsBadChannel(iAn)) hchstatus->Fill(0);
86         if(!cal->IsBadChannel(iAn) && !cal->IsChipBad(ic) && !cal->IsBad() ){
87           hbase->Fill(base);
88           hchstatus->Fill(1);
89           hnoise->Fill(noise);
90           hgain->Fill(gain);
91         }
92       } 
93     }
94     printf("Set Point %d run %d Base = %f Noise =%f Entries = %d\n",iPoint,nrun,hbase->GetMean(),hnoise->GetMean(),(Int_t)hbase->GetEntries());
95     if((Int_t)hbase->GetEntries()==0) continue;
96     gbasevstim->SetPoint(iPoint,(Double_t)nrun,hbase->GetMean());
97     gbasevstim->SetPointError(iPoint,0.,hbase->GetRMS());
98     gnoisevstim->SetPoint(iPoint,(Double_t)nrun,hnoise->GetMean());
99     gnoisevstim->SetPointError(iPoint,0.,hnoise->GetRMS());
100     ggainvstim->SetPoint(iPoint,(Double_t)nrun,hgain->GetMean());
101     ggainvstim->SetPointError(iPoint,0.,hgain->GetRMS());
102     gstatvstim->SetPoint(iPoint,(Double_t)nrun,hchstatus->GetBinContent(2));
103     iPoint++;
104     f->Close();
105   }
106
107   TFile *ofil=new TFile("Calib2009VsTime.root","recreate");
108   gbasevstim->Write();
109   gnoisevstim->Write();
110   ggainvstim->Write();
111   gstatvstim->Write();
112   ofil->Close();
113
114   TCanvas* cbase=new TCanvas("cbase","Baselines");
115   gbasevstim->SetFillColor(kOrange-2);
116   gbasevstim->SetMarkerStyle(20);
117   gbasevstim->Draw("AP3");
118   gbasevstim->Draw("PLXSAME");
119   gbasevstim->SetMinimum(0.);
120   gbasevstim->SetMaximum(70.);  
121   gbasevstim->GetXaxis()->SetTitle("Run number");
122   gbasevstim->GetYaxis()->SetTitle("<Baseline> (ADC counts)");
123   cbase->SaveAs("BaseRun09.gif");
124
125   TCanvas* cnoise=new TCanvas("cnoise","Noise");
126   gnoisevstim->SetFillColor(kOrange-2);
127   gnoisevstim->SetMarkerStyle(20);
128   gnoisevstim->Draw("AP3");
129   gnoisevstim->Draw("PLXSAME");
130   gnoisevstim->SetMinimum(0.);
131   gnoisevstim->SetMaximum(4.);
132   gnoisevstim->GetXaxis()->SetTitle("Run number");
133   gnoisevstim->GetYaxis()->SetTitle("<Noise> (ADC counts)");
134   cnoise->SaveAs("NoiseRun09.gif");
135
136   TCanvas* cgain=new TCanvas("cgain","Gain");
137   ggainvstim->SetFillColor(kOrange-2);
138   ggainvstim->SetMarkerStyle(20);
139   ggainvstim->Draw("AP3");
140   ggainvstim->Draw("PLXSAME");
141   ggainvstim->SetMinimum(0.);
142   ggainvstim->SetMaximum(4.);
143   ggainvstim->GetXaxis()->SetTitle("Run number");
144   ggainvstim->GetYaxis()->SetTitle("<Gain> (ADC/DAC)");
145   cgain->SaveAs("GainRun09.gif");
146
147   TCanvas* cfrac=new TCanvas("cfrac","Good channels");
148   gstatvstim->SetFillColor(kOrange-2);
149   gstatvstim->SetMarkerStyle(20);
150   gstatvstim->Draw("AP3");
151   gstatvstim->Draw("PLXSAME");
152   gstatvstim->SetMinimum(100000.);
153   gstatvstim->SetMaximum(133000.);
154   gstatvstim->GetXaxis()->SetTitle("Run number");
155   gstatvstim->GetYaxis()->SetTitle("Number of good anodes in acquisition");
156   cfrac->SaveAs("GoodAnodesRun09.gif");
157
158 }