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