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