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>
37 #include <TPluginManager.h>
44 #include "TMinuitMinimizer.h"
47 #include <AliRawReaderDate.h>
48 #include <AliRawEventHeaderBase.h>
49 #include <AliZDCRawStream.h>
51 int main(int argc, char **argv) {
53 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
60 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
61 "Minuit", "TMinuitMinimizer(const char *)");
62 TVirtualFitter::SetDefaultFitter("Minuit");
64 //const Char_t* tableSOD[] = {"ALL", "no", "SOD", "all", NULL, NULL};
65 //monitorDeclareTable(const_cast<char**>(tableSOD));
67 char *monitor_table[] = { "ALL", "no", "PHY", "yes", "SOD", "all", NULL };
68 err = monitorDeclareTable(monitor_table);
70 printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(err));
74 int status=0, nphys=0;
75 int const kNModules = 10;
76 int const kNChannels = 24;
77 int const kNScChannels = 32;
78 int const kZDCTDCGeo=4;
80 int itdc=0, iprevtdc=-1, ihittdc=0;
81 float tdcData[6], tdcL0=-999.;
82 for(int ij=0; ij<6; ij++) tdcData[ij]=-999.;
85 Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
86 for(Int_t kl=0; kl<kNModules; kl++){
87 modGeo[kl]=modType[kl]=modNCh[kl]=0;
91 Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
92 Int_t det[2*kNChannels], sec[2*kNChannels];
93 for(Int_t y=0; y<2*kNChannels; y++){
94 adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
98 Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
99 Int_t scDet[kNScChannels], scSec[kNScChannels];
100 for(Int_t y=0; y<kNScChannels; y++){
101 scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
105 Int_t tdcMod[kNScChannels], tdcCh[kNScChannels], tdcSigCode[kNScChannels];
106 Int_t tdcDet[kNScChannels], tdcSec[kNScChannels];
107 for(Int_t y=0; y<kNScChannels; y++){
108 tdcMod[y]=tdcCh[y]=tdcSigCode[y]=tdcDet[y]=tdcSec[y]=-1;
113 for(Int_t it=0; it<6; it++){
114 if(it==0) hTDC[it] = new TH1F("TDCZNC", "TDC ZNC", 200, -200., 200.);
115 else if(it==1) hTDC[it] = new TH1F("TDCZNA", "TDC ZNA", 200, -200., 200.);
116 else if(it==2) hTDC[it] = new TH1F("TDCZPC", "TDC ZPC", 200, -200., 200.);
117 else if(it==3) hTDC[it] = new TH1F("TDCZPA", "TDC ZPA", 200, -200., 200.);
118 else if(it==4) hTDC[it] = new TH1F("TDCZEM1","TDC ZEM1",200, -200., 200.);
119 else if(it==5) hTDC[it] = new TH1F("TDCZEM2","TDC ZEM2",200, -200., 200.);
122 /* log start of process */
123 printf("\n ZDC MAPPING program started\n");
125 /* check that we got some arguments = list of files */
127 printf("Wrong number of arguments\n");
131 FILE *mapFile4Shuttle;
133 /* read the data files */
137 status=monitorSetDataSource( argv[n] );
139 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
143 /* declare monitoring program */
144 status=monitorDeclareMp( __FILE__ );
146 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
150 monitorSetNoWaitNetworkTimeout(1000);
154 //Bool_t sodRead = kFALSE;
157 /* loop on events (infinite) */
160 struct eventHeaderStruct *event;
161 eventTypeType eventT;
163 /* check shutdown condition */
164 if (daqDA_checkShutdown()) {break;}
167 status=monitorGetEventDynamic((void **)&event);
168 if(status==MON_ERR_EOF){
169 printf ("End of File detected\n");
170 break; /* end of monitoring file has been reached */
173 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
177 /* retry if got no event */
178 if(event==NULL) continue;
180 // Initalize raw-data reading and decoding
181 AliRawReader *reader = new AliRawReaderDate((void*)event);
183 reader->Select("ZDC");
184 // --- Reading event header
185 //UInt_t evtype = reader->GetType();
186 //printf("\t ZDCMAPPINGda -> run # %d\n",reader->GetRunNumber());
188 AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
191 /* use event - here, just write event id to result file */
192 eventT=event->eventType;
194 if(eventT==START_OF_DATA){
196 rawStreamZDC->SetSODReading(kTRUE);
198 // --------------------------------------------------------
199 // --- Writing ascii data file for the Shuttle preprocessor
200 mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
201 //if(rawStreamZDC->Next()){
202 while((rawStreamZDC->Next())){
203 if(rawStreamZDC->IsHeaderMapping()){ // mapping header
205 modGeo[iMod] = rawStreamZDC->GetADCModule();
206 modType[iMod] = rawStreamZDC->GetModType();
207 modNCh[iMod] = rawStreamZDC->GetADCNChannels();
209 if(rawStreamZDC->IsChMapping()){
210 if(modType[iMod]==1){ // ADC mapping ----------------------
211 adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
212 adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
213 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
214 det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
215 sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
218 else if(modType[iMod]==2){ // VME scaler mapping -------------
219 scMod[iScCh] = rawStreamZDC->GetScalerModFromMap(iScCh);
220 scCh[iScCh] = rawStreamZDC->GetScalerChFromMap(iScCh);
221 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
222 scDet[iScCh] = rawStreamZDC->GetScDetectorFromMap(iScCh);
223 scSec[iScCh] = rawStreamZDC->GetScTowerFromMap(iScCh);
226 else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
227 tdcMod[itdcCh] = rawStreamZDC->GetTDCModFromMap(itdcCh);
228 tdcCh[itdcCh] = rawStreamZDC->GetTDCChFromMap(itdcCh);
229 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
234 // Writing data on output FXS file
235 for(Int_t is=0; is<2*kNChannels; is++){
236 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
237 is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
238 //printf(" Mapping DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
239 // is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
241 for(Int_t is=0; is<kNScChannels; is++){
242 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
243 is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
244 //if(scMod[is]!=-1) printf(" Mapping DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
245 // is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
247 for(Int_t is=0; is<kNScChannels; is++){
248 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
249 is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
250 //if(tdcMod[is]!=-1) printf(" Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
251 // is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
253 for(Int_t is=0; is<kNModules; is++){
254 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
255 modGeo[is],modType[is],modNCh[is]);
256 //printf(" Mapping DA -> Module mapping: geo %d type %d #ch %d\n",
257 // modGeo[is],modType[is],modNCh[is]);
261 fclose(mapFile4Shuttle);
263 else if(eventT==PHYSICS_EVENT){
265 rawStreamZDC->SetSODReading(kTRUE);
267 // ----- Setting ch. mapping -----
268 for(Int_t jk=0; jk<2*kNChannels; jk++){
269 rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
270 rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
271 rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
272 rawStreamZDC->SetMapDet(jk, det[jk]);
273 rawStreamZDC->SetMapTow(jk, sec[jk]);
276 while(rawStreamZDC->Next()){
277 if(rawStreamZDC->GetADCModule()==kZDCTDCGeo && rawStreamZDC->IsZDCTDCDatum()==kTRUE){
279 itdc = rawStreamZDC->GetChannel();
280 if((itdc>=8 && itdc<=13) || itdc==15){
281 if(itdc==iprevtdc) ihittdc++;
284 if(ihittdc<1 && itdc!=15) tdcData[itdc-8] = 0.025*rawStreamZDC->GetZDCTDCDatum();
286 if(itdc==15 && ihittdc<1){
287 tdcL0 = 0.025*rawStreamZDC->GetZDCTDCDatum();
288 for(int ic=0; ic<6; ic++){
289 if(tdcData[ic]!=-999. && tdcL0!=-999.) hTDC[ic]->Fill(tdcData[ic]-tdcL0);
298 }//(if PHYSICS_EVENT)
299 else if(eventT==END_OF_RUN){
300 printf("End Of Run detected\n");
315 printf(" \n # of physics events analyzed = %d\n\n", nphys);
317 /* Analysis of the histograms */
320 fileShuttle = fopen(TDCDATA_FILE,"w");
322 Float_t xUp=0., xLow=0., deltaX=0;
323 Int_t binMax=0, nBinsx=0;
324 Float_t mean[6], sigma[6];
326 for(Int_t k=0; k<6; k++){
327 if(hTDC[k]->GetEntries()!=0){
328 binMax = hTDC[k]->GetMaximumBin();
330 printf("\n WARNING! Something wrong with det %d histo \n\n", k);
334 xUp = (hTDC[k]->GetXaxis())->GetXmax();
335 xLow = (hTDC[k]->GetXaxis())->GetXmin();
337 nBinsx = (hTDC[k]->GetXaxis())->GetNbins();
338 //printf(" xMax = %f\n", xLow+binMax*deltaX/nBinsx);
339 hTDC[k]->Fit("gaus","Q","",xLow+binMax*deltaX/nBinsx*0.6,xLow+binMax*deltaX/nBinsx*1.24);
340 fitfun[k] = hTDC[k]->GetFunction("gaus");
341 mean[k] = (Float_t) (fitfun[k]->GetParameter(1));
342 sigma[k] = (Float_t) (fitfun[k]->GetParameter(2));
343 //printf("\t Mean value from fit = %1.2f\n", mean[k]);
345 fprintf(fileShuttle,"\t%f\t%f\n",mean[k], sigma[k]);
347 else printf(" WARNING! TDC histo %d has no entries and can't be processed!\n",k);
352 TFile *histofile = new TFile(TDCHISTO_FILE,"RECREATE");
354 for(int k=0; k<6; k++) hTDC[k]->Write();
357 for(Int_t j=0; j<6; j++) delete hTDC[j];
359 /* store the result files on FES */
360 // [1] File with mapping
361 status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
363 printf("Failed to export file : %d\n",status);
366 // [2] File with TDC data
367 status = daqDA_FES_storeFile(TDCDATA_FILE, "TDCDATA");
369 printf("Failed to export TDC data file to DAQ FES\n");
372 // [3] File with TDC histos
373 status = daqDA_FES_storeFile(TDCHISTO_FILE, "TDCHISTOS");
375 printf("Failed to export TDC histos file to DAQ FES\n");