ROOT logo
//------------------------------------------------------------------------------
// Implementation of AliPerformanceMatch class. It checks the track matching 
// quality (residuals, pulls) bewteen tracking detectors TPC-ITS and TPC-TRD 
// with TPC as the reference detector. In addition, the ITS and TRD track 
// reconstruction efficiency with respect to TPC is calculated.
//
// The matchinig quality parameters are stored in the THnSparse histograms. 
// The analysis of these histograms to extract reference information is done in 
// the Analyse() function (post-processing).  
//  
// The histograms with reference information can be exported to the ROOT folder.
//
// Author: J.Otwinowski 17/10/2009 
// Changes by M.Knichel 22/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");
  AliPerformanceMatch * compObj = (AliPerformanceMatch*)coutput->FindObject("AliPerformanceMatchTPCTRD");
 
  // analyse comparison data
  compObj->Analyse();

  // the output histograms/graphs will be stored in the folder "folderMatch" 
  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_Match.root","recreate");
  compObj->Write(); // compObj->GetAnalysisFolder()->Write();
  fout.Close();

*/

#include "TCanvas.h"
#include "TH1.h"
#include "TH2.h"
#include "TAxis.h"
#include "TF1.h"

#include "AliPerformanceMatch.h" 
#include "AliESDEvent.h" 
#include "AliESDVertex.h"
#include "AliESDtrack.h"
#include "AliESDfriendTrack.h"
#include "AliESDfriend.h"
#include "AliLog.h" 
#include "AliMCEvent.h" 
#include "AliMCParticle.h" 
#include "AliHeader.h" 
#include "AliGenEventHeader.h" 
#include "AliStack.h" 
#include "AliMCInfoCuts.h" 
#include "AliRecInfoCuts.h" 
#include "AliTracker.h" 
#include "AliTRDtrackV1.h" 
#include "AliTreeDraw.h" 

using namespace std;

ClassImp(AliPerformanceMatch)

Bool_t AliPerformanceMatch::fgMergeTHnSparse = kFALSE;
Bool_t AliPerformanceMatch::fgUseMergeTHnSparse = kFALSE;

//_____________________________________________________________________________
AliPerformanceMatch::AliPerformanceMatch(const Char_t* name, const Char_t* title, Int_t analysisMode, Bool_t hptGenerator):
  AliPerformanceObject(name,title),
  fResolHisto(0),
  fPullHisto(0),
  fTrackingEffHisto(0),
  fTPCConstrain(0),
  fFolderObj(0),
  // Cuts 
  fCutsRC(0),  
  fCutsMC(0),  

  // histogram folder 
  fAnalysisFolder(0),
  fUseHLT(0)
{
  // named constructor	
  // 
  SetAnalysisMode(analysisMode);
  SetHptGenerator(hptGenerator);

  Init();
}

//_____________________________________________________________________________
AliPerformanceMatch::~AliPerformanceMatch()
{
  // destructor
   
  if(fResolHisto) delete fResolHisto; fResolHisto=0;     
  if(fPullHisto)  delete fPullHisto;  fPullHisto=0;
  if(fTrackingEffHisto) delete fTrackingEffHisto; fTrackingEffHisto = 0x0;
  if(fTPCConstrain) delete fTPCConstrain; fTPCConstrain = 0x0;
  
  if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
  if(fFolderObj) delete fFolderObj; fFolderObj=0;
}

