]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCLASERda.cxx
Updates + small fixes
[u/mrichter/AliRoot.git] / ZDC / ZDCLASERda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5 It computes the average event size and populates local "./result.txt" file with the 
6 result.
7
8 The program reports about its processing progress.
9
10 Messages on stdout are exported to DAQ log system.
11
12 DA for ZDC standalone pedestal runs
13
14 Contact: Chiara.Oppedisano@to.infn.it
15 Link: 
16 Run Type: STANDALONE_LASER_RUN
17 DA Type: LDC
18 Number of events needed: no constraint (tipically ~10^3)
19 Input Files: 
20 Output Files: ZDCLaser.dat
21 Trigger Types Used: Standalone Trigger
22
23 */
24 #define PEDDATA_FILE  "ZDCPedestal.dat"
25 #define MAPDATA_FILE  "ZDCChMapping.dat"
26 #define LASDATA_FILE  "ZDCLaserCalib.dat"
27
28 #include <stdio.h>
29 #include <stdlib.h>
30 #include <Riostream.h>
31
32 // DATE
33 #include <event.h>
34 #include <monitor.h>
35 #include <daqDA.h>
36
37 //ROOT
38 #include <TRandom.h>
39 #include <TH1F.h>
40 #include <TF1.h>
41 #include <TFile.h>
42 #include <TFitter.h>
43
44 //AliRoot
45 #include <AliRawReaderDate.h>
46 #include <AliRawEventHeaderBase.h>
47 #include <AliZDCRawStream.h>
48
49
50 /* Main routine
51       Arguments: list of DATE raw data files
52 */
53 int main(int argc, char **argv) {
54   
55 //  TVirtualFitter::SetDefaultFitter("Minuit2");
56
57   int status = 0;
58   int const kNChannels = 24;
59
60   /* log start of process */
61   printf("\n ZDC LASER program started\n");  
62
63   /* check that we got some arguments = list of files */
64   if (argc<2) {
65     printf("Wrong number of arguments\n");
66     return -1;
67   }
68
69   // --- Histograms for LASER runs 
70   //     20 signal channels + 2 reference PTMs
71   //
72   TH1F::AddDirectory(0);
73   // --- Histos for reference PMTs (high gain chains)
74   TH1F *hPMRefChg = new TH1F("hPMRefChg","hPMRefChg", 100,0.,1400.);
75   TH1F *hPMRefAhg = new TH1F("hPMRefAhg","hPMRefAhg", 100,0.,1400.);
76   TH1F *hPMRefClg = new TH1F("hPMRefClg","hPMRefClg", 100,0.,4000.);
77   TH1F *hPMRefAlg = new TH1F("hPMRefAlg","hPMRefAlg", 100,0.,4000.);
78   //
79   // --- Histos for detector PMTs 
80   TH1F *hZNChg[5], *hZPChg[5], *hZNAhg[5], *hZPAhg[5], *hZEMhg[2];
81   TH1F *hZNClg[5], *hZPClg[5], *hZNAlg[5], *hZPAlg[5], *hZEMlg[2];
82   char hnamZNChg[20], hnamZPChg[20], hnamZNAhg[20], hnamZPAhg[20];
83   char hnamZNClg[20], hnamZPClg[20], hnamZNAlg[20], hnamZPAlg[20];
84   char hnamZEMhg[20], hnamZEMlg[20];
85   for(Int_t j=0; j<5; j++){
86     sprintf(hnamZNChg,"ZNChg-tow%d",j);
87     sprintf(hnamZPChg,"ZPChg-tow%d",j);
88     sprintf(hnamZNAhg,"ZNAhg-tow%d",j);
89     sprintf(hnamZPAhg,"ZPAhg-tow%d",j);
90     //
91     hZNChg[j] = new TH1F(hnamZNChg, hnamZNChg, 100, 0., 1400.);
92     hZPChg[j] = new TH1F(hnamZPChg, hnamZPChg, 100, 0., 1400.);
93     hZNAhg[j] = new TH1F(hnamZNAhg, hnamZNAhg, 100, 0., 1400.);
94     hZPAhg[j] = new TH1F(hnamZPAhg, hnamZPAhg, 100, 0., 1400.);
95     //
96     sprintf(hnamZNClg,"ZNClg-tow%d",j);
97     sprintf(hnamZPClg,"ZPClg-tow%d",j);
98     sprintf(hnamZNAlg,"ZNAlg-tow%d",j);
99     sprintf(hnamZPAlg,"ZPAlg-tow%d",j);
100     //
101     hZNClg[j] = new TH1F(hnamZNClg, hnamZNClg, 100, 0., 4000.);
102     hZPClg[j] = new TH1F(hnamZPClg, hnamZPClg, 100, 0., 4000.);
103     hZNAlg[j] = new TH1F(hnamZNAlg, hnamZNAlg, 100, 0., 4000.);
104     hZPAlg[j] = new TH1F(hnamZPAlg, hnamZPAlg, 100, 0., 4000.);
105     //
106     if(j<2){
107       sprintf(hnamZEMhg,"ZEM%dhg",j);
108       sprintf(hnamZEMlg,"ZEM%dlg",j);
109       //
110       hZEMhg[j] = new TH1F(hnamZEMhg, hnamZEMhg, 100, 0., 1400.);      
111       hZEMlg[j] = new TH1F(hnamZEMlg, hnamZEMlg, 100, 0., 4000.);      
112     }
113   }
114
115   /* open result file */
116   FILE *fp=NULL;
117   fp=fopen("./result.txt","a");
118   if (fp==NULL) {
119     printf("Failed to open file\n");
120     return -1;
121   }
122   /* report progress */
123   daqDA_progressReport(10);
124         
125   // *** To analyze LASER events you MUST have a pedestal data file!!!
126   // *** -> check if a pedestal run has been analyzed
127   int read = 0;
128   read = daqDA_DB_getFile(PEDDATA_FILE, PEDDATA_FILE);
129   if(read){
130     printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
131     return -1;
132   }
133   else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
134   
135   FILE *filePed = fopen(PEDDATA_FILE,"r");
136   if (filePed==NULL) {
137     printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
138     return -1;
139   }
140
141   // 144 = 48 in-time + 48 out-of-time + 48 correlations
142   Float_t readValues[2][6*kNChannels];
143   Float_t MeanPedhg[kNChannels], MeanPedlg[kNChannels];
144   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
145   // ***************************************************
146   //   Unless we have a narrow correlation to fit we
147   //    don't fit and store in-time vs. out-of-time
148   //    histograms -> mean pedstal subtracted!!!!!!
149   // ***************************************************
150   //
151   for(int jj=0; jj<6*kNChannels; jj++){
152     for(int ii=0; ii<2; ii++){
153        fscanf(filePed,"%f",&readValues[ii][jj]);
154     }
155     if(jj<kNChannels){
156       MeanPedhg[jj] = readValues[0][jj];
157       //printf("\t MeanPedhg[%d] = %1.1f\n",jj, MeanPedhg[jj]);
158     }
159     else if(jj>=kNChannels && jj<2*kNChannels){
160       MeanPedlg[jj-kNChannels] = readValues[0][jj];
161       //printf("\t MeanPedlg[%d] = %1.1f\n",jj-kNChannels, MeanPedlg[jj-kNChannels]);
162     }
163     else if(jj>4*kNChannels){
164       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
165       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
166     }
167   }
168   
169   FILE *mapFile4Shuttle;
170
171   /* report progress */
172   daqDA_progressReport(20);
173
174
175   /* init some counters */
176   int nevents_physics=0;
177   int nevents_total=0;
178
179   /* read the data files */
180   int n;
181   for(n=1;n<argc;n++) {
182    
183     status=monitorSetDataSource( argv[n] );
184     if (status!=0) {
185       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
186       return -1;
187     }
188
189     /* report progress */
190     /* in this example, indexed on the number of files */
191     daqDA_progressReport(20+70*n/argc);
192
193     /* read the file */
194     for(;;) {
195       struct eventHeaderStruct *event;
196       eventTypeType eventT;
197
198       /* get next event */
199       status=monitorGetEventDynamic((void **)&event);
200       if (status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
201       if (status!=0) {
202         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
203         return -1;
204       }
205
206       /* retry if got no event */
207       if (event==NULL) {
208         break;
209       }
210
211       // Initalize raw-data reading and decoding
212       AliRawReader *reader = new AliRawReaderDate((void*)event);
213       reader->Select("ZDC");
214       // --- Reading event header
215       //UInt_t evtype = reader->GetType();
216       //printf("\n\t ZDCLASERda -> ev. type %d\n",evtype);
217       //printf("\t ZDCLASERda -> run # %d\n",reader->GetRunNumber());
218       //
219       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
220         
221
222       /* use event - here, just write event id to result file */
223       eventT=event->eventType;
224       
225       Int_t ich=0;
226       Int_t adcMod[2*kNChannels], adcCh[2*kNChannels];
227       Int_t sigCode[2*kNChannels], det[2*kNChannels], sec[2*kNChannels];
228
229       if(eventT==START_OF_DATA){
230
231         rawStreamZDC->SetSODReading(kTRUE);
232                 
233         // --------------------------------------------------------
234         // --- Writing ascii data file for the Shuttle preprocessor
235         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
236         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
237         else{
238           while((rawStreamZDC->Next()) && (ich<2*kNChannels)){
239             if(rawStreamZDC->IsChMapping()){
240               adcMod[ich] = rawStreamZDC->GetADCModFromMap(ich);
241               adcCh[ich] = rawStreamZDC->GetADCChFromMap(ich);
242               sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
243               det[ich] = rawStreamZDC->GetDetectorFromMap(ich);
244               sec[ich] = rawStreamZDC->GetTowerFromMap(ich);
245               //
246               fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
247                 ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
248               //
249               //printf("ZDCPEDESTALda.cxx -> %d mod %d ch %d, code %d det %d, sec %d\n",
250               //   ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
251               //
252               ich++;
253             }
254           }
255         }
256         fclose(mapFile4Shuttle);
257       }// SOD event
258
259       if(eventT==PHYSICS_EVENT){
260         // --- Reading data header
261         reader->ReadHeader();
262         const AliRawDataHeader* header = reader->GetDataHeader();
263         if(header) {
264          UChar_t message = header->GetAttributes();
265          if(message & 0x30){ // DEDICATED LASER RUN
266             //printf("\t STANDALONE_LASER_RUN raw data found\n");
267             continue;
268          }
269          else{
270             printf("ZDCLASERda.cxx -> NO STANDALONE_LASER_RUN raw data found\n");
271             return -1;
272          }
273         }
274         else{
275            printf("\t ATTENTION! No Raw Data Header found!!!\n");
276            return -1;
277         }
278
279         rawStreamZDC->SetSODReading(kTRUE);
280
281         if (!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
282         //
283         // ----- Setting ch. mapping -----
284         for(Int_t jk=0; jk<2*kNChannels; jk++){
285           rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
286           rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
287           rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
288           rawStreamZDC->SetMapDet(jk, det[jk]);
289           rawStreamZDC->SetMapTow(jk, sec[jk]);
290         }
291         //
292         while(rawStreamZDC->Next()){
293           Int_t index=-1;
294           Int_t detector = rawStreamZDC->GetSector(0);
295           Int_t sector = rawStreamZDC->GetSector(1);
296           
297           if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
298              && !(rawStreamZDC->IsOverflow()) && detector!=-1){
299             
300             //printf("  IsADCWord %d, IsUnderflow %d, IsOverflow %d\n",
301             //  rawStreamZDC->IsADCDataWord(),rawStreamZDC->IsUnderflow(),rawStreamZDC->IsOverflow());
302  
303             if(sector!=5){ // Physics signals
304               if(detector==1) index = sector;        // *** ZNC
305               else if(detector==2) index = sector+5; // *** ZPC
306               else if(detector==3) index = sector+9; // *** ZEM
307               else if(detector==4) index = sector+12;// *** ZNA
308               else if(detector==5) index = sector+17;// *** ZPA
309             }
310             else{ // Reference PMs
311               index = (detector-1)/3+22;
312             }
313             //
314             if(index==-1) printf("ERROR in ZDCLASERda.cxx -> det %d quad %d res %d index %d ADC %d\n", 
315               detector, sector, rawStreamZDC->GetADCGain(), index, rawStreamZDC->GetADCValue());
316             
317             Float_t Pedestal=0.;
318             if(rawStreamZDC->GetADCGain()==0) Pedestal = MeanPedhg[index];
319             else if(rawStreamZDC->GetADCGain()==1) Pedestal = MeanPedlg[index];
320             //
321             Float_t CorrADC = rawStreamZDC->GetADCValue() - Pedestal;
322             //
323             //printf("\tdet %d sec %d res %d index %d ped %1.0f ADCcorr %1.0f\n", 
324             //  detector, sector, rawStreamZDC->GetADCGain(), index, Pedestal,CorrADC);
325             
326             // **** Detector PMs
327             if(sector!=5){
328               if(rawStreamZDC->GetADCGain()==0){ // --- High gain chain ---
329                 // ---- side C
330                 if(detector==1) hZNChg[sector]->Fill(CorrADC);
331                 else if(detector==2) hZPChg[sector]->Fill(CorrADC);
332                 // ---- side A
333                 else if(detector==4) hZNAhg[sector]->Fill(CorrADC);
334                 else if(detector==5) hZPAhg[sector]->Fill(CorrADC);
335                 // ---- ZEM
336                 else if(detector==3) hZEMhg[sector-1]->Fill(CorrADC);
337               }
338               else if(rawStreamZDC->GetADCGain()==1){ // --- Low gain chain ---
339                 // ---- side C
340                 if(detector==1) hZNClg[sector]->Fill(CorrADC);
341                 else if(detector==2) hZPClg[sector]->Fill(CorrADC);
342                 // ---- side A
343                 else if(detector==4) hZNAlg[sector]->Fill(CorrADC);
344                 else if(detector==5) hZPAlg[sector]->Fill(CorrADC);
345                 // ---- ZEM
346                 else if(detector==3) hZEMlg[sector-1]->Fill(CorrADC);
347               }
348             }
349             // **** Reference PMs
350             else if(sector==5){
351               if(rawStreamZDC->GetADCGain()==0){ // --- High gain chain ---
352                 // ---- PMRef chain side C
353                 if(detector==1) hPMRefChg->Fill(CorrADC);
354                 // ---- PMRef side A
355                 else if(detector==4) hPMRefAhg->Fill(CorrADC);
356               }
357               else if(rawStreamZDC->GetADCGain()==1){ // --- Low gain chain ---
358                 // ---- PMRef chain side C
359                 if(detector==1) hPMRefClg->Fill(CorrADC);
360                 // ---- PMRef side A
361                 else if(detector==4) hPMRefAlg->Fill(CorrADC);
362               }
363             }
364           }//IsADCDataWord()+NOunderflow+NOoverflow
365           //
366          }
367          //
368          nevents_physics++;
369          //
370          delete reader;
371          delete rawStreamZDC;
372
373       }//(if PHYSICS_EVENT) 
374       nevents_total++;
375
376       /* free resources */
377       free(event);
378     
379     }
380   }  
381   
382   /* Analysis of the histograms */
383   //
384   Int_t det[2*kNChannels], quad[2*kNChannels];
385   Int_t maxBin[2*kNChannels], nBin[2*kNChannels];
386   Float_t xMax[2*kNChannels], maxXval[2*kNChannels], xlow[2*kNChannels]; 
387   Float_t mean[2*kNChannels], sigma[2*kNChannels];
388   TF1 *fun[2*kNChannels];
389   
390   // ******** High gain chain ********
391   for(Int_t k=0; k<5; k++){
392     // --- ZNC
393     det[k] = 1;
394     quad[k] = k;
395     maxBin[k] = hZNChg[k]->GetMaximumBin();
396     nBin[k] = (hZNChg[k]->GetXaxis())->GetNbins();
397     xMax[k] = (hZNChg[k]->GetXaxis())->GetXmax();
398     if(nBin[k]!=0) maxXval[k] = maxBin[k]*xMax[k]/nBin[k];
399     if(maxXval[k]-150.<0.) xlow[k]=0.;
400     else xlow[k] = maxXval[k]-150.;
401     hZNChg[k]->Fit("gaus","Q","",xlow[k],maxXval[k]+150.);
402     fun[k] = hZNChg[k]->GetFunction("gaus");
403     mean[k]  = (Float_t) (fun[k]->GetParameter(1));
404     sigma[k] = (Float_t) (fun[k]->GetParameter(2));
405     // --- ZPC
406     det[k+5] = 2;
407     quad[k+5] = k;
408     maxBin[k+5] = hZPChg[k]->GetMaximumBin();
409     nBin[k+5] = (hZPChg[k]->GetXaxis())->GetNbins();
410     xMax[k+5] = (hZPChg[k]->GetXaxis())->GetXmax();
411     if(nBin[k+5]!=0) maxXval[k+5] = maxBin[k+5]*xMax[k+5]/nBin[k+5];
412     if(maxXval[k+5]-150.<0.) xlow[k+5]=0.;
413     else xlow[k+5] = maxXval[k+5]-150.;
414     hZPChg[k]->Fit("gaus","Q","",xlow[k+5],maxXval[k+5]+150.);
415     fun[k+5] = hZPChg[k]->GetFunction("gaus");
416     mean[k+5]  = (Float_t) (fun[k+5]->GetParameter(1));
417     sigma[k+5] = (Float_t) (fun[k+5]->GetParameter(2));
418     // --- ZEM1
419     if(k<2){
420       det[k+10] = 3;
421       quad[k+10] = k+1;
422       maxBin[k+10] = hZEMhg[k]->GetMaximumBin();
423       nBin[k+10] = (hZEMhg[k]->GetXaxis())->GetNbins();
424       xMax[k+10] = (hZEMhg[k]->GetXaxis())->GetXmax();
425       if(nBin[k+10]!=0) maxXval[k+10] = maxBin[k+10]*xMax[k+10]/nBin[k+10];
426       if(maxXval[k+10]-150.<0.) xlow[k+10]=0.;
427       else xlow[k+10] = maxXval[k+10]-150.;
428       hZEMhg[k]->Fit("gaus","Q","",xlow[k+10],maxXval[k+10]+150.);
429       fun[k+10] = hZEMhg[k]->GetFunction("gaus");
430       mean[k+10]  = (Float_t) (fun[k+10]->GetParameter(1));
431       sigma[k+10] = (Float_t) (fun[k+10]->GetParameter(2));
432     }
433     // --- ZNA
434     det[k+12] = 4;
435     quad[k+12] = k;
436     maxBin[k+12] = hZNAhg[k]->GetMaximumBin();
437     nBin[k+12] = (hZNAhg[k]->GetXaxis())->GetNbins();
438     xMax[k+12] = (hZNAhg[k]->GetXaxis())->GetXmax();
439     if(nBin[k+12]!=0) maxXval[k+12] = maxBin[k+12]*xMax[k+12]/nBin[k+12];
440     if(maxXval[k+12]-150.<0.) xlow[k+12]=0.;
441     else xlow[k+12] = maxXval[k+12]-150.;
442     hZNAhg[k]->Fit("gaus","Q","",xlow[k+12],maxXval[k+12]+150.);
443     fun[k+12] = hZNAhg[k]->GetFunction("gaus");
444     mean[k+12]  = (Float_t) (fun[k+12]->GetParameter(1));
445     sigma[k+12] = (Float_t) (fun[k+12]->GetParameter(2));
446     // --- ZPA
447     det[k+17] = 4;
448     quad[k+17] = 5;
449     maxBin[k+17] = hZPAhg[k]->GetMaximumBin();
450     nBin[k+17] = (hZPAhg[k]->GetXaxis())->GetNbins();
451     xMax[k+17] = (hZPAhg[k]->GetXaxis())->GetXmax();
452     if(nBin[k+17]!=0) maxXval[k+17] = maxBin[k+17]*xMax[k+17]/nBin[k+17];
453     if(maxXval[k+17]-150.<0.) xlow[k+17]=0.;
454     else xlow[k+17] = maxXval[k+17]-150.;
455     hZPAhg[k]->Fit("gaus","Q","",xlow[k+17],maxXval[k+17]+150.);
456     fun[k+17] = hZPAhg[k]->GetFunction("gaus");
457     mean[k+17]  = (Float_t) (fun[k+17]->GetParameter(1));
458     sigma[k+17] = (Float_t) (fun[k+17]->GetParameter(2));    
459   }
460   // ~~~~~~~~ PM Ref side C ~~~~~~~~
461   det[22] = 1;
462   quad[22] = 5;
463   maxBin[22] = hPMRefChg->GetMaximumBin();
464   nBin[22] = (hPMRefChg->GetXaxis())->GetNbins();
465   xMax[22] = (hPMRefChg->GetXaxis())->GetXmax();
466   if(nBin[22]!=0) maxXval[22] = maxBin[22]*xMax[22]/nBin[22];
467   if(maxXval[22]-150.<0.) xlow[22]=0.;
468   else xlow[22] = maxXval[22];
469   hPMRefChg->Fit("gaus","Q","",xlow[22],maxXval[22]+150.);
470   fun[22] = hPMRefChg->GetFunction("gaus");
471   mean[22]  = (Float_t) (fun[22]->GetParameter(1));
472   sigma[22] = (Float_t) (fun[22]->GetParameter(2));
473   // ~~~~~~~~ PM Ref side A ~~~~~~~~
474   det[23] = 4;
475   quad[23] = 5;
476   maxBin[23] = hPMRefAhg->GetMaximumBin();
477   nBin[23] = (hPMRefAhg->GetXaxis())->GetNbins();
478   xMax[23] = (hPMRefAhg->GetXaxis())->GetXmax();
479   if(nBin[23]!=0) maxXval[23] = maxBin[23]*xMax[23]/nBin[23];
480   if(maxXval[23]-100.<0.) xlow[23]=0.;
481   else xlow[23] = maxXval[23];
482   hPMRefAhg->Fit("gaus","Q","",xlow[23],maxXval[23]+100.);
483   fun[23] = hPMRefAhg->GetFunction("gaus");
484   mean[23]  = (Float_t) (fun[23]->GetParameter(1));
485   sigma[23] = (Float_t) (fun[23]->GetParameter(2));
486   
487   // ******** Low gain chain ********
488   Int_t kOffset = 24;
489   for(Int_t k=0; k<5; k++){
490     // --- ZNC
491     det[k+kOffset] = 1;
492     quad[k+kOffset] = k;
493     maxBin[k+kOffset] = hZNClg[k]->GetMaximumBin();
494     nBin[k+kOffset] = (hZNClg[k]->GetXaxis())->GetNbins();
495     xMax[k+kOffset] = (hZNClg[k]->GetXaxis())->GetXmax();
496     if(nBin[k+kOffset]!=0) maxXval[k+kOffset] = maxBin[k+kOffset]*xMax[k+kOffset]/nBin[k+kOffset];
497     if(maxXval[k+kOffset]-150.<0.) xlow[k+kOffset]=0.;
498     else xlow[k+kOffset] = maxXval[k+kOffset]-150.;
499     hZNClg[k]->Fit("gaus","Q","",xlow[k+kOffset],maxXval[k+kOffset]+150.);
500     fun[k+kOffset] = hZNClg[k]->GetFunction("gaus");
501     mean[k+kOffset]  = (Float_t) (fun[k+kOffset]->GetParameter(1));
502     sigma[k+kOffset] = (Float_t) (fun[k+kOffset]->GetParameter(2));
503     // --- ZPC
504     det[k+kOffset+5] = 2;
505     quad[k+kOffset+5] = k;
506     maxBin[k+kOffset+5] = hZPClg[k]->GetMaximumBin();
507     nBin[k+kOffset+5] = (hZPClg[k]->GetXaxis())->GetNbins();
508     xMax[k+kOffset+5] = (hZPClg[k]->GetXaxis())->GetXmax();
509     if(nBin[k+kOffset+5]!=0) maxXval[k+kOffset+5] = maxBin[k+kOffset+5]*xMax[k+kOffset+5]/nBin[k+kOffset+5];
510     if(maxXval[k+kOffset+5]-150.<0.) xlow[k+kOffset+5]=0.;
511     else xlow[k+kOffset+5] = maxXval[k+kOffset+5]-150.;
512     hZPClg[k]->Fit("gaus","Q","",xlow[k+kOffset+5],maxXval[k+kOffset+5]+150.);
513     fun[k+kOffset+5] = hZPClg[k]->GetFunction("gaus");
514     mean[k+kOffset+5]  = (Float_t) (fun[k+kOffset+5]->GetParameter(1));
515     sigma[k+kOffset+5] = (Float_t) (fun[k+kOffset+5]->GetParameter(2));
516     // --- ZEM1
517     if(k+kOffset<2){
518       det[k+kOffset+10] = 3;
519       quad[k+kOffset+10] = k+1;
520       maxBin[k+kOffset+10] = hZEMlg[k]->GetMaximumBin();
521       nBin[k+kOffset+10] = (hZEMlg[k]->GetXaxis())->GetNbins();
522       xMax[k+kOffset+10] = (hZEMlg[k]->GetXaxis())->GetXmax();
523       if(nBin[k+kOffset+10]!=0) maxXval[k+kOffset+10] = maxBin[k+kOffset+10]*xMax[k+kOffset+10]/nBin[k+kOffset+10];
524       if(maxXval[k+kOffset+10]-150.<0.) xlow[k+kOffset+10]=0.;
525       else xlow[k+kOffset+10] = maxXval[k+kOffset+10]-150.;
526       hZEMlg[k]->Fit("gaus","Q","",xlow[k+kOffset+10],maxXval[k+kOffset+10]+150.);
527       fun[k+kOffset+10] = hZEMlg[k]->GetFunction("gaus");
528       mean[k+kOffset+10]  = (Float_t) (fun[k+kOffset+10]->GetParameter(1));
529       sigma[k+kOffset+10] = (Float_t) (fun[k+kOffset+10]->GetParameter(2));
530     }
531     // --- ZNA
532     det[k+kOffset+12] = 4;
533     quad[k+kOffset+12] = k;
534     maxBin[k+kOffset+12] = hZNAlg[k]->GetMaximumBin();
535     nBin[k+kOffset+12] = (hZNAlg[k]->GetXaxis())->GetNbins();
536     xMax[k+kOffset+12] = (hZNAlg[k]->GetXaxis())->GetXmax();
537     if(nBin[k+kOffset+12]!=0) maxXval[k+kOffset+12] = maxBin[k+kOffset+12]*xMax[k+kOffset+12]/nBin[k+kOffset+12];
538     if(maxXval[k+kOffset+12]-150.<0.) xlow[k+kOffset+12]=0.;
539     else xlow[k+kOffset+12] = maxXval[k+kOffset+12]-150.;
540     hZNAlg[k]->Fit("gaus","Q","",xlow[k+kOffset+12],maxXval[k+kOffset+12]+150.);
541     fun[k+kOffset+12] = hZNAlg[k]->GetFunction("gaus");
542     mean[k+kOffset+12]  = (Float_t) (fun[k+kOffset+12]->GetParameter(1));
543     sigma[k+kOffset+12] = (Float_t) (fun[k+kOffset+12]->GetParameter(2));
544     // --- ZPA
545     det[k+kOffset+17] = 5;
546     quad[k+kOffset+17] = k;
547     maxBin[k+kOffset+17] = hZPAlg[k]->GetMaximumBin();
548     nBin[k+kOffset+17] = (hZPAlg[k]->GetXaxis())->GetNbins();
549     xMax[k+kOffset+17] = (hZPAlg[k]->GetXaxis())->GetXmax();
550     if(nBin[k+kOffset+17]!=0) maxXval[k+kOffset+17] = maxBin[k+kOffset+17]*xMax[k+kOffset+17]/nBin[k+kOffset+17];
551     if(maxXval[k+kOffset+17]-150.<0.) xlow[k+kOffset+17]=0.;
552     else xlow[k+kOffset+17] = maxXval[k+kOffset+17]-150.;
553     hZPAlg[k]->Fit("gaus","Q","",xlow[k+kOffset+17],maxXval[k+kOffset+17]+150.);
554     fun[k+kOffset+17] = hZPAlg[k]->GetFunction("gaus");
555     mean[k+kOffset+17]  = (Float_t) (fun[k+kOffset+17]->GetParameter(1));
556     sigma[k+kOffset+17] = (Float_t) (fun[k+kOffset+17]->GetParameter(2));    
557   }
558   // ~~~~~~~~ PM Ref side C ~~~~~~~~
559   det[46] = 1;
560   quad[46] = 5;
561   maxBin[46] = hPMRefClg->GetMaximumBin();
562   nBin[46] = (hPMRefClg->GetXaxis())->GetNbins();
563   xMax[46] = (hPMRefClg->GetXaxis())->GetXmax();
564   if(nBin[46]!=0) maxXval[46] = maxBin[46]*xMax[46]/nBin[46];
565   if(maxXval[46]-150.<0.) xlow[46]=0.;
566   else xlow[46] = maxXval[46];
567   hPMRefClg->Fit("gaus","Q","",xlow[46],maxXval[46]+150.);
568   fun[46] = hPMRefClg->GetFunction("gaus");
569   mean[46]  = (Float_t) (fun[46]->GetParameter(1));
570   sigma[46] = (Float_t) (fun[46]->GetParameter(2));
571   // ~~~~~~~~ PM Ref side A ~~~~~~~~
572   det[47] = 4;
573   quad[47] = 5;
574   maxBin[47] = hPMRefAlg->GetMaximumBin();
575   nBin[47] = (hPMRefAlg->GetXaxis())->GetNbins();
576   xMax[47] = (hPMRefAlg->GetXaxis())->GetXmax();
577   if(nBin[47]!=0) maxXval[47] = maxBin[47]*xMax[47]/nBin[47];
578   if(maxXval[47]-100.<0.) xlow[47]=0.;
579   else xlow[47] = maxXval[47];
580   hPMRefAlg->Fit("gaus","Q","",xlow[47],maxXval[47]+100.);
581   fun[47] = hPMRefAlg->GetFunction("gaus");
582   mean[47]  = (Float_t) (fun[47]->GetParameter(1));
583   sigma[47] = (Float_t) (fun[47]->GetParameter(2));
584     
585   FILE *fileShuttle;
586   fileShuttle = fopen(LASDATA_FILE,"w");
587   for(Int_t i=0; i<2*kNChannels; i++){
588     fprintf(fileShuttle,"\t%d\t%d\t%f\t%f\n",det[i],quad[i],mean[i], sigma[i]); 
589   }
590   //                                                   
591   fclose(fileShuttle);
592   //
593   for(Int_t j=0; j<5; j++){
594     delete hZNChg[j];
595     delete hZPChg[j];
596     delete hZNAhg[j];
597     delete hZPAhg[j];
598     delete hZNClg[j];
599     delete hZPClg[j];
600     delete hZNAlg[j];
601     delete hZPAlg[j];
602     if(j<2){
603       delete hZEMhg[j];
604       delete hZEMlg[j];
605     }
606   }
607   delete hPMRefChg;
608   delete hPMRefAhg;
609   delete hPMRefClg;
610   delete hPMRefAlg;
611
612   /* write report */
613   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
614
615   /* close result file */
616   fclose(fp);
617   
618   /* report progress */
619   daqDA_progressReport(90);
620   
621   /* store the result file on FES */
622   status = daqDA_FES_storeFile(MAPDATA_FILE,MAPDATA_FILE);
623   if(status){
624     printf("Failed to export file : %d\n",status);
625     return -1;
626   }
627   //
628   status = daqDA_FES_storeFile(LASDATA_FILE,LASDATA_FILE);
629   if(status){
630     printf("Failed to export file : %d\n",status);
631     return -1;
632   }
633
634   /* report progress */
635   daqDA_progressReport(100);
636
637   return status;
638 }