]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/VZEROda.cxx
automatic histogram scaling for AODs now available
[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 - Reference run number :              
6 - Run Type:   PHYSICS
7 - DA Type:    MON
8 - Number of events needed: >=500
9 - Input Files:  argument list
10 - Output Files: local files  VZERO_Histos.root, V0_Pedestals.dat (Online mapping)
11                 FXS file     V0_Ped_Width_Gain.dat (Offline mapping)
12 - Trigger types used: PHYSICS_EVENT
13 **********************************************************************************/
14
15 /**********************************************************************************
16 *                                                                                 *
17 * VZERO Detector Algorithm used for extracting calibration parameters             *
18 *                                                                                 *
19 * This program connects to the DAQ data source passed as argument.                *
20 * It computes calibration parameters, populates local "./V0_Ped_Width_Gain.dat"   *            
21 * file, exports it to the FES, and stores it into DAQ DB                          *
22 * The program exits when being asked to shut down (daqDA_checkshutdown)           *
23 * or on End of Run event.                                                         *
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 *                                                                                 *
27 ***********************************************************************************/
28
29 // DATE
30 #include "event.h"
31 #include "monitor.h"
32 #include "daqDA.h"
33
34 //AliRoot
35 #include <AliVZERORawStream.h>
36 #include <AliRawReaderDate.h>
37 #include <AliRawReader.h>
38 #include <AliDAQ.h>
39
40 // standard
41 #include <stdio.h>
42 #include <stdlib.h>
43
44 //ROOT
45 #include "TROOT.h"
46 #include "TPluginManager.h"
47 #include <TFile.h>
48 #include <TH1F.h>
49 #include <TMath.h>
50
51 Int_t GetOfflineChannel(Int_t channel);
52
53 /* Main routine --- Arguments: monitoring data source */
54       
55 int main(int argc, char **argv) {
56
57 /* magic line from Cvetan */
58   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
59                     "*",
60                     "TStreamerInfo",
61                     "RIO",
62                     "TStreamerInfo()");
63   int status;
64   if (argc!=2) {
65      printf("Wrong number of arguments\n");
66      return -1;
67   }
68   
69   // Online values (using FEE channel numbering), 
70   // stored into local V0_Pedestals.dat:
71   Double_t ADCmean[128];
72   Double_t ADCsigma[128];
73   Double_t PEDmean[128];
74   Double_t PEDsigma[128];
75   // Offline values(same but ordered as in aliroot for offliners)
76   // stored into V0_Ped_Width_Gain.dat: 
77   Double_t ADCmean_Off[128];
78   Double_t ADCsigma_Off[128];
79   Double_t PEDmean_Off[128];
80   Double_t PEDsigma_Off[128];
81        
82 //___________________________________________________
83 // Get cuts from V00DA.config file
84
85   Int_t    kLowCut;     // = 60;   low cut on signal distribution - to be tuned
86   Int_t    kHighCut;    // = 50;   high cut on pedestal distribution - to be tuned
87   
88   status = daqDA_DB_getFile("V00DA.config","./V00DA.config");
89   if (status) {
90       printf("Failed to get config file (V00DA.config) from DAQ DB, status=%d\n", status);
91       return -1;   
92   }
93   /* open the config file and retrieve cuts */
94   FILE *fpConfig = fopen("V00DA.config","r");
95   fscanf(fpConfig,"%d %d",&kLowCut,&kHighCut);
96   fclose(fpConfig);
97   
98   printf("LowCut on signal = %d ; HighCut on pedestal = %d\n",kLowCut,kHighCut);
99
100 //___________________________________________________
101 // Book HISTOGRAMS - dynamics of p-p collisions -
102       
103   char     ADCname[6]; 
104   char     PEDname[6]; 
105   TH1F     *hADCname[128];
106   TH1F     *hPEDname[128];  
107   
108   char  texte[12];
109   for (Int_t i=0; i<128; i++) {
110        sprintf(ADCname,"hADC%d",i);
111        sprintf(texte,"ADC cell%d",i);
112        hADCname[i]  = new TH1F(ADCname,texte,1024,-0.5, 1023.5);
113        sprintf(PEDname,"hPED%d",i);
114        sprintf(texte,"PED cell%d",i);
115        hPEDname[i]  = new TH1F(PEDname,texte,1024,-0.5, 1023.5);
116   }
117 //___________________________________________________ 
118
119   /* open result file to be exported to FES */
120   FILE *fpLocal=NULL;
121   fpLocal=fopen("./V0_Pedestals.dat","w");
122   if (fpLocal==NULL) {
123       printf("Failed to open local result file\n");
124       return -1;}
125    
126   /* open result file to be exported to FES */
127   FILE *fp=NULL;
128   fp=fopen("./V0_Ped_Width_Gain.dat","w");
129   if (fp==NULL) {
130       printf("Failed to open result file\n");
131       return -1;}
132
133   /* define data source : this is argument 1 */  
134   status=monitorSetDataSource( argv[1] );
135   if (status!=0) {
136      printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
137      return -1;
138   }
139
140   /* declare monitoring program */
141   status=monitorDeclareMp( __FILE__ );
142   if (status!=0) {
143     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
144     return -1;
145   }
146
147   /* define wait event timeout - 1s max */
148   monitorSetNowait();
149   monitorSetNoWaitNetworkTimeout(1000);
150   
151   /* init counters on events */
152   int nevents_physics=0;
153   int nevents_total=0;
154
155   /* loop on events (infinite) */
156   for(;;) {
157       struct eventHeaderStruct *event;
158       eventTypeType eventT;
159
160       /* check shutdown condition */
161       if (daqDA_checkShutdown()) {break;}   
162
163       /* get next event (blocking call until timeout) */
164       status=monitorGetEventDynamic((void **)&event);
165       if (status==MON_ERR_EOF) {
166           printf ("End of File detected\n");
167       break; /* end of monitoring file has been reached */
168       }
169     
170       if (status!=0) {
171            printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
172       break;
173       }
174
175       /* retry if got no event */
176       if (event==NULL) continue;
177                
178       /* decode event */
179       eventT=event->eventType;
180         
181       switch (event->eventType){
182       
183       case START_OF_RUN:
184            break;
185       
186       case END_OF_RUN:
187            printf("End Of Run detected\n");
188            break;
189       
190       case PHYSICS_EVENT:
191            nevents_physics++;
192                          
193            AliRawReader *rawReader = new AliRawReaderDate((void*)event);
194   
195            AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
196            rawStream->Next();   
197            for(Int_t i=0; i<64; i++) {
198                 Int_t nFlag = 0;
199                 for(Int_t j=0; j<21; j++) {
200                    if((rawStream->GetBBFlag(i,j)) || (rawStream->GetBGFlag(i,j))) nFlag++; 
201                 }
202                 if(nFlag == 0){       // Pedestal
203                    for(Int_t j=5;j<16;j++){
204                        Int_t Integrator = rawStream->GetIntegratorFlag(i,j);
205                        Float_t pedestal = (float)(rawStream->GetPedestal(i,j));
206                        hPEDname[i + 64 * Integrator]->Fill(pedestal);
207                    }    
208                 } 
209                 if((rawStream->GetBBFlag(i,10)) || (rawStream->GetBGFlag(i,10))){ // Charge
210                     Int_t Integrator = rawStream->GetIntegratorFlag(i,10);
211                     Float_t charge = (float)(rawStream->GetADC(i));
212                     hADCname[i + 64 * Integrator]->Fill(charge);
213                 }                          
214            }    
215            delete rawStream;
216            rawStream = 0x0;      
217            delete rawReader;
218            rawReader = 0x0;                                                                      
219       } // end of switch on event type 
220         
221       nevents_total++;
222       /* free resources */
223       free(event);
224         
225       /* exit when last event received, no need to wait for TERM signal */
226       if (eventT==END_OF_RUN) {
227            printf("End Of Run event detected\n");
228       break;
229     }
230
231   }  // loop over events
232   
233   printf("%d physics events processed\n",nevents_physics);
234     
235 //___________________________________________________________________________
236 //  Computes mean values, converts FEE channels into Offline AliRoot channels
237 //  and dumps the ordered values into the output text file for SHUTTLE
238         
239   for(Int_t i=0; i<128; i++) {
240       hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
241       PEDmean[i]  = hPEDname[i]->GetMean(); 
242       PEDsigma[i] = hPEDname[i]->GetRMS(); 
243       hADCname[i]->GetXaxis()->SetRange(kLowCut,1024);
244       ADCmean[i]  = hADCname[i]->GetMean();
245       ADCsigma[i] = hADCname[i]->GetRMS(); 
246 //      printf(" i = %d, %.3f %.3f %.3f %.3f\n",i,PEDmean[i],PEDsigma[i],ADCmean[i],ADCsigma[i]);
247       fprintf(fpLocal," %.3f %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],
248                                                ADCmean[i],ADCsigma[i]);
249       if (i < 64) {
250           Int_t j = GetOfflineChannel(i);     
251           PEDmean_Off[j]  = PEDmean[i];
252           PEDsigma_Off[j] = PEDsigma[i];
253           ADCmean_Off[j]  = ADCmean[i];
254           ADCsigma_Off[j] = ADCsigma[i]; }
255       else{
256           Int_t j = GetOfflineChannel(i-64);     
257           PEDmean_Off[j+64]  = PEDmean[i];
258           PEDsigma_Off[j+64] = PEDsigma[i];
259           ADCmean_Off[j+64]  = ADCmean[i];
260           ADCsigma_Off[j+64] = ADCsigma[i];
261       }
262   }
263   
264   for(Int_t j=0; j<128; j++) {
265 //      printf(" j = %d, %.3f %.3f %.3f %.3f\n",j,PEDmean_Off[j],PEDsigma_Off[j],
266 //                                                ADCmean_Off[j],ADCsigma_Off[j]);
267       fprintf(fp," %.3f %.3f %.3f %.3f\n",PEDmean_Off[j],PEDsigma_Off[j],
268                                           ADCmean_Off[j],ADCsigma_Off[j]);                                     
269   }
270    
271 //________________________________________________________________________
272 // Write root file with histos for users further check - just in case - 
273
274   TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
275     
276   for (Int_t i=0; i<128; i++) {
277        hADCname[i]->GetXaxis()->SetRange(0,1024);
278        hADCname[i]->Write(); 
279        hPEDname[i]->Write(); }
280
281   histoFile->Close(); 
282   delete histoFile;
283   
284 //________________________________________________________________________
285    
286   /* close local result file and FXS result file*/
287   fclose(fpLocal);
288   fclose(fp);
289  
290   /* export result file to FES */
291   status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
292   if (status)    {
293       printf("Failed to export file : %d\n",status);
294       return -1; }
295
296   /* store result file into Online DB */
297   status=daqDA_DB_storeFile("./V0_Pedestals.dat","V00da_results");
298   if (status)    {
299       printf("Failed to store file into Online DB: %d\n",status);
300       return -1; }
301       
302   return status;
303 }
304
305  Int_t GetOfflineChannel(Int_t channel) {
306
307 // Channel mapping Online - Offline:
308  
309  Int_t fOfflineChannel[64] = {39, 38, 37, 36, 35, 34, 33, 32, 
310                               47, 46, 45, 44, 43, 42, 41, 40, 
311                               55, 54, 53, 52, 51, 50, 49, 48, 
312                               63, 62, 61, 60, 59, 58, 57, 56,
313                                7,  6,  5,  4,  3,  2,  1,  0, 
314                               15, 14, 13, 12, 11, 10,  9,  8,
315                               23, 22, 21, 20, 19, 18, 17, 16, 
316                               31, 30, 29, 28, 27, 26, 25, 24};
317  return fOfflineChannel[channel];                     
318 }