]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQACheckerBase.cxx
Coverity fixes.
[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 <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)) 
203       continue ; 
204     if (list[specie]->GetEntries() == 0)  
205       test[specie] = 0. ; // nothing to check
206     else {
207       if (!fRefSubDir && !fRefOCDBSubDir)
208         test[specie] = -1 ; // no reference data
209       else {
210         TIter next(list[specie]) ; 
211         TH1 * hdata ;
212         count[specie] = 0 ; 
213         while ( (hdata = static_cast<TH1 *>(next())) ) {
214           if ( hdata->IsA()->InheritsFrom("TH1") ) {
215             if ( hdata->TestBit(AliQAv1::GetExpertBit()) )  // does not perform the test for expert data
216               continue ; 
217             TH1 * href = NULL ; 
218             if (fRefSubDir) 
219               href  = static_cast<TH1*>(fRefSubDir->Get(hdata->GetName())) ;
220             else if (fRefOCDBSubDir[specie])
221               href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hdata->GetName())) ;
222             if (!href) 
223               test[specie] = -1 ; // no reference data ; 
224             else {
225               Double_t rv =  DiffK(hdata, href) ;
226               AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ; 
227               test[specie] += rv ; 
228               count[specie]++ ;
229             }
230           } else
231             AliError("Data type cannot be processed") ;
232           if (count[specie] != 0) 
233             test[specie] /= count[specie] ;
234         }
235       }
236     }
237   }
238 }  
239
240
241 //____________________________________________________________________________ 
242 Double_t AliQACheckerBase::DiffC(const TH1 * href, const TH1 * hin) const
243 {
244   // compares two histograms using the Chi2 test
245   if ( hin->Integral() == 0 ) {
246     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s is empty", hin->GetName())) ; 
247     return 0. ;
248   }
249     
250   return hin->Chi2Test(href) ;  
251 }
252
253 //____________________________________________________________________________ 
254 Double_t AliQACheckerBase::DiffK(const TH1 * href, const TH1 * hin) const
255 {
256   // compares two histograms using the Kolmogorov test
257   if ( hin->Integral() == 0 || href->Integral() == 0) {
258     AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s or its reference is empty", hin->GetName())) ; 
259     return 0. ;
260   }
261     
262   return hin->KolmogorovTest(href) ;  
263 }
264
265   //_____________________________________________________________________________
266 void AliQACheckerBase::GetRefSubDir(const char * det, const char * task, TDirectory *& dirFile, TObjArray **& dirOCDB)     
267
268     // Opens and returns the file with the reference data 
269   dirFile = NULL ; 
270   TString refStorage(AliQAv1::GetQARefStorage()) ;
271   if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
272     AliError(Form("%s is not a valid location for reference data", refStorage.Data())) ; 
273     return ; 
274   } else {
275     AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::GetTaskIndex(task)) ;
276       //    dirOCDB = new TObjArray*[AliRecoParam::kNSpecies] ; 
277     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
278       dirOCDB[specie] = NULL ; 
279       if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
280         continue ; 
281       AliQAv1::SetQARefDataDirName(specie) ;
282       if ( ! manQA->GetLock() ) { 
283         manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ; 
284         manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
285         manQA->SetRun(AliCDBManager::Instance()->GetRun()) ; 
286         manQA->SetLock() ; 
287       }
288       char * detOCDBDir = Form("%s/%s/%s", det, AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
289       AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
290       if (entry) {
291         TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
292         if ( strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
293           AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), det)) ; 
294           listDetQAD = NULL ; 
295           continue ; 
296         } 
297         if ( listDetQAD ) {
298           TIter next(listDetQAD) ;
299           TObjArray * ar ; 
300           while ( (ar = (TObjArray*)next()) ) 
301             dirOCDB[specie] = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", task, AliRecoParam::GetEventSpecieName(specie)))) ;             
302         }
303       }
304     }
305   }
306 }
307
308 //____________________________________________________________________________
309 void AliQACheckerBase::PrintExternParam() 
310 {
311     // Print the list of external parameter list
312   TIter next(fExternParamList) ; 
313   TParameter<double> *pp ; 
314   TString printit("\n") ;
315   while( (pp = (TParameter<double>*)next()) )
316     printit += Form("%s = %f\n", pp->GetName(), pp->GetVal());  
317   AliInfo(Form("%s", printit.Data())) ;
318 }
319   
320 //____________________________________________________________________________
321 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, AliDetectorRecoParam * recoParam) 
322
323         AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
324   
325         Double_t * rv = NULL ;
326   Check(rv, index, recoParam) ;
327         SetQA(index, rv) ;      
328         
329   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
330         
331   if (rv) 
332     delete [] rv ; 
333   Finish() ; 
334 }
335
336 //____________________________________________________________________________
337 void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, TObjArray ** list, AliDetectorRecoParam * recoParam) 
338
339         AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
340   
341         Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
342   Check(rv, index, list, recoParam) ;
343         SetQA(index, rv) ;      
344         
345   AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
346         
347   delete [] rv ; 
348   Finish() ; 
349 }
350
351 //____________________________________________________________________________
352 void AliQACheckerBase::Finish() const 
353 {
354         // wrap up and save QA in proper file
355         AliQAv1::GetQAResultFile() ; 
356         AliQAv1 * qa = AliQAv1::Instance() ; 
357         qa->Write(AliQAv1::GetQAName(), kWriteDelete) ;   
358 }
359
360 //____________________________________________________________________________ 
361 void AliQACheckerBase::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode) 
362 {
363   // makes the QA image for sim and rec
364   Int_t nImages = 0 ;
365   for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
366     if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) ) 
367       continue ;
368     TIter next(list[esIndex]) ;  
369     TH1 * hdata = NULL ; 
370     while ( (hdata=static_cast<TH1 *>(next())) ) {
371       TString cln(hdata->ClassName()) ; 
372       if ( ! cln.Contains("TH") )
373         continue ; 
374       if ( hdata->TestBit(AliQAv1::GetImageBit()) )
375         nImages++; 
376     }
377     break ; 
378   }
379   if ( nImages == 0 ) {
380     AliDebug(AliQAv1::GetQADebugLevel(), Form("No histogram will be plotted for %s %s\n", GetName(), AliQAv1::GetTaskName(task).Data())) ;  
381   } else {
382     AliDebug(AliQAv1::GetQADebugLevel(), Form("%d histograms will be plotted for %s %s\n", nImages, GetName(), AliQAv1::GetTaskName(task).Data())) ;  
383     for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
384       if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0) 
385         continue ;
386       const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)) ; 
387       if ( !fImage[esIndex] ) {
388         fImage[esIndex] = new TCanvas(title, title) ;
389       }
390       fImage[esIndex]->Clear() ; 
391       fImage[esIndex]->SetTitle(title) ; 
392       fImage[esIndex]->cd() ; 
393       TPaveText someText(0.015, 0.015, 0.98, 0.98) ;
394       someText.AddText(title) ;
395       someText.Draw() ; 
396       fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ; 
397       fImage[esIndex]->Clear() ; 
398       Int_t nx = TMath::Nint(TMath::Sqrt(nImages));
399       Int_t ny = nx  ; 
400       if (nx < TMath::Sqrt(nImages))
401         ny++ ; 
402       
403       fImage[esIndex]->Divide(nx, ny) ; 
404       TIter nexthist(list[esIndex]) ; 
405       TH1* hist = NULL ;
406       Int_t npad = 1 ; 
407       fImage[esIndex]->cd(npad) ; 
408       while ( (hist=static_cast<TH1*>(nexthist())) ) {
409         TString cln(hist->ClassName()) ; 
410         if ( ! cln.Contains("TH") )
411           continue ; 
412         if(hist->TestBit(AliQAv1::GetImageBit())) {
413           TString opts = hist->GetDrawOption();
414           if (opts.Contains("logy",TString::kIgnoreCase)) {
415             gPad->SetLogy();
416             opts.ReplaceAll("logy", "");
417           }
418           if (opts.Contains("logx", TString::kIgnoreCase)) {
419             gPad->SetLogx();
420             opts.ReplaceAll("logx", "");
421           }
422           hist->DrawCopy() ; 
423           fImage[esIndex]->cd(++npad) ; 
424         }
425       }
426       fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps") ; 
427     }
428   }  
429 }
430
431 //____________________________________________________________________________
432 void AliQACheckerBase::SetHiLo(Float_t * hiValue, Float_t * lowValue) 
433 {
434   AliDebug(AliQAv1::GetQADebugLevel(), "Previous setting was:") ;
435   if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
436     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", 
437                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
438                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
439                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
440                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
441     AliInfo(Form("%s", text)) ; 
442   }
443   
444   for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
445     fLowTestValue[index]  = lowValue[index] ; 
446     fUpTestValue[index]   = hiValue[index] ; 
447   }
448   AliDebug(AliQAv1::GetQADebugLevel(), "Current setting is:") ;
449   if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
450     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", 
451                               fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
452                               fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
453                               fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
454                               fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;     AliInfo(Form("%s", text)) ; 
455   }
456 }
457
458 //____________________________________________________________________________
459 void AliQACheckerBase::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
460 {
461         // sets the QA according the return value of the Check
462
463   AliQAv1 * qa = AliQAv1::Instance(index) ;
464
465   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
466     if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
467       continue ;
468     if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
469       qa->Set(AliQAv1::kFATAL, specie) ; 
470     } else {
471       if ( value[specie] >= fLowTestValue[AliQAv1::kFATAL] && value[specie] < fUpTestValue[AliQAv1::kFATAL] ) 
472         qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
473       else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
474         qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
475       else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
476         qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
477       else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
478         qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ;   
479     }
480   }
481 }