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