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 notifce   *
 * 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
//  Produces the data needed to calculate the quality assurance for Reconstruction
//  All data must be mergeable objects.
//  Y. Schutz CERN July 2007
//

// --- ROOT system ---
#include <TFile.h>
#include <TTree.h>
#include <TNtupleD.h>
#include <TObjArray.h>

// --- Standard library ---

// --- AliRoot header files ---
#include "AliCDBPath.h"
#include "AliCDBEntry.h"
#include "AliDetectorRecoParam.h"
#include "AliCDBManager.h"

#include "AliLog.h"
#include "AliQADataMakerRec.h"
#include "AliQAManager.h"
#include "AliESDEvent.h"
#include "AliRawReader.h"

ClassImp(AliQADataMakerRec)
             
//____________________________________________________________________________ 
AliQADataMakerRec::AliQADataMakerRec(const char * name, const char * title) : 
AliQADataMaker(name, title), 
  fDigitsQAList(NULL),
  fESDsQAList(NULL), 
  fRawsQAList(NULL), 
  fRecPointsQAList(NULL),
  fCorrNt(NULL), 
  fRecoParam(NULL),
  fRecPointsArray(NULL)
{
  // ctor
  fDetectorDirName = GetName() ; 
}

//____________________________________________________________________________ 
AliQADataMakerRec::AliQADataMakerRec(const AliQADataMakerRec& qadm) :
  AliQADataMaker(qadm.GetName(), qadm.GetTitle()), 
  fDigitsQAList(qadm.fDigitsQAList),
  fESDsQAList(qadm.fESDsQAList),
  fRawsQAList(qadm.fRawsQAList),
  fRecPointsQAList(qadm.fRecPointsQAList),
  fCorrNt(qadm.fCorrNt),  
  fRecoParam(qadm.fRecoParam),
  fRecPointsArray(NULL)
{
  //copy ctor
  fDetectorDirName = GetName() ; 
}

//____________________________________________________________________________ 
AliQADataMakerRec::~AliQADataMakerRec()
{
  //dtor: delete the TObjArray and thei content
  if ( fESDsQAList ) {
    for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
      if ( fESDsQAList[specie] ) {
	fESDsQAList[specie]->Delete() ;     
      }
    }
    delete[] fESDsQAList ;
  }
  if ( fRawsQAList ) {
    for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
      if ( fRawsQAList[specie] ) {
        fRawsQAList[specie]->Delete() ;
      }
    }
    delete[] fRawsQAList ;
  }
  if ( fDigitsQAList ) {
    for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
      if ( fDigitsQAList[specie] ) {
        fDigitsQAList[specie]->Delete() ;
      }
    }
    delete[] fDigitsQAList ; 
  }
  if ( fRecPointsQAList ) {
    for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
      if ( fRecPointsQAList[specie] ) {
        fRecPointsQAList[specie]->Delete() ;
      }
    }
    delete[] fRecPointsQAList ; 
  }
  if (fRecPointsArray) {
    fRecPointsArray->Clear() ; 
    delete fRecPointsArray ; 
  }
}

//__________________________________________________________________
AliQADataMakerRec& AliQADataMakerRec::operator = (const AliQADataMakerRec& qadm )
{
  // Assignment operator.
  this->~AliQADataMakerRec();
  new(this) AliQADataMakerRec(qadm);
  return *this;
}

//____________________________________________________________________________
void AliQADataMakerRec::EndOfCycle() 
{
  // Finishes a cycle of QA for all the tasks
  EndOfCycle(AliQAv1::kRAWS) ; 
  EndOfCycle(AliQAv1::kDIGITSR) ; 
  EndOfCycle(AliQAv1::kRECPOINTS) ; 
  EndOfCycle(AliQAv1::kESDS) ; 
  ResetCycle() ; 
}

