3 This program reads the DAQ data files passed as argument using the monitoring library.
5 The program reports about its processing progress.
7 Messages on stdout are exported to DAQ log system.
9 DA to write mapping for ADC modules and VME scaler
11 Contact: Chiara.Oppedisano@to.infn.it
15 Number of events needed: at least 50%
17 Output Files: ZDCChMapping.dat
18 Trigger Types Used: different trigger types are used
22 #define MAPDATA_FILE "ZDCChMapping.dat"
23 #define TDCDATA_FILE "ZDCTDCCalib.dat"
24 #define TDCHISTO_FILE "ZDCTDCHisto.root"
28 #include <Riostream.h>
38 #include <TPluginManager.h>
45 #include "TMinuitMinimizer.h"
48 #include <AliRawReaderDate.h>
49 #include <AliRawEventHeaderBase.h>
50 #include <AliZDCRawStream.h>
52 int main(int argc, char **argv) {
54 // needed for streamer application
55 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
60 // needed for Minuit plugin
61 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
65 "TMinuitMinimizer(const char*)");
68 TVirtualFitter::SetDefaultFitter("Minuit");
70 char *monitor_table[] = { "ALL", "no", "PHY", "yes", "SOD", "all", NULL };
71 int err = monitorDeclareTable(monitor_table);
73 printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(err));
77 int status=0, nphys=0;
78 int const kNModules = 9;
79 int const kNChannels = 24;
80 int const kNScChannels = 32;
81 int const kZDCTDCGeo = 4;
83 int itdc=0, iprevtdc=-1, ihittdc=0;
84 float tdcData[6], tdcL0=-999.;
85 for(int ij=0; ij<6; ij++) tdcData[ij]=-999.;
88 Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
89 for(Int_t kl=0; kl<kNModules; kl++){
90 modGeo[kl]=modType[kl]=modNCh[kl]=0;
94 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
95 Int_t det[2*kNChannels], sec[2*kNChannels];
96 for(Int_t y=0; y<2*kNChannels; y++){
97 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
101 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
102 Int_t scDet[kNScChannels], scSec[kNScChannels];
103 for(Int_t y=0; y<kNScChannels; y++){
104 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
108 Int_t tdcMod[kNScChannels], tdcCh[kNScChannels], tdcSigCode[kNScChannels];
109 Int_t tdcDet[kNScChannels], tdcSec[kNScChannels];
110 for(Int_t y=0; y<kNScChannels; y++){
111 tdcMod[y]=tdcCh[y]=tdcSigCode[y]=tdcDet[y]=tdcSec[y]=-1;
114 TH1F * hTDC[6]={0x0,0x0,0x0,0x0,0x0,0x0};
115 for(int it=0; it<6; it++){
116 if(it==0) hTDC[it] = new TH1F("TDCZEM1","TDC ZEM1",200, -200., 200.);
117 else if(it==1) hTDC[it] = new TH1F("TDCZEM2","TDC ZEM2",200, -200., 200.);
118 else if(it==2) hTDC[it] = new TH1F("TDCZNC", "TDC ZNC", 200, -200., 200.);
119 else if(it==3) hTDC[it] = new TH1F("TDCZPC", "TDC ZPC", 200, -200., 200.);
120 else if(it==4) hTDC[it] = new TH1F("TDCZNA", "TDC ZNA", 200, -200., 200.);
121 else if(it==5) hTDC[it] = new TH1F("TDCZPA", "TDC ZPA", 200, -200., 200.);
124 /* log start of process */
125 printf("\n ZDC MAPPING program started\n");
126 signal(SIGSEGV, SIG_DFL);
128 /* check that we got some arguments = list of files */
130 printf("Wrong number of arguments\n");
134 FILE *mapFile4Shuttle;
136 /* read the data files */
137 for(int n=1;n<argc;n++){
139 status=monitorSetDataSource( argv[n] );
141 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
145 /* declare monitoring program */
146 status=monitorDeclareMp( __FILE__ );
148 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
152 monitorSetNoWaitNetworkTimeout(1000);
156 //Bool_t sodRead = kFALSE;
159 /* loop on events (infinite) */
162 if(nphys > 50000) break;
164 struct eventHeaderStruct *event;
165 eventTypeType eventT;
167 /* check shutdown condition */
168 if (daqDA_checkShutdown()) {break;}
171 status=monitorGetEventDynamic((void **)&event);
172 if(status==MON_ERR_EOF){
173 printf ("End of File detected\n");
174 break; /* end of monitoring file has been reached */
177 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
181 /* retry if got no event */
182 if(event==NULL) continue;
184 // Initalize raw-data reading and decoding
185 AliRawReader *reader = new AliRawReaderDate((void*)event);
187 reader->Select("ZDC");
188 // --- Reading event header
189 //UInt_t evtype = reader->GetType();
190 //printf("\t ZDCMAPPINGda -> run # %d\n",reader->GetRunNumber());
192 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
195 /* use event - here, just write event id to result file */
196 eventT=event->eventType;
198 if(eventT==START_OF_DATA){
200 rawStreamZDC->SetSODReading(kTRUE);
202 // --------------------------------------------------------
203 // --- Writing ascii data file for the Shuttle preprocessor
204 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
205 //if(rawStreamZDC->Next()){
206 while((rawStreamZDC->Next())){
207 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
209 modGeo[iMod] = rawStreamZDC->GetADCModule();
210 modType[iMod] = rawStreamZDC->GetModType();
211 modNCh[iMod] = rawStreamZDC->GetADCNChannels();
213 if(rawStreamZDC->IsChMapping()){
214 if(modType[iMod]==1){ // ADC mapping ----------------------
215 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
216 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
217 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
218 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
219 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
222 else if(modType[iMod]==2){ // VME scaler mapping -------------
223 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
224 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
225 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
226 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
227 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
230 else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
231 tdcMod[itdcCh] = rawStreamZDC->GetTDCModFromMap(itdcCh);
232 tdcCh[itdcCh] = rawStreamZDC->GetTDCChFromMap(itdcCh);
233 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
238 // Writing data on output FXS file
239 for(Int_t is=0; is<2*kNChannels; is++){
240 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
241 is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
242 //printf(" Mapping DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
243 //is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
245 for(Int_t is=0; is<kNScChannels; is++){
246 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
247 is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
248 //if(scMod[is]!=-1) printf(" Mapping DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
249 //is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
251 for(Int_t is=0; is<kNScChannels; is++){
252 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
253 is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
254 //if(tdcMod[is]!=-1) printf(" Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
255 //is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
257 for(Int_t is=0; is<kNModules; is++){
258 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
259 modGeo[is],modType[is],modNCh[is]);
260 //printf(" Mapping DA -> Module mapping: geo %d type %d #ch %d\n",
261 //modGeo[is],modType[is],modNCh[is]);
265 fclose(mapFile4Shuttle);
267 else if(eventT==PHYSICS_EVENT){
269 for(int ij=0; ij<6; ij++) tdcData[ij]=-999.;
270 rawStreamZDC->SetSODReading(kTRUE);
272 // ----- Setting ch. mapping -----
273 for(int jk=0; jk<2*kNChannels; jk++){
274 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
275 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
276 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
277 rawStreamZDC->SetMapDet(jk, det[jk]);
278 rawStreamZDC->SetMapTow(jk, sec[jk]);
281 while(rawStreamZDC->Next()){
282 if(rawStreamZDC->GetADCModule()!=kZDCTDCGeo) continue; //skipping ADCs, scalers and trigger cards
284 if(rawStreamZDC->GetADCModule()==kZDCTDCGeo && rawStreamZDC->IsZDCTDCDatum()==kTRUE){
286 itdc = rawStreamZDC->GetChannel();
287 if(itdc==iprevtdc) ihittdc++;
291 if(((itdc>=8 && itdc<=13) || itdc==15) && ihittdc==0){
293 tdcData[itdc-8] = 0.025*rawStreamZDC->GetZDCTDCDatum();
295 //printf(" ev.%d **** TDC%d %1.0f ns \n",nphys,itdc, tdcData[itdc-8]);
299 tdcL0 = 0.025*rawStreamZDC->GetZDCTDCDatum();
301 //printf(" ev.%d ++++ TDCL0 %1.0f ns \n",nphys,tdcL0);
303 for(int ic=0; ic<6; ic++){
304 if(tdcData[ic]!=-999. && tdcL0!=-999.){
305 hTDC[ic]->Fill(tdcData[ic]-tdcL0);
306 //printf(" ev.%d -> Filling histo%d: %f ns\n",nphys,ic, tdcData[ic]-tdcL0);
320 }//(if PHYSICS_EVENT)
321 else if(eventT==END_OF_RUN){
322 printf("End Of Run detected\n");
334 printf(" \n # of physics events analyzed = %d\n\n", nphys);
336 /* Analysis of the histograms */
339 fileShuttle = fopen(TDCDATA_FILE,"w");
341 Float_t xUp=0., xLow=0., deltaX=0;
342 Int_t binMax=0, nBinsx=0;
343 Float_t mean[6]={0.,0.,0.,0.,0.,0.}, sigma[6]={0.,0.,0.,0.,0.,0.};
344 TF1 *fitfun[6]={0x0,0x0,0x0,0x0,0x0,0x0};
345 for(int k=0; k<6; k++){
346 if(hTDC[k]->GetEntries()!=0){
347 binMax = hTDC[k]->GetMaximumBin();
349 printf("\n WARNING! Something wrong with det %d histo \n\n", k);
353 xUp = (hTDC[k]->GetXaxis())->GetXmax();
354 xLow = (hTDC[k]->GetXaxis())->GetXmin();
356 nBinsx = (hTDC[k]->GetXaxis())->GetNbins();
357 //printf(" xMax = %f\n", xLow+binMax*deltaX/nBinsx);
358 hTDC[k]->Fit("gaus","Q","",xLow+binMax*deltaX/nBinsx*0.75,xLow+binMax*deltaX/nBinsx*1.25);
359 fitfun[k] = hTDC[k]->GetFunction("gaus");
360 mean[k] = (Float_t) (fitfun[k]->GetParameter(1));
361 sigma[k] = (Float_t) (fitfun[k]->GetParameter(2));
362 //printf("\t Mean value from fit = %1.2f ns\n", mean[k]);
364 fprintf(fileShuttle,"\t%f\t%f\n",mean[k], sigma[k]);
366 else printf(" WARNING! TDC histo %d has no entries and can't be processed!\n",k);
371 TFile *histofile = new TFile(TDCHISTO_FILE,"RECREATE");
373 for(int k=0; k<6; k++) hTDC[k]->Write();
376 for(Int_t j=0; j<6; j++) delete hTDC[j];
378 /* store the result files on FES */
379 // [1] File with mapping
380 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
382 printf("Failed to export file : %d\n",status);
385 // [2] File with TDC data
386 status = daqDA_FES_storeFile(TDCDATA_FILE, "TDCDATA");
388 printf("Failed to export TDC data file to DAQ FES\n");
391 // [3] File with TDC histos
392 status = daqDA_FES_storeFile(TDCHISTO_FILE, "TDCHISTOS");
394 printf("Failed to export TDC histos file to DAQ FES\n");