]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AnalyzeSDDGainAllMod.C
Bug fix. This eliminates a compilation warning. (R. Shahoyan)
[u/mrichter/AliRoot.git] / ITS / AnalyzeSDDGainAllMod.C
CommitLineData
91282711 1#if !defined(__CINT__) || defined(__MAKECINT__)
2#include <TH2F.h>
3#include <TCanvas.h>
4#include <TStopwatch.h>
5#include <TStyle.h>
6#include <TLatex.h>
7#include <TFile.h>
8#include <TGrid.h>
9#include "AliRawReader.h"
10#include "AliRawReaderDate.h"
11#include "AliRawReaderRoot.h"
12#include "AliITSOnlineSDDBase.h"
13#include "AliITSOnlineSDDCMN.h"
14#include "AliITSOnlineSDDTP.h"
15#include "AliITSRawStreamSDD.h"
16#endif
17
18// Macro for the analysis of PULSER runs (equivalent to ITSSDDGAINda.cxx)
19// Two functions named AnalyzeSDDGainAllModules:
20// The first is for analyzing a local raw data file and takes as agrument the file name.
21// The second is for running on ALIEN
22// All DDLs are analyzed, the argument nDDL selects the DDL to be plotted
23// Origin: F. Prino (prino@to.infn.it)
24
25
26void AnalyzeSDDGainAllMod(Char_t *datafil, Int_t nDDL, Int_t firstEv=10, Int_t lastEv=16, Float_t pascalDAC=100){
27
91282711 28 const Int_t kTotDDL=24;
29 const Int_t kModPerDDL=12;
30 const Int_t kSides=2;
31
32
33 TH2F** histo = new TH2F*[kTotDDL*kModPerDDL*kSides];
34 AliITSOnlineSDDTP **anal=new AliITSOnlineSDDTP*[kTotDDL*kModPerDDL*kSides];
35 Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
36
37 Char_t hisnam[20];
38 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
39 for(Int_t imod=0; imod<kModPerDDL;imod++){
40 for(Int_t isid=0;isid<kSides;isid++){
41 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
42 anal[index]=new AliITSOnlineSDDTP(iddl,imod,isid,pascalDAC);
43 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
44 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
45 isFilled[index]=0;
46 }
47 }
48 }
49
50 TCanvas* c0 = new TCanvas("c0","Ev Display",900,900);
51 gStyle->SetPalette(1);
52 Char_t text[50];
53
54 Int_t iev=firstEv;
55 AliRawReader *rd;
56 if(strstr(datafil,".root")!=0){
57 rd=new AliRawReaderRoot(datafil,iev);
58 }else{
59 rd=new AliRawReaderDate(datafil,iev);
60 }
61 TLatex *t0=new TLatex();
62 t0->SetNDC();
63 t0->SetTextSize(0.06);
64 t0->SetTextColor(4);
65
66 do{
67 c0->Clear();
68 c0->Divide(4,6,0.001,0.001);
69 printf("Event # %d\n",iev);
91282711 70 rd->Reset();
71 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
72 for(Int_t imod=0; imod<kModPerDDL;imod++){
73 for(Int_t isid=0;isid<kSides;isid++){
74 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
75 histo[index]->Reset();
76 }
77 }
78 }
79 AliITSRawStreamSDD s(rd);
91282711 80 while(s.Next()){
81 Int_t iDDL=rd->GetDDLID();
82 Int_t iCarlos=s.GetCarlosId();
de075dae 83 if(s.IsCompletedModule()) continue;
84 if(s.IsCompletedDDL()) continue;
85 if(iDDL>=0 && iDDL<kTotDDL){
91282711 86 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel();
87 histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
88 isFilled[index]=1;
89 }
90 }
91 iev++;
92 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
93 for(Int_t imod=0; imod<kModPerDDL;imod++){
94 for(Int_t isid=0;isid<kSides;isid++){
95 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
96 anal[index]->AddEvent(histo[index]);
97 if(iddl==nDDL){
98 Int_t index2=kSides*imod+isid;
99 c0->cd(index2+1);
100 histo[index]->DrawCopy("colz");
101 sprintf(text,"DDL %d channel %d Side %d",nDDL,imod,isid);
102 t0->DrawLatex(0.15,0.92,text);
103 c0->Update();
104 }
105 }
106 }
107 }
108 printf(" --- OK\n");
109 }while(rd->NextEvent()&&iev<=lastEv);
110
de075dae 111 TH1F *htotgain=new TH1F("htotgain","",100,0.2,4.2);
91282711 112 TH1F *htotpeakpos=new TH1F("htotpeakpos","",256,-0.5,255.5);
113 TH1F *hstatus=new TH1F("hstatus","",2,-0.5,1.5);
114
115 TFile *outfil=new TFile("SDDgain-results.root","recreate");
116 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
117 for(Int_t imod=0; imod<kModPerDDL;imod++){
118 for(Int_t isid=0;isid<kSides;isid++){
119 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
120 if(isFilled[index]){
121 anal[index]->ValidateAnodes();
122 anal[index]->WriteToASCII();
123 anal[index]->WriteToROOT(outfil);
124 for(Int_t ian=0; ian<256;ian++){
125 Float_t gain=anal[index]->GetChannelGain(ian);
126 Float_t ppos=anal[index]->GetTimeBinTPPeak(ian);
127 Int_t anstatus=anal[index]->IsAnodeGood(ian);
128 hstatus->Fill(anstatus);
129 htotgain->Fill(gain);
130 htotpeakpos->Fill(ppos);
131 }
132 }
133 }
134 }
135 }
136 outfil->Close();
137
138 // Draw Statistics of baselines and noise
139 TCanvas *call=new TCanvas("call","General stats",700,700);
140 call->Divide(2,2);
141 call->cd(1);
142 htotpeakpos->Draw();
143 htotpeakpos->GetXaxis()->SetTitle("TP peak position (Time Bin)");
144 htotpeakpos->GetXaxis()->SetTitleSize(0.07);
145 htotpeakpos->GetXaxis()->SetTitleOffset(0.6);
146 call->cd(2);
147 htotgain->Draw();
148 htotgain->GetXaxis()->SetTitle("Gain (ADC/DAC)");
149 htotgain->GetXaxis()->SetTitleSize(0.07);
150 htotgain->GetXaxis()->SetTitleOffset(0.6);
151 call->cd(3);
152 hstatus->Draw();
153 hstatus->GetXaxis()->SetTitle("Anode Status (0=bad 1=good)");
154 hstatus->GetXaxis()->SetTitleSize(0.07);
155 hstatus->GetXaxis()->SetTitleOffset(0.6);
156 call->Update();
157 call->SaveAs("GenStatsPulser.gif");
158
159 // Draw baselines and noisegain and TestPulse time bin for all modules
160
161 TH1F** hgain = new TH1F*[kModPerDDL*kSides];
162 TH1F** htptb = new TH1F*[kModPerDDL*kSides];
163
164 TCanvas *c1=new TCanvas("c1","DDL: TP position",900,900);
165 c1->SetBottomMargin(0.14);
166 c1->Divide(4,6,0.001,0.001);
167 TCanvas *c2=new TCanvas("c2","DDL: gain",900,900);
168 c2->SetBottomMargin(0.14);
169 c2->Divide(4,6,0.001,0.001);
170
171 for(Int_t imod=0; imod<kModPerDDL;imod++){
172 for(Int_t isid=0;isid<kSides;isid++){
173 Int_t index1=kSides*(kModPerDDL*nDDL+imod)+isid;
174 Int_t index2=kSides*imod+isid;
175 sprintf(text,"DDL %d channel %d Side %d",nDDL,imod,isid);
176
177 TLatex *t3=new TLatex(0.15,0.92,text);
178 t3->SetNDC();
179 t3->SetTextSize(0.06);
180 t3->SetTextColor(4);
181 sprintf(hisnam,"hgain%ds%d",imod,isid);
182 hgain[index2]=new TH1F(hisnam,"",256,-0.5,255.5);
183 sprintf(hisnam,"htptb%ds%d",imod,isid);
184 htptb[index2]=new TH1F(hisnam,"",256,-0.5,255.5);
185 for(Int_t ian=0;ian<256;ian++){
186 hgain[index2]->SetBinContent(ian+1,anal[index1]->GetChannelGain(ian));
187 htptb[index2]->SetBinContent(ian+1,anal[index1]->GetTimeBinTPPeak(ian));
188 }
189
190 c1->cd(index2+1);
191 htptb[index2]->Draw();
192 // htptb[imod]->SetMinimum(0);
193 // htptb[imod]->SetMaximum(75);
194 htptb[index2]->GetXaxis()->SetTitle("Anode");
195 htptb[index2]->GetYaxis()->SetTitle("TP position (Time Bin)");
196 htptb[index2]->GetXaxis()->SetTitleSize(0.07);
197 htptb[index2]->GetYaxis()->SetTitleSize(0.07);
198 htptb[index2]->GetXaxis()->SetTitleOffset(0.6);
199 htptb[index2]->GetYaxis()->SetTitleOffset(0.7);
200 t3->Draw();
201 c1->Update();
202
203
204 c2->cd(index2+1);
205 hgain[index2]->SetMinimum(0.);
206 hgain[index2]->SetMaximum(4.);
207 hgain[index2]->Draw();
208 hgain[index2]->GetXaxis()->SetTitle("Anode");
209 hgain[index2]->GetYaxis()->SetTitle("Gain");
210 hgain[index2]->GetXaxis()->SetTitleSize(0.07);
211 hgain[index2]->GetYaxis()->SetTitleSize(0.07);
212 hgain[index2]->GetXaxis()->SetTitleOffset(0.6);
213 hgain[index2]->GetYaxis()->SetTitleOffset(0.7);
214 hgain[index2]->SetStats(0);
215 t3->Draw();
216 c2->Update();
217 }
218 }
219
220 c1->SaveAs("TPtimebin.gif");
221 c2->SaveAs("Gain.gif");
222
223}
224
a7996467 225void AnalyzeSDDGainAllMod(Int_t nrun, Int_t n2, Char_t* dir="LHC08d_SDD", Int_t nDDL=0, Int_t firstEv=15, Int_t lastEv=20, Float_t pascalDAC=100){
91282711 226 TGrid::Connect("alien:",0,0,"t");
227 Char_t filnam[200];
a7996467 228 sprintf(filnam,"alien:///alice/data/2008/%s/%09d/raw/08%09d%03d.10.root",dir,nrun,nrun,n2);
91282711 229 printf("Open file %s\n",filnam);
230 AnalyzeSDDGainAllMod(filnam,nDDL,firstEv,lastEv,pascalDAC);
231}