ROOT logo
//------------------------------------------------------------------------------
// Implementation of AliPerformanceTPC class. It keeps information from 
// comparison of reconstructed and MC particle tracks. In addtion, 
// it keeps selection cuts used during comparison. The comparison 
// information is stored in the ROOT histograms. Analysis of these 
// histograms can be done by using Analyse() class function. The result of 
// the analysis (histograms/graphs) are stored in the folder which is
// a data member of AliPerformanceTPC.
//
// Author: J.Otwinowski 04/02/2008 
// Changes by M.Knichel 15/10/2010
//------------------------------------------------------------------------------

/*
 
  // after running comparison task, read the file, and get component
  gROOT->LoadMacro("$ALICE_ROOT/PWGPP/Macros/LoadMyLibs.C");
  LoadMyLibs();

  TFile f("Output.root");
  AliPerformanceTPC * compObj = (AliPerformanceTPC*)coutput->FindObject("AliPerformanceTPC");
 
  // analyse comparison data
  compObj->Analyse();

  // the output histograms/graphs will be stored in the folder "folderTPC" 
  compObj->GetAnalysisFolder()->ls("*");

  // user can save whole comparison object (or only folder with anlysed histograms) 
  // in the seperate output file (e.g.)
  TFile fout("Analysed_TPC.root","recreate");
  compObj->Write(); // compObj->GetAnalysisFolder()->Write();
  fout.Close();

*/

#include "TCanvas.h"
#include "TH1.h"
#include "TH2.h"
#include "TH3.h"
#include "TAxis.h"
#include "TPostScript.h"
#include "TString.h"
#include "TUUID.h"
#include "TTree.h"
#include "TChain.h"
#include "AliTPCPerformanceSummary.h"
#include "TSystem.h"

#include "AliPerformanceTPC.h" 
#include "AliESDEvent.h" 
#include "AliESDVertex.h"
#include "AliESDtrack.h"
#include "AliLog.h" 
#include "AliMCEvent.h" 
#include "AliHeader.h" 
#include "AliGenEventHeader.h" 
#include "AliStack.h" 
#include "AliMCInfoCuts.h" 
#include "AliRecInfoCuts.h" 
#include "AliTracker.h" 
#include "AliTreeDraw.h" 
#include "AliTPCTransform.h" 
#include "AliTPCseed.h" 
#include "AliTPCcalibDB.h" 
#include "AliESDfriend.h" 
#include "AliESDfriendTrack.h" 
#include "AliTPCclusterMI.h" 

using namespace std;

ClassImp(AliPerformanceTPC)

Bool_t AliPerformanceTPC::fgMergeTHnSparse = kFALSE;
Bool_t AliPerformanceTPC::fgUseMergeTHnSparse = kFALSE;


//_____________________________________________________________________________
/*
AliPerformanceTPC::AliPerformanceTPC():
  AliPerformanceObject("AliPerformanceTPC"),
  fTPCClustHisto(0),
  fTPCEventHisto(0),
  fTPCTrackHisto(0),
  fFolderObj(0),

  // Cuts 
  fCutsRC(0),  
  fCutsMC(0),  

  // histogram folder 
  fAnalysisFolder(0),
  
  fUseHLT(kFALSE)

{
  Init();
}
*/

//_____________________________________________________________________________
AliPerformanceTPC::AliPerformanceTPC(const Char_t* name, const Char_t* title,Int_t analysisMode,Bool_t hptGenerator, Int_t run, Bool_t highMult):
  AliPerformanceObject(name,title,run,highMult),
  fTPCClustHisto(0),
  fTPCEventHisto(0),
  fTPCTrackHisto(0),
  fFolderObj(0),

  // Cuts 
  fCutsRC(0),  
  fCutsMC(0),  

  // histogram folder 
  fAnalysisFolder(0),
  
  fUseHLT(kFALSE)

{
  // named constructor	
  // 
  SetAnalysisMode(analysisMode);
  SetHptGenerator(hptGenerator);

  Init();
}


//_____________________________________________________________________________
AliPerformanceTPC::~AliPerformanceTPC()
{
  // destructor
   
  if(fTPCClustHisto) delete fTPCClustHisto; fTPCClustHisto=0;     
  if(fTPCEventHisto) delete fTPCEventHisto; fTPCEventHisto=0;     
  if(fTPCTrackHisto) delete fTPCTrackHisto; fTPCTrackHisto=0;   
  if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
  if(fFolderObj) delete fFolderObj; fFolderObj=0;
}


//_____________________________________________________________________________
void AliPerformanceTPC::Init()
{
  //
  // histogram bining
  //

 
  // set pt bins
  Int_t nPtBins = 50;
  Double_t ptMin = 1.e-2, ptMax = 20.;

  Double_t *binsPt = 0;

  if (IsHptGenerator())  { 
        ptMax = 100.;
  } 
   binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);


  /*
  const Int_t  nCOverPtBins = 80;
  Double_t coverptMin = -10, coverptMax = 10;
  Double_t *binsCOverPtP = 0;
  Double_t *binsCOverPt = new Double_t[nCOverPtBins+1];
  binsCOverPtP = CreateLogAxis(nCOverPtBins/2,0.04,coverptMax-0.04);
  for(Int_t i=0; i < nCOverPtBins/2; i++){
    binsCOverPt[nCOverPtBins - i] = binsCOverPtP[nCOverPtBins/2  - i];
    binsCOverPt[i] = 0 - binsCOverPtP[nCOverPtBins/2  - i];
 }
 */

  /*
  Int_t nPtBins = 31;
  Double_t binsPt[32] = {0.,0.05,0.1,0.15,0.2,0.25,0.3,0.35,0.4,0.45,0.5,0.55,0.6,0.7,0.8,0.9,1.0,1.2,1.4,1.6,1.8,2.0,2.25,2.5,2.75,3.,3.5,4.,5.,6.,8.,10.};
  Double_t ptMin = 0., ptMax = 10.; 

  if(IsHptGenerator() == kTRUE) {
    nPtBins = 100;
    ptMin = 0.; ptMax = 100.; 
  }
  */
  //

  //
  //padRow:phi:TPCSide
  Int_t binsTPCClustHisto[3] =   {160,  144,  2};
  Double_t minTPCClustHisto[3] = {0.,   0.,   0.};
  Double_t maxTPCClustHisto[3] = {160., 2.*TMath::Pi(), 2.};

  fTPCClustHisto = new THnSparseF("fTPCClustHisto","padRow:phi:TPCSide",3,binsTPCClustHisto,minTPCClustHisto,maxTPCClustHisto);
  fTPCClustHisto->GetAxis(0)->SetTitle("padRow");
  fTPCClustHisto->GetAxis(1)->SetTitle("phi (rad)");
  fTPCClustHisto->GetAxis(2)->SetTitle("TPCSide");
  //fTPCClustHisto->Sumw2();

  //padRow:phi:TPCSide:pad:detector:glZ
  /*
  Int_t binsTPCClustHisto[6] =   {160,  144,  2, 128, 72, 50};
  Double_t minTPCClustHisto[6] = {0.,   0.,   0., 0, 0, -250};
  Double_t maxTPCClustHisto[6] = {160., 2.*TMath::Pi(), 2., 128, 72,250};

  fTPCClustHisto = new THnSparseF("fTPCClustHisto","padRow:phi:TPCSide:pad:detector:gZ",6,binsTPCClustHisto,minTPCClustHisto,maxTPCClustHisto);
  fTPCClustHisto->GetAxis(0)->SetTitle("padRow");
  fTPCClustHisto->GetAxis(1)->SetTitle("phi (rad)");
  fTPCClustHisto->GetAxis(2)->SetTitle("TPCSide");
  fTPCClustHisto->GetAxis(3)->SetTitle("pad");
  fTPCClustHisto->GetAxis(4)->SetTitle("detector");
  fTPCClustHisto->GetAxis(5)->SetTitle("glZ (cm)");
  //fTPCClustHisto->Sumw2();
  */
  
  Float_t scaleVxy = 1.0;
  if(fAnalysisMode !=0) scaleVxy = 0.1; 

  Int_t maxMult;
  if (fHighMultiplicity) { maxMult = 4001; scaleVxy = 0.1;} else { maxMult = 151; }
  // Xv:Yv:Zv:mult:multP:multN:vertStatus
  Int_t binsTPCEventHisto[7]=  {100,  100,   100,  maxMult,  maxMult,  maxMult, 2   };
  Double_t minTPCEventHisto[7]={-10.*scaleVxy, -10.*scaleVxy, -30., -0.5,  -0.5,  -0.5, -0.5  };
  Double_t maxTPCEventHisto[7]={ 10.*scaleVxy,  10.*scaleVxy,  30.,  maxMult-0.5,  maxMult-0.5, maxMult-0.5, 1.5 };

  fTPCEventHisto = new THnSparseF("fTPCEventHisto","Xv:Yv:Zv:mult:multP:multN:vertStatus",7,binsTPCEventHisto,minTPCEventHisto,maxTPCEventHisto);
  fTPCEventHisto->GetAxis(0)->SetTitle("Xv (cm)");
  fTPCEventHisto->GetAxis(1)->SetTitle("Yv (cm)");
  fTPCEventHisto->GetAxis(2)->SetTitle("Zv (cm)");
  fTPCEventHisto->GetAxis(3)->SetTitle("mult");
  fTPCEventHisto->GetAxis(4)->SetTitle("multP");
  fTPCEventHisto->GetAxis(5)->SetTitle("multN");
  fTPCEventHisto->GetAxis(6)->SetTitle("vertStatus");
  //fTPCEventHisto->Sumw2();

  Float_t scaleDCA = 1.0;
  if(fAnalysisMode !=0) scaleDCA = 0.1; 
  // nTPCClust:chi2PerTPCClust:nTPCClustFindRatio:DCAr:DCAz:eta:phi:pt:charge:vertStatus
   //Int_t binsTPCTrackHisto[10]=  { 160,  20,  60,  30, 30,  30,   144,             nPtBins,   nCOverPtBins, 2 };
   //Double_t minTPCTrackHisto[10]={ 0.,   0.,  0., -3*scaleDCA, -3.*scaleDCA, -1.5, 0.,             ptMin,   coverptMin, -0.5 };
   //Double_t maxTPCTrackHisto[10]={ 160., 5., 1.2, 3*scaleDCA,  3.*scaleDCA,  1.5, 2.*TMath::Pi(), ptMax,    coverptMax,  1.5 };
   Int_t binsTPCTrackHisto[10]=  { 160,  20,  60,  30, 30,  30,   144,             nPtBins,   3, 2 };
   Double_t minTPCTrackHisto[10]={ 0.,   0.,  0., -3*scaleDCA, -3.*scaleDCA, -1.5, 0.,             ptMin,  -1.5, -0.5 };
   Double_t maxTPCTrackHisto[10]={ 160., 5., 1.2, 3*scaleDCA,  3.*scaleDCA,  1.5, 2.*TMath::Pi(), ptMax,    1.5,  1.5 };
  
  //fTPCTrackHisto = new THnSparseF("fTPCTrackHisto","nClust:chi2PerClust:nClust/nFindableClust:DCAr:DCAz:eta:phi:pt:charge/pt:vertStatus",10,binsTPCTrackHisto,minTPCTrackHisto,maxTPCTrackHisto);
  fTPCTrackHisto = new THnSparseF("fTPCTrackHisto","nClust:chi2PerClust:nClust/nFindableClust:DCAr:DCAz:eta:phi:pt:charge:vertStatus",10,binsTPCTrackHisto,minTPCTrackHisto,maxTPCTrackHisto);
  fTPCTrackHisto->SetBinEdges(7,binsPt);
  //fTPCTrackHisto->SetBinEdges(8,binsCOverPt);
  fTPCTrackHisto->GetAxis(0)->SetTitle("nClust");
  fTPCTrackHisto->GetAxis(1)->SetTitle("chi2PerClust");
  fTPCTrackHisto->GetAxis(2)->SetTitle("nClust/nFindableClust");
  fTPCTrackHisto->GetAxis(3)->SetTitle("DCAr (cm)");
  fTPCTrackHisto->GetAxis(4)->SetTitle("DCAz (cm)");
  fTPCTrackHisto->GetAxis(5)->SetTitle("#eta");
  fTPCTrackHisto->GetAxis(6)->SetTitle("#phi (rad)");
  fTPCTrackHisto->GetAxis(7)->SetTitle("p_{T} (GeV/c)");
  //fTPCTrackHisto->GetAxis(8)->SetTitle("charge/pt");
  fTPCTrackHisto->GetAxis(8)->SetTitle("charge");
  fTPCTrackHisto->GetAxis(9)->SetTitle("vertStatus");
  //fTPCTrackHisto->Sumw2();

  // Init cuts 
  if(!fCutsMC) {
    AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
  }
  if(!fCutsRC) {
    AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object"); 
  }

  // init folder
  fAnalysisFolder = CreateFolder("folderTPC","Analysis Resolution Folder");

  //delete []binsCOverPt;
  
   // save merge status in object
  fMergeTHnSparseObj = fgMergeTHnSparse;
  
}


