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 **************************************************************************/
20 // Base class for detectors quality assurance checkers
21 // Compares Data made by QADataMakers with reference data
22 // Y. Schutz CERN August 2007
25 // --- ROOT system ---
30 #include <TIterator.h>
35 #include <TParameter.h>
36 #include <TPaveText.h>
38 // --- Standard library ---
40 // --- AliRoot header files ---
41 #include "AliCDBEntry.h"
44 #include "AliQAChecker.h"
45 #include "AliQACheckerBase.h"
46 #include "AliQADataMaker.h"
47 #include "AliQAManager.h"
48 #include "AliDetectorRecoParam.h"
50 ClassImp(AliQACheckerBase)
53 //____________________________________________________________________________
54 AliQACheckerBase::AliQACheckerBase(const char * name, const char * title) :
58 fRefOCDBSubDir(new TObjArray*[AliRecoParam::kNSpecies]),
59 fLowTestValue(new Float_t[AliQAv1::kNBIT]),
60 fUpTestValue(new Float_t[AliQAv1::kNBIT]),
61 fImage(new TCanvas*[AliRecoParam::kNSpecies]),
63 fExternParamList(new TList())
66 fLowTestValue[AliQAv1::kINFO] = 0.5 ;
67 fUpTestValue[AliQAv1::kINFO] = 1.0 ;
68 fLowTestValue[AliQAv1::kWARNING] = 0.002 ;
69 fUpTestValue[AliQAv1::kWARNING] = 0.5 ;
70 fLowTestValue[AliQAv1::kERROR] = 0.0 ;
71 fUpTestValue[AliQAv1::kERROR] = 0.002 ;
72 fLowTestValue[AliQAv1::kFATAL] = -1.0 ;
73 fUpTestValue[AliQAv1::kFATAL] = 0.0 ;
75 AliDebug(AliQAv1::GetQADebugLevel(), "Default setting is:") ;
76 if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
77 const Char_t * text= Form(" INFO -> %1.5f < value < %1.5f WARNING -> %1.5f < value <= %1.5f \n ERROR -> %1.5f < value <= %1.5f \n FATAL -> %1.5f <= value < %1.5f \n",
78 fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO],
79 fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING],
80 fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR],
81 fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;
82 AliInfo(Form("%s", text)) ;
85 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
86 fImage[specie] = NULL ;
87 fRefOCDBSubDir[specie] = NULL ;
91 //____________________________________________________________________________
92 AliQACheckerBase::AliQACheckerBase(const AliQACheckerBase& qac) :
93 TNamed(qac.GetName(), qac.GetTitle()),
94 fDataSubDir(qac.fDataSubDir),
95 fRefSubDir(qac.fRefSubDir),
96 fRefOCDBSubDir(qac.fRefOCDBSubDir),
97 fLowTestValue(new Float_t[AliQAv1::kNBIT]),
98 fUpTestValue(new Float_t[AliQAv1::kNBIT]),
99 fImage(new TCanvas*[AliRecoParam::kNSpecies]),
101 fExternParamList(new TList())
104 for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
105 fLowTestValue[index] = qac.fLowTestValue[index] ;
106 fUpTestValue[index] = qac.fUpTestValue[index] ;
108 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
109 fImage[specie] = qac.fImage[specie] ;
110 fRefOCDBSubDir[specie] = qac.fRefOCDBSubDir[specie] ;
112 if (qac.fExternParamList) {
113 TIter next(qac.fExternParamList) ;
114 TParameter<double> * p ;
115 while ( (p = (TParameter<double>*)next()) )
116 fExternParamList->Add(p) ;
120 //____________________________________________________________________________
121 AliQACheckerBase& AliQACheckerBase::operator = (const AliQACheckerBase& qac )
124 this->~AliQACheckerBase();
125 new(this) AliQACheckerBase(qac);
129 //____________________________________________________________________________
130 AliQACheckerBase::~AliQACheckerBase()
132 delete [] fLowTestValue ;
133 delete [] fUpTestValue ;
136 delete[] fRefOCDBSubDir ;
137 AliQAv1::GetQAResultFile()->Close() ;
138 if (fExternParamList) {
139 fExternParamList->Clear() ;
140 delete fExternParamList ;
144 //____________________________________________________________________________
145 void AliQACheckerBase::Check(Double_t * test, AliQAv1::ALITASK_t index, const AliDetectorRecoParam * recoParam)
147 // Performs a basic checking
148 // Compares all the histograms stored in the directory
149 // With reference histograms either in a file of in OCDB
151 TObjArray ** list = new TObjArray *[AliRecoParam::kNSpecies] ;
153 for (specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
154 list[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ;
155 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
158 TList * keyList = fDataSubDir->GetListOfKeys() ;
159 TIter next(keyList) ;
161 while ( (key = static_cast<TKey *>(next())) ) {
162 TDirectory * specieDir = fDataSubDir->GetDirectory(key->GetName()) ;
163 TList * keykeyList = specieDir->GetListOfKeys() ;
164 TIter next2(keykeyList) ;
166 while ( (keykey = static_cast<TKey *>(next2())) ) {
167 TObject * odata = specieDir->Get(keykey->GetName()) ;
168 if ( odata->IsA()->InheritsFrom("TH1") ) {
169 TH1 * hdata = static_cast<TH1*>(odata) ;
170 list[specie]->Add(hdata) ;
171 } else if (!odata->IsA()->InheritsFrom("TDirectory")) // skip the expert directory
172 AliError(Form("%s Is a Classname that cannot be processed", key->GetClassName())) ;
178 Check(test, index, list, recoParam) ;
184 //____________________________________________________________________________
185 void AliQACheckerBase::Check(Double_t * test, AliQAv1::ALITASK_t task, TObjArray ** list, const AliDetectorRecoParam * recoParam)
187 // Performs a basic checking
188 // Compares all the histograms in the list
190 Int_t count[AliRecoParam::kNSpecies] = { 0 };
192 GetRefSubDir(GetName(), AliQAv1::GetTaskName(task), fRefSubDir, fRefOCDBSubDir) ;
193 // SetRefandData(refDir, refOCDBDir) ;
195 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
197 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) continue ;
198 if (list[specie]->GetEntries() == 0)
199 test[specie] = 0. ; // nothing to check
201 if (!fRefSubDir && !fRefOCDBSubDir)
202 test[specie] = -1 ; // no reference data
204 TIter next(list[specie]) ;
207 while ( (hdata = static_cast<TH1 *>(next())) ) {
208 if ( hdata->IsA()->InheritsFrom("TH1") ) {
209 if ( hdata->TestBit(AliQAv1::GetExpertBit()) ) // does not perform the test for expert data
212 // First try to find the ref histo with exact name (with possible trigger clon ending)
213 TString hname = hdata->GetName();
215 if (fRefSubDir) href = static_cast<TH1*>(fRefSubDir->Get(hname.Data())) ;
216 else if (fRefOCDBSubDir[specie]) href = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hname.Data()));
218 if (!href && hdata->TestBit(AliQAv1::GetClonedBit())) { // try to find the histo for the base name (w/o trigger ending
219 int ind = hname.Index(AliQADataMaker::GetTriggerPrefix());
222 if (fRefSubDir) href = static_cast<TH1*>(fRefSubDir->Get(hname.Data())) ;
223 else if (fRefOCDBSubDir[specie]) href = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hname.Data()));
228 test[specie] = -1 ; // no reference data ;
230 Double_t rv = DiffK(hdata, href) ;
231 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ;
236 AliError("Data type cannot be processed") ;
237 if (count[specie] != 0)
238 test[specie] /= count[specie] ;
246 //____________________________________________________________________________
247 void AliQACheckerBase::DeleteImages()
250 for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
251 if ( fImage[esIndex] ) {delete fImage[esIndex]; fImage[esIndex] = 0;}
252 if ( fRefOCDBSubDir[esIndex] ) {delete fRefOCDBSubDir[esIndex]; fRefOCDBSubDir[esIndex] = 0;}
256 //____________________________________________________________________________
257 Double_t AliQACheckerBase::DiffC(const TH1 * href, const TH1 * hin) const
259 // compares two histograms using the Chi2 test
260 if ( hin->Integral() == 0 ) {
261 AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s is empty", hin->GetName())) ;
265 return hin->Chi2Test(href) ;
268 //____________________________________________________________________________
269 Double_t AliQACheckerBase::DiffK(const TH1 * href, const TH1 * hin) const
271 // compares two histograms using the Kolmogorov test
272 if ( hin->Integral() == 0 || href->Integral() == 0) {
273 AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s or its reference is empty", hin->GetName())) ;
277 return hin->KolmogorovTest(href) ;
280 //_____________________________________________________________________________
281 void AliQACheckerBase::GetRefSubDir(const char * det, const char * task, TDirectory *& dirFile, TObjArray **& dirOCDB)
283 // Opens and returns the file with the reference data
285 TString refStorage(AliQAv1::GetQARefStorage()) ;
286 if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
287 AliError(Form("%s is not a valid location for reference data", refStorage.Data())) ;
290 AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::GetTaskIndex(task)) ;
291 // dirOCDB = new TObjArray*[AliRecoParam::kNSpecies] ;
292 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
293 dirOCDB[specie] = NULL ;
294 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
296 AliQAv1::SetQARefDataDirName(specie) ;
297 if ( ! manQA->GetLock() ) {
298 manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ;
299 manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
300 manQA->SetRun(AliCDBManager::Instance()->GetRun()) ;
303 char * detOCDBDir = Form("%s/%s/%s", det, AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
304 AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
306 TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
307 if ( listDetQAD && strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
308 AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), det)) ;
313 TIter next(listDetQAD) ;
314 while ( (TObjArray*)next() )
315 dirOCDB[specie] = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", task, AliRecoParam::GetEventSpecieName(specie)))) ;
322 //____________________________________________________________________________
323 void AliQACheckerBase::PrintExternParam()
325 // Print the list of external parameter list
326 TIter next(fExternParamList) ;
327 TParameter<double> *pp ;
328 TString printit("\n") ;
329 while( (pp = (TParameter<double>*)next()) )
330 printit += Form("%s = %f\n", pp->GetName(), pp->GetVal());
331 AliInfo(Form("%s", printit.Data())) ;
334 //____________________________________________________________________________
335 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, const AliDetectorRecoParam * recoParam)
337 //Run the checker for all kind of species
338 AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ;
340 Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
341 for (int i=AliRecoParam::kNSpecies;i--;) rv[i] = 0.0;
342 Check(rv, index, recoParam) ;
345 AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
351 //____________________________________________________________________________
352 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * recoParam)
354 // RS: perform check for all trigger classes in loop
355 Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
357 TObjArray ** listTrig = new TObjArray *[AliRecoParam::kNSpecies];
359 for (int itc=-1;itc<AliQADataMaker::GetNTrigClasses();itc++) {
361 // RS: fetch the histograms for each specie and for given trigger
362 //AliInfo(Form("Processing %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc)));
364 for (int specie=0;specie<AliRecoParam::kNSpecies;specie++) {
365 listTrig[specie] = 0;
366 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) || !list[specie]) continue;
367 listTrig[specie] = new TObjArray( list[specie]->GetSize() ); // destination for clones of this trigger
368 AliQADataMaker::GetDataOfTrigClass(list[specie],itc, listTrig[specie]);
370 AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc)));
371 Check(rv, index, listTrig, recoParam) ;
373 AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc)));
375 for (int specie=0;specie<AliRecoParam::kNSpecies;specie++) if (listTrig[specie]) delete listTrig[specie]; // clean temporary container
382 //____________________________________________________________________________
383 void AliQACheckerBase::Finish() const
385 // wrap up and save QA in proper file
386 AliQAv1::GetQAResultFile() ;
387 AliQAv1 * qa = AliQAv1::Instance() ;
388 qa->Write(AliQAv1::GetQAName(), kWriteDelete) ;
391 //____________________________________________________________________________
392 void AliQACheckerBase::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode)
394 // makes the QA image for sim and rec
395 TObjArray tmpArr; // array to store flat version of original array (which may contain clones)
397 for (Int_t esIndex = 0; esIndex < AliRecoParam::kNSpecies; esIndex++) {
398 if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0) continue;
400 TIter next(list[esIndex]);
401 TObject* hdata = NULL;
403 while ( (hdata=(next())) ) { // count histos and transfere to flat array
404 if (hdata->InheritsFrom(TH1::Class()) && hdata->TestBit(AliQAv1::GetImageBit()) ) { // histo, not cloned
406 tmpArr.AddLast(hdata);
409 if (!hdata->TestBit(AliQAv1::GetClonedBit())) continue; // not an array of clones, unknown object
410 TIter nextCl((TObjArray*)hdata); // array of histo clones
412 while ((hcl=nextCl())) if (hcl->InheritsFrom(TH1::Class()) && hcl->TestBit(AliQAv1::GetImageBit())) {tmpArr.AddLast(hcl); nImages++;}
415 if ( nImages == 0 ) {
416 AliDebug(AliQAv1::GetQADebugLevel(), Form("No histogram will be plotted for %s %s %s\n", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)));
419 AliDebug(AliQAv1::GetQADebugLevel(), Form("%d histograms will be plotted for %s %s %s\n", nImages, GetName(), AliQAv1::GetTaskName(task).Data(),AliRecoParam::GetEventSpecieName(esIndex)));
421 const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex));
423 if ( !fImage[esIndex] ) fImage[esIndex] = new TCanvas(title, title);
425 fImage[esIndex]->Clear();
426 fImage[esIndex]->SetTitle(title);
427 fImage[esIndex]->cd();
428 TPaveText someText(0.015, 0.015, 0.98, 0.98);
429 someText.AddText(title);
431 fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()));
432 fImage[esIndex]->Clear();
433 Int_t nx = TMath::Nint(TMath::Sqrt(nImages));
435 if (nx < TMath::Sqrt(nImages)) ny++;
437 fImage[esIndex]->Divide(nx, ny);
438 TIter nexthist(&tmpArr);
440 fImage[esIndex]->cd(npad);
442 while ( (histo=(TH1*)nexthist()) ) { // tmpArr is guaranteed to contain only plottable histos, no checks needed
443 TString opts = histo->GetDrawOption();
444 if (opts.Contains("logy",TString::kIgnoreCase)) gPad->SetLogy();
445 if (opts.Contains("logx",TString::kIgnoreCase)) gPad->SetLogx();
447 fImage[esIndex]->cd(++npad);
449 fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps");
453 //____________________________________________________________________________
454 void AliQACheckerBase::SetHiLo(Float_t * hiValue, Float_t * lowValue)
456 AliDebug(AliQAv1::GetQADebugLevel(), "Previous setting was:") ;
457 if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
458 const Char_t * text= Form(" INFO -> %1.5f < value < %1.5f WARNING -> %1.5f < value <= %1.5f \n ERROR -> %1.5f < value <= %1.5f \n FATAL -> %1.5f <= value < %1.5f \n",
459 fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO],
460 fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING],
461 fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR],
462 fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;
463 AliInfo(Form("%s", text)) ;
466 for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
467 fLowTestValue[index] = lowValue[index] ;
468 fUpTestValue[index] = hiValue[index] ;
470 AliDebug(AliQAv1::GetQADebugLevel(), "Current setting is:") ;
471 if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
472 const Char_t * text= Form(" INFO -> %1.5f < value < %1.5f WARNING -> %1.5f < value <= %1.5f \n ERROR -> %1.5f < value <= %1.5f \n FATAL -> %1.5f <= value < %1.5f \n",
473 fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO],
474 fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING],
475 fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR],
476 fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; AliInfo(Form("%s", text)) ;
480 //____________________________________________________________________________
481 void AliQACheckerBase::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
483 // sets the QA according the return value of the Check
485 AliQAv1 * qa = AliQAv1::Instance(index) ;
487 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
488 if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
490 if ( value == NULL ) { // No checker is implemented, set all QA to Fatal
491 qa->Set(AliQAv1::kFATAL, specie) ;
493 if ( value[specie] >= fLowTestValue[AliQAv1::kFATAL] && value[specie] < fUpTestValue[AliQAv1::kFATAL] )
494 qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ;
495 else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR] )
496 qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ;
497 else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING] )
498 qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
499 else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] )
500 qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;