]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQACheckerBase.cxx
MakeImage is now a method of AliCheckerBase, was AliQADataMaker before. This will...
[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 "AliLog.h"
41 #include "AliQAv1.h"
42 #include "AliQAChecker.h"
43 #include "AliQACheckerBase.h"
44 #include "AliQADataMaker.h"
45
46 ClassImp(AliQACheckerBase)
47
48            
49 //____________________________________________________________________________ 
50 AliQACheckerBase::AliQACheckerBase(const char * name, const char * title) : 
51   TNamed(name, title), 
52   fDataSubDir(0x0),
53   fRefSubDir(0x0), 
54   fRefOCDBSubDir(0x0), 
55   fLowTestValue(0x0),
56   fUpTestValue(0x0),
57   fImage(new TCanvas*[AliRecoParam::kNSpecies]), 
58   fPrintImage(kTRUE)
59 {
60   // ctor
61   fLowTestValue = new Float_t[AliQAv1::kNBIT] ; 
62   fUpTestValue  = new Float_t[AliQAv1::kNBIT] ; 
63   fLowTestValue[AliQAv1::kINFO]    =  0.5   ; 
64   fUpTestValue[AliQAv1::kINFO]     = 1.0 ; 
65   fLowTestValue[AliQAv1::kWARNING] =  0.002 ; 
66   fUpTestValue[AliQAv1::kWARNING]  = 0.5 ; 
67   fLowTestValue[AliQAv1::kERROR]   =  0.0   ; 
68   fUpTestValue[AliQAv1::kERROR]    = 0.002 ; 
69   fLowTestValue[AliQAv1::kFATAL]   = -1.0   ; 
70   fUpTestValue[AliQAv1::kFATAL]    = 0.0 ; 
71   
72   AliDebug(AliQAv1::GetQADebugLevel(), "Default setting is:") ;
73   if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
74     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", 
75                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
76                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
77                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
78                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
79     AliInfo(Form("%s", text)) ; 
80   }
81   
82   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
83     fImage[specie] = NULL ; 
84 }
85
86 //____________________________________________________________________________ 
87 AliQACheckerBase::AliQACheckerBase(const AliQACheckerBase& qac) :
88   TNamed(qac.GetName(), qac.GetTitle()),
89   fDataSubDir(qac.fDataSubDir), 
90   fRefSubDir(qac.fRefSubDir), 
91   fRefOCDBSubDir(qac.fRefOCDBSubDir), 
92   fLowTestValue(qac.fLowTestValue),
93   fUpTestValue(qac.fLowTestValue), 
94   fImage(NULL),  
95   fPrintImage(kTRUE)
96 {
97   //copy ctor
98   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
99     fLowTestValue[index]  = qac.fLowTestValue[index] ; 
100     fUpTestValue[index] = qac.fUpTestValue[index] ; 
101   }
102   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
103     fImage[specie] = qac.fImage[specie] ; 
104 }
105
106 //____________________________________________________________________________
107 AliQACheckerBase& AliQACheckerBase::operator = (const AliQACheckerBase& qac )
108 {
109   // Equal operator.
110   this->~AliQACheckerBase();
111   new(this) AliQACheckerBase(qac);
112   return *this;
113 }
114
115 //____________________________________________________________________________ 
116 AliQACheckerBase::~AliQACheckerBase()
117 {
118   delete [] fLowTestValue ; 
119   delete [] fUpTestValue ; 
120   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
121     if ( fImage[esIndex] ) 
122       delete fImage[esIndex] ;
123   }
124   delete[] fImage ; 
125 }
126
127 //____________________________________________________________________________
128 Double_t * AliQACheckerBase::Check(AliQAv1::ALITASK_t /*index*/) 
129 {
130   // Performs a basic checking
131   // Compares all the histograms stored in the directory
132   // With reference histograms either in a file of in OCDB  
133
134         Double_t * test = new Double_t[AliRecoParam::kNSpecies] ;
135         Int_t count[AliRecoParam::kNSpecies]   = { 0 }; 
136
137   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
138     test[specie] = 1.0 ; 
139     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
140       continue ; 
141     if (!fDataSubDir) {
142       test[specie] = 0. ; // nothing to check
143     } else if (!fRefSubDir && !fRefOCDBSubDir) {
144         test[specie] = -1 ; // no reference data
145     } else {
146       TList * keyList = fDataSubDir->GetListOfKeys() ; 
147       TIter next(keyList) ; 
148       TKey * key ;
149       count[specie] = 0 ; 
150       while ( (key = static_cast<TKey *>(next())) ) {
151         TObject * odata = fRefSubDir->Get(key->GetName()) ; 
152         if ( odata->IsA()->InheritsFrom("TH1") ) {
153           TH1 * hdata = static_cast<TH1*>(odata) ;
154           TH1 * href = NULL ; 
155           if (fRefSubDir) 
156             href  = static_cast<TH1*>(fRefSubDir->Get(key->GetName())) ;
157           else if (fRefOCDBSubDir[specie]) {  
158             href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(key->GetName())) ;
159           }
160           if (!href) 
161             test[specie] = -1 ; // no reference data ; 
162           else {
163             Double_t rv =  DiffK(hdata, href) ;
164             AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ; 
165             test[specie] += rv ; 
166             count[specie]++ ; 
167           }
168         } else
169           AliError(Form("%s Is a Classname that cannot be processed", key->GetClassName())) ;
170       }
171       if (count[specie] != 0) 
172         test[specie] /= count[specie] ;
173     }
174   }
175         return test ;
176 }  
177
178 //____________________________________________________________________________
179 Double_t * AliQACheckerBase::Check(AliQAv1::ALITASK_t /*index*/, TObjArray ** list) 
180 {
181   // Performs a basic checking
182   // Compares all the histograms in the list
183
184         Double_t * test = new Double_t[AliRecoParam::kNSpecies] ;
185         Int_t count[AliRecoParam::kNSpecies]   = { 0 }; 
186
187   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
188     test[specie] = 1.0 ; 
189     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) 
190       continue ; 
191     if (list[specie]->GetEntries() == 0)  
192       test[specie] = 0. ; // nothing to check
193     else {
194       if (!fRefSubDir && !fRefOCDBSubDir)
195         test[specie] = -1 ; // no reference data
196       else {
197         TIter next(list[specie]) ; 
198         TH1 * hdata ;
199         count[specie] = 0 ; 
200         while ( (hdata = static_cast<TH1 *>(next())) ) {
201           TString cln(hdata->ClassName()) ; 
202           if ( ! cln.Contains("TH1") )
203             continue ;           
204           if ( hdata) { 
205             if ( hdata->TestBit(AliQAv1::GetExpertBit()) )  // does not perform the test for expert data
206               continue ; 
207             TH1 * href = NULL ; 
208             if (fRefSubDir) 
209               href  = static_cast<TH1*>(fRefSubDir->Get(hdata->GetName())) ;
210             else if (fRefOCDBSubDir[specie])
211               href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hdata->GetName())) ;
212             if (!href) 
213               test[specie] = -1 ; // no reference data ; 
214             else {
215               Double_t rv =  DiffK(hdata, href) ;
216               AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ; 
217               test[specie] += rv ; 
218               count[specie]++ ; 
219             }
220           } 
221           else
222             AliError("Data type cannot be processed") ;
223         }
224         if (count[specie] != 0) 
225           test[specie] /= count[specie] ;
226       }
227     }
228   }
229         return test ;
230 }  
231
232
233 //____________________________________________________________________________ 
234 Double_t AliQACheckerBase::DiffC(const TH1 * href, const TH1 * hin) const
235 {
236   // compares two histograms using the Chi2 test
237   if ( hin->Integral() == 0 ) {
238     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s is empty", hin->GetName())) ; 
239     return 0. ;
240   }
241     
242   return hin->Chi2Test(href) ;  
243 }
244
245 //____________________________________________________________________________ 
246 Double_t AliQACheckerBase::DiffK(const TH1 * href, const TH1 * hin) const
247 {
248   // compares two histograms using the Kolmogorov test
249   if ( hin->Integral() == 0 || href->Integral() == 0) {
250     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s or its reference is empty", hin->GetName())) ; 
251     return 0. ;
252   }
253     
254   return hin->KolmogorovTest(href) ;  
255 }
256
257 //____________________________________________________________________________
258 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, TObjArray ** list) 
259
260         AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
261   
262         Double_t * rv = NULL ;
263   if ( !list) 
264     rv = Check(index) ;
265   else 
266     rv = Check(index, list) ;
267         SetQA(index, rv) ;      
268         
269   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
270         
271   if (rv) 
272     delete [] rv ; 
273   Finish() ; 
274 }
275
276 //____________________________________________________________________________
277 void AliQACheckerBase::Finish() const 
278 {
279         // wrap up and save QA in proper file
280         AliQAv1 * qa = AliQAv1::Instance() ; 
281         //qa->Show() ;
282         AliQAv1::GetQAResultFile()->cd() ; 
283         qa->Write(qa->GetName(), kWriteDelete) ;   
284         AliQAv1::GetQAResultFile()->Close() ; 
285 }
286
287 //____________________________________________________________________________ 
288 void AliQACheckerBase::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode) 
289 {
290   // makes the QA image for sim and rec
291   Int_t nImages = 0 ;
292   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
293     if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) ) 
294       continue ;
295     TIter next(list[esIndex]) ;  
296     TH1 * hdata = NULL ; 
297     while ( (hdata=static_cast<TH1 *>(next())) ) {
298       TString cln(hdata->ClassName()) ; 
299       if ( ! cln.Contains("TH1") )
300         continue ; 
301       if ( hdata->TestBit(AliQAv1::GetImageBit()) )
302         nImages++; 
303     }
304     break ; 
305   }
306   if ( nImages == 0 ) {
307     AliDebug(AliQAv1::GetQADebugLevel(), Form("No histogram will be plotted for %s %s\n", GetName(), AliQAv1::GetTaskName(task).Data())) ;  
308   } else {
309     AliDebug(AliQAv1::GetQADebugLevel(), Form("%d histograms will be plotted for %s %s\n", nImages, GetName(), AliQAv1::GetTaskName(task).Data())) ;  
310     for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
311       if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) ) 
312         continue ;
313       const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)) ; 
314       if ( !fImage[esIndex] ) {
315         fImage[esIndex] = new TCanvas(title, title) ;
316       }
317       fImage[esIndex]->Clear() ; 
318       fImage[esIndex]->SetTitle(title) ; 
319       fImage[esIndex]->cd() ; 
320       TPaveText someText(0.015, 0.015, 0.98, 0.98) ;
321       someText.AddText(title) ;
322       someText.Draw() ; 
323       fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat())) ; 
324       fImage[esIndex]->Clear() ; 
325       Int_t nx = TMath::Sqrt(nImages) ; 
326       Int_t ny = nx  ; 
327       if ( nx < TMath::Sqrt(nImages)) 
328         ny++ ; 
329       fImage[esIndex]->Divide(nx, ny) ; 
330       TIter nexthist(list[esIndex]) ; 
331       TH1* hist = NULL ;
332       Int_t npad = 1 ; 
333       fImage[esIndex]->cd(npad) ; 
334       while ( (hist=static_cast<TH1*>(nexthist())) ) {
335         TString cln(hist->ClassName()) ; 
336         if ( ! cln.Contains("TH1") )
337           continue ; 
338         if(hist->TestBit(AliQAv1::GetImageBit())) {
339           hist->Draw() ; 
340           fImage[esIndex]->cd(++npad) ; 
341         }
342       }
343       fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat())) ; 
344     }
345   }  
346 }
347
348 //____________________________________________________________________________
349 void AliQACheckerBase::SetHiLo(Float_t * hiValue, Float_t * lowValue) 
350 {
351   AliDebug(AliQAv1::GetQADebugLevel(), "Previous setting was:") ;
352   if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
353     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", 
354                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
355                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
356                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
357                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
358     AliInfo(Form("%s", text)) ; 
359   }
360   
361   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
362     fLowTestValue[index]  = lowValue[index] ; 
363     fUpTestValue[index]   = hiValue[index] ; 
364   }
365   AliDebug(AliQAv1::GetQADebugLevel(), "Current setting is:") ;
366   if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
367     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", 
368                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
369                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
370                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
371                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;     AliInfo(Form("%s", text)) ; 
372   }
373 }
374
375 //____________________________________________________________________________
376 void AliQACheckerBase::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
377 {
378         // sets the QA according the return value of the Check
379
380   AliQAv1 * qa = AliQAv1::Instance(index) ;
381   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
382     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
383       qa->Set(AliQAv1::kFATAL, specie) ; 
384     } else {
385       if ( value[specie] >= fLowTestValue[AliQAv1::kFATAL] && value[specie] < fUpTestValue[AliQAv1::kFATAL] ) 
386         qa->Set(AliQAv1::kFATAL, specie) ; 
387       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
388         qa->Set(AliQAv1::kERROR, specie) ; 
389       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
390         qa->Set(AliQAv1::kWARNING, specie) ;
391       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
392         qa->Set(AliQAv1::kINFO, specie) ;       
393     }
394   }
395 }