]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQACheckerBase.cxx
5c35416aa671b1dd85da26850116c7a174a02b3b
[u/mrichter/AliRoot.git] / STEER / AliQACheckerBase.cxx
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
17 /* $Id$ */
18
19 //
20 //  Base class for detectors quality assurance checkers 
21 //  Compares Data made by QADataMakers with reference data
22 //  Y. Schutz CERN August 2007
23 //
24
25 // --- ROOT system ---
26 #include <TCanvas.h>
27 #include <TClass.h>
28 #include <TH1F.h> 
29 #include <TH1I.h> 
30 #include <TIterator.h> 
31 #include <TKey.h> 
32 #include <TFile.h> 
33 #include <TList.h>
34 #include <TNtupleD.h>
35 #include <TPaveText.h>
36
37 // --- Standard library ---
38
39 // --- AliRoot header files ---
40 #include "AliCDBEntry.h"
41 #include "AliLog.h"
42 #include "AliQAv1.h"
43 #include "AliQAChecker.h"
44 #include "AliQACheckerBase.h"
45 #include "AliQADataMaker.h"
46 #include "AliQAManager.h"
47 #include "AliDetectorRecoParam.h"
48
49 ClassImp(AliQACheckerBase)
50
51            
52 //____________________________________________________________________________ 
53 AliQACheckerBase::AliQACheckerBase(const char * name, const char * title) : 
54   TNamed(name, title), 
55   fDataSubDir(0x0),
56   fRefSubDir(0x0), 
57   fRefOCDBSubDir(new TObjArray*[AliRecoParam::kNSpecies]), 
58   fLowTestValue(0x0),
59   fUpTestValue(0x0),
60   fImage(new TCanvas*[AliRecoParam::kNSpecies]), 
61   fPrintImage(kTRUE)
62 {
63   // ctor
64   fLowTestValue = new Float_t[AliQAv1::kNBIT] ; 
65   fUpTestValue  = new Float_t[AliQAv1::kNBIT] ; 
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 ; 
74   
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)) ; 
83   }
84   
85   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
86     fImage[specie] = NULL ; 
87     fRefOCDBSubDir[specie] = NULL ;
88   }
89 }
90
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(qac.fLowTestValue),
98   fUpTestValue(qac.fLowTestValue), 
99   fImage(NULL),  
100   fPrintImage(kTRUE)
101 {
102   //copy ctor
103   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
104     fLowTestValue[index]  = qac.fLowTestValue[index] ; 
105     fUpTestValue[index] = qac.fUpTestValue[index] ; 
106   }
107     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
108       fImage[specie] = qac.fImage[specie] ; 
109       fRefOCDBSubDir[specie] = qac.fRefOCDBSubDir[specie] ; 
110     }
111 }
112
113 //____________________________________________________________________________
114 AliQACheckerBase& AliQACheckerBase::operator = (const AliQACheckerBase& qac )
115 {
116   // Equal operator.
117   this->~AliQACheckerBase();
118   new(this) AliQACheckerBase(qac);
119   return *this;
120 }
121
122 //____________________________________________________________________________ 
123 AliQACheckerBase::~AliQACheckerBase()
124 {
125   delete [] fLowTestValue ; 
126   delete [] fUpTestValue ; 
127   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
128     if ( fImage[esIndex] ) 
129       delete fImage[esIndex] ;
130     if ( fRefOCDBSubDir[esIndex] ) 
131       delete fRefOCDBSubDir[esIndex] ; 
132   }
133   delete[] fImage ; 
134   delete[] fRefOCDBSubDir ; 
135   AliQAv1::GetQAResultFile()->Close() ; 
136 }
137
138 //____________________________________________________________________________
139 Double_t * AliQACheckerBase::Check(AliQAv1::ALITASK_t index, AliDetectorRecoParam * recoParam) 
140 {
141   // Performs a basic checking
142   // Compares all the histograms stored in the directory
143   // With reference histograms either in a file of in OCDB  
144
145   TObjArray ** list = new TObjArray *[AliRecoParam::kNSpecies] ; 
146   Int_t specie ;
147   for (specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
148     list[specie] =  new TObjArray(AliQAv1::GetMaxQAObj()) ; 
149     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
150       continue ; 
151     if (fDataSubDir) {
152       TList * keyList = fDataSubDir->GetListOfKeys() ; 
153       TIter next(keyList) ; 
154       TKey * key ;
155       while ( (key = static_cast<TKey *>(next())) ) {
156         TDirectory * specieDir = fDataSubDir->GetDirectory(key->GetName()) ; 
157         TList * keykeyList = specieDir->GetListOfKeys() ; 
158         TIter next2(keykeyList) ; 
159         TKey * keykey ;
160         while ( (keykey = static_cast<TKey *>(next2())) ) {
161           TObject * odata = specieDir->Get(keykey->GetName()) ; 
162           if ( odata->IsA()->InheritsFrom("TH1") ) {
163             TH1 * hdata = static_cast<TH1*>(odata) ;
164             list[specie]->Add(hdata) ; 
165           } else if (!odata->IsA()->InheritsFrom("TDirectory")) // skip the expert directory
166             AliError(Form("%s Is a Classname that cannot be processed", key->GetClassName())) ;
167         }
168       }
169     }
170   }
171  
172   Double_t * test = Check(index, list, recoParam) ;
173   
174   delete[] list ; 
175     
176   return test ;
177 }  
178
179 //____________________________________________________________________________
180 Double_t * AliQACheckerBase::Check(AliQAv1::ALITASK_t task, TObjArray ** list, AliDetectorRecoParam * /*recoParam*/) 
181 {
182   // Performs a basic checking
183   // Compares all the histograms in the list
184
185         Double_t * test = new Double_t[AliRecoParam::kNSpecies] ;
186         Int_t count[AliRecoParam::kNSpecies]   = { 0 }; 
187
188 //  TDirectory * refDir     = NULL ; 
189 //      TObjArray ** refOCDBDir = NULL  ;       
190   GetRefSubDir(GetName(), AliQAv1::GetTaskName(task), fRefSubDir, fRefOCDBSubDir) ;
191  // SetRefandData(refDir, refOCDBDir) ; 
192   
193   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
194     test[specie] = 1.0 ; 
195     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) 
196       continue ; 
197     if (list[specie]->GetEntries() == 0)  
198       test[specie] = 0. ; // nothing to check
199     else {
200       if (!fRefSubDir && !fRefOCDBSubDir)
201         test[specie] = -1 ; // no reference data
202       else {
203         TIter next(list[specie]) ; 
204         TH1 * hdata ;
205         count[specie] = 0 ; 
206         while ( (hdata = static_cast<TH1 *>(next())) ) {
207           if ( hdata->IsA()->InheritsFrom("TH1") ) {
208             if ( hdata->TestBit(AliQAv1::GetExpertBit()) )  // does not perform the test for expert data
209               continue ; 
210             TH1 * href = NULL ; 
211             if (fRefSubDir) 
212               href  = static_cast<TH1*>(fRefSubDir->Get(hdata->GetName())) ;
213             else if (fRefOCDBSubDir[specie])
214               href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hdata->GetName())) ;
215             if (!href) 
216               test[specie] = -1 ; // no reference data ; 
217             else {
218               Double_t rv =  DiffK(hdata, href) ;
219               AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ; 
220               test[specie] += rv ; 
221               count[specie]++ ;
222             }
223           } else
224             AliError("Data type cannot be processed") ;
225           if (count[specie] != 0) 
226             test[specie] /= count[specie] ;
227         }
228       }
229     }
230   }
231   return test ;
232 }  
233
234
235 //____________________________________________________________________________ 
236 Double_t AliQACheckerBase::DiffC(const TH1 * href, const TH1 * hin) const
237 {
238   // compares two histograms using the Chi2 test
239   if ( hin->Integral() == 0 ) {
240     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s is empty", hin->GetName())) ; 
241     return 0. ;
242   }
243     
244   return hin->Chi2Test(href) ;  
245 }
246
247 //____________________________________________________________________________ 
248 Double_t AliQACheckerBase::DiffK(const TH1 * href, const TH1 * hin) const
249 {
250   // compares two histograms using the Kolmogorov test
251   if ( hin->Integral() == 0 || href->Integral() == 0) {
252     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s or its reference is empty", hin->GetName())) ; 
253     return 0. ;
254   }
255     
256   return hin->KolmogorovTest(href) ;  
257 }
258
259   //_____________________________________________________________________________
260 void AliQACheckerBase::GetRefSubDir(const char * det, const char * task, TDirectory *& dirFile, TObjArray **& dirOCDB)     
261
262     // Opens and returns the file with the reference data 
263   dirFile = NULL ; 
264   TString refStorage(AliQAv1::GetQARefStorage()) ;
265   if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
266     AliError(Form("%s is not a valid location for reference data", refStorage.Data())) ; 
267     return ; 
268   } else {
269     AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::GetTaskIndex(task)) ;
270     dirOCDB = new TObjArray*[AliRecoParam::kNSpecies] ; 
271     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
272       dirOCDB[specie] = NULL ; 
273       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
274         continue ; 
275       AliQAv1::SetQARefDataDirName(specie) ;
276       if ( ! manQA->GetLock() ) { 
277         manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ; 
278         manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
279         manQA->SetRun(AliCDBManager::Instance()->GetRun()) ; 
280         manQA->SetLock() ; 
281       }
282       char * detOCDBDir = Form("%s/%s/%s", det, AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
283       AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
284       if (entry) {
285         TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
286         if ( strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
287           AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), det)) ; 
288           listDetQAD = NULL ; 
289           continue ; 
290         } 
291         if ( listDetQAD ) {
292           TIter next(listDetQAD) ;
293           TObjArray * ar ; 
294           while ( (ar = (TObjArray*)next()) ) 
295             dirOCDB[specie] = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", task, AliRecoParam::GetEventSpecieName(specie)))) ;             
296         }
297       }
298     }
299   }
300 }
301
302 //____________________________________________________________________________
303 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, AliDetectorRecoParam * recoParam) 
304
305         AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
306   
307         Double_t * rv = NULL ;
308   rv = Check(index, recoParam) ;
309         SetQA(index, rv) ;      
310         
311   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
312         
313   if (rv) 
314     delete [] rv ; 
315   Finish() ; 
316 }
317
318 //____________________________________________________________________________
319 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, TObjArray ** list, AliDetectorRecoParam * recoParam) 
320
321         AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
322   
323         Double_t * rv = NULL ;
324   rv = Check(index, list, recoParam) ;
325         SetQA(index, rv) ;      
326         
327   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
328         
329   if (rv) 
330     delete [] rv ; 
331   Finish() ; 
332 }
333
334 //____________________________________________________________________________
335 void AliQACheckerBase::Finish() const 
336 {
337         // wrap up and save QA in proper file
338         AliQAv1::GetQAResultFile() ; 
339         AliQAv1 * qa = AliQAv1::Instance() ; 
340         qa->Write(AliQAv1::GetQAName(), kWriteDelete) ;   
341 }
342
343 //____________________________________________________________________________ 
344 void AliQACheckerBase::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode) 
345 {
346   // makes the QA image for sim and rec
347   Int_t nImages = 0 ;
348   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
349     if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) ) 
350       continue ;
351     TIter next(list[esIndex]) ;  
352     TH1 * hdata = NULL ; 
353     while ( (hdata=static_cast<TH1 *>(next())) ) {
354       TString cln(hdata->ClassName()) ; 
355       if ( ! cln.Contains("TH") )
356         continue ; 
357       if ( hdata->TestBit(AliQAv1::GetImageBit()) )
358         nImages++; 
359     }
360     break ; 
361   }
362   if ( nImages == 0 ) {
363     AliDebug(AliQAv1::GetQADebugLevel(), Form("No histogram will be plotted for %s %s\n", GetName(), AliQAv1::GetTaskName(task).Data())) ;  
364   } else {
365     AliDebug(AliQAv1::GetQADebugLevel(), Form("%d histograms will be plotted for %s %s\n", nImages, GetName(), AliQAv1::GetTaskName(task).Data())) ;  
366     for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
367       if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0) 
368         continue ;
369       const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)) ; 
370       if ( !fImage[esIndex] ) {
371         fImage[esIndex] = new TCanvas(title, title) ;
372       }
373       fImage[esIndex]->Clear() ; 
374       fImage[esIndex]->SetTitle(title) ; 
375       fImage[esIndex]->cd() ; 
376       TPaveText someText(0.015, 0.015, 0.98, 0.98) ;
377       someText.AddText(title) ;
378       someText.Draw() ; 
379       fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat())) ; 
380       fImage[esIndex]->Clear() ; 
381       Int_t nx = TMath::Sqrt(nImages) ; 
382       Int_t ny = nx  ; 
383       while ( nx*ny <= nImages) 
384         ny++ ; 
385       
386       fImage[esIndex]->Divide(nx, ny) ; 
387       TIter nexthist(list[esIndex]) ; 
388       TH1* hist = NULL ;
389       Int_t npad = 1 ; 
390       fImage[esIndex]->cd(npad) ; 
391       while ( (hist=static_cast<TH1*>(nexthist())) ) {
392         TString cln(hist->ClassName()) ; 
393         if ( ! cln.Contains("TH") )
394           continue ; 
395         if(hist->TestBit(AliQAv1::GetImageBit())) {
396           hist->Draw() ; 
397           fImage[esIndex]->cd(++npad) ; 
398         }
399       }
400       fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat())) ; 
401     }
402   }  
403 }
404
405 //____________________________________________________________________________
406 void AliQACheckerBase::SetHiLo(Float_t * hiValue, Float_t * lowValue) 
407 {
408   AliDebug(AliQAv1::GetQADebugLevel(), "Previous setting was:") ;
409   if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
410     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", 
411                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
412                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
413                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
414                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
415     AliInfo(Form("%s", text)) ; 
416   }
417   
418   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
419     fLowTestValue[index]  = lowValue[index] ; 
420     fUpTestValue[index]   = hiValue[index] ; 
421   }
422   AliDebug(AliQAv1::GetQADebugLevel(), "Current setting is:") ;
423   if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
424     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", 
425                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
426                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
427                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
428                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;     AliInfo(Form("%s", text)) ; 
429   }
430 }
431
432 //____________________________________________________________________________
433 void AliQACheckerBase::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
434 {
435         // sets the QA according the return value of the Check
436
437   AliQAv1 * qa = AliQAv1::Instance(index) ;
438
439   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
440     if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
441       continue ;
442     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
443       qa->Set(AliQAv1::kFATAL, specie) ; 
444     } else {
445       if ( value[specie] >= fLowTestValue[AliQAv1::kFATAL] && value[specie] < fUpTestValue[AliQAv1::kFATAL] ) 
446         qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
447       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
448         qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
449       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
450         qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
451       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
452         qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;   
453     }
454   }
455 }