//____________________________________________________________________________
void AliQADataMakerRec::EndOfCycle(AliQAv1::TASKINDEX_t task) 
{
  // Finishes a cycle of QA 
	
    
  TObjArray ** list = NULL ; 
	
  if ( task == AliQAv1::kRAWS )     
    list = fRawsQAList ; 
  else if ( task == AliQAv1::kDIGITSR ) 
    list = fDigitsQAList ; 
  else if ( task == AliQAv1::kRECPOINTS ) 
    list = fRecPointsQAList ; 
  else if ( task == AliQAv1::kESDS )
    list = fESDsQAList ; 

 
  if ( ! list && ! fCorrNt ) 
    return ; 
  //DefaultEndOfDetectorCycle(task) ;
  EndOfDetectorCycle(task, list) ;
  
  if (! AliQAManager::QAManager(AliQAv1::kRECMODE)->IsSaveData()) return; 

  fDetectorDir = fOutput->GetDirectory(GetDetectorDirName()) ; 
  if (!fDetectorDir) fDetectorDir = fOutput->mkdir(GetDetectorDirName()) ; 
  TDirectory * subDir = fDetectorDir->GetDirectory(AliQAv1::GetTaskName(task)) ; 
  if (!subDir) subDir = fDetectorDir->mkdir(AliQAv1::GetTaskName(task)) ;  
  subDir->cd();
  // 
  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) { // skip Default
    if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)) || AliRecoParam::ConvertIndex(specie) == AliRecoParam::kDefault) continue ; 
    TDirectory * eventSpecieDir = subDir->GetDirectory(AliRecoParam::GetEventSpecieName(specie)) ;
    if (!eventSpecieDir) eventSpecieDir = subDir->mkdir(AliRecoParam::GetEventSpecieName(specie)) ; 
    eventSpecieDir->cd();    
    if (list) {
      if (list[specie]) {
        TIter next(list[specie]) ; 
        TObject * obj ; 
        while( (obj = next()) ) {
          if (!obj->TestBit(AliQAv1::GetExpertBit())) obj->Write(); // RS: Note, this can be also TObjArray of clones
        }
        if (WriteExpert()) {
          TDirectory * expertDir = eventSpecieDir->GetDirectory(AliQAv1::GetExpert()) ; 
          if (!expertDir)
            expertDir = eventSpecieDir->mkdir(AliQAv1::GetExpert()) ; 
          expertDir->cd() ;
          next.Reset() ; 
          while( (obj = next()) ) {
            if (!obj->TestBit(AliQAv1::GetExpertBit())) continue; 
            obj->Write();  // RS: Note, this can be also TObjArray of clones
          }      
        }
      }
    } 
    else if ( fCorrNt ) {
      if (fCorrNt[specie] && AliQAv1::GetDetIndex(GetName()) == AliQAv1::kCORR) {
        if (fCorrNt[specie]->GetNvar() != 0) {
          eventSpecieDir->cd() ; 
          fCorrNt[specie]->Write() ; 
        }
      }
      fOutput->Save() ; 
    }
  }
}

//____________________________________________________________________________
void AliQADataMakerRec::Exec(AliQAv1::TASKINDEX_t task, TObject * data) 
{ 
  // creates the quality assurance data for the various tasks (Hits, SDigits, Digits, ESDs)
	
  if ( task == AliQAv1::kRAWS ) {
    AliDebug(AliQAv1::GetQADebugLevel(), "Processing Raws QA") ; 
    AliRawReader * rawReader = static_cast<AliRawReader *>(data) ; 
    if (rawReader) 
      MakeRaws(rawReader) ;
    else
      AliDebug(AliQAv1::GetQADebugLevel(), "Raw data are not processed") ;     
  } else if ( task == AliQAv1::kDIGITSR ) {
    AliDebug(AliQAv1::GetQADebugLevel(), "Processing Digits QA") ; 
    TTree * tree = static_cast<TTree *>(data) ; 
    if (strcmp(tree->ClassName(), "TTree") == 0) {
      MakeDigits(tree) ; 
    } else {
      AliWarning("data are not a TTree") ; 
    }
  } else if ( task == AliQAv1::kRECPOINTS ) {
    AliDebug(AliQAv1::GetQADebugLevel(), "Processing RecPoints QA") ; 
    TTree * tree = static_cast<TTree *>(data) ; 
    if (strcmp(tree->ClassName(), "TTree") == 0) {
      MakeRecPoints(tree) ; 
    } else {
      AliWarning("data are not a TTree") ; 
    }
  } else if ( task == AliQAv1::kESDS ) {
    AliDebug(AliQAv1::GetQADebugLevel(), "Processing ESDs QA") ; 
    AliESDEvent * esd = static_cast<AliESDEvent *>(data) ; 
    if (strcmp(esd->ClassName(), "AliESDEvent") == 0) 
      MakeESDs(esd) ;
    else 
      AliError("Wrong type of esd container") ; 
  }  
}

