]>
Commit | Line | Data |
---|---|---|
1 | /* | |
2 | ||
3 | This program reads the DAQ data files passed as argument using the monitoring library. | |
4 | ||
5 | It computes the average event size and populates local "./result.txt" file with the | |
6 | result. | |
7 | ||
8 | The program reports about its processing progress. | |
9 | ||
10 | Messages on stdout are exported to DAQ log system. | |
11 | ||
12 | DA for ZDC standalone CALIBRATION_MB runs | |
13 | ||
14 | Contact: Chiara.Oppedisano@to.infn.it | |
15 | Link: | |
16 | Run Type: CALIBRATION_MB_RUN | |
17 | DA Type: LDC | |
18 | Number of events needed: at least 10^4 | |
19 | Input Files: ZDCPedestal.dat | |
20 | Output Files: ZDCMBCalib.root, ZDCChMapping.dat | |
21 | Trigger Types Used: Standalone Trigger | |
22 | ||
23 | */ | |
24 | #define PEDDATA_FILE "ZDCPedestal.dat" | |
25 | #define MAPDATA_FILE "ZDCChMapping.dat" | |
26 | #define MBCALIBDATA_FILE "ZDCMBCalib.root" | |
27 | ||
28 | #include <stdio.h> | |
29 | #include <Riostream.h> | |
30 | #include <Riostream.h> | |
31 | ||
32 | // DATE | |
33 | #include <daqDA.h> | |
34 | #include <event.h> | |
35 | #include <monitor.h> | |
36 | ||
37 | //ROOT | |
38 | #include <TROOT.h> | |
39 | #include <TPluginManager.h> | |
40 | #include <TH1F.h> | |
41 | #include <TH2F.h> | |
42 | #include <TProfile.h> | |
43 | #include <TF1.h> | |
44 | #include <TFile.h> | |
45 | #include <TFitter.h> | |
46 | #include "TMinuitMinimizer.h" | |
47 | ||
48 | //AliRoot | |
49 | #include <AliRawReaderDate.h> | |
50 | #include <AliRawEventHeaderBase.h> | |
51 | #include <AliZDCRawStream.h> | |
52 | ||
53 | ||
54 | /* Main routine | |
55 | Arguments: | |
56 | 1- monitoring data source | |
57 | */ | |
58 | int main(int argc, char **argv) { | |
59 | ||
60 | ||
61 | gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo", | |
62 | "*", | |
63 | "TStreamerInfo", | |
64 | "RIO", | |
65 | "TStreamerInfo()"); | |
66 | ||
67 | TMinuitMinimizer m; | |
68 | gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer", | |
69 | "Minuit", "TMinuitMinimizer(const char *)"); | |
70 | TVirtualFitter::SetDefaultFitter("Minuit"); | |
71 | ||
72 | int status = 0; | |
73 | int const kNModules = 10; | |
74 | int const kNChannels = 24; | |
75 | int const kNScChannels = 32; | |
76 | Int_t kFirstADCGeo=0, kLastADCGeo=3; | |
77 | ||
78 | Int_t iMod=-1; | |
79 | Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules]; | |
80 | for(Int_t kl=0; kl<kNModules; kl++){ | |
81 | modGeo[kl]=modType[kl]=modNCh[kl]=0; | |
82 | } | |
83 | ||
84 | Int_t ich=0; | |
85 | Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels]; | |
86 | Int_t det[2*kNChannels], sec[2*kNChannels]; | |
87 | for(Int_t y=0; y<2*kNChannels; y++){ | |
88 | adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0; | |
89 | } | |
90 | ||
91 | Int_t iScCh=0; | |
92 | Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels]; | |
93 | Int_t scDet[kNScChannels], scSec[kNScChannels]; | |
94 | for(Int_t y=0; y<kNScChannels; y++){ | |
95 | scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0; | |
96 | } | |
97 | ||
98 | /* log start of process */ | |
99 | printf("ZDC EMD program started\n"); | |
100 | ||
101 | /* check that we got some arguments = list of files */ | |
102 | if (argc<2) { | |
103 | printf("Wrong number of arguments\n"); | |
104 | return -1; | |
105 | } | |
106 | ||
107 | // --- Preparing histos for EM dissociation spectra | |
108 | // | |
109 | TH2F * hZDCvsZEM = new TH2F("hZDCvsZEM","EZDCTot vs. EZEM",100,0.,5.,100,0.,800.); | |
110 | TH2F * hZDCCvsZEM = new TH2F("hZDCCvsZEM","EZDCC vs. EZEM",100,0.,5.,100,0.,400.); | |
111 | TH2F * hZDCAvsZEM = new TH2F("hZDCAvsZEM","EZDCA vs. EZEM",100,0.,5.,100,0.,400.); | |
112 | ||
113 | /* open result file */ | |
114 | FILE *fp=NULL; | |
115 | fp=fopen("./result.txt","a"); | |
116 | if (fp==NULL) { | |
117 | printf("Failed to open file\n"); | |
118 | return -1; | |
119 | } | |
120 | ||
121 | FILE *mapFile4Shuttle; | |
122 | ||
123 | // *** To analyze CALIBRATION_MB events a pedestal data file is needed!!! | |
124 | // *** -> check if a pedestal run has been analyzed | |
125 | int read = 0; | |
126 | read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE); | |
127 | if(read){ | |
128 | printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n"); | |
129 | return -1; | |
130 | } | |
131 | else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n"); | |
132 | ||
133 | FILE *filePed = fopen(PEDDATA_FILE,"r"); | |
134 | if (filePed==NULL) { | |
135 | printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n"); | |
136 | return -1; | |
137 | } | |
138 | ||
139 | // 144 = 48 in-time + 48 out-of-time + 48 correlations | |
140 | Float_t readValues[2][6*kNChannels]; | |
141 | Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels]; | |
142 | Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels]; | |
143 | // *************************************************** | |
144 | // Unless we have a narrow correlation to fit we | |
145 | // don't fit and store in-time vs. out-of-time | |
146 | // histograms -> mean pedstal subtracted!!!!!! | |
147 | // *************************************************** | |
148 | // | |
149 | for(int jj=0; jj<6*kNChannels; jj++){ | |
150 | for(int ii=0; ii<2; ii++){ | |
151 | fscanf(filePed,"%f",&readValues[ii][jj]); | |
152 | } | |
153 | if(jj<kNChannels){ | |
154 | MeanPedhg[jj] = readValues[0][jj]; | |
155 | //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]); | |
156 | } | |
157 | else if(jj>=kNChannels && jj<2*kNChannels){ | |
158 | MeanPedlg[jj-kNChannels] = readValues[0][jj]; | |
159 | //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]); | |
160 | } | |
161 | else if(jj>4*kNChannels){ | |
162 | CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; | |
163 | CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];; | |
164 | } | |
165 | } | |
166 | ||
167 | /* report progress */ | |
168 | daqDA_progressReport(10); | |
169 | ||
170 | ||
171 | /* init some counters */ | |
172 | int nevents_physics=0; | |
173 | int nevents_total=0; | |
174 | ||
175 | struct eventHeaderStruct *event; | |
176 | eventTypeType eventT; | |
177 | ||
178 | /* read the data files */ | |
179 | int n; | |
180 | for(n=1;n<argc;n++){ | |
181 | ||
182 | status=monitorSetDataSource( argv[n] ); | |
183 | if (status!=0) { | |
184 | printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status)); | |
185 | return -1; | |
186 | } | |
187 | ||
188 | /* report progress */ | |
189 | /* in this example, indexed on the number of files */ | |
190 | daqDA_progressReport(10+80*n/argc); | |
191 | ||
192 | /* read the file */ | |
193 | for(;;){ | |
194 | ||
195 | /* get next event */ | |
196 | status=monitorGetEventDynamic((void **)&event); | |
197 | if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */ | |
198 | if(status!=0) { | |
199 | printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status)); | |
200 | return -1; | |
201 | } | |
202 | ||
203 | /* retry if got no event */ | |
204 | if(event==NULL) { | |
205 | break; | |
206 | } | |
207 | ||
208 | // Initalize raw-data reading and decoding | |
209 | AliRawReader *reader = new AliRawReaderDate((void*)event); | |
210 | reader->Select("ZDC"); | |
211 | // --- Reading event header | |
212 | //UInt_t evtype = reader->GetType(); | |
213 | //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype); | |
214 | //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber()); | |
215 | // | |
216 | AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader); | |
217 | ||
218 | ||
219 | /* use event - here, just write event id to result file */ | |
220 | eventT=event->eventType; | |
221 | ||
222 | if(eventT==START_OF_DATA){ | |
223 | ||
224 | rawStreamZDC->SetSODReading(kTRUE); | |
225 | ||
226 | // -------------------------------------------------------- | |
227 | // --- Writing ascii data file for the Shuttle preprocessor | |
228 | mapFile4Shuttle = fopen(MAPDATA_FILE,"w"); | |
229 | if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n"); | |
230 | else{ | |
231 | while((rawStreamZDC->Next())){ | |
232 | if(rawStreamZDC->IsHeaderMapping()){ // mapping header | |
233 | iMod++; | |
234 | modGeo[iMod] = rawStreamZDC->GetADCModule(); | |
235 | modType[iMod] = rawStreamZDC->GetModType(); | |
236 | modNCh[iMod] = rawStreamZDC->GetADCNChannels(); | |
237 | } | |
238 | if(rawStreamZDC->IsChMapping()){ | |
239 | if(modType[iMod]==1){ // ADC mapping ---------------------- | |
240 | adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich); | |
241 | adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich); | |
242 | sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich); | |
243 | det[ich] = rawStreamZDC->GetDetectorFromMap(ich); | |
244 | sec[ich] = rawStreamZDC->GetTowerFromMap(ich); | |
245 | ich++; | |
246 | } | |
247 | else if(modType[iMod]==2){ //VME scaler mapping -------------------- | |
248 | scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh); | |
249 | scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh); | |
250 | scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh); | |
251 | scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh); | |
252 | scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh); | |
253 | iScCh++; | |
254 | } | |
255 | } | |
256 | } | |
257 | // Writing data on output FXS file | |
258 | for(Int_t is=0; is<2*kNChannels; is++){ | |
259 | fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n", | |
260 | is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]); | |
261 | //printf(" CalibMB DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n", | |
262 | // is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]); | |
263 | } | |
264 | for(Int_t is=0; is<kNScChannels; is++){ | |
265 | fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n", | |
266 | is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]); | |
267 | //printf(" CalibMB DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n", | |
268 | // is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]); | |
269 | } | |
270 | for(Int_t is=0; is<kNModules; is++){ | |
271 | fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n", | |
272 | modGeo[is],modType[is],modNCh[is]); | |
273 | //printf(" CalibMB DA -> Module mapping: geo %d type %d #ch %d\n", | |
274 | // modGeo[is],modType[is],modNCh[is]); | |
275 | } | |
276 | ||
277 | } | |
278 | fclose(mapFile4Shuttle); | |
279 | }// SOD event | |
280 | ||
281 | if(eventT==PHYSICS_EVENT){ | |
282 | // --- Reading data header | |
283 | reader->ReadHeader(); | |
284 | const AliRawDataHeader* header = reader->GetDataHeader(); | |
285 | if(header){ | |
286 | UChar_t message = header->GetAttributes(); | |
287 | if((message & 0xf0) == 0x60){ // DEDICATED CALIBRATION MB RUN | |
288 | //printf("\t STANDALONE_CALIBRAION_MB_RUN raw data found\n"); | |
289 | continue; | |
290 | } | |
291 | else{ | |
292 | printf("\t NO STANDALONE_MB_RUN raw data found\n"); | |
293 | return -1; | |
294 | } | |
295 | } | |
296 | else{ | |
297 | printf("\t ATTENTION! No Raw Data Header found!!!\n"); | |
298 | return -1; | |
299 | } | |
300 | ||
301 | rawStreamZDC->SetSODReading(kTRUE); | |
302 | ||
303 | if (!rawStreamZDC->Next()) printf(" \t No raw data found!! "); | |
304 | // | |
305 | // ----- Setting ch. mapping ----- | |
306 | for(Int_t jk=0; jk<2*kNChannels; jk++){ | |
307 | rawStreamZDC->SetMapADCMod(jk, adcMod[jk]); | |
308 | rawStreamZDC->SetMapADCCh(jk, adcCh[jk]); | |
309 | rawStreamZDC->SetMapADCSig(jk, sigCode[jk]); | |
310 | rawStreamZDC->SetMapDet(jk, det[jk]); | |
311 | rawStreamZDC->SetMapTow(jk, sec[jk]); | |
312 | } | |
313 | // | |
314 | Float_t rawZDC=0., rawZEM=0.; | |
315 | Float_t corrZDC=0., corrZEM=0.; | |
316 | Float_t corrZDCTot=0., corrZDCC=0., corrZDCA=0., corrZEMTot=0.; | |
317 | // | |
318 | while(rawStreamZDC->Next()){ | |
319 | Int_t detector = rawStreamZDC->GetSector(0); | |
320 | Int_t quad = rawStreamZDC->GetSector(1); | |
321 | Int_t index=-1; | |
322 | ||
323 | if( (rawStreamZDC->IsADCDataWord()) && !(rawStreamZDC->IsUnderflow()) | |
324 | && !(rawStreamZDC->IsOverflow()) && (detector!=-1) | |
325 | && ((rawStreamZDC->GetADCGain()) == 0) && // Selecting HIGH RES ch.s | |
326 | (rawStreamZDC->GetADCModule()>=kFirstADCGeo) && (rawStreamZDC->GetADCModule()<=kLastADCGeo)){ | |
327 | ||
328 | //printf(" IsADCWord %d, IsUnderflow %d, IsOverflow %d\n", | |
329 | // rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow()); | |
330 | ||
331 | if(quad!=5){ // Physics signals | |
332 | if(detector==1) index = quad; // *** ZNC | |
333 | else if(detector==2) index = quad+5; // *** ZPC | |
334 | else if(detector==3) index = quad+9; // *** ZEM | |
335 | else if(detector==4) index = quad+12;// *** ZNA | |
336 | else if(detector==5) index = quad+17;// *** ZPA | |
337 | } | |
338 | else continue; | |
339 | // | |
340 | if(index==-1) printf("ERROR in ZDCCALIBMBda.cxx -> det %d quad %d index %d ADC %d\n", | |
341 | detector, quad, index, rawStreamZDC->GetADCValue()); | |
342 | ||
343 | // Mean pedestal subtraction | |
344 | Float_t Pedestal = MeanPedhg[index]; | |
345 | // Pedestal subtraction from correlation with out-of-time signals | |
346 | //Float_t Pedestal = CorrCoeff0[index]+CorrCoeff1[index]*MeanPedOOT[index]; | |
347 | // | |
348 | if(index!=-1 && quad!=5){ | |
349 | // | |
350 | if(detector==3){ | |
351 | rawZEM = (Float_t) rawStreamZDC->GetADCValue(); | |
352 | corrZEM = (rawStreamZDC->GetADCValue()) - Pedestal; | |
353 | corrZEMTot += corrZEM; | |
354 | } | |
355 | else{ | |
356 | rawZDC = (Float_t) rawStreamZDC->GetADCValue(); | |
357 | corrZDC = rawZDC - Pedestal; | |
358 | if(detector==1 || detector==2) corrZDCC += corrZDC; | |
359 | else corrZDCA += corrZDC; | |
360 | corrZDCTot += corrZDC; | |
361 | } | |
362 | } | |
363 | }//IsADCDataWord() | |
364 | ||
365 | } // Next() | |
366 | ||
367 | hZDCvsZEM ->Fill(corrZDCTot, corrZEMTot); | |
368 | hZDCCvsZEM->Fill(corrZDCC, corrZEMTot); | |
369 | hZDCAvsZEM->Fill(corrZDCA, corrZEMTot); | |
370 | ||
371 | nevents_physics++; | |
372 | }//(if PHYSICS_EVENT) | |
373 | ||
374 | nevents_total++; | |
375 | ||
376 | } | |
377 | ||
378 | /* free resources */ | |
379 | free(event); | |
380 | } | |
381 | ||
382 | // | |
383 | TFile *histofile = new TFile(MBCALIBDATA_FILE,"RECREATE"); | |
384 | histofile->cd(); | |
385 | hZDCvsZEM ->Write(); | |
386 | hZDCCvsZEM->Write(); | |
387 | hZDCAvsZEM->Write(); | |
388 | histofile->Close(); | |
389 | ||
390 | if(hZDCvsZEM) delete hZDCvsZEM ; | |
391 | if(hZDCCvsZEM) delete hZDCCvsZEM; | |
392 | if(hZDCAvsZEM) delete hZDCAvsZEM; | |
393 | ||
394 | /* write report */ | |
395 | fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total); | |
396 | ||
397 | /* close result file */ | |
398 | fclose(fp); | |
399 | ||
400 | /* report progress */ | |
401 | daqDA_progressReport(90); | |
402 | ||
403 | /* store the result file on FES */ | |
404 | status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING"); | |
405 | if(status){ | |
406 | printf("Failed to export file : %d\n",status); | |
407 | return -1; | |
408 | } | |
409 | // | |
410 | status = daqDA_FES_storeFile(MBCALIBDATA_FILE, "MBCALIB"); | |
411 | if(status){ | |
412 | printf("Failed to export file : %d\n",status); | |
413 | return -1; | |
414 | } | |
415 | ||
416 | /* report progress */ | |
417 | daqDA_progressReport(100); | |
418 | ||
419 | return status; | |
420 | } |