Fix
[u/mrichter/AliRoot.git] / ZDC / ZDCEMDda.cxx
CommitLineData
f16f149b 1/*
2
442e1b18 3This program reads the DAQ data files passed as argument using the monitoring library.
f16f149b 4
442e1b18 5It computes the average event size and populates local "./result.txt" file with the
6result.
f16f149b 7
442e1b18 8The program reports about its processing progress.
f16f149b 9
10Messages on stdout are exported to DAQ log system.
11
f9641c3b 12DA for ZDC standalone CALIBRATION_EMD runs
442e1b18 13
ea628de7 14Contact: Chiara.Oppedisano@to.infn.it
442e1b18 15Link:
ea628de7 16Run Type: STANDALONE_EMD_RUN
442e1b18 17DA Type: LDC
18Number of events needed: at least ~5*10^3
ea628de7 19Input Files: ZDCPedestal.dat
442e1b18 20Output Files: ZDCEMDCalib.dat, ZDCChMapping.dat
ea628de7 21Trigger Types Used: Standalone Trigger
f16f149b 22
23*/
78beff0d 24#define PEDDATA_FILE "ZDCPedestal.dat"
218f916a 25#define MAPDATA_FILE "ZDCChMapping.dat"
92fdeb0e 26#define ENCALIBDATA_FILE "ZDCEnergyCalib.dat"
f9641c3b 27#define TOWCALIBDATA_FILE "ZDCTowerCalib.dat"
f16f149b 28
29#include <stdio.h>
30#include <Riostream.h>
442e1b18 31#include <Riostream.h>
f16f149b 32
33// DATE
34#include <daqDA.h>
35#include <event.h>
36#include <monitor.h>
37
38//ROOT
f9641c3b 39#include <TROOT.h>
40#include <TPluginManager.h>
f16f149b 41#include <TH1F.h>
42#include <TH2F.h>
43#include <TProfile.h>
44#include <TF1.h>
45#include <TFile.h>
442e1b18 46#include <TFitter.h>
f9641c3b 47#include "TMinuitMinimizer.h"
f16f149b 48
49//AliRoot
50#include <AliRawReaderDate.h>
442e1b18 51#include <AliRawEventHeaderBase.h>
f16f149b 52#include <AliZDCRawStream.h>
53
54
55/* Main routine
56 Arguments:
57 1- monitoring data source
58*/
59int main(int argc, char **argv) {
60
f9641c3b 61
62 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
63 "*",
64 "TStreamerInfo",
65 "RIO",
66 "TStreamerInfo()");
67
68 TMinuitMinimizer m;
69 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
70 "Minuit", "TMinuitMinimizer(const char *)");
71 TVirtualFitter::SetDefaultFitter("Minuit");
442e1b18 72
73 int status = 0;
f9641c3b 74 // No. of ZDC cabled ch.
75 int const kNChannels = 24;
27afc0c8 76 int const kNScChannels = 32;
442e1b18 77
78 /* log start of process */
79 printf("ZDC EMD program started\n");
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
f16f149b 87 //
88 // --- Preparing histos for EM dissociation spectra
89 //
90 TH1F* histoEMDRaw[4];
91 TH1F* histoEMDCorr[4];
92
93 char namhistr[50], namhistc[50];
94 for(Int_t i=0; i<4; i++) {
95 if(i==0){
96 sprintf(namhistr,"ZN%d-EMDRaw",i+1);
97 sprintf(namhistc,"ZN%d-EMDCorr",i+1);
98 }
99 else if(i==1){
100 sprintf(namhistr,"ZP%d-EMDRaw",i);
101 sprintf(namhistc,"ZP%d-EMDCorr",i);
102 }
103 else if(i==2){
104 sprintf(namhistr,"ZN%d-EMDRaw",i);
105 sprintf(namhistc,"ZN%d-EMDCorr",i);
106 }
107 else if(i==3){
108 sprintf(namhistr,"ZP%d-EMDRaw",i-1);
109 sprintf(namhistc,"ZP%d-EMDCorr",i-1);
110 }
111 histoEMDRaw[i] = new TH1F(namhistr,namhistr,100,0.,4000.);
112 histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
113 }
114
f16f149b 115 /* open result file */
116 FILE *fp=NULL;
117 fp=fopen("./result.txt","a");
118 if (fp==NULL) {
119 printf("Failed to open file\n");
120 return -1;
121 }
122
442e1b18 123 FILE *mapFile4Shuttle;
218f916a 124
f9641c3b 125 // *** To analyze EMD events you MUST have a pedestal data file!!!
126 // *** -> check if a pedestal run has been analyzed
78beff0d 127 int read = 0;
218f916a 128 read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
78beff0d 129 if(read){
130 printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
f16f149b 131 return -1;
132 }
78beff0d 133 else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
134
135 FILE *filePed = fopen(PEDDATA_FILE,"r");
136 if (filePed==NULL) {
137 printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
f16f149b 138 return -1;
139 }
140
78beff0d 141 // 144 = 48 in-time + 48 out-of-time + 48 correlations
f9641c3b 142 Float_t readValues[2][3*2*kNChannels];
143 Float_t MeanPed[2*kNChannels];
144 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
78beff0d 145 // ***************************************************
146 // Unless we have a narrow correlation to fit we
147 // don't fit and store in-time vs. out-of-time
148 // histograms -> mean pedstal subtracted!!!!!!
149 // ***************************************************
78beff0d 150 //
f9641c3b 151 for(int jj=0; jj<6*kNChannels; jj++){
78beff0d 152 for(int ii=0; ii<2; ii++){
153 fscanf(filePed,"%f",&readValues[ii][jj]);
154 }
f9641c3b 155 if(jj<kNChannels && jj<2*kNChannels){
78beff0d 156 MeanPed[jj] = readValues[0][jj];
f9641c3b 157 //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
78beff0d 158 }
f9641c3b 159 else if(jj>2*kNChannels && jj>4*kNChannels){
160 CorrCoeff0[jj-4*kNChannels] = readValues[0][jj];
161 CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
78beff0d 162 }
78beff0d 163 }
f16f149b 164
442e1b18 165 /* report progress */
166 daqDA_progressReport(10);
167
168
f16f149b 169 /* init some counters */
170 int nevents_physics=0;
171 int nevents_total=0;
172
442e1b18 173 /* read the data files */
174 int n;
175 for(n=1;n<argc;n++){
176
177 status=monitorSetDataSource( argv[n] );
f16f149b 178 if (status!=0) {
442e1b18 179 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
180 return -1;
f16f149b 181 }
182
442e1b18 183 /* report progress */
184 /* in this example, indexed on the number of files */
185 daqDA_progressReport(10+80*n/argc);
186
187 /* read the file */
188 for(;;){
189 struct eventHeaderStruct *event;
190 eventTypeType eventT;
f16f149b 191
442e1b18 192 /* get next event */
193 status=monitorGetEventDynamic((void **)&event);
194 if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
195 if(status!=0) {
196 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
197 return -1;
198 }
f16f149b 199
442e1b18 200 /* retry if got no event */
201 if(event==NULL) {
202 break;
203 }
204
205 // Initalize raw-data reading and decoding
206 AliRawReader *reader = new AliRawReaderDate((void*)event);
207 reader->Select("ZDC");
208 // --- Reading event header
209 //UInt_t evtype = reader->GetType();
210 //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
211 //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
212 //
213 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
214
215
216 /* use event - here, just write event id to result file */
217 eventT=event->eventType;
218
f9641c3b 219 Int_t ich=0;
27afc0c8 220 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
221 Int_t det[2*kNChannels], sec[2*kNChannels];
222 for(Int_t y=0; y<2*kNChannels; y++){
223 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
224 }
225
226 Int_t iScCh=0;
227 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
228 Int_t scDet[kNScChannels], scSec[kNScChannels];
229 for(Int_t y=0; y<kNScChannels; y++){
230 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
231 }
232 //
233 Int_t modNum=-1, modType=-1;
f9641c3b 234
442e1b18 235 if(eventT==START_OF_DATA){
236
27afc0c8 237 rawStreamZDC->SetSODReading(kTRUE);
238
f9641c3b 239 // --------------------------------------------------------
240 // --- Writing ascii data file for the Shuttle preprocessor
241 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
442e1b18 242 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
243 else{
27afc0c8 244 while((rawStreamZDC->Next())){
245 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
246 modNum = rawStreamZDC->GetADCModule();
247 modType = rawStreamZDC->GetModType();
248 }
249 if(rawStreamZDC->IsChMapping()){
250 if(modType==1){ // ADC mapping ----------------------
251 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
252 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
253 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
254 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
255 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
256 //
257 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
258 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
259 //
260 //printf(" Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
261 // ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
262 //
263 ich++;
264 }
265 else if(modType==2){ //VME scaler mapping --------------------
266 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
267 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
268 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
269 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
270 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
271 //
272 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
273 iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
274 //
275 //printf(" Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
276 // iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
277 //
278 iScCh++;
279 }
442e1b18 280 }
281 }
282 }
442e1b18 283 fclose(mapFile4Shuttle);
27afc0c8 284 }// SOD event
f16f149b 285
286 if(eventT==PHYSICS_EVENT){
78beff0d 287 // --- Reading data header
288 reader->ReadHeader();
f16f149b 289 const AliRawDataHeader* header = reader->GetDataHeader();
290 if(header){
ea628de7 291 UChar_t message = header->GetAttributes();
292 if(message & 0x70){ // DEDICATED EMD RUN
442e1b18 293 //printf("\t STANDALONE_EMD_RUN raw data found\n");
ea628de7 294 continue;
295 }
296 else{
297 printf("\t NO STANDALONE_EMD_RUN raw data found\n");
298 return -1;
299 }
f16f149b 300 }
442e1b18 301 else{
302 printf("\t ATTENTION! No Raw Data Header found!!!\n");
303 return -1;
304 }
7f4bde92 305
306 rawStreamZDC->SetSODReading(kTRUE);
442e1b18 307
f16f149b 308 if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
309 //
442e1b18 310 // ----- Setting ch. mapping -----
f9641c3b 311 for(Int_t jk=0; jk<2*kNChannels; jk++){
442e1b18 312 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
313 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
314 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
315 rawStreamZDC->SetMapDet(jk, det[jk]);
316 rawStreamZDC->SetMapTow(jk, sec[jk]);
317 }
318 //
f16f149b 319 Float_t ZDCRawADC[4], ZDCCorrADC[4], ZDCCorrADCSum[4];
320 for(Int_t g=0; g<4; g++){
321 ZDCCorrADCSum[g] = 0.;
322 ZDCRawADC[g] = 0.;
323 }
324 //
325 while(rawStreamZDC->Next()){
f9641c3b 326 Int_t det = rawStreamZDC->GetSector(0);
327 Int_t quad = rawStreamZDC->GetSector(1);
328
329 if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
330 && !(rawStreamZDC->IsOverflow()) && det!=-1
331 && (rawStreamZDC->GetADCGain() == 1)){ // Selecting LOW RES ch.s
332
333 //printf(" IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
334 // rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
335
336 // Taking LOW RES channels -> channel+kNChannels !!!!
ea628de7 337 Int_t DetIndex=999, PedIndex=999;
f9641c3b 338 if(det != 3 && quad != 5){ // Not ZEM nor PMRef
339 if(det == 1){
340 DetIndex = det-1;
341 PedIndex = quad+kNChannels;
342 }
343 else if(det==2){
344 DetIndex = det-1;
345 PedIndex = quad+5+kNChannels;
346 }
347 else if(det == 4){
348 DetIndex = det-2;
349 PedIndex = quad+12+kNChannels;
350 }
351 else if(det == 5){
352 DetIndex = det-2;
353 PedIndex = quad+17+kNChannels;
354 }
355 //EMD -> LR ADCs
356 if(rawStreamZDC->GetADCGain() == 1 && (DetIndex!=999 || PedIndex!=999)){
357 //
358 ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
359 //
360 // Mean pedestal subtraction
361 Float_t Pedestal = MeanPed[PedIndex];
362 // Pedestal subtraction from correlation with out-of-time signals
363 //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
364 //
365 ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
366 ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
367 //
368 /*printf("\t det %d quad %d res %d pedInd %d detInd %d:"
369 "ADCCorr = %d, ZDCCorrADCSum = %d\n",
370 det,quad,rawStreamZDC->GetADCGain(),PedIndex,DetIndex,
371 (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);
372 */
373 }
374 if(DetIndex==999 || PedIndex==999)
375 printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
376
377 }
f16f149b 378 }//IsADCDataWord()
379 //
380 }
381 //
382 nevents_physics++;
383 //
384 for(Int_t j=0; j<4; j++){
385 histoEMDRaw[j]->Fill(ZDCRawADC[j]);
386 histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
387 }
388 }
389
390 nevents_total++;
391
f16f149b 392 /* free resources */
393 free(event);
394
395 /* exit when last event received, no need to wait for TERM signal */
396 if (eventT==END_OF_RUN) {
397 printf("EOR event detected\n");
398 break;
399 }
442e1b18 400
401 }
f16f149b 402 }
442e1b18 403
f16f149b 404 /* Analysis of the histograms */
405 //
f9641c3b 406 FILE *fileShuttle1 = fopen(ENCALIBDATA_FILE,"w");
f16f149b 407 //
408 Int_t BinMax[4];
409 Float_t YMax[4];
410 Int_t NBinsx[4];
411 Float_t MeanFitVal[4];
412 TF1 *fitfun[4];
413 for(Int_t k=0; k<4; k++){
f9641c3b 414 if(histoEMDCorr[k]->GetEntries() == 0){
415 printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
416 return -1;
417 }
f16f149b 418 BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
419 YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
420 NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
ea628de7 421// printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
f16f149b 422 histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
423 fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
424 MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
c38bc7ef 425 //printf("\n\t Mean Value from gaussian fit = %f\n", MeanFitVal[k]);
f16f149b 426 }
427 //
218f916a 428 Float_t CalibCoeff[6];
429 Float_t icoeff[5];
430 //
f9641c3b 431 for(Int_t j=0; j<6; j++){
218f916a 432 if(j<4){
433 CalibCoeff[j] = MeanFitVal[j];
f9641c3b 434 fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
218f916a 435 }
f9641c3b 436 // ZEM energy calib. coeff. = 1
218f916a 437 else if(j==4 || j==5){
438 CalibCoeff[j] = 1.;
f9641c3b 439 fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
218f916a 440 }
f9641c3b 441 }
442 fclose(fileShuttle1);
443 //
444 FILE *fileShuttle2 = fopen(TOWCALIBDATA_FILE,"w");
445 for(Int_t j=0; j<4; j++){
446 // Note -> For the moment the inter-calibration coeff. are set to 1
447 for(Int_t k=0; k<5; k++){
448 icoeff[k] = 1.;
449 fprintf(fileShuttle2,"\t%f",icoeff[k]);
450 if(k==4) fprintf(fileShuttle2,"\n");
218f916a 451 }
f16f149b 452 }
f9641c3b 453 fclose(fileShuttle2);
f16f149b 454
442e1b18 455 for(Int_t ij=0; ij<4; ij++){
456 delete histoEMDRaw[ij];
457 delete histoEMDCorr[ij];
458 }
459
460 //delete minuitFit;
461 TVirtualFitter::SetFitter(0);
f16f149b 462
463 /* write report */
464 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
465
466 /* close result file */
467 fclose(fp);
442e1b18 468
469 /* report progress */
470 daqDA_progressReport(90);
471
472 /* store the result file on FES */
218f916a 473 status = daqDA_FES_storeFile(MAPDATA_FILE, MAPDATA_FILE);
442e1b18 474 if(status){
475 printf("Failed to export file : %d\n",status);
476 return -1;
477 }
478 //
f9641c3b 479 status = daqDA_FES_storeFile(ENCALIBDATA_FILE, ENCALIBDATA_FILE);
480 if(status){
481 printf("Failed to export file : %d\n",status);
482 return -1;
483 }
484 //
485 status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, TOWCALIBDATA_FILE);
442e1b18 486 if(status){
487 printf("Failed to export file : %d\n",status);
488 return -1;
489 }
f16f149b 490
442e1b18 491 /* report progress */
492 daqDA_progressReport(100);
f16f149b 493
494 return status;
495}