ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/


/* $Id$ */

//
//  Base class for detectors quality assurance checkers 
//  Compares Data made by QADataMakers with reference data
//  Y. Schutz CERN August 2007
//

// --- ROOT system ---
#include <TCanvas.h>
#include <TClass.h>
#include <TH1F.h> 
#include <TH1I.h> 
#include <TIterator.h> 
#include <TKey.h> 
#include <TFile.h> 
#include <TList.h>
#include <TNtupleD.h>
#include <TParameter.h>
#include <TPaveText.h>

// --- Standard library ---

// --- AliRoot header files ---
#include "AliCDBEntry.h"
#include "AliLog.h"
#include "AliQAv1.h"
#include "AliQAChecker.h"
#include "AliQACheckerBase.h"
#include "AliQADataMaker.h"
#include "AliQAManager.h"
#include "AliDetectorRecoParam.h"

ClassImp(AliQACheckerBase)

           
//____________________________________________________________________________ 
AliQACheckerBase::AliQACheckerBase(const char * name, const char * title) : 
  TNamed(name, title), 
  fDataSubDir(0x0),
  fRefSubDir(0x0), 
  fRefOCDBSubDir(new TObjArray*[AliRecoParam::kNSpecies]), 
  fLowTestValue(new Float_t[AliQAv1::kNBIT]),
  fUpTestValue(new Float_t[AliQAv1::kNBIT]),
  fImage(new TCanvas*[AliRecoParam::kNSpecies]), 
  fPrintImage(kTRUE), 
  fExternParamList(new TList())
{
  // ctor
  fLowTestValue[AliQAv1::kINFO]    =  0.5   ; 
  fUpTestValue[AliQAv1::kINFO]     = 1.0 ; 
  fLowTestValue[AliQAv1::kWARNING] =  0.002 ; 
  fUpTestValue[AliQAv1::kWARNING]  = 0.5 ; 
  fLowTestValue[AliQAv1::kERROR]   =  0.0   ; 
  fUpTestValue[AliQAv1::kERROR]    = 0.002 ; 
  fLowTestValue[AliQAv1::kFATAL]   = -1.0   ; 
  fUpTestValue[AliQAv1::kFATAL]    = 0.0 ; 
  
  AliDebug(AliQAv1::GetQADebugLevel(), "Default setting is:") ;
  if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
    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", 
                              fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
                              fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
                              fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
                              fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
    AliInfo(Form("%s", text)) ; 
  }
  
  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
    fImage[specie] = NULL ; 
    fRefOCDBSubDir[specie] = NULL ;
  }
}

//____________________________________________________________________________ 
AliQACheckerBase::~AliQACheckerBase()
{
  delete [] fLowTestValue ; 
  delete [] fUpTestValue ; 
  DeleteImages();  
  delete[] fImage ; 
  delete[] fRefOCDBSubDir ; 
  AliQAv1::GetQAResultFile()->Close() ; 
  if (fExternParamList) {
    fExternParamList->Clear() ; 
    delete fExternParamList ; 
  }
}

//____________________________________________________________________________
void AliQACheckerBase::PrivateCheck(Double_t * test, AliQAv1::ALITASK_t index, const AliDetectorRecoParam * recoParam) 
{
  // Performs a basic checking
  // Compares all the histograms stored in the directory
  // With reference histograms either in a file of in OCDB  

  TObjArray ** list = new TObjArray *[AliRecoParam::kNSpecies] ; 
  Int_t specie ;
  for (specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
    list[specie] =  new TObjArray(AliQAv1::GetMaxQAObj()) ; 
    if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
      continue ; 
    if (fDataSubDir) {
      TList * keyList = fDataSubDir->GetListOfKeys() ; 
      TIter next(keyList) ; 
      TKey * key ;
      while ( (key = static_cast<TKey *>(next())) ) {
        TDirectory * specieDir = fDataSubDir->GetDirectory(key->GetName()) ; 
        TList * keykeyList = specieDir->GetListOfKeys() ; 
        TIter next2(keykeyList) ; 
        TKey * keykey ;
        while ( (keykey = static_cast<TKey *>(next2())) ) {
          TObject * odata = specieDir->Get(keykey->GetName()) ; 
          if ( odata->IsA()->InheritsFrom("TH1") ) {
            TH1 * hdata = static_cast<TH1*>(odata) ;
            list[specie]->Add(hdata) ; 
          } else if (!odata->IsA()->InheritsFrom("TDirectory")) // skip the expert directory
            AliError(Form("%s Is a Classname that cannot be processed", key->GetClassName())) ;
        }
      }
    }
  }
 
  Check(test, index, list, recoParam) ;
  
  delete[] list ; 
    
}  