//_____________________________________________________________________________
void AliPerformanceTPC::ProcessTPC(AliStack* const stack, AliESDtrack *const esdTrack, AliESDEvent *const esdEvent, Bool_t vertStatus)
{
//
// fill TPC QA info
//
  if(!esdEvent) return;
  if(!esdTrack) return;

  if(IsUseTOFBunchCrossing())
    if(esdTrack->GetTOFBunchCrossing(esdEvent->GetMagneticField())!=0)
      return;

  if( IsUseTrackVertex() ) 
  { 
    // Relate TPC inner params to prim. vertex
    const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTracks();
    Double_t x[3]; esdTrack->GetXYZ(x);
    Double_t b[3]; AliTracker::GetBxByBz(x,b);
    //    Bool_t isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
    Bool_t isOK=kFALSE;
    if(fabs(b[2])>0.000001)
     isOK = esdTrack->RelateToVertexTPCBxByBz(vtxESD, b, kVeryBig);
    if(!isOK) return;

    /*
      // JMT -- recaluclate DCA for HLT if not present
      if ( dca[0] == 0. && dca[1] == 0. ) {
        track->GetDZ( vtxESD->GetX(), vtxESD->GetY(), vtxESD->GetZ(), esdEvent->GetMagneticField(), dca );
      }
    */
  }

  // Fill TPC only resolution comparison information 
  const AliExternalTrackParam *track = esdTrack->GetTPCInnerParam();
  if(!track) return;

  Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
  esdTrack->GetImpactParametersTPC(dca,cov);

  Float_t q = esdTrack->Charge();
  Float_t pt = track->Pt();
  Float_t eta = track->Eta();
  Float_t phi = track->Phi();
  Int_t nClust = esdTrack->GetTPCclusters(0);
  Int_t nFindableClust = esdTrack->GetTPCNclsF();

  Float_t chi2PerCluster = 0.;
  if(nClust>0.) chi2PerCluster = esdTrack->GetTPCchi2()/Float_t(nClust);

  Float_t clustPerFindClust = 0.;
  if(nFindableClust>0.) clustPerFindClust = Float_t(nClust)/nFindableClust;
  
  Float_t qpt = 0;
  if( fabs(pt)>0 ) qpt = q/fabs(pt);

  // filter out noise tracks
  if(esdTrack->GetTPCsignal() < 5) return;

  //
  // select primaries
  //
  Double_t dcaToVertex = -1;
  if( fCutsRC->GetDCAToVertex2D() ) 
  {
      dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY() + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ()); 
  }
  if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
  if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
  if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;

  //Double_t vTPCTrackHisto[10] = {nClust,chi2PerCluster,clustPerFindClust,dca[0],dca[1],eta,phi,pt,qpt,vertStatus};
  Double_t vTPCTrackHisto[10] = {static_cast<Double_t>(nClust),static_cast<Double_t>(chi2PerCluster),static_cast<Double_t>(clustPerFindClust),static_cast<Double_t>(dca[0]),static_cast<Double_t>(dca[1]),static_cast<Double_t>(eta),static_cast<Double_t>(phi),static_cast<Double_t>(pt),static_cast<Double_t>(q),static_cast<Double_t>(vertStatus)};
  fTPCTrackHisto->Fill(vTPCTrackHisto); 
 
  //
  // Fill rec vs MC information
  //
  if(!stack) return;

}


