1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliQAManager.cxx 30894 2009-02-05 13:46:48Z schutz $ */
17 ///////////////////////////////////////////////////////////////////////////////
19 // class for running the QA makers //
21 // AliQAManager qas; //
22 // qas.Run(AliQAv1::kRAWS, rawROOTFileName); //
23 // qas.Run(AliQAv1::kHITS); //
24 // qas.Run(AliQAv1::kSDIGITS); //
25 // qas.Run(AliQAv1::kDIGITS); //
26 // qas.Run(AliQAv1::kRECPOINTS); //
27 // qas.Run(AliQAv1::kESDS); //
29 ///////////////////////////////////////////////////////////////////////////////
34 #include <TFileMerger.h>
36 #include <TGridCollection.h>
37 #include <TGridResult.h>
38 #include <TPluginManager.h>
42 #include <TStopwatch.h>
44 #include "AliCDBManager.h"
45 #include "AliCDBEntry.h"
47 #include "AliCDBMetaData.h"
48 #include "AliCodeTimer.h"
49 #include "AliCorrQADataMakerRec.h"
50 #include "AliDetectorRecoParam.h"
51 #include "AliESDEvent.h"
52 #include "AliGeomManager.h"
53 #include "AliGlobalQADataMaker.h"
54 #include "AliHeader.h"
56 #include "AliModule.h"
58 #include "AliQAChecker.h"
59 #include "AliQACheckerBase.h"
60 #include "AliQADataMakerRec.h"
61 #include "AliQADataMakerSim.h"
62 #include "AliQAManager.h"
63 #include "AliRawReaderDate.h"
64 #include "AliRawReaderFile.h"
65 #include "AliRawReaderRoot.h"
67 #include "AliRunLoader.h"
68 #include "AliRunTag.h"
70 ClassImp(AliQAManager)
71 AliQAManager* AliQAManager::fgQAInstance = 0x0;
73 //_____________________________________________________________________________
74 AliQAManager::AliQAManager() :
87 fNumberOfEvents(999999),
91 fRawReaderDelete(kTRUE),
94 fEventSpecie(AliRecoParam::kDefault),
99 fMaxEvents = fNumberOfEvents ;
100 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
101 if (IsSelected(AliQAv1::GetDetName(iDet))) {
102 fLoader[iDet] = NULL ;
103 fQADataMaker[iDet] = NULL ;
104 fQACycles[iDet] = 999999 ;
110 //_____________________________________________________________________________
111 AliQAManager::AliQAManager(AliQAv1::MODE_t mode, const Char_t* gAliceFilename) :
120 fGAliceFileName(gAliceFilename),
123 fMode(AliQAv1::GetModeName(mode)),
124 fNumberOfEvents(999999),
128 fRawReaderDelete(kTRUE),
131 fEventSpecie(AliRecoParam::kDefault),
136 fMaxEvents = fNumberOfEvents ;
137 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
138 if (IsSelected(AliQAv1::GetDetName(iDet))) {
139 fLoader[iDet] = NULL ;
140 fQADataMaker[iDet] = NULL ;
141 fQACycles[iDet] = 999999 ;
147 //_____________________________________________________________________________
148 AliQAManager::AliQAManager(const AliQAManager & qas) :
150 fCurrentEvent(qas.fCurrentEvent),
152 fDetectors(qas.fDetectors),
153 fDetectorsW(qas.fDetectorsW),
157 fGAliceFileName(qas.fGAliceFileName),
158 fFirstEvent(qas.fFirstEvent),
159 fMaxEvents(qas.fMaxEvents),
161 fNumberOfEvents(qas.fNumberOfEvents),
163 fRunNumber(qas.fRunNumber),
165 fRawReaderDelete(kTRUE),
168 fEventSpecie(qas.fEventSpecie),
169 fPrintImage(qas.fPrintImage),
170 fSaveData(qas.fSaveData)
174 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
175 fLoader[iDet] = qas.fLoader[iDet] ;
176 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
177 fQACycles[iDet] = qas.fQACycles[iDet] ;
178 fQAWriteExpert[iDet] = qas.fQAWriteExpert[iDet] ;
182 //_____________________________________________________________________________
183 AliQAManager & AliQAManager::operator = (const AliQAManager & qas)
185 // assignment operator
186 this->~AliQAManager() ;
187 new(this) AliQAManager(qas) ;
191 //_____________________________________________________________________________
192 AliQAManager::~AliQAManager()
195 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
196 if (IsSelected(AliQAv1::GetDetName(iDet))) {
197 fLoader[iDet] = NULL;
198 if (fQADataMaker[iDet]) {
199 (fQADataMaker[iDet])->Finish() ;
200 delete fQADataMaker[iDet] ;
204 if (fRawReaderDelete) {
210 //_____________________________________________________________________________
211 Bool_t AliQAManager::DoIt(const AliQAv1::TASKINDEX_t taskIndex)
213 // Runs all the QA data Maker for every detector
216 // Fill QA data in event loop
217 for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
220 if ( iEvent%10 == 0 )
221 AliDebug(AliQAv1::GetQADebugLevel(), Form("processing event %d", iEvent));
222 if ( taskIndex == AliQAv1::kRAWS ) {
223 if ( !fRawReader->NextEvent() )
225 } else if ( taskIndex == AliQAv1::kESDS ) {
226 if ( fESDTree->GetEntry(iEvent) == 0 )
229 if ( fRunLoader->GetEvent(iEvent) != 0 )
232 // loop over active loaders
235 detList = GetEventInfo()->GetTriggerCluster() ;
236 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
237 if (!detList.IsNull() && !detList.Contains(AliQAv1::GetDetName(iDet)))
239 if (IsSelected(AliQAv1::GetDetName(iDet)) ){
240 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
241 if (!qadm) continue; // This detector doesn't have any QA (for example, HLT)
242 if ( qadm->IsCycleDone() ) {
243 qadm->EndOfCycle(taskIndex) ;
245 TTree * data = NULL ;
246 AliLoader* loader = GetLoader(qadm->GetUniqueID());
248 case AliQAv1::kNULLTASKINDEX :
250 case AliQAv1::kRAWS :
251 qadm->Exec(taskIndex, fRawReader) ;
253 case AliQAv1::kHITS :
256 data = loader->TreeH() ;
258 AliWarning(Form(" Hit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
261 qadm->Exec(taskIndex, data) ;
264 case AliQAv1::kSDIGITS :
267 TString fileName(Form("%s.SDigits.root", AliQAv1::GetDetName(iDet))) ;
268 if (gSystem->FindFile("./", fileName)) {
270 loader->LoadSDigits() ;
271 data = loader->TreeS() ;
273 AliWarning(Form(" SDigit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
276 qadm->Exec(taskIndex, data) ;
281 case AliQAv1::kDIGITS :
283 loader->LoadDigits() ;
284 data = loader->TreeD() ;
286 AliWarning(Form(" Digit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
289 qadm->Exec(taskIndex, data) ;
292 case AliQAv1::kDIGITSR :
294 loader->LoadDigits() ;
295 data = loader->TreeD() ;
297 AliWarning(Form(" Digit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
300 qadm->Exec(taskIndex, data) ;
303 case AliQAv1::kRECPOINTS :
305 loader->LoadRecPoints() ;
306 data = loader->TreeR() ;
308 AliWarning(Form("RecPoints not found for %s", AliQAv1::GetDetName(iDet))) ;
311 qadm->Exec(taskIndex, data) ;
314 case AliQAv1::kTRACKSEGMENTS :
316 case AliQAv1::kRECPARTICLES :
318 case AliQAv1::kESDS :
319 qadm->Exec(taskIndex, fESD) ;
321 case AliQAv1::kNTASKINDEX :
326 Increment(taskIndex) ;
328 // Save QA data for all detectors
332 if ( taskIndex == AliQAv1::kRAWS )
333 fRawReader->RewindEvents() ;
338 //_____________________________________________________________________________
339 Bool_t AliQAManager::Finish(const AliQAv1::TASKINDEX_t taskIndex)
341 // write output to file for all detectors
343 AliQAChecker::Instance()->SetRunNumber(fRunNumber) ;
345 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
346 if (IsSelected(AliQAv1::GetDetName(iDet))) {
347 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
349 qadm->EndOfCycle(taskIndex) ;
355 //_____________________________________________________________________________
356 TObjArray * AliQAManager::GetFromOCDB(AliQAv1::DETECTORINDEX_t det, AliQAv1::TASKINDEX_t task, const Char_t * year) const
358 // Retrieve the list of QA data for a given detector and a given task
359 TObjArray * rv = NULL ;
360 if ( !strlen(AliQAv1::GetQARefStorage()) ) {
361 AliError("No storage defined, use AliQAv1::SetQARefStorage") ;
364 if ( ! IsDefaultStorageSet() ) {
365 TString tmp(AliQAv1::GetQARefDefaultStorage()) ;
368 Instance()->SetDefaultStorage(tmp.Data()) ;
369 Instance()->SetSpecificStorage(Form("%s/*", AliQAv1::GetQAName()), AliQAv1::GetQARefStorage()) ;
371 TString detOCDBDir(Form("%s/%s/%s", AliQAv1::GetQAName(), AliQAv1::GetDetName((Int_t)det), AliQAv1::GetRefOCDBDirName())) ;
372 AliDebug(AliQAv1::GetQADebugLevel(), Form("Retrieving reference data from %s/%s for %s", AliQAv1::GetQARefStorage(), detOCDBDir.Data(), AliQAv1::GetTaskName(task).Data())) ;
373 AliCDBEntry* entry = QAManager()->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
374 TList * listDetQAD = static_cast<TList *>(entry->GetObject()) ;
376 rv = static_cast<TObjArray *>(listDetQAD->FindObject(AliQAv1::GetTaskName(task))) ;
380 //_____________________________________________________________________________
381 TCanvas ** AliQAManager::GetImage(Char_t * detName)
383 // retrieves QA Image for the given detector
384 TCanvas ** rv = NULL ;
385 Int_t detIndex = AliQAv1::GetDetIndex(detName) ;
386 if ( detIndex != AliQAv1::kNULLDET) {
387 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(detIndex) ;
388 rv = qac->GetImage() ;
393 //_____________________________________________________________________________
394 AliLoader * AliQAManager::GetLoader(Int_t iDet)
396 // get the loader for a detector
398 if ( !fRunLoader || iDet == AliQAv1::kCORR || iDet == AliQAv1::kGLOBAL )
401 TString detName = AliQAv1::GetDetName(iDet) ;
402 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
404 return fLoader[iDet] ;
406 // load the QA data maker object
407 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
408 TString loaderName = "Ali" + detName + "Loader" ;
410 AliLoader * loader = NULL ;
411 // first check if a plugin is defined for the quality assurance data maker
412 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
413 // if not, add a plugin for it
414 if (!pluginHandler) {
415 AliDebug(AliQAv1::GetQADebugLevel(), Form("defining plugin for %s", loaderName.Data())) ;
416 TString libs = gSystem->GetLibraries() ;
417 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
418 pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
420 pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
422 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
424 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
425 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
428 fLoader[iDet] = loader ;
432 //_____________________________________________________________________________
433 AliQAv1 * AliQAManager::GetQA(UInt_t run, UInt_t evt)
435 // retrieves the QA object stored in a file named "Run{run}.Event{evt}_1.ESD.tag.root"
436 Char_t * fileName = Form("Run%d.Event%d_1.ESD.tag.root", run, evt) ;
437 TFile * tagFile = TFile::Open(fileName) ;
439 AliError(Form("File %s not found", fileName)) ;
442 TTree * tagTree = static_cast<TTree *>(tagFile->Get("T")) ;
444 AliError(Form("Tree T not found in %s", fileName)) ;
448 AliRunTag * tag = new AliRunTag ;
449 tagTree->SetBranchAddress("AliTAG", &tag) ;
450 tagTree->GetEntry(evt) ;
451 AliQAv1 * qa = AliQAv1::Instance(tag->GetQALength(), tag->GetQAArray(), tag->GetESLength(), tag->GetEventSpecies()) ;
456 //_____________________________________________________________________________
457 AliQADataMaker * AliQAManager::GetQADataMaker(const Int_t iDet)
459 // get the quality assurance data maker for a detector
461 AliQADataMaker * qadm = fQADataMaker[iDet] ;
465 qadm->SetEventSpecie(fEventSpecie) ;
466 if ( qadm->GetRecoParam() )
467 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
468 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
470 } else if (iDet == AliQAv1::kGLOBAL && strcmp(GetMode(), AliQAv1::GetModeName(AliQAv1::kRECMODE)) == 0) { //Global QA
472 qadm = new AliGlobalQADataMaker();
473 qadm->SetName(AliQAv1::GetDetName(iDet));
474 qadm->SetUniqueID(iDet);
475 fQADataMaker[iDet] = qadm;
476 qadm->SetEventSpecie(fEventSpecie) ;
477 if ( qadm->GetRecoParam() )
478 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
479 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
481 } else if (iDet == AliQAv1::kCORR && strcmp(GetMode(), AliQAv1::GetModeName(AliQAv1::kRECMODE)) == 0 ) { //the data maker for correlations among detectors
482 qadm = new AliCorrQADataMakerRec(fQADataMaker) ;
483 qadm->SetName(AliQAv1::GetDetName(iDet));
484 qadm->SetUniqueID(iDet);
485 fQADataMaker[iDet] = qadm;
486 qadm->SetEventSpecie(fEventSpecie) ;
487 if ( qadm->GetRecoParam() )
488 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
489 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
491 } else if ( iDet < AliQAv1::kGLOBAL ) {
492 TString smode(GetMode()) ;
493 if (smode.Contains(AliQAv1::GetModeName(AliQAv1::kQAMODE)))
494 smode = AliQAv1::GetModeName(AliQAv1::kRECMODE) ;
495 // load the QA data maker object
496 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
497 TString detName = AliQAv1::GetDetName(iDet) ;
498 TString qadmName = "Ali" + detName + "QADataMaker" + smode ;
500 // first check if a plugin is defined for the quality assurance data maker
501 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
502 // if not, add a plugin for it
503 if (!pluginHandler) {
504 AliDebug(AliQAv1::GetQADebugLevel(), Form("defining plugin for %s", qadmName.Data())) ;
505 TString libs = gSystem->GetLibraries() ;
506 TString temp(smode) ;
508 if (libs.Contains("lib" + detName + smode + ".so") || (gSystem->Load("lib" + detName + temp.Data() + ".so") >= 0)) {
509 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
511 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
513 pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
515 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
516 qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
519 qadm->SetName(AliQAv1::GetDetName(iDet));
520 qadm->SetUniqueID(iDet);
521 fQADataMaker[iDet] = qadm ;
522 qadm->SetEventSpecie(fEventSpecie) ;
523 if ( qadm->GetRecoParam() )
524 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
525 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
531 //_____________________________________________________________________________
532 void AliQAManager::EndOfCycle(TObjArray * detArray)
534 // End of cycle QADataMakers
536 AliQAChecker::Instance()->SetRunNumber(fRunNumber) ;
539 fakeCanvas.Print(Form("%s%s%d.%s[", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps") ;
540 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
541 if (IsSelected(AliQAv1::GetDetName(iDet))) {
542 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
545 // skip non active detectors
547 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQAv1::GetDetName(iDet))) ;
548 if (!det || !det->IsActive())
551 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(iDet) ;
553 qac->SetPrintImage(fPrintImage) ;
554 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
555 if ( fTasks.Contains(Form("%d", taskIndex)) )
556 qadm->EndOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex))) ;
562 fakeCanvas.Print(Form("%s%s%d.%s]", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps");
565 //_____________________________________________________________________________
566 void AliQAManager::EndOfCycle(TString detectors)
568 // End of cycle QADataMakers
570 AliQAChecker::Instance()->SetRunNumber(fRunNumber) ;
573 fakeCanvas.Print(Form("%s%s%d.%s[", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps") ;
574 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
575 if (IsSelected(AliQAv1::GetDetName(iDet))) {
576 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
579 // skip non active detectors
580 if (!detectors.Contains(AliQAv1::GetDetName(iDet)))
582 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(iDet) ;
584 qac->SetPrintImage(fPrintImage) ;
585 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
586 if ( fTasks.Contains(Form("%d", taskIndex)) )
587 qadm->EndOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex))) ;
593 fakeCanvas.Print(Form("%s%s%d.%s]", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps");
596 //_____________________________________________________________________________
597 AliRecoParam::EventSpecie_t AliQAManager::GetEventSpecieFromESD()
599 AliRecoParam::EventSpecie_t runtype = AliRecoParam::kDefault ;
600 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
601 TFile * esdFile = TFile::Open("AliESDs.root") ;
602 TTree * esdTree = static_cast<TTree *> (esdFile->Get("esdTree")) ;
604 AliError("esdTree not found") ;
606 AliESDEvent * esd = new AliESDEvent() ;
607 esd->ReadFromTree(esdTree) ;
608 esdTree->GetEntry(0) ;
609 runtype = AliRecoParam::Convert(esd->GetEventType()) ;
612 AliError("AliESDs.root not found") ;
617 //_____________________________________________________________________________
618 void AliQAManager::Increment(const AliQAv1::TASKINDEX_t taskIndex)
620 // Increments the cycle counter for all QA Data Makers
621 static AliQAv1::TASKINDEX_t currentTask = AliQAv1::kNTASKINDEX ;
622 if ( (currentTask == taskIndex) && taskIndex != AliQAv1::kNULLTASKINDEX )
625 currentTask = taskIndex ;
626 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
627 if (IsSelected(AliQAv1::GetDetName(iDet))) {
628 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
635 //_____________________________________________________________________________
636 Bool_t AliQAManager::InitQA(const AliQAv1::TASKINDEX_t taskIndex, const Char_t * input )
638 // Initialize the event source and QA data makers
640 fTasks += Form("%d", taskIndex) ;
642 if (taskIndex == AliQAv1::kRAWS) {
644 fRawReader = AliRawReader::Create(input);
648 fRawReaderDelete = kTRUE ;
649 fRawReader->NextEvent() ;
650 fRunNumber = fRawReader->GetRunNumber() ;
652 fRawReader->RewindEvents();
653 fNumberOfEvents = 999999 ;
654 if ( fMaxEvents < 0 )
655 fMaxEvents = fNumberOfEvents ;
656 } else if (taskIndex == AliQAv1::kESDS) {
657 fTasks = AliQAv1::GetTaskName(AliQAv1::kESDS) ;
658 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
659 TFile * esdFile = TFile::Open("AliESDs.root") ;
660 fESDTree = static_cast<TTree *> (esdFile->Get("esdTree")) ;
662 AliError("esdTree not found") ;
665 fESD = new AliESDEvent() ;
666 fESD->ReadFromTree(fESDTree) ;
667 fESDTree->GetEntry(0) ;
668 fRunNumber = fESD->GetRunNumber() ;
669 fNumberOfEvents = fESDTree->GetEntries() ;
670 if ( fMaxEvents < 0 )
671 fMaxEvents = fNumberOfEvents ;
674 AliError("AliESDs.root not found") ;
678 if ( !InitRunLoader() ) {
679 AliWarning("No Run Loader not found") ;
681 fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
682 if ( fMaxEvents < 0 )
683 fMaxEvents = fNumberOfEvents ;
688 TObjArray* detArray = NULL ;
689 if (fRunLoader) // check if RunLoader exists
690 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
691 detArray = fRunLoader->GetAliRun()->Detectors() ;
692 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
695 // Initialize all QA data makers for all detectors
696 fRunNumber = AliCDBManager::Instance()->GetRun() ;
697 if ( ! AliGeomManager::GetGeometry() )
698 AliGeomManager::LoadGeometry() ;
700 InitQADataMaker(fRunNumber, detArray) ; //, fCycleSame, kTRUE, detArray) ;
705 while (timer.CpuTime()<5) {
707 gSystem->ProcessEvents();
709 fakeCanvas.Print(Form("%s%s%d.%s[", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps") ;
714 //_____________________________________________________________________________
715 void AliQAManager::InitQADataMaker(UInt_t run, TObjArray * detArray)
717 // Initializes The QADataMaker for all active detectors and for all active tasks
719 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
720 if (IsSelected(AliQAv1::GetDetName(iDet))) {
721 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
723 AliError(Form("AliQADataMaker not found for %s", AliQAv1::GetDetName(iDet))) ;
724 fDetectorsW.ReplaceAll(AliQAv1::GetDetName(iDet), "") ;
726 if (fQAWriteExpert[iDet])
727 qadm->SetWriteExpert() ;
728 AliDebug(AliQAv1::GetQADebugLevel(), Form("Data Maker found for %s %d", qadm->GetName(), qadm->WriteExpert())) ;
729 // skip non active detectors
731 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQAv1::GetDetName(iDet))) ;
732 if (!det || !det->IsActive())
735 // Set default reco params
736 Bool_t sameCycle = kFALSE ;
737 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
738 if ( fTasks.Contains(Form("%d", taskIndex)) ) {
739 qadm->Init(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
740 qadm->StartOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex)), run, sameCycle) ;
750 //_____________________________________________________________________________
751 Bool_t AliQAManager::InitRunLoader()
753 // get or create the run loader
757 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
758 // load all base libraries to get the loader classes
759 TString libs = gSystem->GetLibraries() ;
760 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
761 if (!IsSelected(AliQAv1::GetDetName(iDet)))
763 TString detName = AliQAv1::GetDetName(iDet) ;
764 if (detName == "HLT")
766 if (libs.Contains("lib" + detName + "base.so"))
768 gSystem->Load("lib" + detName + "base.so");
770 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
772 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
775 fRunLoader->CdGAFile();
776 if (fRunLoader->LoadgAlice() == 0) {
777 gAlice = fRunLoader->GetAliRun();
781 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
785 } else { // galice.root does not exist
786 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
792 fRunLoader->LoadHeader();
793 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
798 //_____________________________________________________________________________
799 Bool_t AliQAManager::IsSelected(const Char_t * det)
801 // check whether detName is contained in detectors
802 // if yes, it is removed from detectors
805 const TString detName(det) ;
806 // always activates Correlation
807 // if ( detName.Contains(AliQAv1::GetDetName(AliQAv1::kCORR)) || detName.Contains(AliQAv1::GetDetName(AliQAv1::kGLOBAL))) {
810 // check if all detectors are selected
811 if (fDetectors.Contains("ALL")) {
814 } else if ((fDetectors.CompareTo(detName) == 0) ||
815 fDetectors.BeginsWith(detName+" ") ||
816 fDetectors.EndsWith(" "+detName) ||
817 fDetectors.Contains(" "+detName+" ")) {
824 //_____________________________________________________________________________
825 Bool_t AliQAManager::Merge(Int_t runNumber, const char *fileName) const
827 // Merge data from all detectors from a given run in one single file
828 // Merge the QA results from all the data chunks in one run
829 // The 'fileName' is name of the output file with merged QA data
830 if ( runNumber == -1)
831 runNumber = fRunNumber ;
832 Bool_t rv = MergeData(runNumber,fileName) ;
833 //rv *= MergeResults(runNumber) ; // not needed for the time being
837 //______________________________________________________________________
838 Bool_t AliQAManager::MergeXML(const Char_t * collectionFile, const Char_t * subFile, const Char_t * outFile)
840 // merges files listed in a xml collection
841 // usage Merge(collection, outputFile))
842 // collection: is a xml collection
846 if ( strstr(collectionFile, ".xml") == 0 ) {
847 AliError("Input collection file must be an \".xml\" file\n") ;
852 TGrid::Connect("alien://");
856 // Open the file collection
857 AliInfoClass(Form("*** Create Collection ***\n*** Wk-Dir = |%s| \n*** Coll = |%s| \n",gSystem->WorkingDirectory(), collectionFile));
859 TGridCollection * collection = (TGridCollection*)gROOT->ProcessLine(Form("TAlienCollection::Open(\"%s\")",collectionFile));
860 TGridResult* result = collection->GetGridResult("", 0, 0);
863 const Char_t * turl ;
864 TFileMerger merger(kFALSE) ;
866 TString tempo(collectionFile) ;
868 tempo.ReplaceAll(".xml", subFile) ;
870 tempo.ReplaceAll(".xml", "_Merged.root") ;
871 outFile = tempo.Data() ;
873 merger.OutputFile(outFile) ;
875 while ( (turl = result->GetKey(index, "turl")) ) {
878 file = Form("%s#%s", turl, subFile) ;
880 file = Form("%s", turl) ;
882 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s\n", file)) ;
883 merger.AddFile(file) ;
890 AliDebug(AliQAv1::GetQADebugLevel(), Form("Files merged into %s\n", outFile)) ;
896 //_____________________________________________________________________________
897 void AliQAManager::MergeCustom() const
899 // Custom Merge of QA data from all detectors for all runs in one single file
900 // search all the run numbers
901 // search all the run numbers
902 gROOT->ProcessLine(".! ls *QA*.root > QAtempo.txt") ;
904 FILE * theQAfiles = fopen("QAtempo.txt", "r") ;
907 TIter nextRun(&srunList) ;
908 TObjString * srun = NULL ;
909 Int_t loRun = 999999999 ;
911 while ( theQAfile.Gets(theQAfiles) ) {
912 Bool_t runExist = kFALSE ;
913 TString srunNew(theQAfile(theQAfile.Index("QA.")+3, theQAfile.Index(".root")-(theQAfile.Index("QA.")+3))) ;
914 Int_t cuRun = srunNew.Atoi() ;
919 while ( (srun = static_cast<TObjString *> (nextRun())) ) {
920 if ( cuRun == (srun->String()).Atoi() ) {
927 srunList.Add(new TObjString(srunNew.Data()));
930 Int_t runNumber = 0 ;
931 TFile mergedFile(Form("Merged.%s.Data.root", AliQAv1::GetQADataFileName()), "RECREATE") ;
932 TH1I * hisRun = new TH1I("hLMR", "List of merged runs", hiRun-loRun+10, loRun, hiRun+10) ;
933 // create the structure into the merged file
934 for (Int_t iDet = 0; iDet < AliQAv1::kNDET ; iDet++) {
935 TDirectory * detDir = mergedFile.mkdir(AliQAv1::GetDetName(iDet)) ;
936 for (Int_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
938 TDirectory * taskDir = gDirectory->mkdir(AliQAv1::GetTaskName(taskIndex)) ;
939 for (Int_t es = 0 ; es < AliRecoParam::kNSpecies ; es++) {
941 TDirectory * esDir = gDirectory->mkdir(AliRecoParam::GetEventSpecieName(es)) ;
943 gDirectory->mkdir(AliQAv1::GetExpert()) ;
947 while ( (srun = static_cast<TObjString *> (nextRun())) ) {
948 runNumber = (srun->String()).Atoi() ;
949 hisRun->Fill(runNumber) ;
950 AliDebug(AliQAv1::GetQADebugLevel(), Form("Merging run number %d", runNumber)) ;
951 // search all QA files for runNumber in the current directory
952 Char_t * fileList[AliQAv1::kNDET] ;
954 for (Int_t iDet = 0; iDet < AliQAv1::kNDET ; iDet++) {
955 Char_t * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQAv1::GetDetName(iDet), AliQAv1::GetQADataFileName(), runNumber));
957 fileList[index++] = file ;
960 AliError("No QA data file found\n") ;
963 for ( Int_t i = 0 ; i < index ; i++) {
964 TFile * inFile = TFile::Open(fileList[i]) ;
965 TList * listOfKeys =inFile->GetListOfKeys() ;
966 TIter nextkey(listOfKeys) ;
968 TString dirName("") ;
969 while ( (obj1 = nextkey()) ) {
970 TDirectory * directoryDet = inFile->GetDirectory(obj1->GetName()) ;
971 if ( directoryDet ) {
972 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s dir = %s", inFile->GetName(), directoryDet->GetName())) ;
973 dirName += Form("%s/", directoryDet->GetName() ) ;
975 TList * listOfTasks = directoryDet->GetListOfKeys() ;
976 TIter nextTask(listOfTasks) ;
978 while ( (obj2 = nextTask()) ) {
979 TDirectory * directoryTask = directoryDet->GetDirectory(obj2->GetName()) ;
980 if ( directoryTask ) {
981 dirName += Form("%s", obj2->GetName()) ;
982 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s", dirName.Data())) ;
983 directoryTask->cd() ;
984 TList * listOfEventSpecie = directoryTask->GetListOfKeys() ;
985 TIter nextEventSpecie(listOfEventSpecie) ;
987 while ( (obj3 = nextEventSpecie()) ) {
988 TDirectory * directoryEventSpecie = directoryTask->GetDirectory(obj3->GetName()) ;
989 if ( directoryEventSpecie ) {
990 dirName += Form("/%s/", obj3->GetName()) ;
991 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s\n", dirName.Data())) ;
992 directoryEventSpecie->cd() ;
993 // histograms are here
994 TDirectory * mergedDirectory = mergedFile.GetDirectory(dirName.Data()) ;
995 TList * listOfData = directoryEventSpecie->GetListOfKeys() ;
996 TIter nextData(listOfData) ;
998 while ( (key = static_cast<TKey *>(nextData())) ) {
999 TString className(key->GetClassName()) ;
1000 if ( className.Contains("TH") || className.Contains("TProfile") ) {
1001 TH1 * histIn = static_cast<TH1*> (key->ReadObj()) ;
1002 TH1 * histOu = static_cast<TH1*> (mergedDirectory->FindObjectAny(histIn->GetName())) ;
1003 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s %p %p\n", key->GetName(), histIn, histOu)) ;
1004 mergedDirectory->cd() ;
1008 histOu->Add(histIn) ;
1009 histOu->Write(histOu->GetName(), kOverwrite) ;
1012 else if ( className.Contains("TDirectoryFile") ) {
1013 TDirectory * dirExpert = directoryEventSpecie->GetDirectory(key->GetName()) ;
1015 TDirectory * mergedDirectoryExpert = mergedDirectory->GetDirectory(dirExpert->GetName()) ;
1016 TList * listOfExpertData = dirExpert->GetListOfKeys() ;
1017 TIter nextExpertData(listOfExpertData) ;
1019 while ( (keykey = static_cast<TKey *>(nextExpertData())) ) {
1020 TString classNameExpert(keykey->GetClassName()) ;
1021 if (classNameExpert.Contains("TH")) {
1022 TH1 * histInExpert = static_cast<TH1*> (keykey->ReadObj()) ;
1023 TH1 * histOuExpert = static_cast<TH1*> (mergedDirectory->FindObjectAny(histInExpert->GetName())) ;
1024 mergedDirectoryExpert->cd() ;
1025 if ( ! histOuExpert ) {
1026 histInExpert->Write() ;
1028 histOuExpert->Add(histInExpert) ;
1029 histOuExpert->Write(histOuExpert->GetName(), kOverwrite) ;
1034 AliError(Form("No merge done for this object %s in %s", key->GetName(), dirName.Data())) ;
1037 dirName.ReplaceAll(Form("/%s/",obj3->GetName()), "") ;
1040 dirName.ReplaceAll(obj2->GetName(), "") ;
1050 mergedFile.Close() ;
1054 //_____________________________________________________________________________
1055 Bool_t AliQAManager::MergeData(const Int_t runNumber, const char *fileName) const
1057 // Merge QA data from all detectors for a given run in one single file
1059 TFileMerger merger(kFALSE) ;
1060 TString outFileName = fileName;
1061 if (outFileName.IsNull()) outFileName.Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName());
1062 merger.OutputFile(outFileName.Data()) ;
1063 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
1064 Char_t * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQAv1::GetDetName(iDet), AliQAv1::GetQADataFileName(), runNumber));
1066 merger.AddFile(file);
1073 //_____________________________________________________________________________
1074 Bool_t AliQAManager::MergeResults(const Int_t runNumber) const
1076 // Merge the QA result from all the data chunks in a run
1077 // to be revised whwn it will be used (see MergeData)
1079 cmd = Form(".! ls %s*.root > tempo.txt", AliQAv1::GetQADataFileName()) ;
1080 gROOT->ProcessLine(cmd.Data()) ;
1081 ifstream in("tempo.txt") ;
1082 const Int_t chunkMax = 100 ;
1083 TString fileList[chunkMax] ;
1088 in >> fileList[index] ;
1091 AliDebug(AliQAv1::GetQADebugLevel(), Form("index = %d file = %s", index, (fileList[index].Data()))) ;
1096 AliError("No QA Result File found") ;
1100 TFileMerger merger ;
1101 TString outFileName ;
1102 if (runNumber != -1)
1103 outFileName = Form("Merged.%s.Result.%d.root",AliQAv1::GetQADataFileName(),runNumber);
1105 outFileName = Form("Merged.%s.Result.root",AliQAv1::GetQADataFileName());
1106 merger.OutputFile(outFileName.Data()) ;
1107 for (Int_t ifile = 0 ; ifile < index ; ifile++) {
1108 TString file = fileList[ifile] ;
1109 merger.AddFile(file) ;
1116 //_____________________________________________________________________________
1117 void AliQAManager::Reset(const Bool_t sameCycle)
1119 // Reset the default data members
1121 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
1122 if (IsSelected(AliQAv1::GetDetName(iDet))) {
1123 AliQADataMaker * qadm = GetQADataMaker(iDet);
1128 if (fRawReaderDelete) {
1133 fCycleSame = sameCycle ;
1137 fNumberOfEvents = 999999 ;
1140 //_____________________________________________________________________________
1141 void AliQAManager::ResetDetectors(AliQAv1::TASKINDEX_t task, AliQAv1::DETECTORINDEX_t det)
1143 //calls ResetDetector of specified or all detectors
1145 UInt_t iDetMax = fgkNDetectors ;
1146 if ( det != AliQAv1::kNULLDET ) {
1151 for (iDet = 0; iDet < iDetMax ; iDet++) {
1152 if (IsSelected(AliQAv1::GetDetName(iDet))) {
1153 AliQADataMaker * qadm = GetQADataMaker(iDet);
1154 qadm->ResetDetector(task);
1159 //_____________________________________________________________________________
1160 AliQAManager * AliQAManager::QAManager(AliQAv1::MODE_t mode, TMap *entryCache, Int_t run)
1162 // returns AliQAManager instance (singleton)
1164 if (!fgQAInstance) {
1165 if ( (mode != AliQAv1::kSIMMODE) && (mode != AliQAv1::kRECMODE) && (mode != AliQAv1::kQAMODE) ) {
1166 AliWarningClass("You must specify kSIMMODE or kRECMODE or kQAMODE") ;
1169 fgQAInstance = new AliQAManager(mode) ;
1171 fgQAInstance->Init();
1173 fgQAInstance->InitFromCache(entryCache,run);
1175 return fgQAInstance;
1178 //_____________________________________________________________________________
1179 AliQAManager * AliQAManager::QAManager(AliQAv1::TASKINDEX_t task)
1181 // returns AliQAManager instance (singleton)
1182 return QAManager(AliQAv1::Mode(task)) ;
1185 //_____________________________________________________________________________
1186 TString AliQAManager::Run(const Char_t * detectors, AliRawReader * rawReader, const Bool_t sameCycle)
1188 //Runs all the QA data Maker for Raws only
1190 fCycleSame = sameCycle ;
1191 fRawReader = rawReader ;
1192 fDetectors = detectors ;
1193 fDetectorsW = detectors ;
1195 AliCDBManager* man = AliCDBManager::Instance() ;
1197 if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
1198 rawReader->NextEvent() ;
1199 man->SetRun(fRawReader->GetRunNumber()) ;
1200 rawReader->RewindEvents() ;
1204 if ( !InitQA(AliQAv1::kRAWS) )
1206 fRawReaderDelete = kFALSE ;
1208 DoIt(AliQAv1::kRAWS) ;
1209 return fDetectorsW ;
1212 //_____________________________________________________________________________
1213 TString AliQAManager::Run(const Char_t * detectors, const Char_t * fileName, const Bool_t sameCycle)
1215 //Runs all the QA data Maker for Raws only
1217 fCycleSame = sameCycle ;
1218 fDetectors = detectors ;
1219 fDetectorsW = detectors ;
1221 AliCDBManager* man = AliCDBManager::Instance() ;
1222 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
1223 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
1225 AliFatal("galice.root file not found in current directory") ;
1229 if ( ! rl->GetAliRun() ) {
1230 AliFatal("AliRun not found in galice.root") ;
1233 man->SetRun(rl->GetHeader()->GetRun());
1239 if ( !InitQA(AliQAv1::kRAWS, fileName) )
1242 DoIt(AliQAv1::kRAWS) ;
1243 return fDetectorsW ;
1246 //_____________________________________________________________________________
1247 TString AliQAManager::Run(const Char_t * detectors, const AliQAv1::TASKINDEX_t taskIndex, Bool_t const sameCycle, const Char_t * fileName )
1249 // Runs all the QA data Maker for every detector
1251 fCycleSame = sameCycle ;
1252 fDetectors = detectors ;
1253 fDetectorsW = detectors ;
1255 AliCDBManager* man = AliCDBManager::Instance() ;
1256 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
1257 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
1259 AliFatal("galice.root file not found in current directory") ;
1263 if ( ! rl->GetAliRun() ) {
1264 AliDebug(AliQAv1::GetQADebugLevel(), "AliRun not found in galice.root") ;
1267 man->SetRun(rl->GetHeader()->GetRun()) ;
1271 if ( taskIndex == AliQAv1::kNULLTASKINDEX) {
1272 for (UInt_t task = 0; task < AliQAv1::kNTASKINDEX; task++) {
1273 if ( fTasks.Contains(Form("%d", task)) ) {
1275 if ( !InitQA(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(task)), fileName) )
1277 DoIt(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(task))) ;
1282 if ( !InitQA(taskIndex, fileName) )
1286 return fDetectorsW ;
1289 //_____________________________________________________________________________
1290 void AliQAManager::RunOneEvent(AliRawReader * rawReader)
1292 //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
1296 if (fTasks.Contains(Form("%d", AliQAv1::kRAWS))){
1298 if ( GetEventInfo())
1299 detList = GetEventInfo()->GetTriggerCluster() ;
1300 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1301 if (!IsSelected(AliQAv1::GetDetName(iDet)) || (!detList.IsNull() && !detList.Contains(AliQAv1::GetDetName(iDet))))
1303 AliQADataMaker *qadm = GetQADataMaker(iDet);
1306 if ( qadm->IsCycleDone() ) {
1307 qadm->EndOfCycle() ;
1309 qadm->SetEventSpecie(fEventSpecie) ;
1310 if ( qadm->GetRecoParam() )
1311 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
1312 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
1313 qadm->Exec(AliQAv1::kRAWS, rawReader) ;
1318 //_____________________________________________________________________________
1319 void AliQAManager::RunOneEvent(AliESDEvent *& esd, AliESDEvent *& hltesd)
1321 //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1323 if (fTasks.Contains(Form("%d", AliQAv1::kESDS))) {
1325 if ( GetEventInfo())
1326 detList = GetEventInfo()->GetTriggerCluster() ;
1327 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1328 if (!IsSelected(AliQAv1::GetDetName(iDet)) || (!detList.IsNull() && !detList.Contains(AliQAv1::GetDetName(iDet))))
1330 AliQADataMaker *qadm = GetQADataMaker(iDet);
1333 qadm->SetEventSpecie(fEventSpecie) ;
1334 if ( qadm->GetRecoParam() )
1335 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
1336 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
1337 if ( qadm->IsCycleDone() ) {
1338 qadm->EndOfCycle() ;
1340 if (iDet == AliQAv1::kHLT) {
1343 esdarray.Add(hltesd);
1344 qadm->Exec(AliQAv1::kESDS, &esdarray);
1346 qadm->Exec(AliQAv1::kESDS, esd) ;
1352 //_____________________________________________________________________________
1353 void AliQAManager::RunOneEventInOneDetector(Int_t det, TTree * tree)
1355 // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1358 if ( GetEventInfo())
1359 detList = GetEventInfo()->GetTriggerCluster() ;
1360 if (!detList.IsNull() && !detList.Contains(AliQAv1::GetDetName(det)))
1363 TString test(tree->GetName()) ;
1364 if (fTasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
1365 if (IsSelected(AliQAv1::GetDetName(det))) {
1366 AliQADataMaker *qadm = GetQADataMaker(det);
1368 qadm->SetEventSpecie(fEventSpecie) ;
1369 if ( qadm->GetRecoParam() ) {
1370 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
1371 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
1373 AliError(Form("%d defined by %s is not an event specie", qadm->GetRecoParam()->GetEventSpecie(), qadm->GetName())) ;
1375 if ( qadm->IsCycleDone() ) {
1376 qadm->EndOfCycle() ;
1378 if (test.Contains("TreeD")) {
1379 qadm->Exec(AliQAv1::kDIGITSR, tree) ;
1380 } else if (test.Contains("TreeR")) {
1381 qadm->Exec(AliQAv1::kRECPOINTS, tree) ;
1388 //_____________________________________________________________________________
1389 Bool_t AliQAManager::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const Char_t * year, const Char_t * detectors) const
1391 // take the locasl QA data merge into a single file and save in OCDB
1393 TString tmp(AliQAv1::GetQARefStorage()) ;
1394 if ( tmp.IsNull() ) {
1395 AliError("No storage defined, use AliQAv1::SetQARefStorage") ;
1398 if ( !(tmp.Contains(AliQAv1::GetLabLocalOCDB()) || tmp.Contains(AliQAv1::GetLabAliEnOCDB())) ) {
1399 AliError(Form("%s is a wrong storage, use %s or %s", AliQAv1::GetQARefStorage(), AliQAv1::GetLabLocalOCDB().Data(), AliQAv1::GetLabAliEnOCDB().Data())) ;
1402 TString sdet(detectors) ;
1405 if ( sdet.Contains("ALL") ) {
1406 rv = Merge(runNumber) ;
1409 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQAv1::GetQADataFileName(), runNumber)) ;
1410 inputFile = TFile::Open(inputFileName.Data()) ;
1411 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ;
1413 for (Int_t index = 0; index < AliQAv1::kNDET; index++) {
1414 if (sdet.Contains(AliQAv1::GetDetName(index))) {
1415 TString inputFileName(Form("%s.%s.%d.root", AliQAv1::GetDetName(index), AliQAv1::GetQADataFileName(), runNumber)) ;
1416 inputFile = TFile::Open(inputFileName.Data()) ;
1417 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ;
1424 //_____________________________________________________________________________
1425 Bool_t AliQAManager::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const Char_t * year, AliRecoParam::EventSpecie_t es) const
1427 // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1429 AliDebug(AliQAv1::GetQADebugLevel(), Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQAv1::GetQARefStorage())) ;
1430 if ( ! IsDefaultStorageSet() ) {
1431 TString tmp( AliQAv1::GetQARefStorage() ) ;
1432 if ( tmp.Contains(AliQAv1::GetLabLocalOCDB()) )
1433 Instance()->SetDefaultStorage(AliQAv1::GetQARefStorage()) ;
1435 TString tmp1(AliQAv1::GetQARefDefaultStorage()) ;
1437 tmp1.Append("?user=alidaq") ;
1438 Instance()->SetDefaultStorage(tmp1.Data()) ;
1441 Instance()->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
1443 Instance()->SetRun(runNumber);
1445 AliCDBMetaData mdr ;
1446 mdr.SetResponsible("yves schutz");
1448 for ( Int_t detIndex = 0 ; detIndex < AliQAv1::kNDET ; detIndex++) {
1449 TDirectory * detDir = inputFile->GetDirectory(AliQAv1::GetDetName(detIndex)) ;
1451 AliDebug(AliQAv1::GetQADebugLevel(), Form("Entering %s", detDir->GetName())) ;
1452 AliQAv1::SetQARefDataDirName(es) ;
1453 TString detOCDBDir(Form("%s/%s/%s", AliQAv1::GetDetName(detIndex), AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName())) ;
1454 AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity()) ;
1455 TList * listDetQAD = new TList() ;
1456 TString listName(Form("%s QA data Reference", AliQAv1::GetDetName(detIndex))) ;
1457 mdr.SetComment(Form("%s QA stuff", AliQAv1::GetDetName(detIndex)));
1458 listDetQAD->SetName(listName) ;
1459 TList * taskList = detDir->GetListOfKeys() ;
1460 TIter nextTask(taskList) ;
1462 while ( (taskKey = static_cast<TKey*>(nextTask())) ) {
1463 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ;
1464 TDirectory * esDir = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ;
1465 AliDebug(AliQAv1::GetQADebugLevel(), Form("Saving %s", esDir->GetName())) ;
1466 TObjArray * listTaskQAD = new TObjArray(100) ;
1467 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1468 listDetQAD->Add(listTaskQAD) ;
1469 TList * histList = esDir->GetListOfKeys() ;
1470 TIter nextHist(histList) ;
1472 while ( (histKey = static_cast<TKey*>(nextHist())) ) {
1473 TObject * odata = esDir->Get(histKey->GetName()) ;
1475 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1477 if ( AliQAv1::GetExpert() == histKey->GetName() ) {
1478 TDirectory * expertDir = esDir->GetDirectory(histKey->GetName()) ;
1479 TList * expertHistList = expertDir->GetListOfKeys() ;
1480 TIter nextExpertHist(expertHistList) ;
1481 TKey * expertHistKey ;
1482 while ( (expertHistKey = static_cast<TKey*>(nextExpertHist())) ) {
1483 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ;
1484 if ( !expertOdata ) {
1485 AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1487 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", expertHistKey->GetName())) ;
1488 if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
1489 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", expertHistKey->GetName())) ;
1490 TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ;
1491 listTaskQAD->Add(hExpertdata) ;
1496 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", histKey->GetName())) ;
1497 if ( odata->IsA()->InheritsFrom("TH1") ) {
1498 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", histKey->GetName())) ;
1499 TH1 * hdata = static_cast<TH1*>(odata) ;
1500 listTaskQAD->Add(hdata) ;
1505 Instance()->Put(listDetQAD, idr, &mdr) ;
1511 //_____________________________________________________________________________
1513 void AliQAManager::SetCheckerExternParam(AliQAv1::DETECTORINDEX_t detIndex, TList * parameterList)
1515 // set the external parameters list for the detector checkers
1516 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(detIndex) ;
1517 qac->SetExternParamlist(parameterList) ;
1518 qac->PrintExternParam() ;
1521 //_____________________________________________________________________________
1522 void AliQAManager::SetEventSpecie(AliRecoParam::EventSpecie_t es)
1524 // set the current event specie and inform AliQAv1 that this event specie has been encountered
1526 AliQAv1::Instance()->SetEventSpecie(es) ;
1529 //_____________________________________________________________________________
1530 void AliQAManager::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par)
1532 // Set custom reconstruction parameters for a given detector
1533 // Single set of parameters for all the events
1534 GetQADataMaker(det)->SetRecoParam(par) ;
1537 //_____________________________________________________________________________
1538 void AliQAManager::SetWriteExpert()
1540 // enable the writing of QA expert data
1541 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1542 if (IsSelected(AliQAv1::GetDetName(iDet)))
1543 fQAWriteExpert[iDet] = kTRUE ;
1547 //_____________________________________________________________________________
1548 void AliQAManager::Destroy() {
1549 // delete AliQAManager instance and
1550 // all associated objects
1553 delete fgQAInstance ;
1554 fgQAInstance = NULL ;
1558 //_____________________________________________________________________________
1559 void AliQAManager::ShowQA() {
1560 // Show the result of the QA checking
1561 // for all detectors
1562 for ( Int_t detIndex = 0 ; detIndex < AliQAv1::kNDET ; detIndex++)
1563 if ( IsSelected(AliQAv1::GetDetName(detIndex)) )
1564 AliQAv1::Instance(AliQAv1::GetDetIndex(AliQAv1::GetDetName(detIndex)))->Show() ;