//____________________________________________________________________________
void AliQACheckerBase::Check(Double_t * test, AliQAv1::ALITASK_t task, TObjArray ** list, const AliDetectorRecoParam * /* recoParam */) 
{
  // Performs a basic checking
  // Compares all the histograms in the list

  Int_t count[AliRecoParam::kNSpecies]   = { 0 }; 

  GetRefSubDir(GetName(), AliQAv1::GetTaskName(task), fRefSubDir, fRefOCDBSubDir) ;
 // SetRefandData(refDir, refOCDBDir) ; 
  
  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
    test[specie] = 1.0 ; 
    if ( !AliQAv1::Instance()->IsEventSpecieSet(specie)) continue ; 
    if (list[specie]->GetEntries() == 0)  
      test[specie] = 0. ; // nothing to check
    else {
      if (!fRefSubDir && !fRefOCDBSubDir)
        test[specie] = -1 ; // no reference data
      else {
        TIter next(list[specie]) ; 
        TH1 * hdata ;
        count[specie] = 0 ; 
        while ( (hdata = static_cast<TH1 *>(next())) ) {
          if ( hdata->IsA()->InheritsFrom("TH1") ) {
            if ( hdata->TestBit(AliQAv1::GetExpertBit()) )  // does not perform the test for expert data
              continue ; 
            // 
	    // First try to find the ref histo with exact name (with possible trigger clon ending)
	    TString hname = hdata->GetName();
	    TH1 * href = NULL ; 
            if (fRefSubDir)                  href  = static_cast<TH1*>(fRefSubDir->Get(hname.Data())) ;
            else if (fRefOCDBSubDir[specie]) href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hname.Data()));
	    //
	    if (!href && hdata->TestBit(AliQAv1::GetClonedBit())) { // try to find the histo for the base name (w/o trigger ending
	      int ind = hname.Index(AliQADataMaker::GetTriggerPrefix());
	      if (ind>0) {
		hname.Resize(ind);
		if (fRefSubDir)                  href  = static_cast<TH1*>(fRefSubDir->Get(hname.Data())) ;
		else if (fRefOCDBSubDir[specie]) href  = static_cast<TH1*>(fRefOCDBSubDir[specie]->FindObject(hname.Data()));
	      }		    
	    }
	    //
            if (!href) 
              test[specie] = -1 ; // no reference data ; 
            else {
              Double_t rv =  DiffK(hdata, href) ;
              AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ->Test = %f", hdata->GetName(), rv)) ; 
              test[specie] += rv ; 
              count[specie]++ ;
            }
          } else
            AliError("Data type cannot be processed") ;
          if (count[specie] != 0) 
            test[specie] /= count[specie] ;
        }
      }
    }
  }
}  


//____________________________________________________________________________ 
void AliQACheckerBase::DeleteImages()
{
  // clean images
  for (Int_t esIndex = 0 ; esIndex < AliRecoParam::kNSpecies ; esIndex++) {
    if ( fImage[esIndex] )          {delete fImage[esIndex];          fImage[esIndex] = 0;}
    if ( fRefOCDBSubDir[esIndex] )  {delete fRefOCDBSubDir[esIndex];  fRefOCDBSubDir[esIndex] = 0;}
  }
}

//____________________________________________________________________________ 
Double_t AliQACheckerBase::DiffC(const TH1 * href, const TH1 * hin) const
{
  // compares two histograms using the Chi2 test
  if ( hin->Integral() == 0 ) {
    AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s is empty", hin->GetName())) ; 
    return 0. ;
  }
    
  return hin->Chi2Test(href) ;  
}

//____________________________________________________________________________ 
Double_t AliQACheckerBase::DiffK(const TH1 * href, const TH1 * hin) const
{
  // compares two histograms using the Kolmogorov test
  if ( hin->Integral() == 0 || href->Integral() == 0) {
    AliDebug(AliQAv1::GetQADebugLevel(), Form("Spectrum %s or its reference is empty", hin->GetName())) ; 
    return 0. ;
  }
    
  return hin->KolmogorovTest(href) ;  
}

  //_____________________________________________________________________________
void AliQACheckerBase::GetRefSubDir(const char * det, const char * task, TDirectory *& dirFile, TObjArray **& dirOCDB)     
{ 
    // Opens and returns the file with the reference data 
  dirFile = NULL ; 
  TString refStorage(AliQAv1::GetQARefStorage()) ;
  if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
    AliError(Form("%s is not a valid location for reference data", refStorage.Data())) ; 
    return ; 
  } else {
    AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::GetTaskIndex(task)) ;
      //    dirOCDB = new TObjArray*[AliRecoParam::kNSpecies] ;	
    for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
      dirOCDB[specie] = NULL ; 
      if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
        continue ; 
      AliQAv1::SetQARefDataDirName(specie) ;
      if ( ! manQA->GetLock() ) { 
        manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ; 
        manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
        manQA->SetRun(AliCDBManager::Instance()->GetRun()) ; 
        manQA->SetLock() ; 
      }
      char * detOCDBDir = Form("%s/%s/%s", det, AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
      AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
      if (entry) {
        TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
        if ( listDetQAD && strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
          AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), det)) ; 
          listDetQAD = NULL ; 
          continue ; 
        } 
        if ( listDetQAD ) {
          TIter next(listDetQAD) ;
          while ( (TObjArray*)next() ) 
            dirOCDB[specie] = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", task, AliRecoParam::GetEventSpecieName(specie)))) ;             
        }
      }
    }
  }
}

