]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCEMDda.cxx
Updates for new ADD + ZDC modules
[u/mrichter/AliRoot.git] / ZDC / ZDCEMDda.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 CALIBRATION_EMD runs
13
14 Contact: Chiara.Oppedisano@to.infn.it
15 Link: 
16 Run Type: STANDALONE_EMD_RUN
17 DA Type: LDC
18 Number of events needed: at least ~5*10^3
19 Input Files: ZDCPedestal.dat
20 Output Files: ZDCEMDCalib.dat, ZDCChMapping.dat
21 Trigger Types Used: Standalone Trigger
22
23 */
24 #define PEDDATA_FILE  "ZDCPedestal.dat"
25 #define MAPDATA_FILE  "ZDCChMapping.dat"
26 #define ENCALIBDATA_FILE   "ZDCEnergyCalib.dat"
27 #define TOWCALIBDATA_FILE  "ZDCTowerCalib.dat"
28
29 #include <stdio.h>
30 #include <Riostream.h>
31 #include <Riostream.h>
32
33 // DATE
34 #include <daqDA.h>
35 #include <event.h>
36 #include <monitor.h>
37
38 //ROOT
39 #include <TROOT.h>
40 #include <TPluginManager.h>
41 #include <TH1F.h>
42 #include <TH2F.h>
43 #include <TProfile.h>
44 #include <TF1.h>
45 #include <TFile.h>
46 #include <TFitter.h>
47 #include "TMinuitMinimizer.h"
48
49 //AliRoot
50 #include <AliRawReaderDate.h>
51 #include <AliRawEventHeaderBase.h>
52 #include <AliZDCRawStream.h>
53
54
55 /* Main routine
56       Arguments: 
57       1- monitoring data source
58 */
59 int main(int argc, char **argv) {
60   
61
62   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
63                                         "*",
64                                         "TStreamerInfo",
65                                         "RIO",
66                                         "TStreamerInfo()"); 
67
68   TMinuitMinimizer m; 
69   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit","TMinuitMinimizer",
70       "Minuit", "TMinuitMinimizer(const char *)");
71   TVirtualFitter::SetDefaultFitter("Minuit");
72
73   int status = 0;
74   // No. of ZDC cabled ch.
75   int const kNModules = 10;
76   int const kNChannels = 24;
77   int const kNScChannels = 32;
78   Int_t kFirstADCGeo=0, kLastADCGeo=3;
79             
80   Int_t iMod=-1;
81   Int_t modGeo[kNModules], modType[kNModules],modNCh[kNModules];
82   for(Int_t kl=0; kl<kNModules; kl++){
83      modGeo[kl]=modType[kl]=modNCh[kl]=0;
84   }
85   
86   Int_t ich=0;
87   Int_t adcMod[2*kNChannels], adcCh[2*kNChannels], sigCode[2*kNChannels];
88   Int_t det[2*kNChannels], sec[2*kNChannels];
89   for(Int_t y=0; y<2*kNChannels; y++){
90     adcMod[y]=adcCh[y]=sigCode[y]=det[y]=sec[y]=0;
91   }
92   
93   Int_t iScCh=0;
94   Int_t scMod[kNScChannels], scCh[kNScChannels], scSigCode[kNScChannels];
95   Int_t scDet[kNScChannels], scSec[kNScChannels];
96   for(Int_t y=0; y<kNScChannels; y++){
97     scMod[y]=scCh[y]=scSigCode[y]=scDet[y]=scSec[y]=0;
98   }
99
100   /* log start of process */
101   printf("ZDC EMD program started\n");  
102
103   /* check that we got some arguments = list of files */
104   if (argc<2) {
105     printf("Wrong number of arguments\n");
106     return -1;
107   }
108   
109   // --- Preparing histos for EM dissociation spectra
110   //
111   TH1F* histoEMDRaw[4];
112   TH1F* histoEMDCorr[4];
113   //
114   char namhistr[50], namhistc[50];
115   for(Int_t i=0; i<4; i++) {
116      if(i==0){
117        sprintf(namhistr,"ZN%d-EMDRaw",i+1);
118        sprintf(namhistc,"ZN%d-EMDCorr",i+1);
119      }
120      else if(i==1){
121        sprintf(namhistr,"ZP%d-EMDRaw",i);
122        sprintf(namhistc,"ZP%d-EMDCorr",i);
123      }
124      else if(i==2){
125        sprintf(namhistr,"ZN%d-EMDRaw",i);
126        sprintf(namhistc,"ZN%d-EMDCorr",i);
127      }
128      else if(i==3){
129        sprintf(namhistr,"ZP%d-EMDRaw",i-1);
130        sprintf(namhistc,"ZP%d-EMDCorr",i-1);
131      }
132      histoEMDRaw[i] = new TH1F(namhistr,namhistr,100,0.,4000.);
133      histoEMDCorr[i] = new TH1F(namhistc,namhistc,100,0.,4000.);
134   }
135
136    // --- Preparing histos for tower inter-calibration
137   //
138   TH1F* histZNCtow[4]; TH1F* histZPCtow[4];
139   TH1F* histZNAtow[4]; TH1F* histZPAtow[4];
140   //
141   char namhistznc[50], namhistzpc[50];
142   char namhistzna[50], namhistzpa[50];
143   for(Int_t i=0; i<4; i++) {
144      sprintf(namhistznc,"ZNC-tow%d",i+1);
145      sprintf(namhistzpc,"ZPC-tow%d",i+1);
146      sprintf(namhistzna,"ZNA-tow%d",i+1);
147      sprintf(namhistzpa,"ZPA-tow%d",i+1);
148      //
149      histZNCtow[i] = new TH1F(namhistznc,namhistznc,100,0.,4000.);
150      histZPCtow[i] = new TH1F(namhistzpc,namhistzpc,100,0.,4000.);
151      histZNAtow[i] = new TH1F(namhistzna,namhistzna,100,0.,4000.);
152      histZPAtow[i] = new TH1F(namhistzpa,namhistzpa,100,0.,4000.);
153   }
154   
155   /* open result file */
156   FILE *fp=NULL;
157   fp=fopen("./result.txt","a");
158   if (fp==NULL) {
159     printf("Failed to open file\n");
160     return -1;
161   }
162   
163   FILE *mapFile4Shuttle;
164
165   // *** To analyze EMD events you MUST have a pedestal data file!!!
166   // *** -> check if a pedestal run has been analyzed
167   int read = 0;
168   read = daqDA_DB_getFile(PEDDATA_FILE,PEDDATA_FILE);
169   if(read){
170     printf("\t ERROR!!! ZDCPedestal.dat file NOT FOUND in DAQ db!!!\n");
171     return -1;
172   }
173   else printf("\t ZDCPedestal.dat file retrieved from DAQ db\n");
174   
175   FILE *filePed = fopen(PEDDATA_FILE,"r");
176   if (filePed==NULL) {
177     printf("\t ERROR!!! Can't open ZDCPedestal.dat file!!!\n");
178     return -1;
179   }
180
181   // 144 = 48 in-time + 48 out-of-time + 48 correlations
182   Float_t readValues[2][3*2*kNChannels];
183   Float_t MeanPed[2*kNChannels];
184   Float_t CorrCoeff0[2*kNChannels], CorrCoeff1[2*kNChannels];
185   // ***************************************************
186   //   Unless we have a narrow correlation to fit we
187   //    don't fit and store in-time vs. out-of-time
188   //    histograms -> mean pedstal subtracted!!!!!!
189   // ***************************************************
190   //
191   for(int jj=0; jj<6*kNChannels; jj++){
192     for(int ii=0; ii<2; ii++){
193        fscanf(filePed,"%f",&readValues[ii][jj]);
194     }
195     if(jj<2*kNChannels){
196       MeanPed[jj] = readValues[0][jj];
197       printf("\t MeanPed[%d] = %1.1f\n",jj, MeanPed[jj]);
198     }
199     else if(jj>2*kNChannels){
200       CorrCoeff0[jj-4*kNChannels] = readValues[0][jj]; 
201       CorrCoeff1[jj-4*kNChannels] = readValues[1][jj];;
202     }
203   }
204
205   /* report progress */
206   daqDA_progressReport(10);
207
208
209   /* init some counters */
210   int nevents_physics=0;
211   int nevents_total=0;
212
213   struct eventHeaderStruct *event;
214   eventTypeType eventT;
215
216   /* read the data files */
217   int n;
218   for(n=1;n<argc;n++){
219    
220     status=monitorSetDataSource( argv[n] );
221     if (status!=0) {
222       printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
223       return -1;
224     }
225
226     /* report progress */
227     /* in this example, indexed on the number of files */
228     daqDA_progressReport(10+80*n/argc);
229
230     /* read the file */
231     for(;;){
232
233       /* get next event */
234       status=monitorGetEventDynamic((void **)&event);
235       if(status==MON_ERR_EOF) break; /* end of monitoring file has been reached */
236       if(status!=0) {
237         printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
238         return -1;
239       }
240
241       /* retry if got no event */
242       if(event==NULL) {
243         break;
244       }
245       
246       // Initalize raw-data reading and decoding
247       AliRawReader *reader = new AliRawReaderDate((void*)event);
248       reader->Select("ZDC");
249       // --- Reading event header
250       //UInt_t evtype = reader->GetType();
251       //printf("\n\t ZDCEMDda -> ev. type %d\n",evtype);
252       //printf("\t ZDCEMDda -> run # %d\n",reader->GetRunNumber());
253       //
254       AliZDCRawStream *rawStreamZDC = new AliZDCRawStream(reader);
255         
256
257       /* use event - here, just write event id to result file */
258       eventT=event->eventType;
259       
260       if(eventT==START_OF_DATA){
261                         
262         rawStreamZDC->SetSODReading(kTRUE);
263         
264         // --------------------------------------------------------
265         // --- Writing ascii data file for the Shuttle preprocessor
266         mapFile4Shuttle = fopen(MAPDATA_FILE,"w");
267         if(!rawStreamZDC->Next()) printf(" \t No raw data found!! \n");
268         else{
269           while((rawStreamZDC->Next())){
270             if(rawStreamZDC->IsHeaderMapping()){ // mapping header
271                iMod++;
272                modGeo[iMod]  = rawStreamZDC->GetADCModule();
273                modType[iMod] = rawStreamZDC->GetModType();
274                modNCh[iMod]  = rawStreamZDC->GetADCNChannels();
275             }
276             if(rawStreamZDC->IsChMapping()){ 
277               if(modType[iMod]==1){ // ADC mapping ----------------------
278                 adcMod[ich]  = rawStreamZDC->GetADCModFromMap(ich);
279                 adcCh[ich]   = rawStreamZDC->GetADCChFromMap(ich);
280                 sigCode[ich] = rawStreamZDC->GetADCSignFromMap(ich);
281                 det[ich]     = rawStreamZDC->GetDetectorFromMap(ich);
282                 sec[ich]     = rawStreamZDC->GetTowerFromMap(ich);
283                 ich++;
284               }
285               else if(modType[iMod]==2){ //VME scaler mapping --------------------
286                 scMod[iScCh]     = rawStreamZDC->GetScalerModFromMap(iScCh);
287                 scCh[iScCh]      = rawStreamZDC->GetScalerChFromMap(iScCh);
288                 scSigCode[iScCh] = rawStreamZDC->GetScalerSignFromMap(iScCh);
289                 scDet[iScCh]     = rawStreamZDC->GetScDetectorFromMap(iScCh);
290                 scSec[iScCh]     = rawStreamZDC->GetScTowerFromMap(iScCh);
291                 iScCh++;
292               }
293             }
294           }
295           // Writing data on output FXS file
296           for(Int_t is=0; is<kNModules; is++){
297              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\n",
298              modGeo[is],modType[is],modNCh[is]);
299              //printf("  EMD DA -> Module mapping: geo %d type %d #ch %d\n",
300              //  modGeo[is],modType[is],modNCh[is]);
301           }
302           for(Int_t is=0; is<2*kNChannels; is++){
303              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
304                is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
305              //printf("  EMD DA -> %d ADC: mod %d ch %d, code %d det %d, sec %d\n",
306              //  is,adcMod[is],adcCh[is],sigCode[is],det[is],sec[is]);
307           }
308           for(Int_t is=0; is<kNScChannels; is++){
309              fprintf(mapFile4Shuttle,"\t%d\t%d\t%d\t%d\t%d\t%d\n",
310                is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
311              //printf("  EMD DA -> %d Scaler: mod %d ch %d, code %d det %d, sec %d\n",
312              //  is,scMod[is],scCh[is],scSigCode[is],scDet[is],scSec[is]);
313           }
314           
315         }
316         fclose(mapFile4Shuttle);
317       }// SOD event
318     
319     if(eventT==PHYSICS_EVENT){
320       // --- Reading data header
321       reader->ReadHeader();
322       const AliRawDataHeader* header = reader->GetDataHeader();
323       if(header){
324          UChar_t message = header->GetAttributes();
325          if((message & 0x70) == 0x70){ // DEDICATED EMD RUN
326             //printf("\t STANDALONE_EMD_RUN raw data found\n");
327             continue;
328          }
329          else{
330             printf("\t NO STANDALONE_EMD_RUN raw data found\n");
331             return -1;
332          }
333       }
334       else{
335          printf("\t ATTENTION! No Raw Data Header found!!!\n");
336          return -1;
337       }
338       
339       rawStreamZDC->SetSODReading(kTRUE);
340
341       if (!rawStreamZDC->Next()) printf(" \t No raw data found!! ");
342       //
343       // ----- Setting ch. mapping -----
344       for(Int_t jk=0; jk<2*kNChannels; jk++){
345         rawStreamZDC->SetMapADCMod(jk, adcMod[jk]);
346         rawStreamZDC->SetMapADCCh(jk, adcCh[jk]);
347         rawStreamZDC->SetMapADCSig(jk, sigCode[jk]);
348         rawStreamZDC->SetMapDet(jk, det[jk]);
349         rawStreamZDC->SetMapTow(jk, sec[jk]);
350       }
351       //
352       Float_t ZDCRawADC[4], ZDCCorrADC[4], ZDCCorrADCSum[4];
353       for(Int_t g=0; g<4; g++){
354            ZDCCorrADCSum[g] = 0.;
355            ZDCRawADC[g] = 0.;
356       }
357       //
358       while(rawStreamZDC->Next()){
359         Int_t det = rawStreamZDC->GetSector(0);
360         Int_t quad = rawStreamZDC->GetSector(1);
361         
362         if(rawStreamZDC->IsADCDataWord() && !(rawStreamZDC->IsUnderflow())
363              && !(rawStreamZDC->IsOverflow()) && det!=-1 && det!=3 
364              && (rawStreamZDC->GetADCGain() == 1 && // Selecting LOW RES ch.s
365              rawStreamZDC->GetADCModule()>=kFirstADCGeo && rawStreamZDC->GetADCModule()<=kLastADCGeo)){
366           
367           // Taking LOW RES channels -> ch.+kNChannels !!!!
368           Int_t DetIndex=999, PedIndex=999;
369           // Not PMRef
370           if(quad!=5){ 
371             if(det == 1){
372               DetIndex = det-1;
373               PedIndex = quad+kNChannels; 
374             }
375             else if(det==2){
376               DetIndex = det-1;
377               PedIndex = quad+5+kNChannels;
378             }
379             else if(det == 4){
380               DetIndex = det-2;
381               PedIndex = quad+12+kNChannels;
382             }
383             else if(det == 5){
384               DetIndex = det-2;
385               PedIndex = quad+17+kNChannels;
386             }
387             // Mean pedestal subtraction 
388             Float_t Pedestal = MeanPed[PedIndex];
389             // Pedestal subtraction from correlation with out-of-time signals
390             //Float_t Pedestal = CorrCoeff0[PedIndex]+CorrCoeff1[PedIndex]*MeanPedOOT[PedIndex];
391             //
392             if(DetIndex!=999 || PedIndex!=999){ 
393               //
394               ZDCRawADC[DetIndex] += (Float_t) rawStreamZDC->GetADCValue();
395               //
396               //
397               ZDCCorrADC[DetIndex] = (rawStreamZDC->GetADCValue()) - Pedestal;
398               ZDCCorrADCSum[DetIndex] += ZDCCorrADC[DetIndex];
399               //
400               /*printf("\t det %d quad %d res %d pedInd %d "
401                  "Pedestal %1.0f -> ADCCorr = %d ZDCCorrADCSum = %d\n", 
402                  det,quad,rawStreamZDC->GetADCGain(),PedIndex,Pedestal, 
403                  (Int_t) ZDCCorrADC[DetIndex],(Int_t) ZDCCorrADCSum[DetIndex]);*/
404                       
405             }
406             // Not common PM 
407             if(quad!=0){
408               Float_t corrADCval = (rawStreamZDC->GetADCValue()) - Pedestal;
409               if(det==1)      histZNCtow[quad-1]->Fill(corrADCval);
410               else if(det==2) histZPCtow[quad-1]->Fill(corrADCval);
411               else if(det==4) histZNAtow[quad-1]->Fill(corrADCval);
412               else if(det==5) histZPAtow[quad-1]->Fill(corrADCval);
413               //
414               //printf("\t det %d tow %d fill histo w. value %1.0f\n", 
415               //  det,quad,corrADCval);
416             }
417
418             if(DetIndex==999 || PedIndex==999) 
419                printf(" WARNING! Detector a/o pedestal index are WRONG!!!\n");
420  
421           }//quad!=5
422         }//IsADCDataWord()
423         
424        }
425        //
426        nevents_physics++;
427        //
428        delete reader;
429        delete rawStreamZDC;
430        //
431        for(Int_t j=0; j<4; j++){
432           histoEMDRaw[j]->Fill(ZDCRawADC[j]);
433           histoEMDCorr[j]->Fill(ZDCCorrADCSum[j]);
434        }
435     }//(if PHYSICS_EVENT)
436       
437     /* exit when last event received, no need to wait for TERM signal */
438     else if(eventT==END_OF_RUN) {
439       printf(" -> EOR event detected\n");
440       break;
441     }
442     
443     nevents_total++;
444
445    }
446
447    /* free resources */
448    free(event);
449   }
450     
451   /* Analysis of the histograms */
452   //
453   FILE *fileShuttle1 = fopen(ENCALIBDATA_FILE,"w");
454   //
455   Int_t BinMax[4]={0,0,0,0};
456   Float_t YMax[4]={0.,0.,0.,0.};
457   Int_t NBinsx[4]={0,0,0,0};
458   Float_t MeanFitVal[4]={0.,0.,0.,0.};
459   TF1 *fitfun[4];
460   for(Int_t k=0; k<4; k++){
461      if(histoEMDCorr[k]->GetEntries() == 0){
462        printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
463        return -1;
464      } 
465      //
466      BinMax[k] = histoEMDCorr[k]->GetMaximumBin();
467      if(BinMax[k]<=6){
468        printf("\n WARNING! Something wrong with det %d histo -> ending DA WITHOUT writing output\n\n", k);
469        return -1;
470      }
471      // 
472      YMax[k] = (histoEMDCorr[k]->GetXaxis())->GetXmax();
473      NBinsx[k] = (histoEMDCorr[k]->GetXaxis())->GetNbins();
474      //printf("\n\t Det%d -> BinMax = %d, ChXMax = %f\n", k+1, BinMax[k], BinMax[k]*YMax[k]/NBinsx[k]);
475      histoEMDCorr[k]->Fit("gaus","Q","",BinMax[k]*YMax[k]/NBinsx[k]*0.7,BinMax[k]*YMax[k]/NBinsx[k]*1.25);
476      fitfun[k] = histoEMDCorr[k]->GetFunction("gaus");
477      MeanFitVal[k] = (Float_t) (fitfun[k]->GetParameter(1));
478      //printf("\n\t Mean Value from gaussian fit = %f\n", MeanFitVal[k]);
479   }
480   //
481   Float_t CalibCoeff[6];     
482   //
483   for(Int_t j=0; j<6; j++){
484      if(j<4){
485        CalibCoeff[j] = MeanFitVal[j];
486        fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
487      }
488      // ZEM energy calib. coeff. = 1
489      else if(j==4 || j==5){
490        CalibCoeff[j] = 1.; 
491        fprintf(fileShuttle1,"\t%f\n",CalibCoeff[j]);
492      }
493   }
494   fclose(fileShuttle1);
495  
496   FILE *fileShuttle2 = fopen(TOWCALIBDATA_FILE,"w");
497   //Float_t meanvalznc[4], meanvalzpc[4], meanvalzna[4], meanvalzpa[4];
498   for(Int_t j=0; j<4; j++){
499      /*if(histZNCtow[j]->GetEntries() == 0){
500        printf("\n WARNING! Empty histos -> ending DA WITHOUT writing output\n\n");
501        return -1;
502      } 
503      meanvalznc[j] = histZNCtow[j]->GetMean();
504      meanvalzpc[j] = histZPCtow[j]->GetMean();
505      meanvalzna[j] = histZNAtow[j]->GetMean();
506      meanvalzpa[j] = histZPAtow[j]->GetMean();*/
507      
508      // Note -> For the moment the inter-calibration coeff. are set to 1 
509      for(Int_t k=0; k<5; k++){  
510        Float_t icoeff = 1.;
511        fprintf(fileShuttle2,"\t%f",icoeff);
512        if(k==5) fprintf(fileShuttle2,"\n");
513      }
514   }
515   //
516   /*if(meanvalznc[1]!=0 && meanvalznc[2]!=0 && meanvalznc[3]!=0 && 
517      meanvalzpc[1]!=0 && meanvalzpc[2]!=0 && meanvalzpc[3]!=0 &&
518      meanvalzna[1]!=0 && meanvalzna[2]!=0 && meanvalzna[3]!=0 &&
519      meanvalzpa[1]!=0 && meanvalzpa[2]!=0 && meanvalzpa[3]!=0){
520     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
521         1.0,meanvalznc[0]/meanvalznc[1],meanvalznc[0]/meanvalznc[2],meanvalznc[0]/meanvalznc[3]);
522     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
523         1.0,meanvalzpc[0]/meanvalzpc[1],meanvalzpc[0]/meanvalzpc[2],meanvalzpc[0]/meanvalzpc[3]);
524     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
525         1.0,meanvalzna[0]/meanvalzna[1],meanvalzpc[0]/meanvalzna[2],meanvalzpc[0]/meanvalzna[3]);
526     fprintf(fileShuttle2,"\t%f\t%f\t%f\t%f\n",
527         1.0,meanvalzpa[0]/meanvalzpa[1],meanvalzpc[0]/meanvalzpa[2],meanvalzpc[0]/meanvalzpa[3]);
528   }
529   else{
530     printf("\n Tower intercalib. coeff. CAN'T be calculated (some mean values are ZERO)!!!\n\n");
531     return -1;
532   }*/
533   fclose(fileShuttle2);
534   
535   for(Int_t ij=0; ij<4; ij++){
536     delete histoEMDRaw[ij];
537     delete histoEMDCorr[ij];
538   }
539   
540   //delete minuitFit;
541   TVirtualFitter::SetFitter(0);
542
543   /* write report */
544   fprintf(fp,"Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
545
546   /* close result file */
547   fclose(fp);
548   
549   /* report progress */
550   daqDA_progressReport(90);
551
552   /* store the result file on FES */
553   status = daqDA_FES_storeFile(MAPDATA_FILE, "MAPPING");
554   if(status){
555     printf("Failed to export file : %d\n",status);
556     return -1;
557   }
558   //
559   status = daqDA_FES_storeFile(ENCALIBDATA_FILE, "EMDENERGYCALIB");
560   if(status){
561     printf("Failed to export file : %d\n",status);
562     return -1;
563   }
564   //
565   status = daqDA_FES_storeFile(TOWCALIBDATA_FILE, "EMDTOWERCALIB");
566   if(status){
567     printf("Failed to export file : %d\n",status);
568     return -1;
569   }
570
571   /* report progress */
572   daqDA_progressReport(100);
573
574   return status;
575 }