//_____________________________________________________________________________
void AliPerformanceTPC::ProcessTPCITS(AliStack* const stack, AliESDtrack *const esdTrack, AliESDEvent* const esdEvent, Bool_t vertStatus)
{
  // Fill comparison information (TPC+ITS) 
  if(!esdTrack) return;
  if(!esdEvent) return;

  if( IsUseTrackVertex() ) 
  { 
    // Relate TPC inner params to prim. vertex
    const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTracks();
    Double_t x[3]; esdTrack->GetXYZ(x);
    Double_t b[3]; AliTracker::GetBxByBz(x,b);
    Bool_t isOK = esdTrack->RelateToVertexBxByBz(vtxESD, b, kVeryBig);
    if(!isOK) return;

    /*
      // JMT -- recaluclate DCA for HLT if not present
      if ( dca[0] == 0. && dca[1] == 0. ) {
        track->GetDZ( vtxESD->GetX(), vtxESD->GetY(), vtxESD->GetZ(), esdEvent->GetMagneticField(), dca );
      }
    */
  }

  Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
  esdTrack->GetImpactParameters(dca,cov);

  if ((esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0) return; // ITS refit
  if ((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return; // TPC refit
  if ((esdTrack->HasPointOnITSLayer(0)==kFALSE)&&(esdTrack->HasPointOnITSLayer(1)==kFALSE)) return; // at least one SPD
  //if (esdTrack->GetITSclusters(0)<fCutsRC->GetMinNClustersITS()) return;  // min. nb. ITS clusters

  Float_t q = esdTrack->Charge();
  Float_t pt = esdTrack->Pt();
  Float_t eta = esdTrack->Eta();
  Float_t phi = esdTrack->Phi();
  Int_t nClust = esdTrack->GetTPCclusters(0);
  Int_t nFindableClust = esdTrack->GetTPCNclsF();

  Float_t chi2PerCluster = 0.;
  if(nClust>0.) chi2PerCluster = esdTrack->GetTPCchi2()/Float_t(nClust);

  Float_t clustPerFindClust = 0.;
  if(nFindableClust>0.) clustPerFindClust = Float_t(nClust)/nFindableClust;
  
  Float_t qpt = 0;
  if( fabs(pt)>0 ) qpt = q/fabs(pt);

  //
  // select primaries
  //
  Double_t dcaToVertex = -1;
  if( fCutsRC->GetDCAToVertex2D() ) 
  {
      dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY() + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ()); 
  }
  if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
  if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
  if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;

  Double_t vTPCTrackHisto[10] = {static_cast<Double_t>(nClust),static_cast<Double_t>(chi2PerCluster),static_cast<Double_t>(clustPerFindClust),static_cast<Double_t>(dca[0]),static_cast<Double_t>(dca[1]),static_cast<Double_t>(eta),static_cast<Double_t>(phi),static_cast<Double_t>(pt),static_cast<Double_t>(q),static_cast<Double_t>(vertStatus)};
  fTPCTrackHisto->Fill(vTPCTrackHisto); 
 
  //
  // Fill rec vs MC information
  //
  if(!stack) return;
}


//_____________________________________________________________________________
void AliPerformanceTPC::ProcessConstrained(AliStack* const /*stack*/, AliESDtrack *const /*esdTrack*/, AliESDEvent* const /*esdEvent*/)
{
  // Fill comparison information (constarained parameters) 
  AliDebug(AliLog::kWarning, "Warning: Not implemented");
}


//_____________________________________________________________________________
void AliPerformanceTPC::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
{
  // Process comparison information 
  //

  if(!esdEvent) 
  {
    Error("Exec","esdEvent not available");
    return;
  }
  AliHeader* header = 0;
  AliGenEventHeader* genHeader = 0;
  AliStack* stack = 0;
  TArrayF vtxMC(3);
  
  if(bUseMC)
  {
    if(!mcEvent) {
      Error("Exec","mcEvent not available");
      return;
    }
    // get MC event header
    header = mcEvent->Header();
    if (!header) {
      Error("Exec","Header not available");
      return;
    }
    // MC particle stack
    stack = mcEvent->Stack();
    if (!stack) {
      Error("Exec","Stack not available");
      return;
    }
    // get MC vertex
    genHeader = header->GenEventHeader();
    if (!genHeader) {
      Error("Exec","Could not retrieve genHeader from Header");
      return;
    }
    genHeader->PrimaryVertex(vtxMC);
  } 
  
  // trigger
  if(!bUseMC &&GetTriggerClass()) {
    Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
    if(!isEventTriggered) return; 
  }

  // get TPC event vertex
  const AliESDVertex *vtxESD = NULL; 
  if(fUseTrackVertex) {
    vtxESD = esdEvent->GetPrimaryVertexTracks();
  } else {
    vtxESD = esdEvent->GetPrimaryVertexTPC();
  }
  if(!vtxESD) return;

  //  events with rec. vertex
  Int_t mult=0; Int_t multP=0; Int_t multN=0;
  
  // changed to take all events but store vertex status
//  if(vtxESD->GetStatus() >0)
//  {
  // store vertex status
  Bool_t vertStatus = vtxESD->GetStatus();
  //  Process ESD events
  for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
  { 
    AliESDtrack *track = esdEvent->GetTrack(iTrack);
    if(!track) continue;

    // if not fUseKinkDaughters don't use tracks with kink index > 0
    if(!fUseKinkDaughters && track->GetKinkIndex(0) > 0) continue;
    
    if(bUseESDfriend && esdFriend && esdFriend->TestSkipBit()==kFALSE && iTrack<esdFriend->GetNumberOfTracks()) 
    {
      AliESDfriendTrack *friendTrack=esdFriend->GetTrack(iTrack);
      if(friendTrack) 
      {
        //
        TObject *calibObject=0;
        AliTPCseed *seed=0;
        for (Int_t j=0;(calibObject=friendTrack->GetCalibObject(j));++j) {
	    if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) {
	    break;
	  }
        }

        // 
	for (Int_t irow=0;irow<159;irow++) {
	if(!seed) continue;
	  
	  AliTPCclusterMI *cluster=seed->GetClusterPointer(irow);
	  if (!cluster) continue;

	     Float_t gclf[3];
	     cluster->GetGlobalXYZ(gclf);

	     //Double_t x[3]={cluster->GetRow(),cluster->GetPad(),cluster->GetTimeBin()};
	     //Int_t i[1]={cluster->GetDetector()};
             //transform->Transform(x,i,0,1);
	     //printf("gx %f gy  %f  gz %f \n", cluster->GetX(), cluster->GetY(),cluster->GetZ());
	     //printf("gclf[0] %f gclf[1]  %f  gclf[2] %f \n", gclf[0], gclf[1],  gclf[2]);
     
             Int_t TPCside; 
	     if(gclf[2]>0.) TPCside=0; // A side 
	     else TPCside=1;

	     //
             //Double_t vTPCClust1[3] = { gclf[0], gclf[1],  TPCside };
             //fTPCClustHisto1->Fill(vTPCClust1);

             //  
	     Double_t phi = TMath::ATan2(gclf[1],gclf[0]);
	     if(phi < 0) phi += 2.*TMath::Pi();
	    
             //Float_t pad = cluster->GetPad();
             //Int_t detector = cluster->GetDetector();
             //Double_t vTPCClust[6] = { irow, phi, TPCside, pad, detector, gclf[2] };
             Double_t vTPCClust[3] = { static_cast<Double_t>(irow), phi, static_cast<Double_t>(TPCside) };
             fTPCClustHisto->Fill(vTPCClust);
        }
      }
    }

    if(GetAnalysisMode() == 0) ProcessTPC(stack,track,esdEvent,vertStatus);
    else if(GetAnalysisMode() == 1) ProcessTPCITS(stack,track,esdEvent,vertStatus);
    else if(GetAnalysisMode() == 2) ProcessConstrained(stack,track,esdEvent);
    else {
      printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
      return;
    }

    // TPC only
    if(!fUseHLT){
      if(GetAnalysisMode() == 0) {
        AliESDtrack *tpcTrack = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
        if(!tpcTrack) continue;
        // track selection
        if( fCutsRC->AcceptTrack(tpcTrack) ) { 
	  mult++;
	  if(tpcTrack->Charge()>0.) multP++;
	  if(tpcTrack->Charge()<0.) multN++;
        }
        if(tpcTrack) delete tpcTrack;
      }
      else {
       // track selection
        if( fCutsRC->AcceptTrack(track) ) { 
	  mult++;
	  if(track->Charge()>0.) multP++;
	  if(track->Charge()<0.) multN++;
        }
      } 
    }
    else {
      if( fCutsRC->AcceptTrack(track) ) { 
	//Printf("Still here for HLT");
	mult++;
	if(track->Charge()>0.) multP++;
	if(track->Charge()<0.) multN++;
      }
    }
  }

  Double_t vTPCEvent[7] = {vtxESD->GetX(),vtxESD->GetY(),vtxESD->GetZ(),static_cast<Double_t>(mult),static_cast<Double_t>(multP),static_cast<Double_t>(multN),static_cast<Double_t>(vtxESD->GetStatus())};
  fTPCEventHisto->Fill(vTPCEvent);
}


//_____________________________________________________________________________
void AliPerformanceTPC::Analyse()
{
    //
    // Analyse comparison information and store output histograms
    // in the folder "folderTPC"
    //
    TH1::AddDirectory(kFALSE);
    TH1::SetDefaultSumw2(kFALSE);
    TObjArray *aFolderObj = new TObjArray;
    //aFolderObj->SetOwner(); // objects are owned by fanalysisFolder
    TString selString;

    //
    // Cluster histograms
    //
    AddProjection(aFolderObj, "clust", fTPCClustHisto, 0, 1, 2);
    
    selString = "all";
    for(Int_t i=0; i <= 2; i++) {
        AddProjection(aFolderObj, "clust", fTPCClustHisto, i, &selString);
    }
    
    //fTPCClustHisto->GetAxis(2)->SetRange(1,1); // A-side
    //selString = "A_side";
    //AddProjection(aFolderObj, fTPCClustHisto, 0, 1, &selString);
    
    //fTPCClustHisto->GetAxis(2)->SetRange(2,2); // C-side
    //selString = "C_side";
    //AddProjection(aFolderObj, fTPCClustHisto, 0, 1, &selString);
    
    //reset range
    fTPCClustHisto->GetAxis(2)->SetRange(1,2); 
    
    //
    // event histograms
    //
    for(Int_t i=0; i<=6; i++) {
      AddProjection(aFolderObj, "event", fTPCEventHisto, i);
    }    
    AddProjection(aFolderObj, "event", fTPCEventHisto, 4, 5);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 1);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 3);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 1, 3);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 2, 3);

    // reconstructed vertex status > 0
    fTPCEventHisto->GetAxis(6)->SetRange(2,2);
    selString = "recVertex";
    for(Int_t i=0; i<=5; i++) {
      AddProjection(aFolderObj, "event", fTPCEventHisto, i, &selString);
    }
    AddProjection(aFolderObj, "event", fTPCEventHisto, 4, 5, &selString);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 1, &selString);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 0, 3, &selString);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 1, 3, &selString);
    AddProjection(aFolderObj, "event", fTPCEventHisto, 2, 3, &selString);

    // reset cuts
    fTPCEventHisto->GetAxis(6)->SetRange(1,2);

    //
    // Track histograms 
    // 
    // all with vertex
    fTPCTrackHisto->GetAxis(8)->SetRangeUser(-1.5,1.5);
    fTPCTrackHisto->GetAxis(9)->SetRangeUser(0.5,1.5);
    selString = "all_recVertex";
    for(Int_t i=0; i <= 9; i++) {
        AddProjection(aFolderObj, "track", fTPCTrackHisto, i, &selString);        
    }

     AddProjection(aFolderObj, "track", fTPCTrackHisto, 5, 8, &selString); 

    for(Int_t i=0; i <= 4; i++) {
        AddProjection(aFolderObj, "track", fTPCTrackHisto, i, 5, 7, &selString);        
    }    



    // Track histograms (pos with vertex)
    fTPCTrackHisto->GetAxis(8)->SetRangeUser(0,1.5);
    selString = "pos_recVertex";
    for(Int_t i=0; i <= 9; i++) {
        AddProjection(aFolderObj, "track", fTPCTrackHisto, i, &selString);
    }
    for(Int_t i=0; i <= 4; i++) { for(Int_t j=5; j <= 5; j++) { for(Int_t k=j+1; k <= 7; k++) {
        AddProjection(aFolderObj, "track", fTPCTrackHisto, i, j, k, &selString);
    }  }  }
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 2, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 2, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 1, 2, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 3, 4, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 5, 6, 7, &selString);
  
    // Track histograms (neg with vertex)
    fTPCTrackHisto->GetAxis(8)->SetRangeUser(-1.5,0);
    selString = "neg_recVertex";
    for(Int_t i=0; i <= 9; i++) {
        AddProjection(aFolderObj, "track", fTPCTrackHisto, i, &selString);
    }
    for(Int_t i=0; i <= 4; i++) { for(Int_t j=5; j <= 5; j++) { for(Int_t k=j+1; k <= 7; k++) {
        AddProjection(aFolderObj, "track", fTPCTrackHisto, i, j, k, &selString);
    }  }  }
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 2, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 1, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 0, 2, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 1, 2, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 3, 4, 5, &selString);
    AddProjection(aFolderObj, "track", fTPCTrackHisto, 5, 6, 7, &selString);

    //restore cuts
    fTPCTrackHisto->GetAxis(8)->SetRangeUser(-1.5,1.5);
    fTPCTrackHisto->GetAxis(9)->SetRangeUser(-0.5,1.5);
  
  
    printf("exportToFolder\n");
    // export objects to analysis folder
    fAnalysisFolder = ExportToFolder(aFolderObj);
    if (fFolderObj) delete fFolderObj;
    fFolderObj = aFolderObj;
    aFolderObj=0;
}