//____________________________________________________________________________
void AliQACheckerBase::PrintExternParam() 
{
    // Print the list of external parameter list
  TIter next(fExternParamList) ; 
  TParameter<double> *pp ; 
  TString printit("\n") ;
  while( (pp = (TParameter<double>*)next()) )
    printit += Form("%s = %f\n", pp->GetName(), pp->GetVal());  
  AliInfo(Form("%s", printit.Data())) ;
}
  
//____________________________________________________________________________
void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, const AliDetectorRecoParam * recoParam) 
{ 
    //Run the checker for all kind of species
  AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s", AliQAv1::GetAliTaskName(index))) ; 
  
  Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
  for (int i=AliRecoParam::kNSpecies;i--;) rv[i] = 0.0;
  PrivateCheck(rv, index, recoParam) ;
  SetQA(index, rv) ; 	
  
  AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s", AliQAv1::GetAliTaskName(index))) ;
  
  delete [] rv ; 
  Finish() ; 
}

//____________________________________________________________________________
void AliQACheckerBase::Run(AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * recoParam) 
{ 
  // RS: perform check for all trigger classes in loop
  Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
  //
  TObjArray ** listTrig = new TObjArray *[AliRecoParam::kNSpecies];
  //
  for (int itc=-1;itc<AliQADataMaker::GetNTrigClasses();itc++) {
    //
    // RS: fetch the histograms for each specie and for given trigger
    //AliInfo(Form("Processing %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc))); 
    
    for (int specie=0;specie<AliRecoParam::kNSpecies;specie++) {
      listTrig[specie] = 0;
      if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) || !list[specie]) continue;
      listTrig[specie] = new TObjArray( list[specie]->GetSize() ); // destination for clones of this trigger
      AliQADataMaker::GetDataOfTrigClass(list[specie],itc, listTrig[specie]);
    }
    AliDebug(AliQAv1::GetQADebugLevel(), Form("Processing %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc))); 
    Check(rv, index, listTrig, recoParam) ;
    SetQA(index, rv) ; 	
    AliDebug(AliQAv1::GetQADebugLevel(), Form("Test result of %s for trigger: %s", AliQAv1::GetAliTaskName(index),AliQADataMaker::GetTrigClassName(itc)));
    //
    for (int specie=0;specie<AliRecoParam::kNSpecies;specie++) if (listTrig[specie]) delete listTrig[specie]; // clean temporary container
  }
  delete [] rv ; 
  delete [] listTrig;
  Finish() ; 
}

//____________________________________________________________________________
void AliQACheckerBase::Finish() const 
{
  // wrap up and save QA in proper file
  AliQAv1::GetQAResultFile() ; 
  AliQAv1 * qa = AliQAv1::Instance() ; 
  qa->Write(AliQAv1::GetQAName(), kWriteDelete) ;   
}

