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 MBCALIBDATA_FILE "ZDCMBCalib.root"
29 #include <Riostream.h>
30 #include <Riostream.h>
39 #include <TPluginManager.h>
46 #include "TMinuitMinimizer.h"
49 #include <AliRawReaderDate.h>
50 #include <AliRawEventHeaderBase.h>
51 #include <AliZDCRawStream.h>
56 1- monitoring data source
58 int main(int argc, char **argv) {
61 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
68 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
69 "Minuit", "TMinuitMinimizer(const char *)");
70 TVirtualFitter::SetDefaultFitter("Minuit");
73 // No. of ZDC cabled ch.
74 int const kNChannels = 24;
75 int const kNScChannels = 32;
76 Int_t kFirstADCGeo=0, kLastADCGeo=3;
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");
87 // --- Preparing histos for EM dissociation spectra
89 TH2F * hZDCvsZEM = new TH2F("hZDCvsZEM","EZDCTot vs. EZEM",100,0.,5.,100,0.,800.);
90 TH2F * hZDCCvsZEM = new TH2F("hZDCCvsZEM","EZDCC vs. EZEM",100,0.,5.,100,0.,400.);
91 TH2F * hZDCAvsZEM = new TH2F("hZDCAvsZEM","EZDCA vs. EZEM",100,0.,5.,100,0.,400.);
93 /* open result file */
95 fp=fopen("./result.txt","a");
97 printf("Failed to open file\n");
101 FILE *mapFile4Shuttle;
103 // *** To analyze CALIBRATION_MB events a pedestal data file is needed!!!
104 // *** -> check if a pedestal run has been analyzed
106 read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
108 printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
111 else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
113 FILE *filePed = fopen(PEDDATA_FILE,"r");
115 printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
119 // 144 = 48 in-time + 48 out-of-time + 48 correlations
120 Float_t readValues[2][6*kNChannels];
121 Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
122 Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
123 // ***************************************************
124 // Unless we have a narrow correlation to fit we
125 // don't fit and store in-time vs. out-of-time
126 // histograms -> mean pedstal subtracted!!!!!!
127 // ***************************************************
129 for(int jj=0; jj<6*kNChannels; jj++){
130 for(int ii=0; ii<2; ii++){
131 fscanf(filePed,"%f",&readValues[ii][jj]);
134 MeanPedhg[jj] = readValues[0][jj];
135 //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
137 else if(jj>=kNChannels && jj<2*kNChannels){
138 MeanPedlg[jj-kNChannels] = readValues[0][jj];
139 //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
141 else if(jj>4*kNChannels){
142 CorrCoeff0[jj-4*kNChannels] = readValues[0][jj];
143 CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
147 /* report progress */
148 daqDA_progressReport(10);
151 /* init some counters */
152 int nevents_physics=0;
155 struct eventHeaderStruct *event;
156 eventTypeType eventT;
158 /* read the data files */
162 status=monitorSetDataSource( argv[n] );
164 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
168 /* report progress */
169 /* in this example, indexed on the number of files */
170 daqDA_progressReport(10+80*n/argc);
176 status=monitorGetEventDynamic((void **)&event);
177 if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
179 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
183 /* retry if got no event */
188 // Initalize raw-data reading and decoding
189 AliRawReader *reader = new AliRawReaderDate((void*)event);
190 reader->Select("ZDC");
191 // --- Reading event header
192 //UInt_t evtype = reader->GetType();
193 //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
194 //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
196 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
199 /* use event - here, just write event id to result file */
200 eventT=event->eventType;
203 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
204 Int_t det[2*kNChannels], sec[2*kNChannels];
205 for(Int_t y=0; y<2*kNChannels; y++){
206 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
210 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
211 Int_t scDet[kNScChannels], scSec[kNScChannels];
212 for(Int_t y=0; y<kNScChannels; y++){
213 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
216 Int_t modNum=-1, modType=-1;
218 if(eventT==START_OF_DATA){
220 rawStreamZDC->SetSODReading(kTRUE);
222 // --------------------------------------------------------
223 // --- Writing ascii data file for the Shuttle preprocessor
224 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
225 if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
227 while((rawStreamZDC->Next())){
228 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
229 modNum = rawStreamZDC->GetADCModule();
230 modType = rawStreamZDC->GetModType();
232 if(rawStreamZDC->IsChMapping()){
233 if(modType==1){ // ADC mapping ----------------------
234 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
235 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
236 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
237 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
238 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
240 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
241 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
243 //printf(" Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
244 // ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
248 else if(modType==2){ //VME scaler mapping --------------------
249 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
250 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
251 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
252 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
253 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
255 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
256 iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
258 //printf(" Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
259 // iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
266 fclose(mapFile4Shuttle);
269 if(eventT==PHYSICS_EVENT){
270 // --- Reading data header
271 reader->ReadHeader();
272 const AliRawDataHeader* header = reader->GetDataHeader();
274 UChar_t message = header->GetAttributes();
275 if((message & 0xf0) == 0x60){ // DEDICATED CALIBRATION MB RUN
276 //printf("\t STANDALONE_CALIBRAION_MB_RUN raw data found\n");
280 printf("\t NO STANDALONE_MB_RUN raw data found\n");
285 printf("\t ATTENTION! No Raw Data Header found!!!\n");
289 rawStreamZDC->SetSODReading(kTRUE);
291 if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
293 // ----- Setting ch. mapping -----
294 for(Int_t jk=0; jk<2*kNChannels; jk++){
295 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
296 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
297 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
298 rawStreamZDC->SetMapDet(jk, det[jk]);
299 rawStreamZDC->SetMapTow(jk, sec[jk]);
302 Float_t rawZDC=0., rawZEM=0.;
303 Float_t corrZDC=0., corrZEM=0.;
304 Float_t corrZDCTot=0., corrZDCC=0., corrZDCA=0., corrZEMTot=0.;
306 while(rawStreamZDC->Next()){
307 Int_t detector = rawStreamZDC->GetSector(0);
308 Int_t quad = rawStreamZDC->GetSector(1);
311 if( (rawStreamZDC->IsADCDataWord()) && !(rawStreamZDC->IsUnderflow())
312 && !(rawStreamZDC->IsOverflow()) && (detector!=-1)
313 && ((rawStreamZDC->GetADCGain()) == 0) && // Selecting HIGH RES ch.s
314 (rawStreamZDC->GetADCModule()>=kFirstADCGeo) && (rawStreamZDC->GetADCModule()<=kLastADCGeo)){
316 //printf(" IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
317 // rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
319 if(quad!=5){ // Physics signals
320 if(detector==1) index = quad; // *** ZNC
321 else if(detector==2) index = quad+5; // *** ZPC
322 else if(detector==3) index = quad+9; // *** ZEM
323 else if(detector==4) index = quad+12;// *** ZNA
324 else if(detector==5) index = quad+17;// *** ZPA
328 if(index==-1) printf("ERROR in ZDCCALIBMBda.cxx -> det %d quad %d index %d ADC %d\n",
329 detector, quad, index, rawStreamZDC->GetADCValue());
331 // Mean pedestal subtraction
332 Float_t Pedestal = MeanPedhg[index];
333 // Pedestal subtraction from correlation with out-of-time signals
334 //Float_t Pedestal = CorrCoeff0[index]+CorrCoeff1[index]*MeanPedOOT[index];
336 if(index!=-1 && quad!=5){
339 rawZEM = (Float_t) rawStreamZDC->GetADCValue();
340 corrZEM = (rawStreamZDC->GetADCValue()) - Pedestal;
341 corrZEMTot += corrZEM;
344 rawZDC = (Float_t) rawStreamZDC->GetADCValue();
345 corrZDC = rawZDC - Pedestal;
346 if(detector==1 || detector==2) corrZDCC += corrZDC;
347 else corrZDCA += corrZDC;
348 corrZDCTot += corrZDC;
355 hZDCvsZEM ->Fill(corrZDCTot, corrZEMTot);
356 hZDCCvsZEM->Fill(corrZDCC, corrZEMTot);
357 hZDCAvsZEM->Fill(corrZDCA, corrZEMTot);
360 }//(if PHYSICS_EVENT)
371 TFile *histofile = new TFile(MBCALIBDATA_FILE,"RECREATE");
378 if(hZDCvsZEM) delete hZDCvsZEM ;
379 if(hZDCCvsZEM) delete hZDCCvsZEM;
380 if(hZDCAvsZEM) delete hZDCAvsZEM;
383 fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
385 /* close result file */
388 /* report progress */
389 daqDA_progressReport(90);
391 /* store the result file on FES */
392 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
394 printf("Failed to export file : %d\n",status);
398 status = daqDA_FES_storeFile(MBCALIBDATA_FILE, "MBCALIB");
400 printf("Failed to export file : %d\n",status);
404 /* report progress */
405 daqDA_progressReport(100);