//_____________________________________________________________________________
TFolder* AliPerformanceTPC::ExportToFolder(TObjArray * array) 
{
  // recreate folder avery time and export objects to new one
  //
  AliPerformanceTPC * comp=this;
  TFolder *folder = comp->GetAnalysisFolder();

  TString name, title;
  TFolder *newFolder = 0;
  Int_t i = 0;
  Int_t size = array->GetSize();

  if(folder) { 
     // get name and title from old folder
     name = folder->GetName();  
     title = folder->GetTitle();  

	 // delete old one
     delete folder;

	 // create new one
     newFolder = CreateFolder(name.Data(),title.Data());
     newFolder->SetOwner();

	 // add objects to folder
     while(i < size) {
	   newFolder->Add(array->At(i));
	   i++;
	 }
  }

return newFolder;
}

//_____________________________________________________________________________
Long64_t AliPerformanceTPC::Merge(TCollection* const list) 
{
  // Merge list of objects (needed by PROOF)

  if (!list)
  return 0;

  if (list->IsEmpty())
  return 1;
  
  Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));

  TIterator* iter = list->MakeIterator();
  TObject* obj = 0;
  TObjArray* objArrayList = 0;
  objArrayList = new TObjArray();

  // collection of generated histograms
  Int_t count=0;
  while((obj = iter->Next()) != 0) 
  {
    AliPerformanceTPC* entry = dynamic_cast<AliPerformanceTPC*>(obj);
    if (entry == 0) continue; 
    if (merge) {
        if ((fTPCClustHisto) && (entry->fTPCClustHisto)) { fTPCClustHisto->Add(entry->fTPCClustHisto); }
        if ((fTPCEventHisto) && (entry->fTPCEventHisto)) { fTPCEventHisto->Add(entry->fTPCEventHisto); }
        if ((fTPCTrackHisto) && (entry->fTPCTrackHisto)) { fTPCTrackHisto->Add(entry->fTPCTrackHisto); }
    }
    // the analysisfolder is only merged if present
    if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }

    count++;
  }
  if (fFolderObj) { fFolderObj->Merge(objArrayList); } 
  // to signal that track histos were not merged: reset
  if (!merge) { fTPCTrackHisto->Reset(); fTPCClustHisto->Reset(); fTPCEventHisto->Reset(); }
  // delete
  if (objArrayList)  delete objArrayList;  objArrayList=0;
return count;
}


//_____________________________________________________________________________
TFolder* AliPerformanceTPC::CreateFolder(TString name, TString title) 
{ 
// create folder for analysed histograms
//
TFolder *folder = 0;
  folder = new TFolder(name.Data(),title.Data());

  return folder;
}

