]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/VZEROda.cxx
2cc1255aabd9d6fea8eefdb47432ca11874f99c6
[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_Ped_Width_Gain.dat
11                 FXS file     V0_Ped_Width_Gain.dat
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 and exports it to the FES.                                                 *
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 /* Main routine --- Arguments: monitoring data source */
52       
53 int main(int argc, char **argv) {
54
55 /* magic line from Cvetan */
56   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
57                     "*",
58                     "TStreamerInfo",
59                     "RIO",
60                     "TStreamerInfo()");
61   int status;
62   if (argc!=2) {
63      printf("Wrong number of arguments\n");
64      return -1;
65   }
66   
67   Double_t ADCmean[128];
68   Double_t ADCsigma[128];
69   Double_t PEDmean[128];
70   Double_t PEDsigma[128];
71      
72 //___________________________________________________
73 // Get cuts from V00DA.config file
74
75   Int_t    kHighCut;    // = 50; high cut on pedestal distribution - to be tuned
76   Int_t    kLowCut;     // = 30; low cut on signal distribution - to be tuned
77
78   status = daqDA_DB_getFile("V00DA.config","./V00DA.config");
79   if (status) {
80       printf("Failed to get config file (V00DA.config) from DAQdetDB, status=%d\n", status);
81       return -1;   
82   }
83   /* open the config file and retrieve cuts */
84   FILE *fpConfig = fopen("V00DA.config","r");
85   fscanf(fpConfig,"%d %d",&kLowCut,&kHighCut);
86   fclose(fpConfig);
87   
88   printf("LowCut on signal = %d ; HighCut on pedestal = %d\n",kLowCut,kHighCut);
89
90 //___________________________________________________
91 // Book HISTOGRAMS - dynamics of p-p collisions -
92       
93   char     ADCname[6]; 
94   char     PEDname[6]; 
95   TH1F     *hADCname[128];
96   TH1F     *hPEDname[128];  
97   
98   char  texte[12];
99   for (Int_t i=0; i<128; i++) {
100        sprintf(ADCname,"hADC%d",i);
101        sprintf(texte,"ADC cell%d",i);
102        hADCname[i]  = new TH1F(ADCname,texte,1024,0,1023);
103        sprintf(PEDname,"hPED%d",i);
104        sprintf(texte,"PED cell%d",i);
105        hPEDname[i]  = new TH1F(PEDname,texte,1024,0,1023);
106   }
107 //___________________________________________________ 
108   
109  
110   /* open result file to be exported to FES */
111   FILE *fp=NULL;
112   fp=fopen("./V0_Ped_Width_Gain.dat","a");
113   if (fp==NULL) {
114       printf("Failed to open result file\n");
115       return -1;}
116
117   /* define data source : this is argument 1 */  
118   status=monitorSetDataSource( argv[1] );
119   if (status!=0) {
120      printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
121      return -1;
122   }
123
124   /* declare monitoring program */
125   status=monitorDeclareMp( __FILE__ );
126   if (status!=0) {
127     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
128     return -1;
129   }
130
131   /* define wait event timeout - 1s max */
132   monitorSetNowait();
133   monitorSetNoWaitNetworkTimeout(1000);
134   
135   /* init counters on events */
136   int nevents_physics=0;
137   int nevents_total=0;
138
139   /* loop on events (infinite) */
140   for(;;) {
141       struct eventHeaderStruct *event;
142       eventTypeType eventT;
143
144       /* check shutdown condition */
145       if (daqDA_checkShutdown()) {break;}   
146
147       /* get next event (blocking call until timeout) */
148       status=monitorGetEventDynamic((void **)&event);
149       if (status==MON_ERR_EOF) {
150           printf ("End of File detected\n");
151       break; /* end of monitoring file has been reached */
152       }
153     
154       if (status!=0) {
155            printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
156       break;
157       }
158
159       /* retry if got no event */
160       if (event==NULL) continue;
161                
162       /* decode event */
163       eventT=event->eventType;
164         
165       switch (event->eventType){
166       
167       case START_OF_RUN:
168            break;
169       
170       case END_OF_RUN:
171            printf("End Of Run detected\n");
172            break;
173       
174       case PHYSICS_EVENT:
175            nevents_physics++;
176                          
177            AliRawReader *rawReader = new AliRawReaderDate((void*)event);
178   
179            AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
180            rawStream->Next();   
181            for(Int_t i=0; i<64; i++) {
182                if(!rawStream->GetIntegratorFlag(i,10))
183                    hADCname[i]->Fill(float(rawStream->GetADC(i)));       // even integrator - fills 0 to 63
184                else 
185                    hADCname[i+64]->Fill(float(rawStream->GetADC(i)));    // odd integrator  - fills 64 to 123
186                for(Int_t j=0; j<21; j++) {
187                    if(j==10) continue;
188                    if(!rawStream->GetIntegratorFlag(i,j))
189                        { hPEDname[i]->Fill(float(rawStream->GetPedestal(i,j))); }     // even integrator
190                    else 
191                        { hPEDname[i+64]->Fill(float(rawStream->GetPedestal(i,j))); }  // odd integrator 
192                }
193            }    
194            delete rawStream;
195            rawStream = 0x0;      
196            delete rawReader;
197            rawReader = 0x0;                                                                      
198       } // end of switch on event type 
199         
200       nevents_total++;
201       /* free resources */
202       free(event);
203         
204       /* exit when last event received, no need to wait for TERM signal */
205       if (eventT==END_OF_RUN) {
206            printf("End Of Run event detected\n");
207       break;
208     }
209
210   }  // loop over events
211   
212   printf("%d physics events processed\n",nevents_physics);
213     
214 //________________________________________________________________________
215 //  Computes mean values, dumps them into the output text file
216         
217   for(Int_t i=0; i<128; i++) {
218       hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
219       PEDmean[i]  = hPEDname[i]->GetMean(); 
220       PEDsigma[i] = hPEDname[i]->GetRMS(); 
221       hADCname[i]->GetXaxis()->SetRange(kLowCut,1023);
222       ADCmean[i]  = hADCname[i]->GetMean();
223       ADCsigma[i] = hADCname[i]->GetRMS(); 
224       fprintf(fp," %.3f %.3f %.3f %.3f\n",PEDmean[i],PEDsigma[i],ADCmean[i],ADCsigma[i]);
225   }
226    
227 //________________________________________________________________________
228 // Write root file with histos for users further check - just in case - 
229
230   TFile *histoFile = new TFile("VZERO_histos.root","RECREATE");
231     
232   for (Int_t i=0; i<128; i++) {
233        hADCname[i]->GetXaxis()->SetRange(0,1023);
234        hADCname[i]->Write(); 
235        hPEDname[i]->Write(); }
236
237   histoFile->Close(); 
238   delete histoFile;
239   
240 //________________________________________________________________________
241    
242
243   /* close result file */
244   fclose(fp);
245  
246   /* export result file to FES */
247   status=daqDA_FES_storeFile("./V0_Ped_Width_Gain.dat","V00da_results");
248   if (status)    {
249       printf("Failed to export file : %d\n",status);
250       return -1; }
251
252   return status;
253 }