ROOT logo
/**************************************************************************
 * Copyright(c) 1998-2007, 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: $ */

/*
  Based on AliPHOSQADataMaker
  Produces the data needed to calculate the quality assurance. 
  All data must be mergeable objects.
  P. Christiansen, Lund, January 2008

  Updated July 2011:
  ==================

  Major changes to accomodate updates of general DQM/QA changes to have per
  trigger histograms (for a given event specie).

  1) One instance of AliTPCdataQA only. (This also solves some old wishes by
  offline team to use less memory because event the 2d arrays for this object
  is not used). This now has a new flag for only keeping DQM info event by
  event! For this reason there is no need for a special DQM reset any more
  between runs!

  2) Fill the histogram for each event. The histograms are no longer filled
  from the AliTPCdataQA but per event.

  3) Use profiles for the RAW info. By adding the profiles event by event we
  get the correct event averages WITHOUT having to normalize in the end!
  Results should therefore also be directly mergable when that feature will
  come. (none of the other histograms are merged).

  This means that from the DQM/QA point of view the TPC DQM is now fully
  standard and should ease future developments.

  Updated June 2010:
  ==================

  The "beautification" of the online DQM histograms have been moved to
  an amore macro.  

  The per event RAW histograms have been modified in AliTPCdataQA and
  the copies have therefore also been modified here.

  The AliTPCdataQA can now be configured a bit from here: time bin
  range (extended default range to 1-1000, event range at start:
  0-100000, 1000 events per bin). (At least the parameters are not
  hardcoded:-)

  Implementation:
  ===============

  For the QA of the RAW data we use the class, AliTPCdataQA, from the
  existing TPC Calibration framework (which is more advanced than the
  standard QA framework) and extract the histograms at the end. The
  Analyse method of the AliTPCdataQA class is called in the method,
  EndOfDetectorCycle, and there also: 1d histogram(s) are projected
  and added to the QA list.
*/

#include "AliTPCQADataMakerRec.h"

// --- ROOT system ---
#include <TClonesArray.h>
#include <TString.h>
#include <TSystem.h>
#include <TBox.h>
#include <TLine.h>
#include <TAxis.h>
#include <TH1.h> 
#include <TProfile.h> 
#include <TProfile2D.h> 

// --- Standard library ---

// --- AliRoot header files ---
#include "AliQAChecker.h"
#include "AliESDEvent.h"
#include "AliESDtrack.h"
#include "AliLog.h"
#include "AliTPCCalPad.h"
#include "AliTPCCalROC.h"
#include "AliTPCClustersRow.h"
#include "AliTPCclusterMI.h"
#include "AliSimDigits.h"
#include <AliDetectorRecoParam.h>


ClassImp(AliTPCQADataMakerRec)

//____________________________________________________________________________ 
AliTPCQADataMakerRec::AliTPCQADataMakerRec() : 
AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTPC), 
		  "TPC Rec Quality Assurance Data Maker"),
fTPCdataQA(NULL),
fRawFirstTimeBin(1),
fRawLastTimeBin(1000)
{
  // ctor
  
  for(Int_t i = 0; i < 6; i++)
    fMapping[i] = 0;
}

//____________________________________________________________________________ 
AliTPCQADataMakerRec::AliTPCQADataMakerRec(const AliTPCQADataMakerRec& qadm) :
  AliQADataMakerRec(),
  fTPCdataQA(NULL),
  fRawFirstTimeBin(qadm.GetRawFirstTimeBin()),
  fRawLastTimeBin(qadm.GetRawLastTimeBin())
{
  //copy ctor 
  // Does not copy the calibration object, instead InitRaws have to be
  // called again
  SetName((const char*)qadm.GetName()) ; 
  SetTitle((const char*)qadm.GetTitle()); 

  for(Int_t i = 0; i < 6; i++)
    fMapping[i] = 0;
}

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

//__________________________________________________________________
AliTPCQADataMakerRec::~AliTPCQADataMakerRec()
{
  // Destructor
  delete fTPCdataQA; 

  for(Int_t i = 0; i < 6; i++) 
    delete fMapping[i];
}
 
//____________________________________________________________________________ 
void AliTPCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
{
  //Detector specific actions at end of cycle
  ResetEventTrigClasses();

  AliQAChecker::Instance()->Run(AliQAv1::kTPC, task, list) ;  
}