//____________________________________________________________________________ 
void AliQACheckerBase::MakeImage( TObjArray ** list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode) 
{
  // makes the QA image for sim and rec
  TObjArray tmpArr;  // array to store flat version of original array (which may contain clones)
  //
  for (Int_t esIndex = 0; esIndex < AliRecoParam::kNSpecies; esIndex++) {
    if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(esIndex)) || list[esIndex]->GetEntries() == 0) continue;
    Int_t nImages = 0;
    TIter next(list[esIndex]);
    TObject* hdata = NULL;
    tmpArr.Clear();
    while ( (hdata=(next())) ) { // count histos and transfere to flat array
      if (hdata->InheritsFrom(TH1::Class()) && hdata->TestBit(AliQAv1::GetImageBit()) ) {  // histo, not cloned
	nImages++; 
	tmpArr.AddLast(hdata); 
	continue;
      }
      if (!hdata->TestBit(AliQAv1::GetClonedBit())) continue;  // not an array of clones, unknown object
      TIter nextCl((TObjArray*)hdata);   // array of histo clones
      TObject* hcl = 0;
      while ((hcl=nextCl())) if (hcl->InheritsFrom(TH1::Class()) && hcl->TestBit(AliQAv1::GetImageBit())) {tmpArr.AddLast(hcl); nImages++;}
    }
    //
    if ( nImages == 0 ) {
      AliDebug(AliQAv1::GetQADebugLevel(), Form("No histogram will be plotted for %s %s %s\n", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)));  
      continue;
    }
    AliDebug(AliQAv1::GetQADebugLevel(), Form("%d histograms will be plotted for %s %s %s\n", nImages, GetName(), AliQAv1::GetTaskName(task).Data(),AliRecoParam::GetEventSpecieName(esIndex)));  
    //        
    const Char_t * title = Form("QA_%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(esIndex)); 
    //
    if ( !fImage[esIndex] ) fImage[esIndex] = new TCanvas(title, title);
    //
    fImage[esIndex]->Clear(); 
    fImage[esIndex]->SetTitle(title); 
    fImage[esIndex]->cd(); 
    TPaveText someText(0.015, 0.015, 0.98, 0.98);
    someText.AddText(title);
    someText.Draw(); 
    fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat())); 
    fImage[esIndex]->Clear(); 
    Int_t nx = TMath::Nint(TMath::Sqrt(nImages));
    Int_t ny = nx; 
    if (nx < TMath::Sqrt(nImages)) ny++; 
    //
    fImage[esIndex]->Divide(nx, ny); 
    TIter nexthist(&tmpArr);
    Int_t npad = 1; 
    fImage[esIndex]->cd(npad); 
    TH1* histo = 0;
    while ( (histo=(TH1*)nexthist()) ) { // tmpArr is guaranteed to contain only plottable histos, no checks needed
      TString opts = histo->GetDrawOption();
      if (opts.Contains("logy",TString::kIgnoreCase)) gPad->SetLogy();
      if (opts.Contains("logx",TString::kIgnoreCase)) gPad->SetLogx();
      histo->DrawCopy(); 
      fImage[esIndex]->cd(++npad); 
    }
    fImage[esIndex]->Print(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), AliQAv1::GetModeName(mode), AliQAChecker::Instance()->GetRunNumber(), AliQAv1::GetImageFileFormat()), "ps"); 
  }
}

//____________________________________________________________________________
void AliQACheckerBase::SetHiLo(Float_t * hiValue, Float_t * lowValue) 
{
  AliDebug(AliQAv1::GetQADebugLevel(), "Previous setting was:") ;
  if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
    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", 
                              fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
                              fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
                              fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
                              fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ; 
    AliInfo(Form("%s", text)) ; 
  }
  
  for (Int_t index = 0 ; index < AliQAv1::kNBIT ; index++) {
    fLowTestValue[index]  = lowValue[index] ; 
    fUpTestValue[index]   = hiValue[index] ; 
  }
  AliDebug(AliQAv1::GetQADebugLevel(), "Current setting is:") ;
  if ( AliDebugLevel()  == AliQAv1::GetQADebugLevel() ) {
    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", 
                              fLowTestValue[AliQAv1::kINFO], fUpTestValue[AliQAv1::kINFO], 
                              fLowTestValue[AliQAv1::kWARNING], fUpTestValue[AliQAv1::kWARNING], 
                              fLowTestValue[AliQAv1::kERROR], fUpTestValue[AliQAv1::kERROR], 
                              fLowTestValue[AliQAv1::kFATAL], fUpTestValue[AliQAv1::kFATAL]) ;     AliInfo(Form("%s", text)) ; 
  }
}

