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