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