]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQAManager.cxx
fixed a problem with writing the expert QA Data. Added a method in simulation and...
[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                 }
99         }       
100   SetWriteExpert() ; 
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     }
133   }
134   SetWriteExpert() ; 
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 %d", qadm->GetName(), qadm->WriteExpert())) ; 
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               // Set default reco params
621         Bool_t sameCycle = kFALSE ; 
622                                 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
623                                         if ( fTasks.Contains(Form("%d", taskIndex)) ) {
624                                                 qadm->Init(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
625             qadm->StartOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), run,  sameCycle) ;
626             sameCycle = kTRUE ;
627                                         }
628                                 }
629                         }
630                 }
631         }
632 }
633
634
635 //_____________________________________________________________________________
636 Bool_t AliQAManager::InitRunLoader()
637 {
638         // get or create the run loader
639         if (fRunLoader) {
640                 fCycleSame = kTRUE ; 
641         } else {
642                 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
643                         // load all base libraries to get the loader classes
644                         TString libs = gSystem->GetLibraries() ;
645                         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
646                                 if (!IsSelected(AliQA::GetDetName(iDet))) 
647                                         continue ; 
648                                 TString detName = AliQA::GetDetName(iDet) ;
649                                 if (detName == "HLT") 
650                                         continue;
651                                 if (libs.Contains("lib" + detName + "base.so")) 
652                                         continue;
653                                 gSystem->Load("lib" + detName + "base.so");
654                         }
655                         fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
656                         if (!fRunLoader) {
657                                 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
658                                 return kFALSE;
659                         }
660                         fRunLoader->CdGAFile();
661                         if (fRunLoader->LoadgAlice() == 0) {
662                                 gAlice = fRunLoader->GetAliRun();
663                         }
664
665                         if (!gAlice) {
666                                 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
667                                 return kFALSE;
668                         }
669
670                 } else {               // galice.root does not exist
671                         AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
672                         return kFALSE;
673                 }
674         }
675
676         if (!fRunNumber) { 
677                 fRunLoader->LoadHeader();
678                 fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
679         }
680         return kTRUE;
681 }
682
683 //_____________________________________________________________________________
684 Bool_t AliQAManager::IsSelected(const char * det) 
685 {
686   // check whether detName is contained in detectors
687         // if yes, it is removed from detectors
688         
689         Bool_t rv = kFALSE;
690         const TString detName(det) ;
691   // always activates Correlation
692   if ( detName.Contains(AliQA::GetDetName(AliQA::kCORR))) {
693     rv = kTRUE ; 
694   } else {
695     // check if all detectors are selected
696     if (fDetectors.Contains("ALL")) {
697       fDetectors = "ALL";
698       rv = kTRUE;
699     } else if ((fDetectors.CompareTo(detName) == 0) ||
700                fDetectors.BeginsWith(detName+" ") ||
701                fDetectors.EndsWith(" "+detName) ||
702                fDetectors.Contains(" "+detName+" ")) {
703       rv = kTRUE;
704     }
705   }
706         return rv ;
707 }
708
709 //_____________________________________________________________________________
710 Bool_t AliQAManager::Merge(Int_t runNumber) const
711 {
712         // Merge data from all detectors from a given run in one single file 
713         // Merge the QA results from all the data chunks in one run 
714  if ( runNumber == -1)
715    runNumber = fRunNumber ; 
716  Bool_t rv = MergeData(runNumber) ; 
717  //rv *= MergeResults(runNumber) ; // not needed for the time being
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 = (TGridCollection*)gROOT->ProcessLine(Form("TAlienCollection::Open(\"%s\")",collectionFile));
746   TGridResult* result = collection->GetGridResult("", 0, 0);
747   
748   Int_t index = 0  ;
749   const char * turl ;
750   TFileMerger merger(kFALSE) ; 
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 QA data from all detectors for a given run in one single file 
786   
787   TFileMerger merger ; 
788   TString outFileName = Form("Merged.%s.Data.root",AliQA::GetQADataFileName()) ;
789   merger.OutputFile(outFileName.Data()) ; 
790   for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
791     char * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQA::GetDetName(iDet), AliQA::GetQADataFileName(), runNumber)); 
792     if (file) 
793       merger.AddFile(file) ; 
794   }
795   merger.Merge() ; 
796         return kTRUE ; 
797 }
798
799 //_____________________________________________________________________________
800 Bool_t AliQAManager::MergeResults(const Int_t runNumber) const
801 {
802         // Merge the QA result from all the data chunks in a run 
803   // to be revised whwn it will be used (see MergeData)
804         TString cmd ;
805         cmd = Form(".! ls %s*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
806         gROOT->ProcessLine(cmd.Data()) ;
807         ifstream in("tempo.txt") ; 
808         const Int_t chunkMax = 100 ;  
809         TString fileList[chunkMax] ;
810         
811         Int_t index = 0 ; 
812         while ( 1 ) {
813                 TString file ; 
814                 in >> fileList[index] ; 
815                 if ( !in.good() ) 
816                         break ; 
817                 AliInfo(Form("index = %d file = %s", index, (fileList[index].Data()))) ; 
818                 index++ ;
819         }
820         
821         if ( index == 0 ) { 
822                 AliError("No QA Result File found") ; 
823                 return kFALSE ; 
824         }
825         
826         TFileMerger merger ; 
827   TString outFileName ;
828   if (runNumber != -1) 
829     outFileName = Form("Merged.%s.Result.%d.root",AliQA::GetQADataFileName(),runNumber); 
830   else 
831     outFileName = Form("Merged.%s.Result.root",AliQA::GetQADataFileName()); 
832         merger.OutputFile(outFileName.Data()) ; 
833         for (Int_t ifile = 0 ; ifile < index ; ifile++) {
834                 TString file = fileList[ifile] ; 
835                 merger.AddFile(file) ; 
836         }
837         merger.Merge() ; 
838         
839         return kTRUE ; 
840 }
841
842 //_____________________________________________________________________________
843 void AliQAManager::Reset(const Bool_t sameCycle)
844 {
845         // Reset the default data members
846
847         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
848                 if (IsSelected(AliQA::GetDetName(iDet))) {
849                         AliQADataMaker * qadm = GetQADataMaker(iDet);
850                         qadm->Reset();
851                 }
852         } 
853         if (fRawReaderDelete) { 
854                 delete fRawReader ;
855                 fRawReader      = NULL ;
856         }
857
858         fCycleSame      = sameCycle ; 
859         fESD            = NULL ; 
860         fESDTree        = NULL ; 
861         //fFirst          = kTRUE ;   
862         fNumberOfEvents = 999999 ;  
863 }
864
865 //_____________________________________________________________________________
866 AliQAManager * AliQAManager::QAManager(const Char_t * mode, TMap *entryCache, Int_t run) 
867 {
868   // returns AliQAManager instance (singleton)
869   
870         if (!fgQAInstance) {
871     fgQAInstance = new AliQAManager(mode) ;  
872     if (!entryCache)
873                   fgQAInstance->Init();
874                 else
875                   fgQAInstance->InitFromCache(entryCache,run);
876   }
877         return fgQAInstance;
878 }
879
880 //_____________________________________________________________________________
881 TString AliQAManager::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle) 
882 {
883         //Runs all the QA data Maker for Raws only
884         
885         fCycleSame       = sameCycle ;
886         fRawReader       = rawReader ;
887         fDetectors       = detectors ; 
888         fDetectorsW      = detectors ;  
889         
890         AliCDBManager* man = AliCDBManager::Instance() ; 
891
892         if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
893                 rawReader->NextEvent() ; 
894                 man->SetRun(fRawReader->GetRunNumber()) ;
895                 rawReader->RewindEvents() ;
896         }       
897         
898         if (!fCycleSame) 
899     if ( !InitQA(AliQA::kRAWS) ) 
900       return kFALSE ; 
901   fRawReaderDelete = kFALSE ; 
902
903         DoIt(AliQA::kRAWS) ; 
904         return  fDetectorsW ;
905 }
906
907 //_____________________________________________________________________________
908 TString AliQAManager::Run(const char * detectors, const char * fileName, const Bool_t sameCycle) 
909 {
910         //Runs all the QA data Maker for Raws only
911
912         fCycleSame       = sameCycle ;
913         fDetectors       = detectors ; 
914         fDetectorsW      = detectors ;  
915         
916         AliCDBManager* man = AliCDBManager::Instance() ; 
917         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
918                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
919                 if ( ! rl ) {
920                         AliFatal("galice.root file not found in current directory") ; 
921                 } else {
922                         rl->CdGAFile() ; 
923                         rl->LoadgAlice() ;
924                         if ( ! rl->GetAliRun() ) {
925                                 AliFatal("AliRun not found in galice.root") ;
926                         } else {
927                                 rl->LoadHeader() ;
928                                 man->SetRun(rl->GetHeader()->GetRun());
929                         }
930                 }
931         }
932         
933         if (!fCycleSame) 
934     if ( !InitQA(AliQA::kRAWS, fileName) ) 
935       return kFALSE ; 
936         
937         DoIt(AliQA::kRAWS) ; 
938         return  fDetectorsW ;
939 }
940
941 //_____________________________________________________________________________
942 TString AliQAManager::Run(const char * detectors, const AliQA::TASKINDEX_t taskIndex, Bool_t const sameCycle, const  char * fileName ) 
943 {
944         // Runs all the QA data Maker for every detector
945         
946         fCycleSame       = sameCycle ;
947         fDetectors       = detectors ; 
948         fDetectorsW      = detectors ;          
949         
950         AliCDBManager* man = AliCDBManager::Instance() ;        
951         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
952                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
953                 if ( ! rl ) {
954                         AliFatal("galice.root file not found in current directory") ; 
955                 } else {
956                         rl->CdGAFile() ; 
957                         rl->LoadgAlice() ;
958                         if ( ! rl->GetAliRun() ) {
959                                 AliInfo("AliRun not found in galice.root") ;
960                         } else {
961                                 rl->LoadHeader() ;
962                                 man->SetRun(rl->GetHeader()->GetRun()) ;
963                         }
964                 }
965         }
966         
967
968         if ( taskIndex == AliQA::kNULLTASKINDEX) { 
969                 for (UInt_t task = 0; task < AliQA::kNTASKINDEX; task++) {
970                         if ( fTasks.Contains(Form("%d", task)) ) {
971         if (!fCycleSame)
972           if ( !InitQA(AliQA::GetTaskIndex(AliQA::GetTaskName(task)), fileName) ) 
973             return kFALSE ;
974         DoIt(AliQA::GetTaskIndex(AliQA::GetTaskName(task))) ;
975                         }
976                 }
977         } else {
978     if (! fCycleSame )
979       if ( !InitQA(taskIndex, fileName) ) 
980         return kFALSE ; 
981       DoIt(taskIndex) ; 
982   }             
983         
984         return fDetectorsW ;
985
986 }
987
988 //_____________________________________________________________________________
989 void AliQAManager::RunOneEvent(AliRawReader * rawReader) 
990 {
991         //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
992   if ( ! rawReader ) 
993     return ; 
994         AliCodeTimerAuto("") ;
995   if (fTasks.Contains(Form("%d", AliQA::kRAWS))){
996     for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
997       if (!IsSelected(AliQA::GetDetName(iDet))) 
998         continue;
999       AliQADataMaker *qadm = GetQADataMaker(iDet);  
1000       if (!qadm) 
1001         continue;
1002       if ( qadm->IsCycleDone() ) {
1003         qadm->EndOfCycle() ;
1004       }
1005       AliCodeTimerStart(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet))); 
1006       qadm->SetEventSpecie(fEventSpecie) ; 
1007                         qadm->Exec(AliQA::kRAWS, rawReader) ;
1008       AliCodeTimerStop(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1009                 }
1010   }
1011 }
1012
1013 //_____________________________________________________________________________
1014 void AliQAManager::RunOneEvent(AliESDEvent *& esd) 
1015 {
1016         //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1017         
1018   AliCodeTimerAuto("") ;
1019   if (fTasks.Contains(Form("%d", AliQA::kESDS))) {
1020     for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1021       if (!IsSelected(AliQA::GetDetName(iDet))) 
1022         continue;
1023       AliQADataMaker *qadm = GetQADataMaker(iDet);  
1024       if (!qadm) 
1025         continue;
1026       if ( qadm->IsCycleDone() ) {
1027         qadm->EndOfCycle() ;
1028       }
1029       AliCodeTimerStart(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1030                         qadm->Exec(AliQA::kESDS, esd) ;
1031       AliCodeTimerStop(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1032                 }
1033         }
1034 }
1035
1036 //_____________________________________________________________________________
1037 void AliQAManager::RunOneEventInOneDetector(Int_t det, TTree * tree) 
1038 {
1039         // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1040         AliCodeTimerAuto("") ;
1041   if (fTasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1042     if (IsSelected(AliQA::GetDetName(det))) {
1043       AliQADataMaker *qadm = GetQADataMaker(det);  
1044       if (qadm) { 
1045         if ( qadm->IsCycleDone() ) {
1046           qadm->EndOfCycle() ;
1047         }
1048         AliCodeTimerStart(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
1049         qadm->Exec(AliQA::kRECPOINTS, tree) ;
1050         AliCodeTimerStop(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
1051       }
1052     }
1053   }
1054 }
1055
1056 //_____________________________________________________________________________
1057 Bool_t AliQAManager::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const char * year, const char * detectors) const
1058 {
1059         // take the locasl QA data merge into a single file and save in OCDB 
1060         Bool_t rv = kTRUE ; 
1061         TString tmp(AliQA::GetQARefStorage()) ; 
1062         if ( tmp.IsNull() ) { 
1063                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
1064                 return kFALSE ; 
1065         }
1066         if ( !(tmp.Contains(AliQA::GetLabLocalOCDB()) || tmp.Contains(AliQA::GetLabAliEnOCDB())) ) {
1067                 AliError(Form("%s is a wrong storage, use %s or %s", AliQA::GetQARefStorage(), AliQA::GetLabLocalOCDB().Data(), AliQA::GetLabAliEnOCDB().Data())) ; 
1068                 return kFALSE ; 
1069         }
1070         TString sdet(detectors) ; 
1071         sdet.ToUpper() ;
1072         TFile * inputFile ; 
1073         if ( sdet.Contains("ALL") ) {
1074                 rv = Merge(runNumber) ; 
1075                 if ( ! rv )
1076                         return kFALSE ; 
1077                 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQA::GetQADataFileName(), runNumber)) ; 
1078                 inputFile = TFile::Open(inputFileName.Data()) ; 
1079                 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ; 
1080         } else {
1081                 for (Int_t index = 0; index < AliQA::kNDET; index++) {
1082                         if (sdet.Contains(AliQA::GetDetName(index))) {
1083                                 TString inputFileName(Form("%s.%s.%d.root", AliQA::GetDetName(index), AliQA::GetQADataFileName(), runNumber)) ; 
1084                                 inputFile = TFile::Open(inputFileName.Data()) ;                         
1085                                 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ; 
1086                         }
1087                 }
1088         }
1089         return rv ; 
1090 }
1091
1092 //_____________________________________________________________________________
1093 Bool_t AliQAManager::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year, AliRecoParam::EventSpecie_t es) const
1094 {
1095         // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1096         Bool_t rv = kTRUE ;
1097         AliInfo(Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQA::GetQARefStorage())) ; 
1098         if ( ! IsDefaultStorageSet() ) {
1099                 TString tmp( AliQA::GetQARefStorage() ) ; 
1100                 if ( tmp.Contains(AliQA::GetLabLocalOCDB()) ) 
1101                         Instance()->SetDefaultStorage(AliQA::GetQARefStorage()) ;
1102                 else {
1103                         TString tmp1(AliQA::GetQARefDefaultStorage()) ; 
1104                         tmp1.Append(year) ; 
1105                         tmp1.Append("?user=alidaq") ; 
1106                         Instance()->SetDefaultStorage(tmp1.Data()) ; 
1107                 }
1108         }
1109         Instance()->SetSpecificStorage("*", AliQA::GetQARefStorage()) ; 
1110         if(GetRun() < 0) 
1111                 Instance()->SetRun(runNumber);
1112
1113         AliCDBMetaData mdr ;
1114         mdr.SetResponsible("yves schutz");
1115
1116         for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++) {
1117                 TDirectory * detDir = inputFile->GetDirectory(AliQA::GetDetName(detIndex)) ; 
1118                 if ( detDir ) {
1119                         AliInfo(Form("Entering %s", detDir->GetName())) ;
1120       AliQA::SetQARefDataDirName(es) ;
1121                         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetDetName(detIndex), AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName())) ; 
1122                         AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity())  ;
1123                         TList * listDetQAD = new TList() ;
1124                         TString listName(Form("%s QA data Reference", AliQA::GetDetName(detIndex))) ; 
1125                         mdr.SetComment(Form("%s QA stuff", AliQA::GetDetName(detIndex)));
1126                         listDetQAD->SetName(listName) ; 
1127                         TList * taskList = detDir->GetListOfKeys() ; 
1128                         TIter nextTask(taskList) ; 
1129                         TKey * taskKey ; 
1130                         while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
1131                                 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ; 
1132         TDirectory * esDir   = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ; 
1133                                 AliInfo(Form("Saving %s", esDir->GetName())) ; 
1134                                 TObjArray * listTaskQAD = new TObjArray(100) ; 
1135                                 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1136                                 listDetQAD->Add(listTaskQAD) ; 
1137                                 TList * histList = esDir->GetListOfKeys() ; 
1138                                 TIter nextHist(histList) ; 
1139                                 TKey * histKey ; 
1140                                 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
1141                                         TObject * odata = esDir->Get(histKey->GetName()) ; 
1142                                         if ( !odata ) {
1143                                                 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1144                                         } else {
1145             if ( AliQA::GetExpert() == histKey->GetName() ) {
1146               TDirectory * expertDir   = esDir->GetDirectory(histKey->GetName()) ; 
1147               TList * expertHistList = expertDir->GetListOfKeys() ; 
1148               TIter nextExpertHist(expertHistList) ; 
1149               TKey * expertHistKey ; 
1150               while ( (expertHistKey = dynamic_cast<TKey*>(nextExpertHist())) ) {
1151                 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ; 
1152                 if ( !expertOdata ) {
1153                   AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1154                 } else {
1155                   AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1156                   if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
1157                     AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1158                     TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ; 
1159                     listTaskQAD->Add(hExpertdata) ; 
1160                   }                  
1161                 }                
1162               }
1163             }
1164                                                 AliInfo(Form("Adding %s", histKey->GetName())) ;
1165                                                 if ( odata->IsA()->InheritsFrom("TH1") ) {
1166                                                         AliInfo(Form("Adding %s", histKey->GetName())) ;
1167                                                         TH1 * hdata = static_cast<TH1*>(odata) ; 
1168                                                         listTaskQAD->Add(hdata) ; 
1169                                                 }
1170                                         }
1171                                 }
1172                         }
1173                         Instance()->Put(listDetQAD, idr, &mdr) ;
1174                 }
1175         }
1176         return rv ; 
1177 }       
1178
1179 //_____________________________________________________________________________
1180 void AliQAManager::SetEventSpecie(AliRecoParam::EventSpecie_t es) 
1181 {
1182   // set the current event specie and inform AliQA that this event specie has been encountered
1183   fEventSpecie = es ;
1184   AliQA::Instance()->SetEventSpecie(es) ; 
1185 }
1186
1187 //_____________________________________________________________________________
1188 void AliQAManager::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par) 
1189 {
1190   // Set custom reconstruction parameters for a given detector
1191   // Single set of parameters for all the events
1192   GetQADataMaker(det)->SetRecoParam(par) ; 
1193 }
1194
1195 //_____________________________________________________________________________
1196 void AliQAManager::SetWriteExpert()
1197 {
1198   // enable the writing of QA expert data
1199   for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1200         if (IsSelected(AliQA::GetDetName(iDet))) 
1201       fQAWriteExpert[iDet] = kTRUE ;
1202   }
1203 }