]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDBASda.cxx
better organization of private and protected methods. Bug fix in AliITSPlaneEff:...
[u/mrichter/AliRoot.git] / ITS / ITSSDDBASda.cxx
CommitLineData
00252d8f 1/*
2- Contact: - prino@to.infn.it
f2c71038 3- Link: - http://www.to.infn.it/~prino/alice/RawData/run11161.date
4- Run Type: - PEDESTAL_RUN
00252d8f 5- DA Type: - LDC
6- Number of events needed: 100
7- Input Files: -
979b5a5f 8- Output Files: - SDDbase_step1_ddl*c*_sid*.data SDDbase_step2_ddl*c*_sid*.data
00252d8f 9- Trigger types used:
10*/
11
12
e44f571c 13//////////////////////////////////////////////////////////////////////////////
14// Detector Algorithm for analysis of SDD baseline runs. //
15// //
16// Produces ASCII and ROOT output files with: //
17// 1. anode quality bit //
18// 1. Baseline values //
19// 2. Raw noise //
20// 3. Common mode coefficients //
21// 4. Noise corrected for common mode //
22// Files are stored locally. //
23// The next DAQ-DA step on Test Pulse run writes files to FXS //
24// //
25// Author: F. Prino (prino@to.infn.it) //
26// //
27//////////////////////////////////////////////////////////////////////////////
28
29
30
31extern "C" {
32#include "daqDA.h"
33}
34
35#include "event.h"
36#include "monitor.h"
37
38#include <stdio.h>
39#include <stdlib.h>
40
41// ROOT includes
42#include <TFile.h>
43#include <TH1F.h>
44#include <TH2F.h>
00252d8f 45#include <TROOT.h>
46#include <TPluginManager.h>
e44f571c 47
48// AliRoot includes
49#include "AliRawReaderDate.h"
50#include "AliITSOnlineSDDBase.h"
51#include "AliITSOnlineSDDCMN.h"
52#include "AliITSRawStreamSDD.h"
53/* Main routine
54 Arguments: list of DATE raw data files
55*/
56int main(int argc, char **argv) {
57
58 int status = 0;
59
00252d8f 60 // line added to solve IO problems
61 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
62 "*",
63 "TStreamerInfo",
64 "RIO",
65 "TStreamerInfo()");
e44f571c 66
67 /* log start of process */
979b5a5f 68 printf("ITS SDD BASELINE+NOISE algorithm program started\n");
e44f571c 69
70
71 /* check that we got some arguments = list of files */
72 if (argc<2) {
73 printf("Wrong number of arguments\n");
74 return -1;
75 }
76
77
767934e8 78 Int_t eqOffset = 256;
79 Int_t DDLrange = 24;
f2c71038 80 Int_t maxNEvents=18; // maximum number of events to be analyzed
979b5a5f 81 const Int_t kTotDDL=24;
82 const Int_t kModPerDDL=12;
83 const Int_t kSides=2;
84 AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides];
85 AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides];
86 TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
e44f571c 87
88 Char_t hisnam[20];
979b5a5f 89 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
90 for(Int_t imod=0; imod<kModPerDDL;imod++){
91 for(Int_t isid=0;isid<kSides;isid++){
92 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
93 base[index]=new AliITSOnlineSDDBase(iddl,imod,isid);
94 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
95 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
96 }
e44f571c 97 }
98 }
99
100 /* report progress */
3083967f 101 daqDA_progressReport(8);
e44f571c 102 Int_t iev;
103 for(Int_t iStep=0;iStep<2;iStep++){
e44f571c 104 printf("Start Analysis Step %d\n",iStep);
105 iev=0;
106 if(iStep==1){
979b5a5f 107 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
108 for(Int_t imod=0; imod<kModPerDDL;imod++){
109 for(Int_t isid=0;isid<kSides;isid++){
110 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
111 corr[index]=new AliITSOnlineSDDCMN(iddl,imod,isid);
112 }
e44f571c 113 }
114 }
115 }
116
117 /* read the data files */
118 int n;
119 for (n=1;n<argc;n++) {
120
121 status=monitorSetDataSource( argv[n] );
122 if (status!=0) {
123 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
124 return -1;
125 }
126
127 /* report progress */
128 /* in this example, indexed on the number of files */
129 // Progress report inside the event loop as well?
3083967f 130 daqDA_progressReport(10+40*iStep*n/argc);
e44f571c 131
132 /* read the file */
133 for(;;) {
134 struct eventHeaderStruct *event;
135 eventTypeType eventT;
136
137 /* get next event */
138 status=monitorGetEventDynamic((void **)&event);
139 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
140 if (status!=0) {
141 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
142 return -1;
143 }
144
145 /* retry if got no event */
146 if (event==NULL) {
147 break;
148 }
149
150 iev++;
3083967f 151 if(iev>maxNEvents) break;
e44f571c 152
153 /* use event - here, just write event id to result file */
154 eventT=event->eventType;
155 switch (event->eventType){
156
157 /* START OF RUN */
158 case START_OF_RUN:
159 break;
160
161 /* END OF RUN */
162 case END_OF_RUN:
163 break;
164
3083967f 165 // case PHYSICS_EVENT: // comment this line for test raw data
166 // break; // comment this line for test raw data
167
e44f571c 168
169 case CALIBRATION_EVENT:
3083967f 170 break; // uncomment this line for test raw data
171 case PHYSICS_EVENT: // uncomment this line for test raw data
e44f571c 172 printf(" event number = %i \n",iev);
173 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
767934e8 174 rawReader->SelectEquipment(17,eqOffset+1,eqOffset+DDLrange);
e44f571c 175
e44f571c 176
177 Int_t evtyp=0;
178 while(rawReader->ReadHeader()){
179 const UInt_t *subev = rawReader->GetSubEventAttributes();
180 if(subev[0]==0 && subev[1]==0 && subev[2]==0) evtyp=1;
181 }
182
183 rawReader->Reset();
979b5a5f 184 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
185 for(Int_t imod=0; imod<kModPerDDL;imod++){
186 for(Int_t isid=0;isid<kSides;isid++){
187 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
188 histo[index]->Reset();
189 }
e44f571c 190 }
191 }
192 AliITSRawStreamSDD s(rawReader);
e44f571c 193
194 while(s.Next()){
979b5a5f 195 Int_t iDDL=rawReader->GetDDLID();
196 Int_t iCarlos=s.GetCarlosId();
197 if(iDDL>=0 && iDDL<kTotDDL && s.IsCompletedModule()==kFALSE){
198 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s.GetChannel();
e44f571c 199 histo[index]->Fill(s.GetCoord2(),s.GetCoord1(),s.GetSignal());
200 }
201 }
202 delete rawReader;
979b5a5f 203 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
204 for(Int_t imod=0; imod<kModPerDDL;imod++){
205 for(Int_t isid=0;isid<kSides;isid++){
206 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
207 if(iStep==0) base[index]->AddEvent(histo[index]);
208 if(iStep==1) corr[index]->AddEvent(histo[index]);
209 }
e44f571c 210 }
211 }
212
213 /* free resources */
214 free(event);
215 }
216 }
e44f571c 217 }
979b5a5f 218
219 if(iStep==0){
220 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
221 for(Int_t imod=0; imod<kModPerDDL;imod++){
222 for(Int_t isid=0;isid<kSides;isid++){
223 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
224 base[index]->ValidateAnodes();
225 base[index]->WriteToASCII();
226 }
e44f571c 227 }
979b5a5f 228 }
229 }
e44f571c 230 }
231
232 /* write report */
233 printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),iev);
234
235 /* report progress */
236 daqDA_progressReport(90);
237
238
239
3083967f 240 TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
979b5a5f 241 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
242 for(Int_t imod=0; imod<kModPerDDL;imod++){
243 for(Int_t isid=0;isid<kSides;isid++){
244 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
245 corr[index]->ValidateAnodes();
246 corr[index]->WriteToASCII();
247 corr[index]->WriteToROOT(fh);
248 }
e44f571c 249 }
250 }
251 fh->Close();
252
253
e44f571c 254 /* report progress */
255 daqDA_progressReport(100);
256
257
258
259 return status;
260}