3 This program reads the DAQ data files passed as argument using the monitoring library.
5 It computes the average event size and populates local "./result.txt" file with the
8 The program reports about its processing progress.
10 Messages on stdout are exported to DAQ log system.
12 DA for ZDC standalone CALIBRATION_EMD runs
14 Contact: Chiara.Oppedisano@to.infn.it
16 Run Type: STANDALONE_EMD_RUN
18 Number of events needed: at least ~5*10^3
19 Input Files: ZDCPedestal.dat
20 Output Files: ZDCEMDCalib.dat, ZDCChMapping.dat
21 Trigger Types Used: Standalone Trigger
24 #define PEDDATA_FILE "ZDCPedestal.dat"
25 #define MAPDATA_FILE "ZDCChMapping.dat"
26 #define ENCALIBDATA_FILE "ZDCEnergyCalib.dat"
27 #define TOWCALIBDATA_FILE "ZDCTowerCalib.dat"
30 #include <Riostream.h>
31 #include <Riostream.h>
40 #include <TPluginManager.h>
47 #include "TMinuitMinimizer.h"
50 #include <AliRawReaderDate.h>
51 #include <AliRawEventHeaderBase.h>
52 #include <AliZDCRawStream.h>
57 1- monitoring data source
59 int main(int argc, char **argv) {
62 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
69 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
70 "Minuit", "TMinuitMinimizer(const char *)");
71 TVirtualFitter::SetDefaultFitter("Minuit");
74 // No. of ZDC cabled ch.
75 int const kNModules = 10;
76 int const kNChannels = 24;
77 int const kNScChannels = 32;
78 Int_t kFirstADCGeo=0, kLastADCGeo=3;
81 Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
82 for(Int_t kl=0; kl<kNModules; kl++){
83 modGeo[kl]=modType[kl]=modNCh[kl]=0;
87 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
88 Int_t det[2*kNChannels], sec[2*kNChannels];
89 for(Int_t y=0; y<2*kNChannels; y++){
90 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
94 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
95 Int_t scDet[kNScChannels], scSec[kNScChannels];
96 for(Int_t y=0; y<kNScChannels; y++){
97 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
100 /* log start of process */
101 printf("ZDC EMD program started\n");
103 /* check that we got some arguments = list of files */
105 printf("Wrong number of arguments\n");
109 // --- Preparing histos for EM dissociation spectra
111 TH1F* histoEMDRaw[4];
112 TH1F* histoEMDCorr[4];
114 char namhistr[50], namhistc[50];
115 for(Int_t i=0; i<4; i++) {
117 sprintf(namhistr,"ZN%d-EMDRaw",i+1);
118 sprintf(namhistc,"ZN%d-EMDCorr",i+1);
121 sprintf(namhistr,"ZP%d-EMDRaw",i);
122 sprintf(namhistc,"ZP%d-EMDCorr",i);
125 sprintf(namhistr,"ZN%d-EMDRaw",i);
126 sprintf(namhistc,"ZN%d-EMDCorr",i);
129 sprintf(namhistr,"ZP%d-EMDRaw",i-1);
130 sprintf(namhistc,"ZP%d-EMDCorr",i-1);
132 histoEMDRaw[i] = new TH1F(namhistr,namhistr,100,0.,4000.);
133 histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
136 // --- Preparing histos for tower inter-calibration
138 TH1F* histZNCtow[4]; TH1F* histZPCtow[4];
139 TH1F* histZNAtow[4]; TH1F* histZPAtow[4];
141 char namhistznc[50], namhistzpc[50];
142 char namhistzna[50], namhistzpa[50];
143 for(Int_t i=0; i<4; i++) {
144 sprintf(namhistznc,"ZNC-tow%d",i+1);
145 sprintf(namhistzpc,"ZPC-tow%d",i+1);
146 sprintf(namhistzna,"ZNA-tow%d",i+1);
147 sprintf(namhistzpa,"ZPA-tow%d",i+1);
149 histZNCtow[i] = new TH1F(namhistznc,namhistznc,100,0.,4000.);
150 histZPCtow[i] = new TH1F(namhistzpc,namhistzpc,100,0.,4000.);
151 histZNAtow[i] = new TH1F(namhistzna,namhistzna,100,0.,4000.);
152 histZPAtow[i] = new TH1F(namhistzpa,namhistzpa,100,0.,4000.);
155 /* open result file */
157 fp=fopen("./result.txt","a");
159 printf("Failed to open file\n");
163 FILE *mapFile4Shuttle;
165 // *** To analyze EMD events you MUST have a pedestal data file!!!
166 // *** -> check if a pedestal run has been analyzed
168 read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
170 printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
173 else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
175 FILE *filePed = fopen(PEDDATA_FILE,"r");
177 printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
181 // 144 = 48 in-time + 48 out-of-time + 48 correlations
182 Float_t readValues[2][3*2*kNChannels];
183 Float_t MeanPed[2*kNChannels];
184 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
185 // ***************************************************
186 // Unless we have a narrow correlation to fit we
187 // don't fit and store in-time vs. out-of-time
188 // histograms -> mean pedstal subtracted!!!!!!
189 // ***************************************************
191 for(int jj=0; jj<6*kNChannels; jj++){
192 for(int ii=0; ii<2; ii++){
193 fscanf(filePed,"%f",&readValues[ii][jj]);
196 MeanPed[jj] = readValues[0][jj];
197 printf("\t MeanPed[%d] = %1.1f\n",jj, MeanPed[jj]);
199 else if(jj>2*kNChannels){
200 CorrCoeff0[jj-4*kNChannels] = readValues[0][jj];
201 CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
205 /* report progress */
206 daqDA_progressReport(10);
209 /* init some counters */
210 int nevents_physics=0;
213 struct eventHeaderStruct *event;
214 eventTypeType eventT;
216 /* read the data files */
220 status=monitorSetDataSource( argv[n] );
222 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
226 /* report progress */
227 /* in this example, indexed on the number of files */
228 daqDA_progressReport(10+80*n/argc);
234 status=monitorGetEventDynamic((void **)&event);
235 if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
237 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
241 /* retry if got no event */
246 // Initalize raw-data reading and decoding
247 AliRawReader *reader = new AliRawReaderDate((void*)event);
248 reader->Select("ZDC");
249 // --- Reading event header
250 //UInt_t evtype = reader->GetType();
251 //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
252 //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
254 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
257 /* use event - here, just write event id to result file */
258 eventT=event->eventType;
260 if(eventT==START_OF_DATA){
262 rawStreamZDC->SetSODReading(kTRUE);
264 // --------------------------------------------------------
265 // --- Writing ascii data file for the Shuttle preprocessor
266 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
267 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
269 while((rawStreamZDC->Next())){
270 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
272 modGeo[iMod] = rawStreamZDC->GetADCModule();
273 modType[iMod] = rawStreamZDC->GetModType();
274 modNCh[iMod] = rawStreamZDC->GetADCNChannels();
276 if(rawStreamZDC->IsChMapping()){
277 if(modType[iMod]==1){ // ADC mapping ----------------------
278 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
279 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
280 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
281 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
282 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
285 else if(modType[iMod]==2){ //VME scaler mapping --------------------
286 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
287 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
288 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
289 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
290 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
295 // Writing data on output FXS file
296 for(Int_t is=0; is<kNModules; is++){
297 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
298 modGeo[is],modType[is],modNCh[is]);
299 //printf(" EMD DA -> Module mapping: geo %d type %d #ch %d\n",
300 // modGeo[is],modType[is],modNCh[is]);
302 for(Int_t is=0; is<2*kNChannels; is++){
303 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
304 is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
305 //printf(" EMD DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
306 // is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
308 for(Int_t is=0; is<kNScChannels; is++){
309 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
310 is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
311 //printf(" EMD DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
312 // is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
316 fclose(mapFile4Shuttle);
319 if(eventT==PHYSICS_EVENT){
320 // --- Reading data header
321 reader->ReadHeader();
322 const AliRawDataHeader* header = reader->GetDataHeader();
324 UChar_t message = header->GetAttributes();
325 if((message & 0x70) == 0x70){ // DEDICATED EMD RUN
326 //printf("\t STANDALONE_EMD_RUN raw data found\n");
330 printf("\t NO STANDALONE_EMD_RUN raw data found\n");
335 printf("\t ATTENTION! No Raw Data Header found!!!\n");
339 rawStreamZDC->SetSODReading(kTRUE);
341 if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
343 // ----- Setting ch. mapping -----
344 for(Int_t jk=0; jk<2*kNChannels; jk++){
345 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
346 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
347 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
348 rawStreamZDC->SetMapDet(jk, det[jk]);
349 rawStreamZDC->SetMapTow(jk, sec[jk]);
352 Float_t ZDCRawADC[4], ZDCCorrADC[4], ZDCCorrADCSum[4];
353 for(Int_t g=0; g<4; g++){
354 ZDCCorrADCSum[g] = 0.;
358 while(rawStreamZDC->Next()){
359 Int_t det = rawStreamZDC->GetSector(0);
360 Int_t quad = rawStreamZDC->GetSector(1);
362 if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
363 && !(rawStreamZDC->IsOverflow()) && det!=-1 && det!=3
364 && (rawStreamZDC->GetADCGain() == 1 && // Selecting LOW RES ch.s
365 rawStreamZDC->GetADCModule()>=kFirstADCGeo && rawStreamZDC->GetADCModule()<=kLastADCGeo)){
367 // Taking LOW RES channels -> ch.+kNChannels !!!!
368 Int_t DetIndex=999, PedIndex=999;
373 PedIndex = quad+kNChannels;
377 PedIndex = quad+5+kNChannels;
381 PedIndex = quad+12+kNChannels;
385 PedIndex = quad+17+kNChannels;
387 // Mean pedestal subtraction
388 Float_t Pedestal = MeanPed[PedIndex];
389 // Pedestal subtraction from correlation with out-of-time signals
390 //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
392 if(DetIndex!=999 || PedIndex!=999){
394 ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
397 ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
398 ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
400 /*printf("\t det %d quad %d res %d pedInd %d "
401 "Pedestal %1.0f -> ADCCorr = %d ZDCCorrADCSum = %d\n",
402 det,quad,rawStreamZDC->GetADCGain(),PedIndex,Pedestal,
403 (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);*/
408 Float_t corrADCval = (rawStreamZDC->GetADCValue()) - Pedestal;
409 if(det==1) histZNCtow[quad-1]->Fill(corrADCval);
410 else if(det==2) histZPCtow[quad-1]->Fill(corrADCval);
411 else if(det==4) histZNAtow[quad-1]->Fill(corrADCval);
412 else if(det==5) histZPAtow[quad-1]->Fill(corrADCval);
414 //printf("\t det %d tow %d fill histo w. value %1.0f\n",
415 // det,quad,corrADCval);
418 if(DetIndex==999 || PedIndex==999)
419 printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
431 for(Int_t j=0; j<4; j++){
432 histoEMDRaw[j]->Fill(ZDCRawADC[j]);
433 histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
435 }//(if PHYSICS_EVENT)
437 /* exit when last event received, no need to wait for TERM signal */
438 else if(eventT==END_OF_RUN) {
439 printf(" -> EOR event detected\n");
451 /* Analysis of the histograms */
453 FILE *fileShuttle1 = fopen(ENCALIBDATA_FILE,"w");
455 Int_t BinMax[4]={0,0,0,0};
456 Float_t YMax[4]={0.,0.,0.,0.};
457 Int_t NBinsx[4]={0,0,0,0};
458 Float_t MeanFitVal[4]={0.,0.,0.,0.};
460 for(Int_t k=0; k<4; k++){
461 if(histoEMDCorr[k]->GetEntries() == 0){
462 printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
466 BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
468 printf("\n WARNING! Something wrong with det %d histo -> ending DA WITHOUT writing output\n\n", k);
472 YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
473 NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
474 //printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
475 histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
476 fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
477 MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
478 //printf("\n\t Mean Value from gaussian fit = %f\n", MeanFitVal[k]);
481 Float_t CalibCoeff[6];
483 for(Int_t j=0; j<6; j++){
485 CalibCoeff[j] = MeanFitVal[j];
486 fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
488 // ZEM energy calib. coeff. = 1
489 else if(j==4 || j==5){
491 fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
494 fclose(fileShuttle1);
496 FILE *fileShuttle2 = fopen(TOWCALIBDATA_FILE,"w");
497 //Float_t meanvalznc[4], meanvalzpc[4], meanvalzna[4], meanvalzpa[4];
498 for(Int_t j=0; j<4; j++){
499 /*if(histZNCtow[j]->GetEntries() == 0){
500 printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
503 meanvalznc[j] = histZNCtow[j]->GetMean();
504 meanvalzpc[j] = histZPCtow[j]->GetMean();
505 meanvalzna[j] = histZNAtow[j]->GetMean();
506 meanvalzpa[j] = histZPAtow[j]->GetMean();*/
508 // Note -> For the moment the inter-calibration coeff. are set to 1
509 for(Int_t k=0; k<5; k++){
511 fprintf(fileShuttle2,"\t%f",icoeff);
512 if(k==5) fprintf(fileShuttle2,"\n");
516 /*if(meanvalznc[1]!=0 && meanvalznc[2]!=0 && meanvalznc[3]!=0 &&
517 meanvalzpc[1]!=0 && meanvalzpc[2]!=0 && meanvalzpc[3]!=0 &&
518 meanvalzna[1]!=0 && meanvalzna[2]!=0 && meanvalzna[3]!=0 &&
519 meanvalzpa[1]!=0 && meanvalzpa[2]!=0 && meanvalzpa[3]!=0){
520 fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
521 1.0,meanvalznc[0]/meanvalznc[1],meanvalznc[0]/meanvalznc[2],meanvalznc[0]/meanvalznc[3]);
522 fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
523 1.0,meanvalzpc[0]/meanvalzpc[1],meanvalzpc[0]/meanvalzpc[2],meanvalzpc[0]/meanvalzpc[3]);
524 fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
525 1.0,meanvalzna[0]/meanvalzna[1],meanvalzpc[0]/meanvalzna[2],meanvalzpc[0]/meanvalzna[3]);
526 fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
527 1.0,meanvalzpa[0]/meanvalzpa[1],meanvalzpc[0]/meanvalzpa[2],meanvalzpc[0]/meanvalzpa[3]);
530 printf("\n Tower intercalib. coeff. CAN'T be calculated (some mean values are ZERO)!!!\n\n");
533 fclose(fileShuttle2);
535 for(Int_t ij=0; ij<4; ij++){
536 delete histoEMDRaw[ij];
537 delete histoEMDCorr[ij];
541 TVirtualFitter::SetFitter(0);
544 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
546 /* close result file */
549 /* report progress */
550 daqDA_progressReport(90);
552 /* store the result file on FES */
553 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
555 printf("Failed to export file : %d\n",status);
559 status = daqDA_FES_storeFile(ENCALIBDATA_FILE, "EMDENERGYCALIB");
561 printf("Failed to export file : %d\n",status);
565 status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, "EMDTOWERCALIB");
567 printf("Failed to export file : %d\n",status);
571 /* report progress */
572 daqDA_progressReport(100);