//____________________________________________________________________________ 
void AliTPCQADataMakerRec::InitESDs()
{
  //create ESDs histograms in ESDs subdir  
  const Bool_t expert   = kTRUE ; 
  const Bool_t image    = kTRUE ; 
  
  TH1F * histESDclusters = 
    new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
	     160, 0, 160);
  histESDclusters->Sumw2();
  Add2ESDsList(histESDclusters, kClusters, !expert, image);

  TH1F * histESDratio = 
    new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
	     100, 0, 1);
  histESDratio->Sumw2();
  Add2ESDsList(histESDratio, kRatio, !expert, image);
  
  TH1F * histESDpt = 
    new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
	     50, 0, 5);
  histESDpt->Sumw2();
  Add2ESDsList(histESDpt, kPt, !expert, image);
  //
  ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
}

//____________________________________________________________________________ 
void AliTPCQADataMakerRec::InitRaws()
{
  //
  // Adding the raw 
  //  

  // Modified: 7/7 - 2008
  // Laurent Aphecetche pointed out that the mapping was read from file
  // for each event, so now we read in the map here and set if for 
  // the raw data qa
  const Bool_t expert   = kTRUE ; 
  const Bool_t saveCorr = kTRUE ; 
  const Bool_t image    = kTRUE ; 
  
  // It might happen that we will be in this method a few times
  // (we create the dataQA at the first call to this method)
  if(!fTPCdataQA) {
    fTPCdataQA = new AliTPCdataQA();
    LoadMaps(); // Load Altro maps
    fTPCdataQA->SetAltroMapping(fMapping); // set Altro mapping
    fTPCdataQA->SetRangeTime(fRawFirstTimeBin, fRawLastTimeBin); // set time bin interval 
    fTPCdataQA->SetIsDQM(kTRUE);
  }

  TProfile * histRawsOccupancyVsSector = 
    new TProfile("hRawsOccupancyVsSector", "Occupancy vs sector; Sector; Occupancy",
	     72, 0, 72);
  histRawsOccupancyVsSector->SetMarkerStyle(20);
  histRawsOccupancyVsSector->SetOption("P");
  histRawsOccupancyVsSector->SetStats(kFALSE);
  if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib )
    Add2RawsList(histRawsOccupancyVsSector, kRawsOccupancyVsSector, expert, image, !saveCorr);
  else
    Add2RawsList(histRawsOccupancyVsSector, kRawsOccupancyVsSector, !expert, image, !saveCorr);
    
  TProfile * histRawsQVsSector = 
    new TProfile("hRawsQVsSector", "<Q> vs sector; Sector; <Q>",
	     72, 0, 72);
  Add2RawsList(histRawsQVsSector, kRawsQVsSector, expert, !image, !saveCorr);

  TProfile * histRawsQmaxVsSector = 
    new TProfile("hRawsQmaxVsSector", "<Qmax> vs sector; Sector; <Qmax>",
	     72, 0, 72);
  histRawsQmaxVsSector->SetMarkerStyle(20);
  histRawsQmaxVsSector->SetOption("P");
  histRawsQmaxVsSector->SetStats(kFALSE);
  if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib )
    Add2RawsList(histRawsQmaxVsSector, kRawsQmaxVsSector, expert, image, !saveCorr);
  else
    Add2RawsList(histRawsQmaxVsSector, kRawsQmaxVsSector, !expert, image, !saveCorr);

  TProfile2D * histRawsOccupancy2dVsSector = 
    new TProfile2D("hRawsOccupancy2dVsSector", "Occupancy vs sector; Sector; Patch",
		   72, 0, 36, 6, 0, 6);
  histRawsOccupancy2dVsSector->SetOption("COLZ");
  histRawsOccupancy2dVsSector->SetStats(kFALSE);
  if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib )
    Add2RawsList(histRawsOccupancy2dVsSector, kRawsOccupancy2dVsSector, expert, image, !saveCorr);
  else
    Add2RawsList(histRawsOccupancy2dVsSector, kRawsOccupancy2dVsSector, !expert, image, !saveCorr);
    
  //
  ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
}

