]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AnalyzeSDDInjectorsAllMod.C
Adding output data type of rootifier component to the constants class.
[u/mrichter/AliRoot.git] / ITS / AnalyzeSDDInjectorsAllMod.C
CommitLineData
91282711 1#if !defined(__CINT__) || defined(__MAKECINT__)
2#include <TH2F.h>
3#include <TCanvas.h>
4#include <TGraphErrors.h>
5#include <TStopwatch.h>
6#include <TStyle.h>
7#include <TLatex.h>
8#include <TFile.h>
9#include <TGrid.h>
10#include <TF1.h>
11#include <TLine.h>
12#include "AliRawReader.h"
13#include "AliRawReaderDate.h"
14#include "AliRawReaderRoot.h"
15#include "AliITSOnlineSDDInjectors.h"
16#include "AliITSRawStreamSDD.h"
17#include "AliITSDDLModuleMapSDD.h"
18#endif
19
20// Macro for the analysis of PULSER runs (equivalent to ITSSDDINJda.cxx)
21// Two functions named AnalyzeSDDInjectorsAllModules:
22// The first is for analyzing a local raw data file and takes as agrument the file name.
23// The second is for running on ALIEN
24// All DDLs are analyzed, the argument nDDL selects the DDL to be plotted
25// Origin: F. Prino (prino@to.infn.it)
26
27
28void AnalyzeSDDInjectorsAllMod(Char_t *datafil, Int_t nDDL, Int_t firstEv=10, Int_t lastEv=15){
5a4e63d5 29
91282711 30 const Int_t kTotDDL=24;
31 const Int_t kModPerDDL=12;
32 const Int_t kSides=2;
33
34 AliITSDDLModuleMapSDD* dmap=new AliITSDDLModuleMapSDD();
35 dmap->SetJun08Map();
36
37 TH2F** histo = new TH2F*[kTotDDL*kModPerDDL*kSides];
38 Int_t nWrittenEv[kTotDDL*kModPerDDL*kSides];
39 TGraphErrors** gvel = new TGraphErrors*[kTotDDL*kModPerDDL*kSides];
40 AliITSOnlineSDDInjectors **anal=new AliITSOnlineSDDInjectors*[kTotDDL*kModPerDDL*kSides];
41
42 Char_t hisnam[20];
43 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
44 for(Int_t imod=0; imod<kModPerDDL;imod++){
45 for(Int_t isid=0;isid<kSides;isid++){
46 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
47 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
48 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
49 anal[index]=new AliITSOnlineSDDInjectors(iddl,imod,isid);
5a4e63d5 50/* Uncomment these lines for analysis of runs with 40 MHz sapling */
51// anal[index]->SetInjLineRange(0,20,50);
52// anal[index]->SetInjLineRange(1,90,160);
53// anal[index]->SetInjLineRange(2,170,240);
54// anal[index]->SetTimeStep(25.);
55/* END of lines to be uncommented */
91282711 56 nWrittenEv[index]=0;
57 }
58 }
59 }
60 TGraph *gvvsmod0=new TGraph(0);
61 TGraph *gvvsmod1=new TGraph(0);
62
63 TCanvas* c0 = new TCanvas("c0","",900,900);
64 gStyle->SetPalette(1);
65 TCanvas* c1 = new TCanvas("c1","",900,900);
66 Char_t text[50];
67
68 Int_t iev=firstEv;
69 AliRawReader *rd;
70 if(strstr(datafil,".root")!=0){
71 rd=new AliRawReaderRoot(datafil,iev);
72 }else{
73 rd=new AliRawReaderDate(datafil,iev);
74 }
75
76 Char_t gname[15];
77 TF1 *funz=new TF1("funz","[0]+[1]*x+[2]*x*x+[3]*x*x*x",0.,255.);
78 TLatex *t0=new TLatex();
79 t0->SetNDC();
80 t0->SetTextSize(0.06);
81 t0->SetTextColor(4);
82
83 do{
84 c0->Clear();
85 c0->Divide(4,6,0.001,0.001);
86 c1->Clear();
87 c1->Divide(4,6,0.001,0.001);
88 printf("Event # %d\n",iev);
89 UInt_t timeSt=rd->GetTimestamp();
91282711 90 rd->Reset();
91 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
92 for(Int_t imod=0; imod<kModPerDDL;imod++){
93 for(Int_t isid=0;isid<kSides;isid++){
94 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
95 histo[index]->Reset();
96 }
97 }
98 }
99
100 AliITSRawStreamSDD s(rd);
91282711 101 while(s.Next()){
102 Int_t iDDL=rd->GetDDLID();
103 Int_t iCarlos=s.GetCarlosId();
104 if(iDDL>=0 && iDDL<kTotDDL && s.IsCompletedModule()==kFALSE){
105 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel();
106 histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
107 }
108 }
109
110 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
111 for(Int_t imod=0; imod<kModPerDDL;imod++){
112 for(Int_t isid=0;isid<kSides;isid++){
113 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
114 anal[index]->AnalyzeEvent(histo[index]);
115 anal[index]->WriteToASCII(iev,timeSt,nWrittenEv[index]);
116 nWrittenEv[index]++;
117 if(iev==firstEv && anal[index]->GetInjPadStatus(16)>=6){
118 Float_t vel=anal[index]->GetDriftSpeed(16);
119 Int_t iMod=dmap->GetModuleNumber(iddl,imod);
120 if(isid==0) gvvsmod0->SetPoint(gvvsmod0->GetN(),(Float_t)iMod,vel);
121 if(isid==1) gvvsmod1->SetPoint(gvvsmod1->GetN(),(Float_t)iMod,vel);
122 }
123 if(iddl==nDDL){
124 Int_t index2=kSides*imod+isid;
125 c0->cd(index2+1);
126 histo[index]->SetMaximum(100.);
127 histo[index]->DrawCopy("colz");
128 sprintf(text,"DDL %d channel %d Side %d",nDDL,imod,isid);
129 t0->DrawLatex(0.15,0.92,text);
130 c0->Update();
131 c1->cd(index2+1);
132 gvel[index]=anal[index]->GetDriftSpeedGraph();
133 gvel[index]->SetMarkerStyle(20);
134 gvel[index]->SetTitle("");
135 sprintf(gname,"gvel%dev%d",index,iev);
136 gvel[index]->SetName(gname);
137 // gvel[index]->SetMinimum(0);
138 //gvel[index]->SetMaximum(200);
139
140 gvel[index]->GetXaxis()->SetLimits(0,256);
141 gvel[index]->GetXaxis()->SetTitle("Anode");
142 gvel[index]->GetYaxis()->SetTitle("Drift vel.");
143 gvel[index]->GetXaxis()->SetTitleSize(0.07);
144 gvel[index]->GetYaxis()->SetTitleSize(0.07);
145 gvel[index]->GetXaxis()->SetTitleOffset(0.6);
146 gvel[index]->GetYaxis()->SetTitleOffset(0.6);
147 if(gvel[index]->GetN()>0) gvel[index]->Draw("AP");
148 Float_t *param=anal[index]->GetDriftSpeedFitParam();
149 funz->SetParameters(param[0],param[1],param[2],param[3]);
150 funz->SetLineColor(2);
151 funz->DrawCopy("LSAME");
152 t0->DrawLatex(0.15,0.92,text);
153 c1->Update();
154 }
155 }
156 }
157 }
158 iev++;
159 printf(" --- OK\n");
160 }while(rd->NextEvent()&&iev<=lastEv);
161
162 TCanvas* c8=new TCanvas("c8");
163 gvvsmod0->SetTitle("");
164 gvvsmod1->SetTitle("");
165
166 gvvsmod0->SetMarkerStyle(20);
167 gvvsmod1->SetMarkerStyle(21);
168 gvvsmod1->SetMarkerColor(2);
169 gvvsmod0->Draw("AP");
170 gvvsmod0->SetMinimum(6.2);
171 gvvsmod0->SetMaximum(7.2);
172 gvvsmod0->GetXaxis()->SetTitle("Module Number");
173 gvvsmod0->GetYaxis()->SetTitle("Vdrift at injector pad 16");
174 gvvsmod1->Draw("PSAME");
175 TLatex* tleft=new TLatex(0.7,0.82,"Side 0");
176 tleft->SetNDC();
177 tleft->SetTextColor(1);
178 tleft->Draw();
179 TLatex* tright=new TLatex(0.7,0.75,"Side 1");
180 tright->SetNDC();
181 tright->SetTextColor(2);
182 tright->Draw();
183
184 TLine *lin=new TLine(323,6.2,323,7.2);
185 lin->SetLineColor(4);
186 lin->Draw();
187 c8->Update();
188}
189
190void AnalyzeSDDInjectorsAllMod(Int_t nrun, Int_t n2, Int_t nDDL=0, Int_t firstEv=10, Int_t lastEv=15){
191 TGrid::Connect("alien:",0,0,"t");
192 Char_t filnam[200];
193 sprintf(filnam,"alien:///alice/data/2008/LHC08c_SDD/%09d/raw/08%09d%03d.10.root",nrun,nrun,n2);
194 printf("Open file %s\n",filnam);
195 AnalyzeSDDInjectorsAllMod(filnam,nDDL,firstEv,lastEv);
196}
197
198
199