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