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(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); //
29 ///////////////////////////////////////////////////////////////////////////////
31 #include <TAlienCollection.h>
34 #include <TFileMerger.h>
36 #include <TGridCollection.h>
37 #include <TGridResult.h>
38 #include <TPluginManager.h>
43 #include "AliCDBManager.h"
44 #include "AliCDBEntry.h"
46 #include "AliCDBMetaData.h"
47 #include "AliCodeTimer.h"
48 #include "AliCorrQADataMakerRec.h"
49 #include "AliDetectorRecoParam.h"
50 #include "AliESDEvent.h"
51 #include "AliGeomManager.h"
52 #include "AliGlobalQADataMaker.h"
53 #include "AliHeader.h"
55 #include "AliModule.h"
57 #include "AliQADataMakerRec.h"
58 #include "AliQADataMakerSim.h"
59 #include "AliQAManager.h"
60 #include "AliRawReaderDate.h"
61 #include "AliRawReaderFile.h"
62 #include "AliRawReaderRoot.h"
64 #include "AliRunLoader.h"
65 #include "AliRunTag.h"
67 ClassImp(AliQAManager)
68 AliQAManager* AliQAManager::fgQAInstance = 0x0;
70 //_____________________________________________________________________________
71 AliQAManager::AliQAManager() :
83 fNumberOfEvents(999999),
87 fRawReaderDelete(kTRUE),
90 fEventSpecie(AliRecoParam::kDefault)
93 fMaxEvents = fNumberOfEvents ;
94 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
95 if (IsSelected(AliQA::GetDetName(iDet))) {
96 fLoader[iDet] = NULL ;
97 fQADataMaker[iDet] = NULL ;
98 fQACycles[iDet] = 999999 ;
99 fQAWriteExpert[iDet] = kTRUE ;
104 //_____________________________________________________________________________
105 AliQAManager::AliQAManager(const Char_t * mode, const Char_t* gAliceFilename) :
113 fGAliceFileName(gAliceFilename),
117 fNumberOfEvents(999999),
121 fRawReaderDelete(kTRUE),
124 fEventSpecie(AliRecoParam::kDefault)
127 fMaxEvents = fNumberOfEvents ;
128 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
129 if (IsSelected(AliQA::GetDetName(iDet))) {
130 fLoader[iDet] = NULL ;
131 fQADataMaker[iDet] = NULL ;
132 fQACycles[iDet] = 999999 ;
133 fQAWriteExpert[iDet] = kTRUE ;
138 //_____________________________________________________________________________
139 AliQAManager::AliQAManager(const AliQAManager & qas) :
141 fCurrentEvent(qas.fCurrentEvent),
143 fDetectors(qas.fDetectors),
144 fDetectorsW(qas.fDetectorsW),
147 fGAliceFileName(qas.fGAliceFileName),
148 fFirstEvent(qas.fFirstEvent),
149 fMaxEvents(qas.fMaxEvents),
151 fNumberOfEvents(qas.fNumberOfEvents),
153 fRunNumber(qas.fRunNumber),
155 fRawReaderDelete(kTRUE),
158 fEventSpecie(qas.fEventSpecie)
161 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
162 fLoader[iDet] = qas.fLoader[iDet] ;
163 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
164 fQACycles[iDet] = qas.fQACycles[iDet] ;
165 fQAWriteExpert[iDet] = qas.fQAWriteExpert[iDet] ;
169 //_____________________________________________________________________________
170 AliQAManager & AliQAManager::operator = (const AliQAManager & qas)
172 // assignment operator
173 this->~AliQAManager() ;
174 new(this) AliQAManager(qas) ;
178 //_____________________________________________________________________________
179 AliQAManager::~AliQAManager()
182 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
183 if (IsSelected(AliQA::GetDetName(iDet))) {
184 fLoader[iDet] = NULL;
185 if (fQADataMaker[iDet]) {
186 (fQADataMaker[iDet])->Finish() ;
187 delete fQADataMaker[iDet] ;
191 if (fRawReaderDelete) {
198 //_____________________________________________________________________________
199 Bool_t AliQAManager::DoIt(const AliQA::TASKINDEX_t taskIndex)
201 // Runs all the QA data Maker for every detector
204 // Fill QA data in event loop
205 for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
208 if ( iEvent%10 == 0 )
209 AliInfo(Form("processing event %d", iEvent));
210 if ( taskIndex == AliQA::kRAWS ) {
211 if ( !fRawReader->NextEvent() )
213 } else if ( taskIndex == AliQA::kESDS ) {
214 if ( fESDTree->GetEntry(iEvent) == 0 )
217 if ( fRunLoader->GetEvent(iEvent) != 0 )
220 // loop over active loaders
221 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
222 if (IsSelected(AliQA::GetDetName(iDet))) {
223 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
224 if (!qadm) continue; // This detector doesn't have any QA (for example, HLT)
225 if ( qadm->IsCycleDone() ) {
226 qadm->EndOfCycle(taskIndex) ;
228 TTree * data = NULL ;
229 AliLoader* loader = GetLoader(qadm->GetUniqueID());
231 case AliQA::kNULLTASKINDEX :
234 qadm->Exec(taskIndex, fRawReader) ;
239 data = loader->TreeH() ;
241 AliWarning(Form(" Hit Tree not found for %s", AliQA::GetDetName(iDet))) ;
245 qadm->Exec(taskIndex, data) ;
247 case AliQA::kSDIGITS :
249 loader->LoadSDigits() ;
250 data = loader->TreeS() ;
252 AliWarning(Form(" SDigit Tree not found for %s", AliQA::GetDetName(iDet))) ;
256 qadm->Exec(taskIndex, data) ;
258 case AliQA::kDIGITS :
260 loader->LoadDigits() ;
261 data = loader->TreeD() ;
263 AliWarning(Form(" Digit Tree not found for %s", AliQA::GetDetName(iDet))) ;
267 qadm->Exec(taskIndex, data) ;
269 case AliQA::kRECPOINTS :
271 loader->LoadRecPoints() ;
272 data = loader->TreeR() ;
274 AliWarning(Form("RecPoints not found for %s", AliQA::GetDetName(iDet))) ;
278 qadm->Exec(taskIndex, data) ;
280 case AliQA::kTRACKSEGMENTS :
282 case AliQA::kRECPARTICLES :
285 qadm->Exec(taskIndex, fESD) ;
287 case AliQA::kNTASKINDEX :
294 // Save QA data for all detectors
295 rv = Finish(taskIndex) ;
297 if ( taskIndex == AliQA::kRAWS )
298 fRawReader->RewindEvents() ;
303 //_____________________________________________________________________________
304 Bool_t AliQAManager::Finish(const AliQA::TASKINDEX_t taskIndex)
306 // write output to file for all detectors
307 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
308 if (IsSelected(AliQA::GetDetName(iDet))) {
309 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
311 qadm->EndOfCycle(taskIndex) ;
317 //_____________________________________________________________________________
318 TObjArray * AliQAManager::GetFromOCDB(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, const char * year) const
320 // Retrieve the list of QA data for a given detector and a given task
321 TObjArray * rv = NULL ;
322 if ( !strlen(AliQA::GetQARefStorage()) ) {
323 AliError("No storage defined, use AliQA::SetQARefStorage") ;
326 if ( ! IsDefaultStorageSet() ) {
327 TString tmp(AliQA::GetQARefDefaultStorage()) ;
330 Instance()->SetDefaultStorage(tmp.Data()) ;
331 Instance()->SetSpecificStorage(Form("%s/*", AliQA::GetQAName()), AliQA::GetQARefStorage()) ;
333 TString detOCDBDir(Form("%s/%s/%s", AliQA::GetQAName(), AliQA::GetDetName((Int_t)det), AliQA::GetRefOCDBDirName())) ;
334 AliInfo(Form("Retrieving reference data from %s/%s for %s", AliQA::GetQARefStorage(), detOCDBDir.Data(), AliQA::GetTaskName(task).Data())) ;
335 AliCDBEntry* entry = QAManager()->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
336 TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
338 rv = dynamic_cast<TObjArray *>(listDetQAD->FindObject(AliQA::GetTaskName(task))) ;
342 //_____________________________________________________________________________
343 AliLoader * AliQAManager::GetLoader(Int_t iDet)
345 // get the loader for a detector
347 if ( !fRunLoader || iDet == AliQA::kCORR)
350 TString detName = AliQA::GetDetName(iDet) ;
351 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
353 return fLoader[iDet] ;
355 // load the QA data maker object
356 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
357 TString loaderName = "Ali" + detName + "Loader" ;
359 AliLoader * loader = NULL ;
360 // first check if a plugin is defined for the quality assurance data maker
361 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
362 // if not, add a plugin for it
363 if (!pluginHandler) {
364 AliDebug(1, Form("defining plugin for %s", loaderName.Data())) ;
365 TString libs = gSystem->GetLibraries() ;
366 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
367 pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
369 pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
371 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
373 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
374 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
377 fLoader[iDet] = loader ;
381 //_____________________________________________________________________________
382 AliQA * AliQAManager::GetQA(UInt_t run, UInt_t evt)
384 // retrieves the QA object stored in a file named "Run{run}.Event{evt}_1.ESD.tag.root"
385 char * fileName = Form("Run%d.Event%d_1.ESD.tag.root", run, evt) ;
386 TFile * tagFile = TFile::Open(fileName) ;
388 AliError(Form("File %s not found", fileName)) ;
391 TTree * tagTree = dynamic_cast<TTree *>(tagFile->Get("T")) ;
393 AliError(Form("Tree T not found in %s", fileName)) ;
397 AliRunTag * tag = new AliRunTag ;
398 tagTree->SetBranchAddress("AliTAG", &tag) ;
399 tagTree->GetEntry(evt) ;
400 AliQA * qa = AliQA::Instance(tag->GetQALength(), tag->GetQA(), tag->GetESLength(), tag->GetEventSpecies()) ;
405 //_____________________________________________________________________________
406 AliQADataMaker * AliQAManager::GetQADataMaker(const Int_t iDet)
408 // get the quality assurance data maker for a detector
410 if (fQADataMaker[iDet]) {
411 fQADataMaker[iDet]->SetEventSpecie(fEventSpecie) ;
412 return fQADataMaker[iDet] ;
415 AliQADataMaker * qadm = NULL ;
417 if (iDet == AliQA::kGLOBAL) { //Global QA
418 qadm = new AliGlobalQADataMaker();
419 qadm->SetName(AliQA::GetDetName(iDet));
420 qadm->SetUniqueID(iDet);
421 fQADataMaker[iDet] = qadm;
422 qadm->SetEventSpecie(fEventSpecie) ;
426 if (iDet == AliQA::kCORR) { //the data maker for correlations among detectors
427 qadm = new AliCorrQADataMakerRec(fQADataMaker) ;
428 qadm->SetName(AliQA::GetDetName(iDet));
429 qadm->SetUniqueID(iDet);
430 fQADataMaker[iDet] = qadm;
431 qadm->SetEventSpecie(fEventSpecie) ;
435 // load the QA data maker object
436 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
437 TString detName = AliQA::GetDetName(iDet) ;
439 if (tmp.Contains("sim"))
440 tmp.ReplaceAll("s", "S") ;
441 else if (tmp.Contains("rec"))
442 tmp.ReplaceAll("r", "R") ;
444 TString qadmName = "Ali" + detName + "QADataMaker" + tmp ;
446 // first check if a plugin is defined for the quality assurance data maker
447 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
448 // if not, add a plugin for it
449 if (!pluginHandler) {
450 AliDebug(1, Form("defining plugin for %s", qadmName.Data())) ;
451 TString libs = gSystem->GetLibraries() ;
452 if (libs.Contains("lib" + detName + fMode + ".so") || (gSystem->Load("lib" + detName + fMode + ".so") >= 0)) {
453 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
455 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
457 pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
459 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
460 qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
463 qadm->SetName(AliQA::GetDetName(iDet));
464 qadm->SetUniqueID(iDet);
465 fQADataMaker[iDet] = qadm ;
466 qadm->SetEventSpecie(fEventSpecie) ;
472 //_____________________________________________________________________________
473 void AliQAManager::EndOfCycle(TObjArray * detArray)
475 // End of cycle QADataMakers
477 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
478 if (IsSelected(AliQA::GetDetName(iDet))) {
479 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
482 // skip non active detectors
484 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
485 if (!det || !det->IsActive())
488 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
489 if ( fTasks.Contains(Form("%d", taskIndex)) )
490 qadm->EndOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex))) ;
497 //_____________________________________________________________________________
498 void AliQAManager::EndOfCycle(TString detectors)
500 // End of cycle QADataMakers
502 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
503 if (IsSelected(AliQA::GetDetName(iDet))) {
504 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
507 // skip non active detectors
508 if (!detectors.Contains(AliQA::GetDetName(iDet)))
510 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
511 if ( fTasks.Contains(Form("%d", taskIndex)) )
512 qadm->EndOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex))) ;
519 //_____________________________________________________________________________
520 void AliQAManager::Increment()
522 // Increments the cycle counter for all QA Data Makers
523 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
524 if (IsSelected(AliQA::GetDetName(iDet))) {
525 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
532 //_____________________________________________________________________________
533 Bool_t AliQAManager::InitQA(const AliQA::TASKINDEX_t taskIndex, const char * input )
535 // Initialize the event source and QA data makers
537 fTasks += Form("%d", taskIndex) ;
539 if (taskIndex == AliQA::kRAWS) {
541 fRawReader = AliRawReader::Create(input);
545 fRawReaderDelete = kTRUE ;
546 fRawReader->NextEvent() ;
547 fRunNumber = fRawReader->GetRunNumber() ;
549 fRawReader->RewindEvents();
550 fNumberOfEvents = 999999 ;
551 if ( fMaxEvents < 0 )
552 fMaxEvents = fNumberOfEvents ;
553 } else if (taskIndex == AliQA::kESDS) {
554 fTasks = AliQA::GetTaskName(AliQA::kESDS) ;
555 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
556 TFile * esdFile = TFile::Open("AliESDs.root") ;
557 fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ;
559 AliError("esdTree not found") ;
562 fESD = new AliESDEvent() ;
563 fESD->ReadFromTree(fESDTree) ;
564 fESDTree->GetEntry(0) ;
565 fRunNumber = fESD->GetRunNumber() ;
566 fNumberOfEvents = fESDTree->GetEntries() ;
567 if ( fMaxEvents < 0 )
568 fMaxEvents = fNumberOfEvents ;
571 AliError("AliESDs.root not found") ;
575 if ( !InitRunLoader() ) {
576 AliWarning("No Run Loader not found") ;
578 fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
579 if ( fMaxEvents < 0 )
580 fMaxEvents = fNumberOfEvents ;
585 TObjArray* detArray = NULL ;
586 if (fRunLoader) // check if RunLoader exists
587 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
588 detArray = fRunLoader->GetAliRun()->Detectors() ;
589 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
592 // Initialize all QA data makers for all detectors
593 fRunNumber = AliCDBManager::Instance()->GetRun() ;
594 if ( ! AliGeomManager::GetGeometry() )
595 AliGeomManager::LoadGeometry() ;
597 InitQADataMaker(fRunNumber, detArray) ; //, fCycleSame, kTRUE, detArray) ;
601 //_____________________________________________________________________________
602 void AliQAManager::InitQADataMaker(UInt_t run, TObjArray * detArray)
604 // Initializes The QADataMaker for all active detectors and for all active tasks
605 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
606 if (IsSelected(AliQA::GetDetName(iDet))) {
607 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
609 AliError(Form("AliQADataMaker not found for %s", AliQA::GetDetName(iDet))) ;
610 fDetectorsW.ReplaceAll(AliQA::GetDetName(iDet), "") ;
612 if (fQAWriteExpert[iDet])
613 qadm->SetWriteExpert() ;
614 AliDebug(1, Form("Data Maker found for %s", qadm->GetName())) ;
615 // skip non active detectors
617 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
618 if (!det || !det->IsActive())
621 if (fQAWriteExpert[iDet]) qadm->SetWriteExpert() ;
622 // Set default reco params
623 Bool_t sameCycle = kFALSE ;
624 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
625 if ( fTasks.Contains(Form("%d", taskIndex)) ) {
626 qadm->Init(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
627 qadm->StartOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), run, sameCycle) ;
637 //_____________________________________________________________________________
638 Bool_t AliQAManager::InitRunLoader()
640 // get or create the run loader
644 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
645 // load all base libraries to get the loader classes
646 TString libs = gSystem->GetLibraries() ;
647 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
648 if (!IsSelected(AliQA::GetDetName(iDet)))
650 TString detName = AliQA::GetDetName(iDet) ;
651 if (detName == "HLT")
653 if (libs.Contains("lib" + detName + "base.so"))
655 gSystem->Load("lib" + detName + "base.so");
657 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
659 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
662 fRunLoader->CdGAFile();
663 if (fRunLoader->LoadgAlice() == 0) {
664 gAlice = fRunLoader->GetAliRun();
668 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
672 } else { // galice.root does not exist
673 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
679 fRunLoader->LoadHeader();
680 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
685 //_____________________________________________________________________________
686 Bool_t AliQAManager::IsSelected(const char * det)
688 // check whether detName is contained in detectors
689 // if yes, it is removed from detectors
692 const TString detName(det) ;
693 // always activates Correlation
694 if ( detName.Contains(AliQA::GetDetName(AliQA::kCORR))) {
697 // check if all detectors are selected
698 if (fDetectors.Contains("ALL")) {
701 } else if ((fDetectors.CompareTo(detName) == 0) ||
702 fDetectors.BeginsWith(detName+" ") ||
703 fDetectors.EndsWith(" "+detName) ||
704 fDetectors.Contains(" "+detName+" ")) {
711 //_____________________________________________________________________________
712 Bool_t AliQAManager::Merge(const Int_t runNumber) const
714 // Merge data from all the cycles from all detectors in one single file per run
715 // Merge the QA results from all the data chunks in one run
716 Bool_t rv = MergeData(runNumber) ;
717 rv *= MergeResults(runNumber) ;
721 //______________________________________________________________________
722 Bool_t AliQAManager::MergeXML(const char * collectionFile, const char * subFile, const char * outFile)
724 // merges files listed in a xml collection
725 // usage Merge(collection, outputFile))
726 // collection: is a xml collection
730 if ( strstr(collectionFile, ".xml") == 0 ) {
731 AliError("Input collection file must be an \".xml\" file\n") ;
736 TGrid::Connect("alien://");
740 // Open the file collection
741 printf("*** Create Collection ***\n");
742 printf("*** Wk-Dir = |%s| \n",gSystem->WorkingDirectory());
743 printf("*** Coll = |%s| \n",collectionFile);
745 TGridCollection * collection = TAlienCollection::Open(collectionFile);
746 TGridResult* result = collection->GetGridResult("", 0, 0);
752 TString tempo(collectionFile) ;
754 tempo.ReplaceAll(".xml", subFile) ;
756 tempo.ReplaceAll(".xml", "_Merged.root") ;
757 outFile = tempo.Data() ;
759 merger.OutputFile(outFile) ;
761 while ( (turl = result->GetKey(index, "turl")) ) {
764 file = Form("%s#%s", turl, subFile) ;
766 file = Form("%s", turl) ;
768 AliInfo(Form("%s\n", file)) ;
769 merger.AddFile(file) ;
776 AliInfo(Form("Files merged into %s\n", outFile)) ;
782 //_____________________________________________________________________________
783 Bool_t AliQAManager::MergeData(const Int_t runNumber) const
785 // Merge all the cycles from all detectors in one single file per run
788 cmd = Form(".! ls *%s*.%d.root > tempo.txt", AliQA::GetQADataFileName(), runNumber) ;
790 cmd = Form(".! ls *%s*.*.root > tempo.txt", AliQA::GetQADataFileName()) ;
791 gROOT->ProcessLine(cmd.Data()) ;
792 ifstream in("tempo.txt") ;
793 const Int_t runMax = 10 ;
794 TString file[AliQA::kNDET*runMax] ;
801 AliInfo(Form("index = %d file = %s", index, (file[index]).Data())) ;
806 AliError(Form("run number %d not found", runNumber)) ;
811 TString outFileName ;
813 outFileName = Form("Merged.%s.Data.%d.root",AliQA::GetQADataFileName(),runNumber);
815 outFileName = Form("Merged.%s.Data.root",AliQA::GetQADataFileName());
816 merger.OutputFile(outFileName.Data()) ;
817 for (Int_t ifile = 0 ; ifile < index-1 ; ifile++) {
818 TString pattern(Form("%s.%d.", AliQA::GetQADataFileName(), runNumber));
819 TString tmp(file[ifile]) ;
820 if (tmp.Contains(pattern)) {
821 merger.AddFile(tmp) ;
828 //_____________________________________________________________________________
829 Bool_t AliQAManager::MergeResults(const Int_t runNumber) const
831 // Merge the QA result from all the data chunks in a run
833 cmd = Form(".! ls %s*.root > tempo.txt", AliQA::GetQADataFileName()) ;
834 gROOT->ProcessLine(cmd.Data()) ;
835 ifstream in("tempo.txt") ;
836 const Int_t chunkMax = 100 ;
837 TString fileList[chunkMax] ;
842 in >> fileList[index] ;
845 AliInfo(Form("index = %d file = %s", index, (fileList[index].Data()))) ;
850 AliError("No QA Result File found") ;
855 TString outFileName ;
857 outFileName = Form("Merged.%s.Result.%d.root",AliQA::GetQADataFileName(),runNumber);
859 outFileName = Form("Merged.%s.Result.root",AliQA::GetQADataFileName());
860 merger.OutputFile(outFileName.Data()) ;
861 for (Int_t ifile = 0 ; ifile < index ; ifile++) {
862 TString file = fileList[ifile] ;
863 merger.AddFile(file) ;
870 //_____________________________________________________________________________
871 void AliQAManager::Reset(const Bool_t sameCycle)
873 // Reset the default data members
875 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
876 if (IsSelected(AliQA::GetDetName(iDet))) {
877 AliQADataMaker * qadm = GetQADataMaker(iDet);
881 if (fRawReaderDelete) {
886 fCycleSame = sameCycle ;
890 fNumberOfEvents = 999999 ;
893 //_____________________________________________________________________________
894 AliQAManager * AliQAManager::QAManager(const Char_t * mode, TMap *entryCache, Int_t run)
896 // returns AliQAManager instance (singleton)
899 fgQAInstance = new AliQAManager(mode) ;
901 fgQAInstance->Init();
903 fgQAInstance->InitFromCache(entryCache,run);
908 //_____________________________________________________________________________
909 TString AliQAManager::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle)
911 //Runs all the QA data Maker for Raws only
913 fCycleSame = sameCycle ;
914 fRawReader = rawReader ;
915 fDetectors = detectors ;
916 fDetectorsW = detectors ;
918 AliCDBManager* man = AliCDBManager::Instance() ;
920 if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
921 rawReader->NextEvent() ;
922 man->SetRun(fRawReader->GetRunNumber()) ;
923 rawReader->RewindEvents() ;
927 if ( !InitQA(AliQA::kRAWS) )
929 fRawReaderDelete = kFALSE ;
935 //_____________________________________________________________________________
936 TString AliQAManager::Run(const char * detectors, const char * fileName, const Bool_t sameCycle)
938 //Runs all the QA data Maker for Raws only
940 fCycleSame = sameCycle ;
941 fDetectors = detectors ;
942 fDetectorsW = detectors ;
944 AliCDBManager* man = AliCDBManager::Instance() ;
945 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
946 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
948 AliFatal("galice.root file not found in current directory") ;
952 if ( ! rl->GetAliRun() ) {
953 AliFatal("AliRun not found in galice.root") ;
956 man->SetRun(rl->GetHeader()->GetRun());
962 if ( !InitQA(AliQA::kRAWS, fileName) )
969 //_____________________________________________________________________________
970 TString AliQAManager::Run(const char * detectors, const AliQA::TASKINDEX_t taskIndex, Bool_t const sameCycle, const char * fileName )
972 // Runs all the QA data Maker for every detector
974 fCycleSame = sameCycle ;
975 fDetectors = detectors ;
976 fDetectorsW = detectors ;
978 AliCDBManager* man = AliCDBManager::Instance() ;
979 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
980 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
982 AliFatal("galice.root file not found in current directory") ;
986 if ( ! rl->GetAliRun() ) {
987 AliInfo("AliRun not found in galice.root") ;
990 man->SetRun(rl->GetHeader()->GetRun()) ;
996 if ( taskIndex == AliQA::kNULLTASKINDEX) {
997 for (UInt_t task = 0; task < AliQA::kNTASKINDEX; task++) {
998 if ( fTasks.Contains(Form("%d", task)) ) {
1000 if ( !InitQA(AliQA::GetTaskIndex(AliQA::GetTaskName(task)), fileName) )
1002 DoIt(AliQA::GetTaskIndex(AliQA::GetTaskName(task))) ;
1007 if ( !InitQA(taskIndex, fileName) )
1012 return fDetectorsW ;
1016 //_____________________________________________________________________________
1017 void AliQAManager::RunOneEvent(AliRawReader * rawReader)
1019 //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
1022 AliCodeTimerAuto("") ;
1023 if (fTasks.Contains(Form("%d", AliQA::kRAWS))){
1024 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1025 if (!IsSelected(AliQA::GetDetName(iDet)))
1027 AliQADataMaker *qadm = GetQADataMaker(iDet);
1030 if ( qadm->IsCycleDone() ) {
1031 qadm->EndOfCycle() ;
1033 AliCodeTimerStart(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1034 qadm->SetEventSpecie(fEventSpecie) ;
1035 qadm->Exec(AliQA::kRAWS, rawReader) ;
1036 AliCodeTimerStop(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1041 //_____________________________________________________________________________
1042 void AliQAManager::RunOneEvent(AliESDEvent *& esd)
1044 //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1046 AliCodeTimerAuto("") ;
1047 if (fTasks.Contains(Form("%d", AliQA::kESDS))) {
1048 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1049 if (!IsSelected(AliQA::GetDetName(iDet)))
1051 AliQADataMaker *qadm = GetQADataMaker(iDet);
1054 if ( qadm->IsCycleDone() ) {
1055 qadm->EndOfCycle() ;
1057 AliCodeTimerStart(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1058 qadm->Exec(AliQA::kESDS, esd) ;
1059 AliCodeTimerStop(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
1064 //_____________________________________________________________________________
1065 void AliQAManager::RunOneEventInOneDetector(Int_t det, TTree * tree)
1067 // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1068 AliCodeTimerAuto("") ;
1069 if (fTasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1070 if (IsSelected(AliQA::GetDetName(det))) {
1071 AliQADataMaker *qadm = GetQADataMaker(det);
1073 if ( qadm->IsCycleDone() ) {
1074 qadm->EndOfCycle() ;
1076 AliCodeTimerStart(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
1077 qadm->Exec(AliQA::kRECPOINTS, tree) ;
1078 AliCodeTimerStop(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
1084 //_____________________________________________________________________________
1085 Bool_t AliQAManager::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const char * year, const char * detectors) const
1087 // take the locasl QA data merge into a single file and save in OCDB
1089 TString tmp(AliQA::GetQARefStorage()) ;
1090 if ( tmp.IsNull() ) {
1091 AliError("No storage defined, use AliQA::SetQARefStorage") ;
1094 if ( !(tmp.Contains(AliQA::GetLabLocalOCDB()) || tmp.Contains(AliQA::GetLabAliEnOCDB())) ) {
1095 AliError(Form("%s is a wrong storage, use %s or %s", AliQA::GetQARefStorage(), AliQA::GetLabLocalOCDB().Data(), AliQA::GetLabAliEnOCDB().Data())) ;
1098 TString sdet(detectors) ;
1101 if ( sdet.Contains("ALL") ) {
1102 rv = Merge(runNumber) ;
1105 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQA::GetQADataFileName(), runNumber)) ;
1106 inputFile = TFile::Open(inputFileName.Data()) ;
1107 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ;
1109 for (Int_t index = 0; index < AliQA::kNDET; index++) {
1110 if (sdet.Contains(AliQA::GetDetName(index))) {
1111 TString inputFileName(Form("%s.%s.%d.root", AliQA::GetDetName(index), AliQA::GetQADataFileName(), runNumber)) ;
1112 inputFile = TFile::Open(inputFileName.Data()) ;
1113 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ;
1120 //_____________________________________________________________________________
1121 Bool_t AliQAManager::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year, AliRecoParam::EventSpecie_t es) const
1123 // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1125 AliInfo(Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQA::GetQARefStorage())) ;
1126 if ( ! IsDefaultStorageSet() ) {
1127 TString tmp( AliQA::GetQARefStorage() ) ;
1128 if ( tmp.Contains(AliQA::GetLabLocalOCDB()) )
1129 Instance()->SetDefaultStorage(AliQA::GetQARefStorage()) ;
1131 TString tmp1(AliQA::GetQARefDefaultStorage()) ;
1133 tmp1.Append("?user=alidaq") ;
1134 Instance()->SetDefaultStorage(tmp1.Data()) ;
1137 Instance()->SetSpecificStorage("*", AliQA::GetQARefStorage()) ;
1139 Instance()->SetRun(runNumber);
1141 AliCDBMetaData mdr ;
1142 mdr.SetResponsible("yves schutz");
1144 for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++) {
1145 TDirectory * detDir = inputFile->GetDirectory(AliQA::GetDetName(detIndex)) ;
1147 AliInfo(Form("Entering %s", detDir->GetName())) ;
1148 AliQA::SetQARefDataDirName(es) ;
1149 TString detOCDBDir(Form("%s/%s/%s", AliQA::GetDetName(detIndex), AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName())) ;
1150 AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity()) ;
1151 TList * listDetQAD = new TList() ;
1152 TString listName(Form("%s QA data Reference", AliQA::GetDetName(detIndex))) ;
1153 mdr.SetComment(Form("%s QA stuff", AliQA::GetDetName(detIndex)));
1154 listDetQAD->SetName(listName) ;
1155 TList * taskList = detDir->GetListOfKeys() ;
1156 TIter nextTask(taskList) ;
1158 while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
1159 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ;
1160 TDirectory * esDir = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ;
1161 AliInfo(Form("Saving %s", esDir->GetName())) ;
1162 TObjArray * listTaskQAD = new TObjArray(100) ;
1163 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1164 listDetQAD->Add(listTaskQAD) ;
1165 TList * histList = esDir->GetListOfKeys() ;
1166 TIter nextHist(histList) ;
1168 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
1169 TObject * odata = esDir->Get(histKey->GetName()) ;
1171 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1173 if ( AliQA::GetExpert() == histKey->GetName() ) {
1174 TDirectory * expertDir = esDir->GetDirectory(histKey->GetName()) ;
1175 TList * expertHistList = expertDir->GetListOfKeys() ;
1176 TIter nextExpertHist(expertHistList) ;
1177 TKey * expertHistKey ;
1178 while ( (expertHistKey = dynamic_cast<TKey*>(nextExpertHist())) ) {
1179 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ;
1180 if ( !expertOdata ) {
1181 AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1183 AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1184 if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
1185 AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1186 TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ;
1187 listTaskQAD->Add(hExpertdata) ;
1192 AliInfo(Form("Adding %s", histKey->GetName())) ;
1193 if ( odata->IsA()->InheritsFrom("TH1") ) {
1194 AliInfo(Form("Adding %s", histKey->GetName())) ;
1195 TH1 * hdata = static_cast<TH1*>(odata) ;
1196 listTaskQAD->Add(hdata) ;
1201 Instance()->Put(listDetQAD, idr, &mdr) ;
1207 //_____________________________________________________________________________
1208 void AliQAManager::SetEventSpecie(AliRecoParam::EventSpecie_t es)
1210 // set the current event specie and inform AliQA that this event specie has been encountered
1212 AliQA::Instance()->SetEventSpecie(es) ;
1215 //_____________________________________________________________________________
1216 void AliQAManager::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par)
1218 // Set custom reconstruction parameters for a given detector
1219 // Single set of parameters for all the events
1220 GetQADataMaker(det)->SetRecoParam(par) ;