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