]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliQA.cxx
Set charge-array to zero if old class version was read.
[u/mrichter/AliRoot.git] / STEER / AliQA.cxx
CommitLineData
8661738e 1
421ab0fb 2/**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
7 * *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
16/* $Id$ */
17
18//////////////////////////////////////////////////////////////////////////////
19//
20// Quality Assurance Object//_________________________________________________________________________
21// Quality Assurance object. The QA status is held in one word per detector,
22// each bit corresponds to a different status.
12581b04 23// bit 0-3 : QA raised during simulation (RAW)
24// bit 4-7 : QA raised during simulation (SIM)
25// bit 8-11 : QA raised during reconstruction (REC)
26// bit 12-15 : QA raised during ESD checking (ESD)
27// bit 16-19 : QA raised during analysis (ANA)
421ab0fb 28// Each of the 4 bits corresponds to a severity level of increasing importance
29// from lower to higher bit (INFO, WARNING, ERROR, FATAL)
30//
31//*-- Yves Schutz CERN, July 2007
32//////////////////////////////////////////////////////////////////////////////
33
34
b09247a2 35#include <cstdlib>
421ab0fb 36// --- ROOT system ---
37#include <TFile.h>
38#include <TSystem.h>
0d13dd63 39#include <TROOT.h>
421ab0fb 40
41// --- Standard library ---
42
43// --- AliRoot header files ---
44#include "AliLog.h"
2e42b4d4 45#include "AliQA.h"
421ab0fb 46
47
2e42b4d4 48ClassImp(AliQA)
f73f556a 49AliQA * AliQA::fgQA = 0x0 ;
50TFile * AliQA::fgQADataFile = 0x0 ;
51TString AliQA::fgQADataFileName = "QA" ; // will transform into Det.QA.run.cycle.root
52TFile * AliQA::fgQARefFile = 0x0 ;
53TString AliQA::fgQARefDirName = "" ;
54TString AliQA::fgQARefFileName = "QA.root" ;
55TFile * AliQA::fgQAResultFile = 0x0 ;
56TString AliQA::fgQAResultDirName = "" ;
57TString AliQA::fgQAResultFileName = "QA.root" ;
58TString AliQA::fgDetNames[] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD",
a2b64fbd 59 "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT", "Global", "CORR"} ;
4fb24ca4 60TString AliQA::fgGRPPath = "GRP/GRP/Data" ;
61TString AliQA::fgRTNames[] = {"UNKNOWN", "AUTO_TEST", "CALIBRATION", "CALIBRATION_PULSER", "CHANNEL_DELAY_TUNING", "COSMIC",
7ce961eb 62 "COSMICS", "DAQ_FO_UNIF_SCAN", "DAQ_GEN_DAC_SCAN", "DAQ_MEAN_TH_SCAN", "DAQ_MIN_TH_SCAN",
63 "DAQ_NOISY_PIX_SCAN", "DAQ_PIX_DELAY_SCAN", "DAQ_UNIFORMITY_SCAN", "DCS_FO_UNIF_SCAN",
64 "DCS_MEAN_TH_SCAN", "DCS_MIN_TH_SCAN", "DCS_PIX_DELAY_SCAN", "DCS_UNIFORMITY_SCAN",
65 "DDL_TEST", "GAIN", "PEDESTAL", "INJECTOR", "LASER", "MONTECARLO", "NOISE", "NOISY_PIX_SCAN",
66 "PHYSICS", "PULSER", "STANDALONE", "STANDALONE_BC", "STANDALONE_CENTRAL", "STANDALONE_COSMIC",
67 "STANDALONE_EMD", "STANDALONE_LASER", "STANDALONE_MB", "STANDALONE_PEDESTAL",
68 "STANDALONE_SEMICENTRAL", "STANDALONE_PULSER" } ;
f73f556a 69TString AliQA::fgTaskNames[] = {"Raws", "Hits", "SDigits", "Digits", "RecPoints", "TrackSegments", "RecParticles", "ESDs"} ;
70const TString AliQA::fkgLabLocalFile = "file://" ;
71const TString AliQA::fkgLabLocalOCDB = "local://" ;
72const TString AliQA::fkgLabAliEnOCDB = "alien://" ;
73const TString AliQA::fkgRefFileName = "QA.root" ;
96d67a8d 74const TString AliQA::fkgQAName = "QA" ;
a2b64fbd 75const TString AliQA::fkgQACorrNtName = "CorrQA" ;
f73f556a 76const TString AliQA::fkgRefOCDBDirName = "Ref" ;
4fb24ca4 77TString AliQA::fkgRefDataDirName = "" ;
13c8b505 78const TString AliQA::fkgQARefOCDBDefault = "alien://folder=/alice/QA/20" ;
421ab0fb 79//____________________________________________________________________________
2e42b4d4 80AliQA::AliQA() :
421ab0fb 81 TNamed("", ""),
9a223722 82 fNdet(kNDET),
83 fQA(new ULong_t[fNdet]),
421ab0fb 84 fDet(kNULLDET),
85 fTask(kNULLTASK)
4edbc5bc 86
421ab0fb 87{
88 // default constructor
89 // beware singleton: not to be used
4ecde5fc 90
9a223722 91 for (Int_t index = 0 ; index < fNdet ; index++)
92 fQA[index] = 0 ;
421ab0fb 93}
94
95//____________________________________________________________________________
2e42b4d4 96AliQA::AliQA(const AliQA& qa) :
421ab0fb 97 TNamed(qa),
9a223722 98 fNdet(qa.fNdet),
421ab0fb 99 fQA(qa.fQA),
100 fDet(qa.fDet),
101 fTask(qa.fTask)
102{
103 // cpy ctor
104}
105
106//_____________________________________________________________________________
2e42b4d4 107AliQA& AliQA::operator = (const AliQA& qa)
421ab0fb 108{
109// assignment operator
110
2e42b4d4 111 this->~AliQA();
112 new(this) AliQA(qa);
421ab0fb 113 return *this;
114}
115
116//_______________________________________________________________
96d67a8d 117AliQA::AliQA(const DETECTORINDEX_t det) :
9a223722 118 TNamed("QA", "Quality Assurance status"),
119 fNdet(kNDET),
120 fQA(new ULong_t[fNdet]),
a5fa6165 121 fDet(det),
4edbc5bc 122 fTask(kNULLTASK)
a5fa6165 123{
124 // constructor to be used
125 if (! CheckRange(det) ) {
126 fDet = kNULLDET ;
127 return ;
a4976ef3 128 }
a5fa6165 129 Int_t index ;
9a223722 130 for (index = 0; index < fNdet; index++)
a5fa6165 131 fQA[index] = 0 ;
132}
133
134//_______________________________________________________________
96d67a8d 135AliQA::AliQA(const ALITASK_t tsk) :
421ab0fb 136 TNamed("QA", "Quality Assurance status"),
9a223722 137 fNdet(kNDET),
138 fQA(new ULong_t[fNdet]),
421ab0fb 139 fDet(kNULLDET),
140 fTask(tsk)
141{
142 // constructor to be used in the AliRoot module (SIM, REC, ESD or ANA)
143 if (! CheckRange(tsk) ) {
144 fTask = kNULLTASK ;
145 return ;
a4976ef3 146 }
a5fa6165 147 Int_t index ;
9a223722 148 for (index = 0; index < fNdet; index++)
a5fa6165 149 fQA[index] = 0 ;
421ab0fb 150}
151
152//____________________________________________________________________________
2e42b4d4 153AliQA::~AliQA()
421ab0fb 154{
155 // dtor
156 delete[] fQA ;
157}
158
46ca5304 159//_______________________________________________________________
160void AliQA::Close()
161{
162 // close the open files
163 if (fgQADataFile)
164 if (fgQADataFile->IsOpen())
165 fgQADataFile->Close() ;
166 if (fgQAResultFile)
167 if (fgQAResultFile->IsOpen())
168 fgQAResultFile->Close() ;
169 if (fgQARefFile)
170 if (fgQARefFile->IsOpen())
171 fgQARefFile->Close() ;
172}
173
4ecde5fc 174//_______________________________________________________________
175const Bool_t AliQA::CheckFatal() const
176{
177 // check if any FATAL status is set
178 Bool_t rv = kFALSE ;
179 Int_t index ;
180 for (index = 0; index < kNDET ; index++)
96d67a8d 181 rv = rv || IsSet(DETECTORINDEX_t(index), fTask, kFATAL) ;
4ecde5fc 182 return rv ;
183}
184
421ab0fb 185//_______________________________________________________________
96d67a8d 186const Bool_t AliQA::CheckRange(DETECTORINDEX_t det) const
421ab0fb 187{
a4976ef3 188 // check if detector is in given detector range: 0-kNDET
421ab0fb 189
a5fa6165 190 Bool_t rv = ( det < 0 || det > kNDET ) ? kFALSE : kTRUE ;
421ab0fb 191 if (!rv)
192 AliFatal(Form("Detector index %d is out of range: 0 <= index <= %d", det, kNDET)) ;
193 return rv ;
194}
195
196//_______________________________________________________________
96d67a8d 197const Bool_t AliQA::CheckRange(ALITASK_t task) const
421ab0fb 198{
199 // check if task is given taskk range: 0:kNTASK
6c18591a 200 Bool_t rv = ( task < kRAW || task > kNTASK ) ? kFALSE : kTRUE ;
421ab0fb 201 if (!rv)
202 AliFatal(Form("Module index %d is out of range: 0 <= index <= %d", task, kNTASK)) ;
203 return rv ;
204}
205
206//_______________________________________________________________
96d67a8d 207const Bool_t AliQA::CheckRange(QABIT_t bit) const
421ab0fb 208{
209 // check if bit is in given bit range: 0-kNBit
210
211 Bool_t rv = ( bit < 0 || bit > kNBIT ) ? kFALSE : kTRUE ;
212 if (!rv)
213 AliFatal(Form("Status bit %d is out of range: 0 <= bit <= %d", bit, kNBIT)) ;
214 return rv ;
215}
216
421ab0fb 217
4ecde5fc 218
a5fa6165 219//_______________________________________________________________
96d67a8d 220const char * AliQA::GetAliTaskName(ALITASK_t tsk)
421ab0fb 221{
0b96c27c 222 // returns the char name corresponding to module index
223 TString tskName ;
224 switch (tsk) {
225 case kNULLTASK:
226 break ;
227 case kRAW:
228 tskName = "RAW" ;
229 break ;
230 case kSIM:
231 tskName = "SIM" ;
232 break ;
233 case kREC:
234 tskName = "REC" ;
235 break ;
236 case kESD:
237 tskName = "ESD" ;
238 break ;
239 case kANA:
240 tskName = "ANA" ;
241 break ;
242 default:
243 tsk = kNULLTASK ;
244 break ;
245 }
246 return tskName.Data() ;
247}
0acdb2bb 248
0b96c27c 249//_______________________________________________________________
250const char * AliQA::GetBitName(QABIT_t bit) const
251{
252 // returns the char name corresponding to bit
253 TString bitName ;
254 switch (bit) {
255 case kNULLBit:
256 break ;
257 case kINFO:
258 bitName = "INFO" ;
259 break ;
260 case kWARNING:
261 bitName = "WARNING" ;
262 break ;
263 case kERROR:
264 bitName = "ERROR" ;
265 break ;
266 case kFATAL:
267 bitName = "FATAL" ;
268 break ;
269 default:
270 bit = kNULLBit ;
271 break ;
272 }
273 return bitName.Data() ;
421ab0fb 274}
275
46ca5304 276//_______________________________________________________________
7c002d48 277const AliQA::DETECTORINDEX_t AliQA::GetDetIndex(const char * name)
96d67a8d 278{
279 // returns the detector index corresponding to a given name
280 TString sname(name) ;
281 DETECTORINDEX_t rv = kNULLDET ;
282 for (Int_t det = 0; det < kNDET ; det++) {
283 if ( GetDetName(det) == sname ) {
284 rv = DETECTORINDEX_t(det) ;
285 break ;
286 }
287 }
288 return rv ;
289}
290
7c002d48 291//_______________________________________________________________
292const char * AliQA::GetDetName(Int_t det)
293{
294 // returns the detector name corresponding to a given index (needed in a loop)
295
296 if ( det >= 0 && det < kNDET)
297 return (fgDetNames[det]).Data() ;
298 else
299 return NULL ;
300}
301
46ca5304 302//_______________________________________________________________
303TFile * AliQA::GetQADataFile(const char * name, const Int_t run, const Int_t cycle)
304{
305 // opens the file to store the detectors Quality Assurance Data Maker results
96d67a8d 306 const char * temp = Form("%s.%s.%d.%d.root", name, fgQADataFileName.Data(), run, cycle) ;
307 TString opt ;
308 if (! fgQADataFile ) {
309 if (gSystem->AccessPathName(temp))
310 opt = "NEW" ;
311 else
312 opt = "UPDATE" ;
313 fgQADataFile = TFile::Open(temp, opt.Data()) ;
314 } else {
315 if ( strcmp(temp, fgQADataFile->GetName()) != 0 ) {
316 fgQADataFile = dynamic_cast<TFile *>(gROOT->FindObject(temp)) ;
317 if ( !fgQADataFile ) {
318 if (gSystem->AccessPathName(temp))
319 opt = "NEW" ;
320 else
321 opt = "UPDATE" ;
322 fgQADataFile = TFile::Open(temp, opt.Data()) ;
323 }
324 }
46ca5304 325 }
96d67a8d 326 return fgQADataFile ;
46ca5304 327}
328
329//_____________________________________________________________________________
330TFile * AliQA::GetQADataFile(const char * fileName)
331{
332 // Open if necessary the Data file and return its pointer
333
334 if (!fgQADataFile)
335 if (!fileName)
336 fileName = AliQA::GetQADataFileName() ;
337 if (!gSystem->AccessPathName(fileName)) {
338 fgQADataFile = TFile::Open(fileName) ;
339 } else {
340 printf("File %s not found", fileName) ;
341 exit(1) ;
342 }
343 return fgQADataFile ;
344}
345
4ecde5fc 346//_______________________________________________________________
347TFile * AliQA::GetQAResultFile()
348{
349 // opens the file to store the Quality Assurance Data Checker results
2ba0b5f5 350 if (fgQAResultFile)
351 fgQAResultFile->Close() ;
352 fgQAResultFile = 0x0 ;
353// if (!fgQAResultFile) {
1c0190ec 354 TString dirName(fgQAResultDirName) ;
4edbc5bc 355 if ( dirName.Contains(fkgLabLocalFile))
356 dirName.ReplaceAll(fkgLabLocalFile, "") ;
1c0190ec 357 TString fileName(dirName + fgQAResultFileName) ;
46ca5304 358 TString opt("") ;
359 if ( !gSystem->AccessPathName(fileName) )
360 opt = "UPDATE" ;
5bd22145 361 else {
1c0190ec 362 if ( gSystem->AccessPathName(dirName) )
363 gSystem->mkdir(dirName) ;
46ca5304 364 opt = "NEW" ;
5bd22145 365 }
46ca5304 366 fgQAResultFile = TFile::Open(fileName, opt) ;
2ba0b5f5 367// }
46ca5304 368
369 return fgQAResultFile ;
421ab0fb 370}
371
51757634 372//_______________________________________________________________
373const TString AliQA::GetRunTypeName(RUNTYPE_t rt)
374{
375 TString rv("Invalid Run Type") ;
376 if ( rt == kNULLTYPE ) {
377 rv = "Known RUN_TYPE are: \n" ;
378 for (Int_t index = 0 ; index < kNTYPE; index++) {
379 rv += Form("%2d -- %s\n", index, fgRTNames[index].Data()) ;
380 }
381 printf("%s", rv.Data()) ;
382 return "" ;
383 }
384 else {
385 if ( rt > kNULLTYPE && rt < kNTYPE )
386 rv = fgRTNames[rt] ;
387 }
388 return rv ;
389}
390
940d8e5f 391//_______________________________________________________________
392const AliQA::TASKINDEX_t AliQA::GetTaskIndex(const char * name)
393{
394 // returns the detector index corresponding to a given name
395 TString sname(name) ;
396 TASKINDEX_t rv ;
75f6ebf0 397 for (Int_t tsk = 0; tsk < kNTASKINDEX ; tsk++) {
940d8e5f 398 if ( GetTaskName(tsk) == sname ) {
399 rv = TASKINDEX_t(tsk) ;
400 break ;
401 }
402 }
403 return rv ;
404}
405
421ab0fb 406//_______________________________________________________________
96d67a8d 407const Bool_t AliQA::IsSet(DETECTORINDEX_t det, ALITASK_t tsk, QABIT_t bit) const
421ab0fb 408{
409 // Checks is the requested bit is set
384c0618 410
421ab0fb 411 CheckRange(det) ;
412 CheckRange(tsk) ;
413 CheckRange(bit) ;
384c0618 414
421ab0fb 415 ULong_t offset = Offset(tsk) ;
416 ULong_t status = GetStatus(det) ;
417 offset+= bit ;
418 status = (status & 1 << offset) != 0 ;
419 return status ;
420}
421
384c0618 422//_______________________________________________________________
423const Bool_t AliQA::IsSetAny(DETECTORINDEX_t det, ALITASK_t tsk) const
424{
425 // Checks is the requested bit is set
426
427 CheckRange(det) ;
428 CheckRange(tsk) ;
429
430 ULong_t offset = Offset(tsk) ;
431 ULong_t status = GetStatus(det) ;
432 UShort_t st = 0 ;
433 for ( Int_t bit = 0 ; bit < kNBIT ; bit++) {
434 offset+= bit ;
435 st += (status & 1 << offset) != 0 ;
436 }
437 if ( st == 0 )
438 return kFALSE ;
439 else
440 return kTRUE ;
441}
442//_______________________________________________________________
443const Bool_t AliQA::IsSetAny(DETECTORINDEX_t det) const
444{
445 // Checks is the requested bit is set
446
447 CheckRange(det) ;
448
449 ULong_t status = GetStatus(det) ;
450 UShort_t st = 0 ;
451 for ( Int_t tsk = 0 ; tsk < kNTASK ; tsk++) {
452 ULong_t offset = Offset(ALITASK_t(tsk)) ;
453 for ( Int_t bit = 0 ; bit < kNBIT ; bit++) {
454 offset+= bit ;
455 st += (status & 1 << offset) != 0 ;
456 }
457 }
458 if ( st == 0 )
459 return kFALSE ;
460 else
461 return kTRUE ;
462}
463
421ab0fb 464//_______________________________________________________________
2e42b4d4 465AliQA * AliQA::Instance()
421ab0fb 466{
467 // Get an instance of the singleton.
468 // Object must have been instantiated with Instance(ALITASK) first
469
470 return fgQA ;
471}
472
473//_______________________________________________________________
96d67a8d 474AliQA * AliQA::Instance(const DETECTORINDEX_t det)
421ab0fb 475{
476 // Get an instance of the singleton. The only authorized way to call the ctor
477
9a223722 478 if ( ! fgQA) {
46ca5304 479 TFile * f = GetQAResultFile() ;
9a223722 480 fgQA = dynamic_cast<AliQA *>(f->Get("QA")) ;
481 if ( ! fgQA )
482 fgQA = new AliQA(det) ;
483 }
421ab0fb 484 fgQA->Set(det) ;
485 return fgQA ;
486}
487
488//_______________________________________________________________
96d67a8d 489AliQA * AliQA::Instance(const ALITASK_t tsk)
421ab0fb 490{
491 // get an instance of the singleton.
492
493 if ( ! fgQA)
494 switch (tsk) {
495 case kNULLTASK:
496 break ;
6c18591a 497 case kRAW:
2e42b4d4 498 fgQA = new AliQA(tsk) ;
6c18591a 499 break ;
500 case kSIM:
2e42b4d4 501 fgQA = new AliQA(tsk) ;
421ab0fb 502 break ;
503 case kREC:
504 printf("fgQA = gAlice->GetQA()") ;
505 break ;
506 case kESD:
507 printf("fgQA = dynamic_cast<AliQA *> (esdFile->Get(\"QA\")") ;
508 break ;
509 case kANA:
510 printf("fgQA = dynamic_cast<AliQA *> (esdFile->Get(\"QA\")") ;
511 break ;
512 case kNTASK:
513 break ;
514 }
515 if (fgQA)
516 fgQA->Set(tsk) ;
517 return fgQA ;
518}
519
520//_______________________________________________________________
96d67a8d 521AliQA * AliQA::Instance(const TASKINDEX_t tsk)
522{
523 // get an instance of the singleton.
524
525 ALITASK_t index = kNULLTASK ;
526
527 if ( tsk == kRAWS )
528 index = kRAW ;
529 else if (tsk < kDIGITS)
530 index = kSIM ;
531 else if (tsk < kRECPARTICLES)
532 index = kREC ;
533 else if (tsk == kESDS)
534 index = kESD ;
535
536 return Instance(index) ;
537}
538
0acdb2bb 539//_______________________________________________________________
540void AliQA::Merge(TCollection * list) {
541 // Merge the QA resuls in the list into this single AliQA object
542
543 for (Int_t det = 0 ; det < kNDET ; det++) {
544 Set(DETECTORINDEX_t(det)) ;
545 for (Int_t task = 0 ; task < kNTASK ; task++) {
546 Set(ALITASK_t(task)) ;
547 for (Int_t bit = 0 ; bit < kNBIT ; bit++) {
548 TIter next(list) ;
549 AliQA * qa ;
550 while ( (qa = (AliQA*)next() ) ) {
551 qa->IsSet(DETECTORINDEX_t(det), ALITASK_t(task), QABIT_t(bit)) ;
552 Set(QABIT_t(bit)) ;
553 } // qa list
554 } // bit
555 } // task
556 } // detector
557}
558
96d67a8d 559//_______________________________________________________________
560const ULong_t AliQA::Offset(ALITASK_t tsk) const
421ab0fb 561{
562 // Calculates the bit offset for a given module (SIM, REC, ESD, ANA)
563
564 CheckRange(tsk) ;
565
566 ULong_t offset = 0 ;
567 switch (tsk) {
568 case kNULLTASK:
569 break ;
6c18591a 570 case kRAW:
421ab0fb 571 offset+= 0 ;
572 break ;
6c18591a 573 case kSIM:
421ab0fb 574 offset+= 4 ;
575 break ;
6c18591a 576 case kREC:
421ab0fb 577 offset+= 8 ;
578 break ;
6c18591a 579 case kESD:
421ab0fb 580 offset+= 12 ;
581 break ;
6c18591a 582 case kANA:
583 offset+= 16 ;
584 break ;
421ab0fb 585 case kNTASK:
586 break ;
587 }
588
589 return offset ;
590}
591
592//_______________________________________________________________
96d67a8d 593void AliQA::Set(QABIT_t bit)
421ab0fb 594{
595 // Set the status bit of the current detector in the current module
596
597 SetStatusBit(fDet, fTask, bit) ;
598}
599
4ecde5fc 600//_____________________________________________________________________________
4edbc5bc 601void AliQA::SetQARefStorage(const char * name)
4ecde5fc 602{
72c25501 603 // Set the root directory where the QA reference data are stored
604
605 fgQARefDirName = name ;
606 if ( fgQARefDirName.Contains(fkgLabLocalFile) )
607 fgQARefFileName = fkgRefFileName ;
608 else if ( fgQARefDirName.Contains(fkgLabLocalOCDB) )
96d67a8d 609 fgQARefFileName = fkgQAName ;
72c25501 610 else if ( fgQARefDirName.Contains(fkgLabAliEnOCDB) )
96d67a8d 611 fgQARefFileName = fkgQAName ;
4ecde5fc 612
4edbc5bc 613 else {
614 printf("ERROR: %s is an invalid storage definition\n", name) ;
615 fgQARefDirName = "" ;
616 fgQARefFileName = "" ;
617 }
72c25501 618 TString tmp(fgQARefDirName) ; // + fgQARefFileName) ;
4edbc5bc 619 printf("AliQA::SetQARefDir: QA references are in %s\n", tmp.Data() ) ;
4ecde5fc 620}
621
4fb24ca4 622//_____________________________________________________________________________
623void AliQA::SetQARefDataDirName(const char * name)
624{
625 // Set the lower level directory name where reference data are found
626 TString test(name) ;
7ce961eb 627 RUNTYPE_t rt = kNULLTYPE ;
628 for (Int_t index = 0; index < kNTYPE; index++) {
629 if (test.CompareTo(fgRTNames[index]) == 0) {
4fb24ca4 630 rt = (RUNTYPE_t) index ;
631 break ;
7ce961eb 632 }
633 }
634
635 if (rt == kNULLTYPE) {
4fb24ca4 636 printf("AliQA::SetQARefDataDirName: %s is an unknown RUN TYPE name\n", name) ;
637 return ;
7ce961eb 638 }
639
640 SetQARefDataDirName(rt) ;
4fb24ca4 641}
642
4ecde5fc 643//_____________________________________________________________________________
644void AliQA::SetQAResultDirName(const char * name)
645{
646 // Set the root directory where to store the QA status object
647
648 fgQAResultDirName.Prepend(name) ;
649 printf("AliQA::SetQAResultDirName: QA results are in %s\n", fgQAResultDirName.Data()) ;
4edbc5bc 650 if ( fgQAResultDirName.Contains(fkgLabLocalFile))
651 fgQAResultDirName.ReplaceAll(fkgLabLocalFile, "") ;
4ecde5fc 652 fgQAResultFileName.Prepend(fgQAResultDirName) ;
653}
654
421ab0fb 655//_______________________________________________________________
96d67a8d 656void AliQA::SetStatusBit(DETECTORINDEX_t det, ALITASK_t tsk, QABIT_t bit)
421ab0fb 657{
658 // Set the status bit for a given detector and a given task
659
660 CheckRange(det) ;
661 CheckRange(tsk) ;
662 CheckRange(bit) ;
663
664 ULong_t offset = Offset(tsk) ;
665 ULong_t status = GetStatus(det) ;
666 offset+= bit ;
667 status = status | 1 << offset ;
668 SetStatus(det, status) ;
669}
670
671//_______________________________________________________________
2e42b4d4 672void AliQA::ShowAll() const
421ab0fb 673{
674 // dispplay the QA status word
675 Int_t index ;
15df7cc8 676 for (index = 0 ; index < kNDET ; index++) {
677 for (Int_t tsk = kRAW ; tsk < kNTASK ; tsk++) {
678 ShowStatus(DETECTORINDEX_t(index), ALITASK_t(tsk)) ;
679 }
680 }
421ab0fb 681}
682
683//_______________________________________________________________
15df7cc8 684void AliQA::ShowStatus(DETECTORINDEX_t det, ALITASK_t tsk) const
421ab0fb 685{
0b96c27c 686 // Prints the full QA status of a given detector
687 CheckRange(det) ;
688 ULong_t status = GetStatus(det) ;
689 ULong_t tskStatus[kNTASK] ;
690 tskStatus[kRAW] = status & 0x0000f ;
691 tskStatus[kSIM] = status & 0x000f0 ;
692 tskStatus[kREC] = status & 0x00f00 ;
693 tskStatus[kESD] = status & 0x0f000 ;
694 tskStatus[kANA] = status & 0xf0000 ;
695
696 AliInfo(Form("====> QA Status for %8s raw =0x%x, sim=0x%x, rec=0x%x, esd=0x%x, ana=0x%x", GetDetName(det).Data(),
697 tskStatus[kRAW], tskStatus[kSIM], tskStatus[kREC], tskStatus[kESD], tskStatus[kANA] )) ;
15df7cc8 698 if (tsk == kNULLTASK) {
699 for (Int_t itsk = kRAW ; itsk < kNTASK ; itsk++) {
700 ShowASCIIStatus(det, ALITASK_t(itsk), tskStatus[itsk]) ;
701 }
702 } else {
703 ShowASCIIStatus(det, tsk, tskStatus[tsk]) ;
0b96c27c 704 }
705}
421ab0fb 706
0b96c27c 707//_______________________________________________________________
708void AliQA::ShowASCIIStatus(DETECTORINDEX_t det, ALITASK_t tsk, const ULong_t status) const
709{
710 // print the QA status in human readable format
711 TString text;
712 for (Int_t bit = kINFO ; bit < kNBIT ; bit++) {
713 if (IsSet(det, tsk, QABIT_t(bit))) {
714 text = GetBitName(QABIT_t(bit)) ;
715 text += " " ;
716 }
717 }
718 if (! text.IsNull())
719 printf(" %8s %4s 0x%4x, Problem signalled: %8s \n", GetDetName(det).Data(), GetAliTaskName(tsk), status, text.Data()) ;
421ab0fb 720}
721
96d67a8d 722//_______________________________________________________________
723void AliQA::UnSet(QABIT_t bit)
724{
725 // UnSet the status bit of the current detector in the current module
726
727 UnSetStatusBit(fDet, fTask, bit) ;
728}
729
730//_______________________________________________________________
731void AliQA::UnSetStatusBit(DETECTORINDEX_t det, ALITASK_t tsk, QABIT_t bit)
732{
733 // UnSet the status bit for a given detector and a given task
734
735 CheckRange(det) ;
736 CheckRange(tsk) ;
737 CheckRange(bit) ;
738
739 ULong_t offset = Offset(tsk) ;
740 ULong_t status = GetStatus(det) ;
741 offset+= bit ;
742 status = status & 0 << offset ;
743 SetStatus(det, status) ;
92a357bf 744}