]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCPEDESTALda.cxx
Updated calibration object
[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 standalon pedestal runs
15 contact: Chiara.Oppedisano@cern.ch
16
17 */
18
19 #include <stdio.h>
20 #include <Riostream.h>
21
22 // DATE
23 #include <daqDA.h>
24 #include <event.h>
25 #include <monitor.h>
26
27 //ROOT
28 #include <TRandom.h>
29 #include <TH1F.h>
30 #include <TH2F.h>
31 #include <TProfile.h>
32 #include <TF1.h>
33 #include <TFile.h>
34
35 //AliRoot
36 #include <AliRawReaderDate.h>
37 #include <AliZDCRawStream.h>
38
39
40 /* Main routine
41       Arguments: 
42       1- monitoring data source
43 */
44 int main(int argc, char **argv) {
45
46   TH1F::AddDirectory(0);
47   // --- Histograms for ADC pedestals 
48   //     [22 signal channels x 2 gain chains + 2 reference PTMs]
49   //
50   TH1F *hPed[44], *hPedOutOfTime[44];
51   TH2F *hPedCorr[44];
52   //
53   char namhist1[50], namhist2[50], namhist3[50];
54   for(Int_t j=0; j<44; j++){
55      if(j<10){
56        sprintf(namhist1,"PedZN1_%d",j);
57        sprintf(namhist2,"PedZN1OutOfTime_%d",j);
58        sprintf(namhist3,"PedCorrZN1_%d",j);
59      }
60      else if(j>=10 && j<20){
61        sprintf(namhist1,"PedZP1_%d",j-10);
62        sprintf(namhist2,"PedZP1OutOfTime_%d",j-10);
63        sprintf(namhist3,"PedCorrZP1_%d",j-10);
64      }
65      else if(j>=20 && j<24){
66        sprintf(namhist1,"PedZEM_%d",j-20);
67        sprintf(namhist2,"PedZEMOutOfTime_%d",j-20);
68        sprintf(namhist3,"PedCorrZEM_%d",j-20);
69      }
70      else if(j>=24 && j<33){
71        sprintf(namhist1,"PedZN2_%d",j-24);
72        sprintf(namhist2,"PedZN2OutOfTime_%d",j-24);
73        sprintf(namhist3,"PedCorrZN2_%d",j-24);
74      }
75      else if(j>=33 && j<43){
76        sprintf(namhist1,"PedZP2_%d",j-33);
77        sprintf(namhist2,"PedZP2OutOfTime_%d",j-33);
78        sprintf(namhist3,"PedCorrZP2_%d",j-33);
79      }
80      hPed[j] = new TH1F(namhist1, namhist1, 100,0., 200.);
81      hPedOutOfTime[j] = new TH1F(namhist2, namhist2, 100,0., 200.);
82      hPedCorr[j] = new TH2F(namhist3,namhist3,100,0.,200.,100,0.,200.);
83   }
84
85   int status;
86   
87   if (argc!=2) {
88     printf("Wrong number of arguments\n");
89     return -1;
90   }
91
92
93   /* open result file */
94   FILE *fp=NULL;
95   fp=fopen("./result.txt","a");
96   if (fp==NULL) {
97     printf("Failed to open file\n");
98     return -1;
99   }
100   
101
102   /* define data source : this is argument 1 */  
103   status = monitorSetDataSource( argv[1] );
104   if(status!=0) {
105     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
106     return -1;
107   }
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
118   /* define wait event timeout - 1s max */
119   monitorSetNowait();
120   monitorSetNoWaitNetworkTimeout(1000);
121   
122
123   /* log start of process */
124   printf("ZDC PEDESTAL monitoring program started\n");  
125   
126   /* init some counters */
127   int nevents_physics=0;
128   int nevents_total=0;
129
130   struct equipmentStruct *equipment;
131   int *eventEnd;
132   int *eventData;
133   int *equipmentEnd;
134   int *equipmentData;
135   int *equipmentID;
136
137   struct eventHeaderStruct *event;
138   eventTypeType eventT;
139   Int_t iev=0;
140   
141   /* main loop (infinite) */
142   for(;;) {
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) {
161       continue;
162     }
163
164     iev++; 
165
166     /* use event - here, just write event id to result file */
167     eventT=event->eventType;
168     
169     if(eventT==PHYSICS_EVENT){
170       //fprintf(fp,"Run #%lu, event size: %lu, BC:%u, Orbit:%u, Period:%u\n",
171       //
172       // Initalize raw-data reading and decoding
173       AliRawReader *reader = new AliRawReaderDate((void*)event);
174       const AliRawDataHeader* header = reader->GetDataHeader();
175       if(header) {
176          UChar_t message = header->GetL1TriggerMessage();
177          if(message & 0x40000){ // DEDICATED PEDESTAL RUN
178             printf("\t L1 message -> PEDESTAL raw data\n");
179             continue;
180          }
181          else{
182             printf("\t L1 message -> NO PEDESTAL raw data found\n");
183             return -1;
184          }
185       }
186       //Commented until we won't have Raw Data Header...
187       /*else{
188          //printf("\t ERROR! No Raw Data Header found!!!\n");
189          //return -1;
190       }*/
191       //
192     
193       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
194       //
195       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
196       Int_t counter=0;
197       Int_t RawADC[44], RawADCoot[44];
198       for(Int_t j=0; j<44; j++){
199          RawADC[j]=0;
200          RawADCoot[j]=0;
201       }
202       while(rawStreamZDC->Next()){
203         Int_t index=-1;
204         if(rawStreamZDC->IsADCDataWord()){
205           if(rawStreamZDC->GetSector(0)==1 || rawStreamZDC->GetSector(0)==2){ // *** ZN1, ZP1
206             index = 10*(rawStreamZDC->GetSector(0)-1)+rawStreamZDC->GetSector(1)+5*rawStreamZDC->GetADCGain();
207           }
208           else if(rawStreamZDC->GetSector(0)==3){ // *** ZEM 
209             index = 10*(rawStreamZDC->GetSector(0)-1)+(rawStreamZDC->GetSector(1)-1)+2*rawStreamZDC->GetADCGain();
210           }
211           else if(rawStreamZDC->GetSector(0)==4 || rawStreamZDC->GetSector(0)==5){ // *** ZN2, ZP2
212             index = 10*(rawStreamZDC->GetSector(0)-2)+rawStreamZDC->GetSector(1)+5*rawStreamZDC->GetADCGain()+4;
213           }
214           if(counter<44){
215             hPed[index]->Fill(rawStreamZDC->GetADCValue()); 
216             RawADC[counter] = rawStreamZDC->GetADCValue();
217           }
218           else{ 
219             hPedOutOfTime[index]->Fill(rawStreamZDC->GetADCValue());
220             RawADCoot[counter-44] = rawStreamZDC->GetADCValue();
221           }
222           counter++;
223          }//IsADCDataWord()
224          //
225          if(counter == 88){ // Last ADC channel
226            for(Int_t k=0; k<44; k++){
227               hPedCorr[k]->Fill(RawADCoot[k], RawADC[k]);
228            }
229          }
230        }
231        //
232        nevents_physics++;
233
234     }
235     
236     nevents_total++;
237
238
239     /* free resources */
240     free(event);
241     
242     /* exit when last event received, no need to wait for TERM signal */
243     if (eventT==END_OF_RUN) {
244       printf("EOR event detected\n");
245       break;
246     }
247   }
248   
249   /* Analysis of the histograms */
250   //
251   FILE *fileShuttle;
252   fileShuttle = fopen("ZDCPedestal.dat","w");
253   //
254   Float_t MeanPed[44], MeanPedWidth[44], 
255         MeanPedOOT[44], MeanPedWidthOOT[44],
256         CorrCoeff0[44], CorrCoeff1[44];
257   // --- Out-of-time pedestals
258   TF1 *ADCfunc[44];
259   for(Int_t i=0; i<44; i++){
260      hPed[i]->Fit("gaus","Q");
261      ADCfunc[i] = hPed[i]->GetFunction("gaus");
262      MeanPed[i] = ADCfunc[i]->GetParameter(1);
263      MeanPedWidth[i] = ADCfunc[i]->GetParameter(2);
264      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i],MeanPedWidth[i]);
265      //printf("\t MeanPed[%d] = %f\n",i, MeanPed[i]);
266   }
267   // --- Out-of-time pedestals
268   TF1 *ADCootfunc[44];
269   for(Int_t i=0; i<44; i++){
270      hPedOutOfTime[i]->Fit("gaus","Q");
271      ADCootfunc[i] = hPedOutOfTime[i]->GetFunction("gaus");
272      MeanPedOOT[i] = ADCootfunc[i]->GetParameter(1);
273      MeanPedWidthOOT[i] = ADCootfunc[i]->GetParameter(2);
274      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i],MeanPedWidthOOT[i]);
275      //printf("\t MeanPedOOT[%d] = %f\n",i, MeanPedOOT[i]);
276   }
277   //
278   // --- Fit of correlations
279   TProfile* hPedCorrProf[44];
280   TF1 *ffunc[44];
281   char namhist4[50];
282   for(int i=0;i<44;i++) {
283      sprintf(namhist4,"ADCvsOOT%d_Prof",i);
284      hPedCorrProf[i] = hPedCorr[i]->ProfileX(namhist4,-1,-1,"S");
285      hPedCorrProf[i]->SetName(namhist4);
286      hPedCorrProf[i]->Fit("pol1","Q");
287      ffunc[i] = hPedCorrProf[i]->GetFunction("pol1");
288      CorrCoeff0[i] = ffunc[i]->GetParameter(0);
289      CorrCoeff1[i] = ffunc[i]->GetParameter(1);
290      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i],CorrCoeff1[i]);
291      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",i, CorrCoeff0[i], i, CorrCoeff1[i]);
292   }    
293   //                                                   
294   fclose(fileShuttle);
295   
296
297   /* write report */
298   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
299
300   /* close result file */
301   fclose(fp);
302
303
304   return status;
305 }