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