]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCMAPPINGda.cxx
I need to change work
[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
30 // DATE
31 #include <event.h>
32 #include <monitor.h>
33 #include <daqDA.h>
34
35 //ROOT
36 #include <TROOT.h>
37 #include <TPluginManager.h>
38 #include <TH1F.h>
39 #include <TH2F.h>
40 #include <TProfile.h>
41 #include <TF1.h>
42 #include <TFile.h>
43 #include <TFitter.h>
44 #include "TMinuitMinimizer.h"
45
46 //AliRoot
47 #include <AliRawReaderDate.h>
48 #include <AliRawEventHeaderBase.h>
49 #include <AliZDCRawStream.h>
50
51 int main(int argc, char **argv) {
52
53   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
54                                         "*",
55                                         "TStreamerInfo",
56                                         "RIO",
57                                         "TStreamerInfo()"); 
58
59   TMinuitMinimizer m; 
60   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
61       "Minuit", "TMinuitMinimizer(const char *)");
62   TVirtualFitter::SetDefaultFitter("Minuit");
63
64   //const Char_t* tableSOD[]  = {"ALL", "no", "SOD", "all", NULL, NULL};
65   //monitorDeclareTable(const_cast<char**>(tableSOD));
66   
67   char *monitor_table[] = { "ALL", "no", "PHY", "yes", "SOD", "all", NULL };
68   int err = monitorDeclareTable(monitor_table);
69   if(err){
70     printf("monitorDeclareTable() failed: %s\n", monitorDecodeError(err));
71     return -1;
72   } 
73   
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;
79   
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.;
83   
84   Int_t iMod=-1;
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;
88   }
89   
90   Int_t ich=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;
95   }
96   
97   Int_t iScCh=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;
102   }
103       
104   Int_t itdcCh=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;
109   }
110  
111   TH1F * hTDC[6];
112   char ntdchist[20];
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.);
120   }
121   
122   /* log start of process */
123   printf("\n ZDC MAPPING program started\n");  
124   signal(SIGSEGV, SIG_DFL);
125
126   /* check that we got some arguments = list of files */
127   if (argc<2) {
128     printf("Wrong number of arguments\n");
129     return -1;
130   }
131   
132   FILE *mapFile4Shuttle;
133
134   /* read the data files */
135   int n;
136   for(n=1;n<argc;n++){
137    
138     status=monitorSetDataSource( argv[n] );
139     if (status!=0) {
140       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
141       return -1;
142     }
143     
144     /* declare monitoring program */
145     status=monitorDeclareMp( __FILE__ );
146     if (status!=0) {
147       printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
148       return -1;
149     }
150     monitorSetNowait();
151     monitorSetNoWaitNetworkTimeout(1000);
152
153
154     int iev = 0;
155     //Bool_t sodRead = kFALSE;
156     //while(!sodRead){
157     
158     /* loop on events (infinite) */
159     for(;;) {
160
161       struct eventHeaderStruct *event;
162       eventTypeType eventT;
163  
164       /* check shutdown condition */
165       if (daqDA_checkShutdown()) {break;}
166
167       /* get next event */
168       status=monitorGetEventDynamic((void **)&event);
169       if(status==MON_ERR_EOF){
170         printf ("End of File detected\n");
171         break; /* end of monitoring file has been reached */
172       }
173       if(status!=0) {
174         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
175         return -1;
176       }
177
178       /* retry if got no event */
179       if(event==NULL) continue;
180       
181       // Initalize raw-data reading and decoding
182       AliRawReader *reader = new AliRawReaderDate((void*)event);
183       reader->Reset();
184       reader->Select("ZDC");
185       // --- Reading event header
186       //UInt_t evtype = reader->GetType();
187       //printf("\t ZDCMAPPINGda -> run # %d\n",reader->GetRunNumber());
188       //
189       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
190         
191         
192       /* use event - here, just write event id to result file */
193       eventT=event->eventType;
194       
195       if(eventT==START_OF_DATA){
196                 
197         rawStreamZDC->SetSODReading(kTRUE);
198         
199         // --------------------------------------------------------
200         // --- Writing ascii data file for the Shuttle preprocessor
201         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
202         //if(rawStreamZDC->Next()){
203           while((rawStreamZDC->Next())){
204             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
205                iMod++;
206                modGeo[iMod]  = rawStreamZDC->GetADCModule();
207                modType[iMod] = rawStreamZDC->GetModType();
208                modNCh[iMod]  = rawStreamZDC->GetADCNChannels();
209             }
210             if(rawStreamZDC->IsChMapping()){ 
211               if(modType[iMod]==1){ // ADC mapping ----------------------
212                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
213                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
214                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
215                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
216                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
217                 ich++;
218               }
219               else if(modType[iMod]==2){ // VME scaler mapping -------------
220                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
221                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
222                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
223                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
224                 scSec[iScCh]     = rawStreamZDC->GetScTowerFromMap(iScCh);
225                 iScCh++;
226               }
227               else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
228                 tdcMod[itdcCh]     = rawStreamZDC->GetTDCModFromMap(itdcCh);
229                 tdcCh[itdcCh]      = rawStreamZDC->GetTDCChFromMap(itdcCh);
230                 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
231                 itdcCh++;
232               }
233             }
234           }
235           // Writing data on output FXS file
236           for(Int_t is=0; is<2*kNChannels; is++){
237              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
238                is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
239              //printf("  Mapping DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
240              //  is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
241           }
242           for(Int_t is=0; is<kNScChannels; is++){
243              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
244                is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
245              //if(scMod[is]!=-1) printf("  Mapping DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
246              //  is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
247           }
248           for(Int_t is=0; is<kNScChannels; is++){
249              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
250                is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
251              //if(tdcMod[is]!=-1) printf("  Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
252              //  is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
253           }
254           for(Int_t is=0; is<kNModules; is++){
255              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
256              modGeo[is],modType[is],modNCh[is]);
257              //printf("  Mapping DA -> Module mapping: geo %d type %d #ch %d\n",
258              //  modGeo[is],modType[is],modNCh[is]);
259           }
260           
261         //} //if (rawstream)
262         fclose(mapFile4Shuttle);
263       }// SOD event
264       else if(eventT==PHYSICS_EVENT){ 
265
266         rawStreamZDC->SetSODReading(kTRUE);
267
268         // ----- Setting ch. mapping -----
269         for(Int_t jk=0; jk<2*kNChannels; jk++){
270           rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
271           rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
272           rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
273           rawStreamZDC->SetMapDet(jk, det[jk]);
274           rawStreamZDC->SetMapTow(jk, sec[jk]);
275         }
276         
277         while(rawStreamZDC->Next()){
278           if(rawStreamZDC->GetADCModule()==kZDCTDCGeo && rawStreamZDC->IsZDCTDCDatum()==kTRUE){
279              //
280              itdc = rawStreamZDC->GetChannel(); 
281              if((itdc>=8 && itdc<=13) || itdc==15){
282                if(itdc==iprevtdc) ihittdc++;
283                else ihittdc=0;
284                iprevtdc=itdc;
285                if(ihittdc<1 && itdc!=15) tdcData[itdc-8] = 0.025*rawStreamZDC->GetZDCTDCDatum();
286                //
287                if(itdc==15 && ihittdc<1){
288                  tdcL0 = 0.025*rawStreamZDC->GetZDCTDCDatum();
289                  for(int ic=0; ic<6; ic++){
290                    if(tdcData[ic]!=-999. && tdcL0!=-999.) hTDC[ic]->Fill(tdcData[ic]-tdcL0);
291                  }
292                }
293              }
294           }
295         }
296         
297         nphys++;
298         
299       }//(if PHYSICS_EVENT) 
300       else if(eventT==END_OF_RUN){
301         printf("End Of Run detected\n");
302         break;
303       }
304       
305       delete rawStreamZDC;
306       rawStreamZDC = 0x0;
307       
308       iev++; 
309
310       /* free resources */
311       free(event);
312  
313     } // event loop    
314       
315   }
316   printf(" \n # of physics events analyzed = %d\n\n", nphys);
317
318   /* Analysis of the histograms */
319   //
320   FILE *fileShuttle;
321   fileShuttle = fopen(TDCDATA_FILE,"w");
322   //
323   Float_t xUp=0., xLow=0., deltaX=0;
324   Int_t binMax=0, nBinsx=0;
325   Float_t mean[6], sigma[6];
326   TF1 *fitfun[6];
327   for(Int_t k=0; k<6; k++){
328     if(hTDC[k]->GetEntries()!=0){
329        binMax = hTDC[k]->GetMaximumBin();
330        if(binMax<=1){
331          printf("\n WARNING! Something wrong with det %d histo \n\n", k);
332          continue;
333        }
334        // 
335        xUp = (hTDC[k]->GetXaxis())->GetXmax();
336        xLow = (hTDC[k]->GetXaxis())->GetXmin();
337        deltaX = xUp-xLow;
338        nBinsx = (hTDC[k]->GetXaxis())->GetNbins();
339        //printf(" xMax = %f\n", xLow+binMax*deltaX/nBinsx);
340        hTDC[k]->Fit("gaus","Q","",xLow+binMax*deltaX/nBinsx*0.6,xLow+binMax*deltaX/nBinsx*1.24);
341        fitfun[k] = hTDC[k]->GetFunction("gaus");
342        mean[k] = (Float_t) (fitfun[k]->GetParameter(1));
343        sigma[k] = (Float_t) (fitfun[k]->GetParameter(2));
344        //printf("\t Mean value from fit = %1.2f\n", mean[k]);
345        //
346        fprintf(fileShuttle,"\t%f\t%f\n",mean[k], sigma[k]);
347      }
348      else printf("  WARNING! TDC histo %d has no entries and can't be processed!\n",k);
349   }
350   //                                                   
351   fclose(fileShuttle);
352   
353   TFile *histofile = new TFile(TDCHISTO_FILE,"RECREATE");
354   histofile->cd();
355   for(int k=0; k<6; k++) hTDC[k]->Write();                                                     
356   histofile->Close();
357   //
358   for(Int_t j=0; j<6; j++) delete hTDC[j];
359   
360   /* store the result files on FES */
361   // [1] File with mapping
362   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
363   if(status){
364     printf("Failed to export file : %d\n",status);
365     return -1;
366   }
367   // [2] File with TDC data
368   status = daqDA_FES_storeFile(TDCDATA_FILE, "TDCDATA");
369   if(status){
370     printf("Failed to export TDC data file to DAQ FES\n");
371     return -1;
372   }
373   // [3] File with TDC histos
374   status = daqDA_FES_storeFile(TDCHISTO_FILE, "TDCHISTOS");
375   if(status){
376     printf("Failed to export TDC histos file to DAQ FES\n");
377     return -1;
378   }
379
380   return status;
381 }