//____________________________________________________________________________ 
void AliTPCQADataMakerRec::InitDigits()
{
  const Bool_t expert   = kTRUE ; 
  const Bool_t image    = kTRUE ; 
  TH1F * histDigitsADC = 
    new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
             1000, 0, 1000);
  histDigitsADC->Sumw2();
  Add2DigitsList(histDigitsADC, kDigitsADC, !expert, image);
  //
  ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
}

//____________________________________________________________________________ 
void AliTPCQADataMakerRec::InitRecPoints()
{
  const Bool_t expert   = kTRUE ; 
  const Bool_t image    = kTRUE ; 
  
  TH1F * histRecPointsQmaxShort = 
    new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
	     100, 0, 300);
  histRecPointsQmaxShort->Sumw2();
  Add2RecPointsList(histRecPointsQmaxShort, kQmaxShort, !expert, image);

  TH1F * histRecPointsQmaxMedium = 
    new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
	     100, 0, 300);
  histRecPointsQmaxMedium->Sumw2();
  Add2RecPointsList(histRecPointsQmaxMedium, kQmaxMedium, !expert, image);

  TH1F * histRecPointsQmaxLong = 
    new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
	     100, 0, 300);
  histRecPointsQmaxLong->Sumw2();
  Add2RecPointsList(histRecPointsQmaxLong, kQmaxLong, !expert, image);

  TH1F * histRecPointsQShort = 
    new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
	     100, 0, 2000);
  histRecPointsQShort->Sumw2();
  Add2RecPointsList(histRecPointsQShort, kQShort, !expert, image);

  TH1F * histRecPointsQMedium = 
    new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
	     100, 0, 2000);
  histRecPointsQMedium->Sumw2();
  Add2RecPointsList(histRecPointsQMedium, kQMedium, !expert, image);

  TH1F * histRecPointsQLong = 
    new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
	     100, 0, 2000);
  histRecPointsQLong->Sumw2();
  Add2RecPointsList(histRecPointsQLong, kQLong, !expert, image);

  TH1F * histRecPointsRow = 
    new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
	     159, 0, 159);
  histRecPointsRow->Sumw2();
  Add2RecPointsList(histRecPointsRow, kRow, !expert, image);
  //
  ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
}

