]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCEMDda.cxx
Just something left from v2...now clean
[u/mrichter/AliRoot.git] / ZDC / ZDCEMDda.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 Contact: Chiara.Oppedisano@to.infn.it
15 Link: /afs/cern.ch/user/c/chiarao/public/RawEMD.date
16 Run Type: STANDALONE_EMD_RUN
17 DA Type: MON
18 Number of events needed: at least ~10^3
19 Input Files: ZDCPedestal.dat
20 Output Files: ZDCEMDCalib.dat
21 Trigger Types Used: Standalone Trigger
22
23 */
24
25 #include <stdio.h>
26 #include <Riostream.h>
27
28 // DATE
29 #include <daqDA.h>
30 #include <event.h>
31 #include <monitor.h>
32
33 //ROOT
34 #include <TRandom.h>
35 #include <TH1F.h>
36 #include <TH2F.h>
37 #include <TProfile.h>
38 #include <TF1.h>
39 #include <TFile.h>
40
41 //AliRoot
42 #include <AliRawReaderDate.h>
43 #include <AliZDCRawStream.h>
44
45
46 /* Main routine
47       Arguments: 
48       1- monitoring data source
49 */
50 int main(int argc, char **argv) {
51   
52   //
53   // --- Preparing histos for EM dissociation spectra
54   //
55   TH1F* histoEMDRaw[4];
56   TH1F* histoEMDCorr[4];
57
58   char namhistr[50], namhistc[50];
59   for(Int_t i=0; i<4; i++) {
60      if(i==0){
61        sprintf(namhistr,"ZN%d-EMDRaw",i+1);
62        sprintf(namhistc,"ZN%d-EMDCorr",i+1);
63      }
64      else if(i==1){
65        sprintf(namhistr,"ZP%d-EMDRaw",i);
66        sprintf(namhistc,"ZP%d-EMDCorr",i);
67      }
68      else if(i==2){
69        sprintf(namhistr,"ZN%d-EMDRaw",i);
70        sprintf(namhistc,"ZN%d-EMDCorr",i);
71      }
72      else if(i==3){
73        sprintf(namhistr,"ZP%d-EMDRaw",i-1);
74        sprintf(namhistc,"ZP%d-EMDCorr",i-1);
75      }
76      histoEMDRaw[i] = new TH1F(namhistr,namhistr,100,0.,4000.);
77      histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
78   }
79
80   int status;
81   
82   if (argc!=2) {
83     printf("Wrong number of arguments\n");
84     return -1;
85   }
86
87
88   /* open result file */
89   FILE *fp=NULL;
90   fp=fopen("./result.txt","a");
91   if (fp==NULL) {
92     printf("Failed to open file\n");
93     return -1;
94   }
95   
96
97   /* define data source : this is argument 1 */  
98   status = monitorSetDataSource( argv[1] );
99   if(status!=0) {
100     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
101     return -1;
102   }
103
104
105   /* declare monitoring program */
106   status = monitorDeclareMp( __FILE__ );
107   if (status!=0) {
108     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
109     return -1;
110   }
111
112
113   /* define wait event timeout - 1s max */
114   monitorSetNowait();
115   monitorSetNoWaitNetworkTimeout(1000);
116   
117
118   /* log start of process */
119   printf("ZDC PEDESTAL monitoring program started\n");  
120   
121   /* init some counters */
122   int nevents_physics=0;
123   int nevents_total=0;
124
125   struct equipmentStruct *equipment;
126   int *eventEnd;
127   int *eventData;
128   int *equipmentEnd;
129   int *equipmentData;
130   int *equipmentID;
131   struct eventHeaderStruct *event;
132   eventTypeType eventT;
133   Int_t iev=0;
134   
135   /* main loop (infinite) */
136   for(;;) {
137   
138     /* check shutdown condition */
139     if (daqDA_checkShutdown()) {break;}
140     
141     /* get next event (blocking call until timeout) */
142     status=monitorGetEventDynamic((void **)&event);
143     if (status==MON_ERR_EOF) {
144       printf ("End of File detected\n");
145       break; /* end of monitoring file has been reached */
146     }
147     
148     if (status!=0) {
149       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
150       break;
151     }
152
153     /* retry if got no event */
154     if (event==NULL) {
155       continue;
156     }
157
158     iev++; 
159
160     /* use event - here, just write event id to result file */
161     eventT=event->eventType;
162     
163     if(eventT==PHYSICS_EVENT){
164       //
165       // *** To analyze EMD events you MUST have a pedestal data file!!!
166       // *** -> check if a pedestal run has been analyzied
167       FILE *filePed=NULL;
168       filePed=fopen("./ZDCPedestal.dat","r");
169       if (filePed==NULL) {
170         printf("\t ERROR!!! You MUST have a ZDCPedestal.dat file!!!\n");
171         return -1;
172       }
173       // 144 = 48 in-time + 48 out-of-time + 48 correlations
174       Float_t readValues[2][144], MeanPed[44], MeanPedWidth[44], 
175               MeanPedOOT[44], MeanPedWidthOOT[44],
176               CorrCoeff0[44], CorrCoeff1[44];
177       //
178       for(int jj=0; jj<144; jj++){
179         for(int ii=0; ii<2; ii++){
180            fscanf(filePed,"%f",&readValues[ii][jj]);
181         }
182         if(jj<48){
183           MeanPed[jj] = readValues[0][jj];
184           MeanPedWidth[jj] = readValues[1][jj];
185         }
186         else if(jj>48 && jj<96){
187           MeanPedOOT[jj-48] = readValues[0][jj];
188           MeanPedWidthOOT[jj-48] = readValues[1][jj];
189         }
190         else if(jj>144){
191           CorrCoeff0[jj-96] = readValues[0][jj]; 
192           CorrCoeff1[jj-96] = readValues[1][jj];;
193         }
194       }
195       //
196       //
197       // Initalize raw-data reading and decoding
198       AliRawReader *reader = new AliRawReaderDate((void*)event);
199       const AliRawDataHeader* header = reader->GetDataHeader();
200       if(header){
201          UChar_t message = header->GetAttributes();
202          if(message & 0x70){ // DEDICATED EMD RUN
203             printf("\t STANDALONE_EMD_RUN raw data found\n");
204             continue;
205          }
206          else{
207             printf("\t NO STANDALONE_EMD_RUN raw data found\n");
208             return -1;
209          }
210       }
211       //Commented until we won't have true Raw Data Header...
212       //else{
213       //   printf("\t ATTENTION! No Raw Data Header found!!!\n");
214       //   return -1;
215       //}
216       //
217       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
218       //
219       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
220       //
221       Float_t ZDCRawADC[4], ZDCCorrADC[4], ZDCCorrADCSum[4];
222       for(Int_t g=0; g<4; g++){
223            ZDCCorrADCSum[g] = 0.;
224            ZDCRawADC[g] = 0.;
225       }
226       //
227       while(rawStreamZDC->Next()){
228         if(rawStreamZDC->IsADCDataWord()){
229           Int_t DetIndex=999, PedIndex=999;
230           if(rawStreamZDC->GetSector(0) == 1 || rawStreamZDC->GetSector(0) == 2){
231             DetIndex = rawStreamZDC->GetSector(0)-1;
232             PedIndex = (rawStreamZDC->GetSector(0)+1)+4*rawStreamZDC->GetSector(1);
233           }
234           else if(rawStreamZDC->GetSector(0) == 4 || rawStreamZDC->GetSector(0) == 5){
235             DetIndex = rawStreamZDC->GetSector(0)-2;
236             PedIndex = (rawStreamZDC->GetSector(0)-2)+4*rawStreamZDC->GetSector(1)+24;
237           }
238           //
239           if(rawStreamZDC->GetADCGain() == 1){ //EMD -> LR ADC
240             //
241             ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
242             // Mean pedestal subtraction 
243             Float_t Pedestal = MeanPed[PedIndex];
244             // Pedestal subtraction from correlation with out-of-time signals
245             //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
246             //
247             ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
248             ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
249             //
250             /*printf("\t det %d quad %d res %d pedInd %d ADCCorr %d ZDCCorrADCSum[%d] = %d\n", 
251                rawStreamZDC->GetSector(0),rawStreamZDC->GetSector(1),
252                rawStreamZDC->GetADCGain(),PedIndex,  
253                (Int_t) (rawStreamZDC->GetADCValue() - Pedestal), DetIndex, 
254                (Int_t) ZDCCorrADCSum[DetIndex]);
255             */
256           }   
257         }//IsADCDataWord()
258          //
259        }
260        //
261        nevents_physics++;
262        //
263        for(Int_t j=0; j<4; j++){
264           histoEMDRaw[j]->Fill(ZDCRawADC[j]);
265           histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
266        }
267     }
268     
269     nevents_total++;
270
271
272     /* free resources */
273     free(event);
274     
275     /* exit when last event received, no need to wait for TERM signal */
276     if (eventT==END_OF_RUN) {
277       printf("EOR event detected\n");
278       break;
279     }
280   }
281   
282   /* Analysis of the histograms */
283   //
284   FILE *fileShuttle;
285   fileShuttle = fopen("ZDCEMDCalib.dat","w");
286   //
287   Int_t BinMax[4];
288   Float_t YMax[4];
289   Int_t NBinsx[4];
290   Float_t MeanFitVal[4];
291   TF1 *fitfun[4];
292   for(Int_t k=0; k<4; k++){
293      BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
294      YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
295      NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
296 //     printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
297      histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
298      fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
299      MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
300      //printf("\n\t Mean Value from gaussian fit = %f\n", MeanFitVal[k]);
301   }
302   //
303    Float_t CalibCoeff[6];
304   /*for(Int_t j=0; j<6; j++){
305      if(j<4) CalibCoeff[j] = MeanFitVal[j];
306      else  CalibCoeff[j] = 1.;
307      fprintf(fileShuttle,"\t%f\n",CalibCoeff[j]);
308   }
309   */
310   // --- For the moment we have sim data only for ZN1!!!
311   for(Int_t j=0; j<6; j++){
312      if(j==0) CalibCoeff[j] = MeanFitVal[j];
313      else if(j>0 && j<4) CalibCoeff[j] = CalibCoeff[0];
314      else  CalibCoeff[j] = 1.;
315      fprintf(fileShuttle,"\t%f\n",CalibCoeff[j]);
316   }
317   //                                                   
318   fclose(fileShuttle);
319   
320
321   /* write report */
322   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
323
324   /* close result file */
325   fclose(fp);
326
327
328   return status;
329 }