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