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