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 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
61 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
62 "Minuit", "TMinuitMinimizer(const char *)");
63 TVirtualFitter::SetDefaultFitter("Minuit");
65 //const Char_t* tableSOD[] = {"ALL", "no", "SOD", "all", NULL, NULL};
66 //monitorDeclareTable(const_cast<char**>(tableSOD));
68 char *monitor_table[] = { "ALL", "no", "PHY", "yes", "SOD", "all", NULL };
69 int err = monitorDeclareTable(monitor_table);
71 printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(err));
75 int status=0, nphys=0;
76 int const kNModules = 10;
77 int const kNChannels = 24;
78 int const kNScChannels = 32;
79 int const kZDCTDCGeo=4;
81 int itdc=0, iprevtdc=-1, ihittdc=0;
82 float tdcData[6], tdcL0=-999.;
83 for(int ij=0; ij<6; ij++) tdcData[ij]=-999.;
86 Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
87 for(Int_t kl=0; kl<kNModules; kl++){
88 modGeo[kl]=modType[kl]=modNCh[kl]=0;
92 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
93 Int_t det[2*kNChannels], sec[2*kNChannels];
94 for(Int_t y=0; y<2*kNChannels; y++){
95 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
99 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
100 Int_t scDet[kNScChannels], scSec[kNScChannels];
101 for(Int_t y=0; y<kNScChannels; y++){
102 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
106 Int_t tdcMod[kNScChannels], tdcCh[kNScChannels], tdcSigCode[kNScChannels];
107 Int_t tdcDet[kNScChannels], tdcSec[kNScChannels];
108 for(Int_t y=0; y<kNScChannels; y++){
109 tdcMod[y]=tdcCh[y]=tdcSigCode[y]=tdcDet[y]=tdcSec[y]=-1;
114 for(Int_t it=0; it<6; it++){
115 if(it==0) hTDC[it] = new TH1F("TDCZNC", "TDC ZNC", 200, -200., 200.);
116 else if(it==1) hTDC[it] = new TH1F("TDCZNA", "TDC ZNA", 200, -200., 200.);
117 else if(it==2) hTDC[it] = new TH1F("TDCZPC", "TDC ZPC", 200, -200., 200.);
118 else if(it==3) hTDC[it] = new TH1F("TDCZPA", "TDC ZPA", 200, -200., 200.);
119 else if(it==4) hTDC[it] = new TH1F("TDCZEM1","TDC ZEM1",200, -200., 200.);
120 else if(it==5) hTDC[it] = new TH1F("TDCZEM2","TDC ZEM2",200, -200., 200.);
123 /* log start of process */
124 printf("\n ZDC MAPPING program started\n");
125 signal(SIGSEGV, SIG_DFL);
127 /* check that we got some arguments = list of files */
129 printf("Wrong number of arguments\n");
133 FILE *mapFile4Shuttle;
135 /* read the data files */
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 struct eventHeaderStruct *event;
163 eventTypeType eventT;
165 /* check shutdown condition */
166 if (daqDA_checkShutdown()) {break;}
169 status=monitorGetEventDynamic((void **)&event);
170 if(status==MON_ERR_EOF){
171 printf ("End of File detected\n");
172 break; /* end of monitoring file has been reached */
175 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
179 /* retry if got no event */
180 if(event==NULL) continue;
182 // Initalize raw-data reading and decoding
183 AliRawReader *reader = new AliRawReaderDate((void*)event);
185 reader->Select("ZDC");
186 // --- Reading event header
187 //UInt_t evtype = reader->GetType();
188 //printf("\t ZDCMAPPINGda -> run # %d\n",reader->GetRunNumber());
190 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
193 /* use event - here, just write event id to result file */
194 eventT=event->eventType;
196 if(eventT==START_OF_DATA){
198 rawStreamZDC->SetSODReading(kTRUE);
200 // --------------------------------------------------------
201 // --- Writing ascii data file for the Shuttle preprocessor
202 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
203 //if(rawStreamZDC->Next()){
204 while((rawStreamZDC->Next())){
205 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
207 modGeo[iMod] = rawStreamZDC->GetADCModule();
208 modType[iMod] = rawStreamZDC->GetModType();
209 modNCh[iMod] = rawStreamZDC->GetADCNChannels();
211 if(rawStreamZDC->IsChMapping()){
212 if(modType[iMod]==1){ // ADC mapping ----------------------
213 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
214 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
215 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
216 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
217 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
220 else if(modType[iMod]==2){ // VME scaler mapping -------------
221 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
222 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
223 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
224 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
225 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
228 else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
229 tdcMod[itdcCh] = rawStreamZDC->GetTDCModFromMap(itdcCh);
230 tdcCh[itdcCh] = rawStreamZDC->GetTDCChFromMap(itdcCh);
231 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
236 // Writing data on output FXS file
237 for(Int_t is=0; is<2*kNChannels; is++){
238 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
239 is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
240 //printf(" Mapping DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
241 // is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
243 for(Int_t is=0; is<kNScChannels; is++){
244 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
245 is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
246 //if(scMod[is]!=-1) printf(" Mapping DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
247 // is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
249 for(Int_t is=0; is<kNScChannels; is++){
250 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
251 is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
252 //if(tdcMod[is]!=-1) printf(" Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
253 // is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
255 for(Int_t is=0; is<kNModules; is++){
256 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
257 modGeo[is],modType[is],modNCh[is]);
258 //printf(" Mapping DA -> Module mapping: geo %d type %d #ch %d\n",
259 // modGeo[is],modType[is],modNCh[is]);
263 fclose(mapFile4Shuttle);
265 else if(eventT==PHYSICS_EVENT){
267 rawStreamZDC->SetSODReading(kTRUE);
269 // ----- Setting ch. mapping -----
270 for(Int_t jk=0; jk<2*kNChannels; jk++){
271 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
272 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
273 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
274 rawStreamZDC->SetMapDet(jk, det[jk]);
275 rawStreamZDC->SetMapTow(jk, sec[jk]);
278 while(rawStreamZDC->Next()){
279 if(eventT==PHYSICS_EVENT && rawStreamZDC->GetADCModule()==kZDCTDCGeo && rawStreamZDC->IsZDCTDCDatum()==kTRUE){
281 itdc = rawStreamZDC->GetChannel();
282 if((itdc>=8 && itdc<=13) || itdc==15){
283 if(itdc==iprevtdc) ihittdc++;
286 if(ihittdc<1 && itdc!=15) tdcData[itdc-8] = 0.025*rawStreamZDC->GetZDCTDCDatum();
288 if(itdc==15 && ihittdc<1){
289 tdcL0 = 0.025*rawStreamZDC->GetZDCTDCDatum();
290 for(int ic=0; ic<6; ic++){
291 if(tdcData[ic]!=-999. && tdcL0!=-999.) hTDC[ic]->Fill(tdcData[ic]-tdcL0);
300 }//(if PHYSICS_EVENT)
301 else if(eventT==END_OF_RUN){
302 printf("End Of Run detected\n");
317 printf(" \n # of physics events analyzed = %d\n\n", nphys);
319 /* Analysis of the histograms */
322 fileShuttle = fopen(TDCDATA_FILE,"w");
324 Float_t xUp=0., xLow=0., deltaX=0;
325 Int_t binMax=0, nBinsx=0;
326 Float_t mean[6], sigma[6];
328 for(Int_t k=0; k<6; k++){
329 if(hTDC[k]->GetEntries()!=0){
330 binMax = hTDC[k]->GetMaximumBin();
332 printf("\n WARNING! Something wrong with det %d histo \n\n", k);
336 xUp = (hTDC[k]->GetXaxis())->GetXmax();
337 xLow = (hTDC[k]->GetXaxis())->GetXmin();
339 nBinsx = (hTDC[k]->GetXaxis())->GetNbins();
340 //printf(" xMax = %f\n", xLow+binMax*deltaX/nBinsx);
341 hTDC[k]->Fit("gaus","Q","",xLow+binMax*deltaX/nBinsx*0.75,xLow+binMax*deltaX/nBinsx*1.25);
342 fitfun[k] = hTDC[k]->GetFunction("gaus");
343 mean[k] = (Float_t) (fitfun[k]->GetParameter(1));
344 sigma[k] = (Float_t) (fitfun[k]->GetParameter(2));
345 //printf("\t Mean value from fit = %1.2f\n", mean[k]);
347 fprintf(fileShuttle,"\t%f\t%f\n",mean[k], sigma[k]);
349 else printf(" WARNING! TDC histo %d has no entries and can't be processed!\n",k);
354 TFile *histofile = new TFile(TDCHISTO_FILE,"RECREATE");
356 for(int k=0; k<6; k++) hTDC[k]->Write();
359 for(Int_t j=0; j<6; j++) delete hTDC[j];
361 /* store the result files on FES */
362 // [1] File with mapping
363 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
365 printf("Failed to export file : %d\n",status);
368 // [2] File with TDC data
369 status = daqDA_FES_storeFile(TDCDATA_FILE, "TDCDATA");
371 printf("Failed to export TDC data file to DAQ FES\n");
374 // [3] File with TDC histos
375 status = daqDA_FES_storeFile(TDCHISTO_FILE, "TDCHISTOS");
377 printf("Failed to export TDC histos file to DAQ FES\n");