]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEER/AliQACheckerBase.cxx
Moving the classes that belong to the following libraries: STEERBase, ESD, CDB, AOD...
[u/mrichter/AliRoot.git] / STEER / 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 <TParameter.h>
36 #include <TPaveText.h>
37
38 // --- Standard library ---
39
40 // --- AliRoot header files ---
41 #include "AliCDBEntry.h"
42 #include "AliLog.h"
43 #include "AliQAv1.h"
44 #include "AliQAChecker.h"
45 #include "AliQACheckerBase.h"
46 #include "AliQADataMaker.h"
47 #include "AliQAManager.h"
48 #include "AliDetectorRecoParam.h"
49
50 ClassImp(AliQACheckerBase)
51
52            
53 //____________________________________________________________________________ 
54 AliQACheckerBase::AliQACheckerBase(const char * name, const char * title) : 
55   TNamed(name, title), 
56   fDataSubDir(0x0),
57   fRefSubDir(0x0), 
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]), 
62   fPrintImage(kTRUE), 
63   fExternParamList(new TList())
64 {
65   // ctor
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(new Float_t[AliQAv1::kNBIT]),
98   fUpTestValue(new Float_t[AliQAv1::kNBIT]), 
99   fImage(new TCanvas*[AliRecoParam::kNSpecies]),  
100   fPrintImage(kTRUE), 
101   fExternParamList(new TList())  
102 {
103   //copy ctor
104   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
105     fLowTestValue[index]  = qac.fLowTestValue[index] ; 
106     fUpTestValue[index] = qac.fUpTestValue[index] ; 
107   }
108     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
109       fImage[specie] = qac.fImage[specie] ; 
110       fRefOCDBSubDir[specie] = qac.fRefOCDBSubDir[specie] ; 
111     }
112   if (qac.fExternParamList) {
113     TIter next(qac.fExternParamList) ; 
114     TParameter<double> * p ; 
115     while ( (p = (TParameter<double>*)next()) )
116       fExternParamList->Add(p) ;
117   }
118 }
119
120 //____________________________________________________________________________
121 AliQACheckerBase& AliQACheckerBase::operator = (const AliQACheckerBase& qac )
122 {
123   // Equal operator.
124   this->~AliQACheckerBase();
125   new(this) AliQACheckerBase(qac);
126   return *this;
127 }
128
129 //____________________________________________________________________________ 
130 AliQACheckerBase::~AliQACheckerBase()
131 {
132   delete [] fLowTestValue ; 
133   delete [] fUpTestValue ; 
134   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
135     if ( fImage[esIndex] ) 
136       delete fImage[esIndex] ;
137     if ( fRefOCDBSubDir[esIndex] ) 
138       delete fRefOCDBSubDir[esIndex] ; 
139   }
140   delete[] fImage ; 
141   delete[] fRefOCDBSubDir ; 
142   AliQAv1::GetQAResultFile()->Close() ; 
143   if (fExternParamList) {
144     fExternParamList->Clear() ; 
145     delete fExternParamList ; 
146   }
147 }
148
149 //____________________________________________________________________________
150 void AliQACheckerBase::Check(Double_t * test, AliQAv1::ALITASK_t index, const AliDetectorRecoParam * recoParam) 
151 {
152   // Performs a basic checking
153   // Compares all the histograms stored in the directory
154   // With reference histograms either in a file of in OCDB  
155
156   TObjArray ** list = new TObjArray *[AliRecoParam::kNSpecies] ; 
157   Int_t specie ;
158   for (specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
159     list[specie] =  new TObjArray(AliQAv1::GetMaxQAObj()) ; 
160     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
161       continue ; 
162     if (fDataSubDir) {
163       TList * keyList = fDataSubDir->GetListOfKeys() ; 
164       TIter next(keyList) ; 
165       TKey * key ;
166       while ( (key = static_cast<TKey *>(next())) ) {
167         TDirectory * specieDir = fDataSubDir->GetDirectory(key->GetName()) ; 
168         TList * keykeyList = specieDir->GetListOfKeys() ; 
169         TIter next2(keykeyList) ; 
170         TKey * keykey ;
171         while ( (keykey = static_cast<TKey *>(next2())) ) {
172           TObject * odata = specieDir->Get(keykey->GetName()) ; 
173           if ( odata->IsA()->InheritsFrom("TH1") ) {
174             TH1 * hdata = static_cast<TH1*>(odata) ;
175             list[specie]->Add(hdata) ; 
176           } else if (!odata->IsA()->InheritsFrom("TDirectory")) // skip the expert directory
177             AliError(Form("%s Is a Classname that cannot be processed", key->GetClassName())) ;
178         }
179       }
180     }
181   }
182  
183   Check(test, index, list, recoParam) ;
184   
185   delete[] list ; 
186     
187 }  
188
189 //____________________________________________________________________________
190 void AliQACheckerBase::Check(Double_t * test, AliQAv1::ALITASK_t task, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/) 
191 {
192   // Performs a basic checking
193   // Compares all the histograms in the list
194
195   Int_t count[AliRecoParam::kNSpecies]   = { 0 }; 
196
197   GetRefSubDir(GetName(), AliQAv1::GetTaskName(task), fRefSubDir, fRefOCDBSubDir) ;
198  // SetRefandData(refDir, refOCDBDir) ; 
199   
200   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
201     test[specie] = 1.0 ; 
202     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) continue ; 
203     if (list[specie]->GetEntries() == 0)  
204       test[specie] = 0. ; // nothing to check
205     else {
206       if (!fRefSubDir && !fRefOCDBSubDir)
207         test[specie] = -1 ; // no reference data
208       else {
209         TIter next(list[specie]) ; 
210         TH1 * hdata ;
211         count[specie] = 0 ; 
212         while ( (hdata = static_cast<TH1 *>(next())) ) {
213           if ( hdata->IsA()->InheritsFrom("TH1") ) {
214             if ( hdata->TestBit(AliQAv1::GetExpertBit()) )  // does not perform the test for expert data
215               continue ; 
216             // 
217             // First try to find the ref histo with exact name (with possible trigger clon ending)
218             TString hname = hdata->GetName();
219             TH1 * href = NULL ; 
220             if (fRefSubDir)                  href  = static_cast<TH1*>(fRefSubDir->Get(hname.Data())) ;
221             else if (fRefOCDBSubDir[specie]) href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hname.Data()));
222             //
223             if (!href && hdata->TestBit(AliQAv1::GetClonedBit())) { // try to find the histo for the base name (w/o trigger ending
224               int ind = hname.Index(AliQADataMaker::GetTriggerPrefix());
225               if (ind>0) {
226                 hname.Resize(ind);
227                 if (fRefSubDir)                  href  = static_cast<TH1*>(fRefSubDir->Get(hname.Data())) ;
228                 else if (fRefOCDBSubDir[specie]) href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hname.Data()));
229               }             
230             }
231             //
232             if (!href) 
233               test[specie] = -1 ; // no reference data ; 
234             else {
235               Double_t rv =  DiffK(hdata, href) ;
236               AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ; 
237               test[specie] += rv ; 
238               count[specie]++ ;
239             }
240           } else
241             AliError("Data type cannot be processed") ;
242           if (count[specie] != 0) 
243             test[specie] /= count[specie] ;
244         }
245       }
246     }
247   }
248 }  
249
250
251 //____________________________________________________________________________ 
252 Double_t AliQACheckerBase::DiffC(const TH1 * href, const TH1 * hin) const
253 {
254   // compares two histograms using the Chi2 test
255   if ( hin->Integral() == 0 ) {
256     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s is empty", hin->GetName())) ; 
257     return 0. ;
258   }
259     
260   return hin->Chi2Test(href) ;  
261 }
262
263 //____________________________________________________________________________ 
264 Double_t AliQACheckerBase::DiffK(const TH1 * href, const TH1 * hin) const
265 {
266   // compares two histograms using the Kolmogorov test
267   if ( hin->Integral() == 0 || href->Integral() == 0) {
268     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s or its reference is empty", hin->GetName())) ; 
269     return 0. ;
270   }
271     
272   return hin->KolmogorovTest(href) ;  
273 }
274
275   //_____________________________________________________________________________
276 void AliQACheckerBase::GetRefSubDir(const char * det, const char * task, TDirectory *& dirFile, TObjArray **& dirOCDB)     
277
278     // Opens and returns the file with the reference data 
279   dirFile = NULL ; 
280   TString refStorage(AliQAv1::GetQARefStorage()) ;
281   if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
282     AliError(Form("%s is not a valid location for reference data", refStorage.Data())) ; 
283     return ; 
284   } else {
285     AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::GetTaskIndex(task)) ;
286       //    dirOCDB = new TObjArray*[AliRecoParam::kNSpecies] ; 
287     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
288       dirOCDB[specie] = NULL ; 
289       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
290         continue ; 
291       AliQAv1::SetQARefDataDirName(specie) ;
292       if ( ! manQA->GetLock() ) { 
293         manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ; 
294         manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
295         manQA->SetRun(AliCDBManager::Instance()->GetRun()) ; 
296         manQA->SetLock() ; 
297       }
298       char * detOCDBDir = Form("%s/%s/%s", det, AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
299       AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
300       if (entry) {
301         TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
302         if ( listDetQAD && strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
303           AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), det)) ; 
304           listDetQAD = NULL ; 
305           continue ; 
306         } 
307         if ( listDetQAD ) {
308           TIter next(listDetQAD) ;
309           while ( (TObjArray*)next() ) 
310             dirOCDB[specie] = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", task, AliRecoParam::GetEventSpecieName(specie)))) ;             
311         }
312       }
313     }
314   }
315 }
316
317 //____________________________________________________________________________
318 void AliQACheckerBase::PrintExternParam() 
319 {
320     // Print the list of external parameter list
321   TIter next(fExternParamList) ; 
322   TParameter<double> *pp ; 
323   TString printit("\n") ;
324   while( (pp = (TParameter<double>*)next()) )
325     printit += Form("%s = %f\n", pp->GetName(), pp->GetVal());  
326   AliInfo(Form("%s", printit.Data())) ;
327 }
328   
329 //____________________________________________________________________________
330 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, AliDetectorRecoParam * recoParam) 
331
332   AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
333   
334   Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
335   Check(rv, index, recoParam) ;
336   SetQA(index, rv) ;    
337   
338   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
339   
340   delete [] rv ; 
341   Finish() ; 
342 }
343
344 //____________________________________________________________________________
345 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, TObjArray ** list, AliDetectorRecoParam * recoParam) 
346
347   // RS: perform check for all trigger classes in loop
348   Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
349   //
350   TObjArray ** listTrig = new TObjArray *[AliRecoParam::kNSpecies];
351   //
352   for (int itc=-1;itc<AliQADataMaker::GetNTrigClasses();itc++) {
353     //
354     // RS: fetch the histograms for each specie and for given trigger
355     //AliInfo(Form("Processing %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc))); 
356     
357     for (int specie=0;specie<AliRecoParam::kNSpecies;specie++) {
358       listTrig[specie] = 0;
359       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) || !list[specie]) continue;
360       listTrig[specie] = new TObjArray( list[specie]->GetSize() ); // destination for clones of this trigger
361       AliQADataMaker::GetDataOfTrigClass(list[specie],itc, listTrig[specie]);
362     }
363     AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc))); 
364     Check(rv, index, listTrig, recoParam) ;
365     SetQA(index, rv) ;  
366     AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc)));
367     //
368     for (int specie=0;specie<AliRecoParam::kNSpecies;specie++) if (listTrig[specie]) delete listTrig[specie]; // clean temporary container
369   }
370   delete [] rv ; 
371   delete [] listTrig;
372   Finish() ; 
373 }
374
375 //____________________________________________________________________________
376 void AliQACheckerBase::Finish() const 
377 {
378   // wrap up and save QA in proper file
379   AliQAv1::GetQAResultFile() ; 
380   AliQAv1 * qa = AliQAv1::Instance() ; 
381   qa->Write(AliQAv1::GetQAName(), kWriteDelete) ;   
382 }
383
384 //____________________________________________________________________________ 
385 void AliQACheckerBase::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode) 
386 {
387   // makes the QA image for sim and rec
388   TObjArray tmpArr;  // array to store flat version of original array (which may contain clones)
389   //
390   for (Int_t esIndex = 0; esIndex < AliRecoParam::kNSpecies; esIndex++) {
391     if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0) continue;
392     Int_t nImages = 0;
393     TIter next(list[esIndex]);
394     TObject* hdata = NULL;
395     tmpArr.Clear();
396     while ( (hdata=(next())) ) { // count histos and transfere to flat array
397       if (hdata->InheritsFrom(TH1::Class()) && hdata->TestBit(AliQAv1::GetImageBit()) ) {  // histo, not cloned
398         nImages++; 
399         tmpArr.AddLast(hdata); 
400         continue;
401       }
402       if (!hdata->TestBit(AliQAv1::GetClonedBit())) continue;  // not an array of clones, unknown object
403       TIter nextCl((TObjArray*)hdata);   // array of histo clones
404       TObject* hcl = 0;
405       while ((hcl=nextCl())) if (hcl->InheritsFrom(TH1::Class()) && hcl->TestBit(AliQAv1::GetImageBit())) {tmpArr.AddLast(hcl); nImages++;}
406     }
407     //
408     if ( nImages == 0 ) {
409       AliDebug(AliQAv1::GetQADebugLevel(), Form("No histogram will be plotted for %s %s %s\n", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)));  
410       continue;
411     }
412     AliDebug(AliQAv1::GetQADebugLevel(), Form("%d histograms will be plotted for %s %s %s\n", nImages, GetName(), AliQAv1::GetTaskName(task).Data(),AliRecoParam::GetEventSpecieName(esIndex)));  
413     //        
414     const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)); 
415     //
416     if ( !fImage[esIndex] ) fImage[esIndex] = new TCanvas(title, title);
417     //
418     fImage[esIndex]->Clear(); 
419     fImage[esIndex]->SetTitle(title); 
420     fImage[esIndex]->cd(); 
421     TPaveText someText(0.015, 0.015, 0.98, 0.98);
422     someText.AddText(title);
423     someText.Draw(); 
424     fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat())); 
425     fImage[esIndex]->Clear(); 
426     Int_t nx = TMath::Nint(TMath::Sqrt(nImages));
427     Int_t ny = nx; 
428     if (nx < TMath::Sqrt(nImages)) ny++; 
429     //
430     fImage[esIndex]->Divide(nx, ny); 
431     TIter nexthist(&tmpArr);
432     Int_t npad = 1; 
433     fImage[esIndex]->cd(npad); 
434     TH1* histo = 0;
435     while ( (histo=(TH1*)nexthist()) ) { // tmpArr is guaranteed to contain only plottable histos, no checks needed
436       TString opts = histo->GetDrawOption();
437       if (opts.Contains("logy",TString::kIgnoreCase)) gPad->SetLogy();
438       if (opts.Contains("logx",TString::kIgnoreCase)) gPad->SetLogx();
439       histo->DrawCopy(); 
440       fImage[esIndex]->cd(++npad); 
441     }
442     fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps"); 
443   }
444 }
445
446 //____________________________________________________________________________
447 void AliQACheckerBase::SetHiLo(Float_t * hiValue, Float_t * lowValue) 
448 {
449   AliDebug(AliQAv1::GetQADebugLevel(), "Previous setting was:") ;
450   if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
451     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", 
452                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
453                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
454                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
455                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
456     AliInfo(Form("%s", text)) ; 
457   }
458   
459   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
460     fLowTestValue[index]  = lowValue[index] ; 
461     fUpTestValue[index]   = hiValue[index] ; 
462   }
463   AliDebug(AliQAv1::GetQADebugLevel(), "Current setting is:") ;
464   if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
465     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", 
466                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
467                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
468                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
469                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;     AliInfo(Form("%s", text)) ; 
470   }
471 }
472
473 //____________________________________________________________________________
474 void AliQACheckerBase::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
475 {
476         // sets the QA according the return value of the Check
477
478   AliQAv1 * qa = AliQAv1::Instance(index) ;
479
480   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
481     if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
482       continue ;
483     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
484       qa->Set(AliQAv1::kFATAL, specie) ; 
485     } else {
486       if ( value[specie] >= fLowTestValue[AliQAv1::kFATAL] && value[specie] < fUpTestValue[AliQAv1::kFATAL] ) 
487         qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
488       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
489         qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
490       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
491         qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
492       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
493         qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;   
494     }
495   }
496 }