//_____________________________________________________________________________
void AliPerformanceMatch::Init(){

  //
  // Make performance histogrms
  //

  // 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);


  // Double_t yMin = -0.02, yMax = 0.02;
  // Double_t zMin = -12.0, zMax = 12.0;
  // if(GetAnalysisMode() == 0 || GetAnalysisMode() == 1) { 
  //   yMin = -100.; yMax = 100.; 
  //   zMin = -100.; zMax = 100.; 
  // }

  //
  //init ITS TPC Mactching
  //
  if(GetAnalysisMode()==1){
    // res_y:res_z:res_phi,res_lambda:res_pt:y:z:eta:phi:pt:isRec
    Int_t binsResolHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
    Double_t minResolHisto[9]={-1.,-1.,-0.03,-0.03,-0.2, 0., -1.5, ptMin,0};
    Double_t maxResolHisto[9]={ 1., 1., 0.03, 0.03, 0.2, 2.*TMath::Pi(), 1.5, ptMax,2};
    
    fResolHisto = new THnSparseF("fResolHisto","res_y:res_z:res_phi:res_lambda:res_pt:phi:eta:pt:isRec",9,binsResolHisto,minResolHisto,maxResolHisto);
    fResolHisto->SetBinEdges(7,binsPt);
    
    fResolHisto->GetAxis(0)->SetTitle("y-y_{ref} (cm)");
    fResolHisto->GetAxis(1)->SetTitle("z-z_{ref} (cm)");
    fResolHisto->GetAxis(2)->SetTitle("#phi-#phi_{ref} (rad)");
    fResolHisto->GetAxis(3)->SetTitle("#lambda-#lambda_{ref} (rad)");
    fResolHisto->GetAxis(4)->SetTitle("(p_{T}/p_{Tref}-1)");
    fResolHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
    fResolHisto->GetAxis(6)->SetTitle("#eta_{ref}");
    fResolHisto->GetAxis(7)->SetTitle("p_{Tref} (GeV/c)");
    fResolHisto->GetAxis(8)->SetTitle("isReconstructed");
    fResolHisto->Sumw2();
    
    //pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec
    Int_t binsPullHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
    Double_t minPullHisto[9]={-5.,-5.,-5.,-5.,-5., 0.,-1.5, ptMin,0};
    Double_t maxPullHisto[9]={ 5., 5., 5., 5., 5., 2.*TMath::Pi(), 1.5, ptMax,2};
    fPullHisto = new THnSparseF("fPullHisto","pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:phi:eta:1pt:isRec",9,binsPullHisto,minPullHisto,maxPullHisto);
    
    fPullHisto->GetAxis(0)->SetTitle("(y-y_{ref})/#sigma");
    fPullHisto->GetAxis(1)->SetTitle("(z-z_{ref})/#sigma");
    fPullHisto->GetAxis(2)->SetTitle("(sin#phi-sin#phi_{ref})/#sigma");
    fPullHisto->GetAxis(3)->SetTitle("(tan#lambda-tan#lambda_{ref})/#sigma");
    fPullHisto->GetAxis(4)->SetTitle("(p_{Tref}/p_{T}-1)/#sigma");
    fPullHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
    fPullHisto->GetAxis(6)->SetTitle("#eta_{ref}");
    fPullHisto->GetAxis(7)->SetTitle("1/p_{Tref} (GeV/c)^{-1}");
    fPullHisto->GetAxis(8)->SetTitle("isReconstructed");
    fPullHisto->Sumw2();
  }
  
  //
  //TPC  ITS  Mactching
  // 
  if(GetAnalysisMode()==0){
    // -> has match:y:z:snp:tgl:phi:pt:ITSclusters
    Int_t binsTrackingEffHisto[5]    = { 2,   90,       nPtBins, 30,    7   };
    Double_t minTrackingEffHisto[5]  = {-0.5, 0.,          ptMin, -1.5, -0.5 };
    Double_t maxTrackingEffHisto[5]  = { 1.5, 2*TMath::Pi(), ptMax, 1.5,  6.5 };
    
    fTrackingEffHisto = new THnSparseF("fTrackingEffHisto","has match:phi:pt:eta:ITSclusters",5,binsTrackingEffHisto,minTrackingEffHisto,maxTrackingEffHisto);
    fTrackingEffHisto->SetBinEdges(2,binsPt);
    fTrackingEffHisto->GetAxis(0)->SetTitle("IsMatching");
    fTrackingEffHisto->GetAxis(1)->SetTitle("phi (rad)");
    fTrackingEffHisto->GetAxis(2)->SetTitle("p_{T}");
    fTrackingEffHisto->GetAxis(3)->SetTitle("eta");
    fTrackingEffHisto->GetAxis(4)->SetTitle("number of ITS clusters");
    fTrackingEffHisto->Sumw2();
  }

  //
  //TPC constraining to vertex
  //
  if(GetAnalysisMode()==2){
    //initilization of fTPCConstrain
    Int_t  binsTPCConstrain[4] = {100, 90,            nPtBins, 30};
    Double_t minTPCConstrain[4] = {-5, 0,             ptMin,   -1.5};
    Double_t maxTPCConstrain[4] = {5,  2*TMath::Pi(), ptMax,  1.5};
    
    fTPCConstrain = new THnSparseF("fTPCConstrain","pull_phi:phi:pt:eta",4, binsTPCConstrain,minTPCConstrain,maxTPCConstrain);
    fTPCConstrain->SetBinEdges(2,binsPt);
    fTPCConstrain->GetAxis(0)->SetTitle("(#phi-#phi_{ref})/#sigma");
    fTPCConstrain->GetAxis(1)->SetTitle("phi (rad)");
    fTPCConstrain->GetAxis(2)->SetTitle("p_{T}");
    fTPCConstrain->GetAxis(3)->SetTitle("eta");
    fTPCConstrain->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("folderMatch","Analysis Matching Folder");
  
   // save merge status in object
  fMergeTHnSparseObj = fgMergeTHnSparse;

  delete [] binsPt;

}

//_____________________________________________________________________________
void AliPerformanceMatch::ProcessITSTPC(Int_t iTrack, AliESDEvent *const esdEvent, AliStack* /*const stack*/, AliESDtrack *const esdTrack)
{
  //
  // addition to standard analysis - check if ITS stand-alone tracks have a match in the TPC
  // Origin: A. Kalwait
  // Modified: J. Otwinowski
  if(!esdEvent) return;
  if(!esdTrack) return;
  //  if(!esdFriendTrack) return;

  // ITS stand alone tracks with SPD layers 
  if(!(esdTrack->GetStatus() & AliESDtrack::kITSpureSA)) return;
  if(!(esdTrack->GetStatus() & AliESDtrack::kITSrefit)) return;
  if(esdTrack->GetNcls(0)<4) return;
   if(!esdTrack->HasPointOnITSLayer(0) && !esdTrack->HasPointOnITSLayer(1)) return;
  
  const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
  //   const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
 AliESDtrack* tpcTrack2 = NULL;
  Bool_t hasMatch = kFALSE;
  for (Int_t jTrack = 0; jTrack < esdEvent->GetNumberOfTracks(); jTrack++) {
    // loop for all those tracks and check if the corresponding TPC track is found
    if (jTrack==iTrack) continue;
    AliESDtrack *trackTPC = esdEvent->GetTrack(jTrack);
    if (!trackTPC) continue;
    if (!trackTPC->GetTPCInnerParam()) continue;
    if(!(trackTPC->GetStatus() & AliESDtrack::kTPCrefit)) continue;
    
    // TPC nClust/track after first tracking pass
    // if(trackTPC->GetTPCNclsIter1()<fCutsRC->GetMinNClustersTPC()) continue;
    tpcTrack2 = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent, jTrack);
    if(!tpcTrack2) continue;
    if(!tpcTrack2->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack2; tpcTrack2=0; continue; } 
    
     if(!fCutsRC->AcceptTrack(tpcTrack2)) { delete tpcTrack2; tpcTrack2=0; continue; }
    // check matching
    if (TMath::Abs(esdTrack->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; tpcTrack2=0; continue; }
    if (TMath::Abs(esdTrack->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; tpcTrack2=0; continue; }
    if (TMath::Abs(esdTrack->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; tpcTrack2=0; continue; }
    
    hasMatch=kTRUE;
    break;
  }
  
  FillHistograms(tpcTrack2,esdTrack,hasMatch);     
  delete tpcTrack2;
}

//_____________________________________________________________________________
void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack)
{
  //
  // Match TPC and ITS min-bias tracks
  // at radius between detectors
  //
  if(!esdTrack) return;
  //   if(!esdFriendTrack) return;
   
  Bool_t isTPC = kFALSE;
  Bool_t isMatch = kFALSE;


  if(esdTrack->Charge()==0) return;
  if(!esdTrack->GetTPCInnerParam()) return;
  if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;

  if(!fCutsRC->AcceptTrack(esdTrack)) return;

  isTPC = kTRUE;
  
  if( (esdTrack->GetStatus()&AliESDtrack::kITSrefit))
    isMatch = kTRUE;
  
  if(isTPC){
    Double_t vecTrackingEff[5] = { static_cast<Double_t>(isMatch),esdTrack->Phi(), esdTrack->Pt(),esdTrack->Eta(),static_cast<Double_t>(esdTrack->GetITSclusters(0)) };
    fTrackingEffHisto->Fill(vecTrackingEff);
  }
}

//_____________________________________________________________________________
/*void AliPerformanceMatch::ProcessTPCTRD(AliStack* , AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
{
  return;
}*/

void AliPerformanceMatch::ProcessTPCConstrain(AliStack* /*const stack*/, AliESDEvent *const esdEvent, AliESDtrack *const esdTrack){
  //
  // Contrain TPC inner track to the vertex
  // then compare to the global tracks
  //
  if(!esdTrack) return;
  if(esdTrack->Charge()==0) return;

  if(!esdTrack->GetTPCInnerParam()) return;
  if(!(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) return;
  if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
  if(!fCutsRC->AcceptTrack(esdTrack)) return;

  Double_t x[3]; esdTrack->GetXYZ(x);
  Double_t b[3]; AliTracker::GetBxByBz(x,b);
  Bool_t isOK = kFALSE;
  const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
  
  AliExternalTrackParam * TPCinner = (AliExternalTrackParam *)esdTrack->GetTPCInnerParam();
  if(!TPCinner) return;
  
  //  isOK = TPCinner->Rotate(esdTrack->GetAlpha());
  //isOK = TPCinner->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());

  AliExternalTrackParam * TPCinnerC = new AliExternalTrackParam(*TPCinner);
  if (TPCinnerC) {
    isOK = TPCinnerC->ConstrainToVertex(vtxESD, b);

    // transform to the track reference frame 
    isOK = TPCinnerC->Rotate(esdTrack->GetAlpha());
    isOK = TPCinnerC->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
  }
  if(!isOK) return;

  Double_t sigmaPhi=0,deltaPhi=0,pullPhi=0;
  deltaPhi = TPCinnerC->GetSnp() - esdTrack->GetSnp();
  sigmaPhi = TMath::Sqrt(esdTrack->GetSigmaSnp2()+TPCinnerC->GetSigmaSnp2());
  if(sigmaPhi!=0)
    pullPhi = deltaPhi/sigmaPhi;

  Double_t vTPCConstrain[4] = {pullPhi,esdTrack->Phi(),esdTrack->Pt(),esdTrack->Eta()};
  fTPCConstrain->Fill(vTPCConstrain);  
  
  if(TPCinnerC)
    delete TPCinnerC;

  return;
}
//_____________________________________________________________________________
void AliPerformanceMatch::FillHistograms(AliESDtrack *const refParam, AliESDtrack *const param, Bool_t isRec) 
{
  //
  // fill performance histograms 
  // (TPC always as reference)
  //

  if(!refParam) return;
  if(!param) return;
  if(!isRec) return;
  
  //
  // Deltas (dy,dz,dphi,dtheta,dpt)
  //
  Float_t delta[5] = {0};
  if(param && isRec) {
    delta[0] = param->GetY()-refParam->GetY();
    delta[1] = param->GetZ()-refParam->GetZ();
    delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
    delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
    if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
  }
  // 
  // Pulls (y,z,snp,tanl,1/pt)
  //
  Float_t sigma[5] = {0};
  Float_t pull[5] = {0};
  if(param && isRec) {
    sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());   
    sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
    sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
    sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
    sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
    if(sigma[0]) pull[0] = delta[0] / sigma[0]; 
    if(sigma[1]) pull[1] = delta[1] / sigma[1]; 
    if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2]; 
    if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3]; 
    if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4]; 
  }

  // Fill histograms
  Double_t vResolHisto[9] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->Phi(),refParam->Eta(),refParam->Pt(),static_cast<Double_t>(isRec)};
  if(fabs(pull[4])<5)
    fResolHisto->Fill(vResolHisto);

  Double_t vPullHisto[9] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->Phi(),refParam->Eta(),refParam->OneOverPt(),static_cast<Double_t>(isRec)};
  if(fabs(pull[4])<5)
    fPullHisto->Fill(vPullHisto);
}

