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