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