]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCPEDESTALda.cxx
616dfde696dcf8770aa72b6a55734d0d9a24c0c2
[u/mrichter/AliRoot.git] / ZDC / ZDCPEDESTALda.cxx
1 /*
2
3 This program reads the DAQ data files passed as argument using the monitoring library.
4
5 The program reports about its processing progress.
6
7 Messages on stdout are exported to DAQ log system.
8
9 DA for ZDC standalone pedestal runs
10
11 Contact: Chiara.Oppedisano@to.infn.it
12 Link: 
13 Run Type: STANDALONE_PEDESTAL_RUN
14 DA Type: LDC
15 Number of events needed: no constraint (tipically ~10^3)
16 Input Files: none
17 Output Files: ZDCPedestal.dat, ZDCChMapping.dat
18 Trigger Types Used: Standalone Trigger
19
20 */
21 #define PEDDATA_FILE  "ZDCPedestal.dat"
22 #define PEDHISTO_FILE "ZDCPedHisto.root"
23 #define MAPDATA_FILE  "ZDCChMapping.dat"
24
25 #include <stdio.h>
26 #include <stdlib.h>
27 #include <Riostream.h>
28
29 // DATE
30 #include <event.h>
31 #include <monitor.h>
32 #include <daqDA.h>
33
34 //ROOT
35 #include <TROOT.h>
36 #include <TPluginManager.h>
37 #include <TH1F.h>
38 #include <TH2F.h>
39 #include <TProfile.h>
40 #include <TF1.h>
41 #include <TFile.h>
42 #include <TFitter.h>
43 #include "TMinuitMinimizer.h"
44
45 //AliRoot
46 #include <AliRawReaderDate.h>
47 #include <AliRawEventHeaderBase.h>
48 #include <AliZDCRawStream.h>
49
50
51 /* Main routine
52       Arguments: list of DATE raw data files
53 */
54 int main(int argc, char **argv) {
55   
56
57   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
58                                         "*",
59                                         "TStreamerInfo",
60                                         "RIO",
61                                         "TStreamerInfo()"); 
62
63   TMinuitMinimizer m; 
64   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
65       "Minuit", "TMinuitMinimizer(const char *)");
66   TVirtualFitter::SetDefaultFitter("Minuit");
67
68   int status = 0;
69   int const kNChannels = 24;
70   int const kNScChannels = 32;
71   Int_t kFirstADCGeo=0, kLastADCGeo=3;
72   
73   Int_t ich=0;
74   Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
75   Int_t det[2*kNChannels], sec[2*kNChannels];
76   for(Int_t y=0; y<2*kNChannels; y++){
77     adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
78   }
79
80   /* log start of process */
81   printf("\n ZDC PEDESTAL program started\n");  
82
83   /* check that we got some arguments = list of files */
84   if (argc<2){
85     printf("Wrong number of arguments\n");
86     return -1;
87   }
88
89   // --- Histograms for ADC pedestals 
90   //     [22 signal channels + 2 reference PTMs]  x 2 gain chains
91   //
92   TH1F::AddDirectory(0);
93   //
94   TH1F *hPedhg[kNChannels], *hPedOutOfTimehg[kNChannels];
95   TH2F *hPedCorrhg[kNChannels];
96   TH1F *hPedlg[kNChannels], *hPedOutOfTimelg[kNChannels];
97   TH2F *hPedCorrlg[kNChannels];
98   //
99   char namhist1hg[50], namhist2hg[50], namhist3hg[50];
100   char namhist1lg[50], namhist2lg[50], namhist3lg[50];
101   for(Int_t j=0; j<kNChannels; j++){
102      if(j<=4){ // ZNC
103        sprintf(namhist1hg,"PedZNChg_%d",j);
104        sprintf(namhist2hg,"PedZNChgOutOfTime_%d",j);
105        sprintf(namhist3hg,"PedCorrZNChg_%d",j);
106        //
107        sprintf(namhist1lg,"PedZNClg_%d",j);
108        sprintf(namhist2lg,"PedZNClgOutOfTime_%d",j);
109        sprintf(namhist3lg,"PedCorrZNClg_%d",j);
110      }
111      else if(j>=5 && j<=9){ // ZPC
112        sprintf(namhist1hg,"PedZPChg_%d",j-5);
113        sprintf(namhist2hg,"PedZPChgOutOfTime_%d",j-5);
114        sprintf(namhist3hg,"PedCorrZPChg_%d",j-5);
115        //
116        sprintf(namhist1lg,"PedZPClg_%d",j-5);
117        sprintf(namhist2lg,"PedZPClgOutOfTime_%d",j-5);
118        sprintf(namhist3lg,"PedCorrZPClg_%d",j-5);       
119      }
120      else if(j==10 || j==11){ // ZEM
121        sprintf(namhist1hg,"PedZEMhg_%d",j-9);
122        sprintf(namhist2hg,"PedZEMhgOutOfTime_%d",j-9);
123        sprintf(namhist3hg,"PedCorrZEMhg_%d",j-9);
124        //
125        sprintf(namhist1lg,"PedZEMlg_%d",j-9);
126        sprintf(namhist2lg,"PedZEMlgOutOfTime_%d",j-9);
127        sprintf(namhist3lg,"PedCorrZEMlg_%d",j-9);
128      }
129      else if(j>=12 && j<=16){ // ZNA
130        sprintf(namhist1hg,"PedZNAhg_%d",j-12);
131        sprintf(namhist2hg,"PedZNAhgOutOfTime_%d",j-12);
132        sprintf(namhist3hg,"PedCorrZNAhg_%d",j-12);
133        //
134        sprintf(namhist1lg,"PedZNAlg_%d",j-12);
135        sprintf(namhist2lg,"PedZNAlgOutOfTime_%d",j-12);
136        sprintf(namhist3lg,"PedCorrZNAlg_%d",j-12);
137      }
138      else if(j>=17 && j<=21){ // ZPA
139        sprintf(namhist1hg,"PedZPAhg_%d",j-17);
140        sprintf(namhist2hg,"PedZPAhgOutOfTime_%d",j-17);
141        sprintf(namhist3hg,"PedCorrZPAhg_%d",j-17);
142        //
143        sprintf(namhist1lg,"PedZPAlg_%d",j-17);
144        sprintf(namhist2lg,"PedZPAlgOutOfTime_%d",j-17);
145        sprintf(namhist3lg,"PedCorrZPAlg_%d",j-17);
146      }
147      else if(j==22 || j==23){ //Reference PMs
148        sprintf(namhist1hg,"PedRefhg_%d",j-22);
149        sprintf(namhist2hg,"PedRefhgOutOfTime_%d",j-22);
150        sprintf(namhist3hg,"PedCorrRefhg_%d",j-22);
151        //
152        sprintf(namhist1lg,"PedReflg_%d",j-22);
153        sprintf(namhist2lg,"PedReflgOutOfTime_%d",j-22);
154        sprintf(namhist3lg,"PedCorrReflg_%d",j-22);
155      }
156      // --- High gain chain histos
157      hPedhg[j] = new TH1F(namhist1hg, namhist1hg, 200, 0., 200.);
158      hPedOutOfTimehg[j] = new TH1F(namhist2hg, namhist2hg, 200, 0., 200.);
159      hPedCorrhg[j] = new TH2F(namhist3hg,namhist3hg,100,0.,200.,100,0.,200.);
160      // --- Low gain chain histos
161      hPedlg[j] = new TH1F(namhist1lg, namhist1lg, 100, 0., 1000.);
162      hPedOutOfTimelg[j] = new TH1F(namhist2lg, namhist2lg, 100, 0., 1000.);
163      hPedCorrlg[j] = new TH2F(namhist3lg,namhist3lg,100,0.,1000.,100,0.,1000.);
164   }
165
166
167   /* open result file */
168   FILE *fp=NULL;
169   fp=fopen("./result.txt","w");
170   if (fp==NULL) {
171     printf("Failed to open file\n");
172     return -1;
173   }
174   
175   FILE *mapFile4Shuttle;
176
177   /* report progress */
178   daqDA_progressReport(10);
179
180
181   /* init some counters */
182   int nevents_physics=0;
183   int nevents_total=0;
184
185   struct eventHeaderStruct *event;
186   eventTypeType eventT;
187
188   /* read the data files */
189   int n;
190   for(n=1;n<argc;n++){
191    
192     status=monitorSetDataSource( argv[n] );
193     if (status!=0) {
194       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
195       return -1;
196     }
197
198     /* report progress */
199     /* in this example, indexed on the number of files */
200     daqDA_progressReport(10+80*n/argc);
201
202     /* read the file */
203     for(;;) {
204
205       /* get next event */
206       status=monitorGetEventDynamic((void **)&event);
207       if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
208       if(status!=0) {
209         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
210         return -1;
211       }
212
213       /* retry if got no event */
214       if(event==NULL) {
215         break;
216       }
217       
218       // Initalize raw-data reading and decoding
219       AliRawReader *reader = new AliRawReaderDate((void*)event);
220       reader->Select("ZDC");
221       // --- Reading event header
222       //UInt_t evtype = reader->GetType();
223       //printf("\n\t ZDCPEDESTALda -> ev. type %d\n",evtype);
224       //printf("\t ZDCPEDESTALda -> run # %d\n",reader->GetRunNumber());
225       //
226       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
227         
228
229       /* use event - here, just write event id to result file */
230       eventT=event->eventType;
231       
232       Int_t iScCh=0;
233       Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
234       Int_t scDet[kNScChannels], scSec[kNScChannels];
235       for(Int_t y=0; y<kNScChannels; y++){
236         scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
237       }
238       //
239       Int_t modNum=-1, modType=-1;
240       
241       if(eventT==START_OF_DATA){
242                 
243         rawStreamZDC->SetSODReading(kTRUE);
244         
245         // --------------------------------------------------------
246         // --- Writing ascii data file for the Shuttle preprocessor
247         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
248         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
249         else{
250           while((rawStreamZDC->Next())){
251             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
252                modNum = rawStreamZDC->GetADCModule();
253                modType = rawStreamZDC->GetModType();
254             }
255             if(rawStreamZDC->IsChMapping()){ 
256               if(modType==1){ // ADC mapping ----------------------
257                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
258                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
259                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
260                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
261                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
262                 //
263                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
264                   ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
265                 //
266                 //printf("  Mapping in DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
267                 //  ich,adcMod[ich],adcCh[ich],sigCode[ich],det[ich],sec[ich]);
268                 //
269                 ich++;
270               }
271               else if(modType==2){ //VME scaler mapping --------------------
272                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
273                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
274                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
275                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
276                 scSec[iScCh]    = rawStreamZDC->GetScTowerFromMap(iScCh);
277                 //
278                 fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
279                   iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
280                 //
281                 //printf("  Mapping in DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
282                 //  iScCh,scMod[iScCh],scCh[iScCh],scSigCode[iScCh],scDet[iScCh],scSec[iScCh]);
283                 //
284                 iScCh++;
285               }
286             }
287           }
288         }
289         fclose(mapFile4Shuttle);
290       }// SOD event
291       
292       else if(eventT==PHYSICS_EVENT){
293         // --- Reading data header
294         reader->ReadHeader();
295         const AliRawDataHeader* header = reader->GetDataHeader();
296         if(header){
297          UChar_t message = header->GetAttributes();
298          if((message & 0x20) == 0x20){ // PEDESTAL RUN
299             //printf("\t STANDALONE_PEDESTAL RUN raw data found\n");
300          }
301          else{
302             printf("ZDCPEDESTALda.cxx -> NO STANDALONE_PEDESTAL RUN raw data found\n");
303             printf("   CDH attributes: %x\n", message);
304             return -1;
305          }
306         }
307         else{
308            printf("\t ATTENTION! No Raw Data Header found!!!\n");
309            return -1;
310         }
311         
312         rawStreamZDC->SetSODReading(kTRUE);
313
314         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n"); 
315         //
316         // ----- Setting ch. mapping -----
317         for(Int_t jk=0; jk<2*kNChannels; jk++){
318           //printf("ZDCPEDESTALDA.cxx ->  ch.%d mod %d, ch %d, code %d det %d, sec %d\n",
319           //    jk,adcMod[jk],adcCh[jk],sigCode[jk],det[jk],sec[jk]);
320           rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
321           rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
322           rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
323           rawStreamZDC->SetMapDet(jk, det[jk]);
324           rawStreamZDC->SetMapTow(jk, sec[jk]);
325         }
326         //
327         Int_t iraw=0;
328         Int_t RawADChg[kNChannels], RawADCoothg[kNChannels];
329         Int_t RawADClg[kNChannels], RawADCootlg[kNChannels];
330         for(Int_t j=0; j<kNChannels; j++){
331            RawADChg[j]=0; RawADCoothg[j]=0;
332            RawADClg[j]=0; RawADCootlg[j]=0;
333         }
334         //
335         while(rawStreamZDC->Next()){
336          Int_t index=-1;
337          Int_t detector = rawStreamZDC->GetSector(0);
338          Int_t sector = rawStreamZDC->GetSector(1);
339          //
340          //printf(" rawData: det %d sec %d  value %d\n", 
341          //     detector, sector,rawStreamZDC->GetADCGain(),rawStreamZDC->GetADCValue() );
342          
343          if((rawStreamZDC->IsADCDataWord()) && (detector!=-1) &&
344             (rawStreamZDC->GetADCModule()>=kFirstADCGeo && rawStreamZDC->GetADCModule()<=kLastADCGeo)){
345           if(sector!=5){ // Physics signals
346             if(detector==1) index = sector; // *** ZNC
347             else if(detector==2) index = sector+5; // *** ZPC
348             else if(detector==3) index = sector+9; // *** ZEM
349             else if(detector==4) index = sector+12; // *** ZNA
350             else if(detector==5) index = sector+17; // *** ZPA
351           }
352           else{ // Reference PMs
353             index = (detector-1)/3+22;
354           }
355           //
356           if(index==-1) printf("ERROR in ZDCPEDESTALda.cxx -> det %d quad %d index %d \n", 
357             detector,sector,index);
358           
359            //
360            if(iraw<2*kNChannels){ // --- In-time pedestals (1st 48 raw data)
361             if(rawStreamZDC->GetADCGain()==0){ 
362               hPedhg[index]->Fill(rawStreamZDC->GetADCValue()); 
363               RawADChg[index] = rawStreamZDC->GetADCValue();
364               //
365               //printf("\t filling histo hPedhg[%d]\n",index);
366             }
367             else{
368               hPedlg[index]->Fill(rawStreamZDC->GetADCValue()); 
369               RawADClg[index] = rawStreamZDC->GetADCValue();
370               //
371               //printf("\t filling histo hPedlg[%d]\n",index);
372             }
373            }
374            else{  // --- Out-of-time pedestals
375             if(rawStreamZDC->GetADCGain()==0){
376               hPedOutOfTimehg[index]->Fill(rawStreamZDC->GetADCValue());
377               RawADCoothg[index] = rawStreamZDC->GetADCValue();
378               //
379               //printf("\t filling histo hPedOutOfTimehg[%d]\n",index);
380             }
381             else{
382               hPedOutOfTimelg[index]->Fill(rawStreamZDC->GetADCValue());
383               RawADCootlg[index] = rawStreamZDC->GetADCValue();
384               //
385               //printf("\t filling histo hPedOutOfTimelg[%d]\n",index);
386             }
387            }
388             iraw++;
389           }//IsADCDataWord()
390           //
391           if(iraw == 4*kNChannels){ // Last ADC channel -> Filling correlation histos
392             for(Int_t k=0; k<kNChannels; k++){
393               hPedCorrhg[k]->Fill(RawADCoothg[k], RawADChg[k]);
394               hPedCorrlg[k]->Fill(RawADCootlg[k], RawADClg[k]);
395             }
396           }
397         }
398         nevents_physics++;
399         //
400         delete reader;
401         delete rawStreamZDC;
402
403       }//(if PHYSICS_EVENT) 
404       
405       /* exit when last event received, no need to wait for TERM signal */
406       else if(eventT==END_OF_RUN) {
407         printf(" -> EOR event detected\n");
408         break;
409       }
410       
411       nevents_total++;
412     
413     }
414           
415     /* free resources */
416     free(event);
417   }  
418   
419   /* Analysis of the histograms */
420   //
421   FILE *fileShuttle;
422   fileShuttle = fopen(PEDDATA_FILE,"w");
423   //
424   Float_t MeanPed[2*kNChannels], MeanPedWidth[2*kNChannels], 
425         MeanPedOOT[2*kNChannels], MeanPedWidthOOT[2*kNChannels];
426   // --- In-time pedestals
427   TF1 *ADCfunchg[kNChannels];
428   for(Int_t i=0; i<kNChannels; i++){
429      if(hPedhg[i]->GetEntries() == 0){
430        printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
431        return -1;
432      } 
433      hPedhg[i]->Fit("gaus","Q");
434      ADCfunchg[i] = hPedhg[i]->GetFunction("gaus");
435      MeanPed[i] = (Double_t) ADCfunchg[i]->GetParameter(1);
436      MeanPedWidth[i] = (Double_t)  ADCfunchg[i]->GetParameter(2);
437      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i],MeanPedWidth[i]);
438      //printf("\t MeanPedhg[%d] = %f\n",i, MeanPed[i]);
439   }
440   TF1 *ADCfunclg[kNChannels];
441   for(Int_t i=0; i<kNChannels; i++){
442      hPedlg[i]->Fit("gaus","Q");
443      ADCfunclg[i] = hPedlg[i]->GetFunction("gaus");
444      MeanPed[i+kNChannels] = (Double_t)  ADCfunclg[i]->GetParameter(1);
445      MeanPedWidth[i+kNChannels] = (Double_t)  ADCfunclg[i]->GetParameter(2);
446      fprintf(fileShuttle,"\t%f\t%f\n",MeanPed[i+kNChannels],MeanPedWidth[i+kNChannels]);
447      //printf("\t MeanPedlg[%d] = %f\n",i+kNChannels, MeanPed[i+kNChannels]);
448   }
449   // --- Out-of-time pedestals
450   TF1 *ADCootfunchg[kNChannels];
451   for(Int_t i=0; i<kNChannels; i++){
452      if(hPedOutOfTimehg[i]->GetEntries() == 0){
453        printf("\n WARNING! Empty histos for out-of-time channels!!!\n\n");
454        return -1;
455      } 
456      hPedOutOfTimehg[i]->Fit("gaus","Q");
457      ADCootfunchg[i] = hPedOutOfTimehg[i]->GetFunction("gaus");
458      MeanPedOOT[i] = (Double_t)  ADCootfunchg[i]->GetParameter(1);
459      MeanPedWidthOOT[i] = (Double_t)  ADCootfunchg[i]->GetParameter(2);
460      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i],MeanPedWidthOOT[i]);
461      //printf("\t MeanPedOOThg[%d] = %f\n",i, MeanPedOOT[i]);
462   }
463   TF1 *ADCootfunclg[kNChannels];
464   for(Int_t i=0; i<kNChannels; i++){
465      hPedOutOfTimelg[i]->Fit("gaus","Q");
466      ADCootfunclg[i] = hPedOutOfTimelg[i]->GetFunction("gaus");
467      MeanPedOOT[i+kNChannels] = (Double_t)  ADCootfunclg[i]->GetParameter(1);
468      MeanPedWidthOOT[i+kNChannels] = (Double_t)  ADCootfunclg[i]->GetParameter(2);
469      fprintf(fileShuttle,"\t%f\t%f\n",MeanPedOOT[i+kNChannels],MeanPedWidthOOT[i+kNChannels]);
470      //printf("\t MeanPedOOTlg[%d] = %f\n",i+kNChannels, MeanPedOOT[i+kNChannels]);
471   }
472   
473   // --- Correlations
474   // NB -> The correlations are NOT fitted since at the moment
475   // (Sptember 2009) they are NOT correlations and the DA would fail!!!
476 /*  Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
477   TProfile *hPedCorrProfhg[kNChannels], *hPedCorrProflg[kNChannels];
478   TF1 *ffunchg[kNChannels], *ffunclg[kNChannels];
479   char namhist4[50];
480   for(int i=0;i<kNChannels;i++) {
481      sprintf(namhist4,"ADCHRvsOOT%d_Prof",i);
482      hPedCorrProfhg[i] = hPedCorrhg[i]->ProfileX(namhist4,-1,-1,"S");
483      hPedCorrProfhg[i]->SetName(namhist4);
484      hPedCorrProfhg[i]->Fit("pol1","Q");
485      ffunchg[i] = hPedCorrProfhg[i]->GetFunction("pol1");
486      CorrCoeff0[i] = (Double_t)  ffunchg[i]->GetParameter(0);
487      CorrCoeff1[i] = (Double_t) ffunchg[i]->GetParameter(1);
488      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i],CorrCoeff1[i]);
489      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",i, CorrCoeff0[i], i, CorrCoeff1[i]);
490   }    
491   for(int i=0;i<kNChannels;i++) {
492      sprintf(namhist4,"ADCLRvsOOT%d_Prof",i);
493      hPedCorrProflg[i] = hPedCorrlg[i]->ProfileX(namhist4,-1,-1,"S");
494      hPedCorrProflg[i]->SetName(namhist4);
495      hPedCorrProflg[i]->Fit("pol1","Q");
496      ffunclg[i] = hPedCorrProflg[i]->GetFunction("pol1");
497      CorrCoeff0[i+kNChannels] =  (Double_t) ffunclg[i]->GetParameter(0);
498      CorrCoeff1[i+kNChannels] =  (Double_t) ffunclg[i]->GetParameter(1);
499      fprintf(fileShuttle,"\t%f\t%f\n",CorrCoeff0[i+kNChannels],CorrCoeff1[i+kNChannels]);
500      //printf("\t CorrCoeff0[%d] = %f, CorrCoeff1[%d] = %f\n",
501      //         i+kNChannels, CorrCoeff0[i+kNChannels], i+kNChannels, CorrCoeff1[i+kNChannels]);
502   }    
503 */
504   //                                                   
505   fclose(fileShuttle);
506   //
507   /* report progress */
508   daqDA_progressReport(80);
509   //
510   TFile *histofile = new TFile(PEDHISTO_FILE,"RECREATE");
511   histofile->cd();
512   for(int k=0; k<kNChannels; k++){
513      hPedhg[k]->Write();
514      hPedOutOfTimehg[k]->Write();
515      hPedCorrhg[k]->Write();
516      hPedlg[k]->Write();
517      hPedOutOfTimelg[k]->Write();
518      hPedCorrlg[k]->Write();
519   }
520   //                                                   
521   histofile->Close();
522   //
523   for(Int_t j=0; j<kNChannels; j++){
524      delete hPedhg[j];
525      delete hPedOutOfTimehg[j];
526      delete hPedCorrhg[j];
527      delete hPedlg[j];
528      delete hPedOutOfTimelg[j];
529      delete hPedCorrlg[j];
530   }
531
532   /* write report */
533   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
534
535   /* close result file */
536   fclose(fp);
537   
538   /* report progress */
539   daqDA_progressReport(90);
540
541   /* store the result files on FES */
542   // [1] File with mapping
543   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
544   if(status){
545     printf("Failed to export mapping data file to DAQ FES\n");
546     return -1;
547   }
548   // [2] File with pedestal data
549   status = daqDA_FES_storeFile(PEDDATA_FILE, "PEDESTALDATA");
550   if(status){
551     printf("Failed to export pedestal data file to DAQ FES\n");
552     return -1;
553   }
554   // [3] File with pedestal histos
555   status = daqDA_FES_storeFile(PEDHISTO_FILE, "PEDESTALHISTOS");
556   if(status){
557     printf("Failed to export pedestal histos file to DAQ FES\n");
558     return -1;
559   }
560   
561   /* store the result files on DB */
562   status = daqDA_DB_storeFile(PEDDATA_FILE, PEDDATA_FILE);  
563   if(status){
564     printf("Failed to store pedestal data file to DAQ DB\n");
565     return -1;
566   }
567
568   /* report progress */
569   daqDA_progressReport(100);
570
571   return status;
572 }