]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCPEDESTALda.cxx
Added missing files
[u/mrichter/AliRoot.git] / ZDC / ZDCPEDESTALda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5
6 The program reports about its processing progress.
7
8 Messages on stdout are exported to DAQ log system.
9
10 DA for ZDC standalone pedestal runs
11
12 Contact: Chiara.Oppedisano@to.infn.it
13 Link: 
14 Run Type: STANDALONE_PEDESTAL_RUN
15 DA Type: LDC
16 Number of events needed: no constraint (tipically ~10^3)
17 Input Files:  
18 Output Files: ZDCPedestal.dat, ZDCChMapping.dat
19 Trigger Types Used: Standalone Trigger
20
21 */
22 #define PEDDATA_FILE  "ZDCPedestal.dat"
23 #define MAPDATA_FILE  "ZDCChMapping.dat"
24
25 #include <stdio.h>
26 #include <stdlib.h>
27 #include <Riostream.h>
28
29 // DATE
30 #include <event.h>
31 #include <monitor.h>
32 #include <daqDA.h>
33
34 //ROOT
35 #include <TRandom.h>
36 #include <TH1F.h>
37 #include <TH2F.h>
38 #include <TProfile.h>
39 #include <TF1.h>
40 #include <TFile.h>
41 #include <TFitter.h>
42
43 //AliRoot
44 #include <AliRawReaderDate.h>
45 #include <AliRawEventHeaderBase.h>
46 #include <AliZDCRawStream.h>
47
48
49 /* Main routine
50       Arguments: list of DATE raw data files
51 */
52 int main(int argc, char **argv) {
53   
54 //  TVirtualFitter::SetDefaultFitter("Minuit2");
55
56   int status = 0;
57   int const kNChannels = 24;
58
59   /* log start of process */
60   printf("\n ZDC PEDESTAL program started\n");  
61
62   /* check that we got some arguments = list of files */
63   if (argc<2) {
64     printf("Wrong number of arguments\n");
65     return -1;
66   }
67
68   // --- Histograms for ADC pedestals 
69   //     [22 signal channels + 2 reference PTMs]  x 2 gain chains
70   //
71   TH1F::AddDirectory(0);
72   //
73   TH1F *hPedhg[kNChannels], *hPedOutOfTimehg[kNChannels];
74   TH2F *hPedCorrhg[kNChannels];
75   TH1F *hPedlg[kNChannels], *hPedOutOfTimelg[kNChannels];
76   TH2F *hPedCorrlg[kNChannels];
77   //
78   char namhist1hg[50], namhist2hg[50], namhist3hg[50];
79   char namhist1lg[50], namhist2lg[50], namhist3lg[50];
80   for(Int_t j=0; j<kNChannels; j++){
81      if(j<=4){ // ZNC
82        sprintf(namhist1hg,"PedZNChg_%d",j);
83        sprintf(namhist2hg,"PedZNChgOutOfTime_%d",j);
84        sprintf(namhist3hg,"PedCorrZNChg_%d",j);
85        //
86        sprintf(namhist1lg,"PedZNClg_%d",j);
87        sprintf(namhist2lg,"PedZNClgOutOfTime_%d",j);
88        sprintf(namhist3lg,"PedCorrZNClg_%d",j);
89      }
90      else if(j>=5 && j<=9){ // ZPC
91        sprintf(namhist1hg,"PedZPChg_%d",j-5);
92        sprintf(namhist2hg,"PedZPChgOutOfTime_%d",j-5);
93        sprintf(namhist3hg,"PedCorrZPChg_%d",j-5);
94        //
95        sprintf(namhist1lg,"PedZPClg_%d",j-5);
96        sprintf(namhist2lg,"PedZPClgOutOfTime_%d",j-5);
97        sprintf(namhist3lg,"PedCorrZPClg_%d",j-5);       
98      }
99      else if(j==10 || j==11){ // ZEM
100        sprintf(namhist1hg,"PedZEMhg_%d",j-10);
101        sprintf(namhist2hg,"PedZEMhgOutOfTime_%d",j-10);
102        sprintf(namhist3hg,"PedCorrZEMhg_%d",j-10);
103        //
104        sprintf(namhist1lg,"PedZEMlg_%d",j-10);
105        sprintf(namhist2lg,"PedZEMlgOutOfTime_%d",j-10);
106        sprintf(namhist3lg,"PedCorrZEMlg_%d",j-10);
107      }
108      else if(j>=12 && j<=16){ // ZNA
109        sprintf(namhist1hg,"PedZNAhg_%d",j-12);
110        sprintf(namhist2hg,"PedZNAhgOutOfTime_%d",j-12);
111        sprintf(namhist3hg,"PedCorrZNAhg_%d",j-12);
112        //
113        sprintf(namhist1lg,"PedZNAlg_%d",j-12);
114        sprintf(namhist2lg,"PedZNAlgOutOfTime_%d",j-12);
115        sprintf(namhist3lg,"PedCorrZNAlg_%d",j-12);
116      }
117      else if(j>=17 && j<=21){ // ZPA
118        sprintf(namhist1hg,"PedZPAhg_%d",j-17);
119        sprintf(namhist2hg,"PedZPAhgOutOfTime_%d",j-17);
120        sprintf(namhist3hg,"PedCorrZPAhg_%d",j-17);
121        //
122        sprintf(namhist1lg,"PedZPAlg_%d",j-17);
123        sprintf(namhist2lg,"PedZPAlgOutOfTime_%d",j-17);
124        sprintf(namhist3lg,"PedCorrZPAlg_%d",j-17);
125      }
126      else if(j==22 || j==24){ //Reference PMs
127        sprintf(namhist1hg,"PedRefhg_%d",j-22);
128        sprintf(namhist2hg,"PedRefhgOutOfTime_%d",j-22);
129        sprintf(namhist3hg,"PedCorrRefhg_%d",j-22);
130        //
131        sprintf(namhist1lg,"PedReflg_%d",j-22);
132        sprintf(namhist2lg,"PedReflgOutOfTime_%d",j-22);
133        sprintf(namhist3lg,"PedCorrReflg_%d",j-22);
134      }
135      // --- High gain chain histos
136      hPedhg[j] = new TH1F(namhist1hg, namhist1hg, 200, 0., 200.);
137      hPedOutOfTimehg[j] = new TH1F(namhist2hg, namhist2hg, 200, 0., 200.);
138      hPedCorrhg[j] = new TH2F(namhist3hg,namhist3hg,100,0.,200.,100,0.,200.);
139      // --- Low gain chain histos
140      hPedlg[j] = new TH1F(namhist1lg, namhist1lg, 100, 0., 1000.);
141      hPedOutOfTimelg[j] = new TH1F(namhist2lg, namhist2lg, 100, 0., 1000.);
142      hPedCorrlg[j] = new TH2F(namhist3lg,namhist3lg,100,0.,1000.,100,0.,1000.);
143   }
144
145
146   /* open result file */
147   FILE *fp=NULL;
148   fp=fopen("./result.txt","w");
149   if (fp==NULL) {
150     printf("Failed to open file\n");
151     return -1;
152   }
153   
154   FILE *mapFile4Shuttle;
155
156   /* report progress */
157   daqDA_progressReport(10);
158
159
160   /* init some counters */
161   int nevents_physics=0;
162   int nevents_total=0;
163
164   /* read the data files */
165   int n;
166   for(n=1;n<argc;n++){
167    
168     status=monitorSetDataSource( argv[n] );
169     if (status!=0) {
170       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
171       return -1;
172     }
173
174     /* report progress */
175     /* in this example, indexed on the number of files */
176     daqDA_progressReport(10+80*n/argc);
177
178     /* read the file */
179     for(;;) {
180       struct eventHeaderStruct *event;
181       eventTypeType eventT;
182
183       /* get next event */
184       status=monitorGetEventDynamic((void **)&event);
185       if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
186       if(status!=0) {
187         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
188         return -1;
189       }
190
191       /* retry if got no event */
192       if(event==NULL) {
193         break;
194       }
195       
196       // Initalize raw-data reading and decoding
197       AliRawReader *reader = new AliRawReaderDate((void*)event);
198       reader->Select("ZDC");
199       // --- Reading event header
200       //UInt_t evtype = reader->GetType();
201       //printf("\n\t ZDCPEDESTALda -> ev. type %d\n",evtype);
202       //printf("\t ZDCPEDESTALda -> run # %d\n",reader->GetRunNumber());
203       //
204       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
205         
206
207       /* use event - here, just write event id to result file */
208       eventT=event->eventType;
209       
210       Int_t ich=0;
211       Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
212       Int_t det[2*kNChannels], sec[2*kNChannels];
213       
214       if(eventT==START_OF_DATA){
215                 
216         rawStreamZDC->SetSODReading(kTRUE);
217         
218         // --------------------------------------------------------
219         // --- Writing ascii data file for the Shuttle preprocessor
220         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
221         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
222         else{
223           while((rawStreamZDC->Next()) && (ich<2*kNChannels)){
224             if(rawStreamZDC->IsChMapping()){
225               adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
226               adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
227               sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
228               det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
229               sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
230               //
231               fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
232                 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
233               //
234               //printf("ZDCPEDESTALda.cxx -> %d mod %d ch %d, code %d det %d, sec %d\n",
235               //   ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
236               //
237               ich++;
238             }
239           }
240         }
241         fclose(mapFile4Shuttle);
242       }// SOD event
243       
244       if(eventT==PHYSICS_EVENT){
245         // --- Reading data header
246         reader->ReadHeader();
247         const AliRawDataHeader* header = reader->GetDataHeader();
248         if(header){
249          UChar_t message = header->GetAttributes();
250          if(message & 0x20){ // PEDESTAL RUN
251             //printf("\t STANDALONE_PEDESTAL RUN raw data found\n");
252          }
253          else{
254             printf("ZDCPEDESTALda.cxx -> NO STANDALONE_PEDESTAL RUN raw data found\n");
255             return -1;
256          }
257         }
258         else{
259            printf("\t ATTENTION! No Raw Data Header found!!!\n");
260            return -1;
261         }
262         
263         rawStreamZDC->SetSODReading(kTRUE);
264
265         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n"); 
266         //
267         // ----- Setting ch. mapping -----
268         for(Int_t jk=0; jk<2*kNChannels; jk++){
269           //printf("ZDCPEDESTALDA.cxx ->  ch.%d mod %d, ch %d, code %d det %d, sec %d\n",
270           //    jk,adcMod[jk],adcCh[jk],sigCode[jk],det[jk],sec[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         Int_t iraw=0;
279         Int_t RawADChg[kNChannels], RawADCoothg[kNChannels];
280         Int_t RawADClg[kNChannels], RawADCootlg[kNChannels];
281         for(Int_t j=0; j<kNChannels; j++){
282            RawADChg[j]=0; RawADCoothg[j]=0;
283            RawADClg[j]=0; RawADCootlg[j]=0;
284         }
285         //
286         while(rawStreamZDC->Next()){
287          Int_t index=-1;
288          Int_t detector = rawStreamZDC->GetSector(0);
289          Int_t sector = rawStreamZDC->GetSector(1);
290          
291          if(rawStreamZDC->IsADCDataWord() && (detector!=-1)){
292           if(sector!=5){ // Physics signals
293             if(detector==1) index = sector; // *** ZNC
294             else if(detector==2) index = sector+5; // *** ZPC
295             else if(detector==3) index = sector+9; // *** ZEM
296             else if(detector==4) index = sector+12; // *** ZNA
297             else if(detector==5) index = sector+17; // *** ZPA
298           }
299           else{ // Reference PMs
300             index = (detector-1)/3+22;
301           }
302           //
303           if(index==-1) printf("ERROR in ZDCPEDESTALda.cxx -> det %d quad %d res %d index %d \n", 
304             detector,sector,rawStreamZDC->GetADCGain(),index);
305           
306            //
307            if(iraw<2*kNChannels){ // --- In-time pedestals (1st 48 raw data)
308             if(rawStreamZDC->GetADCGain()==0){ 
309               hPedhg[index]->Fill(rawStreamZDC->GetADCValue()); 
310               RawADChg[index] = rawStreamZDC->GetADCValue();
311               //
312               //printf("\t filling histo hPedhg[%d]\n",index);
313             }
314             else{
315               hPedlg[index]->Fill(rawStreamZDC->GetADCValue()); 
316               RawADClg[index] = rawStreamZDC->GetADCValue();
317               //
318               //printf("\t filling histo hPedlg[%d]\n",index);
319             }
320            }
321            else{  // --- Out-of-time pedestals
322             if(rawStreamZDC->GetADCGain()==0){
323               hPedOutOfTimehg[index]->Fill(rawStreamZDC->GetADCValue());
324               RawADCoothg[index] = rawStreamZDC->GetADCValue();
325               //
326               //printf("\t filling histo hPedOutOfTimehg[%d]\n",index);
327             }
328             else{
329               hPedOutOfTimelg[index]->Fill(rawStreamZDC->GetADCValue());
330               RawADCootlg[index] = rawStreamZDC->GetADCValue();
331               //
332               //printf("\t filling histo hPedOutOfTimelg[%d]\n",index);
333             }
334            }
335             iraw++;
336           }//IsADCDataWord()
337           //
338           if(iraw == 4*kNChannels){ // Last ADC channel -> Filling correlation histos
339             for(Int_t k=0; k<kNChannels; k++){
340               hPedCorrhg[k]->Fill(RawADCoothg[k], RawADChg[k]);
341               hPedCorrlg[k]->Fill(RawADCootlg[k], RawADClg[k]);
342             }
343           }
344          }
345          //
346          nevents_physics++;
347          //
348          delete reader;
349          delete rawStreamZDC;
350
351       }//(if PHYSICS_EVENT) 
352       nevents_total++;
353
354       /* free resources */
355       free(event);
356     
357     }
358   }  
359   
360   /* Analysis of the histograms */
361   //
362   FILE *fileShuttle;
363   fileShuttle = fopen(PEDDATA_FILE,"w");
364   //
365   Float_t MeanPed[2*kNChannels], MeanPedWidth[2*kNChannels], 
366         MeanPedOOT[2*kNChannels], MeanPedWidthOOT[2*kNChannels];
367   // --- In-time pedestals
368   TF1 *ADCfunchg[kNChannels];
369   for(Int_t i=0; i<kNChannels; i++){
370      hPedhg[i]->Fit("gaus","Q");
371      ADCfunchg[i] = hPedhg[i]->GetFunction("gaus");
372      MeanPed[i] = (Double_t) ADCfunchg[i]->GetParameter(1);
373      MeanPedWidth[i] = (Double_t)  ADCfunchg[i]->GetParameter(2);
374      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i],MeanPedWidth[i]);
375      //printf("\t MeanPedhg[%d] = %f\n",i, MeanPed[i]);
376   }
377   TF1 *ADCfunclg[kNChannels];
378   for(Int_t i=0; i<kNChannels; i++){
379      hPedlg[i]->Fit("gaus","Q");
380      ADCfunclg[i] = hPedlg[i]->GetFunction("gaus");
381      MeanPed[i+kNChannels] = (Double_t)  ADCfunclg[i]->GetParameter(1);
382      MeanPedWidth[i+kNChannels] = (Double_t)  ADCfunclg[i]->GetParameter(2);
383      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i+kNChannels],MeanPedWidth[i+kNChannels]);
384      //printf("\t MeanPedlg[%d] = %f\n",i+kNChannels, MeanPed[i+kNChannels]);
385   }
386   // --- Out-of-time pedestals
387   TF1 *ADCootfunchg[kNChannels];
388   for(Int_t i=0; i<kNChannels; i++){
389      hPedOutOfTimehg[i]->Fit("gaus","Q");
390      ADCootfunchg[i] = hPedOutOfTimehg[i]->GetFunction("gaus");
391      MeanPedOOT[i] = (Double_t)  ADCootfunchg[i]->GetParameter(1);
392      MeanPedWidthOOT[i] = (Double_t)  ADCootfunchg[i]->GetParameter(2);
393      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i],MeanPedWidthOOT[i]);
394      //printf("\t MeanPedOOThg[%d] = %f\n",i, MeanPedOOT[i]);
395   }
396   TF1 *ADCootfunclg[kNChannels];
397   for(Int_t i=0; i<kNChannels; i++){
398      hPedOutOfTimelg[i]->Fit("gaus","Q");
399      ADCootfunclg[i] = hPedOutOfTimelg[i]->GetFunction("gaus");
400      MeanPedOOT[i+kNChannels] = (Double_t)  ADCootfunclg[i]->GetParameter(1);
401      MeanPedWidthOOT[i+kNChannels] = (Double_t)  ADCootfunclg[i]->GetParameter(2);
402      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i+kNChannels],MeanPedWidthOOT[i+kNChannels]);
403      //printf("\t MeanPedOOTlg[%d] = %f\n",i+kNChannels, MeanPedOOT[i+kNChannels]);
404   }
405   
406   // --- Correlations
407
408   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
409   TProfile *hPedCorrProfhg[kNChannels], *hPedCorrProflg[kNChannels];
410   TF1 *ffunchg[kNChannels], *ffunclg[kNChannels];
411   char namhist4[50];
412   for(int i=0;i<kNChannels;i++) {
413      sprintf(namhist4,"ADCHRvsOOT%d_Prof",i);
414      hPedCorrProfhg[i] = hPedCorrhg[i]->ProfileX(namhist4,-1,-1,"S");
415      hPedCorrProfhg[i]->SetName(namhist4);
416      hPedCorrProfhg[i]->Fit("pol1","Q");
417      ffunchg[i] = hPedCorrProfhg[i]->GetFunction("pol1");
418      CorrCoeff0[i] = (Double_t)  ffunchg[i]->GetParameter(0);
419      CorrCoeff1[i] = (Double_t) ffunchg[i]->GetParameter(1);
420      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i],CorrCoeff1[i]);
421      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",i, CorrCoeff0[i], i, CorrCoeff1[i]);
422   }    
423   for(int i=0;i<kNChannels;i++) {
424      sprintf(namhist4,"ADCLRvsOOT%d_Prof",i);
425      hPedCorrProflg[i] = hPedCorrlg[i]->ProfileX(namhist4,-1,-1,"S");
426      hPedCorrProflg[i]->SetName(namhist4);
427      hPedCorrProflg[i]->Fit("pol1","Q");
428      ffunclg[i] = hPedCorrProflg[i]->GetFunction("pol1");
429      CorrCoeff0[i+kNChannels] =  (Double_t) ffunclg[i]->GetParameter(0);
430      CorrCoeff1[i+kNChannels] =  (Double_t) ffunclg[i]->GetParameter(1);
431      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i+kNChannels],CorrCoeff1[i+kNChannels]);
432      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",
433      //         i+kNChannels, CorrCoeff0[i+kNChannels], i+kNChannels, CorrCoeff1[i+kNChannels]);
434   }    
435
436   //                                                   
437   fclose(fileShuttle);
438   //
439   for(Int_t j=0; j<kNChannels; j++){
440      delete hPedhg[j];
441      delete hPedOutOfTimehg[j];
442      delete hPedCorrhg[j];
443      delete hPedlg[j];
444      delete hPedOutOfTimelg[j];
445      delete hPedCorrlg[j];
446   }
447
448   /* write report */
449   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
450
451   /* close result file */
452   fclose(fp);
453   
454   /* report progress */
455   daqDA_progressReport(90);
456
457   /* store the result files on FES */
458   // [1] File with mapping
459   status = daqDA_FES_storeFile(MAPDATA_FILE,MAPDATA_FILE);
460   if(status){
461     printf("Failed to export file %d to DAQ FES\n",status);
462     return -1;
463   }
464   // [2] File with pedestal data
465   status = daqDA_FES_storeFile(PEDDATA_FILE,PEDDATA_FILE);
466   if(status){
467     printf("Failed to export file %d to DAQ FES\n",status);
468     return -1;
469   }
470   
471   /* store the result files on DB */
472   status = daqDA_DB_storeFile(PEDDATA_FILE,PEDDATA_FILE);  
473   if(status){
474     printf("Failed to export file %d to DAQ DB\n",status);
475     return -1;
476   }
477
478   /* report progress */
479   daqDA_progressReport(100);
480
481   return status;
482 }