]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONTRKGAINda.cxx
Change Mult binning scheme
[u/mrichter/AliRoot.git] / MUON / MUONTRKGAINda.cxx
1 /*
2  Contact: Jean-Luc Charvet <jean-luc.charvet@cern.ch>
3  Link: http://aliceinfo.cern.ch/static/Offline/dimuon/muon_html/README_mchda.html
4  Reference Runs: (station 3)
5  Index          Run
6  1              109303
7  2              109304
8  3              109305
9  4              109306
10  5              109307
11  6              109308
12  7              109309
13  8              109310
14  9              109311
15  10             109312
16  11             109313
17  Run Type: CALIBRATION
18  DA Type: LDC
19  Number of events needed: 400 events for each calibration run (11)
20  Input Files:  mutrkcalibvalues and config_ldc-MTRK-S3-0 in path : /afs/cern.ch/user/j/jcharvet/public/DA_validation
21  Output Files: local dir (not persistent) -> MUONTRKGAINda.par   FXS -> run<#>_MCH_<ldc>_GAINS
22  Trigger types used:
23  */
24
25 /**************************************************************************
26  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
27  *                                                                        *
28  * Author: The ALICE Off-line Project.                                    *
29  * Contributors are mentioned in the code where appropriate.              *
30  *                                                                        *
31  * Permission to use, copy, modify and distribute this software and its   *
32  * documentation strictly for non-commercial purposes is hereby granted   *
33  * without fee, provided that the above copyright notice appears in all   *
34  * copies and that both the copyright notice and this permission notice   *
35  * appear in the supporting documentation. The authors make no claims     *
36  * about the suitability of this software for any purpose. It is          *
37  * provided "as is" without express or implied warranty.                  *
38  **************************************************************************/
39
40 /* $Id$ */
41
42 /*
43  -------------------------------------------------------------------------
44  2012-02-29 New version: MUONTRKGAINda.cxx,v 1.7
45  -------------------------------------------------------------------------
46  
47  Version for MUONTRKGAINda MUON tracking
48  (A. Baldisseri, J.-L. Charvet)
49  
50  
51  Rem:  AliMUON2DMap stores all channels, even those which are not connected
52  if pedMean == -1, channel not connected to a pad  
53  
54  
55  */
56 extern "C" {
57 #include <daqDA.h>
58 }
59
60 #include "event.h"
61 #include "monitor.h"
62
63 #include <Riostream.h>
64 #include <stdio.h>
65 #include <stdlib.h>
66 #include <sstream>
67 #include <math.h> 
68
69 //AliRoot
70 #include "AliMUONRawStreamTrackerHP.h"
71 #include "AliRawReader.h"
72 #include "AliMUONVStore.h"
73 #include "AliMUON2DMap.h"
74 #include "AliMUONCalibParamND.h"
75 #include "AliMpIntPair.h"
76 #include "AliMpConstants.h"
77 #include "AliRawDataErrorLog.h"
78 #include "AliMUONTrackerIO.h"
79 #include "AliLog.h"
80 #include "AliMUONDspHeader.h"
81 #include "AliDAQ.h"
82
83 //ROOT
84 #include "TFile.h"
85 #include "TSystem.h"
86 #include "TTree.h"
87 #include "TH1F.h"
88 #include "TString.h"
89 #include "TStopwatch.h"
90 #include "TMath.h"
91 #include "TTimeStamp.h"
92 #include "TGraphErrors.h"
93 #include "TF1.h"
94 #include "TROOT.h"
95 #include "TPluginManager.h"
96 #include "TFitter.h"
97 #include "TObjString.h"
98 #include "THashTable.h"
99 #include <THashList.h>
100 //
101 //AMORE
102 //
103 #ifdef ALI_AMORE
104 #include <AmoreDA.h>
105 #endif
106
107 #include "AliMUONPedestal.h"
108 #include "AliMUONGain.h"
109 #include "AliMUONErrorCounter.h"
110
111 // main routine
112 int main(Int_t argc, const char** argv) 
113 {
114   Int_t status=0 , status1=0 ;
115   TStopwatch timers;
116   timers.Start(kTRUE); 
117   
118   const char* prefixDA = "MUONTRKGAINda"; // program prefix
119   const char* prefixLDC = getenv("DATE_ROLE_NAME"); // LDC name
120   if(prefixLDC == NULL)  prefixLDC ="MCH" ;
121   printf("%s : -------- Begin execution : %s --------  \n",prefixLDC,prefixDA); 
122
123   // const char* prefixDA = "MUONTRKGAINda"; // program prefix
124   // printf(" ######## Begin execution : %s ######## \n\n",prefixDA); 
125   
126   TString inputFile;
127   // decode the input line
128   if (argc!=2) { printf("Wrong number of arguments\n");  return -1; }
129   
130   inputFile = argv[1];
131   
132   // needed for streamer application
133   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
134                                         "*",
135                                         "TStreamerInfo",
136                                         "RIO",
137                                         "TStreamerInfo()"); 
138   // needed for Minuit plugin
139   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
140                                         "Minuit",
141                                         "TMinuitMinimizer",
142                                         "Minuit",
143                                         "TMinuitMinimizer(const char*)");
144   Int_t skipEvents = 0;
145   Int_t maxEvents  = 1000000;
146   Int_t maxDateEvents  = 1000000;
147   Int_t errorDetail  = 0;
148   
149   Int_t nDateEvents = 0;
150   Int_t nDateRejEvents = 0;
151   Int_t nGlitchErrors= 0;
152   Int_t nParityErrors= 0;
153   Int_t nPaddingErrors= 0;
154   Int_t nTokenlostErrors= 0;
155   Int_t nEventsRecovered = 0;
156   Int_t nEvents = 0;
157   Int_t nEvthreshold = 50; //below this nb_evt the mean of the charge is not calculated and forced to 4085 (sigma)
158   Int_t statusDA = 0 ; // DA return code 
159   
160   TString logOutputFile;
161   Char_t* detail;
162   Char_t flatFile[256]="";
163   TString shuttleFile;
164   
165   UInt_t runNumber   = 0;
166   ofstream filcout;
167   Int_t nIndex = -1; 
168   
169   // For DA Gain
170   Int_t nEntries = daqDA_ECS_getTotalIteration(); // usually = 11 = Nb of calibration runs
171   Int_t nInit=0;  // = 0 all DAC values ; = 1 DAC=0 excluded (default=1)
172   Int_t nbpf1=4;  // nb of points for linear fit (default=6) 
173   Int_t printLevel  = 0;  // printout (default=0, =1 =>.ped , => .peak & .param)
174   Int_t plotLevel  = 1;  // plotout (default=1 => tree , =2 tree+Tgraph+fit)
175   Int_t nbev=0; 
176   Int_t injCharge; 
177   
178   // Reading current iteration
179   nIndex = daqDA_ECS_getCurrentIteration();
180   if(nIndex<0 || nIndex>nEntries) {printf("\n Failed: nIndex = %d\n",nIndex); return -1 ;}
181   
182   //Gain object
183   AliMUONGain* muonGain = new AliMUONGain();
184   muonGain->SetprefixDA(prefixDA);
185   muonGain->SetprefixLDC(prefixLDC);
186   muonGain->SetAliRootDataFileName(); // MUONTRKGAINda_data.root
187   muonGain->SetStatusDA(statusDA);
188  
189   // Output log file initialisations
190   sprintf(flatFile,"%s.log",prefixDA);
191   logOutputFile=flatFile;
192   AliLog::SetStreamOutput(&filcout); // Print details on logfile
193   filcout.open(logOutputFile.Data());
194   filcout<<"//=================================================" << endl;
195   filcout<<"//" << prefixLDC << "       " << prefixDA  << endl;
196   filcout<<"//=================================================" << endl;
197   filcout<<"//  * Date  : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
198
199   muonGain->SetAlifilcout(&filcout);
200   cout<<prefixLDC << " :  Date: " << muonGain->GetDate()->AsString("l") << "\n" << endl;
201
202
203
204   UShort_t manuId;  
205   UChar_t channelId;
206   UShort_t charge;
207   
208   // DAC values stored in array vDAC (reading dbfile in DETDB)
209   //   Int_t vDAC[11] = {0,200,400,800,1200,1600,2000,2500,3000,3500,4000}
210   Int_t nConfig = 1; // flag to read or not configuration ascii file in detDB
211   Int_t vDAC[11]; // DAC values
212   Char_t dbfile[256]="";
213   Int_t nEvthres;
214   Char_t line[80];
215   sprintf(dbfile,"mutrkcalibvalues");
216   status=daqDA_DB_getFile(dbfile,dbfile);
217   if(status) {printf(" Failed  : input file %s is missing, status = %d\n",dbfile,status); return -1; } 
218   
219   ifstream filein(dbfile,ios::in); Int_t k=0, kk;
220   while (k<nEntries ) { filein >> kk >> vDAC[k] ; k++; }
221   injCharge=vDAC[nIndex-1];
222   
223   filein >> line >> nInit ; cout << "mutrkcalibvalues: " << line << nInit << " " ; // = 0 all DAC values fitted ; = 1 DAC=0 excluded (default=1)
224   filein >> line >> nbpf1; cout << line << nbpf1 << " " ; // nb of points for linear fit (default=6) 
225   filein >> line >> printLevel;  cout << line << printLevel << " " ; // printout (default=0, =1 =>.ped /run, =2 => .peak & .param)
226   filein >> line >> plotLevel;   cout << line << plotLevel << " " ; // plotout (default=1 => tree , =2 tree+Tgraph+fit)
227   filein >> line >> nConfig ; cout << line << nConfig << " " ; //nConfig (default=1 => read config in DetDB, otherwise =0)
228   filein >> line >> nEvthres ; if(nEvthres !=0)nEvthreshold=nEvthres;  cout << line << nEvthreshold << " " ; // (default = 0 <=> 50) below nEvthreshold calibration not performed 
229   filein >> line >> nbev ;  if(nbev !=0){maxEvents=nbev; cout << line << maxEvents << " " ;} // Nb of events to read  (default = 0 => reading all events)
230   cout << endl;
231   
232   muonGain->SetAliPrintLevel(printLevel);
233   muonGain->SetAliPlotLevel(plotLevel);
234   muonGain->SetconfigDA(nConfig);
235   muonGain->SetnEvthreshold(nEvthreshold);
236   
237   if(nConfig)
238   {
239     // Laod configuration ascii file from DetDB
240     sprintf(dbfile,"config_%s",getenv("DATE_ROLE_NAME"));
241     status=daqDA_DB_getFile(dbfile,dbfile);
242     if(status) {printf(" Failed  : Configuration file %s is missing, status = %d\n",dbfile,status); return -1; }
243     //    else printf(" *** Copy ascii config file: %s from DetDB to working directory and reading ...*** \n",dbfile);
244     muonGain->LoadConfig(dbfile);  
245   } 
246   
247   // Rawdeader, RawStreamHP
248   AliRawReader* rawReader = AliRawReader::Create(inputFile.Data());
249   AliMUONRawStreamTrackerHP* rawStream  = new AliMUONRawStreamTrackerHP(rawReader);    
250   //  rawStream->DisableWarnings();
251   rawStream->EnabbleErrorLogger();
252   switch (errorDetail)
253     {
254     case 0: rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kLowErrorDetail); break;/// Logs minimal information in the error messages.
255     case 1: rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kMediumErrorDetail); break;/// Logs a medium level of detail in the error messages.
256     case 2: rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kHighErrorDetail); break;/// Logs a medium level of detail in the error messages.
257     default: rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kMediumErrorDetail); break;
258     }
259   
260   printf("\n%s : Reading data from file %s\n",prefixLDC,inputFile.Data());
261
262   Int_t tabTokenError[20][14];
263   for ( Int_t i=0 ; i<20 ; i++) { for ( Int_t j=0 ; j<14 ; j++) { tabTokenError[i][j]=0;}       }
264   
265   while (rawReader->NextEvent())
266   {
267     Int_t eventType = rawReader->GetType();
268     runNumber = rawReader->GetRunNumber();
269      if(nDateEvents==0)  { filcout<<"//  ---->  RUN = " << runNumber << "\n" << endl;}
270     if (nDateEvents >= maxDateEvents) break;
271     if (nEvents >= maxEvents) break;
272     if (nDateEvents>0 &&  nDateEvents % 100 == 0)       
273         cout<< prefixLDC << " :  DATE events = " << nDateEvents << "   Used events = " << nEvents << endl;
274     
275     // check shutdown condition 
276     if (daqDA_checkShutdown()) 
277       break;
278     //Skip events
279     while (skipEvents)
280     {
281       rawReader->NextEvent();
282       skipEvents--;
283     }  
284     nDateEvents++;
285     if (eventType != PHYSICS_EVENT)
286       continue; // for the moment
287     
288     const char* detail = "";
289     // First lopp over DDL's to find good events
290     // Error counters per event (counters in the decoding lib are for each DDL)
291     Bool_t eventIsErrorMessage = kFALSE;
292     int eventGlitchErrors = 0;
293     int eventParityErrors = 0;
294     int eventPaddingErrors = 0;
295     int eventTokenlostErrors = 0;
296     rawStream->First();
297     do
298       {
299         if (rawStream->IsErrorMessage()) eventIsErrorMessage = kTRUE;
300         eventGlitchErrors += rawStream->GetGlitchErrors();
301         eventParityErrors += rawStream->GetParityErrors();
302         eventPaddingErrors += rawStream->GetPaddingErrors();
303         eventTokenlostErrors += rawStream->GetTokenLostErrors();
304         if (rawStream->GetTokenLostErrors())
305           {
306             nTokenlostErrors++;
307             const AliMUONRawStreamTrackerHP::AliBlockHeader*      blkHeader  = 0x0;
308             const AliMUONRawStreamTrackerHP::AliDspHeader*        dspHeader  = 0x0;
309             Int_t nBlock = rawStream->GetBlockCount();
310             for(Int_t iBlock = 0; iBlock < nBlock ;iBlock++)
311               {
312                 blkHeader = rawStream->GetBlockHeader(iBlock);
313                 //                printf("Block %d Total length %d\n",iBlock,blkHeader->GetTotalLength());
314                 Int_t nDsp = rawStream->GetDspCount(iBlock);
315                 //                printf("Block %d DSP %d\n",iBlock,nDsp);                
316                 for(Int_t iDsp = 0; iDsp < nDsp ;iDsp++)
317                   {
318                     dspHeader =  blkHeader->GetDspHeader(iDsp);
319                     //                printf("Dsp %d Add %X\n",iDsp,dspHeader);
320                     if (dspHeader->GetErrorWord())
321                       {
322                         Int_t ddl = rawStream->GetDDL()  ; 
323                         //       Int_t ddl = AliDAQ::DdlID("MUONTRK", rawStream->GetDDL()) - 2560 ; // format 2560 + ddl
324                         Int_t frt = (dspHeader->GetErrorWord() & 0xFFFF0000) >> 16 ; // 4*4bits right shift
325                         tabTokenError[ddl][frt]++;
326                         //       printf(" DDL %d error word %X %d %d\n",ddl,dspHeader->GetErrorWord(),frt,tabTokenError[8][4]);
327                       }
328                       
329                   }
330               }
331           }
332       } while(rawStream->NextDDL()); 
333     
334     AliMUONRawStreamTrackerHP::AliBusPatch* busPatch;
335     if (!eventIsErrorMessage) 
336     {
337       // Good events (no error) -> compute pedestal for all channels
338       rawStream->First(); 
339       nEvents++;
340       muonGain->SetAliNCurrentEvents(nEvents);
341       while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
342             {
343               for(int i = 0; i < busPatch->GetLength(); ++i)
344         {
345           busPatch->GetData(i, manuId, channelId, charge);
346           muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
347         }
348             }
349     }
350     else
351       {
352         // Events with errors
353         if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
354           {
355             filcout << " ----------- Date Event recovered = " << nDateEvents <<  " ----------------" << endl;
356             // Recover parity errors -> compute pedestal for all good buspatches
357             if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
358                                         ATTR_ORBIT_BC )) 
359               {
360                 filcout <<"Event recovered -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
361                         <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
362                         <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                              
363               } 
364             else 
365               {
366                 filcout <<"Event recovered -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
367                         <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
368                         <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
369               }
370             rawStream->First();
371             nEvents++;
372             muonGain->SetAliNCurrentEvents(nEvents);
373             while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
374               {
375                 // Check the buspatch -> if error not use it in the pedestal calculation
376                 int errorCount = 0;
377                 for(int i = 0; i < busPatch->GetLength(); ++i)
378                   {
379                     if (!busPatch->IsParityOk(i)) errorCount++;
380                   }
381                 if (!errorCount) 
382                   {
383                     // Good buspatch
384                     for(int i = 0; i < busPatch->GetLength(); ++i)
385                       {
386                         busPatch->GetData(i, manuId, channelId, charge);
387                         muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
388                       }
389                   }
390                 else
391                   {
392                     char bpname[256];
393                     AliMUONErrorCounter* errorCounter;
394                     // Bad buspatch -> not used (just print)
395                     filcout<<"bpId "<<busPatch->GetBusPatchId()<<" words "<<busPatch->GetLength()
396                            <<" parity errors "<<errorCount<<endl;
397                     // Number of events where this buspatch is missing
398                     sprintf(bpname,"bp%d",busPatch->GetBusPatchId());                                           
399                     if (!(errorCounter = (AliMUONErrorCounter*) (muonGain->GetErrorBuspatchTable()->FindObject(bpname))))
400                       {
401                         // New buspatch
402                         errorCounter = new AliMUONErrorCounter(busPatch->GetBusPatchId());
403                         errorCounter->SetName(bpname);
404                         muonGain->GetErrorBuspatchTable()->Add(errorCounter);
405                       }
406                     else
407                       {
408                         // Existing buspatch
409                         errorCounter->Increment();
410                       } 
411                     // errorCounter->Print();                                           
412                   } // end of if (!errorCount)
413               } // end of while( (busPatch = (AliMUONRawStreamTrackerHP ...
414             //        nEvents++;
415             nEventsRecovered++;
416           } //end of if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
417         else
418           {
419             // Fatal errors reject the event
420             detail = Form(" ----------- Date Event rejected = %d  ----------------",nDateEvents);
421             nDateRejEvents++;
422             filcout << detail << endl;
423             if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
424                                         ATTR_ORBIT_BC )) 
425               {
426                 filcout <<"Event rejected -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
427                         <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
428                         <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                              
429               } 
430             else 
431               {
432                 filcout <<"Event rejected -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
433                         <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
434                         <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
435           
436               }
437           } // end of if (!rawStream->GetGlitchErrors() && !rawStream->GetPaddingErrors() ...
438         filcout<<"Number of errors : Glitch "<<eventGlitchErrors
439                <<" Parity "<<eventParityErrors
440                <<" Padding "<<eventPaddingErrors
441                <<" Token lost "<<eventTokenlostErrors<<endl;
442         filcout<<endl;                  
443       } // end of if (!rawStream->IsErrorMessage())
444     
445     if (eventGlitchErrors)  nGlitchErrors++;
446     if (eventParityErrors)  nParityErrors++;
447     if (eventPaddingErrors) nPaddingErrors++;
448     //if (eventTokenlostErrors) nTokenlostErrors++;
449
450   } // while (rawReader->NextEvent())
451   delete rawReader;
452   delete rawStream;
453   
454   // process and store mean peak values in .root file
455   sprintf(flatFile,"%s.par",prefixDA);
456   if(shuttleFile.IsNull())shuttleFile=flatFile;
457   injCharge=vDAC[nIndex-1];
458   muonGain->SetAliIndex(nIndex); // fIndex 
459   muonGain->SetAliInjCharge(injCharge);
460   muonGain->SetAliNEvents(nEvents);
461   muonGain->SetAliRunNumber(runNumber);
462   muonGain->MakePedStoreForGain(shuttleFile);
463   
464    // writing some counters
465   detail=Form("\n%s : Nb of DATE events           = %d",prefixLDC,nDateEvents) ;                             cout << detail; filcout << detail ;
466   detail=Form("\n%s : Nb of Glitch errors         = %d",prefixLDC,nGlitchErrors) ;                           cout << detail; filcout << detail ;
467   detail=Form("\n%s : Nb of Parity errors         = %d",prefixLDC,nParityErrors) ;                           cout << detail; filcout << detail ;
468   detail=Form("\n%s : Nb of Token lost errors     = %d",prefixLDC,nTokenlostErrors) ;                        cout << detail; filcout << detail ;
469   detail=Form("\n%s : Nb of Rejected DATE events  = %d",prefixLDC,nDateRejEvents) ;                          cout << detail; filcout << detail ;
470   detail=Form("\n%s : Nb of recovered events      = %d",prefixLDC,nEventsRecovered) ;                        cout << detail; filcout << detail ;
471   detail=Form("\n%s : Nb of events without errors = %d",prefixLDC,nEvents-nEventsRecovered) ;                cout << detail; filcout << detail ;
472   detail=Form("\n%s : Nb of used events           = %d (threshold= %d)\n\n",prefixLDC,nEvents,nEvthreshold); cout << detail; filcout << detail ;
473  
474   // Writing Token Error table
475   if(nTokenlostErrors)
476     {
477       char* detail=Form("\nWarning: Token Lost occurence \n");
478       printf("%s",detail);
479       filcout <<  detail ;
480       for ( Int_t i=0 ; i<20 ; i++) 
481         { 
482           for ( Int_t j=4 ; j<14 ; j++) 
483             { 
484               if(tabTokenError[i][j]>0)
485                 {
486                   Int_t tab=tabTokenError[i][j];
487                   Int_t frt = j/2-1;
488                   Int_t station = i/4 +1;
489                   if( j % 2 == 0)detail=Form(" in DDL= %d (station %d) and FRT%d ( Up ) => %d Token errors (address = 0x%X0000)",2560+i,station,frt,tab,j);
490                   else detail=Form(" in DDL= %d (station %d) and FRT%d (Down) => %d Token errors (address = 0x%X0000)",2560+i,station,frt,tab,j);
491                   printf("%s\n",detail);
492                   filcout <<  detail << endl;
493                 }
494             }
495         }
496     }
497
498   // Computing gain 
499   if(nIndex==nEntries)
500   {
501     muonGain->SetAliInit(nInit); // fnInit
502     muonGain->SetAliEntries(nEntries); // fnEntries
503     muonGain->SetAliNbpf1(nbpf1); // fnbpf1
504     muonGain->MakeGainStore(shuttleFile);
505     status = muonGain->GetStatusDA()  ; 
506   }
507   
508   // ouput files
509   detail=Form("%s : Root data file             : %s\n",prefixLDC,muonGain->GetRootDataFileName()); filcout << detail ;  // cout << detail;
510   detail=Form("%s : Output logfile             : %s\n",prefixLDC,logOutputFile.Data()); filcout << detail ;   // cout << detail;
511   detail=Form("%s : Gain Histo file            : %s\n",prefixLDC,muonGain->GetHistoFileName()); filcout << detail ; //  cout << detail; 
512   detail=Form("%s : Gain file (to SHUTTLE)     : %s\n",prefixLDC,shuttleFile.Data()); filcout << detail ;  //  cout << detail;
513
514  // Copying files to local DB folder defined by DAQ_DETDB_LOCAL
515   Char_t *dir;
516   unsigned int nLastVersions=50;
517   dir= getenv("DAQ_DETDB_LOCAL");
518   if(dir != NULL)  {
519     unsigned int nLastVersions=50;
520     printf("\n%s : ---  Local DataBase: %s (Max= %d) ---\n",prefixLDC,dir,nLastVersions);
521     status1 = daqDA_localDB_storeFile(logOutputFile.Data(),nLastVersions);
522
523     if(nIndex==nEntries)
524       {
525         status1 = daqDA_localDB_storeFile(muonGain->GetRootDataFileName(),nLastVersions);
526         status1 = daqDA_localDB_storeFile(muonGain->GetHistoFileName(),nLastVersions);
527         status1 = daqDA_localDB_storeFile(shuttleFile.Data(),nLastVersions);    //   if(status1)printf(" Store file : %s   status = %d\n",shuttleFile.Data(),status1);
528       }  
529   }    
530   filcout.close();
531   
532    // Transferring pedestal file to FES  (be sure that env variable DAQDALIB_PATH is set)
533   cout << endl; 
534   status1 = daqDA_FES_storeFile(shuttleFile.Data(),"GAINS");
535   if (status1) { detail=Form("%s: !!! ERROR: Failed to export calibration file : %s to FES \n",prefixLDC,shuttleFile.Data()); 
536     printf("%s",detail); filcout << detail ; status= -1; }
537   //  else { detail=Form("%s : ----  STORE calibration FILE in FES : OK ---- \n",prefixLDC); printf("%s",detail); filcout << detail ;}
538
539   if(nIndex==nEntries)
540     {
541 #ifdef ALI_AMORE  
542       std::ifstream in(shuttleFile.Data());
543       ostringstream stringout;
544       char line[1024];
545       while ( in.getline(line,1024) )
546         stringout << line << "\n";  
547       in.close();
548           
549       amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
550       TObjString gaindata(stringout.str().c_str());
551      Int_t amoreStatus = amoreDA.Send("Gains",&gaindata);
552       if ( amoreStatus )
553         {cout << prefixLDC << " :  !!! ERROR: Failed to write Gains in the AMORE database : " << amoreStatus << endl ; status=-1 ;}
554       else 
555         cout << prefixLDC << " : amoreDA.Send(Gains) ok" << endl;  
556 #else
557       cout << prefixLDC << " : Warning: MCH DA not compiled with AMORE support" << endl;
558 #endif
559     }
560
561   if(!status)printf("\n%s : -------- End execution : %s -------- (status= %d) \n",prefixLDC,prefixDA,status);
562   else { printf("\n%s : -------- %s ending in ERROR !!!! -------- (status= %d)  \n",prefixLDC,prefixDA,status);}
563
564   timers.Stop();
565   printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());
566   return status;
567 }