]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliQAManager.cxx
Warnings and Coverity fixes (Plamen)
[u/mrichter/AliRoot.git] / STEER / AliQAManager.cxx
CommitLineData
2e331c8b 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// //
bf76b847 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); //
2e331c8b 28// //
29///////////////////////////////////////////////////////////////////////////////
30
fec0891b 31#include <TCanvas.h>
2e331c8b 32#include <TKey.h>
33#include <TFile.h>
34#include <TFileMerger.h>
149a2367 35#include <TGrid.h>
36#include <TGridCollection.h>
37#include <TGridResult.h>
2e331c8b 38#include <TPluginManager.h>
39#include <TROOT.h>
40#include <TString.h>
41#include <TSystem.h>
7c9ff472 42#include <TStopwatch.h>
2e331c8b 43
44#include "AliCDBManager.h"
45#include "AliCDBEntry.h"
46#include "AliCDBId.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"
55#include "AliLog.h"
56#include "AliModule.h"
4e25ac79 57#include "AliQAv1.h"
634696f5 58#include "AliQAChecker.h"
59#include "AliQACheckerBase.h"
2e331c8b 60#include "AliQADataMakerRec.h"
61#include "AliQADataMakerSim.h"
62#include "AliQAManager.h"
63#include "AliRawReaderDate.h"
64#include "AliRawReaderFile.h"
65#include "AliRawReaderRoot.h"
66#include "AliRun.h"
67#include "AliRunLoader.h"
68#include "AliRunTag.h"
69
70ClassImp(AliQAManager)
71AliQAManager* AliQAManager::fgQAInstance = 0x0;
72
73//_____________________________________________________________________________
74AliQAManager::AliQAManager() :
75 AliCDBManager(),
76 fCurrentEvent(0),
77 fCycleSame(kFALSE),
78 fDetectors("ALL"),
79 fDetectorsW("ALL"),
80 fESD(NULL),
81 fESDTree(NULL),
82 fGAliceFileName(""),
83 fFirstEvent(0),
84 fMaxEvents(0),
85 fMode(""),
86 fNumberOfEvents(999999),
87 fRecoParam(),
88 fRunNumber(0),
89 fRawReader(NULL),
90 fRawReaderDelete(kTRUE),
91 fRunLoader(NULL),
bc8761a0 92 fTasks(""),
fec0891b 93 fEventSpecie(AliRecoParam::kDefault),
50dee02c 94 fPrintImage(kTRUE),
95 fSaveData(kTRUE)
2e331c8b 96{
97 // default ctor
98 fMaxEvents = fNumberOfEvents ;
99 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 100 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 101 fLoader[iDet] = NULL ;
102 fQADataMaker[iDet] = NULL ;
103 fQACycles[iDet] = 999999 ;
2e331c8b 104 }
105 }
75373542 106 SetWriteExpert() ;
2e331c8b 107}
108
109//_____________________________________________________________________________
634696f5 110AliQAManager::AliQAManager(AliQAv1::MODE_t mode, const Char_t* gAliceFilename) :
2e331c8b 111 AliCDBManager(),
112 fCurrentEvent(0),
113 fCycleSame(kFALSE),
114 fDetectors("ALL"),
115 fDetectorsW("ALL"),
116 fESD(NULL),
117 fESDTree(NULL),
118 fGAliceFileName(gAliceFilename),
119 fFirstEvent(0),
120 fMaxEvents(0),
634696f5 121 fMode(AliQAv1::GetModeName(mode)),
2e331c8b 122 fNumberOfEvents(999999),
123 fRecoParam(),
124 fRunNumber(0),
125 fRawReader(NULL),
126 fRawReaderDelete(kTRUE),
127 fRunLoader(NULL),
bc8761a0 128 fTasks(""),
fec0891b 129 fEventSpecie(AliRecoParam::kDefault),
50dee02c 130 fPrintImage(kTRUE),
131 fSaveData(kTRUE)
2e331c8b 132{
133 // default ctor
134 fMaxEvents = fNumberOfEvents ;
135 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 136 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 137 fLoader[iDet] = NULL ;
138 fQADataMaker[iDet] = NULL ;
139 fQACycles[iDet] = 999999 ;
65660bf3 140 }
141 }
75373542 142 SetWriteExpert() ;
2e331c8b 143}
144
145//_____________________________________________________________________________
146AliQAManager::AliQAManager(const AliQAManager & qas) :
147 AliCDBManager(),
148 fCurrentEvent(qas.fCurrentEvent),
149 fCycleSame(kFALSE),
150 fDetectors(qas.fDetectors),
151 fDetectorsW(qas.fDetectorsW),
152 fESD(NULL),
153 fESDTree(NULL),
154 fGAliceFileName(qas.fGAliceFileName),
155 fFirstEvent(qas.fFirstEvent),
156 fMaxEvents(qas.fMaxEvents),
157 fMode(qas.fMode),
158 fNumberOfEvents(qas.fNumberOfEvents),
159 fRecoParam(),
160 fRunNumber(qas.fRunNumber),
161 fRawReader(NULL),
162 fRawReaderDelete(kTRUE),
163 fRunLoader(NULL),
bc8761a0 164 fTasks(qas.fTasks),
fec0891b 165 fEventSpecie(qas.fEventSpecie),
50dee02c 166 fPrintImage(qas.fPrintImage),
167 fSaveData(qas.fSaveData)
fec0891b 168
2e331c8b 169{
170 // cpy ctor
171 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
172 fLoader[iDet] = qas.fLoader[iDet] ;
173 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
174 fQACycles[iDet] = qas.fQACycles[iDet] ;
175 fQAWriteExpert[iDet] = qas.fQAWriteExpert[iDet] ;
176 }
177}
178
179//_____________________________________________________________________________
180AliQAManager & AliQAManager::operator = (const AliQAManager & qas)
181{
182 // assignment operator
183 this->~AliQAManager() ;
184 new(this) AliQAManager(qas) ;
185 return *this ;
186}
187
188//_____________________________________________________________________________
189AliQAManager::~AliQAManager()
190{
191 // dtor
192 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 193 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 194 fLoader[iDet] = NULL;
195 if (fQADataMaker[iDet]) {
196 (fQADataMaker[iDet])->Finish() ;
197 delete fQADataMaker[iDet] ;
198 }
199 }
200 }
201 if (fRawReaderDelete) {
202 fRunLoader = NULL ;
203 delete fRawReader ;
204 fRawReader = NULL ;
205 }
206}
2e331c8b 207//_____________________________________________________________________________
4e25ac79 208Bool_t AliQAManager::DoIt(const AliQAv1::TASKINDEX_t taskIndex)
2e331c8b 209{
210 // Runs all the QA data Maker for every detector
211
212 Bool_t rv = kFALSE ;
213 // Fill QA data in event loop
214 for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
215 fCurrentEvent++ ;
216 // Get the event
217 if ( iEvent%10 == 0 )
5379c4a3 218 AliDebug(AliQAv1::GetQADebugLevel(), Form("processing event %d", iEvent));
4e25ac79 219 if ( taskIndex == AliQAv1::kRAWS ) {
2e331c8b 220 if ( !fRawReader->NextEvent() )
221 break ;
4e25ac79 222 } else if ( taskIndex == AliQAv1::kESDS ) {
2e331c8b 223 if ( fESDTree->GetEntry(iEvent) == 0 )
224 break ;
225 } else {
226 if ( fRunLoader->GetEvent(iEvent) != 0 )
227 break ;
228 }
229 // loop over active loaders
230 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 231 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 232 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
233 if (!qadm) continue; // This detector doesn't have any QA (for example, HLT)
234 if ( qadm->IsCycleDone() ) {
235 qadm->EndOfCycle(taskIndex) ;
236 }
237 TTree * data = NULL ;
238 AliLoader* loader = GetLoader(qadm->GetUniqueID());
239 switch (taskIndex) {
4e25ac79 240 case AliQAv1::kNULLTASKINDEX :
2e331c8b 241 break ;
4e25ac79 242 case AliQAv1::kRAWS :
2e331c8b 243 qadm->Exec(taskIndex, fRawReader) ;
244 break ;
4e25ac79 245 case AliQAv1::kHITS :
2e331c8b 246 if( loader ) {
247 loader->LoadHits() ;
248 data = loader->TreeH() ;
249 if ( ! data ) {
4e25ac79 250 AliWarning(Form(" Hit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
2e331c8b 251 break ;
252 }
eca4fa66 253 qadm->Exec(taskIndex, data) ;
2e331c8b 254 }
2e331c8b 255 break ;
4e25ac79 256 case AliQAv1::kSDIGITS :
95144eea 257 {
258
259 TString fileName(Form("%s.SDigits.root", AliQAv1::GetDetName(iDet))) ;
260 if (gSystem->FindFile("./", fileName)) {
261 if( loader ) {
262 loader->LoadSDigits() ;
263 data = loader->TreeS() ;
264 if ( ! data ) {
265 AliWarning(Form(" SDigit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
266 break ;
267 }
268 qadm->Exec(taskIndex, data) ;
269 }
2e331c8b 270 }
95144eea 271 }
2e331c8b 272 break;
4e25ac79 273 case AliQAv1::kDIGITS :
fec0891b 274 if( loader ) {
275 loader->LoadDigits() ;
276 data = loader->TreeD() ;
277 if ( ! data ) {
278 AliWarning(Form(" Digit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
279 break ;
280 }
eca4fa66 281 qadm->Exec(taskIndex, data) ;
fec0891b 282 }
eca4fa66 283 break;
fec0891b 284 case AliQAv1::kDIGITSR :
285 if( loader ) {
286 loader->LoadDigits() ;
287 data = loader->TreeD() ;
288 if ( ! data ) {
289 AliWarning(Form(" Digit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
290 break ;
291 }
eca4fa66 292 qadm->Exec(taskIndex, data) ;
fec0891b 293 }
2e331c8b 294 break;
4e25ac79 295 case AliQAv1::kRECPOINTS :
2e331c8b 296 if( loader ) {
297 loader->LoadRecPoints() ;
298 data = loader->TreeR() ;
299 if (!data) {
4e25ac79 300 AliWarning(Form("RecPoints not found for %s", AliQAv1::GetDetName(iDet))) ;
2e331c8b 301 break ;
302 }
eca4fa66 303 qadm->Exec(taskIndex, data) ;
2e331c8b 304 }
2e331c8b 305 break;
4e25ac79 306 case AliQAv1::kTRACKSEGMENTS :
2e331c8b 307 break;
4e25ac79 308 case AliQAv1::kRECPARTICLES :
2e331c8b 309 break;
4e25ac79 310 case AliQAv1::kESDS :
2e331c8b 311 qadm->Exec(taskIndex, fESD) ;
312 break;
4e25ac79 313 case AliQAv1::kNTASKINDEX :
2e331c8b 314 break;
315 } //task switch
316 }
317 } // detector loop
5e303886 318 Increment(taskIndex) ;
2e331c8b 319 } // event loop
320 // Save QA data for all detectors
5e303886 321
1c00473b 322 EndOfCycle() ;
2e331c8b 323
4e25ac79 324 if ( taskIndex == AliQAv1::kRAWS )
2e331c8b 325 fRawReader->RewindEvents() ;
326
327 return rv ;
328}
329
330//_____________________________________________________________________________
4e25ac79 331Bool_t AliQAManager::Finish(const AliQAv1::TASKINDEX_t taskIndex)
2e331c8b 332{
333 // write output to file for all detectors
8fa875cb 334
335 AliQAChecker::Instance()->SetRunNumber(fRunNumber) ;
336
2e331c8b 337 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 338 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 339 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
340 if (qadm)
341 qadm->EndOfCycle(taskIndex) ;
342 }
343 }
344 return kTRUE ;
345}
346
347//_____________________________________________________________________________
fec0891b 348TObjArray * AliQAManager::GetFromOCDB(AliQAv1::DETECTORINDEX_t det, AliQAv1::TASKINDEX_t task, const Char_t * year) const
2e331c8b 349{
350 // Retrieve the list of QA data for a given detector and a given task
351 TObjArray * rv = NULL ;
4e25ac79 352 if ( !strlen(AliQAv1::GetQARefStorage()) ) {
353 AliError("No storage defined, use AliQAv1::SetQARefStorage") ;
2e331c8b 354 return NULL ;
355 }
356 if ( ! IsDefaultStorageSet() ) {
4e25ac79 357 TString tmp(AliQAv1::GetQARefDefaultStorage()) ;
2e331c8b 358 tmp.Append(year) ;
359 tmp.Append("/") ;
360 Instance()->SetDefaultStorage(tmp.Data()) ;
4e25ac79 361 Instance()->SetSpecificStorage(Form("%s/*", AliQAv1::GetQAName()), AliQAv1::GetQARefStorage()) ;
2e331c8b 362 }
4e25ac79 363 TString detOCDBDir(Form("%s/%s/%s", AliQAv1::GetQAName(), AliQAv1::GetDetName((Int_t)det), AliQAv1::GetRefOCDBDirName())) ;
5379c4a3 364 AliDebug(AliQAv1::GetQADebugLevel(), Form("Retrieving reference data from %s/%s for %s", AliQAv1::GetQARefStorage(), detOCDBDir.Data(), AliQAv1::GetTaskName(task).Data())) ;
2e331c8b 365 AliCDBEntry* entry = QAManager()->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
eca4fa66 366 TList * listDetQAD = static_cast<TList *>(entry->GetObject()) ;
2e331c8b 367 if ( listDetQAD )
eca4fa66 368 rv = static_cast<TObjArray *>(listDetQAD->FindObject(AliQAv1::GetTaskName(task))) ;
2e331c8b 369 return rv ;
370}
371
fec0891b 372//_____________________________________________________________________________
373TCanvas ** AliQAManager::GetImage(Char_t * detName)
374{
375 // retrieves QA Image for the given detector
376 TCanvas ** rv = NULL ;
377 Int_t detIndex = AliQAv1::GetDetIndex(detName) ;
a5b7dd2a 378 if ( detIndex != AliQAv1::kNULLDET) {
379 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(detIndex) ;
380 rv = qac->GetImage() ;
381 }
fec0891b 382 return rv ;
383}
384
2e331c8b 385//_____________________________________________________________________________
386AliLoader * AliQAManager::GetLoader(Int_t iDet)
387{
388 // get the loader for a detector
389
bf76b847 390 if ( !fRunLoader || iDet == AliQAv1::kCORR || iDet == AliQAv1::kGLOBAL )
2e331c8b 391 return NULL ;
392
4e25ac79 393 TString detName = AliQAv1::GetDetName(iDet) ;
2e331c8b 394 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
395 if (fLoader[iDet])
396 return fLoader[iDet] ;
397
398 // load the QA data maker object
399 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
400 TString loaderName = "Ali" + detName + "Loader" ;
401
402 AliLoader * loader = NULL ;
403 // first check if a plugin is defined for the quality assurance data maker
404 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
405 // if not, add a plugin for it
406 if (!pluginHandler) {
5379c4a3 407 AliDebug(AliQAv1::GetQADebugLevel(), Form("defining plugin for %s", loaderName.Data())) ;
2e331c8b 408 TString libs = gSystem->GetLibraries() ;
409 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
410 pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
411 } else {
412 pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
413 }
414 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
415 }
416 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
417 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
418 }
419 if (loader)
420 fLoader[iDet] = loader ;
421 return loader ;
422}
423
424//_____________________________________________________________________________
4e25ac79 425AliQAv1 * AliQAManager::GetQA(UInt_t run, UInt_t evt)
2e331c8b 426{
427// retrieves the QA object stored in a file named "Run{run}.Event{evt}_1.ESD.tag.root"
fec0891b 428 Char_t * fileName = Form("Run%d.Event%d_1.ESD.tag.root", run, evt) ;
2e331c8b 429 TFile * tagFile = TFile::Open(fileName) ;
430 if ( !tagFile ) {
431 AliError(Form("File %s not found", fileName)) ;
432 return NULL ;
433 }
eca4fa66 434 TTree * tagTree = static_cast<TTree *>(tagFile->Get("T")) ;
2e331c8b 435 if ( !tagTree ) {
436 AliError(Form("Tree T not found in %s", fileName)) ;
437 tagFile->Close() ;
438 return NULL ;
439 }
440 AliRunTag * tag = new AliRunTag ;
441 tagTree->SetBranchAddress("AliTAG", &tag) ;
442 tagTree->GetEntry(evt) ;
4e25ac79 443 AliQAv1 * qa = AliQAv1::Instance(tag->GetQALength(), tag->GetQAArray(), tag->GetESLength(), tag->GetEventSpecies()) ;
2e331c8b 444 tagFile->Close() ;
445 return qa ;
446}
447
448//_____________________________________________________________________________
449AliQADataMaker * AliQAManager::GetQADataMaker(const Int_t iDet)
450{
451 // get the quality assurance data maker for a detector
452
eca4fa66 453 AliQADataMaker * qadm = fQADataMaker[iDet] ;
454
455 if (qadm) {
456
5e232cd6 457 qadm->SetEventSpecie(fEventSpecie) ;
bc8761a0 458 if ( qadm->GetRecoParam() )
eca4fa66 459 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
5e232cd6 460 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
2e331c8b 461
5aae1dab 462 } else if (iDet == AliQAv1::kGLOBAL && strcmp(GetMode(), AliQAv1::GetModeName(AliQAv1::kRECMODE)) == 0) { //Global QA
eca4fa66 463
464 qadm = new AliGlobalQADataMaker();
4e25ac79 465 qadm->SetName(AliQAv1::GetDetName(iDet));
2e331c8b 466 qadm->SetUniqueID(iDet);
467 fQADataMaker[iDet] = qadm;
5e232cd6 468 qadm->SetEventSpecie(fEventSpecie) ;
bc8761a0 469 if ( qadm->GetRecoParam() )
5e232cd6 470 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
471 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
2e331c8b 472
5aae1dab 473 } else if (iDet == AliQAv1::kCORR && strcmp(GetMode(), AliQAv1::GetModeName(AliQAv1::kRECMODE)) == 0 ) { //the data maker for correlations among detectors
eca4fa66 474 qadm = new AliCorrQADataMakerRec(fQADataMaker) ;
4e25ac79 475 qadm->SetName(AliQAv1::GetDetName(iDet));
2e331c8b 476 qadm->SetUniqueID(iDet);
eca4fa66 477 fQADataMaker[iDet] = qadm;
5e232cd6 478 qadm->SetEventSpecie(fEventSpecie) ;
bc8761a0 479 if ( qadm->GetRecoParam() )
5e232cd6 480 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
481 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
2e331c8b 482
6090cc85 483 } else if ( iDet < AliQAv1::kGLOBAL ) {
5aae1dab 484 TString smode(GetMode()) ;
485 if (smode.Contains(AliQAv1::GetModeName(AliQAv1::kQAMODE)))
486 smode = AliQAv1::GetModeName(AliQAv1::kRECMODE) ;
eca4fa66 487 // load the QA data maker object
488 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
489 TString detName = AliQAv1::GetDetName(iDet) ;
5aae1dab 490 TString qadmName = "Ali" + detName + "QADataMaker" + smode ;
eca4fa66 491
492 // first check if a plugin is defined for the quality assurance data maker
493 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
494 // if not, add a plugin for it
495 if (!pluginHandler) {
496 AliDebug(AliQAv1::GetQADebugLevel(), Form("defining plugin for %s", qadmName.Data())) ;
497 TString libs = gSystem->GetLibraries() ;
5aae1dab 498 TString temp(smode) ;
bf76b847 499 temp.ToLower() ;
5aae1dab 500 if (libs.Contains("lib" + detName + smode + ".so") || (gSystem->Load("lib" + detName + temp.Data() + ".so") >= 0)) {
eca4fa66 501 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
502 } else {
503 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
504 }
505 pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
506 }
507 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
508 qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
509 }
510 if (qadm) {
511 qadm->SetName(AliQAv1::GetDetName(iDet));
512 qadm->SetUniqueID(iDet);
513 fQADataMaker[iDet] = qadm ;
514 qadm->SetEventSpecie(fEventSpecie) ;
515 if ( qadm->GetRecoParam() )
516 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
517 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
518 }
519 }
2e331c8b 520 return qadm ;
521}
522
523//_____________________________________________________________________________
524void AliQAManager::EndOfCycle(TObjArray * detArray)
525{
526 // End of cycle QADataMakers
527
8fa875cb 528 AliQAChecker::Instance()->SetRunNumber(fRunNumber) ;
cd2a0b54 529 TCanvas fakeCanvas ;
1c00473b 530
1d46dcf4 531 fakeCanvas.Print(Form("%s%s%d.%s[", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps") ;
2e331c8b 532 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 533 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 534 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
535 if (!qadm)
536 continue ;
537 // skip non active detectors
538 if (detArray) {
4e25ac79 539 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQAv1::GetDetName(iDet))) ;
2e331c8b 540 if (!det || !det->IsActive())
541 continue ;
542 }
1c00473b 543 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(iDet) ;
6090cc85 544 if (qac)
545 qac->SetPrintImage(fPrintImage) ;
d338b679 546 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
547 if ( fTasks.Contains(Form("%d", taskIndex)) )
548 qadm->EndOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex))) ;
50dee02c 549 }
2e331c8b 550 qadm->Finish();
551 }
552 }
cd2a0b54 553 if (fPrintImage)
554 fakeCanvas.Print(Form("%s%s%d.%s]", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps");
2e331c8b 555}
556
557//_____________________________________________________________________________
558void AliQAManager::EndOfCycle(TString detectors)
559{
560 // End of cycle QADataMakers
561
634696f5 562 AliQAChecker::Instance()->SetRunNumber(fRunNumber) ;
cd2a0b54 563 TCanvas fakeCanvas ;
564 if (fPrintImage)
1d46dcf4 565 fakeCanvas.Print(Form("%s%s%d.%s[", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps") ;
fec0891b 566 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 567 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 568 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
569 if (!qadm)
570 continue ;
571 // skip non active detectors
4e25ac79 572 if (!detectors.Contains(AliQAv1::GetDetName(iDet)))
2e331c8b 573 continue ;
1c00473b 574 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(iDet) ;
6090cc85 575 if (qac)
576 qac->SetPrintImage(fPrintImage) ;
d338b679 577 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
578 if ( fTasks.Contains(Form("%d", taskIndex)) )
579 qadm->EndOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex))) ;
50dee02c 580 }
2e331c8b 581 qadm->Finish();
582 }
583 }
cd2a0b54 584 if (fPrintImage)
585 fakeCanvas.Print(Form("%s%s%d.%s]", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps");
2e331c8b 586}
587
d34d00a9 588//_____________________________________________________________________________
589AliRecoParam::EventSpecie_t AliQAManager::GetEventSpecieFromESD()
590{
591 AliRecoParam::EventSpecie_t runtype = AliRecoParam::kDefault ;
592 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
593 TFile * esdFile = TFile::Open("AliESDs.root") ;
594 TTree * esdTree = static_cast<TTree *> (esdFile->Get("esdTree")) ;
595 if ( !esdTree ) {
596 AliError("esdTree not found") ;
597 } else {
598 AliESDEvent * esd = new AliESDEvent() ;
599 esd->ReadFromTree(esdTree) ;
600 esdTree->GetEntry(0) ;
bad8f87c 601 runtype = AliRecoParam::Convert(esd->GetEventType()) ;
d34d00a9 602 }
603 } else {
604 AliError("AliESDs.root not found") ;
605 }
606 return runtype ;
607}
608
2e331c8b 609//_____________________________________________________________________________
5e303886 610void AliQAManager::Increment(const AliQAv1::TASKINDEX_t taskIndex)
2e331c8b 611{
612 // Increments the cycle counter for all QA Data Makers
5e303886 613 static AliQAv1::TASKINDEX_t currentTask = AliQAv1::kNTASKINDEX ;
de402dc0 614 if ( (currentTask == taskIndex) && taskIndex != AliQAv1::kNULLTASKINDEX )
5e303886 615 return ;
616 else
617 currentTask = taskIndex ;
2e331c8b 618 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 619 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 620 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
621 if (qadm)
622 qadm->Increment() ;
623 }
624 }
625}
626
627//_____________________________________________________________________________
fec0891b 628Bool_t AliQAManager::InitQA(const AliQAv1::TASKINDEX_t taskIndex, const Char_t * input )
2e331c8b 629{
630 // Initialize the event source and QA data makers
631
632 fTasks += Form("%d", taskIndex) ;
633
4e25ac79 634 if (taskIndex == AliQAv1::kRAWS) {
2e331c8b 635 if (!fRawReader) {
636 fRawReader = AliRawReader::Create(input);
637 }
638 if ( ! fRawReader )
639 return kFALSE ;
640 fRawReaderDelete = kTRUE ;
641 fRawReader->NextEvent() ;
642 fRunNumber = fRawReader->GetRunNumber() ;
643 SetRun(fRunNumber) ;
644 fRawReader->RewindEvents();
645 fNumberOfEvents = 999999 ;
646 if ( fMaxEvents < 0 )
647 fMaxEvents = fNumberOfEvents ;
4e25ac79 648 } else if (taskIndex == AliQAv1::kESDS) {
649 fTasks = AliQAv1::GetTaskName(AliQAv1::kESDS) ;
2e331c8b 650 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
651 TFile * esdFile = TFile::Open("AliESDs.root") ;
eca4fa66 652 fESDTree = static_cast<TTree *> (esdFile->Get("esdTree")) ;
2e331c8b 653 if ( !fESDTree ) {
654 AliError("esdTree not found") ;
655 return kFALSE ;
656 } else {
657 fESD = new AliESDEvent() ;
658 fESD->ReadFromTree(fESDTree) ;
659 fESDTree->GetEntry(0) ;
660 fRunNumber = fESD->GetRunNumber() ;
661 fNumberOfEvents = fESDTree->GetEntries() ;
662 if ( fMaxEvents < 0 )
663 fMaxEvents = fNumberOfEvents ;
664 }
665 } else {
666 AliError("AliESDs.root not found") ;
667 return kFALSE ;
668 }
669 } else {
670 if ( !InitRunLoader() ) {
671 AliWarning("No Run Loader not found") ;
672 } else {
673 fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
674 if ( fMaxEvents < 0 )
675 fMaxEvents = fNumberOfEvents ;
676 }
677 }
678
679 // Get Detectors
680 TObjArray* detArray = NULL ;
681 if (fRunLoader) // check if RunLoader exists
682 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
683 detArray = fRunLoader->GetAliRun()->Detectors() ;
684 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
685 }
686
687 // Initialize all QA data makers for all detectors
688 fRunNumber = AliCDBManager::Instance()->GetRun() ;
689 if ( ! AliGeomManager::GetGeometry() )
690 AliGeomManager::LoadGeometry() ;
691
692 InitQADataMaker(fRunNumber, detArray) ; //, fCycleSame, kTRUE, detArray) ;
fec0891b 693 if (fPrintImage) {
694 TCanvas fakeCanvas ;
7c9ff472 695 TStopwatch timer ;
696 timer.Start() ;
697 while (timer.CpuTime()<5) {
698 timer.Continue();
699 gSystem->ProcessEvents();
700 }
1d46dcf4 701 fakeCanvas.Print(Form("%s%s%d.%s[", AliQAv1::GetImageFileName(), GetMode(), fRunNumber, AliQAv1::GetImageFileFormat()), "ps") ;
fec0891b 702 }
2e331c8b 703 return kTRUE ;
704}
705
706//_____________________________________________________________________________
707void AliQAManager::InitQADataMaker(UInt_t run, TObjArray * detArray)
708{
709 // Initializes The QADataMaker for all active detectors and for all active tasks
fec0891b 710 fRunNumber = run ;
2e331c8b 711 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 712 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 713 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
714 if (!qadm) {
4e25ac79 715 AliError(Form("AliQADataMaker not found for %s", AliQAv1::GetDetName(iDet))) ;
716 fDetectorsW.ReplaceAll(AliQAv1::GetDetName(iDet), "") ;
2e331c8b 717 } else {
718 if (fQAWriteExpert[iDet])
719 qadm->SetWriteExpert() ;
5379c4a3 720 AliDebug(AliQAv1::GetQADebugLevel(), Form("Data Maker found for %s %d", qadm->GetName(), qadm->WriteExpert())) ;
2e331c8b 721 // skip non active detectors
722 if (detArray) {
4e25ac79 723 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQAv1::GetDetName(iDet))) ;
2e331c8b 724 if (!det || !det->IsActive())
725 continue ;
726 }
2e331c8b 727 // Set default reco params
728 Bool_t sameCycle = kFALSE ;
4e25ac79 729 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
2e331c8b 730 if ( fTasks.Contains(Form("%d", taskIndex)) ) {
4e25ac79 731 qadm->Init(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
732 qadm->StartOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex)), run, sameCycle) ;
2e331c8b 733 sameCycle = kTRUE ;
734 }
735 }
736 }
737 }
738 }
739}
740
741
742//_____________________________________________________________________________
743Bool_t AliQAManager::InitRunLoader()
744{
745 // get or create the run loader
746 if (fRunLoader) {
747 fCycleSame = kTRUE ;
748 } else {
749 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
750 // load all base libraries to get the loader classes
751 TString libs = gSystem->GetLibraries() ;
752 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 753 if (!IsSelected(AliQAv1::GetDetName(iDet)))
2e331c8b 754 continue ;
4e25ac79 755 TString detName = AliQAv1::GetDetName(iDet) ;
2e331c8b 756 if (detName == "HLT")
757 continue;
758 if (libs.Contains("lib" + detName + "base.so"))
759 continue;
760 gSystem->Load("lib" + detName + "base.so");
761 }
762 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
763 if (!fRunLoader) {
764 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
765 return kFALSE;
766 }
767 fRunLoader->CdGAFile();
768 if (fRunLoader->LoadgAlice() == 0) {
769 gAlice = fRunLoader->GetAliRun();
770 }
771
772 if (!gAlice) {
773 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
774 return kFALSE;
775 }
776
777 } else { // galice.root does not exist
778 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
779 return kFALSE;
780 }
781 }
782
783 if (!fRunNumber) {
784 fRunLoader->LoadHeader();
785 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
786 }
787 return kTRUE;
788}
789
790//_____________________________________________________________________________
fec0891b 791Bool_t AliQAManager::IsSelected(const Char_t * det)
2e331c8b 792{
793 // check whether detName is contained in detectors
794 // if yes, it is removed from detectors
795
796 Bool_t rv = kFALSE;
797 const TString detName(det) ;
798 // always activates Correlation
d282dbf5 799 if ( detName.Contains(AliQAv1::GetDetName(AliQAv1::kCORR)) || detName.Contains(AliQAv1::GetDetName(AliQAv1::kGLOBAL))) {
2e331c8b 800 rv = kTRUE ;
801 } else {
802 // check if all detectors are selected
803 if (fDetectors.Contains("ALL")) {
804 fDetectors = "ALL";
805 rv = kTRUE;
806 } else if ((fDetectors.CompareTo(detName) == 0) ||
807 fDetectors.BeginsWith(detName+" ") ||
808 fDetectors.EndsWith(" "+detName) ||
809 fDetectors.Contains(" "+detName+" ")) {
810 rv = kTRUE;
811 }
812 }
813 return rv ;
814}
815
816//_____________________________________________________________________________
87da0921 817Bool_t AliQAManager::Merge(Int_t runNumber, const char *fileName) const
2e331c8b 818{
bededc30 819 // Merge data from all detectors from a given run in one single file
87da0921 820 // Merge the QA results from all the data chunks in one run
821 // The 'fileName' is name of the output file with merged QA data
bededc30 822 if ( runNumber == -1)
823 runNumber = fRunNumber ;
87da0921 824 Bool_t rv = MergeData(runNumber,fileName) ;
bededc30 825 //rv *= MergeResults(runNumber) ; // not needed for the time being
2e331c8b 826 return rv ;
827}
828
149a2367 829//______________________________________________________________________
fec0891b 830Bool_t AliQAManager::MergeXML(const Char_t * collectionFile, const Char_t * subFile, const Char_t * outFile)
149a2367 831{
832 // merges files listed in a xml collection
833 // usage Merge(collection, outputFile))
834 // collection: is a xml collection
835
836 Bool_t rv = kFALSE ;
837
838 if ( strstr(collectionFile, ".xml") == 0 ) {
839 AliError("Input collection file must be an \".xml\" file\n") ;
840 return kFALSE ;
841 }
842
843 if ( !gGrid )
844 TGrid::Connect("alien://");
845 if ( !gGrid )
846 return kFALSE ;
847
848 // Open the file collection
5e232cd6 849 AliInfoClass(Form("*** Create Collection ***\n*** Wk-Dir = |%s| \n*** Coll = |%s| \n",gSystem->WorkingDirectory(), collectionFile));
149a2367 850
bededc30 851 TGridCollection * collection = (TGridCollection*)gROOT->ProcessLine(Form("TAlienCollection::Open(\"%s\")",collectionFile));
149a2367 852 TGridResult* result = collection->GetGridResult("", 0, 0);
853
854 Int_t index = 0 ;
fec0891b 855 const Char_t * turl ;
adae62c9 856 TFileMerger merger(kFALSE) ;
149a2367 857 if (!outFile) {
858 TString tempo(collectionFile) ;
859 if ( subFile)
860 tempo.ReplaceAll(".xml", subFile) ;
861 else
862 tempo.ReplaceAll(".xml", "_Merged.root") ;
863 outFile = tempo.Data() ;
864 }
865 merger.OutputFile(outFile) ;
866
867 while ( (turl = result->GetKey(index, "turl")) ) {
fec0891b 868 Char_t * file ;
149a2367 869 if ( subFile )
870 file = Form("%s#%s", turl, subFile) ;
871 else
872 file = Form("%s", turl) ;
873
5379c4a3 874 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s\n", file)) ;
149a2367 875 merger.AddFile(file) ;
876 index++ ;
877 }
878
879 if (index)
880 merger.Merge() ;
881
5379c4a3 882 AliDebug(AliQAv1::GetQADebugLevel(), Form("Files merged into %s\n", outFile)) ;
149a2367 883
884 rv = kFALSE;
885 return rv ;
886}
887
69a898e3 888//_____________________________________________________________________________
889void AliQAManager::MergeCustom() const
890{
891 // Custom Merge of QA data from all detectors for all runs in one single file
892 // search all the run numbers
893 // search all the run numbers
894 gROOT->ProcessLine(".! ls *QA*.root > QAtempo.txt") ;
895 TString QAfile ;
896 FILE * QAfiles = fopen("QAtempo.txt", "r") ;
897 Int_t index = 0 ;
898 TList srunList ;
899 TIter nextRun(&srunList) ;
900 TObjString * srun = NULL ;
901 Int_t loRun = 999999999 ;
902 Int_t hiRun = 0 ;
903 while ( QAfile.Gets(QAfiles) ) {
904 Bool_t runExist = kFALSE ;
905 TString srunNew(QAfile(QAfile.Index("QA.")+3, QAfile.Index(".root")-(QAfile.Index("QA.")+3))) ;
906 Int_t cuRun = srunNew.Atoi() ;
907 if (cuRun < loRun)
908 loRun = cuRun ;
909 if (cuRun > hiRun)
910 hiRun = cuRun ;
eca4fa66 911 while ( (srun = static_cast<TObjString *> (nextRun())) ) {
69a898e3 912 if ( cuRun == (srun->String()).Atoi() ) {
913 runExist = kTRUE ;
914 break ;
915 }
916 }
917 nextRun.Reset() ;
918 if ( ! runExist )
919 srunList.Add(new TObjString(srunNew.Data()));
920 }
921 nextRun.Reset() ;
922 Int_t runNumber = 0 ;
4e25ac79 923 TFile mergedFile(Form("Merged.%s.Data.root", AliQAv1::GetQADataFileName()), "RECREATE") ;
69a898e3 924 TH1I * hisRun = new TH1I("hLMR", "List of merged runs", hiRun-loRun+10, loRun, hiRun+10) ;
925 // create the structure into the merged file
4e25ac79 926 for (Int_t iDet = 0; iDet < AliQAv1::kNDET ; iDet++) {
927 TDirectory * detDir = mergedFile.mkdir(AliQAv1::GetDetName(iDet)) ;
928 for (Int_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
69a898e3 929 detDir->cd() ;
4e25ac79 930 TDirectory * taskDir = gDirectory->mkdir(AliQAv1::GetTaskName(taskIndex)) ;
69a898e3 931 for (Int_t es = 0 ; es < AliRecoParam::kNSpecies ; es++) {
932 taskDir->cd() ;
933 TDirectory * esDir = gDirectory->mkdir(AliRecoParam::GetEventSpecieName(es)) ;
934 esDir->cd() ;
4e25ac79 935 gDirectory->mkdir(AliQAv1::GetExpert()) ;
69a898e3 936 }
937 }
938 }
eca4fa66 939 while ( (srun = static_cast<TObjString *> (nextRun())) ) {
69a898e3 940 runNumber = (srun->String()).Atoi() ;
941 hisRun->Fill(runNumber) ;
5379c4a3 942 AliDebug(AliQAv1::GetQADebugLevel(), Form("Merging run number %d", runNumber)) ;
69a898e3 943 // search all QA files for runNumber in the current directory
fec0891b 944 Char_t * fileList[AliQAv1::kNDET] ;
69a898e3 945 index = 0 ;
4e25ac79 946 for (Int_t iDet = 0; iDet < AliQAv1::kNDET ; iDet++) {
fec0891b 947 Char_t * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQAv1::GetDetName(iDet), AliQAv1::GetQADataFileName(), runNumber));
69a898e3 948 if (file)
949 fileList[index++] = file ;
950 }
951 if ( index == 0 ) {
952 AliError("No QA data file found\n") ;
953 return ;
954 }
955 for ( Int_t i = 0 ; i < index ; i++) {
956 TFile * inFile = TFile::Open(fileList[i]) ;
957 TList * listOfKeys =inFile->GetListOfKeys() ;
958 TIter nextkey(listOfKeys) ;
959 TObject * obj1 ;
960 TString dirName("") ;
961 while ( (obj1 = nextkey()) ) {
962 TDirectory * directoryDet = inFile->GetDirectory(obj1->GetName()) ;
963 if ( directoryDet ) {
5379c4a3 964 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s dir = %s", inFile->GetName(), directoryDet->GetName())) ;
69a898e3 965 dirName += Form("%s/", directoryDet->GetName() ) ;
966 directoryDet->cd() ;
967 TList * listOfTasks = directoryDet->GetListOfKeys() ;
968 TIter nextTask(listOfTasks) ;
969 TObject * obj2 ;
970 while ( (obj2 = nextTask()) ) {
971 TDirectory * directoryTask = directoryDet->GetDirectory(obj2->GetName()) ;
972 if ( directoryTask ) {
973 dirName += Form("%s", obj2->GetName()) ;
5379c4a3 974 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s", dirName.Data())) ;
69a898e3 975 directoryTask->cd() ;
976 TList * listOfEventSpecie = directoryTask->GetListOfKeys() ;
977 TIter nextEventSpecie(listOfEventSpecie) ;
978 TObject * obj3 ;
979 while ( (obj3 = nextEventSpecie()) ) {
980 TDirectory * directoryEventSpecie = directoryTask->GetDirectory(obj3->GetName()) ;
981 if ( directoryEventSpecie ) {
982 dirName += Form("/%s/", obj3->GetName()) ;
5379c4a3 983 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s\n", dirName.Data())) ;
69a898e3 984 directoryEventSpecie->cd() ;
985 // histograms are here
986 TDirectory * mergedDirectory = mergedFile.GetDirectory(dirName.Data()) ;
987 TList * listOfData = directoryEventSpecie->GetListOfKeys() ;
988 TIter nextData(listOfData) ;
989 TKey * key ;
eca4fa66 990 while ( (key = static_cast<TKey *>(nextData())) ) {
69a898e3 991 TString className(key->GetClassName()) ;
992 if ( className.Contains("TH") || className.Contains("TProfile") ) {
eca4fa66 993 TH1 * histIn = static_cast<TH1*> (key->ReadObj()) ;
994 TH1 * histOu = static_cast<TH1*> (mergedDirectory->FindObjectAny(histIn->GetName())) ;
65b25288 995 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s %p %p\n", key->GetName(), histIn, histOu)) ;
69a898e3 996 mergedDirectory->cd() ;
997 if ( ! histOu ) {
998 histIn->Write() ;
999 } else {
1000 histOu->Add(histIn) ;
1001 histOu->Write(histOu->GetName(), kOverwrite) ;
1002 }
1003 }
1004 else if ( className.Contains("TDirectoryFile") ) {
1005 TDirectory * dirExpert = directoryEventSpecie->GetDirectory(key->GetName()) ;
1006 dirExpert->cd() ;
1007 TDirectory * mergedDirectoryExpert = mergedDirectory->GetDirectory(dirExpert->GetName()) ;
1008 TList * listOfExpertData = dirExpert->GetListOfKeys() ;
1009 TIter nextExpertData(listOfExpertData) ;
1010 TKey * keykey ;
eca4fa66 1011 while ( (keykey = static_cast<TKey *>(nextExpertData())) ) {
69a898e3 1012 TString classNameExpert(keykey->GetClassName()) ;
1013 if (classNameExpert.Contains("TH")) {
eca4fa66 1014 TH1 * histInExpert = static_cast<TH1*> (keykey->ReadObj()) ;
1015 TH1 * histOuExpert = static_cast<TH1*> (mergedDirectory->FindObjectAny(histInExpert->GetName())) ;
69a898e3 1016 mergedDirectoryExpert->cd() ;
1017 if ( ! histOuExpert ) {
1018 histInExpert->Write() ;
1019 } else {
1020 histOuExpert->Add(histInExpert) ;
1021 histOuExpert->Write(histOuExpert->GetName(), kOverwrite) ;
1022 }
1023 }
1024 }
1025 } else {
1026 AliError(Form("No merge done for this object %s in %s", key->GetName(), dirName.Data())) ;
1027 }
1028 }
1029 dirName.ReplaceAll(Form("/%s/",obj3->GetName()), "") ;
1030 }
1031 }
1032 dirName.ReplaceAll(obj2->GetName(), "") ;
1033 }
1034 }
1035 }
1036 }
1037 inFile->Close() ;
1038 }
1039 }
1040 mergedFile.cd() ;
1041 hisRun->Write() ;
1042 mergedFile.Close() ;
1043 srunList.Delete() ;
1044}
1045
2e331c8b 1046//_____________________________________________________________________________
87da0921 1047Bool_t AliQAManager::MergeData(const Int_t runNumber, const char *fileName) const
2e331c8b 1048{
bededc30 1049 // Merge QA data from all detectors for a given run in one single file
1050
87da0921 1051 TFileMerger merger(kFALSE) ;
1052 TString outFileName = fileName;
1053 if (outFileName.IsNull()) outFileName.Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName());
2e331c8b 1054 merger.OutputFile(outFileName.Data()) ;
bededc30 1055 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
fec0891b 1056 Char_t * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQAv1::GetDetName(iDet), AliQAv1::GetQADataFileName(), runNumber));
bededc30 1057 if (file)
3136db6f 1058 merger.AddFile(file);
1059 delete[] file;
bededc30 1060 }
2e331c8b 1061 merger.Merge() ;
1062 return kTRUE ;
1063}
1064
1065//_____________________________________________________________________________
1066Bool_t AliQAManager::MergeResults(const Int_t runNumber) const
1067{
1068 // Merge the QA result from all the data chunks in a run
bededc30 1069 // to be revised whwn it will be used (see MergeData)
2e331c8b 1070 TString cmd ;
4e25ac79 1071 cmd = Form(".! ls %s*.root > tempo.txt", AliQAv1::GetQADataFileName()) ;
2e331c8b 1072 gROOT->ProcessLine(cmd.Data()) ;
1073 ifstream in("tempo.txt") ;
1074 const Int_t chunkMax = 100 ;
1075 TString fileList[chunkMax] ;
1076
1077 Int_t index = 0 ;
1078 while ( 1 ) {
1079 TString file ;
1080 in >> fileList[index] ;
1081 if ( !in.good() )
1082 break ;
5379c4a3 1083 AliDebug(AliQAv1::GetQADebugLevel(), Form("index = %d file = %s", index, (fileList[index].Data()))) ;
2e331c8b 1084 index++ ;
1085 }
1086
1087 if ( index == 0 ) {
1088 AliError("No QA Result File found") ;
1089 return kFALSE ;
1090 }
1091
1092 TFileMerger merger ;
fc07289e 1093 TString outFileName ;
1094 if (runNumber != -1)
4e25ac79 1095 outFileName = Form("Merged.%s.Result.%d.root",AliQAv1::GetQADataFileName(),runNumber);
fc07289e 1096 else
4e25ac79 1097 outFileName = Form("Merged.%s.Result.root",AliQAv1::GetQADataFileName());
2e331c8b 1098 merger.OutputFile(outFileName.Data()) ;
1099 for (Int_t ifile = 0 ; ifile < index ; ifile++) {
1100 TString file = fileList[ifile] ;
1101 merger.AddFile(file) ;
1102 }
1103 merger.Merge() ;
1104
1105 return kTRUE ;
1106}
1107
1108//_____________________________________________________________________________
1109void AliQAManager::Reset(const Bool_t sameCycle)
1110{
1111 // Reset the default data members
1112
1113 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
4e25ac79 1114 if (IsSelected(AliQAv1::GetDetName(iDet))) {
2e331c8b 1115 AliQADataMaker * qadm = GetQADataMaker(iDet);
5aae1dab 1116 if (qadm)
1117 qadm->Reset();
2e331c8b 1118 }
1119 }
1120 if (fRawReaderDelete) {
1121 delete fRawReader ;
1122 fRawReader = NULL ;
1123 }
1124
1125 fCycleSame = sameCycle ;
1126 fESD = NULL ;
1127 fESDTree = NULL ;
1128 //fFirst = kTRUE ;
1129 fNumberOfEvents = 999999 ;
1130}
1131
9ac91920 1132//_____________________________________________________________________________
1133void AliQAManager::ResetDetectors(AliQAv1::TASKINDEX_t task, AliQAv1::DETECTORINDEX_t det)
1134{
1135 //calls ResetDetector of specified or all detectors
1136 UInt_t iDet = 0 ;
1137 UInt_t iDetMax = fgkNDetectors ;
1138 if ( det != AliQAv1::kNULLDET ) {
1139 iDet = det ;
1140 iDetMax = det+1 ;
1141 }
1142
1143 for (iDet = 0; iDet < iDetMax ; iDet++) {
1144 if (IsSelected(AliQAv1::GetDetName(iDet))) {
1145 AliQADataMaker * qadm = GetQADataMaker(iDet);
1146 qadm->ResetDetector(task);
1147 }
1148 }
1149}
1150
2e331c8b 1151//_____________________________________________________________________________
634696f5 1152AliQAManager * AliQAManager::QAManager(AliQAv1::MODE_t mode, TMap *entryCache, Int_t run)
2e331c8b 1153{
1154 // returns AliQAManager instance (singleton)
1155
1156 if (!fgQAInstance) {
5aae1dab 1157 if ( (mode != AliQAv1::kSIMMODE) && (mode != AliQAv1::kRECMODE) && (mode != AliQAv1::kQAMODE) ) {
1158 AliWarningClass("You must specify kSIMMODE or kRECMODE or kQAMODE") ;
bf76b847 1159 return NULL ;
1160 }
2e331c8b 1161 fgQAInstance = new AliQAManager(mode) ;
1162 if (!entryCache)
1163 fgQAInstance->Init();
1164 else
1165 fgQAInstance->InitFromCache(entryCache,run);
1166 }
1167 return fgQAInstance;
1168}
1169
5e303886 1170//_____________________________________________________________________________
1171AliQAManager * AliQAManager::QAManager(AliQAv1::TASKINDEX_t task)
1172{
1173 // returns AliQAManager instance (singleton)
634696f5 1174 return QAManager(AliQAv1::Mode(task)) ;
5e303886 1175}
1176
2e331c8b 1177//_____________________________________________________________________________
fec0891b 1178TString AliQAManager::Run(const Char_t * detectors, AliRawReader * rawReader, const Bool_t sameCycle)
2e331c8b 1179{
1180 //Runs all the QA data Maker for Raws only
1181
1182 fCycleSame = sameCycle ;
1183 fRawReader = rawReader ;
1184 fDetectors = detectors ;
1185 fDetectorsW = detectors ;
1186
1187 AliCDBManager* man = AliCDBManager::Instance() ;
1188
1189 if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
1190 rawReader->NextEvent() ;
1191 man->SetRun(fRawReader->GetRunNumber()) ;
1192 rawReader->RewindEvents() ;
1193 }
1194
1195 if (!fCycleSame)
4e25ac79 1196 if ( !InitQA(AliQAv1::kRAWS) )
9b4aee57 1197 return "" ;
2e331c8b 1198 fRawReaderDelete = kFALSE ;
1199
4e25ac79 1200 DoIt(AliQAv1::kRAWS) ;
2e331c8b 1201 return fDetectorsW ;
1202}
1203
1204//_____________________________________________________________________________
fec0891b 1205TString AliQAManager::Run(const Char_t * detectors, const Char_t * fileName, const Bool_t sameCycle)
2e331c8b 1206{
1207 //Runs all the QA data Maker for Raws only
1208
1209 fCycleSame = sameCycle ;
1210 fDetectors = detectors ;
1211 fDetectorsW = detectors ;
1212
1213 AliCDBManager* man = AliCDBManager::Instance() ;
1214 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
1215 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
1216 if ( ! rl ) {
1217 AliFatal("galice.root file not found in current directory") ;
1218 } else {
1219 rl->CdGAFile() ;
1220 rl->LoadgAlice() ;
1221 if ( ! rl->GetAliRun() ) {
1222 AliFatal("AliRun not found in galice.root") ;
1223 } else {
1224 rl->LoadHeader() ;
1225 man->SetRun(rl->GetHeader()->GetRun());
1226 }
1227 }
1228 }
1229
1230 if (!fCycleSame)
4e25ac79 1231 if ( !InitQA(AliQAv1::kRAWS, fileName) )
9b4aee57 1232 return "" ;
2e331c8b 1233
4e25ac79 1234 DoIt(AliQAv1::kRAWS) ;
2e331c8b 1235 return fDetectorsW ;
1236}
1237
1238//_____________________________________________________________________________
fec0891b 1239TString AliQAManager::Run(const Char_t * detectors, const AliQAv1::TASKINDEX_t taskIndex, Bool_t const sameCycle, const Char_t * fileName )
2e331c8b 1240{
1241 // Runs all the QA data Maker for every detector
1242
1243 fCycleSame = sameCycle ;
1244 fDetectors = detectors ;
1245 fDetectorsW = detectors ;
1246
1247 AliCDBManager* man = AliCDBManager::Instance() ;
1248 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
1249 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
1250 if ( ! rl ) {
1251 AliFatal("galice.root file not found in current directory") ;
1252 } else {
1253 rl->CdGAFile() ;
1254 rl->LoadgAlice() ;
1255 if ( ! rl->GetAliRun() ) {
5379c4a3 1256 AliDebug(AliQAv1::GetQADebugLevel(), "AliRun not found in galice.root") ;
2e331c8b 1257 } else {
1258 rl->LoadHeader() ;
1259 man->SetRun(rl->GetHeader()->GetRun()) ;
1260 }
1261 }
1262 }
fec0891b 1263 if ( taskIndex == AliQAv1::kNULLTASKINDEX) {
4e25ac79 1264 for (UInt_t task = 0; task < AliQAv1::kNTASKINDEX; task++) {
2e331c8b 1265 if ( fTasks.Contains(Form("%d", task)) ) {
1266 if (!fCycleSame)
4e25ac79 1267 if ( !InitQA(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(task)), fileName) )
9b4aee57 1268 return "" ;
4e25ac79 1269 DoIt(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(task))) ;
2e331c8b 1270 }
1271 }
1272 } else {
1273 if (! fCycleSame )
1274 if ( !InitQA(taskIndex, fileName) )
9b4aee57 1275 return "" ;
2e331c8b 1276 DoIt(taskIndex) ;
1277 }
2e331c8b 1278 return fDetectorsW ;
2e331c8b 1279}
1280
1281//_____________________________________________________________________________
1282void AliQAManager::RunOneEvent(AliRawReader * rawReader)
1283{
1284 //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
1285 if ( ! rawReader )
1286 return ;
0d105c36 1287 if (fTasks.Contains(Form("%d", AliQAv1::kRAWS))){
2e331c8b 1288 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 1289 if (!IsSelected(AliQAv1::GetDetName(iDet)))
2e331c8b 1290 continue;
1291 AliQADataMaker *qadm = GetQADataMaker(iDet);
1292 if (!qadm)
1293 continue;
1294 if ( qadm->IsCycleDone() ) {
1295 qadm->EndOfCycle() ;
1296 }
5e232cd6 1297 qadm->SetEventSpecie(fEventSpecie) ;
13a2370e 1298 if ( qadm->GetRecoParam() )
5e232cd6 1299 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
1300 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
4e25ac79 1301 qadm->Exec(AliQAv1::kRAWS, rawReader) ;
2e331c8b 1302 }
1303 }
1304}
1305
1306//_____________________________________________________________________________
aeb8fc30 1307void AliQAManager::RunOneEvent(AliESDEvent *& esd, AliESDEvent *& hltesd)
2e331c8b 1308{
1309 //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1310
0d105c36 1311 if (fTasks.Contains(Form("%d", AliQAv1::kESDS))) {
2e331c8b 1312 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 1313 if (!IsSelected(AliQAv1::GetDetName(iDet)))
2e331c8b 1314 continue;
1315 AliQADataMaker *qadm = GetQADataMaker(iDet);
1316 if (!qadm)
1317 continue;
5e232cd6 1318 qadm->SetEventSpecie(fEventSpecie) ;
82a30871 1319 if ( qadm->GetRecoParam() )
5e232cd6 1320 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
1321 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
2e331c8b 1322 if ( qadm->IsCycleDone() ) {
1323 qadm->EndOfCycle() ;
1324 }
aeb8fc30 1325 if (iDet == AliQAv1::kHLT) {
1326 TObjArray esdarray;
1327 esdarray.Add(esd);
1328 esdarray.Add(hltesd);
1329 qadm->Exec(AliQAv1::kESDS, &esdarray);
1330 } else {
1331 qadm->Exec(AliQAv1::kESDS, esd) ;
1332 }
2e331c8b 1333 }
1334 }
1335}
1336
1337//_____________________________________________________________________________
1338void AliQAManager::RunOneEventInOneDetector(Int_t det, TTree * tree)
1339{
1340 // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
44ed7a66 1341
1342 TString test(tree->GetName()) ;
0d105c36 1343 if (fTasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
4e25ac79 1344 if (IsSelected(AliQAv1::GetDetName(det))) {
2e331c8b 1345 AliQADataMaker *qadm = GetQADataMaker(det);
1346 if (qadm) {
5e232cd6 1347 qadm->SetEventSpecie(fEventSpecie) ;
31c31326 1348 if ( qadm->GetRecoParam() ) {
5e232cd6 1349 if ( AliRecoParam::Convert(qadm->GetRecoParam()->GetEventSpecie()) != AliRecoParam::kDefault)
31c31326 1350 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
1351 else
5c5ea798 1352 AliError(Form("%d defined by %s is not an event specie", qadm->GetRecoParam()->GetEventSpecie(), qadm->GetName())) ;
31c31326 1353 }
2e331c8b 1354 if ( qadm->IsCycleDone() ) {
1355 qadm->EndOfCycle() ;
1356 }
44ed7a66 1357 if (test.Contains("TreeD")) {
1358 qadm->Exec(AliQAv1::kDIGITSR, tree) ;
44ed7a66 1359 } else if (test.Contains("TreeR")) {
1360 qadm->Exec(AliQAv1::kRECPOINTS, tree) ;
44ed7a66 1361 }
2e331c8b 1362 }
1363 }
1364 }
1365}
1366
1367//_____________________________________________________________________________
fec0891b 1368Bool_t AliQAManager::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const Char_t * year, const Char_t * detectors) const
2e331c8b 1369{
1370 // take the locasl QA data merge into a single file and save in OCDB
1371 Bool_t rv = kTRUE ;
4e25ac79 1372 TString tmp(AliQAv1::GetQARefStorage()) ;
2e331c8b 1373 if ( tmp.IsNull() ) {
4e25ac79 1374 AliError("No storage defined, use AliQAv1::SetQARefStorage") ;
2e331c8b 1375 return kFALSE ;
1376 }
4e25ac79 1377 if ( !(tmp.Contains(AliQAv1::GetLabLocalOCDB()) || tmp.Contains(AliQAv1::GetLabAliEnOCDB())) ) {
1378 AliError(Form("%s is a wrong storage, use %s or %s", AliQAv1::GetQARefStorage(), AliQAv1::GetLabLocalOCDB().Data(), AliQAv1::GetLabAliEnOCDB().Data())) ;
2e331c8b 1379 return kFALSE ;
1380 }
1381 TString sdet(detectors) ;
1382 sdet.ToUpper() ;
1383 TFile * inputFile ;
1384 if ( sdet.Contains("ALL") ) {
1385 rv = Merge(runNumber) ;
1386 if ( ! rv )
1387 return kFALSE ;
4e25ac79 1388 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQAv1::GetQADataFileName(), runNumber)) ;
2e331c8b 1389 inputFile = TFile::Open(inputFileName.Data()) ;
1390 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ;
1391 } else {
4e25ac79 1392 for (Int_t index = 0; index < AliQAv1::kNDET; index++) {
1393 if (sdet.Contains(AliQAv1::GetDetName(index))) {
1394 TString inputFileName(Form("%s.%s.%d.root", AliQAv1::GetDetName(index), AliQAv1::GetQADataFileName(), runNumber)) ;
2e331c8b 1395 inputFile = TFile::Open(inputFileName.Data()) ;
1396 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ;
1397 }
1398 }
1399 }
1400 return rv ;
1401}
1402
1403//_____________________________________________________________________________
fec0891b 1404Bool_t AliQAManager::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const Char_t * year, AliRecoParam::EventSpecie_t es) const
2e331c8b 1405{
1406 // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1407 Bool_t rv = kTRUE ;
5379c4a3 1408 AliDebug(AliQAv1::GetQADebugLevel(), Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQAv1::GetQARefStorage())) ;
2e331c8b 1409 if ( ! IsDefaultStorageSet() ) {
4e25ac79 1410 TString tmp( AliQAv1::GetQARefStorage() ) ;
1411 if ( tmp.Contains(AliQAv1::GetLabLocalOCDB()) )
1412 Instance()->SetDefaultStorage(AliQAv1::GetQARefStorage()) ;
2e331c8b 1413 else {
4e25ac79 1414 TString tmp1(AliQAv1::GetQARefDefaultStorage()) ;
2e331c8b 1415 tmp1.Append(year) ;
1416 tmp1.Append("?user=alidaq") ;
1417 Instance()->SetDefaultStorage(tmp1.Data()) ;
1418 }
1419 }
4e25ac79 1420 Instance()->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
2e331c8b 1421 if(GetRun() < 0)
1422 Instance()->SetRun(runNumber);
1423
1424 AliCDBMetaData mdr ;
1425 mdr.SetResponsible("yves schutz");
1426
4e25ac79 1427 for ( Int_t detIndex = 0 ; detIndex < AliQAv1::kNDET ; detIndex++) {
1428 TDirectory * detDir = inputFile->GetDirectory(AliQAv1::GetDetName(detIndex)) ;
2e331c8b 1429 if ( detDir ) {
5379c4a3 1430 AliDebug(AliQAv1::GetQADebugLevel(), Form("Entering %s", detDir->GetName())) ;
4e25ac79 1431 AliQAv1::SetQARefDataDirName(es) ;
1432 TString detOCDBDir(Form("%s/%s/%s", AliQAv1::GetDetName(detIndex), AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName())) ;
2e331c8b 1433 AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity()) ;
1434 TList * listDetQAD = new TList() ;
4e25ac79 1435 TString listName(Form("%s QA data Reference", AliQAv1::GetDetName(detIndex))) ;
1436 mdr.SetComment(Form("%s QA stuff", AliQAv1::GetDetName(detIndex)));
2e331c8b 1437 listDetQAD->SetName(listName) ;
1438 TList * taskList = detDir->GetListOfKeys() ;
1439 TIter nextTask(taskList) ;
1440 TKey * taskKey ;
eca4fa66 1441 while ( (taskKey = static_cast<TKey*>(nextTask())) ) {
2e331c8b 1442 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ;
1443 TDirectory * esDir = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ;
5379c4a3 1444 AliDebug(AliQAv1::GetQADebugLevel(), Form("Saving %s", esDir->GetName())) ;
2e331c8b 1445 TObjArray * listTaskQAD = new TObjArray(100) ;
1446 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1447 listDetQAD->Add(listTaskQAD) ;
1448 TList * histList = esDir->GetListOfKeys() ;
1449 TIter nextHist(histList) ;
1450 TKey * histKey ;
eca4fa66 1451 while ( (histKey = static_cast<TKey*>(nextHist())) ) {
2e331c8b 1452 TObject * odata = esDir->Get(histKey->GetName()) ;
1453 if ( !odata ) {
1454 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1455 } else {
4e25ac79 1456 if ( AliQAv1::GetExpert() == histKey->GetName() ) {
2e331c8b 1457 TDirectory * expertDir = esDir->GetDirectory(histKey->GetName()) ;
1458 TList * expertHistList = expertDir->GetListOfKeys() ;
1459 TIter nextExpertHist(expertHistList) ;
1460 TKey * expertHistKey ;
eca4fa66 1461 while ( (expertHistKey = static_cast<TKey*>(nextExpertHist())) ) {
2e331c8b 1462 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ;
1463 if ( !expertOdata ) {
1464 AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1465 } else {
5379c4a3 1466 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", expertHistKey->GetName())) ;
2e331c8b 1467 if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
5379c4a3 1468 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", expertHistKey->GetName())) ;
2e331c8b 1469 TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ;
1470 listTaskQAD->Add(hExpertdata) ;
1471 }
1472 }
1473 }
1474 }
5379c4a3 1475 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", histKey->GetName())) ;
2e331c8b 1476 if ( odata->IsA()->InheritsFrom("TH1") ) {
5379c4a3 1477 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", histKey->GetName())) ;
2e331c8b 1478 TH1 * hdata = static_cast<TH1*>(odata) ;
1479 listTaskQAD->Add(hdata) ;
1480 }
1481 }
1482 }
1483 }
1484 Instance()->Put(listDetQAD, idr, &mdr) ;
1485 }
1486 }
1487 return rv ;
1488}
1489
7c9ff472 1490//_____________________________________________________________________________
1491
1492void AliQAManager::SetCheckerExternParam(AliQAv1::DETECTORINDEX_t detIndex, TList * parameterList)
1493{
1494 // set the external parameters list for the detector checkers
1495 AliQACheckerBase * qac = AliQAChecker::Instance()->GetDetQAChecker(detIndex) ;
1496 qac->SetExternParamlist(parameterList) ;
1497 qac->PrintExternParam() ;
1498}
1499
2e331c8b 1500//_____________________________________________________________________________
1501void AliQAManager::SetEventSpecie(AliRecoParam::EventSpecie_t es)
1502{
4e25ac79 1503 // set the current event specie and inform AliQAv1 that this event specie has been encountered
bc8761a0 1504 fEventSpecie = es ;
4e25ac79 1505 AliQAv1::Instance()->SetEventSpecie(es) ;
2e331c8b 1506}
1507
1508//_____________________________________________________________________________
1509void AliQAManager::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par)
1510{
1511 // Set custom reconstruction parameters for a given detector
1512 // Single set of parameters for all the events
1513 GetQADataMaker(det)->SetRecoParam(par) ;
1514}
1515
75373542 1516//_____________________________________________________________________________
1517void AliQAManager::SetWriteExpert()
1518{
1519 // enable the writing of QA expert data
1520 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
4e25ac79 1521 if (IsSelected(AliQAv1::GetDetName(iDet)))
75373542 1522 fQAWriteExpert[iDet] = kTRUE ;
1523 }
1524}
7cade814 1525
1526//_____________________________________________________________________________
1527void AliQAManager::Destroy() {
bf76b847 1528 // delete AliQAManager instance and
1529 // all associated objects
7cade814 1530
1531 if (fgQAInstance) {
eca4fa66 1532 delete fgQAInstance ;
1533 fgQAInstance = NULL ;
7cade814 1534 }
1535}
1536
bf76b847 1537//_____________________________________________________________________________
1538void AliQAManager::ShowQA() {
1539 // Show the result of the QA checking
1540 // for all detectors
1541 for ( Int_t detIndex = 0 ; detIndex < AliQAv1::kNDET ; detIndex++)
1542 AliQAv1::Instance(AliQAv1::GetDetIndex(AliQAv1::GetDetName(detIndex)))->Show() ;
1543}