//____________________________________________________________________________
void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
{
  // make QA data from ESDs
 
  const Int_t nESDTracks = esd->GetNumberOfTracks();
  Int_t nTPCtracks = 0; 
  for(Int_t i = 0; i < nESDTracks; i++) {
    
    AliESDtrack * track = esd->GetTrack(i);
    
    if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
      continue;
    
    nTPCtracks++;
    
    Int_t nTPCclusters         = track->GetTPCNcls();
    Int_t nTPCclustersFindable = track->GetTPCNclsF();
    if ( nTPCclustersFindable<=0) continue;
    FillESDsData(kClusters,nTPCclusters);
    FillESDsData(kRatio,Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
    FillESDsData(kPt,track->Pt()); 
  }
  //
  IncEvCountCycleESDs();
  IncEvCountTotalESDs();
  //
}

//____________________________________________________________________________
void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
{
  //
  // To make QA for the RAW data we use the TPC Calibration framework 
  // to handle the data and then in the end extract the data
  //
  
  GetRawsData(0); // dummy call to init raw data
  rawReader->Reset() ; 
  if (! fTPCdataQA ) {

    AliError("No TPC data QA (no call to InitRaws?)!!!!") ; 
  } else {  

    if(fTPCdataQA->GetIsDQM() == kFALSE)
      AliError("Data QA has to be initialized as DQM!!!!") ; 

    // Fill profile data
    fTPCdataQA->ResetProfiles();
    
    if(fTPCdataQA->ProcessEvent(rawReader)) { // means that TPC data was processed  

      fTPCdataQA->FillOccupancyProfile();
      
      // Fill histograms    
      TObjArray *arrRW = GetMatchingRawsData(kRawsOccupancyVsSector); // all kRawsOccupancyVsSector clones matching to triggers
      for (int ih=arrRW->GetEntriesFast();ih--;) {
	TProfile* hRawsOccupancyVsSector = dynamic_cast<TProfile*>(arrRW->At(ih));
	if (hRawsOccupancyVsSector) hRawsOccupancyVsSector->Add(fTPCdataQA->GetHistOccVsSector());
      }
      arrRW = GetMatchingRawsData(kRawsOccupancy2dVsSector);
      for (int ih=arrRW->GetEntriesFast();ih--;) {
	TProfile2D* hRawsOccupancy2dVsSector = dynamic_cast<TProfile2D*>(arrRW->At(ih));
	if (hRawsOccupancy2dVsSector) hRawsOccupancy2dVsSector->Add(fTPCdataQA->GetHistOcc2dVsSector());
      }
      arrRW = GetMatchingRawsData(kRawsQVsSector);
      for (int ih=arrRW->GetEntriesFast();ih--;) {
	TProfile* hRawsQVsSector = dynamic_cast<TProfile*>(arrRW->At(ih));
	if (hRawsQVsSector) hRawsQVsSector->Add(fTPCdataQA->GetHistQVsSector());
      }
      arrRW = GetMatchingRawsData(kRawsQmaxVsSector);
      for (int ih=arrRW->GetEntriesFast();ih--;) {
	TProfile* hRawsQmaxVsSector = dynamic_cast<TProfile*>(arrRW->At(ih));
	if (hRawsQmaxVsSector) hRawsQmaxVsSector->Add(fTPCdataQA->GetHistQmaxVsSector());
      }
      //
      IncEvCountCycleRaws();
      IncEvCountTotalRaws();
      //
    }
  }
}

//____________________________________________________________________________
void AliTPCQADataMakerRec::MakeDigits(TTree* digitTree)
{
 
  TBranch* branch = digitTree->GetBranch("Segment");
  AliSimDigits* digArray = 0;
  branch->SetAddress(&digArray);
  
  Int_t nEntries = Int_t(digitTree->GetEntries());
  
  for (Int_t n = 0; n < nEntries; n++) {
    
    digitTree->GetEvent(n);
    
    if (digArray->First())
      do {
        Float_t dig = digArray->CurrentDigit();
        
        FillDigitsData(kDigitsADC,dig);
      } while (digArray->Next());    
  }
  //
  IncEvCountCycleDigits();
  IncEvCountTotalDigits();
  //
}

//____________________________________________________________________________
void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
{
  
  AliTPCClustersRow* clrow = 0x0;
  TBranch* branch = recTree->GetBranch("Segment");  
  branch->SetAddress(&clrow);
  TClonesArray * clarray = 0x0;

  const Int_t nEntries = Int_t(recTree->GetEntries());
  for (Int_t i = 0; i < nEntries; i++) {
    
    branch->GetEntry(i);

    clarray = clrow->GetArray();

    if (!clarray) continue;

    const Int_t nClusters = clarray->GetEntriesFast();
    for (Int_t icl=0; icl < nClusters; icl++){
      
      AliTPCclusterMI* cluster = 
	(AliTPCclusterMI*)clarray->At(icl);
      
      Float_t Qmax = cluster->GetMax();
      Float_t Q    = cluster->GetQ();
      Int_t   row  = cluster->GetRow();

      if(cluster->GetDetector()<36) { // IROC (short pads)

	FillRecPointsData(kQmaxShort,Qmax);
	FillRecPointsData(kQShort,Q);
      } else { // OROC (medium and long pads)
	row += 63;
	if(cluster->GetRow()<64) { // medium pads

	  FillRecPointsData(kQmaxMedium,Qmax);
	  FillRecPointsData(kQMedium,Q);
	} else { // long pads

	  FillRecPointsData(kQmaxLong,Qmax);
	  FillRecPointsData(kQLong,Q);
	}
      }
      
      FillRecPointsData(kRow,row);
    } // end loop over clusters
  } // end loop over tree
  //
  IncEvCountCycleRecPoints();
  IncEvCountTotalRecPoints();
  //
}

//____________________________________________________________________________
void AliTPCQADataMakerRec::LoadMaps()
{
  TString path = gSystem->Getenv("ALICE_ROOT");
  path += "/TPC/mapping/Patch";

  for(Int_t i = 0; i < 6; i++) {

    if(fMapping[i]!=0) // mapping already loaded
      continue;
    TString path2 = path;
    path2 += i;
    path2 += ".data";
    fMapping[i] = new AliTPCAltroMapping(path2.Data());
  }
}

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