]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONTRKPEDda.cxx
Adding data corruption recovery procedure for magnetic field = 0 runs. (http://savann...
[u/mrichter/AliRoot.git] / MUON / MUONTRKPEDda.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 Run: 109302 (station 3 only)
5   Run Type: PEDESTAL
6   DA Type: LDC
7   Number of events needed: 400 events for pedestal run
8   Input Files: mutrkpedvalues and config_ldc-MTRK-S3-0 in path : /afs/cern.ch/user/j/jcharvet/public/DA_validation 
9   Output Files: local dir (not persistent) -> MUONTRKPEDda.ped  FXS -> run<#>_MCH_<ldc>_PEDESTALS
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         2010-04-18 New version: MUONTRKPEDda.cxx,v 1.6
33         -------------------------------------------------------------------------
34
35         Version for MUONTRKPEDda 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 */
45 extern "C" {
46 #include <daqDA.h>
47 }
48
49 #include "event.h"
50 #include "monitor.h"
51
52 #include <Riostream.h>
53 #include <stdio.h>
54 #include <stdlib.h>
55 #include <sstream>
56 #include <math.h> 
57
58 //AliRoot
59 #include "AliMUONRawStreamTrackerHP.h"
60 #include "AliRawReader.h"
61 #include "AliMUONVStore.h"
62 #include "AliMUON2DMap.h"
63 #include "AliMUONCalibParamND.h"
64 #include "AliMpIntPair.h"
65 #include "AliMpConstants.h"
66 #include "AliRawDataErrorLog.h"
67 #include "AliMUONTrackerIO.h"
68 #include "AliLog.h"
69 #include "AliMUONDspHeader.h"
70 #include "AliDAQ.h"
71
72 //ROOT
73 #include "TFile.h"
74 #include "TSystem.h"
75 #include "TTree.h"
76 #include "TH1F.h"
77 #include "TString.h"
78 #include "TStopwatch.h"
79 #include "TMath.h"
80 #include "TTimeStamp.h"
81 #include "TGraphErrors.h"
82 #include "TF1.h"
83 #include "TROOT.h"
84 #include "TPluginManager.h"
85 #include "TFitter.h"
86 #include "TObjString.h"
87 #include "THashTable.h"
88 #include <THashList.h>
89 //
90 //AMORE
91 //
92 #ifdef ALI_AMORE
93 #include <AmoreDA.h>
94 #endif
95
96 #include "AliMUONPedestal.h"
97 #include "AliMUONErrorCounter.h"
98
99  
100 // main routine
101 int main(Int_t argc, const char **argv) 
102 {
103   Int_t status=0;
104   TStopwatch timers;
105   timers.Start(kTRUE); 
106
107   const char* prefixDA = "MUONTRKPEDda"; // program prefix
108   printf(" ######## Begin execution : %s ######## \n\n",prefixDA); 
109
110   // needed for streamer application
111   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
112                                         "*",
113                                         "TStreamerInfo",
114                                         "RIO",
115                                         "TStreamerInfo()"); 
116   // needed for Minuit plugin
117   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
118                                         "Minuit",
119                                         "TMinuitMinimizer",
120                                         "Minuit",
121                                         "TMinuitMinimizer(const char*)");
122 //   cout << argv[0];
123  
124   Int_t skipEvents = 0;
125   Int_t maxEvents  = 1000000;
126   Int_t maxDateEvents  = 1000000;
127   TString inputFile;
128
129   Int_t  nDateEvents = 0;
130   Int_t nGlitchErrors= 0;
131   Int_t nParityErrors= 0;
132   Int_t nPaddingErrors= 0;
133   Int_t nTokenlostErrors= 0;
134   Int_t recoverParityErrors = 1;
135
136   TString logOutputFile;
137
138   Char_t flatFile[256]="";
139   TString shuttleFile;
140
141   Int_t nEventsRecovered = 0;
142   Int_t nEvents = 0;
143   UInt_t runNumber   = 0;
144   Int_t nConfig = 1;
145   Int_t nEvthreshold = 10; //below this nb_evt pedestal are not calculated and forced to 4085 (sigma)
146   ofstream filcout;
147
148   // decode the input line
149   for (Int_t i = 1; i < argc; i++) // argument 0 is the executable name
150     {
151       const char* arg = argv[i];
152
153       arg = argv[i];
154       if (arg[0] != '-') 
155         {
156           // If only one argument and no "-" => DA calling from ECS
157           if (argc == 2)
158             {
159         inputFile=argv[i];
160             }
161           continue;
162         }
163       switch (arg[1])
164         {
165         case 'f' : 
166           i++;
167       inputFile=argv[i];
168           nConfig=0;
169           break;
170         case 'a' : 
171           i++;
172           shuttleFile = argv[i];
173           break;
174         case 's' :
175           i++; 
176           skipEvents=atoi(argv[i]);
177           break;
178         case 'm' :
179           i++; 
180           sscanf(argv[i],"%d",&maxDateEvents);
181           break;
182         case 'n' :
183           i++; 
184           sscanf(argv[i],"%d",&maxEvents);
185           break;
186         case 'h' :
187           i++;
188           printf("\n******************* %s usage **********************",argv[0]);
189           printf("\nOnline (called from ECS) : %s <raw data file> (no inline options)\n",argv[0]);
190           printf("\n%s can be used locally only with options (without DiMuon configuration file)",argv[0]);
191           printf("\n%s -options, the available options are :",argv[0]);
192           printf("\n-h help                    (this screen)");
193           printf("\n");
194           printf("\n Input");
195           printf("\n-f <raw data file>         (default = %s)",inputFile.Data()); 
196           printf("\n");
197           printf("\n Output");
198           printf("\n-a <Flat ASCII file>       (default = %s)",shuttleFile.Data()); 
199           printf("\n");
200           printf("\n Options");
201           printf("\n-m <max date events>       (default = %d)",maxDateEvents);
202           printf("\n-s <skip events>           (default = %d)",skipEvents);
203           printf("\n-n <max events>            (default = %d)",maxEvents);
204
205           printf("\n\n");
206           exit(-1);
207         default :
208           printf("%s : bad argument %s (please check %s -h)\n",argv[0],argv[i],argv[0]);
209           argc = 2; exit(-1); // exit if error
210         } // end of switch  
211     } // end of for i  
212
213   // decoding the events
214
215   UShort_t manuId;  
216   UChar_t channelId;
217   UShort_t charge;
218
219   //Pedestal object
220   AliMUONPedestal* muonPedestal = new AliMUONPedestal();
221   muonPedestal->SetprefixDA(prefixDA);
222
223   Char_t dbfile[256]="";
224   // nConfig=1 : Reading configuration (or not) status via "mutrkpedvalues" file located in DetDB
225   if(nConfig)
226     { 
227       sprintf(dbfile,"mutrkpedvalues");
228       status=daqDA_DB_getFile(dbfile,dbfile);
229       if(status) {printf(" !!! Failed  : input file %s is missing, status = %d\n",dbfile,status); return -1; } 
230       ifstream filein(dbfile,ios::in);
231       filein >> nConfig;
232       //      filein >> nEvthreshold;
233     }
234   else  printf(" ***  Config= %d: no configuration ascii file is used \n",nConfig); 
235   muonPedestal->SetconfigDA(nConfig);
236   muonPedestal->SetnEvthreshold(nEvthreshold);
237
238   // nConfig=1: configuration ascii file config_$DATE_ROLE_NAME read from DetDB
239   if(nConfig)
240     {
241       sprintf(dbfile,"config_%s",getenv("DATE_ROLE_NAME"));
242       status=daqDA_DB_getFile(dbfile,dbfile);
243       if(status) {printf(" !!! Failed  : Configuration file %s is missing, status = %d\n",dbfile,status); return -1; }
244       //      else printf(" *** Copy ascii config file: %s from DetDB to working directory and reading ...*** \n",dbfile);
245       muonPedestal->LoadConfig(dbfile);  
246     } 
247
248   // Rawdeader, RawStreamHP
249   AliRawReader* rawReader = AliRawReader::Create(inputFile);
250   AliMUONRawStreamTrackerHP* rawStream  = new AliMUONRawStreamTrackerHP(rawReader);    
251 //      rawStream->DisableWarnings();
252   rawStream->EnabbleErrorLogger();
253   //
254   // kLowErrorDetail,     /// Logs minimal information in the error messages.
255   // kMediumErrorDetail,  /// Logs a medium level of detail in the error messages.
256   // kHighErrorDetail     /// Logs maximum information in the error messages.
257   //  rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kLowErrorDetail);
258      rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kMediumErrorDetail);
259   //   rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kHighErrorDetail);
260
261   printf("\n %s : Reading data from file %s\n",prefixDA,inputFile.Data());
262
263   Int_t tabTokenError[20][14];
264   for ( Int_t i=0 ; i<20 ; i++) { for ( Int_t j=0 ; j<14 ; j++) { tabTokenError[i][j]=0;}       }
265
266   while (rawReader->NextEvent())
267     {
268       if (nDateEvents >= maxDateEvents) break;
269       if (nEvents >= maxEvents) break;
270       if (nDateEvents>0 &&  nDateEvents % 100 == 0)     
271         cout<<"Cumulated:  DATE events = " << nDateEvents << "   Used events = " << nEvents << endl;
272
273       // check shutdown condition 
274       if (daqDA_checkShutdown()) 
275         break;
276       //Skip events
277       while (skipEvents)
278         {
279           rawReader->NextEvent();
280           skipEvents--;
281         }
282       Int_t eventType = rawReader->GetType();
283       runNumber = rawReader->GetRunNumber();
284
285       // Output log file initialisations
286       if(nDateEvents==0)
287         {
288           sprintf(flatFile,"%s.log",prefixDA);
289           logOutputFile=flatFile;
290                 AliLog::SetStreamOutput(&filcout); // Print details on logfile
291           filcout.open(logOutputFile.Data());
292           filcout<<"//=================================================" << endl;
293           filcout<<"//       " << prefixDA << " for run = " << runNumber << endl;
294           filcout<<"//=================================================" << endl;
295           filcout<<"//   * Date          : " << muonPedestal->GetDate()->AsString("l") << "\n" << endl;
296           cout<<"\n ********  " << prefixDA << " for run = " << runNumber << " ********\n" << endl;
297           cout<<" * Date          : " << muonPedestal->GetDate()->AsString("l") << "\n" << endl;
298
299         }
300
301       muonPedestal->SetAlifilcout(&filcout);
302
303       nDateEvents++;
304       if (eventType != PHYSICS_EVENT)
305         continue; // for the moment
306
307       const char* detail = "";
308       // First lopp over DDL's to find good events
309       // Error counters per event (counters in the decoding lib are for each DDL)
310       Bool_t eventIsErrorMessage = kFALSE;
311       int eventGlitchErrors = 0;
312       int eventParityErrors = 0;
313       int eventPaddingErrors = 0;
314       int eventTokenlostErrors = 0;
315       rawStream->First();
316       do
317         {
318           if (rawStream->IsErrorMessage()) eventIsErrorMessage = kTRUE;
319           eventGlitchErrors += rawStream->GetGlitchErrors();
320           eventParityErrors += rawStream->GetParityErrors();
321           eventPaddingErrors += rawStream->GetPaddingErrors();
322           eventTokenlostErrors += rawStream->GetTokenLostErrors();
323           if (rawStream->GetTokenLostErrors())
324             {
325               nTokenlostErrors++;
326               const AliMUONRawStreamTrackerHP::AliBlockHeader*      blkHeader  = 0x0;
327               const AliMUONRawStreamTrackerHP::AliDspHeader*        dspHeader  = 0x0;
328               Int_t nBlock = rawStream->GetBlockCount();
329               for(Int_t iBlock = 0; iBlock < nBlock ;iBlock++)
330                 {
331                   blkHeader = rawStream->GetBlockHeader(iBlock);
332                   //              printf("Block %d Total length %d\n",iBlock,blkHeader->GetTotalLength());
333                   Int_t nDsp = rawStream->GetDspCount(iBlock);
334                   //              printf("Block %d DSP %d\n",iBlock,nDsp);                
335                   for(Int_t iDsp = 0; iDsp < nDsp ;iDsp++)
336                     {
337                       dspHeader =  blkHeader->GetDspHeader(iDsp);
338                       //                      printf("Dsp %d Add %X\n",iDsp,dspHeader);
339                       if (dspHeader->GetErrorWord())
340                         {
341                           Int_t ddl = rawStream->GetDDL()  ; 
342                           //     Int_t ddl = AliDAQ::DdlID("MUONTRK", rawStream->GetDDL()) - 2560 ; // format 2560 + ddl
343                           Int_t frt = (dspHeader->GetErrorWord() & 0xFFFF0000) >> 16 ; // 4*4bits right shift
344                           tabTokenError[ddl][frt]++;
345                           //     printf(" DDL %d error word %X %d %d\n",ddl,dspHeader->GetErrorWord(),frt,tabTokenError[8][4]);
346                         }
347                       
348                     }
349                 }
350             }
351         } while(rawStream->NextDDL()); 
352
353       AliMUONRawStreamTrackerHP::AliBusPatch* busPatch;
354       if (!eventIsErrorMessage) 
355         {
356           // Good events (no error) -> compute pedestal for all channels
357           rawStream->First(); 
358           nEvents++;
359           muonPedestal->SetAliNCurrentEvents(nEvents);
360           while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
361             {
362               for(int i = 0; i < busPatch->GetLength(); ++i)
363                 {
364                   busPatch->GetData(i, manuId, channelId, charge);
365                   muonPedestal->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
366                 }
367             }
368         }
369       else
370         {
371           // Events with errors
372           if (recoverParityErrors && eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
373             {
374               filcout << " ----------- Date Event recovered = " << nDateEvents <<  " ----------------" << endl;
375               // Recover parity errors -> compute pedestal for all good buspatches
376               if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
377                                           ATTR_ORBIT_BC )) 
378                 {
379                   filcout <<"Event recovered -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
380                               <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
381                               <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                                
382                 } 
383               else 
384                 {
385                   filcout <<"Event recovered -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
386                               <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
387                               <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
388                 }
389               rawStream->First();
390               nEvents++;
391               muonPedestal->SetAliNCurrentEvents(nEvents);
392               while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
393                 {
394                   // Check the buspatch -> if error not use it in the pedestal calculation
395                   int errorCount = 0;
396                   for(int i = 0; i < busPatch->GetLength(); ++i)
397                     {
398                       if (!busPatch->IsParityOk(i)) errorCount++;
399                     }
400                   if (!errorCount) 
401                     {
402                       // Good buspatch
403                       for(int i = 0; i < busPatch->GetLength(); ++i)
404                         {
405                           busPatch->GetData(i, manuId, channelId, charge);
406                           muonPedestal->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
407                         }
408                     }
409                   else
410                     {
411                       AliMUONErrorCounter* errorCounter;
412                       // Bad buspatch -> not used (just print)
413                       filcout<<"bpId "<<busPatch->GetBusPatchId()<<" words "<<busPatch->GetLength()
414                                  <<" parity errors "<<errorCount<<endl;
415                       // Number of events where this buspatch is missing
416                       if (!(errorCounter = (AliMUONErrorCounter*) (muonPedestal->GetErrorBuspatchTable()->FindObject(busPatch->GetBusPatchId()))))
417                         {
418                           // New buspatch
419                           errorCounter = new AliMUONErrorCounter(busPatch->GetBusPatchId());
420                           muonPedestal->GetErrorBuspatchTable()->Add(errorCounter);
421                         }
422                       else
423                         {
424                           // Existing buspatch
425                           errorCounter->Increment();
426                         }       
427                       // errorCounter->Print();                                         
428                     } // end of if (!errorCount)
429                 } // end of while( (busPatch = (AliMUONRawStreamTrackerHP ...
430               //              nEvents++;
431               nEventsRecovered++;
432             } //end of if (recoverParityErrors && eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
433           else
434             {
435               // Fatal errors reject the event
436               detail = Form(" ----------- Date Event rejected = %d  ----------------",nDateEvents);
437               filcout << detail << endl;
438               if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
439                                           ATTR_ORBIT_BC )) 
440                 {
441                   filcout <<"Event rejected -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
442                               <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
443                               <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                                
444                 } 
445               else 
446                 {
447                   filcout <<"Event rejected -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
448                               <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
449                               <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
450
451                 }
452             } // end of if (!rawStream->GetGlitchErrors() && !rawStream->GetPaddingErrors() ...
453           filcout<<"Number of errors : Glitch "<<eventGlitchErrors
454                      <<" Parity "<<eventParityErrors
455                      <<" Padding "<<eventPaddingErrors
456                      <<" Token lost "<<eventTokenlostErrors<<endl;
457           filcout<<endl;                        
458         } // end of if (!rawStream->IsErrorMessage())
459
460       if (eventGlitchErrors)  nGlitchErrors++;
461       if (eventParityErrors)  nParityErrors++;
462       if (eventPaddingErrors) nPaddingErrors++;
463       //      if (eventTokenlostErrors) nTokenlostErrors++;
464       //      muonPedestal->SetAliNCurrentEvents(nEvents);
465
466     } // while (rawReader->NextEvent())
467   delete rawReader;
468   delete rawStream;
469
470   sprintf(flatFile,"%s.ped",prefixDA);
471   if(shuttleFile.IsNull())shuttleFile=flatFile;
472   muonPedestal->SetAliNEvents(nEvents);
473   muonPedestal->SetAliRunNumber(runNumber);
474   
475   muonPedestal->Finalize();  
476   muonPedestal->MakeControlHistos();  
477   if (!shuttleFile.IsNull())  
478   {
479     ofstream out(shuttleFile.Data());  
480     muonPedestal->MakeASCIIoutput(out);
481     out.close();
482 #ifdef ALI_AMORE
483   //
484   //Send objects to the AMORE DB
485   //
486     ostringstream stringout;
487     muonPedestal->MakeASCIIoutput(stringout);
488     
489     amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
490     TObjString peddata(stringout.str().c_str());
491     Int_t amoreStatus = amoreDA.Send("Pedestals",&peddata);
492     if ( amoreStatus )
493       cout << "Warning: Failed to write Pedestals in the AMORE database : " << amoreStatus << endl;
494     else 
495       cout << "amoreDA.Send(Pedestals) ok" << endl;  
496 #else
497     cout << "Warning: MCH DA not compiled with AMORE support" << endl;
498 #endif
499     
500   }
501
502   // writing some counters
503   cout << endl;
504   cout << prefixDA << " : Nb of DATE events           = " << nDateEvents    << endl;
505   cout << prefixDA << " : Nb of Glitch errors         = "   << nGlitchErrors  << endl;
506   cout << prefixDA << " : Nb of Parity errors         = "   << nParityErrors  << endl;
507   cout << prefixDA << " : Nb of Padding errors        = "   << nPaddingErrors << endl;  
508   cout << prefixDA << " : Nb of Token lost errors     = "   << nTokenlostErrors << 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 Token lost errors     = "   << nTokenlostErrors << endl;
519   filcout << prefixDA << " : Nb of events recovered      = "   << nEventsRecovered<< endl;      
520   filcout << prefixDA << " : Nb of events without errors = "   << nEvents-nEventsRecovered<< endl;
521   filcout << prefixDA << " : Nb of events used           = "   << nEvents        << endl;
522
523   // ouput files
524   cout << endl;
525   cout << prefixDA << " : Output logfile         : " << logOutputFile  << endl;
526   cout << prefixDA << " : Pedestal Histo file    : " << muonPedestal->GetHistoFileName() << endl;
527   cout << prefixDA << " : Ped. file (to SHUTTLE) : " << shuttleFile << endl;   
528
529   // Writing Token Error table
530   if(nTokenlostErrors)
531     {
532       char* detail=Form("\nWarning: Token Lost occurence \n");
533       printf("%s",detail);
534       filcout <<  detail ;
535       for ( Int_t i=0 ; i<20 ; i++) 
536         { 
537           for ( Int_t j=4 ; j<14 ; j++) 
538             { 
539               if(tabTokenError[i][j]>0)
540                 {
541                   Int_t tab=tabTokenError[i][j];
542                   Int_t frt=j/2-1;
543                   Int_t station = i/4 +1;
544                   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);
545                   else detail=Form(" in DDL= %d (station %d) and FRT%d (Down) => %d Token errors (address = 0x%X0000)",2560+i,station,frt,tab,j);
546                   printf("%s\n",detail);
547                   filcout <<  detail << endl;
548                 }
549             }
550         }
551     }
552
553   filcout << endl;
554   filcout << prefixDA << " : Output logfile         : " << logOutputFile  << endl;
555   filcout << prefixDA << " : Pedestal Histo file    : " << muonPedestal->GetHistoFileName() << endl;
556   filcout << prefixDA << " : Ped. file (to SHUTTLE) : " << shuttleFile << endl;   
557
558  // Copying files to local DB folder defined by DAQ_DETDB_LOCAL
559   Char_t *dir;
560   dir= getenv("DAQ_DETDB_LOCAL");
561   unsigned int nLastVersions=50;
562   cout << "\n ***  Local DataBase: " << dir << " (Max= " << nLastVersions << ") ***" << endl;
563   status = daqDA_localDB_storeFile(muonPedestal->GetHistoFileName(),nLastVersions);
564   if(status)printf(" Store file : %s   status = %d\n",muonPedestal->GetHistoFileName(),status);
565   status = daqDA_localDB_storeFile(shuttleFile.Data(),nLastVersions);
566   if(status)printf(" Store file : %s   status = %d\n",shuttleFile.Data(),status);
567   status = daqDA_localDB_storeFile(logOutputFile.Data(),nLastVersions);
568   if(status)printf(" Store file : %s   status = %d\n",logOutputFile.Data(),status);
569
570   filcout.close();
571
572   // Transferring to FES  (be sure that env variable DAQDALIB_PATH is set)
573   printf("\n *****  STORE Pedestal FILE in FES ****** \n");
574   status = daqDA_FES_storeFile(shuttleFile.Data(),"PEDESTALS");
575   if (status) { printf(" !!! Failed to export file : %s , status = %d\n",shuttleFile.Data(),status); return -1; }
576   //  else printf(" %s successfully exported to FES  \n",shuttleFile.Data());
577
578   // Transferring to FES  (be sure that env variable DAQDALIB_PATH is set)
579   printf("\n *****  STORE Configuration FILE in FES ****** \n");
580   status = daqDA_FES_storeFile(dbfile,"CONFIG");
581   if (status) { printf(" !!! Failed to export file : %s , status = %d\n",dbfile,status); return -1; }
582   //  else printf(" %s successfully exported to FES  \n",dbfile);
583
584
585   printf("\n ######## End execution : %s ######## \n",prefixDA); 
586   timers.Stop();
587   printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());
588   return status;
589 }