]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCMAPPINGda.cxx
Adding missing include file
[u/mrichter/AliRoot.git] / ZDC / ZDCMAPPINGda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5 The program reports about its processing progress.
6
7 Messages on stdout are exported to DAQ log system.
8
9 DA to write mapping for ADC modules and VME scaler
10
11 Contact: Chiara.Oppedisano@to.infn.it
12 Link: 
13 Run Type: PHYSICS
14 DA Type: MON
15 Number of events needed: at least 50% 
16 Input Files:  none
17 Output Files: ZDCChMapping.dat
18 Trigger Types Used: different trigger types are used
19
20 */
21
22 #define MAPDATA_FILE  "ZDCChMapping.dat"
23 #define TDCDATA_FILE  "ZDCTDCCalib.dat"
24 #define TDCHISTO_FILE "ZDCTDCHisto.root"
25
26 #include <stdio.h>
27 #include <stdlib.h>
28 #include <Riostream.h>
29 #include <signal.h>
30
31 // DATE
32 #include <event.h>
33 #include <monitor.h>
34 #include <daqDA.h>
35
36 //ROOT
37 #include <TROOT.h>
38 #include <TPluginManager.h>
39 #include <TH1F.h>
40 #include <TH2F.h>
41 #include <TProfile.h>
42 #include <TF1.h>
43 #include <TFile.h>
44 #include <TFitter.h>
45 #include "TMinuitMinimizer.h"
46
47 //AliRoot
48 #include <AliRawReaderDate.h>
49 #include <AliRawEventHeaderBase.h>
50 #include <AliZDCRawStream.h>
51
52 int main(int argc, char **argv) {
53
54   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
55                                         "*",
56                                         "TStreamerInfo",
57                                         "RIO",
58                                         "TStreamerInfo()"); 
59
60   TMinuitMinimizer m; 
61   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
62       "Minuit", "TMinuitMinimizer(const char *)");
63   TVirtualFitter::SetDefaultFitter("Minuit");
64
65   //const Char_t* tableSOD[]  = {"ALL", "no", "SOD", "all", NULL, NULL};
66   //monitorDeclareTable(const_cast<char**>(tableSOD));
67   
68   char *monitor_table[] = { "ALL", "no", "PHY", "yes", "SOD", "all", NULL };
69   int err = monitorDeclareTable(monitor_table);
70   if(err){
71     printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(err));
72     return -1;
73   } 
74   
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;
80   
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.;
84   
85   Int_t iMod=-1;
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;
89   }
90   
91   Int_t ich=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;
96   }
97   
98   Int_t iScCh=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;
103   }
104       
105   Int_t itdcCh=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;
110   }
111  
112   TH1F * hTDC[6];
113   char ntdchist[20];
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.);
121   }
122   
123   /* log start of process */
124   printf("\n ZDC MAPPING program started\n");  
125   signal(SIGSEGV, SIG_DFL);
126
127   /* check that we got some arguments = list of files */
128   if (argc<2) {
129     printf("Wrong number of arguments\n");
130     return -1;
131   }
132   
133   FILE *mapFile4Shuttle;
134
135   /* read the data files */
136   int n;
137   for(n=1;n<argc;n++){
138    
139     status=monitorSetDataSource( argv[n] );
140     if (status!=0) {
141       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
142       return -1;
143     }
144     
145     /* declare monitoring program */
146     status=monitorDeclareMp( __FILE__ );
147     if (status!=0) {
148       printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
149       return -1;
150     }
151     monitorSetNowait();
152     monitorSetNoWaitNetworkTimeout(1000);
153
154
155     int iev = 0;
156     //Bool_t sodRead = kFALSE;
157     //while(!sodRead){
158     
159     /* loop on events (infinite) */
160     for(;;) {
161
162       struct eventHeaderStruct *event;
163       eventTypeType eventT;
164  
165       /* check shutdown condition */
166       if (daqDA_checkShutdown()) {break;}
167
168       /* get next event */
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 */
173       }
174       if(status!=0) {
175         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
176         return -1;
177       }
178
179       /* retry if got no event */
180       if(event==NULL) continue;
181       
182       // Initalize raw-data reading and decoding
183       AliRawReader *reader = new AliRawReaderDate((void*)event);
184       reader->Reset();
185       reader->Select("ZDC");
186       // --- Reading event header
187       //UInt_t evtype = reader->GetType();
188       //printf("\t ZDCMAPPINGda -> run # %d\n",reader->GetRunNumber());
189       //
190       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
191         
192         
193       /* use event - here, just write event id to result file */
194       eventT=event->eventType;
195       
196       if(eventT==START_OF_DATA){
197                 
198         rawStreamZDC->SetSODReading(kTRUE);
199         
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
206                iMod++;
207                modGeo[iMod]  = rawStreamZDC->GetADCModule();
208                modType[iMod] = rawStreamZDC->GetModType();
209                modNCh[iMod]  = rawStreamZDC->GetADCNChannels();
210             }
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);
218                 ich++;
219               }
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);
226                 iScCh++;
227               }
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);
232                 itdcCh++;
233               }
234             }
235           }
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]);
242           }
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]);
248           }
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]);
254           }
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]);
260           }
261           
262         //} //if (rawstream)
263         fclose(mapFile4Shuttle);
264       }// SOD event
265       else if(eventT==PHYSICS_EVENT){ 
266
267         rawStreamZDC->SetSODReading(kTRUE);
268
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]);
276         }
277         
278         while(rawStreamZDC->Next()){
279           if(eventT==PHYSICS_EVENT && rawStreamZDC->GetADCModule()==kZDCTDCGeo && rawStreamZDC->IsZDCTDCDatum()==kTRUE){
280              //
281              itdc = rawStreamZDC->GetChannel(); 
282              if((itdc>=8 && itdc<=13) || itdc==15){
283                if(itdc==iprevtdc) ihittdc++;
284                else ihittdc=0;
285                iprevtdc=itdc;
286                if(ihittdc<1 && itdc!=15) tdcData[itdc-8] = 0.025*rawStreamZDC->GetZDCTDCDatum();
287                //
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);
292                  }
293                }
294              }
295           }
296         }
297         
298         nphys++;
299         
300       }//(if PHYSICS_EVENT) 
301       else if(eventT==END_OF_RUN){
302         printf("End Of Run detected\n");
303         break;
304       }
305       
306       delete rawStreamZDC;
307       rawStreamZDC = 0x0;
308       
309       iev++; 
310
311       /* free resources */
312       free(event);
313  
314     } // event loop    
315       
316   }
317   printf(" \n # of physics events analyzed = %d\n\n", nphys);
318
319   /* Analysis of the histograms */
320   //
321   FILE *fileShuttle;
322   fileShuttle = fopen(TDCDATA_FILE,"w");
323   //
324   Float_t xUp=0., xLow=0., deltaX=0;
325   Int_t binMax=0, nBinsx=0;
326   Float_t mean[6], sigma[6];
327   TF1 *fitfun[6];
328   for(Int_t k=0; k<6; k++){
329     if(hTDC[k]->GetEntries()!=0){
330        binMax = hTDC[k]->GetMaximumBin();
331        if(binMax<=1){
332          printf("\n WARNING! Something wrong with det %d histo \n\n", k);
333          continue;
334        }
335        // 
336        xUp = (hTDC[k]->GetXaxis())->GetXmax();
337        xLow = (hTDC[k]->GetXaxis())->GetXmin();
338        deltaX = xUp-xLow;
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]);
346        //
347        fprintf(fileShuttle,"\t%f\t%f\n",mean[k], sigma[k]);
348      }
349      else printf("  WARNING! TDC histo %d has no entries and can't be processed!\n",k);
350   }
351   //                                                   
352   fclose(fileShuttle);
353   
354   TFile *histofile = new TFile(TDCHISTO_FILE,"RECREATE");
355   histofile->cd();
356   for(int k=0; k<6; k++) hTDC[k]->Write();                                                     
357   histofile->Close();
358   //
359   for(Int_t j=0; j<6; j++) delete hTDC[j];
360   
361   /* store the result files on FES */
362   // [1] File with mapping
363   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
364   if(status){
365     printf("Failed to export file : %d\n",status);
366     return -1;
367   }
368   // [2] File with TDC data
369   status = daqDA_FES_storeFile(TDCDATA_FILE, "TDCDATA");
370   if(status){
371     printf("Failed to export TDC data file to DAQ FES\n");
372     return -1;
373   }
374   // [3] File with TDC histos
375   status = daqDA_FES_storeFile(TDCHISTO_FILE, "TDCHISTOS");
376   if(status){
377     printf("Failed to export TDC histos file to DAQ FES\n");
378     return -1;
379   }
380
381   return status;
382 }