]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCPEDESTALda.cxx
Added description of MUONStatusMap.C macro
[u/mrichter/AliRoot.git] / ZDC / ZDCPEDESTALda.cxx
1 /*
2
3 DAcase2.c
4
5 This program connects to the DAQ data source passed as argument
6 and populates local "./result.txt" file with the ids of events received
7 during the run.
8
9 The program exits when being asked to shut down (daqDA_checkshutdown)
10 or End of Run event.
11
12 Messages on stdout are exported to DAQ log system.
13
14 DA for ZDC standalone pedestal runs
15
16 Contact: Chiara.Oppedisano@to.infn.it
17 Link: /afs/cern.ch/user/c/chiarao/public/RawPed.date
18 Run Type: STANDALONE_PEDESTAL_RUN
19 DA Type: MON
20 Number of events needed: no constraint (tipically ~10^3)
21 Input Files: 
22 Output Files: ZDCPedestal.dat
23 Trigger Types Used: Standalone Trigger
24
25 */
26
27 #include <stdio.h>
28 #include <Riostream.h>
29
30 // DATE
31 #include <daqDA.h>
32 #include <event.h>
33 #include <monitor.h>
34
35 //ROOT
36 #include <TRandom.h>
37 #include <TH1F.h>
38 #include <TH2F.h>
39 #include <TProfile.h>
40 #include <TF1.h>
41 #include <TFile.h>
42
43 //AliRoot
44 #include <AliRawReaderDate.h>
45 #include <AliZDCRawStream.h>
46
47
48 /* Main routine
49       Arguments: 
50       1- monitoring data source
51 */
52 int main(int argc, char **argv) {
53
54   TH1F::AddDirectory(0);
55   // --- Histograms for ADC pedestals 
56   //     [22 signal channels + 2 reference PTMs]  x 2 gain chains
57   //
58   int const kNChannels = 24;
59   TH1F *hPedhg[kNChannels], *hPedOutOfTimehg[kNChannels];
60   TH2F *hPedCorrhg[kNChannels];
61   TH1F *hPedlg[kNChannels], *hPedOutOfTimelg[kNChannels];
62   TH2F *hPedCorrlg[kNChannels];
63   //
64   char namhist1hg[50], namhist2hg[50], namhist3hg[50];
65   char namhist1lg[50], namhist2lg[50], namhist3lg[50];
66   for(Int_t j=0; j<kNChannels; j++){
67      if(j<5){ // ZN1
68        sprintf(namhist1hg,"PedZN1hg_%d",j);
69        sprintf(namhist2hg,"PedZN1hgOutOfTime_%d",j);
70        sprintf(namhist3hg,"PedCorrZN1hg_%d",j);
71        //
72        sprintf(namhist1lg,"PedZN1lg_%d",j);
73        sprintf(namhist2lg,"PedZN1lgOutOfTime_%d",j);
74        sprintf(namhist3lg,"PedCorrZN1lg_%d",j);
75      }
76      else if(j>=5 && j<10){ // ZP1
77        sprintf(namhist1hg,"PedZP1hg_%d",j-5);
78        sprintf(namhist2hg,"PedZP1hgOutOfTime_%d",j-5);
79        sprintf(namhist3hg,"PedCorrZP1hg_%d",j-5);
80        //
81        sprintf(namhist1lg,"PedZP1lg_%d",j-5);
82        sprintf(namhist2lg,"PedZP1lgOutOfTime_%d",j-5);
83        sprintf(namhist3lg,"PedCorrZP1lg_%d",j-5);       
84      }
85      else if(j>=10 && j<12){ // ZEM
86        sprintf(namhist1hg,"PedZEMhg_%d",j-10);
87        sprintf(namhist2hg,"PedZEMhgOutOfTime_%d",j-10);
88        sprintf(namhist3hg,"PedCorrZEMhg_%d",j-10);
89        //
90        sprintf(namhist1lg,"PedZEMlg_%d",j-10);
91        sprintf(namhist2lg,"PedZEMlgOutOfTime_%d",j-10);
92        sprintf(namhist3lg,"PedCorrZEMlg_%d",j-10);
93      }
94      else if(j>=12 && j<17){ // ZN2
95        sprintf(namhist1hg,"PedZN2hg_%d",j-12);
96        sprintf(namhist2hg,"PedZN2hgOutOfTime_%d",j-12);
97        sprintf(namhist3hg,"PedCorrZN2hg_%d",j-12);
98        //
99        sprintf(namhist1lg,"PedZN2lg_%d",j-12);
100        sprintf(namhist2lg,"PedZN2lgOutOfTime_%d",j-12);
101        sprintf(namhist3lg,"PedCorrZN2lg_%d",j-12);
102      }
103      else if(j>=17 && j<22){ // ZP2
104        sprintf(namhist1hg,"PedZP2hg_%d",j-17);
105        sprintf(namhist2hg,"PedZP2hgOutOfTime_%d",j-17);
106        sprintf(namhist3hg,"PedCorrZP2hg_%d",j-17);
107        //
108        sprintf(namhist1lg,"PedZP2lg_%d",j-17);
109        sprintf(namhist2lg,"PedZP2lgOutOfTime_%d",j-17);
110        sprintf(namhist3lg,"PedCorrZP2lg_%d",j-17);
111      }
112      else if(j>=22 && j<24){ //Reference PMs
113        sprintf(namhist1hg,"PedRefhg_%d",j-22);
114        sprintf(namhist2hg,"PedRefhgOutOfTime_%d",j-22);
115        sprintf(namhist3hg,"PedCorrRefhg_%d",j-22);
116        //
117        sprintf(namhist1lg,"PedReflg_%d",j-22);
118        sprintf(namhist2lg,"PedReflgOutOfTime_%d",j-22);
119        sprintf(namhist3lg,"PedCorrReflg_%d",j-22);
120      }
121      // --- High gain chain histos
122      hPedhg[j] = new TH1F(namhist1hg, namhist1hg, 100,0., 200.);
123      hPedOutOfTimehg[j] = new TH1F(namhist2hg, namhist2hg, 100,0., 200.);
124      hPedCorrhg[j] = new TH2F(namhist3hg,namhist3hg,100,0.,200.,100,0.,200.);
125      // --- Low gain chain histos
126      hPedlg[j] = new TH1F(namhist1lg, namhist1lg, 100,0., 600.);
127      hPedOutOfTimelg[j] = new TH1F(namhist2lg, namhist2lg, 100,0., 600.);
128      hPedCorrlg[j] = new TH2F(namhist3lg,namhist3lg,100,0.,600.,100,0.,600.);
129   }
130
131   int status;
132   
133   if (argc!=2) {
134     printf("Wrong number of arguments\n");
135     return -1;
136   }
137
138
139   /* open result file */
140   FILE *fp=NULL;
141   fp=fopen("./result.txt","a");
142   if (fp==NULL) {
143     printf("Failed to open file\n");
144     return -1;
145   }
146   
147
148   /* define data source : this is argument 1 */  
149   status = monitorSetDataSource( argv[1] );
150   if(status!=0) {
151     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
152     return -1;
153   }
154
155
156   /* declare monitoring program */
157   status = monitorDeclareMp( __FILE__ );
158   if (status!=0) {
159     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
160     return -1;
161   }
162
163
164   /* define wait event timeout - 1s max */
165   monitorSetNowait();
166   monitorSetNoWaitNetworkTimeout(1000);
167   
168
169   /* log start of process */
170   printf("ZDC PEDESTAL monitoring program started\n");  
171   
172   /* init some counters */
173   int nevents_physics=0;
174   int nevents_total=0;
175
176   struct equipmentStruct *equipment;
177   int *eventEnd;
178   int *eventData;
179   int *equipmentEnd;
180   int *equipmentData;
181   int *equipmentID;
182
183   struct eventHeaderStruct *event;
184   eventTypeType eventT;
185   Int_t iev=0;
186   
187   /* main loop (infinite) */
188   for(;;) {
189   
190     /* check shutdown condition */
191     if (daqDA_checkShutdown()) {break;}
192     
193     /* get next event (blocking call until timeout) */
194     status=monitorGetEventDynamic((void **)&event);
195     if (status==MON_ERR_EOF) {
196       printf ("End of File detected\n");
197       break; /* end of monitoring file has been reached */
198     }
199     
200     if (status!=0) {
201       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
202       break;
203     }
204
205     /* retry if got no event */
206     if (event==NULL) {
207       continue;
208     }
209
210     iev++; 
211
212     /* use event - here, just write event id to result file */
213     eventT=event->eventType;
214     
215     if(eventT==PHYSICS_EVENT){
216       //fprintf(fp,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
217       //
218       // Initalize raw-data reading and decoding
219       AliRawReader *reader = new AliRawReaderDate((void*)event);
220       const AliRawDataHeader* header = reader->GetDataHeader();
221       if(header) {
222          UChar_t message = header->GetAttributes();
223          if(message & 0x20){ // DEDICATED PEDESTAL RUN
224             printf("\t STANDALONE_PEDESTAL_RUN raw data found\n");
225             continue;
226          }
227          else{
228             printf("\t NO STANDALONE_PEDESTAL_RUN raw data found\n");
229             return -1;
230          }
231       }
232       //Commented until we won't have true Raw Data Header...
233       //else{
234       //   printf("\t ATTENTION! No Raw Data Header found!!!\n");
235         // return -1;
236       //}
237       //
238       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
239       //
240       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
241       Int_t counter=0;
242       Int_t RawADChg[kNChannels], RawADCoothg[kNChannels];
243       Int_t RawADClg[kNChannels], RawADCootlg[kNChannels];
244       for(Int_t j=0; j<kNChannels; j++){
245          RawADChg[j]=0; RawADCoothg[j]=0;
246          RawADClg[j]=0; RawADCootlg[j]=0;
247       }
248       while(rawStreamZDC->Next()){
249         Int_t index=-1;
250         if(rawStreamZDC->IsADCDataWord()){
251           if(rawStreamZDC->GetSector(1)!=5){ // Physics signals
252             if(rawStreamZDC->GetSector(0)==1) index = rawStreamZDC->GetSector(1); // *** ZN1
253             else if(rawStreamZDC->GetSector(0)==2) index = rawStreamZDC->GetSector(1)+5; // *** ZP1
254             else if(rawStreamZDC->GetSector(0)==3) index = rawStreamZDC->GetSector(1)+9; // *** ZEM
255             else if(rawStreamZDC->GetSector(0)==4) index = rawStreamZDC->GetSector(1)+12; // *** ZN2
256             else if(rawStreamZDC->GetSector(0)==5) index = rawStreamZDC->GetSector(1)+17; // *** ZP2
257           }
258           else{ // Reference PMs
259             index = (rawStreamZDC->GetSector(0)-1)/3+22;
260           }
261           //
262           /*printf("\t counter %d index %d det %d quad %d res %d ADC %d\n", counter, index,
263                 rawStreamZDC->GetSector(0), rawStreamZDC->GetSector(1), 
264                 rawStreamZDC->GetADCGain(), rawStreamZDC->GetADCValue());
265           */
266           //
267           if(counter<2*kNChannels){ // --- In-time pedestals (1st 48 raw data)
268             if(rawStreamZDC->GetADCGain()==0){ 
269               hPedhg[index]->Fill(rawStreamZDC->GetADCValue()); 
270               RawADChg[index] = rawStreamZDC->GetADCValue();
271             }
272             else{
273               hPedlg[index]->Fill(rawStreamZDC->GetADCValue()); 
274               RawADClg[index] = rawStreamZDC->GetADCValue();
275             }
276           }
277           else{  // --- Out-of-time pedestals
278             if(rawStreamZDC->GetADCGain()==0){
279               hPedOutOfTimehg[index]->Fill(rawStreamZDC->GetADCValue());
280               RawADCoothg[index] = rawStreamZDC->GetADCValue();
281             }
282             else{
283               hPedOutOfTimelg[index]->Fill(rawStreamZDC->GetADCValue());
284               RawADCootlg[index] = rawStreamZDC->GetADCValue();
285             }
286           }
287           counter++;
288          }//IsADCDataWord()
289          //
290          if(counter == 4*kNChannels){ // Last ADC channel -> Filling correlation histos
291            for(Int_t k=0; k<kNChannels; k++){
292             hPedCorrhg[k]->Fill(RawADCoothg[k], RawADChg[k]);
293             hPedCorrlg[k]->Fill(RawADCootlg[k], RawADClg[k]);
294            }
295          }
296        }
297        //
298        nevents_physics++;
299
300     }
301     
302     nevents_total++;
303
304
305     /* free resources */
306     free(event);
307     
308     /* exit when last event received, no need to wait for TERM signal */
309     if (eventT==END_OF_RUN) {
310       printf("EOR event detected\n");
311       break;
312     }
313   }
314   
315   /* Analysis of the histograms */
316   //
317   FILE *fileShuttle;
318   fileShuttle = fopen("ZDCPedestal.dat","w");
319   //
320   Float_t MeanPed[2*kNChannels], MeanPedWidth[2*kNChannels], 
321         MeanPedOOT[2*kNChannels], MeanPedWidthOOT[2*kNChannels],
322         CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
323   // --- In-time pedestals
324   TF1 *ADCfunchg[kNChannels];
325   for(Int_t i=0; i<kNChannels; i++){
326      hPedhg[i]->Fit("gaus","Q");
327      ADCfunchg[i] = hPedhg[i]->GetFunction("gaus");
328      MeanPed[i] = ADCfunchg[i]->GetParameter(1);
329      MeanPedWidth[i] = ADCfunchg[i]->GetParameter(2);
330      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i],MeanPedWidth[i]);
331      //printf("\t MeanPed[%d] = %f\n",i, MeanPed[i]);
332   }
333   TF1 *ADCfunclg[kNChannels];
334   for(Int_t i=0; i<kNChannels; i++){
335      hPedlg[i]->Fit("gaus","Q");
336      ADCfunclg[i] = hPedlg[i]->GetFunction("gaus");
337      MeanPed[i+kNChannels] = ADCfunclg[i]->GetParameter(1);
338      MeanPedWidth[i+kNChannels] = ADCfunclg[i]->GetParameter(2);
339      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i+kNChannels],MeanPedWidth[i+kNChannels]);
340      //printf("\t MeanPed[%d] = %f\n",i+kNChannels, MeanPed[i+kNChannels]);
341   }
342   // --- Out-of-time pedestals
343   TF1 *ADCootfunchg[kNChannels];
344   for(Int_t i=0; i<kNChannels; i++){
345      hPedOutOfTimehg[i]->Fit("gaus","Q");
346      ADCootfunchg[i] = hPedOutOfTimehg[i]->GetFunction("gaus");
347      MeanPedOOT[i] = ADCootfunchg[i]->GetParameter(1);
348      MeanPedWidthOOT[i] = ADCootfunchg[i]->GetParameter(2);
349      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i],MeanPedWidthOOT[i]);
350      //printf("\t MeanPedOOT[%d] = %f\n",i, MeanPedOOT[i]);
351   }
352   TF1 *ADCootfunclg[kNChannels];
353   for(Int_t i=0; i<kNChannels; i++){
354      hPedOutOfTimelg[i]->Fit("gaus","Q");
355      ADCootfunclg[i] = hPedOutOfTimelg[i]->GetFunction("gaus");
356      MeanPedOOT[i+kNChannels] = ADCootfunclg[i]->GetParameter(1);
357      MeanPedWidthOOT[i+kNChannels] = ADCootfunclg[i]->GetParameter(2);
358      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i+kNChannels],MeanPedWidthOOT[i+kNChannels]);
359      //printf("\t MeanPedOOT[%d] = %f\n",i+kNChannels, MeanPedOOT[i+kNChannels]);
360   }
361   //
362   // --- Correlations
363   TProfile *hPedCorrProfhg[kNChannels], *hPedCorrProflg[kNChannels];
364   TF1 *ffunchg[kNChannels], *ffunclg[kNChannels];
365   char namhist4[50];
366   for(int i=0;i<kNChannels;i++) {
367      sprintf(namhist4,"ADCHRvsOOT%d_Prof",i);
368      hPedCorrProfhg[i] = hPedCorrhg[i]->ProfileX(namhist4,-1,-1,"S");
369      hPedCorrProfhg[i]->SetName(namhist4);
370      hPedCorrProfhg[i]->Fit("pol1","Q");
371      ffunchg[i] = hPedCorrProfhg[i]->GetFunction("pol1");
372      CorrCoeff0[i] = ffunchg[i]->GetParameter(0);
373      CorrCoeff1[i] = ffunchg[i]->GetParameter(1);
374      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i],CorrCoeff1[i]);
375      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",i, CorrCoeff0[i], i, CorrCoeff1[i]);
376   }    
377   for(int i=0;i<kNChannels;i++) {
378      sprintf(namhist4,"ADCLRvsOOT%d_Prof",i);
379      hPedCorrProflg[i] = hPedCorrlg[i]->ProfileX(namhist4,-1,-1,"S");
380      hPedCorrProflg[i]->SetName(namhist4);
381      hPedCorrProflg[i]->Fit("pol1","Q");
382      ffunclg[i] = hPedCorrProflg[i]->GetFunction("pol1");
383      CorrCoeff0[i+kNChannels] = ffunclg[i]->GetParameter(0);
384      CorrCoeff1[i+kNChannels] = ffunclg[i]->GetParameter(1);
385      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i+kNChannels],CorrCoeff1[i+kNChannels]);
386      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",
387      //         i+kNChannels, CorrCoeff0[i+kNChannels], i+kNChannels, CorrCoeff1[i+kNChannels]);
388   }    
389   //                                                   
390   fclose(fileShuttle);
391   //
392   
393
394   /* write report */
395   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
396
397   /* close result file */
398   fclose(fp);
399
400
401   return status;
402 }