ROOT logo
/****************************************************************************
 *           Very important, delicate and rather obscure macro.             *
 *                                                                          *
 *               Creates list of "trackable" tracks,                        *
 *             calculates efficiency, resolutions etc.                      *
 *     There is a possibility to run this macro over several events.         *
 *                                                                          *
 *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
 * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de        *
 ****************************************************************************/

#if !defined(__CINT__) || defined(__MAKECINT__)
  #include <TMath.h>
  #include <TError.h>
  #include <Riostream.h>
  #include <TH1F.h>
  #include <TH2F.h>
  #include <TTree.h>
  #include <TParticle.h>
  #include <TCanvas.h>
  #include <TLine.h>
  #include <TText.h>
  #include <TBenchmark.h>
  #include <TStyle.h>
  #include <TFile.h>
  #include <TROOT.h>

  #include "AliStack.h"
  #include "AliHeader.h"
  #include "AliTrackReference.h"
  #include "AliRunLoader.h"
  #include "AliRun.h"
  #include "AliESDEvent.h"
  #include "AliESDtrack.h"

  #include "AliSimDigits.h"
  #include "AliTPC.h"
  #include "AliTPCParamSR.h"
  #include "AliTPCClustersArray.h"
  #include "AliTPCClustersRow.h"
  #include "AliTPCcluster.h"
  #include "AliTPCLoader.h"

  #include "AliCDBManager.h"
  #include "AliTPCcalibDB.h"
#endif

Int_t GoodTracksTPC(const Char_t *dir=".");

extern TBenchmark *gBenchmark;
extern TROOT *gROOT;

static Int_t allgood=0;
static Int_t allselected=0;
static Int_t allfound=0;

