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