]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQAManager.cxx
Revert change 31524: libRAliEn is not linked and TAlienCollection cannot be used
[u/mrichter/AliRoot.git] / STEER / AliQAManager.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: AliQAManager.cxx 30894 2009-02-05 13:46:48Z schutz $ */
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 // class for running the QA makers                                           //
20 //                                                                           //
21 //   AliQAManager 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 <TGrid.h>
35 #include <TGridCollection.h>
36 #include <TGridResult.h>
37 #include <TPluginManager.h>
38 #include <TROOT.h>
39 #include <TString.h>
40 #include <TSystem.h>
41
42 #include "AliCDBManager.h"
43 #include "AliCDBEntry.h"
44 #include "AliCDBId.h"
45 #include "AliCDBMetaData.h"
46 #include "AliCodeTimer.h"
47 #include "AliCorrQADataMakerRec.h"
48 #include "AliDetectorRecoParam.h"
49 #include "AliESDEvent.h"
50 #include "AliGeomManager.h"
51 #include "AliGlobalQADataMaker.h"
52 #include "AliHeader.h"
53 #include "AliLog.h"
54 #include "AliModule.h"
55 #include "AliQA.h"
56 #include "AliQADataMakerRec.h"
57 #include "AliQADataMakerSim.h"
58 #include "AliQAManager.h" 
59 #include "AliRawReaderDate.h"
60 #include "AliRawReaderFile.h"
61 #include "AliRawReaderRoot.h"
62 #include "AliRun.h"
63 #include "AliRunLoader.h"
64 #include "AliRunTag.h"
65
66 ClassImp(AliQAManager) 
67 AliQAManager* AliQAManager::fgQAInstance = 0x0;
68
69 //_____________________________________________________________________________
70 AliQAManager::AliQAManager() :
71   AliCDBManager(), 
72   fCurrentEvent(0),   
73   fCycleSame(kFALSE),
74   fDetectors("ALL"), 
75   fDetectorsW("ALL"), 
76   fESD(NULL), 
77   fESDTree(NULL),
78   fGAliceFileName(""), 
79   fFirstEvent(0),        
80   fMaxEvents(0),   
81   fMode(""), 
82   fNumberOfEvents(999999), 
83   fRecoParam(),
84   fRunNumber(0), 
85   fRawReader(NULL), 
86   fRawReaderDelete(kTRUE), 
87   fRunLoader(NULL), 
88   fTasks(""), 
89   fEventSpecie(AliRecoParam::kDefault)
90 {
91         // default ctor
92         fMaxEvents = fNumberOfEvents ; 
93         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
94                 if (IsSelected(AliQA::GetDetName(iDet))) {
95                         fLoader[iDet]      = NULL ;
96                         fQADataMaker[iDet] = NULL ;
97                         fQACycles[iDet]    = 999999 ;
98       fQAWriteExpert[iDet] = kTRUE ;
99                 }
100         }       
101 }
102
103 //_____________________________________________________________________________
104 AliQAManager::AliQAManager(const Char_t * mode, const Char_t* gAliceFilename) :
105   AliCDBManager(), 
106   fCurrentEvent(0),  
107         fCycleSame(kFALSE),
108         fDetectors("ALL"), 
109         fDetectorsW("ALL"), 
110         fESD(NULL), 
111         fESDTree(NULL),
112         fGAliceFileName(gAliceFilename), 
113         fFirstEvent(0),        
114         fMaxEvents(0),   
115   fMode(mode), 
116         fNumberOfEvents(999999), 
117   fRecoParam(),
118         fRunNumber(0), 
119         fRawReader(NULL), 
120         fRawReaderDelete(kTRUE), 
121         fRunLoader(NULL), 
122   fTasks(""), 
123   fEventSpecie(AliRecoParam::kDefault)
124 {
125         // default ctor
126         fMaxEvents = fNumberOfEvents ; 
127         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
128                 if (IsSelected(AliQA::GetDetName(iDet))) {
129                         fLoader[iDet]      = NULL ;
130                         fQADataMaker[iDet] = NULL ;
131                         fQACycles[iDet]    = 999999 ;
132       fQAWriteExpert[iDet] = kTRUE ;
133     }
134   }
135 }
136
137 //_____________________________________________________________________________
138 AliQAManager::AliQAManager(const AliQAManager & qas) : 
139   AliCDBManager(), 
140         fCurrentEvent(qas.fCurrentEvent),  
141         fCycleSame(kFALSE),
142         fDetectors(qas.fDetectors), 
143         fDetectorsW(qas.fDetectorsW), 
144         fESD(NULL), 
145         fESDTree(NULL), 
146         fGAliceFileName(qas.fGAliceFileName), 
147         fFirstEvent(qas.fFirstEvent),        
148         fMaxEvents(qas.fMaxEvents),    
149         fMode(qas.fMode), 
150         fNumberOfEvents(qas.fNumberOfEvents), 
151         fRecoParam(),           
152         fRunNumber(qas.fRunNumber), 
153         fRawReader(NULL), 
154         fRawReaderDelete(kTRUE), 
155         fRunLoader(NULL), 
156   fTasks(qas.fTasks),
157   fEventSpecie(qas.fEventSpecie)
158 {
159         // cpy ctor
160         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
161                 fLoader[iDet]         = qas.fLoader[iDet] ;
162                 fQADataMaker[iDet]    = qas.fQADataMaker[iDet] ;
163                 fQACycles[iDet]       = qas.fQACycles[iDet] ;   
164     fQAWriteExpert[iDet] = qas.fQAWriteExpert[iDet] ;
165   }
166 }
167
168 //_____________________________________________________________________________
169 AliQAManager & AliQAManager::operator = (const AliQAManager & qas) 
170 {
171         // assignment operator
172   this->~AliQAManager() ;
173   new(this) AliQAManager(qas) ;
174   return *this ;
175 }
176
177 //_____________________________________________________________________________
178 AliQAManager::~AliQAManager() 
179 {
180         // dtor
181         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
182                 if (IsSelected(AliQA::GetDetName(iDet))) {
183                   fLoader[iDet] = NULL;
184                   if (fQADataMaker[iDet]) {
185                           (fQADataMaker[iDet])->Finish() ; 
186                                 delete fQADataMaker[iDet] ;
187                   }
188                 }
189         }
190         if (fRawReaderDelete) { 
191                 fRunLoader = NULL ;
192                 delete fRawReader ;
193                 fRawReader = NULL ;
194         }
195 }
196
197 //_____________________________________________________________________________
198 Bool_t AliQAManager::DoIt(const AliQA::TASKINDEX_t taskIndex)
199 {
200         // Runs all the QA data Maker for every detector
201                 
202         Bool_t rv = kFALSE ;
203     // Fill QA data in event loop 
204         for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
205                 fCurrentEvent++ ; 
206                 // Get the event
207                 if ( iEvent%10 == 0  ) 
208                         AliInfo(Form("processing event %d", iEvent));
209                 if ( taskIndex == AliQA::kRAWS ) {
210                         if ( !fRawReader->NextEvent() )
211                                 break ;
212                 } else if ( taskIndex == AliQA::kESDS ) {
213                         if ( fESDTree->GetEntry(iEvent) == 0 )
214                                 break ;
215                 } else {
216                         if ( fRunLoader->GetEvent(iEvent) != 0 )
217                                 break ;
218                 }
219                 // loop  over active loaders
220                 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
221                         if (IsSelected(AliQA::GetDetName(iDet))) {
222                                 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
223                                 if (!qadm) continue; // This detector doesn't have any QA (for example, HLT)
224                                 if ( qadm->IsCycleDone() ) {
225           qadm->EndOfCycle(taskIndex) ;
226                                 }
227                                 TTree * data = NULL ; 
228                                 AliLoader* loader = GetLoader(qadm->GetUniqueID());
229                                 switch (taskIndex) {
230                                         case AliQA::kNULLTASKINDEX : 
231                                                 break ; 
232                                         case AliQA::kRAWS :
233                                                 qadm->Exec(taskIndex, fRawReader) ; 
234                                                 break ; 
235                                         case AliQA::kHITS :
236             if( loader ) {
237               loader->LoadHits() ; 
238               data = loader->TreeH() ; 
239               if ( ! data ) {
240                 AliWarning(Form(" Hit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
241                 break ; 
242               } 
243             } 
244             qadm->Exec(taskIndex, data) ;
245                                                 break ;
246                                                 case AliQA::kSDIGITS :
247             if( loader ) {      
248               loader->LoadSDigits() ; 
249               data = loader->TreeS() ; 
250               if ( ! data ) {
251                 AliWarning(Form(" SDigit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
252                 break ; 
253               } 
254             }
255             qadm->Exec(taskIndex, data) ; 
256                                                 break; 
257                                                 case AliQA::kDIGITS :
258             if( loader ) {      
259               loader->LoadDigits() ; 
260               data = loader->TreeD() ; 
261               if ( ! data ) {
262                 AliWarning(Form(" Digit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
263                 break ; 
264               } 
265             }
266             qadm->Exec(taskIndex, data) ;
267                                                 break; 
268                                                 case AliQA::kRECPOINTS :
269             if( loader ) {      
270               loader->LoadRecPoints() ; 
271               data = loader->TreeR() ; 
272               if (!data) {
273                 AliWarning(Form("RecPoints not found for %s", AliQA::GetDetName(iDet))) ; 
274                 break ; 
275               } 
276             }
277             qadm->Exec(taskIndex, data) ; 
278             break; 
279                                                 case AliQA::kTRACKSEGMENTS :
280                                                 break; 
281                                                 case AliQA::kRECPARTICLES :
282                                                 break; 
283                                                 case AliQA::kESDS :
284                                                 qadm->Exec(taskIndex, fESD) ;
285                                                 break; 
286                                                 case AliQA::kNTASKINDEX :
287                                                 break; 
288                                 } //task switch
289                         }
290                 } // detector loop
291     Increment() ; 
292         } // event loop 
293         // Save QA data for all detectors
294         rv = Finish(taskIndex) ;
295         
296         if ( taskIndex == AliQA::kRAWS ) 
297                 fRawReader->RewindEvents() ;
298
299         return rv ; 
300 }
301
302 //_____________________________________________________________________________
303 Bool_t AliQAManager::Finish(const AliQA::TASKINDEX_t taskIndex) 
304 {
305         // write output to file for all detectors
306         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
307                 if (IsSelected(AliQA::GetDetName(iDet))) {
308                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
309                         if (qadm) 
310         qadm->EndOfCycle(taskIndex) ;
311                 }
312         }
313         return kTRUE ; 
314 }
315
316 //_____________________________________________________________________________
317 TObjArray * AliQAManager::GetFromOCDB(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, const char * year) const 
318 {
319         // Retrieve the list of QA data for a given detector and a given task 
320         TObjArray * rv = NULL ;
321         if ( !strlen(AliQA::GetQARefStorage()) ) { 
322                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
323                 return NULL ; 
324         }       
325         if ( ! IsDefaultStorageSet() ) {
326                 TString tmp(AliQA::GetQARefDefaultStorage()) ; 
327                 tmp.Append(year) ; 
328                 tmp.Append("/") ; 
329                 Instance()->SetDefaultStorage(tmp.Data()) ;             
330                 Instance()->SetSpecificStorage(Form("%s/*", AliQA::GetQAName()), AliQA::GetQARefStorage()) ;
331         }
332         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetQAName(), AliQA::GetDetName((Int_t)det), AliQA::GetRefOCDBDirName())) ; 
333         AliInfo(Form("Retrieving reference data from %s/%s for %s", AliQA::GetQARefStorage(), detOCDBDir.Data(), AliQA::GetTaskName(task).Data())) ; 
334         AliCDBEntry* entry = QAManager()->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
335         TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
336         if ( listDetQAD ) 
337                 rv = dynamic_cast<TObjArray *>(listDetQAD->FindObject(AliQA::GetTaskName(task))) ; 
338         return rv ; 
339 }
340
341 //_____________________________________________________________________________
342 AliLoader * AliQAManager::GetLoader(Int_t iDet)
343 {
344         // get the loader for a detector
345
346         if ( !fRunLoader || iDet == AliQA::kCORR) 
347                 return NULL ; 
348         
349         TString detName = AliQA::GetDetName(iDet) ;
350     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
351         if (fLoader[iDet]) 
352                 return fLoader[iDet] ;
353         
354         // load the QA data maker object
355         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
356         TString loaderName = "Ali" + detName + "Loader" ;
357
358         AliLoader * loader = NULL ;
359         // first check if a plugin is defined for the quality assurance data maker
360         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
361         // if not, add a plugin for it
362         if (!pluginHandler) {
363                 AliDebug(1, Form("defining plugin for %s", loaderName.Data())) ;
364                 TString libs = gSystem->GetLibraries() ;
365                 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
366                         pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
367                 } else {
368                         pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
369                 }
370                 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
371         }
372         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
373                 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
374         }
375         if (loader) 
376                 fLoader[iDet] = loader ;
377         return loader ;
378 }
379
380 //_____________________________________________________________________________
381 AliQA * AliQAManager::GetQA(UInt_t run, UInt_t evt) 
382 {
383 // retrieves the QA object stored in a file named "Run{run}.Event{evt}_1.ESD.tag.root"  
384   char * fileName = Form("Run%d.Event%d_1.ESD.tag.root", run, evt) ; 
385   TFile * tagFile = TFile::Open(fileName) ;
386   if ( !tagFile ) {
387     AliError(Form("File %s not found", fileName)) ;
388     return NULL ; 
389   }
390   TTree * tagTree = dynamic_cast<TTree *>(tagFile->Get("T")) ; 
391   if ( !tagTree ) {
392     AliError(Form("Tree T not found in %s", fileName)) ; 
393     tagFile->Close() ; 
394     return NULL ; 
395   }
396   AliRunTag * tag = new AliRunTag ; 
397   tagTree->SetBranchAddress("AliTAG", &tag) ; 
398   tagTree->GetEntry(evt) ; 
399   AliQA * qa = AliQA::Instance(tag->GetQALength(), tag->GetQA(), tag->GetESLength(), tag->GetEventSpecies()) ; 
400   tagFile->Close() ; 
401   return qa ; 
402 }
403
404 //_____________________________________________________________________________
405 AliQADataMaker * AliQAManager::GetQADataMaker(const Int_t iDet) 
406 {
407         // get the quality assurance data maker for a detector
408         
409         if (fQADataMaker[iDet]) {
410     fQADataMaker[iDet]->SetEventSpecie(fEventSpecie) ; 
411                 return fQADataMaker[iDet] ;
412   }
413         
414         AliQADataMaker * qadm = NULL ;
415         
416         if (iDet == AliQA::kGLOBAL) { //Global QA
417                 qadm = new AliGlobalQADataMaker();
418                 qadm->SetName(AliQA::GetDetName(iDet));
419                 qadm->SetUniqueID(iDet);
420                 fQADataMaker[iDet] = qadm;
421     qadm->SetEventSpecie(fEventSpecie) ; 
422                 return qadm;
423         }
424
425         if (iDet == AliQA::kCORR) { //the data maker for correlations among detectors
426     qadm = new AliCorrQADataMakerRec(fQADataMaker) ; 
427                 qadm->SetName(AliQA::GetDetName(iDet));
428                 qadm->SetUniqueID(iDet);
429                 fQADataMaker[iDet] = qadm;
430     qadm->SetEventSpecie(fEventSpecie) ; 
431                 return qadm;
432   }
433
434         // load the QA data maker object
435         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
436         TString detName = AliQA::GetDetName(iDet) ;
437         TString tmp(fMode) ; 
438         if (tmp.Contains("sim")) 
439                 tmp.ReplaceAll("s", "S") ; 
440         else if (tmp.Contains("rec")) 
441                 tmp.ReplaceAll("r", "R") ; 
442
443         TString qadmName = "Ali" + detName + "QADataMaker" + tmp ;
444
445         // first check if a plugin is defined for the quality assurance data maker
446         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
447         // if not, add a plugin for it
448         if (!pluginHandler) {
449                 AliDebug(1, Form("defining plugin for %s", qadmName.Data())) ;
450                 TString libs = gSystem->GetLibraries() ;
451                 if (libs.Contains("lib" + detName + fMode + ".so") || (gSystem->Load("lib" + detName + fMode + ".so") >= 0)) {
452                         pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
453                 } else {
454                         pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
455                 }
456                 pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
457         }
458         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
459                 qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
460         }
461         if (qadm) {
462                 qadm->SetName(AliQA::GetDetName(iDet));
463                 qadm->SetUniqueID(iDet);
464                 fQADataMaker[iDet] = qadm ;
465     qadm->SetEventSpecie(fEventSpecie) ; 
466         }
467
468   return qadm ;
469 }
470
471 //_____________________________________________________________________________
472 void  AliQAManager::EndOfCycle(TObjArray * detArray) 
473 {
474         // End of cycle QADataMakers 
475         
476         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
477                 if (IsSelected(AliQA::GetDetName(iDet))) {
478                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
479                         if (!qadm) 
480                                 continue ;      
481                         // skip non active detectors
482                         if (detArray) {
483                                 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
484                                 if (!det || !det->IsActive())  
485                                         continue ;
486                         }
487                         for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
488                                 if ( fTasks.Contains(Form("%d", taskIndex)) ) 
489                                         qadm->EndOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex))) ;
490                         }
491                         qadm->Finish();
492                 }
493         }
494 }
495
496 //_____________________________________________________________________________
497 void  AliQAManager::EndOfCycle(TString detectors) 
498 {
499         // End of cycle QADataMakers 
500         
501         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
502                 if (IsSelected(AliQA::GetDetName(iDet))) {
503                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
504                         if (!qadm) 
505                                 continue ;      
506                         // skip non active detectors
507       if (!detectors.Contains(AliQA::GetDetName(iDet))) 
508         continue ;
509                 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
510                                 if ( fTasks.Contains(Form("%d", taskIndex)) ) 
511                                         qadm->EndOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex))) ;
512                         }
513                         qadm->Finish();
514                 }
515         }
516 }
517
518 //_____________________________________________________________________________
519 void AliQAManager::Increment()
520 {
521   // Increments the cycle counter for all QA Data Makers
522         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
523                 if (IsSelected(AliQA::GetDetName(iDet))) {
524                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
525                         if (qadm) 
526         qadm->Increment() ;
527     }
528   }
529 }
530   
531 //_____________________________________________________________________________
532 Bool_t AliQAManager::InitQA(const AliQA::TASKINDEX_t taskIndex, const  char * input )
533 {
534         // Initialize the event source and QA data makers
535         
536         fTasks += Form("%d", taskIndex) ; 
537
538         if (taskIndex == AliQA::kRAWS) { 
539                 if (!fRawReader) {
540                         fRawReader = AliRawReader::Create(input);
541                 }
542                 if ( ! fRawReader ) 
543                         return kFALSE ; 
544                 fRawReaderDelete = kTRUE ; 
545                 fRawReader->NextEvent() ; 
546                 fRunNumber = fRawReader->GetRunNumber() ; 
547                 SetRun(fRunNumber) ; 
548                 fRawReader->RewindEvents();
549                 fNumberOfEvents = 999999 ;
550                 if ( fMaxEvents < 0 ) 
551                         fMaxEvents = fNumberOfEvents ; 
552                 } else if (taskIndex == AliQA::kESDS) {
553                         fTasks = AliQA::GetTaskName(AliQA::kESDS) ; 
554       if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
555         TFile * esdFile = TFile::Open("AliESDs.root") ;
556         fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ; 
557         if ( !fESDTree ) {
558           AliError("esdTree not found") ; 
559           return kFALSE ; 
560         } else {
561           fESD     = new AliESDEvent() ;
562           fESD->ReadFromTree(fESDTree) ;
563           fESDTree->GetEntry(0) ; 
564           fRunNumber = fESD->GetRunNumber() ; 
565           fNumberOfEvents = fESDTree->GetEntries() ;
566           if ( fMaxEvents < 0 ) 
567             fMaxEvents = fNumberOfEvents ; 
568         }
569       } else {
570         AliError("AliESDs.root not found") ; 
571         return kFALSE ; 
572       }                 
573     } else {
574       if ( !InitRunLoader() ) { 
575         AliWarning("No Run Loader not found") ; 
576       } else {
577         fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
578         if ( fMaxEvents < 0 ) 
579           fMaxEvents = fNumberOfEvents ; 
580       }
581     }
582
583   // Get Detectors 
584   TObjArray* detArray = NULL ; 
585         if (fRunLoader) // check if RunLoader exists 
586                 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
587                         detArray = fRunLoader->GetAliRun()->Detectors() ;
588                         fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
589                 }
590
591         // Initialize all QA data makers for all detectors
592         fRunNumber = AliCDBManager::Instance()->GetRun() ; 
593         if ( !  AliGeomManager::GetGeometry() ) 
594                 AliGeomManager::LoadGeometry() ; 
595         
596         InitQADataMaker(fRunNumber, detArray) ; //, fCycleSame, kTRUE, detArray) ; 
597         return kTRUE ; 
598 }
599
600 //_____________________________________________________________________________
601 void  AliQAManager::InitQADataMaker(UInt_t run, TObjArray * detArray) 
602 {
603         // Initializes The QADataMaker for all active detectors and for all active tasks 
604         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
605                 if (IsSelected(AliQA::GetDetName(iDet))) {
606                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
607                         if (!qadm) {
608                                 AliError(Form("AliQADataMaker not found for %s", AliQA::GetDetName(iDet))) ; 
609                                 fDetectorsW.ReplaceAll(AliQA::GetDetName(iDet), "") ; 
610                         } else {
611         if (fQAWriteExpert[iDet])
612           qadm->SetWriteExpert() ; 
613                                 AliDebug(1, Form("Data Maker found for %s", qadm->GetName())) ; 
614                                 // skip non active detectors
615                                 if (detArray) {
616                                         AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
617                                         if (!det || !det->IsActive())  
618                                                 continue ;
619                                 }
620                                 if (fQAWriteExpert[iDet]) qadm->SetWriteExpert() ; 
621               // Set default reco params
622         Bool_t sameCycle = kFALSE ; 
623                                 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
624                                         if ( fTasks.Contains(Form("%d", taskIndex)) ) {
625                                                 qadm->Init(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
626             qadm->StartOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), run,  sameCycle) ;
627             sameCycle = kTRUE ;
628                                         }
629                                 }
630                         }
631                 }
632         }
633 }
634
635
636 //_____________________________________________________________________________
637 Bool_t AliQAManager::InitRunLoader()
638 {
639         // get or create the run loader
640         if (fRunLoader) {
641                 fCycleSame = kTRUE ; 
642         } else {
643                 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
644                         // load all base libraries to get the loader classes
645                         TString libs = gSystem->GetLibraries() ;
646                         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
647                                 if (!IsSelected(AliQA::GetDetName(iDet))) 
648                                         continue ; 
649                                 TString detName = AliQA::GetDetName(iDet) ;
650                                 if (detName == "HLT") 
651                                         continue;
652                                 if (libs.Contains("lib" + detName + "base.so")) 
653                                         continue;
654                                 gSystem->Load("lib" + detName + "base.so");
655                         }
656                         fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
657                         if (!fRunLoader) {
658                                 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
659                                 return kFALSE;
660                         }
661                         fRunLoader->CdGAFile();
662                         if (fRunLoader->LoadgAlice() == 0) {
663                                 gAlice = fRunLoader->GetAliRun();
664                         }
665
666                         if (!gAlice) {
667                                 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
668                                 return kFALSE;
669                         }
670
671                 } else {               // galice.root does not exist
672                         AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
673                         return kFALSE;
674                 }
675         }
676
677         if (!fRunNumber) { 
678                 fRunLoader->LoadHeader();
679                 fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
680         }
681         return kTRUE;
682 }
683
684 //_____________________________________________________________________________
685 Bool_t AliQAManager::IsSelected(const char * det) 
686 {
687   // check whether detName is contained in detectors
688         // if yes, it is removed from detectors
689         
690         Bool_t rv = kFALSE;
691         const TString detName(det) ;
692   // always activates Correlation
693   if ( detName.Contains(AliQA::GetDetName(AliQA::kCORR))) {
694     rv = kTRUE ; 
695   } else {
696     // check if all detectors are selected
697     if (fDetectors.Contains("ALL")) {
698       fDetectors = "ALL";
699       rv = kTRUE;
700     } else if ((fDetectors.CompareTo(detName) == 0) ||
701                fDetectors.BeginsWith(detName+" ") ||
702                fDetectors.EndsWith(" "+detName) ||
703                fDetectors.Contains(" "+detName+" ")) {
704       rv = kTRUE;
705     }
706   }
707         return rv ;
708 }
709
710 //_____________________________________________________________________________
711 Bool_t AliQAManager::Merge(const Int_t runNumber) const
712 {
713         // Merge data from all the cycles from all detectors in one single file per run
714         // Merge the QA results from all the data chunks in one run 
715  Bool_t rv = MergeData(runNumber) ; 
716  rv *= MergeResults(runNumber) ;
717  return rv ; 
718 }
719         
720 //______________________________________________________________________
721 Bool_t AliQAManager::MergeXML(const char * collectionFile, const char * subFile, const char * outFile) 
722 {
723   // merges files listed in a xml collection 
724   // usage Merge(collection, outputFile))
725   //              collection: is a xml collection  
726   
727   Bool_t rv = kFALSE ; 
728   
729   if ( strstr(collectionFile, ".xml") == 0 ) {
730     AliError("Input collection file must be an \".xml\" file\n") ; 
731     return kFALSE ; 
732   }
733     
734  if ( !gGrid ) 
735    TGrid::Connect("alien://"); 
736  if ( !gGrid ) 
737    return kFALSE ; 
738  
739   // Open the file collection 
740   printf("*** Create Collection       ***\n");
741   printf("***  Wk-Dir = |%s|             \n",gSystem->WorkingDirectory());
742   printf("***  Coll   = |%s|             \n",collectionFile);                   
743   
744   TGridCollection * collection = (TGridCollection*)gROOT->ProcessLine(Form("TAlienCollection::Open(\"%s\", 0)",collectionFile));
745   TGridResult* result = collection->GetGridResult("", 0, 0);
746   
747   Int_t index = 0  ;
748   const char * turl ;
749   TFileMerger merger ; 
750   if (!outFile) {
751     TString tempo(collectionFile) ; 
752     if ( subFile) 
753       tempo.ReplaceAll(".xml", subFile) ; 
754     else 
755       tempo.ReplaceAll(".xml", "_Merged.root") ; 
756     outFile = tempo.Data() ; 
757   }
758   merger.OutputFile(outFile) ; 
759   
760   while ( (turl = result->GetKey(index, "turl")) ) {
761     char * file ;
762     if ( subFile )
763       file = Form("%s#%s", turl, subFile) ; 
764     else 
765       file = Form("%s", turl) ; 
766     
767     AliInfo(Form("%s\n", file)) ; 
768     merger.AddFile(file) ; 
769     index++ ;  
770   }
771   
772   if (index) 
773     merger.Merge() ; 
774   
775   AliInfo(Form("Files merged into %s\n", outFile)) ;
776   
777   rv = kFALSE;
778   return rv ;
779 }
780
781 //_____________________________________________________________________________
782 Bool_t AliQAManager::MergeData(const Int_t runNumber) const
783 {
784         // Merge all the cycles from all detectors in one single file per run
785         TString cmd ;
786         if (runNumber != -1) 
787                 cmd = Form(".! ls *%s*.%d.root > tempo.txt", AliQA::GetQADataFileName(), runNumber) ; 
788         else 
789                 cmd = Form(".! ls *%s*.*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
790         gROOT->ProcessLine(cmd.Data()) ;
791         ifstream in("tempo.txt") ; 
792         const Int_t runMax = 10 ;  
793         TString file[AliQA::kNDET*runMax] ;
794         
795         Int_t index = 0 ; 
796         while ( 1 ) {
797                 in >> file[index] ; 
798                 if ( !in.good() ) 
799                         break ; 
800                 AliInfo(Form("index = %d file = %s", index, (file[index]).Data())) ; 
801                 index++ ;
802         }
803         
804         if ( index == 0 ) { 
805                 AliError(Form("run number %d not found", runNumber)) ; 
806                 return kFALSE ; 
807         }
808
809   TFileMerger merger ; 
810   TString outFileName ;
811   if (runNumber != -1) 
812     outFileName = Form("Merged.%s.Data.%d.root",AliQA::GetQADataFileName(),runNumber); 
813   else 
814     outFileName = Form("Merged.%s.Data.root",AliQA::GetQADataFileName()); 
815   merger.OutputFile(outFileName.Data()) ; 
816   for (Int_t ifile = 0 ; ifile < index-1 ; ifile++) {
817     TString pattern(Form("%s.%d.", AliQA::GetQADataFileName(), runNumber)); 
818     TString tmp(file[ifile]) ; 
819     if (tmp.Contains(pattern)) {
820       merger.AddFile(tmp) ; 
821     }
822         }
823   merger.Merge() ; 
824         return kTRUE ; 
825 }
826
827 //_____________________________________________________________________________
828 Bool_t AliQAManager::MergeResults(const Int_t runNumber) const
829 {
830         // Merge the QA result from all the data chunks in a run 
831         TString cmd ;
832         cmd = Form(".! ls %s*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
833         gROOT->ProcessLine(cmd.Data()) ;
834         ifstream in("tempo.txt") ; 
835         const Int_t chunkMax = 100 ;  
836         TString fileList[chunkMax] ;
837         
838         Int_t index = 0 ; 
839         while ( 1 ) {
840                 TString file ; 
841                 in >> fileList[index] ; 
842                 if ( !in.good() ) 
843                         break ; 
844                 AliInfo(Form("index = %d file = %s", index, (fileList[index].Data()))) ; 
845                 index++ ;
846         }
847         
848         if ( index == 0 ) { 
849                 AliError("No QA Result File found") ; 
850                 return kFALSE ; 
851         }
852         
853         TFileMerger merger ; 
854   TString outFileName ;
855   if (runNumber != -1) 
856     outFileName = Form("Merged.%s.Result.%d.root",AliQA::GetQADataFileName(),runNumber); 
857   else 
858     outFileName = Form("Merged.%s.Result.root",AliQA::GetQADataFileName()); 
859         merger.OutputFile(outFileName.Data()) ; 
860         for (Int_t ifile = 0 ; ifile < index ; ifile++) {
861                 TString file = fileList[ifile] ; 
862                 merger.AddFile(file) ; 
863         }
864         merger.Merge() ; 
865         
866         return kTRUE ; 
867 }
868
869 //_____________________________________________________________________________
870 void AliQAManager::Reset(const Bool_t sameCycle)
871 {
872         // Reset the default data members
873
874         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
875                 if (IsSelected(AliQA::GetDetName(iDet))) {
876                         AliQADataMaker * qadm = GetQADataMaker(iDet);
877                         qadm->Reset();
878                 }
879         } 
880         if (fRawReaderDelete) { 
881                 delete fRawReader ;
882                 fRawReader      = NULL ;
883         }
884
885         fCycleSame      = sameCycle ; 
886         fESD            = NULL ; 
887         fESDTree        = NULL ; 
888         //fFirst          = kTRUE ;   
889         fNumberOfEvents = 999999 ;  
890 }
891
892 //_____________________________________________________________________________
893 AliQAManager * AliQAManager::QAManager(const Char_t * mode, TMap *entryCache, Int_t run) 
894 {
895   // returns AliQAManager instance (singleton)
896   
897         if (!fgQAInstance) {
898     fgQAInstance = new AliQAManager(mode) ;  
899     if (!entryCache)
900                   fgQAInstance->Init();
901                 else
902                   fgQAInstance->InitFromCache(entryCache,run);
903   }
904         return fgQAInstance;
905 }
906
907 //_____________________________________________________________________________
908 TString AliQAManager::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle) 
909 {
910         //Runs all the QA data Maker for Raws only
911         
912         fCycleSame       = sameCycle ;
913         fRawReader       = rawReader ;
914         fDetectors       = detectors ; 
915         fDetectorsW      = detectors ;  
916         
917         AliCDBManager* man = AliCDBManager::Instance() ; 
918
919         if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
920                 rawReader->NextEvent() ; 
921                 man->SetRun(fRawReader->GetRunNumber()) ;
922                 rawReader->RewindEvents() ;
923         }       
924         
925         if (!fCycleSame) 
926     if ( !InitQA(AliQA::kRAWS) ) 
927       return kFALSE ; 
928   fRawReaderDelete = kFALSE ; 
929
930         DoIt(AliQA::kRAWS) ; 
931         return  fDetectorsW ;
932 }
933
934 //_____________________________________________________________________________
935 TString AliQAManager::Run(const char * detectors, const char * fileName, const Bool_t sameCycle) 
936 {
937         //Runs all the QA data Maker for Raws only
938
939         fCycleSame       = sameCycle ;
940         fDetectors       = detectors ; 
941         fDetectorsW      = detectors ;  
942         
943         AliCDBManager* man = AliCDBManager::Instance() ; 
944         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
945                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
946                 if ( ! rl ) {
947                         AliFatal("galice.root file not found in current directory") ; 
948                 } else {
949                         rl->CdGAFile() ; 
950                         rl->LoadgAlice() ;
951                         if ( ! rl->GetAliRun() ) {
952                                 AliFatal("AliRun not found in galice.root") ;
953                         } else {
954                                 rl->LoadHeader() ;
955                                 man->SetRun(rl->GetHeader()->GetRun());
956                         }
957                 }
958         }
959         
960         if (!fCycleSame) 
961     if ( !InitQA(AliQA::kRAWS, fileName) ) 
962       return kFALSE ; 
963         
964         DoIt(AliQA::kRAWS) ; 
965         return  fDetectorsW ;
966 }
967
968 //_____________________________________________________________________________
969 TString AliQAManager::Run(const char * detectors, const AliQA::TASKINDEX_t taskIndex, Bool_t const sameCycle, const  char * fileName ) 
970 {
971         // Runs all the QA data Maker for every detector
972         
973         fCycleSame       = sameCycle ;
974         fDetectors       = detectors ; 
975         fDetectorsW      = detectors ;          
976         
977         AliCDBManager* man = AliCDBManager::Instance() ;        
978         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
979                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
980                 if ( ! rl ) {
981                         AliFatal("galice.root file not found in current directory") ; 
982                 } else {
983                         rl->CdGAFile() ; 
984                         rl->LoadgAlice() ;
985                         if ( ! rl->GetAliRun() ) {
986                                 AliInfo("AliRun not found in galice.root") ;
987                         } else {
988                                 rl->LoadHeader() ;
989                                 man->SetRun(rl->GetHeader()->GetRun()) ;
990                         }
991                 }
992         }
993         
994
995         if ( taskIndex == AliQA::kNULLTASKINDEX) { 
996                 for (UInt_t task = 0; task < AliQA::kNTASKINDEX; task++) {
997                         if ( fTasks.Contains(Form("%d", task)) ) {
998         if (!fCycleSame)
999           if ( !InitQA(AliQA::GetTaskIndex(AliQA::GetTaskName(task)), fileName) ) 
1000             return kFALSE ;
1001         DoIt(AliQA::GetTaskIndex(AliQA::GetTaskName(task))) ;
1002                         }
1003                 }
1004         } else {
1005     if (! fCycleSame )
1006       if ( !InitQA(taskIndex, fileName) ) 
1007         return kFALSE ; 
1008       DoIt(taskIndex) ; 
1009   }             
1010         
1011         return fDetectorsW ;
1012
1013 }
1014
1015 //_____________________________________________________________________________
1016 void AliQAManager::RunOneEvent(AliRawReader * rawReader) 
1017 {
1018         //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
1019   if ( ! rawReader ) 
1020     return ; 
1021         AliCodeTimerAuto("") ;
1022   if (fTasks.Contains(Form("%d", AliQA::kRAWS))){
1023     for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1024       if (!IsSelected(AliQA::GetDetName(iDet))) 
1025         continue;
1026       AliQADataMaker *qadm = GetQADataMaker(iDet);  
1027       if (!qadm) 
1028         continue;
1029       if ( qadm->IsCycleDone() ) {
1030         qadm->EndOfCycle() ;
1031       }
1032       AliCodeTimerStart(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet))); 
1033       qadm->SetEventSpecie(fEventSpecie) ; 
1034                         qadm->Exec(AliQA::kRAWS, rawReader) ;
1035       AliCodeTimerStop(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1036                 }
1037   }
1038 }
1039
1040 //_____________________________________________________________________________
1041 void AliQAManager::RunOneEvent(AliESDEvent *& esd) 
1042 {
1043         //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1044         
1045   AliCodeTimerAuto("") ;
1046   if (fTasks.Contains(Form("%d", AliQA::kESDS))) {
1047     for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1048       if (!IsSelected(AliQA::GetDetName(iDet))) 
1049         continue;
1050       AliQADataMaker *qadm = GetQADataMaker(iDet);  
1051       if (!qadm) 
1052         continue;
1053       if ( qadm->IsCycleDone() ) {
1054         qadm->EndOfCycle() ;
1055       }
1056       AliCodeTimerStart(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1057                         qadm->Exec(AliQA::kESDS, esd) ;
1058       AliCodeTimerStop(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1059                 }
1060         }
1061 }
1062
1063 //_____________________________________________________________________________
1064 void AliQAManager::RunOneEventInOneDetector(Int_t det, TTree * tree) 
1065 {
1066         // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1067         AliCodeTimerAuto("") ;
1068   if (fTasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1069     if (IsSelected(AliQA::GetDetName(det))) {
1070       AliQADataMaker *qadm = GetQADataMaker(det);  
1071       if (qadm) { 
1072         if ( qadm->IsCycleDone() ) {
1073           qadm->EndOfCycle() ;
1074         }
1075         AliCodeTimerStart(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
1076         qadm->Exec(AliQA::kRECPOINTS, tree) ;
1077         AliCodeTimerStop(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
1078       }
1079     }
1080   }
1081 }
1082
1083 //_____________________________________________________________________________
1084 Bool_t AliQAManager::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const char * year, const char * detectors) const
1085 {
1086         // take the locasl QA data merge into a single file and save in OCDB 
1087         Bool_t rv = kTRUE ; 
1088         TString tmp(AliQA::GetQARefStorage()) ; 
1089         if ( tmp.IsNull() ) { 
1090                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
1091                 return kFALSE ; 
1092         }
1093         if ( !(tmp.Contains(AliQA::GetLabLocalOCDB()) || tmp.Contains(AliQA::GetLabAliEnOCDB())) ) {
1094                 AliError(Form("%s is a wrong storage, use %s or %s", AliQA::GetQARefStorage(), AliQA::GetLabLocalOCDB().Data(), AliQA::GetLabAliEnOCDB().Data())) ; 
1095                 return kFALSE ; 
1096         }
1097         TString sdet(detectors) ; 
1098         sdet.ToUpper() ;
1099         TFile * inputFile ; 
1100         if ( sdet.Contains("ALL") ) {
1101                 rv = Merge(runNumber) ; 
1102                 if ( ! rv )
1103                         return kFALSE ; 
1104                 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQA::GetQADataFileName(), runNumber)) ; 
1105                 inputFile = TFile::Open(inputFileName.Data()) ; 
1106                 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ; 
1107         } else {
1108                 for (Int_t index = 0; index < AliQA::kNDET; index++) {
1109                         if (sdet.Contains(AliQA::GetDetName(index))) {
1110                                 TString inputFileName(Form("%s.%s.%d.root", AliQA::GetDetName(index), AliQA::GetQADataFileName(), runNumber)) ; 
1111                                 inputFile = TFile::Open(inputFileName.Data()) ;                         
1112                                 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ; 
1113                         }
1114                 }
1115         }
1116         return rv ; 
1117 }
1118
1119 //_____________________________________________________________________________
1120 Bool_t AliQAManager::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year, AliRecoParam::EventSpecie_t es) const
1121 {
1122         // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1123         Bool_t rv = kTRUE ;
1124         AliInfo(Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQA::GetQARefStorage())) ; 
1125         if ( ! IsDefaultStorageSet() ) {
1126                 TString tmp( AliQA::GetQARefStorage() ) ; 
1127                 if ( tmp.Contains(AliQA::GetLabLocalOCDB()) ) 
1128                         Instance()->SetDefaultStorage(AliQA::GetQARefStorage()) ;
1129                 else {
1130                         TString tmp1(AliQA::GetQARefDefaultStorage()) ; 
1131                         tmp1.Append(year) ; 
1132                         tmp1.Append("?user=alidaq") ; 
1133                         Instance()->SetDefaultStorage(tmp1.Data()) ; 
1134                 }
1135         }
1136         Instance()->SetSpecificStorage("*", AliQA::GetQARefStorage()) ; 
1137         if(GetRun() < 0) 
1138                 Instance()->SetRun(runNumber);
1139
1140         AliCDBMetaData mdr ;
1141         mdr.SetResponsible("yves schutz");
1142
1143         for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++) {
1144                 TDirectory * detDir = inputFile->GetDirectory(AliQA::GetDetName(detIndex)) ; 
1145                 if ( detDir ) {
1146                         AliInfo(Form("Entering %s", detDir->GetName())) ;
1147       AliQA::SetQARefDataDirName(es) ;
1148                         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetDetName(detIndex), AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName())) ; 
1149                         AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity())  ;
1150                         TList * listDetQAD = new TList() ;
1151                         TString listName(Form("%s QA data Reference", AliQA::GetDetName(detIndex))) ; 
1152                         mdr.SetComment(Form("%s QA stuff", AliQA::GetDetName(detIndex)));
1153                         listDetQAD->SetName(listName) ; 
1154                         TList * taskList = detDir->GetListOfKeys() ; 
1155                         TIter nextTask(taskList) ; 
1156                         TKey * taskKey ; 
1157                         while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
1158                                 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ; 
1159         TDirectory * esDir   = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ; 
1160                                 AliInfo(Form("Saving %s", esDir->GetName())) ; 
1161                                 TObjArray * listTaskQAD = new TObjArray(100) ; 
1162                                 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1163                                 listDetQAD->Add(listTaskQAD) ; 
1164                                 TList * histList = esDir->GetListOfKeys() ; 
1165                                 TIter nextHist(histList) ; 
1166                                 TKey * histKey ; 
1167                                 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
1168                                         TObject * odata = esDir->Get(histKey->GetName()) ; 
1169                                         if ( !odata ) {
1170                                                 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1171                                         } else {
1172             if ( AliQA::GetExpert() == histKey->GetName() ) {
1173               TDirectory * expertDir   = esDir->GetDirectory(histKey->GetName()) ; 
1174               TList * expertHistList = expertDir->GetListOfKeys() ; 
1175               TIter nextExpertHist(expertHistList) ; 
1176               TKey * expertHistKey ; 
1177               while ( (expertHistKey = dynamic_cast<TKey*>(nextExpertHist())) ) {
1178                 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ; 
1179                 if ( !expertOdata ) {
1180                   AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1181                 } else {
1182                   AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1183                   if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
1184                     AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1185                     TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ; 
1186                     listTaskQAD->Add(hExpertdata) ; 
1187                   }                  
1188                 }                
1189               }
1190             }
1191                                                 AliInfo(Form("Adding %s", histKey->GetName())) ;
1192                                                 if ( odata->IsA()->InheritsFrom("TH1") ) {
1193                                                         AliInfo(Form("Adding %s", histKey->GetName())) ;
1194                                                         TH1 * hdata = static_cast<TH1*>(odata) ; 
1195                                                         listTaskQAD->Add(hdata) ; 
1196                                                 }
1197                                         }
1198                                 }
1199                         }
1200                         Instance()->Put(listDetQAD, idr, &mdr) ;
1201                 }
1202         }
1203         return rv ; 
1204 }       
1205
1206 //_____________________________________________________________________________
1207 void AliQAManager::SetEventSpecie(AliRecoParam::EventSpecie_t es) 
1208 {
1209   // set the current event specie and inform AliQA that this event specie has been encountered
1210   fEventSpecie = es ;
1211   AliQA::Instance()->SetEventSpecie(es) ; 
1212 }
1213
1214 //_____________________________________________________________________________
1215 void AliQAManager::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par) 
1216 {
1217   // Set custom reconstruction parameters for a given detector
1218   // Single set of parameters for all the events
1219   GetQADataMaker(det)->SetRecoParam(par) ; 
1220 }
1221
1222