Int_t AliTPCComparison
(Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
   gBenchmark->Start("AliTPCComparison");

   ::Info("AliTPCComparison.C","Doing comparison...");



   TH1F *hp=(TH1F*)gROOT->FindObject("hp");
   if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.); 
   hp->SetFillColor(4);

   TH1F *hl=(TH1F*)gROOT->FindObject("hl");
   if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
   hl->SetFillColor(4);

   TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
   if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
   hpt->SetFillColor(2);
 
   TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
   if (!hmpt) 
      hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60); 
   hmpt->SetFillColor(6);



   TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
   if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
    
   TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
   if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);

   TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
   if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);

   TH1F *hg=(TH1F*)gROOT->FindObject("hg");
   if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
   hg->SetLineColor(4); hg->SetLineWidth(2);

   TH1F *hf=(TH1F*)gROOT->FindObject("hf");
   if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
   hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);

   TH1F *he=(TH1F*)gROOT->FindObject("he");
   if (!he) 
      he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);

   TH2F *hep=(TH2F*)gROOT->FindObject("hep");
   if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
   hep->SetMarkerStyle(8);
   hep->SetMarkerSize(0.4);


   Char_t fname[100];
   sprintf(fname,"%s/GoodTracksTPC.root",dir);

   TFile *refFile=TFile::Open(fname,"old");
   if (!refFile || !refFile->IsOpen()) {
     ::Info("AliTPCComparison.C","Marking good tracks (will take a while)...");
     if (GoodTracksTPC(dir)) {
        ::Error("AliTPCComparison.C","Can't generate the reference file !");
        return 1;
     }
   }
   refFile=TFile::Open(fname,"old");
   if (!refFile || !refFile->IsOpen()) {
     ::Error("AliTPCComparison.C","Can't open the reference file !");
     return 2;
   }   
  
   TTree *tpcTree=(TTree*)refFile->Get("tpcTree");
   if (!tpcTree) {
     ::Error("AliTPCComparison.C","Can't get the reference tree !");
     return 3;
   }
   TBranch *branch=tpcTree->GetBranch("TPC");
   if (!branch) {
     ::Error("AliTPCComparison.C","Can't get the TPC branch !");
     return 4;
   }
   TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
   branch->SetAddress(&refs);


   sprintf(fname,"%s/AliESDs.root",dir);
   TFile *ef=TFile::Open(fname);
   if ((!ef)||(!ef->IsOpen())) {
      sprintf(fname,"%s/AliESDtpc.root",dir);
      ef=TFile::Open(fname);
      if ((!ef)||(!ef->IsOpen())) {
         ::Error("AliTPCComparison.C","Can't open AliESDtpc.root !");
         return 5;
      }
   }
   AliESDEvent* event = new AliESDEvent();
   TTree* esdTree = (TTree*) ef->Get("esdTree");
   if (!esdTree) {
      ::Error("AliTPCComparison.C", "no ESD tree found");
      return 6;
   }
   event->ReadFromTree(esdTree);


   //******* Loop over events *********
   Int_t e=0;
   while (esdTree->GetEvent(e)) {
      cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";

      Int_t nentr=event->GetNumberOfTracks();
      allfound+=nentr;

      if (tpcTree->GetEvent(e++)==0) {
	 cerr<<"No reconstructable tracks !\n";
         continue;
      }

      Int_t ngood=refs->GetEntriesFast(); 
      allgood+=ngood;

      const Int_t MAX=15000;
    //MI change
      Int_t track_notfound[MAX], itrack_notfound=0;
      Int_t track_fake[MAX], itrack_fake=0;
      Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0;

      Int_t i;
      for (Int_t k=0; k<ngood; k++) {
	AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k); 
        Int_t lab=ref->Label(), tlab=-1;
        Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());

        if (ptg<1e-33) continue; // for those not crossing 0 pad row
  
        if (ptg<ptcutl) continue;
        if (ptg>ptcuth) continue;

        allselected++;

        hgood->Fill(ptg);

        AliESDtrack *track=0;
        for (i=0; i<nentr; i++) {
          track=event->GetTrack(i);
          tlab=track->GetTPCLabel();
          if (lab==TMath::Abs(tlab)) break;
        }
        if (i==nentr) {
	   track_notfound[itrack_notfound++]=lab;
           continue;
        }
      
        //MI change  - addition
        Int_t micount=0;
        Int_t mi;
        AliESDtrack * mitrack;
        for (mi=0; mi<nentr; mi++) {
	  mitrack=event->GetTrack(mi);          
	  if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
        }
        if (micount>1) {
	  track_multifound[itrack_multifound]=lab;
	  track_multifound_n[itrack_multifound]=micount;
	  itrack_multifound++;
        }
        if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
        if (lab==tlab) hfound->Fill(ptg);
        else { 
	  track_fake[itrack_fake++]=lab;
	  hfake->Fill(ptg); 
        }
      
        Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
        Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
        if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
        if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
        Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
        Float_t lam=TMath::ATan2(pxpypz[2],pt); 
        Float_t pt_1=1/pt;

        Int_t pdg=(Int_t)ref->GetLength();  //this is particle's PDG !

        if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons
           hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
        } else {
           Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
           hp->Fill((phi - phig)*1000.);

           Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
           hl->Fill((lam - lamg)*1000.);

           hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
        }

        Float_t mom=pt/TMath::Cos(lam);
        Float_t dedx=track->GetTPCsignal();
        hep->Fill(mom,dedx,1.);
        if (TMath::Abs(pdg)==211) //pions
	   if (mom>0.4 && mom<0.5) {
              he->Fill(dedx,1.);
           }
      }

      cout<<"\nList of Not found tracks :\n";
      for ( i = 0; i< itrack_notfound; i++){
        cout<<track_notfound[i]<<"\t";
        if ((i%5)==4) cout<<"\n";
      }
      cout<<"\nList of fake  tracks :\n";
      for ( i = 0; i< itrack_fake; i++){
        cout<<track_fake[i]<<"\t";
        if ((i%5)==4) cout<<"\n";
      }
      cout<<"\nList of multiple found tracks :\n";
      for ( i=0; i<itrack_multifound; i++) {
          cout<<"id.   "<<track_multifound[i]
              <<"     found - "<<track_multifound_n[i]<<"times\n";
      }

      cout<<"Number of found tracks ="<<nentr<<endl;
      cout<<"Number of \"good\" tracks ="<<ngood<<endl;

      refs->Clear();
  }// ***** End of the loop over events

   delete event;
   delete esdTree;
   ef->Close();

   delete tpcTree;
   refFile->Close();

   Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
   if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
   cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
   cout<<"Total number of found tracks ="<<allfound<<endl;
   cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
   cout<<endl;

   gStyle->SetOptStat(111110);
   gStyle->SetOptFit(1);

   TCanvas *c1=new TCanvas("c1","",0,0,700,850);

   Int_t minc=33; 

   TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
   p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
   hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); 
   if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();

   TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
   p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
   hl->SetXTitle("(mrad)");
   if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();

   TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
   p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
   hpt->SetXTitle("(%)");
   if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();

   TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
   p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
   hmpt->SetXTitle("(%)");
   if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();

   TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); 
   p5->SetFillColor(41); p5->SetFrameFillColor(10);
   hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
   hg->Divide(hfound,hgood,1,1.,"b");
   hf->Divide(hfake,hgood,1,1.,"b");
   hg->SetMaximum(1.4);
   hg->SetYTitle("Tracking efficiency");
   hg->SetXTitle("Pt (GeV/c)");
   hg->Draw();

   TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
   line1->Draw("same");
   TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
   line2->Draw("same");

   hf->SetFillColor(1);
   hf->SetFillStyle(3013);
   hf->SetLineColor(2);
   hf->SetLineWidth(2);
   hf->Draw("histsame");
   TText *text = new TText(0.461176,0.248448,"Fake tracks");
   text->SetTextSize(0.05);
   text->Draw();
   text = new TText(0.453919,1.11408,"Good tracks");
   text->SetTextSize(0.05);
   text->Draw();



   TCanvas *c2=new TCanvas("c2","",320,32,530,590);

   TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
   p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); 
   he->SetFillColor(2); he->SetFillStyle(3005);  
   he->SetXTitle("Arbitrary Units"); 
   if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();

   TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); 
   p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
   hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); 
   hep->Draw(); c1->cd();

   TFile fc("AliTPCComparison.root","RECREATE");
   c1->Write();
   c2->Write();
   fc.Close();

   gBenchmark->Stop("AliTPCComparison");
   gBenchmark->Show("AliTPCComparison");

   return 0;
}


