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