//____________________________________________________________________________
void AliQACheckerBase::SetQA(AliQAv1::ALITASK_t index, Double_t * value) const
{
	// sets the QA according the return value of the Check

  AliQAv1 * qa = AliQAv1::Instance(index) ;

  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
    if (! qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
      continue ;
    if (  value == NULL ) { // No checker is implemented, set all QA to Fatal
      qa->Set(AliQAv1::kFATAL, specie) ; 
    } else {
      if ( value[specie] >= fLowTestValue[AliQAv1::kFATAL] && value[specie] < fUpTestValue[AliQAv1::kFATAL] ) 
        qa->Set(AliQAv1::kFATAL, AliRecoParam::ConvertIndex(specie)) ; 
      else if ( value[specie] > fLowTestValue[AliQAv1::kERROR] && value[specie] <= fUpTestValue[AliQAv1::kERROR]  )
        qa->Set(AliQAv1::kERROR, AliRecoParam::ConvertIndex(specie)) ; 
      else if ( value[specie] > fLowTestValue[AliQAv1::kWARNING] && value[specie] <= fUpTestValue[AliQAv1::kWARNING]  )
        qa->Set(AliQAv1::kWARNING, AliRecoParam::ConvertIndex(specie)) ;
      else if ( value[specie] > fLowTestValue[AliQAv1::kINFO] && value[specie] <= fUpTestValue[AliQAv1::kINFO] ) 
        qa->Set(AliQAv1::kINFO, AliRecoParam::ConvertIndex(specie)) ; 	
    }
  }
}
 AliQACheckerBase.cxx:1
 AliQACheckerBase.cxx:2
 AliQACheckerBase.cxx:3
 AliQACheckerBase.cxx:4
 AliQACheckerBase.cxx:5
 AliQACheckerBase.cxx:6
 AliQACheckerBase.cxx:7
 AliQACheckerBase.cxx:8
 AliQACheckerBase.cxx:9
 AliQACheckerBase.cxx:10
 AliQACheckerBase.cxx:11
 AliQACheckerBase.cxx:12
 AliQACheckerBase.cxx:13
 AliQACheckerBase.cxx:14
 AliQACheckerBase.cxx:15
 AliQACheckerBase.cxx:16
 AliQACheckerBase.cxx:17
 AliQACheckerBase.cxx:18
 AliQACheckerBase.cxx:19
 AliQACheckerBase.cxx:20
 AliQACheckerBase.cxx:21
 AliQACheckerBase.cxx:22
 AliQACheckerBase.cxx:23
 AliQACheckerBase.cxx:24
 AliQACheckerBase.cxx:25
 AliQACheckerBase.cxx:26
 AliQACheckerBase.cxx:27
 AliQACheckerBase.cxx:28
 AliQACheckerBase.cxx:29
 AliQACheckerBase.cxx:30
 AliQACheckerBase.cxx:31
 AliQACheckerBase.cxx:32
 AliQACheckerBase.cxx:33
 AliQACheckerBase.cxx:34
 AliQACheckerBase.cxx:35
 AliQACheckerBase.cxx:36
 AliQACheckerBase.cxx:37
 AliQACheckerBase.cxx:38
 AliQACheckerBase.cxx:39
 AliQACheckerBase.cxx:40
 AliQACheckerBase.cxx:41
 AliQACheckerBase.cxx:42
 AliQACheckerBase.cxx:43
 AliQACheckerBase.cxx:44
 AliQACheckerBase.cxx:45
 AliQACheckerBase.cxx:46
 AliQACheckerBase.cxx:47
 AliQACheckerBase.cxx:48
 AliQACheckerBase.cxx:49
 AliQACheckerBase.cxx:50
 AliQACheckerBase.cxx:51
 AliQACheckerBase.cxx:52
 AliQACheckerBase.cxx:53
 AliQACheckerBase.cxx:54
 AliQACheckerBase.cxx:55
 AliQACheckerBase.cxx:56
 AliQACheckerBase.cxx:57
 AliQACheckerBase.cxx:58
 AliQACheckerBase.cxx:59
 AliQACheckerBase.cxx:60
 AliQACheckerBase.cxx:61
 AliQACheckerBase.cxx:62
 AliQACheckerBase.cxx:63
 AliQACheckerBase.cxx:64
 AliQACheckerBase.cxx:65
 AliQACheckerBase.cxx:66
 AliQACheckerBase.cxx:67
 AliQACheckerBase.cxx:68
 AliQACheckerBase.cxx:69
 AliQACheckerBase.cxx:70
 AliQACheckerBase.cxx:71
 AliQACheckerBase.cxx:72
 AliQACheckerBase.cxx:73
 AliQACheckerBase.cxx:74
 AliQACheckerBase.cxx:75
 AliQACheckerBase.cxx:76
 AliQACheckerBase.cxx:77
 AliQACheckerBase.cxx:78
 AliQACheckerBase.cxx:79
 AliQACheckerBase.cxx:80
 AliQACheckerBase.cxx:81
 AliQACheckerBase.cxx:82
 AliQACheckerBase.cxx:83
 AliQACheckerBase.cxx:84
 AliQACheckerBase.cxx:85
 AliQACheckerBase.cxx:86
 AliQACheckerBase.cxx:87
 AliQACheckerBase.cxx:88
 AliQACheckerBase.cxx:89
 AliQACheckerBase.cxx:90
 AliQACheckerBase.cxx:91
 AliQACheckerBase.cxx:92
 AliQACheckerBase.cxx:93
 AliQACheckerBase.cxx:94
 AliQACheckerBase.cxx:95
 AliQACheckerBase.cxx:96
 AliQACheckerBase.cxx:97
 AliQACheckerBase.cxx:98
 AliQACheckerBase.cxx:99
 AliQACheckerBase.cxx:100
 AliQACheckerBase.cxx:101
 AliQACheckerBase.cxx:102
 AliQACheckerBase.cxx:103
 AliQACheckerBase.cxx:104
 AliQACheckerBase.cxx:105
 AliQACheckerBase.cxx:106
 AliQACheckerBase.cxx:107
 AliQACheckerBase.cxx:108
 AliQACheckerBase.cxx:109
 AliQACheckerBase.cxx:110
 AliQACheckerBase.cxx:111
 AliQACheckerBase.cxx:112
 AliQACheckerBase.cxx:113
 AliQACheckerBase.cxx:114
 AliQACheckerBase.cxx:115
 AliQACheckerBase.cxx:116
 AliQACheckerBase.cxx:117
 AliQACheckerBase.cxx:118
 AliQACheckerBase.cxx:119
 AliQACheckerBase.cxx:120
 AliQACheckerBase.cxx:121
 AliQACheckerBase.cxx:122
 AliQACheckerBase.cxx:123
 AliQACheckerBase.cxx:124
 AliQACheckerBase.cxx:125
 AliQACheckerBase.cxx:126
 AliQACheckerBase.cxx:127
 AliQACheckerBase.cxx:128
 AliQACheckerBase.cxx:129
 AliQACheckerBase.cxx:130
 AliQACheckerBase.cxx:131
 AliQACheckerBase.cxx:132
 AliQACheckerBase.cxx:133
 AliQACheckerBase.cxx:134
 AliQACheckerBase.cxx:135
 AliQACheckerBase.cxx:136
 AliQACheckerBase.cxx:137
 AliQACheckerBase.cxx:138
 AliQACheckerBase.cxx:139
 AliQACheckerBase.cxx:140
 AliQACheckerBase.cxx:141
 AliQACheckerBase.cxx:142
 AliQACheckerBase.cxx:143
 AliQACheckerBase.cxx:144
 AliQACheckerBase.cxx:145
 AliQACheckerBase.cxx:146
 AliQACheckerBase.cxx:147
 AliQACheckerBase.cxx:148
 AliQACheckerBase.cxx:149
 AliQACheckerBase.cxx:150
 AliQACheckerBase.cxx:151
 AliQACheckerBase.cxx:152
 AliQACheckerBase.cxx:153
 AliQACheckerBase.cxx:154
 AliQACheckerBase.cxx:155
 AliQACheckerBase.cxx:156
 AliQACheckerBase.cxx:157
 AliQACheckerBase.cxx:158
 AliQACheckerBase.cxx:159
 AliQACheckerBase.cxx:160
 AliQACheckerBase.cxx:161
 AliQACheckerBase.cxx:162
 AliQACheckerBase.cxx:163
 AliQACheckerBase.cxx:164
 AliQACheckerBase.cxx:165
 AliQACheckerBase.cxx:166
 AliQACheckerBase.cxx:167
 AliQACheckerBase.cxx:168
 AliQACheckerBase.cxx:169
 AliQACheckerBase.cxx:170
 AliQACheckerBase.cxx:171
 AliQACheckerBase.cxx:172
 AliQACheckerBase.cxx:173
 AliQACheckerBase.cxx:174
 AliQACheckerBase.cxx:175
 AliQACheckerBase.cxx:176
 AliQACheckerBase.cxx:177
 AliQACheckerBase.cxx:178
 AliQACheckerBase.cxx:179
 AliQACheckerBase.cxx:180
 AliQACheckerBase.cxx:181
 AliQACheckerBase.cxx:182
 AliQACheckerBase.cxx:183
 AliQACheckerBase.cxx:184
 AliQACheckerBase.cxx:185
 AliQACheckerBase.cxx:186
 AliQACheckerBase.cxx:187
 AliQACheckerBase.cxx:188
 AliQACheckerBase.cxx:189
 AliQACheckerBase.cxx:190
 AliQACheckerBase.cxx:191
 AliQACheckerBase.cxx:192
 AliQACheckerBase.cxx:193
 AliQACheckerBase.cxx:194
 AliQACheckerBase.cxx:195
 AliQACheckerBase.cxx:196
 AliQACheckerBase.cxx:197
 AliQACheckerBase.cxx:198
 AliQACheckerBase.cxx:199
 AliQACheckerBase.cxx:200
 AliQACheckerBase.cxx:201
 AliQACheckerBase.cxx:202
 AliQACheckerBase.cxx:203
 AliQACheckerBase.cxx:204
 AliQACheckerBase.cxx:205
 AliQACheckerBase.cxx:206
 AliQACheckerBase.cxx:207
 AliQACheckerBase.cxx:208
 AliQACheckerBase.cxx:209
 AliQACheckerBase.cxx:210
 AliQACheckerBase.cxx:211
 AliQACheckerBase.cxx:212
 AliQACheckerBase.cxx:213
 AliQACheckerBase.cxx:214
 AliQACheckerBase.cxx:215
 AliQACheckerBase.cxx:216
 AliQACheckerBase.cxx:217
 AliQACheckerBase.cxx:218
 AliQACheckerBase.cxx:219
 AliQACheckerBase.cxx:220
 AliQACheckerBase.cxx:221
 AliQACheckerBase.cxx:222
 AliQACheckerBase.cxx:223
 AliQACheckerBase.cxx:224
 AliQACheckerBase.cxx:225
 AliQACheckerBase.cxx:226
 AliQACheckerBase.cxx:227
 AliQACheckerBase.cxx:228
 AliQACheckerBase.cxx:229
 AliQACheckerBase.cxx:230
 AliQACheckerBase.cxx:231
 AliQACheckerBase.cxx:232
 AliQACheckerBase.cxx:233
 AliQACheckerBase.cxx:234
 AliQACheckerBase.cxx:235
 AliQACheckerBase.cxx:236
 AliQACheckerBase.cxx:237
 AliQACheckerBase.cxx:238
 AliQACheckerBase.cxx:239
 AliQACheckerBase.cxx:240
 AliQACheckerBase.cxx:241
 AliQACheckerBase.cxx:242
 AliQACheckerBase.cxx:243
 AliQACheckerBase.cxx:244
 AliQACheckerBase.cxx:245
 AliQACheckerBase.cxx:246
 AliQACheckerBase.cxx:247
 AliQACheckerBase.cxx:248
 AliQACheckerBase.cxx:249
 AliQACheckerBase.cxx:250
 AliQACheckerBase.cxx:251
 AliQACheckerBase.cxx:252
 AliQACheckerBase.cxx:253
 AliQACheckerBase.cxx:254
 AliQACheckerBase.cxx:255
 AliQACheckerBase.cxx:256
 AliQACheckerBase.cxx:257
 AliQACheckerBase.cxx:258
 AliQACheckerBase.cxx:259
 AliQACheckerBase.cxx:260
 AliQACheckerBase.cxx:261
 AliQACheckerBase.cxx:262
 AliQACheckerBase.cxx:263
 AliQACheckerBase.cxx:264
 AliQACheckerBase.cxx:265
 AliQACheckerBase.cxx:266
 AliQACheckerBase.cxx:267
 AliQACheckerBase.cxx:268
 AliQACheckerBase.cxx:269
 AliQACheckerBase.cxx:270
 AliQACheckerBase.cxx:271
 AliQACheckerBase.cxx:272
 AliQACheckerBase.cxx:273
 AliQACheckerBase.cxx:274
 AliQACheckerBase.cxx:275
 AliQACheckerBase.cxx:276
 AliQACheckerBase.cxx:277
 AliQACheckerBase.cxx:278
 AliQACheckerBase.cxx:279
 AliQACheckerBase.cxx:280
 AliQACheckerBase.cxx:281
 AliQACheckerBase.cxx:282
 AliQACheckerBase.cxx:283
 AliQACheckerBase.cxx:284
 AliQACheckerBase.cxx:285
 AliQACheckerBase.cxx:286
 AliQACheckerBase.cxx:287
 AliQACheckerBase.cxx:288
 AliQACheckerBase.cxx:289
 AliQACheckerBase.cxx:290
 AliQACheckerBase.cxx:291
 AliQACheckerBase.cxx:292
 AliQACheckerBase.cxx:293
 AliQACheckerBase.cxx:294
 AliQACheckerBase.cxx:295
 AliQACheckerBase.cxx:296
 AliQACheckerBase.cxx:297
 AliQACheckerBase.cxx:298
 AliQACheckerBase.cxx:299
 AliQACheckerBase.cxx:300
 AliQACheckerBase.cxx:301
 AliQACheckerBase.cxx:302
 AliQACheckerBase.cxx:303
 AliQACheckerBase.cxx:304
 AliQACheckerBase.cxx:305
 AliQACheckerBase.cxx:306
 AliQACheckerBase.cxx:307
 AliQACheckerBase.cxx:308
 AliQACheckerBase.cxx:309
 AliQACheckerBase.cxx:310
 AliQACheckerBase.cxx:311
 AliQACheckerBase.cxx:312
 AliQACheckerBase.cxx:313
 AliQACheckerBase.cxx:314
 AliQACheckerBase.cxx:315
 AliQACheckerBase.cxx:316
 AliQACheckerBase.cxx:317
 AliQACheckerBase.cxx:318
 AliQACheckerBase.cxx:319
 AliQACheckerBase.cxx:320
 AliQACheckerBase.cxx:321
 AliQACheckerBase.cxx:322
 AliQACheckerBase.cxx:323
 AliQACheckerBase.cxx:324
 AliQACheckerBase.cxx:325
 AliQACheckerBase.cxx:326
 AliQACheckerBase.cxx:327
 AliQACheckerBase.cxx:328
 AliQACheckerBase.cxx:329
 AliQACheckerBase.cxx:330
 AliQACheckerBase.cxx:331
 AliQACheckerBase.cxx:332
 AliQACheckerBase.cxx:333
 AliQACheckerBase.cxx:334
 AliQACheckerBase.cxx:335
 AliQACheckerBase.cxx:336
 AliQACheckerBase.cxx:337
 AliQACheckerBase.cxx:338
 AliQACheckerBase.cxx:339
 AliQACheckerBase.cxx:340
 AliQACheckerBase.cxx:341
 AliQACheckerBase.cxx:342
 AliQACheckerBase.cxx:343
 AliQACheckerBase.cxx:344
 AliQACheckerBase.cxx:345
 AliQACheckerBase.cxx:346
 AliQACheckerBase.cxx:347
 AliQACheckerBase.cxx:348
 AliQACheckerBase.cxx:349
 AliQACheckerBase.cxx:350
 AliQACheckerBase.cxx:351
 AliQACheckerBase.cxx:352
 AliQACheckerBase.cxx:353
 AliQACheckerBase.cxx:354
 AliQACheckerBase.cxx:355
 AliQACheckerBase.cxx:356
 AliQACheckerBase.cxx:357
 AliQACheckerBase.cxx:358
 AliQACheckerBase.cxx:359
 AliQACheckerBase.cxx:360
 AliQACheckerBase.cxx:361
 AliQACheckerBase.cxx:362
 AliQACheckerBase.cxx:363
 AliQACheckerBase.cxx:364
 AliQACheckerBase.cxx:365
 AliQACheckerBase.cxx:366
 AliQACheckerBase.cxx:367
 AliQACheckerBase.cxx:368
 AliQACheckerBase.cxx:369
 AliQACheckerBase.cxx:370
 AliQACheckerBase.cxx:371
 AliQACheckerBase.cxx:372
 AliQACheckerBase.cxx:373
 AliQACheckerBase.cxx:374
 AliQACheckerBase.cxx:375
 AliQACheckerBase.cxx:376
 AliQACheckerBase.cxx:377
 AliQACheckerBase.cxx:378
 AliQACheckerBase.cxx:379
 AliQACheckerBase.cxx:380
 AliQACheckerBase.cxx:381
 AliQACheckerBase.cxx:382
 AliQACheckerBase.cxx:383
 AliQACheckerBase.cxx:384
 AliQACheckerBase.cxx:385
 AliQACheckerBase.cxx:386
 AliQACheckerBase.cxx:387
 AliQACheckerBase.cxx:388
 AliQACheckerBase.cxx:389
 AliQACheckerBase.cxx:390
 AliQACheckerBase.cxx:391
 AliQACheckerBase.cxx:392
 AliQACheckerBase.cxx:393
 AliQACheckerBase.cxx:394
 AliQACheckerBase.cxx:395
 AliQACheckerBase.cxx:396
 AliQACheckerBase.cxx:397
 AliQACheckerBase.cxx:398
 AliQACheckerBase.cxx:399
 AliQACheckerBase.cxx:400
 AliQACheckerBase.cxx:401
 AliQACheckerBase.cxx:402
 AliQACheckerBase.cxx:403
 AliQACheckerBase.cxx:404
 AliQACheckerBase.cxx:405
 AliQACheckerBase.cxx:406
 AliQACheckerBase.cxx:407
 AliQACheckerBase.cxx:408
 AliQACheckerBase.cxx:409
 AliQACheckerBase.cxx:410
 AliQACheckerBase.cxx:411
 AliQACheckerBase.cxx:412
 AliQACheckerBase.cxx:413
 AliQACheckerBase.cxx:414
 AliQACheckerBase.cxx:415
 AliQACheckerBase.cxx:416
 AliQACheckerBase.cxx:417
 AliQACheckerBase.cxx:418
 AliQACheckerBase.cxx:419
 AliQACheckerBase.cxx:420
 AliQACheckerBase.cxx:421
 AliQACheckerBase.cxx:422
 AliQACheckerBase.cxx:423
 AliQACheckerBase.cxx:424
 AliQACheckerBase.cxx:425
 AliQACheckerBase.cxx:426
 AliQACheckerBase.cxx:427
 AliQACheckerBase.cxx:428
 AliQACheckerBase.cxx:429
 AliQACheckerBase.cxx:430
 AliQACheckerBase.cxx:431
 AliQACheckerBase.cxx:432
 AliQACheckerBase.cxx:433
 AliQACheckerBase.cxx:434
 AliQACheckerBase.cxx:435
 AliQACheckerBase.cxx:436
 AliQACheckerBase.cxx:437
 AliQACheckerBase.cxx:438
 AliQACheckerBase.cxx:439
 AliQACheckerBase.cxx:440
 AliQACheckerBase.cxx:441
 AliQACheckerBase.cxx:442
 AliQACheckerBase.cxx:443
 AliQACheckerBase.cxx:444
 AliQACheckerBase.cxx:445
 AliQACheckerBase.cxx:446
 AliQACheckerBase.cxx:447
 AliQACheckerBase.cxx:448
 AliQACheckerBase.cxx:449
 AliQACheckerBase.cxx:450
 AliQACheckerBase.cxx:451
 AliQACheckerBase.cxx:452
 AliQACheckerBase.cxx:453
 AliQACheckerBase.cxx:454
 AliQACheckerBase.cxx:455
 AliQACheckerBase.cxx:456
 AliQACheckerBase.cxx:457
 AliQACheckerBase.cxx:458
 AliQACheckerBase.cxx:459
 AliQACheckerBase.cxx:460
 AliQACheckerBase.cxx:461
 AliQACheckerBase.cxx:462
 AliQACheckerBase.cxx:463
 AliQACheckerBase.cxx:464
 AliQACheckerBase.cxx:465