]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITSSDDBASda.cxx
Fixing bug #59311
[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
e5b2f7f2 6- Number of events needed: > 10
00252d8f 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>
e5b2f7f2 48#include <TObjArray.h>
49#include <TObjString.h>
50#include <TDatime.h>
e44f571c 51
52// AliRoot includes
53#include "AliRawReaderDate.h"
54#include "AliITSOnlineSDDBase.h"
55#include "AliITSOnlineSDDCMN.h"
56#include "AliITSRawStreamSDD.h"
a8cd2ab2 57#include "AliITSRawStreamSDDCompressed.h"
e5b2f7f2 58
59#ifdef ALI_AMORE
60#include <AmoreDA.h>
61#endif
62
e44f571c 63/* Main routine
64 Arguments: list of DATE raw data files
65*/
66int main(int argc, char **argv) {
67
68 int status = 0;
69
00252d8f 70 // line added to solve IO problems
71 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
72 "*",
73 "TStreamerInfo",
74 "RIO",
75 "TStreamerInfo()");
e44f571c 76
77 /* log start of process */
979b5a5f 78 printf("ITS SDD BASELINE+NOISE algorithm program started\n");
e44f571c 79
80
81 /* check that we got some arguments = list of files */
82 if (argc<2) {
83 printf("Wrong number of arguments\n");
84 return -1;
85 }
86
87
5065707f 88 Int_t maxNEvents=10; // maximum number of events to be analyzed
979b5a5f 89 const Int_t kTotDDL=24;
90 const Int_t kModPerDDL=12;
91 const Int_t kSides=2;
040d0012 92 UInt_t amSamplFreq=40;
93
5065707f 94 gSystem->Exec("rm -f SDDbase_*.data");
456ec396 95 gSystem->Exec("rm -f SDDbase_step2_LDC.tar");
42e868a7 96
979b5a5f 97 AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides];
98 AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides];
99 TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides];
456ec396 100 Bool_t isFilled[kTotDDL*kModPerDDL*kSides];
040d0012 101 Bool_t writtenoutput=kFALSE;
e44f571c 102
103 Char_t hisnam[20];
979b5a5f 104 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
105 for(Int_t imod=0; imod<kModPerDDL;imod++){
106 for(Int_t isid=0;isid<kSides;isid++){
107 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
108 base[index]=new AliITSOnlineSDDBase(iddl,imod,isid);
109 sprintf(hisnam,"h%02dc%02ds%d",iddl,imod,isid);
110 histo[index]=new TH2F(hisnam,"",256,-0.5,255.5,256,-0.5,255.5);
456ec396 111 isFilled[index]=0;
979b5a5f 112 }
e44f571c 113 }
114 }
115
116 /* report progress */
3083967f 117 daqDA_progressReport(8);
e44f571c 118 Int_t iev;
5065707f 119 Int_t ievPed;
e44f571c 120 for(Int_t iStep=0;iStep<2;iStep++){
e44f571c 121 printf("Start Analysis Step %d\n",iStep);
122 iev=0;
5065707f 123 ievPed=0;
e44f571c 124 if(iStep==1){
979b5a5f 125 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
126 for(Int_t imod=0; imod<kModPerDDL;imod++){
127 for(Int_t isid=0;isid<kSides;isid++){
128 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
129 corr[index]=new AliITSOnlineSDDCMN(iddl,imod,isid);
456ec396 130 isFilled[index]=0;
979b5a5f 131 }
e44f571c 132 }
133 }
134 }
135
136 /* read the data files */
137 int n;
138 for (n=1;n<argc;n++) {
139
140 status=monitorSetDataSource( argv[n] );
141 if (status!=0) {
142 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
143 return -1;
144 }
145
146 /* report progress */
147 /* in this example, indexed on the number of files */
148 // Progress report inside the event loop as well?
3083967f 149 daqDA_progressReport(10+40*iStep*n/argc);
e44f571c 150
151 /* read the file */
152 for(;;) {
153 struct eventHeaderStruct *event;
154 eventTypeType eventT;
155
156 /* get next event */
157 status=monitorGetEventDynamic((void **)&event);
158 if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
159 if (status!=0) {
160 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
161 return -1;
162 }
163
164 /* retry if got no event */
165 if (event==NULL) {
166 break;
167 }
168
169 iev++;
e44f571c 170
171 /* use event - here, just write event id to result file */
172 eventT=event->eventType;
173 switch (event->eventType){
174
175 /* START OF RUN */
176 case START_OF_RUN:
177 break;
178
179 /* END OF RUN */
180 case END_OF_RUN:
181 break;
182
3083967f 183 // case PHYSICS_EVENT: // comment this line for test raw data
184 // break; // comment this line for test raw data
185
e44f571c 186
187 case CALIBRATION_EVENT:
3083967f 188 break; // uncomment this line for test raw data
189 case PHYSICS_EVENT: // uncomment this line for test raw data
e44f571c 190 printf(" event number = %i \n",iev);
5065707f 191 ievPed++;
e44f571c 192 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
040d0012 193 rawReader->Reset();
194 UChar_t cdhAttr=AliITSRawStreamSDD::ReadBlockAttributes(rawReader);
195 amSamplFreq=AliITSRawStreamSDD::ReadAMSamplFreqFromCDH(cdhAttr);
196 AliITSRawStream* s=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader,cdhAttr);
197 if(!writtenoutput){
198 printf("Use %s raw stream, sampling frequency %d MHz\n",s->ClassName(),amSamplFreq);
199 writtenoutput=kTRUE;
e44f571c 200 }
201
040d0012 202
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 histo[index]->Reset();
208 }
e44f571c 209 }
210 }
e44f571c 211
a8cd2ab2 212 while(s->Next()){
979b5a5f 213 Int_t iDDL=rawReader->GetDDLID();
a8cd2ab2 214 Int_t iCarlos=s->GetCarlosId();
215 if(s->IsCompletedModule()) continue;
216 if(s->IsCompletedDDL()) continue;
de075dae 217 if(iDDL>=0 && iDDL<kTotDDL){
a8cd2ab2 218 Int_t index=kSides*(kModPerDDL*iDDL+iCarlos)+s->GetChannel();
219 histo[index]->Fill(s->GetCoord2(),s->GetCoord1(),s->GetSignal());
456ec396 220 isFilled[index]=1;
e44f571c 221 }
222 }
a8cd2ab2 223 delete s;
e44f571c 224 delete rawReader;
979b5a5f 225 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
226 for(Int_t imod=0; imod<kModPerDDL;imod++){
227 for(Int_t isid=0;isid<kSides;isid++){
228 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
040d0012 229 if(iStep==0){
230 if(amSamplFreq==20) base[index]->SetLastGoodTB(126);
231 else base[index]->SetLastGoodTB(254);
232 base[index]->AddEvent(histo[index]);
233 }
234 if(iStep==1){
235 if(amSamplFreq==20) corr[index]->SetLastGoodTB(126);
236 else corr[index]->SetLastGoodTB(254);
237 corr[index]->AddEvent(histo[index]);
238 }
979b5a5f 239 }
e44f571c 240 }
241 }
242
243 /* free resources */
244 free(event);
245 }
5065707f 246 if(ievPed>=maxNEvents) break;
e44f571c 247 }
e44f571c 248 }
979b5a5f 249
250 if(iStep==0){
251 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
252 for(Int_t imod=0; imod<kModPerDDL;imod++){
253 for(Int_t isid=0;isid<kSides;isid++){
254 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
255 base[index]->ValidateAnodes();
256 base[index]->WriteToASCII();
257 }
e44f571c 258 }
5065707f 259 }
979b5a5f 260 }
e44f571c 261 }
262
263 /* write report */
e5b2f7f2 264 TDatime time;
265 TObjString timeinfo(Form("%02d%02d%02d%02d%02d%02d",time.GetYear()-2000,time.GetMonth(),time.GetDay(),time.GetHour(),time.GetMinute(),time.GetSecond()));
266 printf("Run #%s, received %d calibration events, time %s\n",getenv("DATE_RUN_NUMBER"),ievPed,timeinfo.GetString().Data());
e44f571c 267
268 /* report progress */
269 daqDA_progressReport(90);
270
e5b2f7f2 271 TObjArray* basHistos=new TObjArray();
272 TObjArray* noiseHistos=new TObjArray();
273 TObjArray* cmnHistos=new TObjArray();
274 TObjArray* corrnHistos=new TObjArray();
275 TObjArray* statusHistos=new TObjArray();
e44f571c 276
277
456ec396 278 Char_t filnam[100],command[150];
979b5a5f 279 for(Int_t iddl=0; iddl<kTotDDL;iddl++){
280 for(Int_t imod=0; imod<kModPerDDL;imod++){
281 for(Int_t isid=0;isid<kSides;isid++){
282 Int_t index=kSides*(kModPerDDL*iddl+imod)+isid;
283 corr[index]->ValidateAnodes();
284 corr[index]->WriteToASCII();
456ec396 285 if(isFilled[index]){
e5b2f7f2 286 basHistos->AddLast(corr[index]->GetBaselineAnodeHisto());
287 noiseHistos->AddLast(corr[index]->GetRawNoiseAnodeHisto());
288 cmnHistos->AddLast(corr[index]->GetCMNCoefAnodeHisto());
289 corrnHistos->AddLast(corr[index]->GetCorrNoiseAnodeHisto());
290 statusHistos->AddLast(corr[index]->GetStatusAnodeHisto());
456ec396 291 sprintf(filnam,"SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
292 sprintf(command,"tar -rf SDDbase_step2_LDC.tar %s",filnam);
293 gSystem->Exec(command);
294 }
979b5a5f 295 }
e44f571c 296 }
297 }
e44f571c 298
e5b2f7f2 299#ifdef ALI_AMORE
300 amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
9ce3c336 301 Int_t statusamore =0;
7f601ef8 302 statusamore += amoreDA.Send("TimeInfoPedestal",&timeinfo);
9ce3c336 303 statusamore += amoreDA.Send("Baselines",basHistos);
304 statusamore += amoreDA.Send("RawNoise",noiseHistos);
305 statusamore += amoreDA.Send("CommonMode",cmnHistos);
306 statusamore += amoreDA.Send("CorrectedNoise",corrnHistos);
307 statusamore += amoreDA.Send("NoisyChannels",statusHistos);
308 if ( statusamore )
e5b2f7f2 309 printf("Warning: Failed to write Arrays in the AMORE database\n");
310 else
311 printf("amoreDA.Send() OK\n");
312#else
313 printf("Warning: SDDBAS DA not compiled with AMORE support\n");
314#endif
315
316
317 TFile *fh=new TFile("SDDbaseHistos.root","RECREATE");
318 basHistos->Write();
319 noiseHistos->Write();
320 cmnHistos->Write();
321 corrnHistos->Write();
322 statusHistos->Write();
323 fh->Close();
e44f571c 324
e44f571c 325 /* report progress */
326 daqDA_progressReport(100);
327
328
329
330 return status;
331}