//____________________________________________________________________________ 
TObjArray **  AliQADataMakerRec::Init(AliQAv1::TASKINDEX_t task, Int_t cycles)
{
  // general intialisation
  InitRecoParams() ;
  TObjArray ** rv = NULL ; 
  
  if (cycles > 0)
    SetCycle(cycles) ;  
	
  if ( task == AliQAv1::kRAWS ) {
    if (! fRawsQAList ) { 
      fRawsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
      for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
        fRawsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ;	 
        fRawsQAList[specie]->SetName(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(specie))) ;
      }
    }
    rv = fRawsQAList ;
  } else if ( task == AliQAv1::kDIGITSR ) {
    if ( ! fDigitsQAList ) {
      fDigitsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
      for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
        fDigitsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ; 
        fDigitsQAList[specie]->SetName(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(specie))) ; 
      }
    }
    rv = fDigitsQAList ;
  } else if ( task == AliQAv1::kRECPOINTS ) {
    if ( ! fRecPointsQAList ) {
      fRecPointsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
      for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
        fRecPointsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ; 
        fRecPointsQAList[specie]->SetName(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(specie))) ; 
      }
    }
    rv = fRecPointsQAList ;
  } else if ( task == AliQAv1::kESDS ) {
    if ( ! fESDsQAList ) {
      fESDsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
      for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
        fESDsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ;
        fESDsQAList[specie]->SetName(Form("%s_%s", GetName(), AliQAv1::GetTaskName(task).Data())); //, AliRecoParam::GetEventSpecieName(specie))) ; 
      }
    }
    rv = fESDsQAList ;
  }
  return rv ; 
}

//____________________________________________________________________________ 
void AliQADataMakerRec::Init(AliQAv1::TASKINDEX_t task, TObjArray ** list, Int_t run, Int_t cycles)
{
  // Intialisation by passing the list of QA data booked elsewhere
  
  InitRecoParams() ;
  fRun = run ;
  if (cycles > 0)
    SetCycle(cycles) ;  
	
  if ( task == AliQAv1::kRAWS ) {
    fRawsQAList = list ;	 
  } else if ( task == AliQAv1::kDIGITSR ) {
    fDigitsQAList = list ; 
  } else if ( task == AliQAv1::kRECPOINTS ) {
    fRecPointsQAList = list ; 
  } else if ( task == AliQAv1::kESDS ) {
    fESDsQAList = list ; 
  }
}

//____________________________________________________________________________
void AliQADataMakerRec::InitRecoParams() 
{
  // Get the recoparam form the OCDB 
  if (!fRecoParam) {
    AliDebug(AliQAv1::GetQADebugLevel(), Form("Loading reconstruction parameter objects for detector %s", GetName()));
    AliCDBPath path(GetName(),"Calib","RecoParam");
    AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());    
    if(!entry) {
      fRecoParam = NULL ; 
      AliDebug(AliQAv1::GetQADebugLevel(), Form("Couldn't find RecoParam entry in OCDB for detector %s",GetName()));
    }
    else {
      //      entry->SetOwner(kTRUE);
      TObject * recoParamObj = entry->GetObject() ; 
      if ( strcmp(recoParamObj->ClassName(), "TObjArray") == 0 ) {
        // The detector has only one set of reco parameters
        AliDebug(AliQAv1::GetQADebugLevel(), Form("Array of reconstruction parameters found for detector %s",GetName()));
        TObjArray *recoParamArray = static_cast<TObjArray*>(recoParamObj) ;
	//        recoParamArray->SetOwner(kTRUE);
        for (Int_t iRP=0; iRP<recoParamArray->GetEntriesFast(); iRP++) {
          fRecoParam = static_cast<AliDetectorRecoParam*>(recoParamArray->At(iRP)) ;
          if (!fRecoParam) 
            break ; 
          else if (fRecoParam->IsDefault()) 
            break ; 
        }
      }
      else if (recoParamObj->InheritsFrom("AliDetectorRecoParam")) {
        // The detector has only one set of reco parameters
        // Registering it in AliRecoParam
        AliDebug(AliQAv1::GetQADebugLevel(), Form("Single set of reconstruction parameters found for detector %s",GetName()));
        fRecoParam = static_cast<AliDetectorRecoParam*>(recoParamObj) ;
        static_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
      } else { 
        AliError(Form("No valid RecoParam object found in the OCDB for detector %s",GetName()));
      }
    }
  }
}