//_____________________________________________________________________________
TTree* AliPerformanceTPC::CreateSummary()
{
    // implementaion removed, switched back to use AliPerformanceSummary (now called in AliPerformanceTask)
    return 0;
}

 AliPerformanceTPC.cxx:1
 AliPerformanceTPC.cxx:2
 AliPerformanceTPC.cxx:3
 AliPerformanceTPC.cxx:4
 AliPerformanceTPC.cxx:5
 AliPerformanceTPC.cxx:6
 AliPerformanceTPC.cxx:7
 AliPerformanceTPC.cxx:8
 AliPerformanceTPC.cxx:9
 AliPerformanceTPC.cxx:10
 AliPerformanceTPC.cxx:11
 AliPerformanceTPC.cxx:12
 AliPerformanceTPC.cxx:13
 AliPerformanceTPC.cxx:14
 AliPerformanceTPC.cxx:15
 AliPerformanceTPC.cxx:16
 AliPerformanceTPC.cxx:17
 AliPerformanceTPC.cxx:18
 AliPerformanceTPC.cxx:19
 AliPerformanceTPC.cxx:20
 AliPerformanceTPC.cxx:21
 AliPerformanceTPC.cxx:22
 AliPerformanceTPC.cxx:23
 AliPerformanceTPC.cxx:24
 AliPerformanceTPC.cxx:25
 AliPerformanceTPC.cxx:26
 AliPerformanceTPC.cxx:27
 AliPerformanceTPC.cxx:28
 AliPerformanceTPC.cxx:29
 AliPerformanceTPC.cxx:30
 AliPerformanceTPC.cxx:31
 AliPerformanceTPC.cxx:32
 AliPerformanceTPC.cxx:33
 AliPerformanceTPC.cxx:34
 AliPerformanceTPC.cxx:35
 AliPerformanceTPC.cxx:36
 AliPerformanceTPC.cxx:37
 AliPerformanceTPC.cxx:38
 AliPerformanceTPC.cxx:39
 AliPerformanceTPC.cxx:40
 AliPerformanceTPC.cxx:41
 AliPerformanceTPC.cxx:42
 AliPerformanceTPC.cxx:43
 AliPerformanceTPC.cxx:44
 AliPerformanceTPC.cxx:45
 AliPerformanceTPC.cxx:46
 AliPerformanceTPC.cxx:47
 AliPerformanceTPC.cxx:48
 AliPerformanceTPC.cxx:49
 AliPerformanceTPC.cxx:50
 AliPerformanceTPC.cxx:51
 AliPerformanceTPC.cxx:52
 AliPerformanceTPC.cxx:53
 AliPerformanceTPC.cxx:54
 AliPerformanceTPC.cxx:55
 AliPerformanceTPC.cxx:56
 AliPerformanceTPC.cxx:57
 AliPerformanceTPC.cxx:58
 AliPerformanceTPC.cxx:59
 AliPerformanceTPC.cxx:60
 AliPerformanceTPC.cxx:61
 AliPerformanceTPC.cxx:62
 AliPerformanceTPC.cxx:63
 AliPerformanceTPC.cxx:64
 AliPerformanceTPC.cxx:65
 AliPerformanceTPC.cxx:66
 AliPerformanceTPC.cxx:67
 AliPerformanceTPC.cxx:68
 AliPerformanceTPC.cxx:69
 AliPerformanceTPC.cxx:70
 AliPerformanceTPC.cxx:71
 AliPerformanceTPC.cxx:72
 AliPerformanceTPC.cxx:73
 AliPerformanceTPC.cxx:74
 AliPerformanceTPC.cxx:75
 AliPerformanceTPC.cxx:76
 AliPerformanceTPC.cxx:77
 AliPerformanceTPC.cxx:78
 AliPerformanceTPC.cxx:79
 AliPerformanceTPC.cxx:80
 AliPerformanceTPC.cxx:81
 AliPerformanceTPC.cxx:82
 AliPerformanceTPC.cxx:83
 AliPerformanceTPC.cxx:84
 AliPerformanceTPC.cxx:85
 AliPerformanceTPC.cxx:86
 AliPerformanceTPC.cxx:87
 AliPerformanceTPC.cxx:88
 AliPerformanceTPC.cxx:89
 AliPerformanceTPC.cxx:90
 AliPerformanceTPC.cxx:91
 AliPerformanceTPC.cxx:92
 AliPerformanceTPC.cxx:93
 AliPerformanceTPC.cxx:94
 AliPerformanceTPC.cxx:95
 AliPerformanceTPC.cxx:96
 AliPerformanceTPC.cxx:97
 AliPerformanceTPC.cxx:98
 AliPerformanceTPC.cxx:99
 AliPerformanceTPC.cxx:100
 AliPerformanceTPC.cxx:101
 AliPerformanceTPC.cxx:102
 AliPerformanceTPC.cxx:103
 AliPerformanceTPC.cxx:104
 AliPerformanceTPC.cxx:105
 AliPerformanceTPC.cxx:106
 AliPerformanceTPC.cxx:107
 AliPerformanceTPC.cxx:108
 AliPerformanceTPC.cxx:109
 AliPerformanceTPC.cxx:110
 AliPerformanceTPC.cxx:111
 AliPerformanceTPC.cxx:112
 AliPerformanceTPC.cxx:113
 AliPerformanceTPC.cxx:114
 AliPerformanceTPC.cxx:115
 AliPerformanceTPC.cxx:116
 AliPerformanceTPC.cxx:117
 AliPerformanceTPC.cxx:118
 AliPerformanceTPC.cxx:119
 AliPerformanceTPC.cxx:120
 AliPerformanceTPC.cxx:121
 AliPerformanceTPC.cxx:122
 AliPerformanceTPC.cxx:123
 AliPerformanceTPC.cxx:124
 AliPerformanceTPC.cxx:125
 AliPerformanceTPC.cxx:126
 AliPerformanceTPC.cxx:127
 AliPerformanceTPC.cxx:128
 AliPerformanceTPC.cxx:129
 AliPerformanceTPC.cxx:130
 AliPerformanceTPC.cxx:131
 AliPerformanceTPC.cxx:132
 AliPerformanceTPC.cxx:133
 AliPerformanceTPC.cxx:134
 AliPerformanceTPC.cxx:135
 AliPerformanceTPC.cxx:136
 AliPerformanceTPC.cxx:137
 AliPerformanceTPC.cxx:138
 AliPerformanceTPC.cxx:139
 AliPerformanceTPC.cxx:140
 AliPerformanceTPC.cxx:141
 AliPerformanceTPC.cxx:142
 AliPerformanceTPC.cxx:143
 AliPerformanceTPC.cxx:144
 AliPerformanceTPC.cxx:145
 AliPerformanceTPC.cxx:146
 AliPerformanceTPC.cxx:147
 AliPerformanceTPC.cxx:148
 AliPerformanceTPC.cxx:149
 AliPerformanceTPC.cxx:150
 AliPerformanceTPC.cxx:151
 AliPerformanceTPC.cxx:152
 AliPerformanceTPC.cxx:153
 AliPerformanceTPC.cxx:154
 AliPerformanceTPC.cxx:155
 AliPerformanceTPC.cxx:156
 AliPerformanceTPC.cxx:157
 AliPerformanceTPC.cxx:158
 AliPerformanceTPC.cxx:159
 AliPerformanceTPC.cxx:160
 AliPerformanceTPC.cxx:161
 AliPerformanceTPC.cxx:162
 AliPerformanceTPC.cxx:163
 AliPerformanceTPC.cxx:164
 AliPerformanceTPC.cxx:165
 AliPerformanceTPC.cxx:166
 AliPerformanceTPC.cxx:167
 AliPerformanceTPC.cxx:168
 AliPerformanceTPC.cxx:169
 AliPerformanceTPC.cxx:170
 AliPerformanceTPC.cxx:171
 AliPerformanceTPC.cxx:172
 AliPerformanceTPC.cxx:173
 AliPerformanceTPC.cxx:174
 AliPerformanceTPC.cxx:175
 AliPerformanceTPC.cxx:176
 AliPerformanceTPC.cxx:177
 AliPerformanceTPC.cxx:178
 AliPerformanceTPC.cxx:179
 AliPerformanceTPC.cxx:180
 AliPerformanceTPC.cxx:181
 AliPerformanceTPC.cxx:182
 AliPerformanceTPC.cxx:183
 AliPerformanceTPC.cxx:184
 AliPerformanceTPC.cxx:185
 AliPerformanceTPC.cxx:186
 AliPerformanceTPC.cxx:187
 AliPerformanceTPC.cxx:188
 AliPerformanceTPC.cxx:189
 AliPerformanceTPC.cxx:190
 AliPerformanceTPC.cxx:191
 AliPerformanceTPC.cxx:192
 AliPerformanceTPC.cxx:193
 AliPerformanceTPC.cxx:194
 AliPerformanceTPC.cxx:195
 AliPerformanceTPC.cxx:196
 AliPerformanceTPC.cxx:197
 AliPerformanceTPC.cxx:198
 AliPerformanceTPC.cxx:199
 AliPerformanceTPC.cxx:200
 AliPerformanceTPC.cxx:201
 AliPerformanceTPC.cxx:202
 AliPerformanceTPC.cxx:203
 AliPerformanceTPC.cxx:204
 AliPerformanceTPC.cxx:205
 AliPerformanceTPC.cxx:206
 AliPerformanceTPC.cxx:207
 AliPerformanceTPC.cxx:208
 AliPerformanceTPC.cxx:209
 AliPerformanceTPC.cxx:210
 AliPerformanceTPC.cxx:211
 AliPerformanceTPC.cxx:212
 AliPerformanceTPC.cxx:213
 AliPerformanceTPC.cxx:214
 AliPerformanceTPC.cxx:215
 AliPerformanceTPC.cxx:216
 AliPerformanceTPC.cxx:217
 AliPerformanceTPC.cxx:218
 AliPerformanceTPC.cxx:219
 AliPerformanceTPC.cxx:220
 AliPerformanceTPC.cxx:221
 AliPerformanceTPC.cxx:222
 AliPerformanceTPC.cxx:223
 AliPerformanceTPC.cxx:224
 AliPerformanceTPC.cxx:225
 AliPerformanceTPC.cxx:226
 AliPerformanceTPC.cxx:227
 AliPerformanceTPC.cxx:228
 AliPerformanceTPC.cxx:229
 AliPerformanceTPC.cxx:230
 AliPerformanceTPC.cxx:231
 AliPerformanceTPC.cxx:232
 AliPerformanceTPC.cxx:233
 AliPerformanceTPC.cxx:234
 AliPerformanceTPC.cxx:235
 AliPerformanceTPC.cxx:236
 AliPerformanceTPC.cxx:237
 AliPerformanceTPC.cxx:238
 AliPerformanceTPC.cxx:239
 AliPerformanceTPC.cxx:240
 AliPerformanceTPC.cxx:241
 AliPerformanceTPC.cxx:242
 AliPerformanceTPC.cxx:243
 AliPerformanceTPC.cxx:244
 AliPerformanceTPC.cxx:245
 AliPerformanceTPC.cxx:246
 AliPerformanceTPC.cxx:247
 AliPerformanceTPC.cxx:248
 AliPerformanceTPC.cxx:249
 AliPerformanceTPC.cxx:250
 AliPerformanceTPC.cxx:251
 AliPerformanceTPC.cxx:252
 AliPerformanceTPC.cxx:253
 AliPerformanceTPC.cxx:254
 AliPerformanceTPC.cxx:255
 AliPerformanceTPC.cxx:256
 AliPerformanceTPC.cxx:257
 AliPerformanceTPC.cxx:258
 AliPerformanceTPC.cxx:259
 AliPerformanceTPC.cxx:260
 AliPerformanceTPC.cxx:261
 AliPerformanceTPC.cxx:262
 AliPerformanceTPC.cxx:263
 AliPerformanceTPC.cxx:264
 AliPerformanceTPC.cxx:265
 AliPerformanceTPC.cxx:266
 AliPerformanceTPC.cxx:267
 AliPerformanceTPC.cxx:268
 AliPerformanceTPC.cxx:269
 AliPerformanceTPC.cxx:270
 AliPerformanceTPC.cxx:271
 AliPerformanceTPC.cxx:272
 AliPerformanceTPC.cxx:273
 AliPerformanceTPC.cxx:274
 AliPerformanceTPC.cxx:275
 AliPerformanceTPC.cxx:276
 AliPerformanceTPC.cxx:277
 AliPerformanceTPC.cxx:278
 AliPerformanceTPC.cxx:279
 AliPerformanceTPC.cxx:280
 AliPerformanceTPC.cxx:281
 AliPerformanceTPC.cxx:282
 AliPerformanceTPC.cxx:283
 AliPerformanceTPC.cxx:284
 AliPerformanceTPC.cxx:285
 AliPerformanceTPC.cxx:286
 AliPerformanceTPC.cxx:287
 AliPerformanceTPC.cxx:288
 AliPerformanceTPC.cxx:289
 AliPerformanceTPC.cxx:290
 AliPerformanceTPC.cxx:291
 AliPerformanceTPC.cxx:292
 AliPerformanceTPC.cxx:293
 AliPerformanceTPC.cxx:294
 AliPerformanceTPC.cxx:295
 AliPerformanceTPC.cxx:296
 AliPerformanceTPC.cxx:297
 AliPerformanceTPC.cxx:298
 AliPerformanceTPC.cxx:299
 AliPerformanceTPC.cxx:300
 AliPerformanceTPC.cxx:301
 AliPerformanceTPC.cxx:302
 AliPerformanceTPC.cxx:303
 AliPerformanceTPC.cxx:304
 AliPerformanceTPC.cxx:305
 AliPerformanceTPC.cxx:306
 AliPerformanceTPC.cxx:307
 AliPerformanceTPC.cxx:308
 AliPerformanceTPC.cxx:309
 AliPerformanceTPC.cxx:310
 AliPerformanceTPC.cxx:311
 AliPerformanceTPC.cxx:312
 AliPerformanceTPC.cxx:313
 AliPerformanceTPC.cxx:314
 AliPerformanceTPC.cxx:315
 AliPerformanceTPC.cxx:316
 AliPerformanceTPC.cxx:317
 AliPerformanceTPC.cxx:318
 AliPerformanceTPC.cxx:319
 AliPerformanceTPC.cxx:320
 AliPerformanceTPC.cxx:321
 AliPerformanceTPC.cxx:322
 AliPerformanceTPC.cxx:323
 AliPerformanceTPC.cxx:324
 AliPerformanceTPC.cxx:325
 AliPerformanceTPC.cxx:326
 AliPerformanceTPC.cxx:327
 AliPerformanceTPC.cxx:328
 AliPerformanceTPC.cxx:329
 AliPerformanceTPC.cxx:330
 AliPerformanceTPC.cxx:331
 AliPerformanceTPC.cxx:332
 AliPerformanceTPC.cxx:333
 AliPerformanceTPC.cxx:334
 AliPerformanceTPC.cxx:335
 AliPerformanceTPC.cxx:336
 AliPerformanceTPC.cxx:337
 AliPerformanceTPC.cxx:338
 AliPerformanceTPC.cxx:339
 AliPerformanceTPC.cxx:340
 AliPerformanceTPC.cxx:341
 AliPerformanceTPC.cxx:342
 AliPerformanceTPC.cxx:343
 AliPerformanceTPC.cxx:344
 AliPerformanceTPC.cxx:345
 AliPerformanceTPC.cxx:346
 AliPerformanceTPC.cxx:347
 AliPerformanceTPC.cxx:348
 AliPerformanceTPC.cxx:349
 AliPerformanceTPC.cxx:350
 AliPerformanceTPC.cxx:351
 AliPerformanceTPC.cxx:352
 AliPerformanceTPC.cxx:353
 AliPerformanceTPC.cxx:354
 AliPerformanceTPC.cxx:355
 AliPerformanceTPC.cxx:356
 AliPerformanceTPC.cxx:357
 AliPerformanceTPC.cxx:358
 AliPerformanceTPC.cxx:359
 AliPerformanceTPC.cxx:360
 AliPerformanceTPC.cxx:361
 AliPerformanceTPC.cxx:362
 AliPerformanceTPC.cxx:363
 AliPerformanceTPC.cxx:364
 AliPerformanceTPC.cxx:365
 AliPerformanceTPC.cxx:366
 AliPerformanceTPC.cxx:367
 AliPerformanceTPC.cxx:368
 AliPerformanceTPC.cxx:369
 AliPerformanceTPC.cxx:370
 AliPerformanceTPC.cxx:371
 AliPerformanceTPC.cxx:372
 AliPerformanceTPC.cxx:373
 AliPerformanceTPC.cxx:374
 AliPerformanceTPC.cxx:375
 AliPerformanceTPC.cxx:376
 AliPerformanceTPC.cxx:377
 AliPerformanceTPC.cxx:378
 AliPerformanceTPC.cxx:379
 AliPerformanceTPC.cxx:380
 AliPerformanceTPC.cxx:381
 AliPerformanceTPC.cxx:382
 AliPerformanceTPC.cxx:383
 AliPerformanceTPC.cxx:384
 AliPerformanceTPC.cxx:385
 AliPerformanceTPC.cxx:386
 AliPerformanceTPC.cxx:387
 AliPerformanceTPC.cxx:388
 AliPerformanceTPC.cxx:389
 AliPerformanceTPC.cxx:390
 AliPerformanceTPC.cxx:391
 AliPerformanceTPC.cxx:392
 AliPerformanceTPC.cxx:393
 AliPerformanceTPC.cxx:394
 AliPerformanceTPC.cxx:395
 AliPerformanceTPC.cxx:396
 AliPerformanceTPC.cxx:397
 AliPerformanceTPC.cxx:398
 AliPerformanceTPC.cxx:399
 AliPerformanceTPC.cxx:400
 AliPerformanceTPC.cxx:401
 AliPerformanceTPC.cxx:402
 AliPerformanceTPC.cxx:403
 AliPerformanceTPC.cxx:404
 AliPerformanceTPC.cxx:405
 AliPerformanceTPC.cxx:406
 AliPerformanceTPC.cxx:407
 AliPerformanceTPC.cxx:408
 AliPerformanceTPC.cxx:409
 AliPerformanceTPC.cxx:410
 AliPerformanceTPC.cxx:411
 AliPerformanceTPC.cxx:412
 AliPerformanceTPC.cxx:413
 AliPerformanceTPC.cxx:414
 AliPerformanceTPC.cxx:415
 AliPerformanceTPC.cxx:416
 AliPerformanceTPC.cxx:417
 AliPerformanceTPC.cxx:418
 AliPerformanceTPC.cxx:419
 AliPerformanceTPC.cxx:420
 AliPerformanceTPC.cxx:421
 AliPerformanceTPC.cxx:422
 AliPerformanceTPC.cxx:423
 AliPerformanceTPC.cxx:424
 AliPerformanceTPC.cxx:425
 AliPerformanceTPC.cxx:426
 AliPerformanceTPC.cxx:427
 AliPerformanceTPC.cxx:428
 AliPerformanceTPC.cxx:429
 AliPerformanceTPC.cxx:430
 AliPerformanceTPC.cxx:431
 AliPerformanceTPC.cxx:432
 AliPerformanceTPC.cxx:433
 AliPerformanceTPC.cxx:434
 AliPerformanceTPC.cxx:435
 AliPerformanceTPC.cxx:436
 AliPerformanceTPC.cxx:437
 AliPerformanceTPC.cxx:438
 AliPerformanceTPC.cxx:439
 AliPerformanceTPC.cxx:440
 AliPerformanceTPC.cxx:441
 AliPerformanceTPC.cxx:442
 AliPerformanceTPC.cxx:443
 AliPerformanceTPC.cxx:444
 AliPerformanceTPC.cxx:445
 AliPerformanceTPC.cxx:446
 AliPerformanceTPC.cxx:447
 AliPerformanceTPC.cxx:448
 AliPerformanceTPC.cxx:449
 AliPerformanceTPC.cxx:450
 AliPerformanceTPC.cxx:451
 AliPerformanceTPC.cxx:452
 AliPerformanceTPC.cxx:453
 AliPerformanceTPC.cxx:454
 AliPerformanceTPC.cxx:455
 AliPerformanceTPC.cxx:456
 AliPerformanceTPC.cxx:457
 AliPerformanceTPC.cxx:458
 AliPerformanceTPC.cxx:459
 AliPerformanceTPC.cxx:460
 AliPerformanceTPC.cxx:461
 AliPerformanceTPC.cxx:462
 AliPerformanceTPC.cxx:463
 AliPerformanceTPC.cxx:464
 AliPerformanceTPC.cxx:465
 AliPerformanceTPC.cxx:466
 AliPerformanceTPC.cxx:467
 AliPerformanceTPC.cxx:468
 AliPerformanceTPC.cxx:469
 AliPerformanceTPC.cxx:470
 AliPerformanceTPC.cxx:471
 AliPerformanceTPC.cxx:472
 AliPerformanceTPC.cxx:473
 AliPerformanceTPC.cxx:474
 AliPerformanceTPC.cxx:475
 AliPerformanceTPC.cxx:476
 AliPerformanceTPC.cxx:477
 AliPerformanceTPC.cxx:478
 AliPerformanceTPC.cxx:479
 AliPerformanceTPC.cxx:480
 AliPerformanceTPC.cxx:481
 AliPerformanceTPC.cxx:482
 AliPerformanceTPC.cxx:483
 AliPerformanceTPC.cxx:484
 AliPerformanceTPC.cxx:485
 AliPerformanceTPC.cxx:486
 AliPerformanceTPC.cxx:487
 AliPerformanceTPC.cxx:488
 AliPerformanceTPC.cxx:489
 AliPerformanceTPC.cxx:490
 AliPerformanceTPC.cxx:491
 AliPerformanceTPC.cxx:492
 AliPerformanceTPC.cxx:493
 AliPerformanceTPC.cxx:494
 AliPerformanceTPC.cxx:495
 AliPerformanceTPC.cxx:496
 AliPerformanceTPC.cxx:497
 AliPerformanceTPC.cxx:498
 AliPerformanceTPC.cxx:499
 AliPerformanceTPC.cxx:500
 AliPerformanceTPC.cxx:501
 AliPerformanceTPC.cxx:502
 AliPerformanceTPC.cxx:503
 AliPerformanceTPC.cxx:504
 AliPerformanceTPC.cxx:505
 AliPerformanceTPC.cxx:506
 AliPerformanceTPC.cxx:507
 AliPerformanceTPC.cxx:508
 AliPerformanceTPC.cxx:509
 AliPerformanceTPC.cxx:510
 AliPerformanceTPC.cxx:511
 AliPerformanceTPC.cxx:512
 AliPerformanceTPC.cxx:513
 AliPerformanceTPC.cxx:514
 AliPerformanceTPC.cxx:515
 AliPerformanceTPC.cxx:516
 AliPerformanceTPC.cxx:517
 AliPerformanceTPC.cxx:518
 AliPerformanceTPC.cxx:519
 AliPerformanceTPC.cxx:520
 AliPerformanceTPC.cxx:521
 AliPerformanceTPC.cxx:522
 AliPerformanceTPC.cxx:523
 AliPerformanceTPC.cxx:524
 AliPerformanceTPC.cxx:525
 AliPerformanceTPC.cxx:526
 AliPerformanceTPC.cxx:527
 AliPerformanceTPC.cxx:528
 AliPerformanceTPC.cxx:529
 AliPerformanceTPC.cxx:530
 AliPerformanceTPC.cxx:531
 AliPerformanceTPC.cxx:532
 AliPerformanceTPC.cxx:533
 AliPerformanceTPC.cxx:534
 AliPerformanceTPC.cxx:535
 AliPerformanceTPC.cxx:536
 AliPerformanceTPC.cxx:537
 AliPerformanceTPC.cxx:538
 AliPerformanceTPC.cxx:539
 AliPerformanceTPC.cxx:540
 AliPerformanceTPC.cxx:541
 AliPerformanceTPC.cxx:542
 AliPerformanceTPC.cxx:543
 AliPerformanceTPC.cxx:544
 AliPerformanceTPC.cxx:545
 AliPerformanceTPC.cxx:546
 AliPerformanceTPC.cxx:547
 AliPerformanceTPC.cxx:548
 AliPerformanceTPC.cxx:549
 AliPerformanceTPC.cxx:550
 AliPerformanceTPC.cxx:551
 AliPerformanceTPC.cxx:552
 AliPerformanceTPC.cxx:553
 AliPerformanceTPC.cxx:554
 AliPerformanceTPC.cxx:555
 AliPerformanceTPC.cxx:556
 AliPerformanceTPC.cxx:557
 AliPerformanceTPC.cxx:558
 AliPerformanceTPC.cxx:559
 AliPerformanceTPC.cxx:560
 AliPerformanceTPC.cxx:561
 AliPerformanceTPC.cxx:562
 AliPerformanceTPC.cxx:563
 AliPerformanceTPC.cxx:564
 AliPerformanceTPC.cxx:565
 AliPerformanceTPC.cxx:566
 AliPerformanceTPC.cxx:567
 AliPerformanceTPC.cxx:568
 AliPerformanceTPC.cxx:569
 AliPerformanceTPC.cxx:570
 AliPerformanceTPC.cxx:571
 AliPerformanceTPC.cxx:572
 AliPerformanceTPC.cxx:573
 AliPerformanceTPC.cxx:574
 AliPerformanceTPC.cxx:575
 AliPerformanceTPC.cxx:576
 AliPerformanceTPC.cxx:577
 AliPerformanceTPC.cxx:578
 AliPerformanceTPC.cxx:579
 AliPerformanceTPC.cxx:580
 AliPerformanceTPC.cxx:581
 AliPerformanceTPC.cxx:582
 AliPerformanceTPC.cxx:583
 AliPerformanceTPC.cxx:584
 AliPerformanceTPC.cxx:585
 AliPerformanceTPC.cxx:586
 AliPerformanceTPC.cxx:587
 AliPerformanceTPC.cxx:588
 AliPerformanceTPC.cxx:589
 AliPerformanceTPC.cxx:590
 AliPerformanceTPC.cxx:591
 AliPerformanceTPC.cxx:592
 AliPerformanceTPC.cxx:593
 AliPerformanceTPC.cxx:594
 AliPerformanceTPC.cxx:595
 AliPerformanceTPC.cxx:596
 AliPerformanceTPC.cxx:597
 AliPerformanceTPC.cxx:598
 AliPerformanceTPC.cxx:599
 AliPerformanceTPC.cxx:600
 AliPerformanceTPC.cxx:601
 AliPerformanceTPC.cxx:602
 AliPerformanceTPC.cxx:603
 AliPerformanceTPC.cxx:604
 AliPerformanceTPC.cxx:605
 AliPerformanceTPC.cxx:606
 AliPerformanceTPC.cxx:607
 AliPerformanceTPC.cxx:608
 AliPerformanceTPC.cxx:609
 AliPerformanceTPC.cxx:610
 AliPerformanceTPC.cxx:611
 AliPerformanceTPC.cxx:612
 AliPerformanceTPC.cxx:613
 AliPerformanceTPC.cxx:614
 AliPerformanceTPC.cxx:615
 AliPerformanceTPC.cxx:616
 AliPerformanceTPC.cxx:617
 AliPerformanceTPC.cxx:618
 AliPerformanceTPC.cxx:619
 AliPerformanceTPC.cxx:620
 AliPerformanceTPC.cxx:621
 AliPerformanceTPC.cxx:622
 AliPerformanceTPC.cxx:623
 AliPerformanceTPC.cxx:624
 AliPerformanceTPC.cxx:625
 AliPerformanceTPC.cxx:626
 AliPerformanceTPC.cxx:627
 AliPerformanceTPC.cxx:628
 AliPerformanceTPC.cxx:629
 AliPerformanceTPC.cxx:630
 AliPerformanceTPC.cxx:631
 AliPerformanceTPC.cxx:632
 AliPerformanceTPC.cxx:633
 AliPerformanceTPC.cxx:634
 AliPerformanceTPC.cxx:635
 AliPerformanceTPC.cxx:636
 AliPerformanceTPC.cxx:637
 AliPerformanceTPC.cxx:638
 AliPerformanceTPC.cxx:639
 AliPerformanceTPC.cxx:640
 AliPerformanceTPC.cxx:641
 AliPerformanceTPC.cxx:642
 AliPerformanceTPC.cxx:643
 AliPerformanceTPC.cxx:644
 AliPerformanceTPC.cxx:645
 AliPerformanceTPC.cxx:646
 AliPerformanceTPC.cxx:647
 AliPerformanceTPC.cxx:648
 AliPerformanceTPC.cxx:649
 AliPerformanceTPC.cxx:650
 AliPerformanceTPC.cxx:651
 AliPerformanceTPC.cxx:652
 AliPerformanceTPC.cxx:653
 AliPerformanceTPC.cxx:654
 AliPerformanceTPC.cxx:655
 AliPerformanceTPC.cxx:656
 AliPerformanceTPC.cxx:657
 AliPerformanceTPC.cxx:658
 AliPerformanceTPC.cxx:659
 AliPerformanceTPC.cxx:660
 AliPerformanceTPC.cxx:661
 AliPerformanceTPC.cxx:662
 AliPerformanceTPC.cxx:663
 AliPerformanceTPC.cxx:664
 AliPerformanceTPC.cxx:665
 AliPerformanceTPC.cxx:666
 AliPerformanceTPC.cxx:667
 AliPerformanceTPC.cxx:668
 AliPerformanceTPC.cxx:669
 AliPerformanceTPC.cxx:670
 AliPerformanceTPC.cxx:671
 AliPerformanceTPC.cxx:672
 AliPerformanceTPC.cxx:673
 AliPerformanceTPC.cxx:674
 AliPerformanceTPC.cxx:675
 AliPerformanceTPC.cxx:676
 AliPerformanceTPC.cxx:677
 AliPerformanceTPC.cxx:678
 AliPerformanceTPC.cxx:679
 AliPerformanceTPC.cxx:680
 AliPerformanceTPC.cxx:681
 AliPerformanceTPC.cxx:682
 AliPerformanceTPC.cxx:683
 AliPerformanceTPC.cxx:684
 AliPerformanceTPC.cxx:685
 AliPerformanceTPC.cxx:686
 AliPerformanceTPC.cxx:687
 AliPerformanceTPC.cxx:688
 AliPerformanceTPC.cxx:689
 AliPerformanceTPC.cxx:690
 AliPerformanceTPC.cxx:691
 AliPerformanceTPC.cxx:692
 AliPerformanceTPC.cxx:693
 AliPerformanceTPC.cxx:694
 AliPerformanceTPC.cxx:695
 AliPerformanceTPC.cxx:696
 AliPerformanceTPC.cxx:697
 AliPerformanceTPC.cxx:698
 AliPerformanceTPC.cxx:699
 AliPerformanceTPC.cxx:700
 AliPerformanceTPC.cxx:701
 AliPerformanceTPC.cxx:702
 AliPerformanceTPC.cxx:703
 AliPerformanceTPC.cxx:704
 AliPerformanceTPC.cxx:705
 AliPerformanceTPC.cxx:706
 AliPerformanceTPC.cxx:707
 AliPerformanceTPC.cxx:708
 AliPerformanceTPC.cxx:709
 AliPerformanceTPC.cxx:710
 AliPerformanceTPC.cxx:711
 AliPerformanceTPC.cxx:712
 AliPerformanceTPC.cxx:713
 AliPerformanceTPC.cxx:714
 AliPerformanceTPC.cxx:715
 AliPerformanceTPC.cxx:716
 AliPerformanceTPC.cxx:717
 AliPerformanceTPC.cxx:718
 AliPerformanceTPC.cxx:719
 AliPerformanceTPC.cxx:720
 AliPerformanceTPC.cxx:721
 AliPerformanceTPC.cxx:722
 AliPerformanceTPC.cxx:723
 AliPerformanceTPC.cxx:724
 AliPerformanceTPC.cxx:725
 AliPerformanceTPC.cxx:726
 AliPerformanceTPC.cxx:727
 AliPerformanceTPC.cxx:728
 AliPerformanceTPC.cxx:729
 AliPerformanceTPC.cxx:730
 AliPerformanceTPC.cxx:731
 AliPerformanceTPC.cxx:732
 AliPerformanceTPC.cxx:733
 AliPerformanceTPC.cxx:734
 AliPerformanceTPC.cxx:735
 AliPerformanceTPC.cxx:736
 AliPerformanceTPC.cxx:737
 AliPerformanceTPC.cxx:738
 AliPerformanceTPC.cxx:739
 AliPerformanceTPC.cxx:740
 AliPerformanceTPC.cxx:741
 AliPerformanceTPC.cxx:742
 AliPerformanceTPC.cxx:743
 AliPerformanceTPC.cxx:744
 AliPerformanceTPC.cxx:745
 AliPerformanceTPC.cxx:746
 AliPerformanceTPC.cxx:747
 AliPerformanceTPC.cxx:748
 AliPerformanceTPC.cxx:749
 AliPerformanceTPC.cxx:750
 AliPerformanceTPC.cxx:751
 AliPerformanceTPC.cxx:752
 AliPerformanceTPC.cxx:753
 AliPerformanceTPC.cxx:754
 AliPerformanceTPC.cxx:755
 AliPerformanceTPC.cxx:756
 AliPerformanceTPC.cxx:757
 AliPerformanceTPC.cxx:758
 AliPerformanceTPC.cxx:759
 AliPerformanceTPC.cxx:760
 AliPerformanceTPC.cxx:761
 AliPerformanceTPC.cxx:762
 AliPerformanceTPC.cxx:763
 AliPerformanceTPC.cxx:764
 AliPerformanceTPC.cxx:765
 AliPerformanceTPC.cxx:766
 AliPerformanceTPC.cxx:767
 AliPerformanceTPC.cxx:768
 AliPerformanceTPC.cxx:769
 AliPerformanceTPC.cxx:770
 AliPerformanceTPC.cxx:771
 AliPerformanceTPC.cxx:772
 AliPerformanceTPC.cxx:773
 AliPerformanceTPC.cxx:774
 AliPerformanceTPC.cxx:775
 AliPerformanceTPC.cxx:776
 AliPerformanceTPC.cxx:777
 AliPerformanceTPC.cxx:778
 AliPerformanceTPC.cxx:779
 AliPerformanceTPC.cxx:780
 AliPerformanceTPC.cxx:781
 AliPerformanceTPC.cxx:782
 AliPerformanceTPC.cxx:783
 AliPerformanceTPC.cxx:784
 AliPerformanceTPC.cxx:785
 AliPerformanceTPC.cxx:786
 AliPerformanceTPC.cxx:787
 AliPerformanceTPC.cxx:788
 AliPerformanceTPC.cxx:789
 AliPerformanceTPC.cxx:790
 AliPerformanceTPC.cxx:791
 AliPerformanceTPC.cxx:792
 AliPerformanceTPC.cxx:793
 AliPerformanceTPC.cxx:794
 AliPerformanceTPC.cxx:795
 AliPerformanceTPC.cxx:796
 AliPerformanceTPC.cxx:797
 AliPerformanceTPC.cxx:798
 AliPerformanceTPC.cxx:799
 AliPerformanceTPC.cxx:800
 AliPerformanceTPC.cxx:801
 AliPerformanceTPC.cxx:802
 AliPerformanceTPC.cxx:803
 AliPerformanceTPC.cxx:804
 AliPerformanceTPC.cxx:805
 AliPerformanceTPC.cxx:806
 AliPerformanceTPC.cxx:807
 AliPerformanceTPC.cxx:808
 AliPerformanceTPC.cxx:809
 AliPerformanceTPC.cxx:810
 AliPerformanceTPC.cxx:811
 AliPerformanceTPC.cxx:812
 AliPerformanceTPC.cxx:813
 AliPerformanceTPC.cxx:814
 AliPerformanceTPC.cxx:815
 AliPerformanceTPC.cxx:816
 AliPerformanceTPC.cxx:817
 AliPerformanceTPC.cxx:818
 AliPerformanceTPC.cxx:819
 AliPerformanceTPC.cxx:820
 AliPerformanceTPC.cxx:821
 AliPerformanceTPC.cxx:822
 AliPerformanceTPC.cxx:823
 AliPerformanceTPC.cxx:824
 AliPerformanceTPC.cxx:825
 AliPerformanceTPC.cxx:826
 AliPerformanceTPC.cxx:827
 AliPerformanceTPC.cxx:828
 AliPerformanceTPC.cxx:829
 AliPerformanceTPC.cxx:830
 AliPerformanceTPC.cxx:831
 AliPerformanceTPC.cxx:832