]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCCALIBMBda.cxx
UPdate to digitize correctly with RELDIS
[u/mrichter/AliRoot.git] / ZDC / ZDCCALIBMBda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5 It computes the average event size and populates local "./result.txt" file with the 
6 result.
7
8 The program reports about its processing progress.
9
10 Messages on stdout are exported to DAQ log system.
11
12 DA for ZDC standalone CALIBRATION_MB runs
13
14 Contact: Chiara.Oppedisano@to.infn.it
15 Link: 
16 Run Type: CALIBRATION_MB_RUN
17 DA Type: LDC
18 Number of events needed: at least 10^4
19 Input Files: ZDCPedestal.dat
20 Output Files: ZDCMBCalib.root, ZDCChMapping.dat
21 Trigger Types Used: Standalone Trigger
22
23 */
24 #define PEDDATA_FILE  "ZDCPedestal.dat"
25 #define MAPDATA_FILE  "ZDCChMapping.dat"
26 #define ENCALIBDATA_FILE   "ZDCEnergyCalib.dat"
27 #define MBCALIBDATA_FILE   "ZDCMBCalib.root"
28
29 #include <stdio.h>
30 #include <Riostream.h>
31 #include <Riostream.h>
32
33 // DATE
34 #include <daqDA.h>
35 #include <event.h>
36 #include <monitor.h>
37
38 //ROOT
39 #include <TROOT.h>
40 #include <TPluginManager.h>
41 #include <TH1F.h>
42 #include <TH2F.h>
43 #include <TProfile.h>
44 #include <TF1.h>
45 #include <TFile.h>
46 #include <TFitter.h>
47 #include "TMinuitMinimizer.h"
48
49 //AliRoot
50 #include <AliRawReaderDate.h>
51 #include <AliRawEventHeaderBase.h>
52 #include <AliZDCRawStream.h>
53
54
55 /* Main routine
56       Arguments: 
57       1- monitoring data source
58 */
59 int main(int argc, char **argv) {
60   
61
62   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
63                                         "*",
64                                         "TStreamerInfo",
65                                         "RIO",
66                                         "TStreamerInfo()"); 
67
68   TMinuitMinimizer m; 
69   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
70       "Minuit", "TMinuitMinimizer(const char *)");
71   TVirtualFitter::SetDefaultFitter("Minuit");
72
73   int status = 0;
74   int const kNModules = 10;
75   int const kNChannels = 24;
76   int const kNScChannels = 32;
77   Int_t kFirstADCGeo=0, kLastADCGeo=1;
78       
79   Int_t iMod=-1;
80   Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
81   for(Int_t kl=0; kl<kNModules; kl++){
82      modGeo[kl]=modType[kl]=modNCh[kl]=0;
83   }
84   
85   Int_t ich=0;
86   Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
87   Int_t det[2*kNChannels], sec[2*kNChannels];
88   for(Int_t y=0; y<2*kNChannels; y++){
89     adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
90   }
91   
92   Int_t iScCh=0;
93   Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
94   Int_t scDet[kNScChannels], scSec[kNScChannels];
95   for(Int_t y=0; y<kNScChannels; y++){
96     scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
97   }
98       
99   Int_t itdcCh=0;
100   Int_t tdcMod[kNScChannels], tdcCh[kNScChannels], tdcSigCode[kNScChannels];
101   Int_t tdcDet[kNScChannels], tdcSec[kNScChannels];
102   for(Int_t y=0; y<kNScChannels; y++){
103     tdcMod[y]=tdcCh[y]=tdcSigCode[y]=tdcDet[y]=tdcSec[y]=-1;
104   }
105
106   /* log start of process */
107   printf("ZDC EMD program started\n");  
108
109   /* check that we got some arguments = list of files */
110   if (argc<2) {
111     printf("Wrong number of arguments\n");
112     return -1;
113   }
114   
115   // --- Preparing histos for EM dissociation spectra
116   //
117   TH2F * hZDCvsZEM  = new TH2F("hZDCvsZEM","EZDCTot vs. EZEM",100,0.,8.,100,0.,800.);
118   TH2F * hZDCCvsZEM = new TH2F("hZDCCvsZEM","EZDCC vs. EZEM",100,0.,8.,100,0.,400.);
119   TH2F * hZDCAvsZEM = new TH2F("hZDCAvsZEM","EZDCA vs. EZEM",100,0.,8.,100,0.,400.);
120   
121   /* open result file */
122   FILE *fp=NULL;
123   fp=fopen("./result.txt","a");
124   if (fp==NULL) {
125     printf("Failed to open file\n");
126     return -1;
127   }
128   
129   FILE *mapFile4Shuttle;
130
131   // *** To analyze CALIBRATION_MB events a pedestal data file is needed!!!
132   // *** -> check if a pedestal run has been analyzed
133   int read = 0;
134   read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
135   if(read){
136     printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
137     return -1;
138   }
139   else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
140   
141   FILE *filePed = fopen(PEDDATA_FILE,"r");
142   if (filePed==NULL) {
143     printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
144     return -1;
145   }
146
147   // 144 = 48 in-time + 48 out-of-time + 48 correlations
148   Float_t readValues[2][6*kNChannels];
149   Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
150   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
151   // ***************************************************
152   //   Unless we have a narrow correlation to fit we
153   //    don't fit and store in-time vs. out-of-time
154   //    histograms -> mean pedstal subtracted!!!!!!
155   // ***************************************************
156   //
157   for(int jj=0; jj<6*kNChannels; jj++){
158     for(int ii=0; ii<2; ii++){
159        fscanf(filePed,"%f",&readValues[ii][jj]);
160     }
161     if(jj<kNChannels){
162       MeanPedhg[jj] = readValues[0][jj];
163       //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
164     }
165     else if(jj>=kNChannels && jj<2*kNChannels){
166       MeanPedlg[jj-kNChannels] = readValues[0][jj];
167       //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
168     }
169     else if(jj>4*kNChannels){
170       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
171       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
172     }
173   }
174
175   /* report progress */
176   daqDA_progressReport(10);
177   
178   FILE *fileEnCal = fopen(ENCALIBDATA_FILE,"r");
179   if(fileEnCal==NULL) {
180     printf("\t ERROR!!! Can't open ZDCEnergyCalib.dat file!!!\n");
181     return -1;
182   }
183   
184   Float_t calibCoeff[6];
185   for(Int_t irow=0; irow<6; irow++){
186      fscanf(fileEnCal,"%f",&calibCoeff[irow]);
187   }
188
189   /* init some counters */
190   int nevents_physics=0;
191   int nevents_total=0;
192
193   struct eventHeaderStruct *event;
194   eventTypeType eventT;
195
196   /* read the data files */
197   int n;
198   for(n=1;n<argc;n++){
199    
200     status=monitorSetDataSource( argv[n] );
201     if (status!=0) {
202       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
203       return -1;
204     }
205
206     /* report progress */
207     /* in this example, indexed on the number of files */
208     daqDA_progressReport(10+80*n/argc);
209
210     /* read the file */
211     for(;;){
212
213       /* get next event */
214       status=monitorGetEventDynamic((void **)&event);
215       if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
216       if(status!=0) {
217         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
218         return -1;
219       }
220
221       /* retry if got no event */
222       if(event==NULL) {
223         break;
224       }
225       
226       // Initalize raw-data reading and decoding
227       AliRawReader *reader = new AliRawReaderDate((void*)event);
228       reader->Select("ZDC");
229       // --- Reading event header
230       //UInt_t evtype = reader->GetType();
231       //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
232       //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
233       //
234       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
235         
236
237       /* use event - here, just write event id to result file */
238       eventT=event->eventType;
239       
240       if(eventT==START_OF_DATA){
241                         
242         rawStreamZDC->SetSODReading(kTRUE);
243         
244         // --------------------------------------------------------
245         // --- Writing ascii data file for the Shuttle preprocessor
246         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
247         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
248         else{
249           while((rawStreamZDC->Next())){
250             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
251                iMod++;
252                modGeo[iMod]  = rawStreamZDC->GetADCModule();
253                modType[iMod] = rawStreamZDC->GetModType();
254                modNCh[iMod]  = rawStreamZDC->GetADCNChannels();
255             }
256             if(rawStreamZDC->IsChMapping()){ 
257               if(modType[iMod]==1){ // ADC mapping ----------------------
258                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
259                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
260                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
261                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
262                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
263                 ich++;
264               }
265               else if(modType[iMod]==2){ //VME scaler mapping --------------------
266                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
267                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
268                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
269                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
270                 scSec[iScCh]     = rawStreamZDC->GetScTowerFromMap(iScCh);
271                 iScCh++;
272               }
273               else if(modType[iMod]==6 && modGeo[iMod]==4){ // ZDC TDC mapping --------------------
274                 tdcMod[itdcCh]     = rawStreamZDC->GetTDCModFromMap(itdcCh);
275                 tdcCh[itdcCh]      = rawStreamZDC->GetTDCChFromMap(itdcCh);
276                 tdcSigCode[itdcCh] = rawStreamZDC->GetTDCSignFromMap(itdcCh);
277                 itdcCh++;
278               }
279             }
280           }
281           // Writing data on output FXS file
282           for(Int_t is=0; is<2*kNChannels; is++){
283              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
284                is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
285              //printf("  CalibMB DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
286              //  is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
287           }
288           for(Int_t is=0; is<kNScChannels; is++){
289              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
290                is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
291              //printf("  CalibMB DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
292              //  is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
293           }
294           for(Int_t is=0; is<kNScChannels; is++){
295              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\n",
296                is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
297              //if(tdcMod[is]!=-1) printf("  Mapping DA -> %d TDC: mod %d ch %d, code %d\n",
298              //  is,tdcMod[is],tdcCh[is],tdcSigCode[is]);
299           }
300           for(Int_t is=0; is<kNModules; is++){
301              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
302              modGeo[is],modType[is],modNCh[is]);
303              //printf("  CalibMB DA -> Module mapping: geo %d type %d #ch %d\n",
304              //  modGeo[is],modType[is],modNCh[is]);
305           }
306           
307         }
308         fclose(mapFile4Shuttle);
309       }// SOD event
310     
311     if(eventT==PHYSICS_EVENT){
312       // --- Reading data header
313       reader->ReadHeader();
314       const AliRawDataHeader* header = reader->GetDataHeader();
315       if(header){
316          UChar_t message = header->GetAttributes();
317          if((message & 0xf0) == 0x60){ // DEDICATED CALIBRATION MB RUN
318             //printf("\t CALIBRATION_MB_RUN raw data found\n");
319             continue;
320          }
321          else{
322             //printf("\t NO CALIBRAION_MB_RUN raw data found\n");
323             //return -1;
324          }
325       }
326       else{
327          printf("\t ATTENTION! No Raw Data Header found!!!\n");
328          return -1;
329       }
330       
331       rawStreamZDC->SetSODReading(kTRUE);
332
333       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
334       //
335       // ----- Setting ch. mapping -----
336       for(Int_t jk=0; jk<2*kNChannels; jk++){
337         rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
338         rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
339         rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
340         rawStreamZDC->SetMapDet(jk, det[jk]);
341         rawStreamZDC->SetMapTow(jk, sec[jk]);
342       }
343       //
344       Float_t rawZNC=0., rawZPC=0., rawZNA=0., rawZPA=0., rawZEM1=0., rawZEM2=0.;
345       Float_t corrZNC=0., corrZPC=0., corrZNA=0., corrZPA=0., corrZEM1=0.,corrZEM2=0.;
346       Float_t calZNC=0., calZPC=0., calZNA=0., calZPA=0., calZEM1=0., calZEM2=0.;
347       Float_t calZDCTot=0., calZDCC=0., calZDCA=0., calZEM=0.;
348       //
349       while(rawStreamZDC->Next()){
350         Int_t detector = rawStreamZDC->GetSector(0);
351         Int_t quad = rawStreamZDC->GetSector(1);
352         Int_t index=-1;
353         
354         if( (rawStreamZDC->IsADCDataWord()) && !(rawStreamZDC->IsUnderflow())
355              && !(rawStreamZDC->IsOverflow()) && (detector!=-1)
356              && ((rawStreamZDC->GetADCGain()) == 0) && // Selecting HIGH RES ch.s
357              (rawStreamZDC->GetADCModule()>=kFirstADCGeo) && (rawStreamZDC->GetADCModule()<=kLastADCGeo)){
358
359             //printf("  IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
360             //  rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
361           
362             if(quad!=5){ // Physics signals
363               if(detector==1)      index = quad;   // *** ZNC
364               else if(detector==2) index = quad+5; // *** ZPC
365               else if(detector==3) index = quad+9; // *** ZEM
366               else if(detector==4) index = quad+12;// *** ZNA
367               else if(detector==5) index = quad+17;// *** ZPA
368             }
369             else continue;
370             //
371             if(index==-1) printf("ERROR in ZDCCALIBMBda.cxx -> det %d quad %d index %d ADC %d\n", 
372               detector, quad, index, rawStreamZDC->GetADCValue());
373             
374             // Mean pedestal subtraction 
375             Float_t Pedestal = MeanPedhg[index];
376             // Pedestal subtraction from correlation with out-of-time signals
377             //Float_t Pedestal = CorrCoeff0[index]+CorrCoeff1[index]*MeanPedOOT[index];
378             //
379             if(index!=-1 && quad!=5){ 
380               //
381               if(detector==1){
382                rawZNC  = (Float_t) rawStreamZDC->GetADCValue();
383                corrZNC = rawZNC - Pedestal;
384               }
385               else if(detector==2){
386                rawZPC  = (Float_t) rawStreamZDC->GetADCValue();
387                corrZPC = rawZPC - Pedestal;
388               }
389               else if(detector==3){
390                if(quad==1){
391                  rawZEM1  = (Float_t) rawStreamZDC->GetADCValue();
392                  corrZEM1 = (rawStreamZDC->GetADCValue()) - Pedestal;
393                }
394                else{
395                  rawZEM2  = (Float_t) rawStreamZDC->GetADCValue();
396                  corrZEM2 = (rawStreamZDC->GetADCValue()) - Pedestal;
397                }
398               }
399               else if(detector==4){
400                rawZNA  = (Float_t) rawStreamZDC->GetADCValue();
401                corrZNA = rawZNA - Pedestal;
402               }
403               else if(detector==5){
404                rawZPA  = (Float_t) rawStreamZDC->GetADCValue();
405                corrZPA = rawZPA - Pedestal;
406               }
407             }
408         }//IsADCDataWord()
409         
410        } // Next()
411        
412        calZNC = calibCoeff[0]*corrZNC;
413        calZPC = calibCoeff[1]*corrZPC;
414        calZNA = calibCoeff[2]*corrZNA;
415        calZPA = calibCoeff[3]*corrZPA;
416        calZEM1 = calibCoeff[4]*corrZEM1;
417        calZEM2 = calibCoeff[4]*corrZEM2;
418        calZDCTot = calZNC+calZPC+calZPA+calZNA;
419        calZDCC = calZNC+calZPC;
420        calZDCA = calZNA+calZPA;
421        calZEM = calZEM1+calZEM2;
422        //
423        hZDCvsZEM ->Fill(calZDCTot/1000, calZEM/1000);
424        hZDCCvsZEM->Fill(calZDCC/1000, calZEM/1000);
425        hZDCAvsZEM->Fill(calZDCA/1000, calZEM/1000);
426
427        nevents_physics++;
428     }//(if PHYSICS_EVENT)
429     
430     nevents_total++;
431
432    }
433
434    /* free resources */
435    free(event);
436   }
437     
438   //
439   TFile *histofile = new TFile(MBCALIBDATA_FILE,"RECREATE");
440   histofile->cd();
441   hZDCvsZEM ->Write();
442   hZDCCvsZEM->Write();
443   hZDCAvsZEM->Write();
444   histofile->Close();
445   
446   if(hZDCvsZEM)  delete hZDCvsZEM ;
447   if(hZDCCvsZEM) delete hZDCCvsZEM;
448   if(hZDCAvsZEM) delete hZDCAvsZEM;
449   
450   /* write report */
451   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
452
453   /* close result file */
454   fclose(fp);
455   
456   /* report progress */
457   daqDA_progressReport(90);
458
459   /* store the result file on FES */
460   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
461   if(status){
462     printf("Failed to export file : %d\n",status);
463     return -1;
464   }
465   //
466   status = daqDA_FES_storeFile(MBCALIBDATA_FILE, "MBCALIB");
467   if(status){
468     printf("Failed to export file : %d\n",status);
469     return -1;
470   }
471  
472   /* report progress */
473   daqDA_progressReport(100);
474
475   return status;
476 }