]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/PlotCalibSDDVsTime.C
Updated DA from Brigitte. Steerable clock range for pedestal evaluation - default...
[u/mrichter/AliRoot.git] / ITS / PlotCalibSDDVsTime.C
CommitLineData
8d024204 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
22void 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,dum;
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,&dum,&dum);
64 if(nrun<85639 && nrun2> 85639) continue;
65 if(nrun<firstRun) continue;
66 if(nrun>lastRun) continue;
67 sprintf(filnamalien,"alien://%s",filnam);
68 printf("Open file: %s\n",filnam);
69 TFile *f= TFile::Open(filnamalien);
70 AliCDBEntry *ent=(AliCDBEntry*)f->Get("AliCDBEntry");
71 TObjArray *calSDD = (TObjArray *)ent->GetObject();
72 printf("Run %d Entries in array=%d Entries in Histo=%d\n",nrun,calSDD->GetEntriesFast(),int(hbase->GetEntries()));
73
74
75 AliITSCalibrationSDD *cal;
76 for(Int_t iMod=0; iMod<260;iMod++){
77 cal=(AliITSCalibrationSDD*)calSDD->At(iMod);
78 if(cal==0) continue;
79 for(Int_t iAn=0; iAn<512; iAn++){
80 Int_t ic=cal->GetChip(iAn);
81 Float_t base=cal->GetBaseline(iAn);
82 Float_t noise=cal->GetNoiseAfterElectronics(iAn);
83 Float_t gain=cal->GetChannelGain(iAn);
84 if(cal->IsBadChannel(iAn)) hchstatus->Fill(0);
85 if(!cal->IsBadChannel(iAn) && !cal->IsChipBad(ic) && !cal->IsBad() ){
86 hbase->Fill(base);
87 hchstatus->Fill(1);
88 hnoise->Fill(noise);
89 hgain->Fill(gain);
90 }
91 }
92 }
93 printf("Set Point %d run %d Base = %f Noise =%f Entries = %d\n",iPoint,nrun,hbase->GetMean(),hnoise->GetMean(),(Int_t)hbase->GetEntries());
94 if((Int_t)hbase->GetEntries()==0) continue;
95 gbasevstim->SetPoint(iPoint,(Double_t)nrun,hbase->GetMean());
96 gbasevstim->SetPointError(iPoint,0.,hbase->GetRMS());
97 gnoisevstim->SetPoint(iPoint,(Double_t)nrun,hnoise->GetMean());
98 gnoisevstim->SetPointError(iPoint,0.,hnoise->GetRMS());
99 ggainvstim->SetPoint(iPoint,(Double_t)nrun,hgain->GetMean());
100 ggainvstim->SetPointError(iPoint,0.,hgain->GetRMS());
101 gstatvstim->SetPoint(iPoint,(Double_t)nrun,hchstatus->GetBinContent(2));
102 iPoint++;
103 f->Close();
104 }
105
106 TFile *ofil=new TFile("Calib2009VsTime.root","recreate");
107 gbasevstim->Write();
108 gnoisevstim->Write();
109 ggainvstim->Write();
110 gstatvstim->Write();
111 ofil->Close();
112
113 TCanvas* cbase=new TCanvas("cbase","Baselines");
114 gbasevstim->SetFillColor(kOrange-2);
115 gbasevstim->SetMarkerStyle(20);
116 gbasevstim->Draw("AP3");
117 gbasevstim->Draw("PLXSAME");
118 gbasevstim->SetMinimum(0.);
119 gbasevstim->SetMaximum(70.);
120 gbasevstim->GetXaxis()->SetTitle("Run number");
121 gbasevstim->GetYaxis()->SetTitle("<Baseline> (ADC counts)");
122 cbase->SaveAs("BaseRun09.gif");
123
124 TCanvas* cnoise=new TCanvas("cnoise","Noise");
125 gnoisevstim->SetFillColor(kOrange-2);
126 gnoisevstim->SetMarkerStyle(20);
127 gnoisevstim->Draw("AP3");
128 gnoisevstim->Draw("PLXSAME");
129 gnoisevstim->SetMinimum(0.);
130 gnoisevstim->SetMaximum(4.);
131 gnoisevstim->GetXaxis()->SetTitle("Run number");
132 gnoisevstim->GetYaxis()->SetTitle("<Noise> (ADC counts)");
133 cnoise->SaveAs("NoiseRun09.gif");
134
135 TCanvas* cgain=new TCanvas("cgain","Gain");
136 ggainvstim->SetFillColor(kOrange-2);
137 ggainvstim->SetMarkerStyle(20);
138 ggainvstim->Draw("AP3");
139 ggainvstim->Draw("PLXSAME");
140 ggainvstim->SetMinimum(0.);
141 ggainvstim->SetMaximum(4.);
142 ggainvstim->GetXaxis()->SetTitle("Run number");
143 ggainvstim->GetYaxis()->SetTitle("<Gain> (ADC/DAC)");
144 cgain->SaveAs("GainRun09.gif");
145
146 TCanvas* cfrac=new TCanvas("cfrac","Good channels");
147 gstatvstim->SetFillColor(kOrange-2);
148 gstatvstim->SetMarkerStyle(20);
149 gstatvstim->Draw("AP3");
150 gstatvstim->Draw("PLXSAME");
151 gstatvstim->SetMinimum(100000.);
152 gstatvstim->SetMaximum(133000.);
153 gstatvstim->GetXaxis()->SetTitle("Run number");
154 gstatvstim->GetYaxis()->SetTitle("Number of good anodes in acquisition");
155 cfrac->SaveAs("GoodAnodesRun09.gif");
156
157}