]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AD/DA/AD0PEDESTALda.cxx
CMake: Retrieve Git information
[u/mrichter/AliRoot.git] / AD / DA / AD0PEDESTALda.cxx
1 /*********************************************************************************
2 - Contact:    Michal Broz     Michal.Broz@cern.ch
3 - Link:       http
4 - Raw data test file :          
5 - Reference run number :              
6 - Run Type:   PEDESTAL
7 - DA Type:    LDC
8 - Number of events needed: >=1000
9 - Input Files:  argument list
10 - Output Files: local files  AD0_Histos.root, AD0_Pedestals_On.dat (Online mapping)
11                 FXS file     AD0_Pedestals_Off.dat (Offline mapping)
12 - Trigger types used: PHYSICS_EVENT
13 **********************************************************************************/
14
15 /**********************************************************************************
16 *                                                                                 *
17 * AD 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 "./AD0_Pedestals_On.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 32 channels instead of 16 as expected for AD0 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 <AliADConst.h>
36 #include <AliADRawStream.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 /* Main routine --- Arguments: monitoring data source */
53       
54 int main(int argc, char **argv) {
55
56 /* magic line from Cvetan */
57   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo","*","TStreamerInfo","RIO","TStreamerInfo()");
58   
59   int status;
60   if (argc!=2) {
61      printf("Wrong number of arguments\n");
62      return -1;
63   }
64   
65   // Online values (using FEE channel numbering), 
66   // stored into local AD0_Pedestals_On.dat:
67   Double_t adcMean[32];
68   Double_t adcSigma[32];
69   Double_t pedMean[32];
70   Double_t pedSigma[32];
71   // Offline values(same but ordered as in aliroot for offliners)
72   // stored into AD0_Pedestals_Off.dat: 
73   Double_t adcMeanOff[32];
74   Double_t adcSigmaOff[32];
75   Double_t pedMeanOff[32];
76   Double_t pedSigmaOff[32];
77        
78 //___________________________________________________
79 // Get cuts from AD0_Pedestal_DA.config file
80
81   Int_t    kClockMin;   // = 16;   LHC Clock Min for pedestal calculation
82   Int_t    kClockMax;   // = 19;   LHC Clock Max for pedestal calculation
83   Int_t    kLowCut;     // = 60;   low cut on signal distribution - to be tuned
84   Int_t    kHighCut;    // = 50;   high cut on pedestal distribution - to be tuned
85   Int_t    kClockMinRef;   // = 16;   LHC Clock Min for Flag checking
86   Int_t    kClockMaxRef;   // = 20;   LHC Clock Max for Flag checking
87   Float_t  kChi2Max;            // = 1.  Maximum chi2 
88
89   status = daqDA_DB_getFile("AD0_Pedestal_DA.config","./AD0_Pedestal_DA.config");
90   if (status) {
91       printf("Failed to get Config file (AD0_Pedestal_DA.config) from DAQ DB, status=%d\n", status);
92       printf("Take default values of parameters for pedestal calculation \n");
93       kClockMin  =  0; 
94       kClockMax  =  20; 
95       kLowCut    =  60;   
96       kHighCut   =  50;  
97       kClockMinRef  =  0; 
98       kClockMaxRef  =  20; 
99       kChi2Max          =  100.;
100   } else {
101       /* open the config file and retrieve cuts */
102       FILE *fpConfig = fopen("AD0_Pedestal_DA.config","r");
103       int res = fscanf(fpConfig,"%d %d %d %d %d %d %f",&kClockMin,&kClockMax,&kLowCut,&kHighCut,&kClockMinRef,&kClockMaxRef,&kChi2Max);
104       if(res!=7) {
105             printf("Failed to get values from Config file (AD0_Pedestal_DA.config): wrong file format - 7 integers are expected - \n");
106             kClockMin  =  0; 
107         kClockMax  =  20; 
108         kLowCut    =  60;   
109         kHighCut   =  50; 
110         kClockMinRef  =  0; 
111         kClockMaxRef  =  20; 
112         kChi2Max          =  100.;
113       }
114       fclose(fpConfig);
115   }
116   
117   printf("LHC Clock Min for pedestal calculation = %d; LHC Clock Max for pedestal calculation = %d; LowCut on signal = %d ; HighCut on pedestal = %d\n",
118           kClockMin, kClockMax, kLowCut, kHighCut);
119
120 //___________________________________________________
121 // Book HISTOGRAMS - dynamics of p-p collisions -
122       
123   char     adcName[6]; 
124   char     pedName[6]; 
125   TH1F     *hADCname[32];
126   TH1F     *hPEDname[32];  
127   
128   char  texte[12];
129   for (Int_t i=0; i<16; i++) {
130         for(Int_t j=0; j<2; j++){
131                 sprintf(adcName,"hADC%d%d",i,j);
132                 sprintf(texte,"ADC cell%d int%d",i,j);
133                 hADCname[i + 16*j]  = new TH1F(adcName,texte,1024,-0.5, 1023.5);
134                 sprintf(pedName,"hPED%d%d",i,j);
135                 sprintf(texte,"PED cell%d int%d",i,j);
136                 hPEDname[i + 16*j]  = new TH1F(pedName,texte,1024,-0.5, 1023.5);
137                 }
138        }
139 //___________________________________________________ 
140
141   /* open result file to be exported to FES */
142   FILE *fpLocal=NULL;
143   fpLocal=fopen("./AD0_Pedestals_On.dat","w");
144   if (fpLocal==NULL) {
145       printf("Failed to open local result file\n");
146       return -1;}
147    
148   /* open result file to be exported to FES */
149   FILE *fp=NULL;
150   fp=fopen("./AD0_Pedestals_Off.dat","w");
151   if (fp==NULL) {
152       printf("Failed to open result file\n");
153       return -1;}
154
155   /* define data source : this is argument 1 */  
156   status=monitorSetDataSource( argv[1] );
157   if (status!=0) {
158      printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
159      return -1;
160   }
161
162   /* declare monitoring program */
163   status=monitorDeclareMp( __FILE__ );
164   if (status!=0) {
165     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
166     return -1;
167   }
168
169   /* define wait event timeout - 1s max */
170   monitorSetNowait();
171   monitorSetNoWaitNetworkTimeout(1000);
172   
173   /* init counters on events */
174   int neventsPhysics=0;
175   int neventsTotal=0;
176
177   /* loop on events (infinite) */
178   for(;;) {
179       struct eventHeaderStruct *event;
180       eventTypeType eventT;
181
182       /* check shutdown condition */
183       if (daqDA_checkShutdown()) {break;}   
184
185       /* get next event (blocking call until timeout) */
186       status=monitorGetEventDynamic((void **)&event);
187       if (status==MON_ERR_EOF) {
188           printf ("End of File detected\n");
189       break; /* end of monitoring file has been reached */
190       }
191     
192       if (status!=0) {
193            printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
194       break;
195       }
196
197       /* retry if got no event */
198       if (event==NULL) continue;
199                
200       /* decode event */
201       eventT=event->eventType;
202         
203       switch (event->eventType){
204       
205       case START_OF_RUN:
206            break;
207       
208       case END_OF_RUN:
209            printf("End Of Run detected\n");
210            break;
211       
212       case PHYSICS_EVENT:
213            neventsPhysics++;
214                          
215            AliRawReader *rawReader = new AliRawReaderDate((void*)event);
216   
217            AliADRawStream* rawStream  = new AliADRawStream(rawReader); 
218            if (rawStream->Next()) {     
219            for(Int_t i=0; i<16; i++) {
220                         Int_t nFlag = 0;
221                         for(Int_t j=kClockMinRef; j <= kClockMaxRef; j++) {  // Check flags on clock range used for pedestal calculation
222                                 if((rawStream->GetBBFlag(i,j)) || (rawStream->GetBGFlag(i,j))) nFlag++; 
223                         }
224                         if(nFlag == 0){       // Fill 16*2 pedestal histograms  - 2 integrators -
225                                 Float_t sum[2] = {0.,0.};
226                                 Float_t sumwi[2] = {0.,0.};
227                                 for(Int_t j=kClockMin;j <= kClockMax;j++){
228                                 Int_t integrator = rawStream->GetIntegratorFlag(i,j);
229                                 Float_t pedestal = (float)(rawStream->GetPedestal(i,j));
230                                 sum[integrator] += pedestal;
231                                         sumwi[integrator] += 1.;
232                                 }       
233                                 Float_t mean[2] =  {0.,0.};
234                                 Float_t chi2[2] =  {0.,0.};
235                                 
236                                 for(int ii=0;ii<2;ii++) if(sumwi[ii]>1.e-6) mean[ii] = sum[ii] / sumwi[ii];
237
238                                 for(Int_t j=kClockMin;j <= kClockMax;j++){
239                                 Int_t integrator = rawStream->GetIntegratorFlag(i,j);
240                                 Float_t pedestal = (float)(rawStream->GetPedestal(i,j));
241                                         chi2[integrator] += (mean[integrator] - pedestal) * (mean[integrator] - pedestal);
242                                 }       
243                                 if(chi2[0]<kChi2Max && chi2[1]<kChi2Max) {
244                                         for(int ii=0;ii<2;ii++){
245                                                 if(mean[ii] >1.e-6) hPEDname[i + 16*ii]->Fill(mean[ii]);
246                                         }
247                                 }
248
249                         } 
250                         if((rawStream->GetBBFlag(i,10)) || (rawStream->GetBGFlag(i,10))){ // Charge
251                         Int_t integrator = rawStream->GetIntegratorFlag(i,10);
252                         Float_t charge = (float)(rawStream->GetADC(i));   // Fill 16*2 ADCmax histograms 
253                         hADCname[i + 16 * integrator]->Fill(charge);
254                         }                          
255         } // End loop over channels
256            }   // End : if rawstream
257            delete rawStream;
258            rawStream = 0x0;      
259            delete rawReader;
260            rawReader = 0x0;                                                                      
261       } // end of switch on event type 
262         
263       neventsTotal++;
264       /* free resources */
265       free(event);
266         
267       /* exit when last event received, no need to wait for TERM signal */
268       if (eventT==END_OF_RUN) {
269            printf("End Of Run event detected\n");
270       break;
271     }
272
273   }  // loop over events
274   
275   printf("%d physics events processed\n",neventsPhysics);
276     
277 //___________________________________________________________________________
278 //  Computes mean values, converts FEE channels into Offline AliRoot channels
279 //  and dumps the ordered values into the output text file for SHUTTLE
280         
281   for(Int_t i=0; i<32; i++) {
282       hPEDname[i]->GetXaxis()->SetRange(0,kHighCut);
283       pedMean[i]  = hPEDname[i]->GetMean(); 
284       pedSigma[i] = hPEDname[i]->GetRMS(); 
285       hADCname[i]->GetXaxis()->SetRange(kLowCut,1024);
286       adcMean[i]  = hADCname[i]->GetMean();
287       adcSigma[i] = hADCname[i]->GetRMS(); 
288 //      printf(" i = %d, %.3f %.3f %.3f %.3f\n",i,pedMean[i],pedSigma[i],adcMean[i],adcSigma[i]);
289       fprintf(fpLocal," %.3f %.3f %.3f %.3f\n",pedMean[i],pedSigma[i],adcMean[i],adcSigma[i]);
290       
291       if (i < 16) {
292           Int_t j = kOfflineChannel[i];;     
293           pedMeanOff[j]  = pedMean[i];
294           pedSigmaOff[j] = pedSigma[i];
295           adcMeanOff[j]  = adcMean[i];
296           adcSigmaOff[j] = adcSigma[i]; 
297           }
298       else{
299           Int_t j = kOfflineChannel[i-16];     
300           pedMeanOff[j+16]  = pedMean[i];
301           pedSigmaOff[j+16] = pedSigma[i];
302           adcMeanOff[j+16]  = adcMean[i];
303           adcSigmaOff[j+16] = adcSigma[i];
304       }
305   }
306   
307   for(Int_t j=0; j<32; j++) {
308 //      printf(" j = %d, %.3f %.3f %.3f %.3f\n",j,pedMeanOff[j],pedSigmaOff[j],
309 //                                                adcMeanOff[j],adcSigmaOff[j]);
310       fprintf(fp," %.3f %.3f %.3f %.3f\n",pedMeanOff[j],pedSigmaOff[j],
311                                           adcMeanOff[j],adcSigmaOff[j]);                                       
312   }
313    
314 //________________________________________________________________________
315 // Write root file with histos for users further check - just in case - 
316
317   TFile *histoFile = new TFile("AD0_Pedestals.root","RECREATE");
318     
319   for (Int_t i=0; i<32; i++) {
320        hADCname[i]->GetXaxis()->SetRange(0,1024);
321        hADCname[i]->Write(); 
322        hPEDname[i]->Write(); }
323
324   histoFile->Close(); 
325   delete histoFile;
326   
327 //________________________________________________________________________
328    
329   /* close local result file and FXS result file*/
330   fclose(fpLocal);
331   fclose(fp);
332  
333   /* export result file to FES */
334   status=daqDA_FES_storeFile("./AD0_Pedestals_Off.dat","AD0da_results");
335   if (status)    {
336       printf("Failed to export file : %d\n",status);
337       return -1; }
338
339   /* store result file into Online DB */
340   status=daqDA_DB_storeFile("./AD0_Pedestals_On.dat","AD0da_results");
341   if (status)    {
342       printf("Failed to store file into Online DB: %d\n",status);
343       return -1; }
344       
345   return status;
346 }
347