//_____________________________________________________________________________
void AliPerformanceMatch::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);
  } 
  
  // use ESD friends
  /*  if(bUseESDfriend) {
    if(!esdFriend) {
      Error("Exec","esdFriend not available");
      return;
    }
    }*/

  // trigger
  if(!bUseMC && GetTriggerClass()) {
    Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
    if(!isEventTriggered) return; 
  }

  // get TPC event vertex
  const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
  if(vtxESD && (vtxESD->GetStatus()<=0)) return;

  //  Process events
  for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
  { 
    AliESDtrack *track = esdEvent->GetTrack(iTrack);
    if(!track) continue;

    /*AliESDfriendTrack *friendTrack=0;
    if(bUseESDfriend) {
      friendTrack=esdFriend->GetTrack(iTrack);
      if(!friendTrack) continue;
      }*/

    if(GetAnalysisMode() == 0){
      if(!IsUseTOFBunchCrossing()){
	ProcessTPCITS(stack,track);
      }
      else
	if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
	  ProcessTPCITS(stack,track);
	}
    }
    /* else if(GetAnalysisMode() == 2) ProcessTPCTRD(stack,track,friendTrack);*/
    else if(GetAnalysisMode() == 1) {ProcessITSTPC(iTrack,esdEvent,stack,track);}
    else if(GetAnalysisMode() == 2){
      if(!IsUseTOFBunchCrossing()){
	ProcessTPCConstrain(stack,esdEvent,track);
      }
      else
	if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
	  ProcessTPCConstrain(stack,esdEvent,track);
	}
    }
    else {
      printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
      return;
    }
  }

}