Int_t GoodTracksTPC(const Char_t *dir) {
   Char_t fname[100];
   sprintf(fname,"%s/galice.root",dir);

   AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
   if (!rl) {
      ::Error("GoodTracksTPC","Can't start session !");
      return 1;
   }

   rl->LoadgAlice();
   rl->LoadHeader();
   rl->LoadKinematics();
   rl->LoadTrackRefs();

   AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
   if (tpcl == 0x0) {
      ::Error("GoodTracksTPC","Can not find TPCLoader !");
      delete rl;
      return 2;
   }
      
   AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
   Int_t ver = TPC->IsVersion(); 
   cout<<"TPC version "<<ver<<" has been found !\n";
   if (ver==1) tpcl->LoadRecPoints();
   else if (ver==2) tpcl->LoadDigits();
   else {
      ::Error("GoodTracksTPC","Wrong TPC version !");
      delete rl;
      return 3;
   } 


   AliCDBManager *man=AliCDBManager::Instance();
   man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
   man->SetRun(0);
   AliTPCParamSR *digp=
   (AliTPCParamSR*)(AliTPCcalibDB::Instance()->GetParameters());
   if (!digp) { 
     ::Error("AliTPCComparison.C","TPC parameters have not been found !");
     delete rl;
     return 4; 
   }

   Int_t nrow_up=digp->GetNRowUp();
   Int_t nrows=digp->GetNRowLow()+nrow_up;
   Int_t zero=digp->GetZeroSup();
   Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
   Int_t good_number=Int_t(0.4*nrows);

   Int_t nev=rl->GetNumberOfEvents();
   ::Info("GoodTracksTPC","Number of events : %d\n",nev);  


   sprintf(fname,"%s/GoodTracksTPC.root",dir);
   TFile *file=TFile::Open(fname,"recreate");

   TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
   TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks");
   tpcTree.Branch("TPC",&refs);

   //********  Loop over generated events 
   for (Int_t e=0; e<nev; e++) {
     Int_t nt=0;
     refs->Clear();

     Int_t i;

     rl->GetEvent(e);  file->cd();

     Int_t np = rl->GetHeader()->GetNtrack();
     cout<<"Event "<<e<<" Number of particles: "<<np<<endl;

     //******** Fill the "good" masks
     Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0;
     switch (ver) {
     case 1:
        {
        AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
        ca.Setup(digp);
        ca.SetClusterType("AliTPCcluster");
        ca.ConnectTree(tpcl->TreeR());
        Int_t nrows=Int_t(ca.GetTree()->GetEntries());
        for (Int_t n=0; n<nrows; n++) {
          AliSegmentID *s=ca.LoadEntry(n);
          Int_t sec,row;
          digp->AdjustSectorRow(s->GetID(),sec,row);
          AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
          Int_t ncl=clrow.GetArray()->GetEntriesFast();
          while (ncl--) {
              AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
              Int_t lab=c->GetLabel(0);
              if (lab<0) continue; //noise cluster
              lab=TMath::Abs(lab);

              if (sec>=digp->GetNInnerSector())
                 if (row==nrow_up-1) good[lab]|=0x4000;
              if (sec>=digp->GetNInnerSector())
                 if (row==nrow_up-1-gap) good[lab]|=0x1000;

              if (sec>=digp->GetNInnerSector())
                 if (row==nrow_up-1-shift) good[lab]|=0x2000;
              if (sec>=digp->GetNInnerSector())
                 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;

              good[lab]++;
          }
          ca.ClearRow(sec,row);
        }
        delete pca;
        }
        break;
     case 2:
        {
        TTree *TD=tpcl->TreeD();
      
        AliSimDigits da, *digits=&da;
        TD->GetBranch("Segment")->SetAddress(&digits);

        Int_t *count = new Int_t[np];
        Int_t i;
        for (i=0; i<np; i++) count[i]=0;

        Int_t sectors_by_rows=(Int_t)TD->GetEntries();
        for (i=0; i<sectors_by_rows; i++) {
          if (!TD->GetEvent(i)) continue;
          Int_t sec,row;
          digp->AdjustSectorRow(digits->GetID(),sec,row);
	  //cerr<<sec<<' '<<row<<'\r';
          digits->First();
          do { //Many thanks to J.Chudoba who noticed this
              Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
              Short_t dig = digits->GetDigit(it,ip);
              Int_t idx0=digits->GetTrackID(it,ip,0); 
              Int_t idx1=digits->GetTrackID(it,ip,1);
              Int_t idx2=digits->GetTrackID(it,ip,2);
              if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
              if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
              if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
          } while (digits->Next());
          for (Int_t j=0; j<np; j++) {
              if (count[j]>1) {
                 if (sec>=digp->GetNInnerSector())
		   if (row==nrow_up-1    ) good[j]|=0x4000;
                 if (sec>=digp->GetNInnerSector())
		   if (row==nrow_up-1-gap) good[j]|=0x1000;

                 if (sec>=digp->GetNInnerSector())
		   if (row==nrow_up-1-shift) good[j]|=0x2000;
                 if (sec>=digp->GetNInnerSector())
		   if (row==nrow_up-1-gap-shift) good[j]|=0x800;
                 good[j]++;
              }
              count[j]=0;
          }
        }
        delete[] count;
        }
        break;
     }


     //****** select tracks which are "good" enough
     AliStack* stack = rl->Stack();
     for (i=0; i<np; i++) {
        if ((good[i]&0x5000) != 0x5000)
        if ((good[i]&0x2800) != 0x2800) continue;
        if ((good[i]&0x7FF ) < good_number) continue;

        TParticle *p = (TParticle*)stack->Particle(i);
        if (p == 0x0) {
         cerr<<"Can not get particle "<<i<<endl;
         continue;
        }
        if (p->Pt()<0.100) continue;
        if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;

        Double_t vx=p->Vx(),vy=p->Vy(),vz=p->Vz();
        if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
        if (TMath::Abs(vz) > 50.) continue;

        AliTrackReference *ref=new((*refs)[nt]) AliTrackReference();

        ref->SetLabel(i);
        Int_t pdg=p->GetPdgCode();
        ref->SetLength(pdg);  //This will the particle's PDG !
        ref->SetMomentum(0.,0.,0.);  
        ref->SetPosition(0.,0.,0.);

        nt++;
     }

     //**** check if there is also information at the entrance of the TPC
     TTree *TR=rl->TreeTR();
     TBranch *branch=TR->GetBranch("TrackReferences");
     if (branch==0) {
        ::Error("GoodTracksTPC","No track references !");
        delete rl;
        return 5;
     }
     TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy;
     branch->SetAddress(&tpcRefs);
     Int_t nr=(Int_t)TR->GetEntries();
     for (Int_t r=0; r<nr; r++) {
         //cerr<<r<<' '<<nr<<'\r';
         TR->GetEvent(r);

	 Int_t nref = tpcRefs->GetEntriesFast();
         if (!nref) continue;
         AliTrackReference *tpcRef= 0x0;	 
	 for (Int_t iref=0; iref<nref; ++iref) {
	   tpcRef = (AliTrackReference*)tpcRefs->UncheckedAt(iref);
	   if (tpcRef->DetectorId() == AliTrackReference::kTPC) break;
	   tpcRef = 0x0;
	 }

	 if (!tpcRef) continue;

         Int_t j;
	 AliTrackReference *ref=0;
         for (j=0; j<nt; j++) {
	   ref=(AliTrackReference *)refs->UncheckedAt(j);
           if (ref->Label()==tpcRef->Label()) break;
           ref=0;
	 }  
         if (ref==0) continue;

         ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz());
         ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z());

	 tpcRefs->Clear();
     }

     tpcTree.Fill();

     delete[] good;
   } //****** end of the loop over generated events

   
   tpcTree.Write();
   file->Close();

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