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