//_____________________________________________________________________________
TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
  // Create resolution histograms
 
  TH1F *hisr, *hism;
  if (!gPad) new TCanvas;
  hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
  if (type) return hism;
  else 
    return hisr;
}

//_____________________________________________________________________________
void AliPerformanceMatch::Analyse() {
  // Analyse comparison information and store output histograms
  // in the folder "folderMatch"
  //
  TString selString;
  /*
  TH1::AddDirectory(kFALSE);
  TH1F *h=0;
  TH1F *h2=0;
  TH2F *h2D=0;
  */
  TObjArray *aFolderObj = new TObjArray;

  // write results in the folder 
  // TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
  // c->cd();

  // char name[256];
  // char title[256];

  if(GetAnalysisMode()==1) { 

  fResolHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
  fPullHisto->GetAxis(8)->SetRangeUser(1.0,2.0);  // only reconstructed
  for(Int_t i=0; i<5; i++) 
    {
      for(Int_t j=5; j<8; j++) 
	{
	  //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
	  if(j!=6) fResolHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
	  else fResolHisto->GetAxis(6)->SetRangeUser(-1.5,1.49);
	  fResolHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
	  
	  selString = "resol";
	  AddProjection(aFolderObj, "match", fResolHisto, i, j, &selString);
	  
	  
	  if(j!=6) fPullHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
	  else  fPullHisto->GetAxis(6)->SetRangeUser(-1.5,1.49); // eta window
	  fPullHisto->GetAxis(7)->SetRangeUser(0.01,10.);  // pt threshold
	  selString = "pull";
	  AddProjection(aFolderObj, "match", fPullHisto, i, j, &selString);
	  
	}
    }
  }
  // 
  // TPC efficiency wrt ITS
  //
  if(GetAnalysisMode()==0) { 
    selString = "trackingeff";
    AddProjection(aFolderObj, "match", fTrackingEffHisto, 0, &selString);
    
    for(Int_t i=1; i<5; i++) 
      {
	// all ITS standalone tracks
	fTrackingEffHisto->GetAxis(0)->SetRange(1,fTrackingEffHisto->GetAxis(0)->GetNbins());
	selString = "trackingeff_all";
	AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
	
      // TPC tracks which has matching with TPC
	fTrackingEffHisto->GetAxis(0)->SetRange(2,2);
	selString = "trackingeff_tpc";
	AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
      }
  }

  //
  //TPC constrained to vertex
  //
  if(GetAnalysisMode()==2) { 
    selString = "tpc";
    //    for(Int_t i=1; i<4; i++)
    AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 2, &selString);
    AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 3, &selString);
    AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 2, 3, &selString);
  }
  
  printf("exportToFolder\n");
  fAnalysisFolder = ExportToFolder(aFolderObj);
  
  // delete only TObjArray
  if(fFolderObj) delete fFolderObj;
  fFolderObj = aFolderObj;  
  aFolderObj=0;  
  
  }


//_____________________________________________________________________________
TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array) 
{
  // recreate folder avery time and export objects to new one
  //
  AliPerformanceMatch * 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;
}
 
//_____________________________________________________________________________
TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) { 
// create folder for analysed histograms
//
TFolder *folder = 0;
  folder = new TFolder(name.Data(),title.Data());

  return folder;
}

//_____________________________________________________________________________
Long64_t AliPerformanceMatch::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) 
  {
    AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
    if (entry == 0) continue; 
    if (merge) {
        if ((fResolHisto) && (entry->fResolHisto)) { fResolHisto->Add(entry->fResolHisto); }
        if ((fPullHisto) && (entry->fPullHisto)) { fPullHisto->Add(entry->fPullHisto); }
        if ((fTrackingEffHisto) && (entry->fTrackingEffHisto)) { fTrackingEffHisto->Add(entry->fTrackingEffHisto); }

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