]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDBASda.cxx
Read data format and sampling frequency from CDH in SDD DAs (F. Prino)
[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>
42e868a7 43#include <TSystem.h>
e44f571c 44#include <TH1F.h>
45#include <TH2F.h>
00252d8f 46#include <TROOT.h>
47#include <TPluginManager.h>
e44f571c 48
49// AliRoot includes
50#include "AliRawReaderDate.h"
51#include "AliITSOnlineSDDBase.h"
52#include "AliITSOnlineSDDCMN.h"
53#include "AliITSRawStreamSDD.h"
a8cd2ab2 54#include "AliITSRawStreamSDDCompressed.h"
e44f571c 55/* Main routine
56 Arguments: list of DATE raw data files
57*/
58int main(int argc, char **argv) {
59
60 int status = 0;
61
00252d8f 62 // line added to solve IO problems
63 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
64 "*",
65 "TStreamerInfo",
66 "RIO",
67 "TStreamerInfo()");
e44f571c 68
69 /* log start of process */
979b5a5f 70 printf("ITS SDD BASELINE+NOISE algorithm program started\n");
e44f571c 71
72
73 /* check that we got some arguments = list of files */
74 if (argc<2) {
75 printf("Wrong number of arguments\n");
76 return -1;
77 }
78
79
5065707f 80 Int_t maxNEvents=10; // 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;
040d0012 84 UInt_t amSamplFreq=40;
85
5065707f 86 gSystem->Exec("rm -f SDDbase_*.data");
456ec396 87 gSystem->Exec("rm -f SDDbase_step2_LDC.tar");
42e868a7 88
979b5a5f 89 AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides];
90 AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides];
91 TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
456ec396 92 Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
040d0012 93 Bool_t writtenoutput=kFALSE;
e44f571c 94
95 Char_t hisnam[20];
979b5a5f 96 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
97 for(Int_t imod=0; imod<kModPerDDL;imod++){
98 for(Int_t isid=0;isid<kSides;isid++){
99 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
100 base[index]=new AliITSOnlineSDDBase(iddl,imod,isid);
101 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
102 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
456ec396 103 isFilled[index]=0;
979b5a5f 104 }
e44f571c 105 }
106 }
107
108 /* report progress */
3083967f 109 daqDA_progressReport(8);
e44f571c 110 Int_t iev;
5065707f 111 Int_t ievPed;
e44f571c 112 for(Int_t iStep=0;iStep<2;iStep++){
e44f571c 113 printf("Start Analysis Step %d\n",iStep);
114 iev=0;
5065707f 115 ievPed=0;
e44f571c 116 if(iStep==1){
979b5a5f 117 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
118 for(Int_t imod=0; imod<kModPerDDL;imod++){
119 for(Int_t isid=0;isid<kSides;isid++){
120 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
121 corr[index]=new AliITSOnlineSDDCMN(iddl,imod,isid);
456ec396 122 isFilled[index]=0;
979b5a5f 123 }
e44f571c 124 }
125 }
126 }
127
128 /* read the data files */
129 int n;
130 for (n=1;n<argc;n++) {
131
132 status=monitorSetDataSource( argv[n] );
133 if (status!=0) {
134 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
135 return -1;
136 }
137
138 /* report progress */
139 /* in this example, indexed on the number of files */
140 // Progress report inside the event loop as well?
3083967f 141 daqDA_progressReport(10+40*iStep*n/argc);
e44f571c 142
143 /* read the file */
144 for(;;) {
145 struct eventHeaderStruct *event;
146 eventTypeType eventT;
147
148 /* get next event */
149 status=monitorGetEventDynamic((void **)&event);
150 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
151 if (status!=0) {
152 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
153 return -1;
154 }
155
156 /* retry if got no event */
157 if (event==NULL) {
158 break;
159 }
160
161 iev++;
e44f571c 162
163 /* use event - here, just write event id to result file */
164 eventT=event->eventType;
165 switch (event->eventType){
166
167 /* START OF RUN */
168 case START_OF_RUN:
169 break;
170
171 /* END OF RUN */
172 case END_OF_RUN:
173 break;
174
3083967f 175 // case PHYSICS_EVENT: // comment this line for test raw data
176 // break; // comment this line for test raw data
177
e44f571c 178
179 case CALIBRATION_EVENT:
3083967f 180 break; // uncomment this line for test raw data
181 case PHYSICS_EVENT: // uncomment this line for test raw data
e44f571c 182 printf(" event number = %i \n",iev);
5065707f 183 ievPed++;
e44f571c 184 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
040d0012 185 rawReader->Reset();
186 UChar_t cdhAttr=AliITSRawStreamSDD::ReadBlockAttributes(rawReader);
187 amSamplFreq=AliITSRawStreamSDD::ReadAMSamplFreqFromCDH(cdhAttr);
188 AliITSRawStream* s=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader,cdhAttr);
189 if(!writtenoutput){
190 printf("Use %s raw stream, sampling frequency %d MHz\n",s->ClassName(),amSamplFreq);
191 writtenoutput=kTRUE;
e44f571c 192 }
193
040d0012 194
979b5a5f 195 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
196 for(Int_t imod=0; imod<kModPerDDL;imod++){
197 for(Int_t isid=0;isid<kSides;isid++){
198 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
199 histo[index]->Reset();
200 }
e44f571c 201 }
202 }
e44f571c 203
a8cd2ab2 204 while(s->Next()){
979b5a5f 205 Int_t iDDL=rawReader->GetDDLID();
a8cd2ab2 206 Int_t iCarlos=s->GetCarlosId();
207 if(s->IsCompletedModule()) continue;
208 if(s->IsCompletedDDL()) continue;
de075dae 209 if(iDDL>=0 && iDDL<kTotDDL){
a8cd2ab2 210 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s->GetChannel();
211 histo[index]->Fill(s->GetCoord2(),s->GetCoord1(),s->GetSignal());
456ec396 212 isFilled[index]=1;
e44f571c 213 }
214 }
a8cd2ab2 215 delete s;
e44f571c 216 delete rawReader;
979b5a5f 217 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
218 for(Int_t imod=0; imod<kModPerDDL;imod++){
219 for(Int_t isid=0;isid<kSides;isid++){
220 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
040d0012 221 if(iStep==0){
222 if(amSamplFreq==20) base[index]->SetLastGoodTB(126);
223 else base[index]->SetLastGoodTB(254);
224 base[index]->AddEvent(histo[index]);
225 }
226 if(iStep==1){
227 if(amSamplFreq==20) corr[index]->SetLastGoodTB(126);
228 else corr[index]->SetLastGoodTB(254);
229 corr[index]->AddEvent(histo[index]);
230 }
979b5a5f 231 }
e44f571c 232 }
233 }
234
235 /* free resources */
236 free(event);
237 }
5065707f 238 if(ievPed>=maxNEvents) break;
e44f571c 239 }
e44f571c 240 }
979b5a5f 241
242 if(iStep==0){
243 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
244 for(Int_t imod=0; imod<kModPerDDL;imod++){
245 for(Int_t isid=0;isid<kSides;isid++){
246 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
247 base[index]->ValidateAnodes();
248 base[index]->WriteToASCII();
249 }
e44f571c 250 }
5065707f 251 }
979b5a5f 252 }
e44f571c 253 }
254
255 /* write report */
5065707f 256 printf("Run #%s, received %d calibration events\n",getenv("DATE_RUN_NUMBER"),ievPed);
e44f571c 257
258 /* report progress */
259 daqDA_progressReport(90);
260
261
262
456ec396 263 Char_t filnam[100],command[150];
3083967f 264 TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
979b5a5f 265 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
266 for(Int_t imod=0; imod<kModPerDDL;imod++){
267 for(Int_t isid=0;isid<kSides;isid++){
268 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
269 corr[index]->ValidateAnodes();
270 corr[index]->WriteToASCII();
456ec396 271 if(isFilled[index]){
272 corr[index]->WriteToROOT(fh);
273 sprintf(filnam,"SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
274 sprintf(command,"tar -rf SDDbase_step2_LDC.tar %s",filnam);
275 gSystem->Exec(command);
276 }
979b5a5f 277 }
e44f571c 278 }
279 }
280 fh->Close();
281
282
e44f571c 283 /* report progress */
284 daqDA_progressReport(100);
285
286
287
288 return status;
289}