]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONTRKGAINda.cxx
AliCaloPID: Correct matching rejection in case of recalculation in the analysis,...
[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  2010-04-18 New version: MUONTRKGAINda.cxx,v 1.6
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;
115   TStopwatch timers;
116   timers.Start(kTRUE); 
117   
118   const char* prefixDA = "MUONTRKGAINda"; // program prefix
119   printf(" ######## Begin execution : %s ######## \n\n",prefixDA); 
120   
121   TString inputFile;
122   // decode the input line
123   if (argc!=2) { printf("Wrong number of arguments\n");  return -1; }
124   
125   inputFile = argv[1];
126   
127   // needed for streamer application
128   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
129                                         "*",
130                                         "TStreamerInfo",
131                                         "RIO",
132                                         "TStreamerInfo()"); 
133   // needed for Minuit plugin
134   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
135                                         "Minuit",
136                                         "TMinuitMinimizer",
137                                         "Minuit",
138                                         "TMinuitMinimizer(const char*)");
139   Int_t skipEvents = 0;
140   Int_t maxEvents  = 1000000;
141   Int_t maxDateEvents  = 1000000;
142   
143   Int_t nDateEvents = 0;
144   Int_t nGlitchErrors= 0;
145   Int_t nParityErrors= 0;
146   Int_t nPaddingErrors= 0;
147   Int_t nTokenlostErrors= 0;
148   Int_t nEventsRecovered = 0;
149   Int_t nEvents = 0;
150   Int_t nEvthreshold = 10; //below this nb_evt the mean of the charge is not calculated and forced to 4085 (sigma)
151   
152   TString logOutputFile;
153   Char_t flatFile[256]="";
154   TString shuttleFile;
155   
156   UInt_t runNumber   = 0;
157   ofstream filcout;
158   Int_t nIndex = -1; 
159   
160   // For DA Gain
161   Int_t nEntries = daqDA_ECS_getTotalIteration(); // usually = 11 = Nb of calibration runs
162   Int_t nInit=1;  // = 0 all DAC values ; = 1 DAC=0 excluded (default=1)
163   Int_t nbpf1=6;  // nb of points for linear fit (default=6) 
164   Int_t printLevel  = 0;  // printout (default=0, =1 =>.ped , => .peak & .param)
165   Int_t plotLevel  = 1;  // plotout (default=1 => tree , =2 tree+Tgraph+fit)
166   Int_t nbev=0; 
167   Int_t injCharge; 
168   
169   // Reading current iteration
170   nIndex = daqDA_ECS_getCurrentIteration();
171   if(nIndex<0 || nIndex>nEntries) {printf("\n Failed: nIndex = %d\n",nIndex); return -1 ;}
172   
173   //Gain object
174   AliMUONGain* muonGain = new AliMUONGain();
175   muonGain->SetprefixDA(prefixDA);
176   muonGain->SetAliRootDataFileName(); // MUONTRKGAINda_data.root
177   
178   UShort_t manuId;  
179   UChar_t channelId;
180   UShort_t charge;
181   
182   // DAC values stored in array vDAC (reading dbfile in DETDB)
183   //   Int_t vDAC[11] = {0,200,400,800,1200,1600,2000,2500,3000,3500,4000}
184   Int_t nConfig = 1; // flag to read or not configuration ascii file in detDB
185   Int_t vDAC[11]; // DAC values
186   Char_t dbfile[256]="";
187   sprintf(dbfile,"mutrkcalibvalues");
188   status=daqDA_DB_getFile(dbfile,dbfile);
189   if(status) {printf(" Failed  : input file %s is missing, status = %d\n",dbfile,status); return -1; } 
190   
191   ifstream filein(dbfile,ios::in); Int_t k=0, kk;
192   while (k<nEntries ) { filein >> kk >> vDAC[k] ; k++; }
193   injCharge=vDAC[nIndex-1];
194   
195   filein >> nInit; // = 0 all DAC values fitted ; = 1 DAC=0 excluded (default=1)
196   filein >> nbpf1; // nb of points for linear fit (default=6) 
197   filein >> printLevel;  // printout (default=0, =1 =>.ped /run, =2 => .peak & .param)
198   filein >> plotLevel;  // plotout (default=1 => tree , =2 tree+Tgraph+fit)
199   filein >> nConfig;  //nConfig (default=1 => read config in DetDB, otherwise =0)
200   filein >> nbev;  // Nb of events to read  (default = 0 => reading all events)
201   if(nbev>0)maxEvents=nbev;
202   
203   //  printf(" *** Copy: %s from DetDB to working directory  ***      Config= %d\n",dbfile,nConfig);
204   printf(" Input parameters:  nInit= %d   Nb linear pts= %d   Print level= %d   Plot Level= %d    nConfig= %d",nInit,nbpf1,printLevel,plotLevel,nConfig);
205   if(nbev==0)printf("\n");
206   else printf("  Nb_max evt = %d\n",maxEvents);
207   
208   muonGain->SetAliPrintLevel(printLevel);
209   muonGain->SetAliPlotLevel(plotLevel);
210   muonGain->SetconfigDA(nConfig);
211   muonGain->SetnEvthreshold(nEvthreshold);
212   
213   if(nConfig)
214   {
215     // Laod configuration ascii file from DetDB
216     sprintf(dbfile,"config_%s",getenv("DATE_ROLE_NAME"));
217     status=daqDA_DB_getFile(dbfile,dbfile);
218     if(status) {printf(" Failed  : Configuration file %s is missing, status = %d\n",dbfile,status); return -1; }
219     //    else printf(" *** Copy ascii config file: %s from DetDB to working directory and reading ...*** \n",dbfile);
220     muonGain->LoadConfig(dbfile);  
221   } 
222   
223   // Rawdeader, RawStreamHP
224   AliRawReader* rawReader = AliRawReader::Create(inputFile.Data());
225   AliMUONRawStreamTrackerHP* rawStream  = new AliMUONRawStreamTrackerHP(rawReader);    
226   //  rawStream->DisableWarnings();
227   rawStream->EnabbleErrorLogger();
228   //
229   // kLowErrorDetail,     /// Logs minimal information in the error messages.
230   // kMediumErrorDetail,  /// Logs a medium level of detail in the error messages.
231   // kHighErrorDetail     /// Logs maximum information in the error messages.
232   //  rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kLowErrorDetail);
233      rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kMediumErrorDetail);
234   //   rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kHighErrorDetail);
235   
236   cout << "\n" << prefixDA << " : Reading data from file " << inputFile.Data()  << endl;
237
238   Int_t tabTokenError[20][14];
239   for ( Int_t i=0 ; i<20 ; i++) { for ( Int_t j=0 ; j<14 ; j++) { tabTokenError[i][j]=0;}       }
240   
241   while (rawReader->NextEvent())
242   {
243     if (nDateEvents >= maxDateEvents) break;
244     if (nEvents >= maxEvents) break;
245     if (nDateEvents>0 &&  nDateEvents % 100 == 0)       
246       cout<<"Cumulated:  DATE events = " << nDateEvents << "   Used events = " << nEvents << endl;
247     
248     // check shutdown condition 
249     if (daqDA_checkShutdown()) 
250       break;
251     //Skip events
252     while (skipEvents)
253     {
254       rawReader->NextEvent();
255       skipEvents--;
256     }  
257     Int_t eventType = rawReader->GetType();
258     runNumber = rawReader->GetRunNumber();
259     
260     // Output log file initialisations
261     if(nDateEvents==0)
262     {
263       sprintf(flatFile,"%s.log",prefixDA);
264       logOutputFile=flatFile;
265                 AliLog::SetStreamOutput(&filcout); // Print details on logfile      
266       filcout.open(logOutputFile.Data());
267       filcout<<"//=================================================" << endl;
268       filcout<<"//       " << prefixDA << " for run = " << runNumber << "  (DAC=" << injCharge << ")" << endl;
269       filcout<<"//=================================================" << endl;
270       filcout<<"//   * Date          : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
271       
272       cout<<"\n ********  " << prefixDA << " for run = " << runNumber << "  (Index= " << nIndex << "/" << nEntries << "  DAC=" << injCharge << ") ********\n" << endl;
273       cout<<" * Date : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
274     }
275     
276     muonGain->SetAlifilcout(&filcout);
277     
278     nDateEvents++;
279     if (eventType != PHYSICS_EVENT)
280       continue; // for the moment
281     
282     const char* detail = "";
283     // First lopp over DDL's to find good events
284     // Error counters per event (counters in the decoding lib are for each DDL)
285     Bool_t eventIsErrorMessage = kFALSE;
286     int eventGlitchErrors = 0;
287     int eventParityErrors = 0;
288     int eventPaddingErrors = 0;
289     int eventTokenlostErrors = 0;
290     rawStream->First();
291     do
292       {
293         if (rawStream->IsErrorMessage()) eventIsErrorMessage = kTRUE;
294         eventGlitchErrors += rawStream->GetGlitchErrors();
295         eventParityErrors += rawStream->GetParityErrors();
296         eventPaddingErrors += rawStream->GetPaddingErrors();
297         eventTokenlostErrors += rawStream->GetTokenLostErrors();
298         if (rawStream->GetTokenLostErrors())
299           {
300             nTokenlostErrors++;
301             const AliMUONRawStreamTrackerHP::AliBlockHeader*      blkHeader  = 0x0;
302             const AliMUONRawStreamTrackerHP::AliDspHeader*        dspHeader  = 0x0;
303             Int_t nBlock = rawStream->GetBlockCount();
304             for(Int_t iBlock = 0; iBlock < nBlock ;iBlock++)
305               {
306                 blkHeader = rawStream->GetBlockHeader(iBlock);
307                 //                printf("Block %d Total length %d\n",iBlock,blkHeader->GetTotalLength());
308                 Int_t nDsp = rawStream->GetDspCount(iBlock);
309                 //                printf("Block %d DSP %d\n",iBlock,nDsp);                
310                 for(Int_t iDsp = 0; iDsp < nDsp ;iDsp++)
311                   {
312                     dspHeader =  blkHeader->GetDspHeader(iDsp);
313                     //                printf("Dsp %d Add %X\n",iDsp,dspHeader);
314                     if (dspHeader->GetErrorWord())
315                       {
316                         Int_t ddl = rawStream->GetDDL()  ; 
317                         //       Int_t ddl = AliDAQ::DdlID("MUONTRK", rawStream->GetDDL()) - 2560 ; // format 2560 + ddl
318                         Int_t frt = (dspHeader->GetErrorWord() & 0xFFFF0000) >> 16 ; // 4*4bits right shift
319                         tabTokenError[ddl][frt]++;
320                         //       printf(" DDL %d error word %X %d %d\n",ddl,dspHeader->GetErrorWord(),frt,tabTokenError[8][4]);
321                       }
322                       
323                   }
324               }
325           }
326       } while(rawStream->NextDDL()); 
327     
328     AliMUONRawStreamTrackerHP::AliBusPatch* busPatch;
329     if (!eventIsErrorMessage) 
330     {
331       // Good events (no error) -> compute pedestal for all channels
332       rawStream->First(); 
333       nEvents++;
334       muonGain->SetAliNCurrentEvents(nEvents);
335       while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
336             {
337               for(int i = 0; i < busPatch->GetLength(); ++i)
338         {
339           busPatch->GetData(i, manuId, channelId, charge);
340           muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
341         }
342             }
343     }
344     else
345       {
346         // Events with errors
347         if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
348           {
349             filcout << " ----------- Date Event recovered = " << nDateEvents <<  " ----------------" << endl;
350             // Recover parity errors -> compute pedestal for all good buspatches
351             if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
352                                         ATTR_ORBIT_BC )) 
353               {
354                 filcout <<"Event recovered -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
355                         <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
356                         <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                              
357               } 
358             else 
359               {
360                 filcout <<"Event recovered -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
361                         <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
362                         <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
363               }
364             rawStream->First();
365             nEvents++;
366             muonGain->SetAliNCurrentEvents(nEvents);
367             while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
368               {
369                 // Check the buspatch -> if error not use it in the pedestal calculation
370                 int errorCount = 0;
371                 for(int i = 0; i < busPatch->GetLength(); ++i)
372                   {
373                     if (!busPatch->IsParityOk(i)) errorCount++;
374                   }
375                 if (!errorCount) 
376                   {
377                     // Good buspatch
378                     for(int i = 0; i < busPatch->GetLength(); ++i)
379                       {
380                         busPatch->GetData(i, manuId, channelId, charge);
381                         muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
382                       }
383                   }
384                 else
385                   {
386                     char bpname[256];
387                     AliMUONErrorCounter* errorCounter;
388                     // Bad buspatch -> not used (just print)
389                     filcout<<"bpId "<<busPatch->GetBusPatchId()<<" words "<<busPatch->GetLength()
390                            <<" parity errors "<<errorCount<<endl;
391                     // Number of events where this buspatch is missing
392                     sprintf(bpname,"bp%d",busPatch->GetBusPatchId());                                           
393                     if (!(errorCounter = (AliMUONErrorCounter*) (muonGain->GetErrorBuspatchTable()->FindObject(bpname))))
394                       {
395                         // New buspatch
396                         errorCounter = new AliMUONErrorCounter(busPatch->GetBusPatchId());
397                         errorCounter->SetName(bpname);
398                         muonGain->GetErrorBuspatchTable()->Add(errorCounter);
399                       }
400                     else
401                       {
402                         // Existing buspatch
403                         errorCounter->Increment();
404                       } 
405                     // errorCounter->Print();                                           
406                   } // end of if (!errorCount)
407               } // end of while( (busPatch = (AliMUONRawStreamTrackerHP ...
408             //        nEvents++;
409             nEventsRecovered++;
410           } //end of if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
411         else
412           {
413             // Fatal errors reject the event
414             detail = Form(" ----------- Date Event rejected = %d  ----------------",nDateEvents);
415             filcout << detail << endl;
416             if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
417                                         ATTR_ORBIT_BC )) 
418               {
419                 filcout <<"Event rejected -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
420                         <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
421                         <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                              
422               } 
423             else 
424               {
425                 filcout <<"Event rejected -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
426                         <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
427                         <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
428           
429               }
430           } // end of if (!rawStream->GetGlitchErrors() && !rawStream->GetPaddingErrors() ...
431         filcout<<"Number of errors : Glitch "<<eventGlitchErrors
432                <<" Parity "<<eventParityErrors
433                <<" Padding "<<eventPaddingErrors
434                <<" Token lost "<<eventTokenlostErrors<<endl;
435         filcout<<endl;                  
436       } // end of if (!rawStream->IsErrorMessage())
437     
438     if (eventGlitchErrors)  nGlitchErrors++;
439     if (eventParityErrors)  nParityErrors++;
440     if (eventPaddingErrors) nPaddingErrors++;
441     //if (eventTokenlostErrors) nTokenlostErrors++;
442
443   } // while (rawReader->NextEvent())
444   delete rawReader;
445   delete rawStream;
446   
447   // process and store mean peak values in .root file
448   sprintf(flatFile,"%s.par",prefixDA);
449   if(shuttleFile.IsNull())shuttleFile=flatFile;
450   injCharge=vDAC[nIndex-1];
451   muonGain->SetAliIndex(nIndex); // fIndex 
452   muonGain->SetAliInjCharge(injCharge);
453   muonGain->SetAliNEvents(nEvents);
454   muonGain->SetAliRunNumber(runNumber);
455   muonGain->MakePedStoreForGain(shuttleFile);
456   
457   
458   // writing some counters
459   cout << endl;
460   cout << prefixDA << " : Nb of DATE events           = " << nDateEvents    << endl;
461   cout << prefixDA << " : Nb of Glitch errors         = "   << nGlitchErrors  << endl;
462   cout << prefixDA << " : Nb of Parity errors         = "   << nParityErrors  << endl;
463   cout << prefixDA << " : Nb of Padding errors        = "   << nPaddingErrors << endl;          
464   cout << prefixDA << " : Nb of Token lost errors     = "   << nTokenlostErrors << endl;
465   cout << prefixDA << " : Nb of events recovered      = "   << nEventsRecovered<< endl;
466   cout << prefixDA << " : Nb of events without errors = "   << nEvents-nEventsRecovered<< endl;
467   cout << prefixDA << " : Nb of events used           = "   << nEvents        << endl;
468   
469   filcout << endl;
470   filcout << prefixDA << " : Nb of DATE events           = " << nDateEvents    << endl;
471   filcout << prefixDA << " : Nb of Glitch errors         = "   << nGlitchErrors << endl;
472   filcout << prefixDA << " : Nb of Parity errors         = "   << nParityErrors << endl;
473   filcout << prefixDA << " : Nb of Padding errors        = "   << nPaddingErrors << endl;
474   filcout << prefixDA << " : Nb of Token lost errors     = "   << nTokenlostErrors << endl;
475   filcout << prefixDA << " : Nb of events recovered      = "   << nEventsRecovered<< endl;      
476   filcout << prefixDA << " : Nb of events without errors = "   << nEvents-nEventsRecovered<< endl;
477   filcout << prefixDA << " : Nb of events used           = "   << nEvents        << endl;
478   
479   // Writing Token Error table
480   if(nTokenlostErrors)
481     {
482       char* detail=Form("\nWarning: Token Lost occurence \n");
483       printf("%s",detail);
484       filcout <<  detail ;
485       for ( Int_t i=0 ; i<20 ; i++) 
486         { 
487           for ( Int_t j=4 ; j<14 ; j++) 
488             { 
489               if(tabTokenError[i][j]>0)
490                 {
491                   Int_t tab=tabTokenError[i][j];
492                   Int_t frt = j/2-1;
493                   Int_t station = i/4 +1;
494                   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);
495                   else detail=Form(" in DDL= %d (station %d) and FRT%d (Down) => %d Token errors (address = 0x%X0000)",2560+i,station,frt,tab,j);
496                   printf("%s\n",detail);
497                   filcout <<  detail << endl;
498                 }
499             }
500         }
501     }
502
503   // Computing gain 
504   if(nIndex==nEntries)
505   {
506     muonGain->SetAliInit(nInit); // fnInit
507     muonGain->SetAliEntries(nEntries); // fnEntries
508     muonGain->SetAliNbpf1(nbpf1); // fnbpf1
509     muonGain->MakeGainStore(shuttleFile);
510 #ifdef ALI_AMORE  
511     std::ifstream in(shuttleFile.Data());
512     ostringstream stringout;
513     char line[1024];
514     while ( in.getline(line,1024) )
515       stringout << line << "\n";  
516     in.close();
517           
518     amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
519     TObjString gaindata(stringout.str().c_str());
520     status = amoreDA.Send("Gains",&gaindata);
521     if ( status )
522       cout << "Warning: Failed to write Pedestals in the AMORE database : " << status << endl;
523     else 
524       cout << "amoreDA.Send(Gains) ok" << endl;  
525 #else
526     cout << "Warning: MCH DA not compiled with AMORE support" << endl;
527 #endif
528   }
529   
530   // ouput files
531   filcout << endl;
532   filcout << prefixDA << " : Root data file         : " << muonGain->GetRootDataFileName() << endl;
533   filcout << prefixDA << " : Output logfile         : " << logOutputFile  << endl;
534   filcout << prefixDA << " : Gain Histo file        : " << muonGain->GetHistoFileName() << endl;
535   filcout << prefixDA << " : Gain file (to SHUTTLE) : " << shuttleFile << endl;
536   
537   //     Copying files to local DB folder defined by DAQ_DETDB_LOCAL
538   Char_t *dir;
539   dir= getenv("DAQ_DETDB_LOCAL");
540   unsigned int nLastVersions = 50;
541   printf("\n ***  Local DataBase: %s  (Max= %d) ***\n",dir,nLastVersions);
542   status = daqDA_localDB_storeFile(logOutputFile.Data(),nLastVersions);
543   if(status)printf(" Store file : %s   status = %d\n",logOutputFile.Data(),status);
544   if(nIndex==nEntries)
545   {
546     status = daqDA_localDB_storeFile(muonGain->GetRootDataFileName(),nLastVersions);
547     if(status)printf(" Store file : %s   status = %d\n",muonGain->GetRootDataFileName(),status);
548     status = daqDA_localDB_storeFile(muonGain->GetHistoFileName(),nLastVersions);
549     if(status)printf(" Store file : %s   status = %d\n",muonGain->GetHistoFileName(),status);
550     status = daqDA_localDB_storeFile(shuttleFile.Data(),nLastVersions);
551     if(status)printf(" Store file : %s   status = %d\n",shuttleFile.Data(),status);
552   }      
553   
554   
555   // ouput files
556   cout << endl;
557   cout << prefixDA << " : Root data file         : " << muonGain->GetRootDataFileName() << endl;
558   cout << prefixDA << " : Output logfile         : " << logOutputFile  << endl;
559   cout << prefixDA << " : Gain Histo file        : " << muonGain->GetHistoFileName() << endl;
560   cout << prefixDA << " : Gain file (to SHUTTLE) : " << shuttleFile << endl;   
561   
562   filcout.close();
563   
564   // Transferring to calibration file to  FES
565   // be sure that env variable DAQDALIB_PATH is set in script file
566   //       gSystem->Setenv("DAQDALIB_PATH", "$DATE_SITE/infoLogger");
567   printf("\n *****  STORE calibration FILE to FES ****** \n");
568   status = daqDA_FES_storeFile(shuttleFile.Data(),"GAINS");
569   if (status) { printf(" Failed to export file : %s , status = %d\n",shuttleFile.Data(),status); return -1; }
570   //  else printf(" %s successfully exported to FES  \n",shuttleFile.Data());
571   
572   printf("\n ######## End execution : %s ######## \n",prefixDA); 
573   timers.Stop();
574   printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());
575   return status;
576 }
577