]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCEMDda.cxx
Error corrected
[u/mrichter/AliRoot.git] / ZDC / ZDCEMDda.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_EMD runs
13
14 Contact: Chiara.Oppedisano@to.infn.it
15 Link: 
16 Run Type: STANDALONE_EMD_RUN
17 DA Type: LDC
18 Number of events needed: at least ~5*10^3
19 Input Files: ZDCPedestal.dat
20 Output Files: ZDCEMDCalib.dat, 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 TOWCALIBDATA_FILE  "ZDCTowerCalib.dat"
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   // No. of ZDC cabled ch.
75   int const kNChannels = 24;
76   int const kNScChannels = 32;
77   Int_t kFirstADCGeo=0, kLastADCGeo=3;
78       
79   Int_t ich=0;
80   Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
81   Int_t det[2*kNChannels], sec[2*kNChannels];
82   for(Int_t y=0; y<2*kNChannels; y++){
83     adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
84   }
85
86   /* log start of process */
87   printf("ZDC EMD program started\n");  
88
89   /* check that we got some arguments = list of files */
90   if (argc<2) {
91     printf("Wrong number of arguments\n");
92     return -1;
93   }
94   
95   // --- Preparing histos for EM dissociation spectra
96   //
97   TH1F* histoEMDRaw[4];
98   TH1F* histoEMDCorr[4];
99   //
100   char namhistr[50], namhistc[50];
101   for(Int_t i=0; i<4; i++) {
102      if(i==0){
103        sprintf(namhistr,"ZN%d-EMDRaw",i+1);
104        sprintf(namhistc,"ZN%d-EMDCorr",i+1);
105      }
106      else if(i==1){
107        sprintf(namhistr,"ZP%d-EMDRaw",i);
108        sprintf(namhistc,"ZP%d-EMDCorr",i);
109      }
110      else if(i==2){
111        sprintf(namhistr,"ZN%d-EMDRaw",i);
112        sprintf(namhistc,"ZN%d-EMDCorr",i);
113      }
114      else if(i==3){
115        sprintf(namhistr,"ZP%d-EMDRaw",i-1);
116        sprintf(namhistc,"ZP%d-EMDCorr",i-1);
117      }
118      histoEMDRaw[i] = new TH1F(namhistr,namhistr,100,0.,4000.);
119      histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
120   }
121
122    // --- Preparing histos for tower inter-calibration
123   //
124   TH1F* histZNCtow[4]; TH1F* histZPCtow[4];
125   TH1F* histZNAtow[4]; TH1F* histZPAtow[4];
126   //
127   char namhistznc[50], namhistzpc[50];
128   char namhistzna[50], namhistzpa[50];
129   for(Int_t i=0; i<4; i++) {
130      sprintf(namhistznc,"ZNC-tow%d",i+1);
131      sprintf(namhistzpc,"ZPC-tow%d",i+1);
132      sprintf(namhistzna,"ZNA-tow%d",i+1);
133      sprintf(namhistzpa,"ZPA-tow%d",i+1);
134      //
135      histZNCtow[i] = new TH1F(namhistznc,namhistznc,100,0.,4000.);
136      histZPCtow[i] = new TH1F(namhistzpc,namhistzpc,100,0.,4000.);
137      histZNAtow[i] = new TH1F(namhistzna,namhistzna,100,0.,4000.);
138      histZPAtow[i] = new TH1F(namhistzpa,namhistzpa,100,0.,4000.);
139   }
140   
141   /* open result file */
142   FILE *fp=NULL;
143   fp=fopen("./result.txt","a");
144   if (fp==NULL) {
145     printf("Failed to open file\n");
146     return -1;
147   }
148   
149   FILE *mapFile4Shuttle;
150
151   // *** To analyze EMD events you MUST have a pedestal data file!!!
152   // *** -> check if a pedestal run has been analyzed
153   int read = 0;
154   read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
155   if(read){
156     printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
157     return -1;
158   }
159   else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
160   
161   FILE *filePed = fopen(PEDDATA_FILE,"r");
162   if (filePed==NULL) {
163     printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
164     return -1;
165   }
166
167   // 144 = 48 in-time + 48 out-of-time + 48 correlations
168   Float_t readValues[2][3*2*kNChannels];
169   Float_t MeanPed[2*kNChannels];
170   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
171   // ***************************************************
172   //   Unless we have a narrow correlation to fit we
173   //    don't fit and store in-time vs. out-of-time
174   //    histograms -> mean pedstal subtracted!!!!!!
175   // ***************************************************
176   //
177   for(int jj=0; jj<6*kNChannels; jj++){
178     for(int ii=0; ii<2; ii++){
179        fscanf(filePed,"%f",&readValues[ii][jj]);
180     }
181     if(jj<2*kNChannels){
182       MeanPed[jj] = readValues[0][jj];
183       printf("\t MeanPed[%d] = %1.1f\n",jj, MeanPed[jj]);
184     }
185     else if(jj>2*kNChannels){
186       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
187       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
188     }
189   }
190
191   /* report progress */
192   daqDA_progressReport(10);
193
194
195   /* init some counters */
196   int nevents_physics=0;
197   int nevents_total=0;
198
199   struct eventHeaderStruct *event;
200   eventTypeType eventT;
201
202   /* read the data files */
203   int n;
204   for(n=1;n<argc;n++){
205    
206     status=monitorSetDataSource( argv[n] );
207     if (status!=0) {
208       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
209       return -1;
210     }
211
212     /* report progress */
213     /* in this example, indexed on the number of files */
214     daqDA_progressReport(10+80*n/argc);
215
216     /* read the file */
217     for(;;){
218
219       /* get next event */
220       status=monitorGetEventDynamic((void **)&event);
221       if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
222       if(status!=0) {
223         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
224         return -1;
225       }
226
227       /* retry if got no event */
228       if(event==NULL) {
229         break;
230       }
231       
232       // Initalize raw-data reading and decoding
233       AliRawReader *reader = new AliRawReaderDate((void*)event);
234       reader->Select("ZDC");
235       // --- Reading event header
236       //UInt_t evtype = reader->GetType();
237       //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
238       //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
239       //
240       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
241         
242
243       /* use event - here, just write event id to result file */
244       eventT=event->eventType;
245       
246       Int_t iScCh=0;
247       Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
248       Int_t scDet[kNScChannels], scSec[kNScChannels];
249       for(Int_t y=0; y<kNScChannels; y++){
250         scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
251       }
252       //
253       Int_t modNum=-1, modType=-1;
254       
255       if(eventT==START_OF_DATA){
256                 
257         rawStreamZDC->SetSODReading(kTRUE);
258         
259         // --------------------------------------------------------
260         // --- Writing ascii data file for the Shuttle preprocessor
261         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
262         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
263         else{
264           while((rawStreamZDC->Next())){
265             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
266                modNum = rawStreamZDC->GetADCModule();
267                modType = rawStreamZDC->GetModType();
268             }
269             if(rawStreamZDC->IsChMapping()){ 
270               if(modType==1){ // ADC mapping ----------------------
271                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
272                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
273                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
274                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
275                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
276                 //
277                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
278                   ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
279                 //
280                 //printf("  Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
281                 //  ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
282                 //
283                 ich++;
284               }
285               else if(modType==2){ //VME scaler mapping --------------------
286                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
287                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
288                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
289                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
290                 scSec[iScCh]    = rawStreamZDC->GetScTowerFromMap(iScCh);
291                 //
292                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
293                   iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
294                 //
295                 //printf("  Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
296                 //  iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
297                 //
298                 iScCh++;
299               }
300             }
301           }
302         }
303         fclose(mapFile4Shuttle);
304       }// SOD event
305     
306     if(eventT==PHYSICS_EVENT){
307       // --- Reading data header
308       reader->ReadHeader();
309       const AliRawDataHeader* header = reader->GetDataHeader();
310       if(header){
311          UChar_t message = header->GetAttributes();
312          if((message & 0x70) == 0x70){ // DEDICATED EMD RUN
313             //printf("\t STANDALONE_EMD_RUN raw data found\n");
314             continue;
315          }
316          else{
317             printf("\t NO STANDALONE_EMD_RUN raw data found\n");
318             return -1;
319          }
320       }
321       else{
322          printf("\t ATTENTION! No Raw Data Header found!!!\n");
323          return -1;
324       }
325       
326       rawStreamZDC->SetSODReading(kTRUE);
327
328       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
329       //
330       // ----- Setting ch. mapping -----
331       for(Int_t jk=0; jk<2*kNChannels; jk++){
332         rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
333         rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
334         rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
335         rawStreamZDC->SetMapDet(jk, det[jk]);
336         rawStreamZDC->SetMapTow(jk, sec[jk]);
337       }
338       //
339       Float_t ZDCRawADC[4], ZDCCorrADC[4], ZDCCorrADCSum[4];
340       for(Int_t g=0; g<4; g++){
341            ZDCCorrADCSum[g] = 0.;
342            ZDCRawADC[g] = 0.;
343       }
344       //
345       while(rawStreamZDC->Next()){
346         Int_t det = rawStreamZDC->GetSector(0);
347         Int_t quad = rawStreamZDC->GetSector(1);
348         
349         if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
350              && !(rawStreamZDC->IsOverflow()) && det!=-1 && det!=3 
351              && (rawStreamZDC->GetADCGain() == 1 && // Selecting LOW RES ch.s
352              rawStreamZDC->GetADCModule()>=kFirstADCGeo && rawStreamZDC->GetADCModule()<=kLastADCGeo)){
353           
354           // Taking LOW RES channels -> ch.+kNChannels !!!!
355           Int_t DetIndex=999, PedIndex=999;
356           // Not PMRef
357           if(quad!=5){ 
358             if(det == 1){
359               DetIndex = det-1;
360               PedIndex = quad+kNChannels; 
361             }
362             else if(det==2){
363               DetIndex = det-1;
364               PedIndex = quad+5+kNChannels;
365             }
366             else if(det == 4){
367               DetIndex = det-2;
368               PedIndex = quad+12+kNChannels;
369             }
370             else if(det == 5){
371               DetIndex = det-2;
372               PedIndex = quad+17+kNChannels;
373             }
374             // Mean pedestal subtraction 
375             Float_t Pedestal = MeanPed[PedIndex];
376             // Pedestal subtraction from correlation with out-of-time signals
377             //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
378             //
379             if(DetIndex!=999 || PedIndex!=999){ 
380               //
381               ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
382               //
383               //
384               ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
385               ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
386               //
387               /*printf("\t det %d quad %d res %d pedInd %d "
388                  "Pedestal %1.0f -> ADCCorr = %d ZDCCorrADCSum = %d\n", 
389                  det,quad,rawStreamZDC->GetADCGain(),PedIndex,Pedestal, 
390                  (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);*/
391                       
392             }
393             // Not common PM 
394             if(quad!=0){
395               Float_t corrADCval = (rawStreamZDC->GetADCValue()) - Pedestal;
396               if(det==1)      histZNCtow[quad-1]->Fill(corrADCval);
397               else if(det==2) histZPCtow[quad-1]->Fill(corrADCval);
398               else if(det==4) histZNAtow[quad-1]->Fill(corrADCval);
399               else if(det==5) histZPAtow[quad-1]->Fill(corrADCval);
400               //
401               //printf("\t det %d tow %d fill histo w. value %1.0f\n", 
402               //  det,quad,corrADCval);
403             }
404
405             if(DetIndex==999 || PedIndex==999) 
406                printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
407  
408           }//quad!=5
409         }//IsADCDataWord()
410         
411        }
412        //
413        nevents_physics++;
414        //
415        delete reader;
416        delete rawStreamZDC;
417        //
418        for(Int_t j=0; j<4; j++){
419           histoEMDRaw[j]->Fill(ZDCRawADC[j]);
420           histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
421        }
422     }//(if PHYSICS_EVENT)
423       
424     /* exit when last event received, no need to wait for TERM signal */
425     else if(eventT==END_OF_RUN) {
426       printf(" -> EOR event detected\n");
427       break;
428     }
429     
430     nevents_total++;
431
432    }
433
434    /* free resources */
435    free(event);
436   }
437     
438   /* Analysis of the histograms */
439   //
440   FILE *fileShuttle1 = fopen(ENCALIBDATA_FILE,"w");
441   //
442   Int_t BinMax[4]={0,0,0,0};
443   Float_t YMax[4]={0.,0.,0.,0.};
444   Int_t NBinsx[4]={0,0,0,0};
445   Float_t MeanFitVal[4]={0.,0.,0.,0.};
446   TF1 *fitfun[4];
447   for(Int_t k=0; k<4; k++){
448      if(histoEMDCorr[k]->GetEntries() == 0){
449        printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
450        return -1;
451      } 
452      //
453      BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
454      if(BinMax[k]<=6){
455        printf("\n WARNING! Something wrong with det %d histo -> ending DA WITHOUT writing output\n\n", k);
456        return -1;
457      }
458      // 
459      YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
460      NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
461      //printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
462      histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
463      fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
464      MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
465      //printf("\n\t Mean Value from gaussian fit = %f\n", MeanFitVal[k]);
466   }
467   //
468   Float_t CalibCoeff[6];     
469   //
470   for(Int_t j=0; j<6; j++){
471      if(j<4){
472        CalibCoeff[j] = MeanFitVal[j];
473        fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
474      }
475      // ZEM energy calib. coeff. = 1
476      else if(j==4 || j==5){
477        CalibCoeff[j] = 1.; 
478        fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
479      }
480   }
481   fclose(fileShuttle1);
482  
483   FILE *fileShuttle2 = fopen(TOWCALIBDATA_FILE,"w");
484   //Float_t meanvalznc[4], meanvalzpc[4], meanvalzna[4], meanvalzpa[4];
485   for(Int_t j=0; j<4; j++){
486      /*if(histZNCtow[j]->GetEntries() == 0){
487        printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
488        return -1;
489      } 
490      meanvalznc[j] = histZNCtow[j]->GetMean();
491      meanvalzpc[j] = histZPCtow[j]->GetMean();
492      meanvalzna[j] = histZNAtow[j]->GetMean();
493      meanvalzpa[j] = histZPAtow[j]->GetMean();*/
494      
495      // Note -> For the moment the inter-calibration coeff. are set to 1 
496      for(Int_t k=0; k<5; k++){  
497        Float_t icoeff = 1.;
498        fprintf(fileShuttle2,"\t%f",icoeff);
499        if(k==5) fprintf(fileShuttle2,"\n");
500      }
501   }
502   //
503   /*if(meanvalznc[1]!=0 && meanvalznc[2]!=0 && meanvalznc[3]!=0 && 
504      meanvalzpc[1]!=0 && meanvalzpc[2]!=0 && meanvalzpc[3]!=0 &&
505      meanvalzna[1]!=0 && meanvalzna[2]!=0 && meanvalzna[3]!=0 &&
506      meanvalzpa[1]!=0 && meanvalzpa[2]!=0 && meanvalzpa[3]!=0){
507     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
508         1.0,meanvalznc[0]/meanvalznc[1],meanvalznc[0]/meanvalznc[2],meanvalznc[0]/meanvalznc[3]);
509     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
510         1.0,meanvalzpc[0]/meanvalzpc[1],meanvalzpc[0]/meanvalzpc[2],meanvalzpc[0]/meanvalzpc[3]);
511     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
512         1.0,meanvalzna[0]/meanvalzna[1],meanvalzpc[0]/meanvalzna[2],meanvalzpc[0]/meanvalzna[3]);
513     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
514         1.0,meanvalzpa[0]/meanvalzpa[1],meanvalzpc[0]/meanvalzpa[2],meanvalzpc[0]/meanvalzpa[3]);
515   }
516   else{
517     printf("\n Tower intercalib. coeff. CAN'T be calculated (some mean values are ZERO)!!!\n\n");
518     return -1;
519   }*/
520   fclose(fileShuttle2);
521   
522   for(Int_t ij=0; ij<4; ij++){
523     delete histoEMDRaw[ij];
524     delete histoEMDCorr[ij];
525   }
526   
527   //delete minuitFit;
528   TVirtualFitter::SetFitter(0);
529
530   /* write report */
531   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
532
533   /* close result file */
534   fclose(fp);
535   
536   /* report progress */
537   daqDA_progressReport(90);
538
539   /* store the result file on FES */
540   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
541   if(status){
542     printf("Failed to export file : %d\n",status);
543     return -1;
544   }
545   //
546   status = daqDA_FES_storeFile(ENCALIBDATA_FILE, "EMDENERGYCALIB");
547   if(status){
548     printf("Failed to export file : %d\n",status);
549     return -1;
550   }
551   //
552   status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, "EMDTOWERCALIB");
553   if(status){
554     printf("Failed to export file : %d\n",status);
555     return -1;
556   }
557
558   /* report progress */
559   daqDA_progressReport(100);
560
561   return status;
562 }