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 kNChannels = 24;
76 int const kNScChannels = 32;
78 /* log start of process */
79 printf("ZDC EMD program started\n");
81 /* check that we got some arguments = list of files */
83 printf("Wrong number of arguments\n");
88 // --- Preparing histos for EM dissociation spectra
91 TH1F* histoEMDCorr[4];
93 char namhistr[50], namhistc[50];
94 for(Int_t i=0; i<4; i++) {
96 sprintf(namhistr,"ZN%d-EMDRaw",i+1);
97 sprintf(namhistc,"ZN%d-EMDCorr",i+1);
100 sprintf(namhistr,"ZP%d-EMDRaw",i);
101 sprintf(namhistc,"ZP%d-EMDCorr",i);
104 sprintf(namhistr,"ZN%d-EMDRaw",i);
105 sprintf(namhistc,"ZN%d-EMDCorr",i);
108 sprintf(namhistr,"ZP%d-EMDRaw",i-1);
109 sprintf(namhistc,"ZP%d-EMDCorr",i-1);
111 histoEMDRaw[i] = new TH1F(namhistr,namhistr,100,0.,4000.);
112 histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
115 /* open result file */
117 fp=fopen("./result.txt","a");
119 printf("Failed to open file\n");
123 FILE *mapFile4Shuttle;
125 // *** To analyze EMD events you MUST have a pedestal data file!!!
126 // *** -> check if a pedestal run has been analyzed
128 read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
130 printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
133 else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
135 FILE *filePed = fopen(PEDDATA_FILE,"r");
137 printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
141 // 144 = 48 in-time + 48 out-of-time + 48 correlations
142 Float_t readValues[2][3*2*kNChannels];
143 Float_t MeanPed[2*kNChannels];
144 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
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 // ***************************************************
151 for(int jj=0; jj<6*kNChannels; jj++){
152 for(int ii=0; ii<2; ii++){
153 fscanf(filePed,"%f",&readValues[ii][jj]);
155 if(jj<kNChannels && jj<2*kNChannels){
156 MeanPed[jj] = readValues[0][jj];
157 //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
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];;
165 /* report progress */
166 daqDA_progressReport(10);
169 /* init some counters */
170 int nevents_physics=0;
173 struct eventHeaderStruct *event;
174 eventTypeType eventT;
176 /* read the data files */
180 status=monitorSetDataSource( argv[n] );
182 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
186 /* report progress */
187 /* in this example, indexed on the number of files */
188 daqDA_progressReport(10+80*n/argc);
194 status=monitorGetEventDynamic((void **)&event);
195 if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
197 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
201 /* retry if got no event */
206 // Initalize raw-data reading and decoding
207 AliRawReader *reader = new AliRawReaderDate((void*)event);
208 reader->Select("ZDC");
209 // --- Reading event header
210 //UInt_t evtype = reader->GetType();
211 //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
212 //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
214 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
217 /* use event - here, just write event id to result file */
218 eventT=event->eventType;
221 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
222 Int_t det[2*kNChannels], sec[2*kNChannels];
223 for(Int_t y=0; y<2*kNChannels; y++){
224 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
228 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
229 Int_t scDet[kNScChannels], scSec[kNScChannels];
230 for(Int_t y=0; y<kNScChannels; y++){
231 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
234 Int_t modNum=-1, modType=-1;
236 if(eventT==START_OF_DATA){
238 rawStreamZDC->SetSODReading(kTRUE);
240 // --------------------------------------------------------
241 // --- Writing ascii data file for the Shuttle preprocessor
242 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
243 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
245 while((rawStreamZDC->Next())){
246 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
247 modNum = rawStreamZDC->GetADCModule();
248 modType = rawStreamZDC->GetModType();
250 if(rawStreamZDC->IsChMapping()){
251 if(modType==1){ // ADC mapping ----------------------
252 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
253 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
254 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
255 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
256 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
258 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
259 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
261 //printf(" Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
262 // ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
266 else if(modType==2){ //VME scaler mapping --------------------
267 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
268 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
269 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
270 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
271 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
273 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
274 iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
276 //printf(" Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
277 // iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
284 fclose(mapFile4Shuttle);
287 if(eventT==PHYSICS_EVENT){
288 // --- Reading data header
289 reader->ReadHeader();
290 const AliRawDataHeader* header = reader->GetDataHeader();
292 UChar_t message = header->GetAttributes();
293 if(message & 0x70){ // DEDICATED EMD RUN
294 //printf("\t STANDALONE_EMD_RUN raw data found\n");
298 printf("\t NO STANDALONE_EMD_RUN raw data found\n");
303 printf("\t ATTENTION! No Raw Data Header found!!!\n");
307 rawStreamZDC->SetSODReading(kTRUE);
309 if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
311 // ----- Setting ch. mapping -----
312 for(Int_t jk=0; jk<2*kNChannels; jk++){
313 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
314 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
315 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
316 rawStreamZDC->SetMapDet(jk, det[jk]);
317 rawStreamZDC->SetMapTow(jk, sec[jk]);
320 Float_t ZDCRawADC[4], ZDCCorrADC[4], ZDCCorrADCSum[4];
321 for(Int_t g=0; g<4; g++){
322 ZDCCorrADCSum[g] = 0.;
326 while(rawStreamZDC->Next()){
327 Int_t det = rawStreamZDC->GetSector(0);
328 Int_t quad = rawStreamZDC->GetSector(1);
330 if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
331 && !(rawStreamZDC->IsOverflow()) && det!=-1
332 && (rawStreamZDC->GetADCGain() == 1)){ // Selecting LOW RES ch.s
334 //printf(" IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
335 // rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
337 // Taking LOW RES channels -> channel+kNChannels !!!!
338 Int_t DetIndex=999, PedIndex=999;
339 if(det != 3 && quad != 5){ // Not ZEM nor PMRef
342 PedIndex = quad+kNChannels;
346 PedIndex = quad+5+kNChannels;
350 PedIndex = quad+12+kNChannels;
354 PedIndex = quad+17+kNChannels;
357 if(rawStreamZDC->GetADCGain() == 1 && (DetIndex!=999 || PedIndex!=999)){
359 ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
361 // Mean pedestal subtraction
362 Float_t Pedestal = MeanPed[PedIndex];
363 // Pedestal subtraction from correlation with out-of-time signals
364 //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
366 ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
367 ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
369 /*printf("\t det %d quad %d res %d pedInd %d detInd %d:"
370 "ADCCorr = %d, ZDCCorrADCSum = %d\n",
371 det,quad,rawStreamZDC->GetADCGain(),PedIndex,DetIndex,
372 (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);
375 if(DetIndex==999 || PedIndex==999)
376 printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
385 for(Int_t j=0; j<4; j++){
386 histoEMDRaw[j]->Fill(ZDCRawADC[j]);
387 histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
389 }//(if PHYSICS_EVENT)
399 /* Analysis of the histograms */
401 FILE *fileShuttle1 = fopen(ENCALIBDATA_FILE,"w");
406 Float_t MeanFitVal[4];
408 for(Int_t k=0; k<4; k++){
409 if(histoEMDCorr[k]->GetEntries() == 0){
410 printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
413 BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
414 YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
415 NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
416 // printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
417 histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
418 fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
419 MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
420 //printf("\n\t Mean Value from gaussian fit = %f\n", MeanFitVal[k]);
423 Float_t CalibCoeff[6];
426 for(Int_t j=0; j<6; j++){
428 CalibCoeff[j] = MeanFitVal[j];
429 fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
431 // ZEM energy calib. coeff. = 1
432 else if(j==4 || j==5){
434 fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
437 fclose(fileShuttle1);
439 FILE *fileShuttle2 = fopen(TOWCALIBDATA_FILE,"w");
440 for(Int_t j=0; j<4; j++){
441 // Note -> For the moment the inter-calibration coeff. are set to 1
442 for(Int_t k=0; k<5; k++){
444 fprintf(fileShuttle2,"\t%f",icoeff[k]);
445 if(k==4) fprintf(fileShuttle2,"\n");
448 fclose(fileShuttle2);
450 for(Int_t ij=0; ij<4; ij++){
451 delete histoEMDRaw[ij];
452 delete histoEMDCorr[ij];
456 TVirtualFitter::SetFitter(0);
459 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
461 /* close result file */
464 /* report progress */
465 daqDA_progressReport(90);
467 /* store the result file on FES */
468 status = daqDA_FES_storeFile(MAPDATA_FILE, MAPDATA_FILE);
470 printf("Failed to export file : %d\n",status);
474 status = daqDA_FES_storeFile(ENCALIBDATA_FILE, ENCALIBDATA_FILE);
476 printf("Failed to export file : %d\n",status);
480 status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, TOWCALIBDATA_FILE);
482 printf("Failed to export file : %d\n",status);
486 /* report progress */
487 daqDA_progressReport(100);