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