]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/VZEROda.cxx
Bug fix (Alex)
[u/mrichter/AliRoot.git] / VZERO / VZEROda.cxx
1 /********************************************************************************
2 *                                                                               *
3 * VZERO Detector Algorithm used for extracting calibration parameters           *
4 *                                                                               *
5 * This program reads the DDL data file passed as argument using the monitoring  *
6 * library.                                                                      *
7 * It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat" *           
8 * file and exports it to the FES                                                *
9 *                                                                               *
10 * The program reports about its processing progress.                            *
11 *                                                                               *
12 *********************************************************************************/
13
14 // DATE
15 #include "event.h"
16 #include "monitor.h"
17 #include "daqDA.h"
18
19 //AliRoot
20 #include <AliVZERORawStream.h>
21 #include <AliRawReaderDate.h>
22 #include <AliRawReader.h>
23 #include <AliDAQ.h>
24
25 // standard
26 #include <stdio.h>
27 #include <stdlib.h>
28
29 //ROOT
30
31 #include <TFile.h>
32 #include <TH1F.h>
33 #include <TMath.h>
34
35
36 /* Main routine --- Arguments: list of DATE raw data files */
37       
38 int main(int argc, char **argv) {
39
40   int status;
41
42   printf(" argc = %d, argv = %s \n",argc, &(**argv));
43
44   Int_t    kHighCut = 50; // high cut on pedestal distribution - to be tuned
45   Int_t    kLowCut  = 30; // low cut on signal distribution - to be tuned
46   Double_t ADCmean[64];
47   Double_t PEDmean[64];
48   Double_t PEDsigma[64];
49   
50 //___________________________________________________
51 // Book HISTOGRAMS - dynamics of p-p collisions -
52       
53   char     ADCname[6]; 
54   char     PEDname[6]; 
55   TH1F     *hADCname[64];
56   TH1F     *hPEDname[64];  
57   
58   char  texte[12];
59   for (Int_t i=0; i<64; i++) {
60        sprintf(ADCname,"hADC%d",i);
61        sprintf(texte,"ADC cell%d",i);
62        hADCname[i]  = new TH1F(ADCname,texte,1024,0,2047);
63        sprintf(PEDname,"hPED%d",i);
64        sprintf(texte,"PED cell%d",i);
65        hPEDname[i]  = new TH1F(PEDname,texte,1024,0,2047);}
66 //___________________________________________________ 
67   
68   
69   /* log start of process */
70   printf("VZERO DA program started\n");  
71
72   /* check that we got some arguments = list of files */
73   if (argc<2)   {
74       printf("Wrong number of arguments\n");
75       return -1;}
76
77   /* open result file to be exported to FES */
78   FILE *fp=NULL;
79   fp=fopen("./V0_Ped_Width_Gain.dat","a");
80   if (fp==NULL) {
81       printf("Failed to open result file\n");
82       return -1;}
83
84   /* open log file to inform user */
85   FILE *flog=NULL;
86   flog=fopen("./V00log.txt","a");
87   if (flog==NULL) {
88       printf("Failed to open log file\n");
89       return -1;  }
90     
91   /* report progress */
92   daqDA_progressReport(10);
93
94
95   /* init counters on events */
96   int nevents_physics=0;
97   int nevents_total=0;
98
99   /* read the data files, considering n files */
100   
101   int n;
102   
103   for (n=1;n<argc;n++) {
104    
105     status=monitorSetDataSource( argv[n] );
106     if (status!=0) {
107         printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
108         return -1; }
109
110     /* report progress - here indexed on the number of files */
111     daqDA_progressReport(10+80*n/argc);
112
113     /* read the data file */
114     for(;;) {
115         struct eventHeaderStruct *event;
116         eventTypeType eventT;
117
118         /* get next event */
119         status=monitorGetEventDynamic((void **)&event);
120         if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
121         if (status!=0) {
122             printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
123             return -1; }
124
125         /* retry if got no event */
126         if (event==NULL) break;
127         
128         /* decode event */
129         eventT=event->eventType;
130         
131         switch (event->eventType){
132       
133         case START_OF_RUN:
134              break;
135       
136         case END_OF_RUN:
137              printf("End Of Run detected\n");
138              break;
139       
140         case PHYSICS_EVENT:
141              nevents_physics++;
142         
143              fprintf(flog,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
144                  (unsigned long)event->eventRunNb,
145                  (unsigned long)event->eventSize,
146                  EVENT_ID_GET_BUNCH_CROSSING(event->eventId),
147                  EVENT_ID_GET_ORBIT(event->eventId),
148                  EVENT_ID_GET_PERIOD(event->eventId) );
149                  
150              AliRawReader *rawReader = new AliRawReaderDate((void*)event);
151              rawReader->RequireHeader(kFALSE);   
152              AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader);    
153              rawStream->Next(); 
154              for(Int_t i=0; i<64; i++) {
155                 hADCname[i]->Fill(rawStream->GetADC(i)); 
156                 for(Int_t j=0; j<21; j++) {
157                     if(j==10) continue;
158                     hPEDname[i]->Fill(rawStream->GetPedestal(i,j)); }
159              }    
160              delete rawStream;
161              rawStream = 0x0;      
162              delete rawReader;
163              rawReader = 0x0;                                                                    
164         } // end of switch on event type 
165         
166         nevents_total++;
167         /* free resources */
168         free(event);
169
170     }  // loop over events
171     
172    }  // loop over data files
173 //________________________________________________________________________
174 //  Computes mean values, dumps them into the output text file
175         
176   for(Int_t i=0; i<64; i++) {
177       hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
178       PEDmean[i]  = hPEDname[i]->GetMean(); 
179       PEDsigma[i] = hPEDname[i]->GetRMS(); 
180       hADCname[i]->GetXaxis()->SetRange(kLowCut,2047);
181       ADCmean[i] = hADCname[i]->GetMean() ; 
182       fprintf(fp," %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i]);
183   } 
184
185 //________________________________________________________________________
186 // Write root file with histos for users further check - just in case - 
187
188   TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
189     
190   for (Int_t i=0; i<64; i++) {
191        hADCname[i]->Write(); 
192        hPEDname[i]->Write(); }
193
194   histoFile->Close(); 
195   
196 //________________________________________________________________________
197    
198   /* write report */
199   fprintf(flog,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
200   printf("Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
201   
202   /* close result and log files */
203   fclose(fp);
204   fclose(flog); 
205   
206   /* report progress */
207   daqDA_progressReport(90);
208
209
210   /* export result file to FES */
211   status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
212   if (status)    {
213       printf("Failed to export file : %d\n",status);
214       return -1; }
215
216   /* report progress */
217   daqDA_progressReport(100);
218   
219   return status;
220 }