]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSteer.cxx
Raw-reader is instantiated via the newly create factory method
[u/mrichter/AliRoot.git] / STEER / AliQADataMakerSteer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 // class for running the QA makers                                           //
20 //                                                                           //
21 //   AliQADataMakerSteer qas;                                                //
22 //   qas.Run(AliQA::kRAWS, rawROOTFileName);                                 //
23 //   qas.Run(AliQA::kHITS);                                                  //
24 //   qas.Run(AliQA::kSDIGITS);                                               //
25 //   qas.Run(AliQA::kDIGITS);                                                //
26 //   qas.Run(AliQA::kRECPOINTS);                                             //
27 //   qas.Run(AliQA::kESDS);                                                  //
28 //                                                                           //
29 ///////////////////////////////////////////////////////////////////////////////
30
31 #include <TKey.h>
32 #include <TFile.h>
33 #include <TFileMerger.h>
34 #include <TPluginManager.h>
35 #include <TROOT.h>
36 #include <TString.h>
37 #include <TSystem.h>
38
39 #include "AliCDBManager.h"
40 #include "AliCDBEntry.h"
41 #include "AliCDBId.h"
42 #include "AliCDBMetaData.h"
43 #include "AliESDEvent.h"
44 #include "AliHeader.h"
45 #include "AliLog.h"
46 #include "AliModule.h"
47 #include "AliQA.h"
48 #include "AliQADataMakerRec.h"
49 #include "AliQADataMakerSim.h"
50 #include "AliQADataMakerSteer.h" 
51 #include "AliRawReaderDate.h"
52 #include "AliRawReaderFile.h"
53 #include "AliRawReaderRoot.h"
54 #include "AliRun.h"
55 #include "AliRunLoader.h"
56
57 ClassImp(AliQADataMakerSteer) 
58
59 //_____________________________________________________________________________
60 AliQADataMakerSteer::AliQADataMakerSteer(const char* gAliceFilename, const char * name, const char * title) :
61         TNamed(name, title), 
62         fCurrentEvent(0),  
63         fCycleSame(kFALSE),
64         fDetectors("ALL"), 
65         fDetectorsW("ALL"), 
66         fESD(NULL), 
67         fESDTree(NULL),
68         fFirst(kTRUE),  
69         fGAliceFileName(gAliceFilename), 
70     fMaxEvents(0),        
71         fNumberOfEvents(999999), 
72     fRunNumber(0), 
73         fRawReader(NULL), 
74         fRawReaderDelete(kTRUE), 
75         fRunLoader(NULL),
76         fQADataMakers()
77
78 {
79         // default ctor
80         fMaxEvents = fNumberOfEvents ; 
81         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
82                 if (IsSelected(AliQA::GetDetName(iDet))) {
83                         fLoader[iDet]      = NULL ;
84                         fQADataMaker[iDet] = NULL ;
85                         fQACycles[iDet]    = 999999 ;
86                 }
87         }       
88 }
89
90 //_____________________________________________________________________________
91 AliQADataMakerSteer::AliQADataMakerSteer(const AliQADataMakerSteer & qas) : 
92         TNamed(qas), 
93         fCurrentEvent(qas.fCurrentEvent),  
94         fCycleSame(kFALSE),
95         fDetectors(qas.fDetectors), 
96         fDetectorsW(qas.fDetectorsW), 
97         fESD(NULL), 
98         fESDTree(NULL), 
99         fFirst(qas.fFirst),  
100         fGAliceFileName(qas.fGAliceFileName), 
101     fMaxEvents(qas.fMaxEvents),        
102         fNumberOfEvents(qas.fNumberOfEvents), 
103     fRunNumber(qas.fRunNumber), 
104         fRawReader(NULL), 
105         fRawReaderDelete(kTRUE), 
106         fRunLoader(NULL),
107         fQADataMakers()
108 {
109         // cpy ctor
110         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
111                 fLoader[iDet]      = qas.fLoader[iDet] ;
112                 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
113                 fQACycles[iDet]    = qas.fQACycles[iDet] ;      
114         }
115 }
116
117 //_____________________________________________________________________________
118 AliQADataMakerSteer & AliQADataMakerSteer::operator = (const AliQADataMakerSteer & qas) 
119 {
120         // assignment operator
121   this->~AliQADataMakerSteer() ;
122   new(this) AliQADataMakerSteer(qas) ;
123   return *this ;
124 }
125
126 //_____________________________________________________________________________
127 AliQADataMakerSteer::~AliQADataMakerSteer() 
128 {
129         // dtor
130   for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
131           if (IsSelected(AliQA::GetDetName(iDet))) {
132                   fLoader[iDet] = NULL;
133                   if (fQADataMaker[iDet]) {
134                           (fQADataMaker[iDet])->Finish() ; 
135                           delete fQADataMaker[iDet] ;
136                           fQADataMaker[iDet] = NULL ;
137                   }
138           }
139   }
140
141   if (fRawReaderDelete) { 
142         fRunLoader = NULL ;
143         delete fRawReader ;
144         fRawReader = NULL ;
145   }
146 }
147
148 //_____________________________________________________________________________
149 Bool_t AliQADataMakerSteer::DoIt(const AliQA::TASKINDEX_t taskIndex, const char * mode)
150 {
151         // Runs all the QA data Maker for every detector
152
153         Bool_t rv = kFALSE ;
154     // Fill QA data in event loop 
155         for (UInt_t iEvent = 0 ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
156                 fCurrentEvent++ ; 
157                 // Get the event
158                 if ( iEvent%10 == 0  ) 
159                         AliInfo(Form("processing event %d", iEvent));
160                 if ( taskIndex == AliQA::kRAWS ) {
161                         if ( !fRawReader->NextEvent() )
162                                 break ;
163                 } else if ( taskIndex == AliQA::kESDS ) {
164                         if ( fESDTree->GetEntry(iEvent) == 0 )
165                                 break ;
166                 } else {
167                         if ( fRunLoader->GetEvent(iEvent) != 0 )
168                                 break ;
169                 }
170                 // loop  over active loaders
171                 for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
172                         AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i)); //GetQADataMaker(iDet, mode) ;
173                         if ( qadm->IsCycleDone() ) {
174                                 qadm->EndOfCycle(AliQA::kRAWS) ;
175                                 qadm->StartOfCycle(AliQA::kRAWS) ;
176                         }
177                         TTree * data ; 
178                         Int_t iDet = qadm->GetUniqueID();
179                         AliLoader* loader = GetLoader(iDet);
180                         switch (taskIndex) {
181                                 case AliQA::kRAWS :
182                                         qadm->Exec(taskIndex, fRawReader) ; 
183                                         break ; 
184                                 case AliQA::kHITS :
185                                         loader->LoadHits() ; 
186                                         data = loader->TreeH() ; 
187                                         if ( ! data ) {
188                                                 AliWarning(Form(" Hit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
189                                         } else {
190                                                 qadm->Exec(taskIndex, data) ;
191                                         } 
192                                         break ;
193                                         case AliQA::kSDIGITS :
194                                         loader->LoadSDigits() ; 
195                                         data = loader->TreeS() ; 
196                                         if ( ! data ) {
197                                                 AliWarning(Form(" SDigit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
198                                         } else {
199                                                 qadm->Exec(taskIndex, data) ; 
200                                         }
201                                         break; 
202                                         case AliQA::kDIGITS :
203                                         loader->LoadDigits() ; 
204                                         data = loader->TreeD() ; 
205                                         if ( ! data ) {
206                                                 AliWarning(Form(" Digit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
207                                         } else {
208                                                 qadm->Exec(taskIndex, data) ;
209                                         }
210                                         break; 
211                                         case AliQA::kRECPOINTS :
212                                         loader->LoadRecPoints() ; 
213                                         data = loader->TreeR() ; 
214                                         if (!data) {
215                                                 AliWarning(Form("RecPoints not found for %s", AliQA::GetDetName(iDet))) ; 
216                                         } else {
217                                                 qadm->Exec(taskIndex, data) ; 
218                                         }
219                                         break; 
220                                         case AliQA::kTRACKSEGMENTS :
221                                         break; 
222                                         case AliQA::kRECPARTICLES :
223                                         break; 
224                                         case AliQA::kESDS :
225                                         qadm->Exec(taskIndex, fESD) ;
226                                         break; 
227                                         case AliQA::kNTASKINDEX :
228                                         break; 
229                         } //task switch
230                         //qadm->Increment() ; 
231                 } // detector loop
232         } // event loop 
233 //      // Save QA data for all detectors
234         rv = Finish(taskIndex, mode) ;
235
236         return rv ; 
237 }
238
239 //_____________________________________________________________________________
240 Bool_t AliQADataMakerSteer::Finish(const AliQA::TASKINDEX_t taskIndex, const char * /*mode*/) 
241 {
242         // write output to file for all detectors
243         for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
244                 AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i));
245                 qadm->EndOfCycle(taskIndex) ; 
246         }    
247         return kTRUE ; 
248 }
249
250 //_____________________________________________________________________________
251 TObjArray * AliQADataMakerSteer::GetFromOCDB(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, const char * year) const 
252 {
253         // Retrieve the list of QA data for a given detector and a given task 
254         TObjArray * rv = NULL ;
255         TString tmp(AliQA::GetQARefStorage()) ; 
256         if ( tmp.IsNull() ) { 
257                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
258                 return NULL ; 
259         }       
260         AliCDBManager* man = AliCDBManager::Instance() ; 
261         if ( ! man->IsDefaultStorageSet() ) {
262                 TString tmp(AliQA::GetQARefDefaultStorage()) ; 
263                 tmp.Append(year) ; 
264                 tmp.Append("/") ; 
265                 man->SetDefaultStorage(tmp.Data()) ;            
266                 man->SetSpecificStorage(Form("%s/*", AliQA::GetQAName()), AliQA::GetQARefStorage()) ;
267         }
268         char detOCDBDir[10] ; 
269         sprintf(detOCDBDir, "%s/%s/%s", AliQA::GetQAName(), AliQA::GetDetName((Int_t)det), AliQA::GetRefOCDBDirName()) ; 
270         AliInfo(Form("Retrieving reference data from %s/%s for %s", AliQA::GetQARefStorage(), detOCDBDir, AliQA::GetTaskName(task).Data())) ; 
271         AliCDBEntry* entry = man->Get(detOCDBDir, 0) ; //FIXME 0 --> Run Number
272         TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
273         if ( listDetQAD ) 
274                 rv = dynamic_cast<TObjArray *>(listDetQAD->FindObject(AliQA::GetTaskName(task))) ; 
275         return rv ; 
276 }
277
278 //_____________________________________________________________________________
279 AliLoader * AliQADataMakerSteer::GetLoader(Int_t iDet)
280 {
281         // get the loader for a detector
282
283         if ( !fRunLoader ) 
284                 return NULL ; 
285         
286         TString detName = AliQA::GetDetName(iDet) ;
287     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
288         if (fLoader[iDet]) 
289                 return fLoader[iDet] ;
290         
291         // load the QA data maker object
292         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
293         TString loaderName = "Ali" + detName + "Loader" ;
294
295         AliLoader * loader = NULL ;
296         // first check if a plugin is defined for the quality assurance data maker
297         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
298         // if not, add a plugin for it
299         if (!pluginHandler) {
300                 AliDebug(1, Form("defining plugin for %s", loaderName.Data())) ;
301                 TString libs = gSystem->GetLibraries() ;
302                 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
303                         pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
304                 } else {
305                         pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
306                 }
307                 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
308         }
309         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
310                 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
311         }
312         if (loader) 
313                 fLoader[iDet] = loader ;
314         return loader ;
315 }
316
317 //_____________________________________________________________________________
318 AliQADataMaker * AliQADataMakerSteer::GetQADataMaker(const Int_t iDet, const char * mode )
319 {
320         // get the quality assurance data maker for a detector
321         
322         if (fQADataMaker[iDet]) 
323                 return fQADataMaker[iDet] ;
324         
325         AliQADataMaker * qadm = NULL ;
326         
327         if (fFirst) {
328                 // load the QA data maker object
329                 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
330                 TString detName = AliQA::GetDetName(iDet) ;
331                 TString tmp(mode) ; 
332                 if (tmp.Contains("sim")) 
333                         tmp.ReplaceAll("s", "S") ; 
334                 else if (tmp.Contains("rec")) 
335                         tmp.ReplaceAll("r", "R") ; 
336                 TString qadmName = "Ali" + detName + "QADataMaker" + tmp ;
337
338                 // first check if a plugin is defined for the quality assurance data maker
339                 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
340                 // if not, add a plugin for it
341                 if (!pluginHandler) {
342                         AliDebug(1, Form("defining plugin for %s", qadmName.Data())) ;
343                         TString libs = gSystem->GetLibraries() ;
344                         if (libs.Contains("lib" + detName + mode + ".so") || (gSystem->Load("lib" + detName + mode + ".so") >= 0)) {
345                                 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
346                         } else {
347                                 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
348                         }
349                         pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
350                 }
351                 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
352                         qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
353                 }
354                 if (qadm) 
355                         fQADataMaker[iDet] = qadm ;
356         }
357         return qadm ;
358 }
359
360 //_____________________________________________________________________________
361 Bool_t AliQADataMakerSteer::Init(const AliQA::TASKINDEX_t taskIndex, const char * mode, const  char * input )
362 {
363         // Initialize the event source and QA data makers
364         
365         if (taskIndex == AliQA::kRAWS) { 
366                 if (!fRawReader) {
367                         fRawReader = AliRawReader::Create(input);
368                 }
369             if ( ! fRawReader ) 
370                         return kFALSE ; 
371                 fRawReaderDelete = kTRUE ; 
372                 fRawReader->NextEvent() ; 
373                 fRunNumber = fRawReader->GetRunNumber() ; 
374                 AliCDBManager::Instance()->SetRun(fRunNumber) ; 
375                 fRawReader->RewindEvents();
376                 fNumberOfEvents = 999999 ;
377                 if ( fMaxEvents < 0 ) 
378                         fMaxEvents = fNumberOfEvents ; 
379                 } else if (taskIndex == AliQA::kESDS) {
380                 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
381                         TFile * esdFile = TFile::Open("AliESDs.root") ;
382                         fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ; 
383                         if ( !fESDTree ) {
384                                 AliError("esdTree not found") ; 
385                                 return kFALSE ; 
386                         } else {
387                                 fESD     = new AliESDEvent() ;
388                                 fESD->ReadFromTree(fESDTree) ;
389                                 fESDTree->GetEntry(0) ; 
390                                 fRunNumber = fESD->GetRunNumber() ; 
391                                 fNumberOfEvents = fESDTree->GetEntries() ;
392                                 if ( fMaxEvents < 0 ) 
393                                         fMaxEvents = fNumberOfEvents ; 
394                         }
395                 } else {
396                         AliError("AliESDs.root not found") ; 
397                         return kFALSE ; 
398                 }                       
399         } else {
400                 if ( !InitRunLoader() ) { 
401                         AliWarning("No Run Loader not found") ; 
402                 } else {
403                         fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
404                         if ( fMaxEvents < 0 ) 
405                                 fMaxEvents = fNumberOfEvents ; 
406
407                 }
408         }
409
410         // Get Detectors 
411         TObjArray* detArray = NULL ; 
412         if (fRunLoader) // check if RunLoader exists 
413                 if ( fRunLoader->GetAliRun() ) // check if AliRun exists in gAlice.root
414                         detArray = fRunLoader->GetAliRun()->Detectors() ;
415         // Build array of QA data makers for all detectors
416         fQADataMakers.Clear();
417         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
418                 if (IsSelected(AliQA::GetDetName(iDet))) {
419                         AliQADataMaker * qadm = GetQADataMaker(iDet, mode) ;
420                         if (!qadm) {
421                                 AliError(Form("AliQADataMaker not found for %s", AliQA::GetDetName(iDet))) ; 
422                                 fDetectorsW.ReplaceAll(AliQA::GetDetName(iDet), "") ; 
423                         } else {
424                                 AliDebug(1, Form("Data Maker found for %s", qadm->GetName())) ; 
425                                 // skip non active detectors
426                                 if (detArray) {
427                                         AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
428                                         if (!det || !det->IsActive())  
429                                                 continue ;
430                                 }
431                                 qadm->SetName(AliQA::GetDetName(iDet));
432                                 qadm->SetUniqueID(iDet);
433                                 fQADataMakers.Add(qadm);
434                         }
435                 }
436         } 
437         // Initialize all QA data makers for all detectors
438         for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
439                 AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i));
440                 qadm->Init(taskIndex, fRunNumber, GetQACycles(qadm->GetUniqueID())) ;
441                 qadm->StartOfCycle(taskIndex, fCycleSame) ;
442         } 
443         fFirst = kFALSE ;
444         return kTRUE ; 
445 }
446
447 //_____________________________________________________________________________
448 Bool_t AliQADataMakerSteer::InitRunLoader()
449 {
450         // get or create the run loader
451         if (fRunLoader) {
452                 fCycleSame = kTRUE ; 
453         } else {
454                 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
455                         // load all base libraries to get the loader classes
456                         TString libs = gSystem->GetLibraries() ;
457                         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
458                                 if (!IsSelected(AliQA::GetDetName(iDet))) 
459                                         continue ; 
460                                 TString detName = AliQA::GetDetName(iDet) ;
461                                 if (detName == "HLT") 
462                                         continue;
463                                 if (libs.Contains("lib" + detName + "base.so")) 
464                                         continue;
465                                 gSystem->Load("lib" + detName + "base.so");
466                         }
467                         fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
468                         if (!fRunLoader) {
469                                 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
470                                 return kFALSE;
471                         }
472                         fRunLoader->CdGAFile();
473                         if (fRunLoader->LoadgAlice() == 0) {
474                                 gAlice = fRunLoader->GetAliRun();
475                         }
476
477                         if (!gAlice) {
478                                 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
479                                 return kFALSE;
480                         }
481
482                 } else {               // galice.root does not exist
483                         AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
484                         return kFALSE;
485                 }
486         }
487
488         if (!fRunNumber) { 
489                 fRunLoader->LoadHeader();
490                 fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
491         }
492         return kTRUE;
493 }
494
495 //_____________________________________________________________________________
496 Bool_t AliQADataMakerSteer::IsSelected(const char * det) 
497 {
498         // check whether detName is contained in detectors
499         // if yes, it is removed from detectors
500         
501         Bool_t rv = kFALSE;
502         const TString detName(det) ;
503         // check if all detectors are selected
504         if (fDetectors.Contains("ALL")) {
505                 fDetectors = "ALL";
506                 rv = kTRUE;
507         } else if ((fDetectors.CompareTo(detName) == 0) ||
508                            fDetectors.BeginsWith(detName+" ") ||
509                            fDetectors.EndsWith(" "+detName) ||
510                            fDetectors.Contains(" "+detName+" ")) {
511                 //              fDetectors.ReplaceAll(detName, "");
512                 rv = kTRUE;
513         }
514         
515         // clean up the detectors string
516         //      while (fDetectors.Contains("  ")) 
517         //              fDetectors.ReplaceAll("  ", " ");
518         //      while (fDetectors.BeginsWith(" ")) 
519         //              fDetectors.Remove(0, 1);
520         //      while (fDetectors.EndsWith(" ")) 
521         //              fDetectors.Remove(fDetectors.Length()-1, 1);
522         
523         return rv ;
524 }
525
526 //_____________________________________________________________________________
527 Bool_t AliQADataMakerSteer::Merge(const Int_t runNumber) const
528 {
529         // Merge all the cycles from all detectors in one single file per run
530         char cmd[80] ;
531         if ( runNumber == -1 )
532                 sprintf(cmd, ".! ls *%s*.*.*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
533         else 
534                 sprintf(cmd, ".! ls *%s*.%d.*.root > tempo.txt", AliQA::GetQADataFileName(), runNumber) ; 
535         gROOT->ProcessLine(cmd) ;
536         ifstream in("tempo.txt") ; 
537         const Int_t runMax = 10 ;  
538         TString file[AliQA::kNDET*runMax] ;
539         Int_t run[AliQA::kNDET*runMax] ;
540         
541         Int_t index = 0 ; 
542         while ( 1 ) {
543                 in >> file[index] ; 
544                 if ( !in.good() ) 
545                         break ; 
546                 AliInfo(Form("index = %d file = %s", index, (file[index]).Data())) ; 
547                 index++ ;
548         }
549         
550         if ( index == 0 ) { 
551                 AliError(Form("run number %d not found", runNumber)) ; 
552                 return kFALSE ; 
553         }
554         
555         Int_t runIndex    = 0 ;  
556         Int_t runIndexMax = 0 ; 
557         char stmp[10] ; 
558         sprintf(stmp, ".%s.", AliQA::GetQADataFileName()) ; 
559         for (Int_t ifile = 0 ; ifile < index ; ifile++) {
560                 TString tmp(file[ifile]) ; 
561                 tmp.ReplaceAll(".root", "") ; 
562                 TString det = tmp(0, tmp.Index(".")) ; 
563                 tmp.Remove(0, tmp.Index(stmp)+4) ; 
564                 TString ttmp = tmp(0, tmp.Index(".")) ; 
565                 Int_t newRun = ttmp.Atoi() ;
566                 for (Int_t irun = 0; irun <= runIndexMax; irun++) {
567                         if (newRun == run[irun]) 
568                                 break ; 
569                         run[runIndex] = newRun ; 
570                         runIndex++ ; 
571                 }
572                 runIndexMax = runIndex ; 
573                 ttmp = tmp(tmp.Index(".")+1, tmp.Length()) ; 
574                 Int_t cycle = ttmp.Atoi() ;  
575                 AliDebug(1, Form("%s : det = %s run = %d cycle = %d \n", file[ifile].Data(), det.Data(), newRun, cycle)) ; 
576         }
577         for (Int_t irun = 0 ; irun < runIndexMax ; irun++) {
578                 TFileMerger merger ; 
579                 char outFileName[20] ; 
580                 sprintf(outFileName, "Merged.%s.%d.root", AliQA::GetQADataFileName(), run[irun]) ; 
581                 merger.OutputFile(outFileName) ; 
582                 for (Int_t ifile = 0 ; ifile < index-1 ; ifile++) {
583                         char pattern[100] ; 
584                         sprintf(pattern, "%s.%d.", AliQA::GetQADataFileName(), run[irun]) ; 
585                         TString tmp(file[ifile]) ; 
586                         if (tmp.Contains(pattern)) {
587                                 merger.AddFile(tmp) ; 
588                         }
589                 }
590                 merger.Merge() ; 
591         }
592         
593         return kTRUE ; 
594 }
595
596 //_____________________________________________________________________________
597 void AliQADataMakerSteer::Reset(const Bool_t sameCycle)
598 {
599         // Reset the default data members
600
601         for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
602                 AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i));
603                 qadm->Reset(sameCycle);
604         } 
605         if (fRawReaderDelete) { 
606                 delete fRawReader ;
607                 fRawReader      = NULL ;
608         }
609
610         fCycleSame      = sameCycle ; 
611         fESD            = NULL ; 
612         fESDTree        = NULL ; 
613         fFirst          = kTRUE ;   
614         fNumberOfEvents = 999999 ;  
615 }
616
617 //_____________________________________________________________________________
618 TString AliQADataMakerSteer::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle) 
619 {
620         //Runs all the QA data Maker for Raws only
621         
622         fCycleSame       = sameCycle ;
623         fRawReader       = rawReader ;
624         fDetectors       = detectors ; 
625         fDetectorsW      = detectors ;  
626         
627         if ( !Init(AliQA::kRAWS, "rec") ) 
628                 return kFALSE ; 
629         fRawReaderDelete = kFALSE ; 
630
631         DoIt(AliQA::kRAWS, "rec") ; 
632         return  fDetectorsW ;
633 }
634
635 //_____________________________________________________________________________
636 TString AliQADataMakerSteer::Run(const char * detectors, const char * fileName, const Bool_t sameCycle) 
637 {
638         //Runs all the QA data Maker for Raws only
639
640         fCycleSame       = sameCycle ;
641         fDetectors       = detectors ; 
642         fDetectorsW      = detectors ;  
643         
644         if ( !Init(AliQA::kRAWS, "rec", fileName) ) 
645                 return kFALSE ; 
646
647         DoIt(AliQA::kRAWS, "rec") ; 
648         return  fDetectorsW ;
649 }
650
651 //_____________________________________________________________________________
652 TString AliQADataMakerSteer::Run(const char * detectors, const AliQA::TASKINDEX_t taskIndex, Bool_t const sameCycle, const  char * fileName ) 
653 {
654         // Runs all the QA data Maker for every detector
655         
656         fCycleSame       = sameCycle ;
657         fDetectors       = detectors ; 
658         fDetectorsW      = detectors ;          
659         
660         TString mode ; 
661         if ( (taskIndex == AliQA::kHITS) || (taskIndex == AliQA::kSDIGITS) || (taskIndex == AliQA::kDIGITS) ) 
662                 mode = "sim" ; 
663         else if ( (taskIndex == AliQA::kRAWS) || (taskIndex == AliQA::kRECPOINTS) || (taskIndex == AliQA::kESDS) )
664                 mode = "rec" ; 
665         else {
666                 AliError(Form("%s not implemented", AliQA::GetTaskName(taskIndex).Data())) ; 
667                 return "" ;
668         }
669
670         if ( !Init(taskIndex, mode.Data(), fileName) ) 
671                 return kFALSE ; 
672
673         DoIt(taskIndex, mode.Data()) ;
674         
675         return fDetectorsW ;
676
677 }
678
679 //_____________________________________________________________________________
680 Bool_t AliQADataMakerSteer::Save2OCDB(const Int_t runNumber, const char * year, const Int_t cycleNumber, const char * detectors) const
681 {
682         // take the locasl QA data merge into a single file and save in OCDB 
683         Bool_t rv = kTRUE ; 
684         TString tmp(AliQA::GetQARefStorage()) ; 
685         if ( tmp.IsNull() ) { 
686                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
687                 return kFALSE ; 
688         }
689         if ( !(tmp.Contains(AliQA::GetLabLocalOCDB()) || tmp.Contains(AliQA::GetLabAliEnOCDB())) ) {
690                 AliError(Form("%s is a wrong storage, use %s or %s", AliQA::GetQARefStorage(), AliQA::GetLabLocalOCDB().Data(), AliQA::GetLabAliEnOCDB().Data())) ; 
691                 return kFALSE ; 
692         }
693         TString sdet(detectors) ; 
694         sdet.ToUpper() ;
695         TFile * inputFile ; 
696         if ( sdet.Contains("ALL") ) {
697                 rv = Merge(runNumber) ; 
698                 if ( ! rv )
699                         return kFALSE ; 
700                 char inputFileName[20] ; 
701                 sprintf(inputFileName, "Merged.%s.%d.root", AliQA::GetQADataFileName(), runNumber) ; 
702                 inputFile = TFile::Open(inputFileName) ; 
703                 rv = SaveIt2OCDB(runNumber, inputFile, year) ; 
704         } else {
705                 for (Int_t index = 0; index < AliQA::kNDET; index++) {
706                         if (sdet.Contains(AliQA::GetDetName(index))) {
707                                 char inputFileName[20] ; 
708                                 sprintf(inputFileName, "%s.%s.%d.%d.root", AliQA::GetDetName(index), AliQA::GetQADataFileName(), runNumber, cycleNumber) ; 
709                                 inputFile = TFile::Open(inputFileName) ;                        
710                                 rv *= SaveIt2OCDB(runNumber, inputFile, year) ; 
711                         }
712                 }
713         }
714         return rv ; 
715 }
716
717 //_____________________________________________________________________________
718 Bool_t AliQADataMakerSteer::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year) const
719 {
720         // reads the TH1 from file and adds it to appropriate list before saving to OCDB
721         Bool_t rv = kTRUE ;
722         AliInfo(Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQA::GetQARefStorage())) ; 
723         AliCDBManager* man = AliCDBManager::Instance() ; 
724         if ( ! man->IsDefaultStorageSet() ) {
725                 TString tmp(AliQA::GetQARefDefaultStorage()) ; 
726                 tmp.Append(year) ; 
727                 tmp.Append("?user=alidaq") ; 
728                 man->SetDefaultStorage(tmp.Data()) ; 
729                 man->SetSpecificStorage("*", AliQA::GetQARefStorage()) ; 
730         }
731         if(man->GetRun() < 0) 
732                 man->SetRun(runNumber);
733
734         AliCDBMetaData mdr ;
735         mdr.SetResponsible("yves schutz");
736
737         for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++) {
738                 TDirectory * detDir = inputFile->GetDirectory(AliQA::GetDetName(detIndex)) ; 
739                 if ( detDir ) {
740                         AliInfo(Form("Entering %s", detDir->GetName())) ;
741                         char detOCDBDir[20] ;
742                         sprintf(detOCDBDir, "%s/%s/%s", AliQA::GetDetName(detIndex), AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName()) ; 
743                         AliCDBId idr(detOCDBDir, runNumber, AliCDBRunRange::Infinity())  ;
744                         TList * listDetQAD = new TList() ;
745                         char listName[20] ; 
746                         sprintf(listName, "%s QA data Reference", AliQA::GetDetName(detIndex)) ; 
747                         mdr.SetComment("HMPID QA stuff");
748                         listDetQAD->SetName(listName) ; 
749                         TList * taskList = detDir->GetListOfKeys() ; 
750                         TIter nextTask(taskList) ; 
751                         TKey * taskKey ; 
752                         while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
753                                 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ; 
754                                 AliInfo(Form("Saving %s", taskDir->GetName())) ; 
755                                 TObjArray * listTaskQAD = new TObjArray(100) ; 
756                                 listTaskQAD->SetName(taskKey->GetName()) ;
757                                 listDetQAD->Add(listTaskQAD) ; 
758                                 TList * histList = taskDir->GetListOfKeys() ; 
759                                 TIter nextHist(histList) ; 
760                                 TKey * histKey ; 
761                                 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
762                                         TObject * odata = taskDir->Get(histKey->GetName()) ; 
763                                         if ( odata->IsA()->InheritsFrom("TH1") ) {
764                                                 AliInfo(Form("Adding %s", histKey->GetName())) ;
765                                                 TH1 * hdata = static_cast<TH1*>(odata) ; 
766                                                 listTaskQAD->Add(hdata) ; 
767                                         }
768                                 }
769                         }
770                         man->Put(listDetQAD, idr, &mdr) ;
771                 }
772         }
773         return rv ; 
774 }       
775