]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCMAPPINGda.cxx
Corrected DA
[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   int status = 0;
68   int const kNModules = 10;
69   int const kNChannels = 24;
70   int const kNScChannels = 32;  
71   int const kZDCTDCGeo=4;
72   
73   int itdc=0, iprevtdc=-1, ihittdc=0;
74   float tdcData[6], tdcL0=-999.;        
75   for(int ij=0; ij<6; ij++) tdcData[ij]=-999.;
76   
77   Int_t iMod=-1;
78   Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
79   for(Int_t kl=0; kl<kNModules; kl++){
80      modGeo[kl]=modType[kl]=modNCh[kl]=0;
81   }
82   
83   Int_t ich=0;
84   Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
85   Int_t det[2*kNChannels], sec[2*kNChannels];
86   for(Int_t y=0; y<2*kNChannels; y++){
87     adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
88   }
89   
90   Int_t iScCh=0;
91   Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
92   Int_t scDet[kNScChannels], scSec[kNScChannels];
93   for(Int_t y=0; y<kNScChannels; y++){
94     scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
95   }
96       
97   Int_t itdcCh=0;
98   Int_t tdcMod[kNScChannels], tdcCh[kNScChannels], tdcSigCode[kNScChannels];
99   Int_t tdcDet[kNScChannels], tdcSec[kNScChannels];
100   for(Int_t y=0; y<kNScChannels; y++){
101     tdcMod[y]=tdcCh[y]=tdcSigCode[y]=tdcDet[y]=tdcSec[y]=-1;
102   }
103  
104   TH1F * hTDC[6];
105   char ntdchist[20];
106   for(Int_t it=0; it<6; it++){
107     if(it==0)      hTDC[it] = new TH1F("TDCZNC", "TDC ZNC", 200, 150., 350.);
108     else if(it==1) hTDC[it] = new TH1F("TDCZNA", "TDC ZNA", 200, 150., 350.);
109     else if(it==2) hTDC[it] = new TH1F("TDCZPC", "TDC ZPC", 200, 150., 350.);
110     else if(it==3) hTDC[it] = new TH1F("TDCZPA", "TDC ZPA", 200, 150., 350.);
111     else if(it==4) hTDC[it] = new TH1F("TDCZEM1","TDC ZEM1",200, 150., 350.);
112     else if(it==5) hTDC[it] = new TH1F("TDCZEM2","TDC ZEM2",200, 150., 350.);
113   }
114   
115   /* log start of process */
116   printf("\n ZDC MAPPING program started\n");  
117
118   /* check that we got some arguments = list of files */
119   if (argc<2) {
120     printf("Wrong number of arguments\n");
121     return -1;
122   }
123   
124   FILE *mapFile4Shuttle;
125
126   /* read the data files */
127   int n;
128   for(n=1;n<argc;n++){
129    
130     status=monitorSetDataSource( argv[n] );
131     if (status!=0) {
132       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
133       return -1;
134     }
135     
136     /* declare monitoring program */
137     status=monitorDeclareMp( __FILE__ );
138     if (status!=0) {
139       printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
140       return -1;
141     }
142     monitorSetNowait();
143     monitorSetNoWaitNetworkTimeout(1000);
144
145     struct eventHeaderStruct *event;
146     eventTypeType eventT;
147
148     Int_t iev = 0;
149     Bool_t sodRead = kFALSE;
150     while(!sodRead){
151  
152       /* check shutdown condition */
153       if (daqDA_checkShutdown()) {break;}
154
155       /* get next event */
156       status=monitorGetEventDynamic((void **)&event);
157       if(status==MON_ERR_EOF){
158         printf ("End of File detected\n");
159         break; /* end of monitoring file has been reached */
160       }
161       if(status!=0) {
162         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
163         return -1;
164       }
165
166       /* retry if got no event */
167       if(event==NULL) continue;
168       
169       // Initalize raw-data reading and decoding
170       AliRawReader *reader = new AliRawReaderDate((void*)event);
171       reader->Select("ZDC");
172       // --- Reading event header
173       //UInt_t evtype = reader->GetType();
174       //printf("\t ZDCMAPPINGda -> run # %d\n",reader->GetRunNumber());
175       //
176       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
177         
178         
179       /* use event - here, just write event id to result file */
180       eventT=event->eventType;
181       
182       if(eventT==START_OF_DATA){
183                 
184         rawStreamZDC->SetSODReading(kTRUE);
185         
186         // --------------------------------------------------------
187         // --- Writing ascii data file for the Shuttle preprocessor
188         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
189         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
190         else{
191           while((rawStreamZDC->Next())){
192             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
193                iMod++;
194                modGeo[iMod]  = rawStreamZDC->GetADCModule();
195                modType[iMod] = rawStreamZDC->GetModType();
196                modNCh[iMod]  = rawStreamZDC->GetADCNChannels();
197             }
198             if(rawStreamZDC->IsChMapping()){ 
199               if(modType[iMod]==1){ // ADC mapping ----------------------
200                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
201                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
202                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
203                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
204                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
205                 ich++;
206               }
207               else if(modType[iMod]==2){ // VME scaler mapping -------------
208                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
209                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
210                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
211                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
212                 scSec[iScCh]     = rawStreamZDC->GetScTowerFromMap(iScCh);
213                 iScCh++;
214               }
215               else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
216                 tdcMod[itdcCh]     = rawStreamZDC->GetTDCModFromMap(itdcCh);
217                 tdcCh[itdcCh]      = rawStreamZDC->GetTDCChFromMap(itdcCh);
218                 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
219                 itdcCh++;
220               }
221             }
222           }
223           // Writing data on output FXS file
224           for(Int_t is=0; is<2*kNChannels; is++){
225              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
226                is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
227              //printf("  Mapping DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
228              //  is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
229           }
230           for(Int_t is=0; is<kNScChannels; is++){
231              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
232                is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
233              //if(scMod[is]!=-1) printf("  Mapping DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
234              //  is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
235           }
236           for(Int_t is=0; is<kNScChannels; is++){
237              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
238                is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
239              //if(tdcMod[is]!=-1) printf("  Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
240              //  is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
241           }
242           for(Int_t is=0; is<kNModules; is++){
243              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
244              modGeo[is],modType[is],modNCh[is]);
245              //printf("  Mapping DA -> Module mapping: geo %d type %d #ch %d\n",
246              //  modGeo[is],modType[is],modNCh[is]);
247           }
248           
249         }
250         fclose(mapFile4Shuttle);
251       }// SOD event
252       else if(eventT==PHYSICS_EVENT){ 
253
254         rawStreamZDC->SetSODReading(kTRUE);
255
256         // ----- Setting ch. mapping -----
257         for(Int_t jk=0; jk<2*kNChannels; jk++){
258           rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
259           rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
260           rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
261           rawStreamZDC->SetMapDet(jk, det[jk]);
262           rawStreamZDC->SetMapTow(jk, sec[jk]);
263         }
264         
265         while(rawStreamZDC->Next()){
266           if(rawStreamZDC->GetADCModule()==kZDCTDCGeo && rawStreamZDC->IsZDCTDCDatum()==kTRUE){
267              //
268              itdc = rawStreamZDC->GetChannel(); 
269              if((itdc>=8 && itdc<=13) || itdc==15){
270                if(itdc==iprevtdc) ihittdc++;
271                else ihittdc=0;
272                iprevtdc=itdc;
273                if(ihittdc<1) tdcData[itdc-8] = 0.025*rawStreamZDC->GetZDCTDCDatum();
274                //
275                if(itdc==15 && ihittdc<1){
276                  tdcL0 = 0.025*rawStreamZDC->GetZDCTDCDatum();
277                  for(int ic=0; ic<6; ic++) if(tdcData[ic]!=-999. && tdcL0!=-999.) hTDC[ic]->Fill(tdcData[ic]-tdcL0);
278                }
279              }
280           }
281         }
282       }//(if PHYSICS_EVENT) 
283       
284       iev++; 
285
286       /* free resources */
287       free(event);
288     }    
289       
290   }
291
292   /* Analysis of the histograms */
293   //
294   FILE *fileShuttle;
295   fileShuttle = fopen(TDCDATA_FILE,"w");
296   //
297   Float_t xUp=0., xLow=0., deltaX=0;
298   Int_t binMax=0, nBinsx=0;
299   Float_t mean[6], sigma[6];
300   TF1 *fitfun[6];
301   for(Int_t k=0; k<6; k++){
302     if(hTDC[k]->GetEntries()!=0){
303        binMax = hTDC[k]->GetMaximumBin();
304        if(binMax<=1){
305          printf("\n WARNING! Something wrong with det %d histo \n\n", k);
306          continue;
307        }
308        // 
309        xUp = (hTDC[k]->GetXaxis())->GetXmax();
310        xLow = (hTDC[k]->GetXaxis())->GetXmin();
311        deltaX = xUp-xLow;
312        nBinsx = (hTDC[k]->GetXaxis())->GetNbins();
313        //printf(" xMax = %f\n", xLow+binMax*deltaX/nBinsx);
314        hTDC[k]->Fit("gaus","Q","",xLow+binMax*deltaX/nBinsx*0.6,xLow+binMax*deltaX/nBinsx*1.24);
315        fitfun[k] = hTDC[k]->GetFunction("gaus");
316        mean[k] = (Float_t) (fitfun[k]->GetParameter(1));
317        sigma[k] = (Float_t) (fitfun[k]->GetParameter(2));
318        //printf("\t Mean value from fit = %1.2f\n", mean[k]);
319        //
320        fprintf(fileShuttle,"\t%f\t%f\n",mean[k], sigma[k]);
321      }
322      else printf("  WARNING! TDC histo %d has no entries and can't be processed!\n",k);
323   }
324   //                                                   
325   fclose(fileShuttle);
326   
327   TFile *histofile = new TFile(TDCHISTO_FILE,"RECREATE");
328   histofile->cd();
329   for(int k=0; k<6; k++) hTDC[k]->Write();                                                     
330   histofile->Close();
331   //
332   for(Int_t j=0; j<6; j++) delete hTDC[j];
333   
334   /* store the result files on FES */
335   // [1] File with mapping
336   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
337   if(status){
338     printf("Failed to export file : %d\n",status);
339     return -1;
340   }
341   // [2] File with TDC data
342   status = daqDA_FES_storeFile(TDCDATA_FILE, "TDCDATA");
343   if(status){
344     printf("Failed to export pedestal data file to DAQ FES\n");
345     return -1;
346   }
347   // [3] File with TDC histos
348   status = daqDA_FES_storeFile(TDCHISTO_FILE, "TDCHISTOS");
349   if(status){
350     printf("Failed to export pedestal histos file to DAQ FES\n");
351     return -1;
352   }
353
354   return status;
355 }