//____________________________________________________________________________ 
void AliQADataMakerRec::ResetDetector(AliQAv1::TASKINDEX_t task)
{
  // default reset that resets all the QA objects.
  // to be overloaded by detectors, if necessary

  TObjArray ** list = NULL ; 
  if ( task == AliQAv1::kRAWS ) {
    list = fRawsQAList ;	 
  } else if ( task == AliQAv1::kDIGITSR ) {
    list = fDigitsQAList ; 
  } else if ( task == AliQAv1::kRECPOINTS ) {
    list = fRecPointsQAList ; 
  } else if ( task == AliQAv1::kESDS ) {
    list = fESDsQAList ; 
  }
  //list was not initialized, skip
  if (!list) return; 
  
  for (int spec = 0; spec < AliRecoParam::kNSpecies; spec++) {
    if (!AliQAv1::Instance()->IsEventSpecieSet(AliRecoParam::ConvertIndex(spec))) continue;
    TIter next(list[spec]) ; 
    TObject *obj = NULL; 
    while ( (obj = next()) ) {
      if (obj->TestBit(AliQAv1::GetClonedBit())) { // this is array of cloned histos
	TObjArray* arr = (TObjArray*)obj;
	for (int ih=arr->GetEntriesFast();ih--;) {
	  TH1* histo = (TH1*)arr->At(ih); 
	  if (!histo) continue;
	  histo->Reset("ICE");
	  histo->ResetStats();
	}
      }
      else {
	((TH1*)obj)->Reset("ICE");
	((TH1*)obj)->ResetStats();
      }
    }
    ResetEvCountCycle(AliRecoParam::ConvertIndex(spec));
    ResetEvCountTotal(AliRecoParam::ConvertIndex(spec));
  }
}

//____________________________________________________________________________
void AliQADataMakerRec::StartOfCycle(Int_t run) 
{
  // Finishes a cycle of QA for all the tasks
  Bool_t samecycle = kFALSE ; 
  StartOfCycle(AliQAv1::kRAWS,      run, samecycle) ;
  samecycle = kTRUE ; 
  StartOfCycle(AliQAv1::kDIGITSR,   run, samecycle) ; 
  StartOfCycle(AliQAv1::kRECPOINTS, run, samecycle) ; 
  StartOfCycle(AliQAv1::kESDS,      run, samecycle) ;
}

//____________________________________________________________________________
void AliQADataMakerRec::StartOfCycle(AliQAv1::TASKINDEX_t task, Int_t run, const Bool_t sameCycle) 
{ 
  // Finishes a cycle of QA data acquistion
  if ( run > 0 ) fRun = run ; 
  if ( !sameCycle || fCurrentCycle == -1) {
    ResetCycle() ;
    if (fOutput) 
      fOutput->Close() ; 
    if (AliQAManager::QAManager(AliQAv1::kRECMODE)->IsSaveData())
      fOutput = AliQAv1::GetQADataFile(GetName(), fRun) ; 	
  }	
  if (AliQAManager::QAManager(AliQAv1::kRECMODE)->IsSaveData())
    AliDebug(AliQAv1::GetQADebugLevel(), Form(" Run %d Cycle %d task %s file %s", 
					    fRun, fCurrentCycle, AliQAv1::GetTaskName(task).Data(), fOutput->GetName() )) ;
  else
    AliDebug(AliQAv1::GetQADebugLevel(), Form(" Run %d Cycle %d task %s not saved", 
                                              fRun, fCurrentCycle, AliQAv1::GetTaskName(task).Data() )) ;    

  //	fDetectorDir = fOutput->GetDirectory(GetDetectorDirName()) ; 
  //	if (!fDetectorDir)
  //		fDetectorDir = fOutput->mkdir(GetDetectorDirName()) ; 
  //  
  //	TDirectory * subDir = fDetectorDir->GetDirectory(AliQAv1::GetTaskName(task)) ; 
  //	if (!subDir)
  //		subDir = fDetectorDir->mkdir(AliQAv1::GetTaskName(task)) ;  
  //  
  //  for ( Int_t specie = AliRecoParam::kDefault ; specie < AliRecoParam::kNSpecies ; specie++ ) {
  //    TDirectory * eventSpecieDir = subDir->GetDirectory(AliRecoParam::GetEventSpecieName(specie)) ; 
  //    if (!eventSpecieDir) 
  //      eventSpecieDir = subDir->mkdir(AliRecoParam::GetEventSpecieName(specie)) ; 
  //    TDirectory * expertDir = eventSpecieDir->GetDirectory(AliQAv1::GetExpert()) ; 
  //    if (!expertDir)
  //      expertDir = eventSpecieDir->mkdir(AliQAv1::GetExpert()) ; 
  //  } 
  StartOfDetectorCycle();
  ResetEvCountCycle();
}

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