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_MB runs
14 Contact: Chiara.Oppedisano@to.infn.it
16 Run Type: CALIBRATION_MB_RUN
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
24 #define PEDDATA_FILE "ZDCPedestal.dat"
25 #define MAPDATA_FILE "ZDCChMapping.dat"
26 #define ENCALIBDATA_FILE "ZDCEnergyCalib.dat"
27 #define MBCALIBDATA_FILE "ZDCMBCalib.root"
31 #include <Riostream.h>
32 #include <Riostream.h>
41 #include <TPluginManager.h>
48 #include "TMinuitMinimizer.h"
51 #include <AliRawReaderDate.h>
52 #include <AliRawEventHeaderBase.h>
53 #include <AliZDCRawStream.h>
58 1- monitoring data source
60 int main(int argc, char **argv) {
63 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
70 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
71 "Minuit", "TMinuitMinimizer(const char *)");
72 TVirtualFitter::SetDefaultFitter("Minuit");
75 int const kNModules = 9;
76 int const kNChannels = 24;
77 int const kNScChannels = 32;
78 Int_t kFirstADCGeo=0, kLastADCGeo=1;
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;
101 Int_t tdcMod[kNScChannels], tdcCh[kNScChannels], tdcSigCode[kNScChannels];
102 Int_t tdcDet[kNScChannels], tdcSec[kNScChannels];
103 for(Int_t y=0; y<kNScChannels; y++){
104 tdcMod[y]=tdcCh[y]=tdcSigCode[y]=tdcDet[y]=tdcSec[y]=-1;
107 /* log start of process */
108 printf("ZDC EMD program started\n");
110 /* check that we got some arguments = list of files */
112 printf("Wrong number of arguments\n");
116 // --- Preparing histos for EM dissociation spectra
118 TH2F * hZDCvsZEM = new TH2F("hZDCvsZEM","EZDCTot vs. EZEM",100,0.,8.,100,0.,800.);
119 TH2F * hZDCCvsZEM = new TH2F("hZDCCvsZEM","EZDCC vs. EZEM",100,0.,8.,100,0.,400.);
120 TH2F * hZDCAvsZEM = new TH2F("hZDCAvsZEM","EZDCA vs. EZEM",100,0.,8.,100,0.,400.);
122 /* open result file */
124 fp=fopen("./result.txt","a");
126 printf("Failed to open file\n");
130 FILE *mapFile4Shuttle;
132 // *** To analyze CALIBRATION_MB events a pedestal data file is needed!!!
133 // *** -> check if a pedestal run has been analyzed
135 read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
137 printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
140 else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
142 FILE *filePed = fopen(PEDDATA_FILE,"r");
144 printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
148 // 144 = 48 in-time + 48 out-of-time + 48 correlations
149 Float_t readValues[2][6*kNChannels];
150 Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
151 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
152 // ***************************************************
153 // Unless we have a narrow correlation to fit we
154 // don't fit and store in-time vs. out-of-time
155 // histograms -> mean pedstal subtracted!!!!!!
156 // ***************************************************
158 for(int jj=0; jj<6*kNChannels; jj++){
159 for(int ii=0; ii<2; ii++){
160 fscanf(filePed,"%f",&readValues[ii][jj]);
163 MeanPedhg[jj] = readValues[0][jj];
164 //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
166 else if(jj>=kNChannels && jj<2*kNChannels){
167 MeanPedlg[jj-kNChannels] = readValues[0][jj];
168 //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
170 else if(jj>4*kNChannels){
171 CorrCoeff0[jj-4*kNChannels] = readValues[0][jj];
172 CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
176 /* report progress */
177 daqDA_progressReport(10);
179 FILE *fileEnCal = fopen(ENCALIBDATA_FILE,"r");
180 if(fileEnCal==NULL) {
181 printf("\t ERROR!!! Can't open ZDCEnergyCalib.dat file!!!\n");
185 Float_t calibCoeff[6];
186 for(Int_t irow=0; irow<6; irow++){
187 fscanf(fileEnCal,"%f",&calibCoeff[irow]);
190 /* init some counters */
191 int nevents_physics=0;
194 struct eventHeaderStruct *event;
195 eventTypeType eventT;
197 /* read the data files */
201 status=monitorSetDataSource( argv[n] );
203 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
207 /* report progress */
208 /* in this example, indexed on the number of files */
209 daqDA_progressReport(10+80*n/argc);
215 status=monitorGetEventDynamic((void **)&event);
216 if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
218 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
222 /* retry if got no event */
227 // Initalize raw-data reading and decoding
228 AliRawReader *reader = new AliRawReaderDate((void*)event);
229 reader->Select("ZDC");
230 // --- Reading event header
231 //UInt_t evtype = reader->GetType();
232 //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
233 //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
235 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
238 /* use event - here, just write event id to result file */
239 eventT=event->eventType;
241 if(eventT==START_OF_DATA){
243 rawStreamZDC->SetSODReading(kTRUE);
245 // --------------------------------------------------------
246 // --- Writing ascii data file for the Shuttle preprocessor
247 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
248 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
250 while((rawStreamZDC->Next())){
251 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
253 modGeo[iMod] = rawStreamZDC->GetADCModule();
254 modType[iMod] = rawStreamZDC->GetModType();
255 modNCh[iMod] = rawStreamZDC->GetADCNChannels();
257 if(rawStreamZDC->IsChMapping()){
258 if(modType[iMod]==1){ // ADC mapping ----------------------
259 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
260 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
261 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
262 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
263 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
266 else if(modType[iMod]==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);
274 else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
275 tdcMod[itdcCh] = rawStreamZDC->GetTDCModFromMap(itdcCh);
276 tdcCh[itdcCh] = rawStreamZDC->GetTDCChFromMap(itdcCh);
277 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
282 // Writing data on output FXS file
283 for(Int_t is=0; is<2*kNChannels; is++){
284 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
285 is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
286 //printf(" CalibMB DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
287 // is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
289 for(Int_t is=0; is<kNScChannels; is++){
290 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
291 is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
292 //printf(" CalibMB DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
293 // is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
295 for(Int_t is=0; is<kNScChannels; is++){
296 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
297 is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
298 //if(tdcMod[is]!=-1) printf(" Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
299 // is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
301 for(Int_t is=0; is<kNModules; is++){
302 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
303 modGeo[is],modType[is],modNCh[is]);
304 //printf(" CalibMB DA -> Module mapping: geo %d type %d #ch %d\n",
305 // modGeo[is],modType[is],modNCh[is]);
309 fclose(mapFile4Shuttle);
312 if(eventT==PHYSICS_EVENT){
313 // --- Reading data header
314 reader->ReadHeader();
315 const AliRawDataHeader* header = reader->GetDataHeader();
317 UChar_t message = header->GetAttributes();
318 if((message & 0xf0) == 0x60){ // DEDICATED CALIBRATION MB RUN
319 //printf("\t CALIBRATION_MB_RUN raw data found\n");
323 //printf("\t NO CALIBRAION_MB_RUN raw data found\n");
328 printf("\t ATTENTION! No Raw Data Header found!!!\n");
332 rawStreamZDC->SetSODReading(kTRUE);
334 if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
336 // ----- Setting ch. mapping -----
337 for(Int_t jk=0; jk<2*kNChannels; jk++){
338 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
339 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
340 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
341 rawStreamZDC->SetMapDet(jk, det[jk]);
342 rawStreamZDC->SetMapTow(jk, sec[jk]);
345 Float_t rawZNC=0., rawZPC=0., rawZNA=0., rawZPA=0., rawZEM1=0., rawZEM2=0.;
346 Float_t corrZNC=0., corrZPC=0., corrZNA=0., corrZPA=0., corrZEM1=0.,corrZEM2=0.;
347 Float_t calZNC=0., calZPC=0., calZNA=0., calZPA=0., calZEM1=0., calZEM2=0.;
348 Float_t calZDCTot=0., calZDCC=0., calZDCA=0., calZEM=0.;
350 while(rawStreamZDC->Next()){
351 Int_t detector = rawStreamZDC->GetSector(0);
352 Int_t quad = rawStreamZDC->GetSector(1);
355 if( (rawStreamZDC->IsADCDataWord()) && !(rawStreamZDC->IsUnderflow())
356 && !(rawStreamZDC->IsOverflow()) && (detector!=-1)
357 && ((rawStreamZDC->GetADCGain()) == 0) && // Selecting HIGH RES ch.s
358 (rawStreamZDC->GetADCModule()>=kFirstADCGeo) && (rawStreamZDC->GetADCModule()<=kLastADCGeo)){
360 //printf(" IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
361 // rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
363 if(quad!=5){ // Physics signals
364 if(detector==1) index = quad; // *** ZNC
365 else if(detector==2) index = quad+5; // *** ZPC
366 else if(detector==3) index = quad+9; // *** ZEM
367 else if(detector==4) index = quad+12;// *** ZNA
368 else if(detector==5) index = quad+17;// *** ZPA
372 if(index==-1) printf("ERROR in ZDCCALIBMBda.cxx -> det %d quad %d index %d ADC %d\n",
373 detector, quad, index, rawStreamZDC->GetADCValue());
375 // Mean pedestal subtraction
376 Float_t Pedestal = MeanPedhg[index];
377 // Pedestal subtraction from correlation with out-of-time signals
378 //Float_t Pedestal = CorrCoeff0[index]+CorrCoeff1[index]*MeanPedOOT[index];
380 if(index!=-1 && quad!=5){
383 rawZNC = (Float_t) rawStreamZDC->GetADCValue();
384 corrZNC = rawZNC - Pedestal;
386 else if(detector==2){
387 rawZPC = (Float_t) rawStreamZDC->GetADCValue();
388 corrZPC = rawZPC - Pedestal;
390 else if(detector==3){
392 rawZEM1 = (Float_t) rawStreamZDC->GetADCValue();
393 corrZEM1 = (rawStreamZDC->GetADCValue()) - Pedestal;
396 rawZEM2 = (Float_t) rawStreamZDC->GetADCValue();
397 corrZEM2 = (rawStreamZDC->GetADCValue()) - Pedestal;
400 else if(detector==4){
401 rawZNA = (Float_t) rawStreamZDC->GetADCValue();
402 corrZNA = rawZNA - Pedestal;
404 else if(detector==5){
405 rawZPA = (Float_t) rawStreamZDC->GetADCValue();
406 corrZPA = rawZPA - Pedestal;
413 calZNC = calibCoeff[0]*corrZNC;
414 calZPC = calibCoeff[1]*corrZPC;
415 calZNA = calibCoeff[2]*corrZNA;
416 calZPA = calibCoeff[3]*corrZPA;
417 calZEM1 = calibCoeff[4]*corrZEM1;
418 calZEM2 = calibCoeff[4]*corrZEM2;
419 calZDCTot = calZNC+calZPC+calZPA+calZNA;
420 calZDCC = calZNC+calZPC;
421 calZDCA = calZNA+calZPA;
422 calZEM = calZEM1+calZEM2;
424 hZDCvsZEM ->Fill(calZDCTot/1000, calZEM/1000);
425 hZDCCvsZEM->Fill(calZDCC/1000, calZEM/1000);
426 hZDCAvsZEM->Fill(calZDCA/1000, calZEM/1000);
429 }//(if PHYSICS_EVENT)
440 TFile *histofile = new TFile(MBCALIBDATA_FILE,"RECREATE");
447 if(hZDCvsZEM) delete hZDCvsZEM ;
448 if(hZDCCvsZEM) delete hZDCCvsZEM;
449 if(hZDCAvsZEM) delete hZDCAvsZEM;
452 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
454 /* close result file */
457 /* report progress */
458 daqDA_progressReport(90);
460 /* store the result file on FES */
461 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
463 printf("Failed to export file : %d\n",status);
467 status = daqDA_FES_storeFile(MBCALIBDATA_FILE, "MBCALIB");
469 printf("Failed to export file : %d\n",status);
473 /* report progress */
474 daqDA_progressReport(100);