ROOT logo
/****************************************************************************
 *           Very important, delicate and rather obscure macro.             *
 *                    ("a la" AliTPCComparison.C)                           *
 *                                                                          *
 *               Creates list of "trackable" tracks,                        *
 *             calculates efficiency, resolutions etc.                      *
 *         (To get the list of the "trackable" tracks one should            *
 *             first run the AliTPCComparison.C macro)                      *
 *                                                                          *
 *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
 ****************************************************************************/

#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 "AliHeader.h"
  #include "AliTrackReference.h"
  #include "AliRunLoader.h"
  #include "AliRun.h"
  #include "AliESDEvent.h"
  #include "AliESDtrack.h"
#endif

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

extern AliRun *gAlice;
extern TBenchmark *gBenchmark;
extern TROOT *gROOT;

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

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

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

   TH1F *hp=(TH1F*)gROOT->FindObject("hp");
   if (!hp) {
      hp=new TH1F("hp","PHI resolution",50,-70.,70.); 
      hp->SetFillColor(4);
      hp->SetXTitle("(mrad)"); 
   }
   TH1F *hl=(TH1F*)gROOT->FindObject("hl");
   if (!hl) {
      hl=new TH1F("hl","LAMBDA resolution",50,-70,70);
      hl->SetFillColor(4);
      hl->SetXTitle("(mrad)");
   }
   TH1F *hc=(TH1F*)gROOT->FindObject("hc");
   if (!hc) {
      hc=new TH1F("hc","Number of the assigned clusters",25,110,135);
      hc->SetLineColor(2);
   }
   TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
   if (!hpt) {
      hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
      hpt->SetFillColor(2);
      hpt->SetXTitle("(%)");
   }
   TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
   if (!hmpt) {
      hmpt=new TH1F("hmpt","Y and Z resolution",30,-30,30); 
      hmpt->SetFillColor(6);
      hmpt->SetXTitle("(mm)");
   }
   TH1F *hz=(TH1F*)gROOT->FindObject("hz");
   if (!hz) hz=new TH1F("hz","Z resolution",30,-30,30); 



   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.,1000.);

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


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

   TFile *refFile=TFile::Open(fname,"old");
   if (!refFile || !refFile->IsOpen()) {
   ::Info("AliTRDComparisonV2.C","Marking good tracks (will take a while)...");
     if (GoodTracksTRD(dir)) {
        ::Error("AliTRDComparisonV2.C","Can't generate the reference file !");
        return 1;
     }
   }
   refFile=TFile::Open(fname,"old");
   if (!refFile || !refFile->IsOpen()) {
     ::Error("AliTRDComparisonV2.C","Can't open the reference file !");
     return 1;
   }   
  
   TTree *trdTree=(TTree*)refFile->Get("trdTree");
   if (!trdTree) {
     ::Error("AliTRDComparisonV2.C","Can't get the reference tree !");
     return 2;
   }
   TBranch *branch=trdTree->GetBranch("TRD");
   if (!branch) {
     ::Error("AliTRDComparisonV2.C","Can't get the TRD branch !");
     return 3;
   }
   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/AliESDtrd.root",dir);
      ef=TFile::Open(fname);
      if ((!ef)||(!ef->IsOpen())) {
         ::Error("AliTRDComparisonV2.C","Can't open AliESDtrd.root !");
         return 4;
      }
   }
   AliESDEvent* event = new AliESDEvent();
   TTree* esdTree = (TTree*) ef->Get("esdTree");
   if (!esdTree) {
      ::Error("AliTRDComparisonV2.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 (trdTree->GetEvent(e++)==0) {
        cerr<<"No reconstructable tracks !\n";
        continue;
     }

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

     const Int_t MAX=15000;
     Int_t notf[MAX], nnotf=0;
     Int_t fake[MAX], nfake=0;
     Int_t mult[MAX], numb[MAX], nmult=0;
     Int_t k;
     for (k=0; k<ngood; k++) {
	AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k); 
        Int_t lab=ref->Label(), tlab=-1;
        Float_t ptg=ref->Pt();

        if (ptg<ptcutl) continue;
        if (ptg>ptcuth) continue;

        allselected++;

        hgood->Fill(ptg);

        AliESDtrack *esd=0;
        Int_t cnt=0;
        for (Int_t i=0; i<nentr; i++) {
           AliESDtrack *t=event->GetTrack(i);
	   UInt_t status=t->GetStatus();

           if ((status&AliESDtrack::kTRDout)==0) continue;

           Int_t lbl=t->GetTRDLabel();
           if (lab==TMath::Abs(lbl)) {
	      if (cnt==0) {esd=t; tlab=lbl;}
              cnt++;
           }
        }
        if (cnt==0) {
           notf[nnotf++]=lab;
           continue;
        } else if (cnt>1){
           mult[nmult]=lab;
           numb[nmult]=cnt; nmult++;        
        }

        if (lab==tlab) hfound->Fill(ptg);
        else {
          fake[nfake++]=lab;
          hfake->Fill(ptg); 
        }

        AliExternalTrackParam out(*esd->GetOuterParam());

        Double_t bz=event->GetMagneticField();
        Double_t xg = ref->LocalX();
        out.PropagateTo(xg,bz);

        Float_t phi=TMath::ASin(out.GetSnp()) + out.GetAlpha();
        if (phi < 0) phi+=2*TMath::Pi();
        if (phi >=2*TMath::Pi()) phi-=2*TMath::Pi();
        Float_t phig=ref->Phi();
        hp->Fill((phi - phig)*1000.);

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

        hc->Fill(esd->GetTRDclusters(0));

        Float_t pt_1=out.OneOverPt();
        hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);

        Float_t y=out.GetY();
        Float_t yg=ref->LocalY();
        hmpt->Fill(10*(y - yg));

        Float_t z=out.GetZ();
        Float_t zg=ref->Z();
        hz->Fill(10*(z - zg));

        Float_t mom=1./(pt_1*TMath::Cos(lam));
        Float_t dedx=esd->GetTRDsignal();
        hep->Fill(mom,dedx,1.);

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

        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 (k=0; k<nnotf; k++){
       cout<<notf[k]<<"\t";
       if ((k%9)==8) cout<<"\n";
     }
     cout<<"\n\nList of fake  tracks :\n";
     for (k=0; k<nfake; k++){
       cout<<fake[k]<<"\t";
       if ((k%9)==8) cout<<"\n";
     }
     cout<<"\n\nList of multiple found tracks :\n";
     for (k=0; k<nmult; k++) {
         cout<<"id.   "<<mult[k]
             <<"     found - "<<numb[k]<<"times\n";
     }
     cout<<endl;

     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 trdTree;
   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); 
   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);
   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); 
   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);
   if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); 
   hz->Draw("same"); 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.2,1.0,6.1,1.0); line1->SetLineStyle(4);
   line1->Draw("same");
   TLine *line2 = new TLine(0.2,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("AliTRDComparisonV2.root","RECREATE");
   c1->Write();
   c2->Write();
   fc.Close();

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

   return 0;
}



Int_t GoodTracksTRD(const Char_t *dir) {
   if (gAlice) { 
       delete gAlice->GetRunLoader();
       delete gAlice;//if everything was OK here it is already NULL
       gAlice = 0x0;
   }

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

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

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

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

   sprintf(fname,"%s/GoodTracksTPC.root",dir);
   TFile *tpcFile=TFile::Open(fname);
   if ((!tpcFile)||(!tpcFile->IsOpen())) {
       ::Error("GoodTracksTRD","Can't open the GoodTracksTPC.root !");
       delete rl;
       return 5; 
   }
   TClonesArray dum("AliTrackReference",1000), *tpcRefs=&dum;
   TTree *tpcTree=(TTree*)tpcFile->Get("tpcTree");
   if (!tpcTree) {
       ::Error("GoodTracksTRD","Can't get the TPC reference tree !");
       delete rl;
       return 6;
   }
   TBranch *tpcBranch=tpcTree->GetBranch("TPC");
   if (!tpcBranch) {
      ::Error("GoodTracksTRD","Can't get the TPC reference branch !");
      delete rl;
      return 7;
   }
   tpcBranch->SetAddress(&tpcRefs);

   sprintf(fname,"%s/GoodTracksTRD.root",dir);
   TFile *trdFile=TFile::Open(fname,"recreate");
   TClonesArray dummy2("AliTrackReference",1000), *trdRefs=&dummy2;
   TTree trdTree("trdTree","Info about the reconstructable TRD tracks");
   trdTree.Branch("TRD",&trdRefs);

   //********  Loop over generated events 
   for (Int_t e=0; e<nev; e++) {
     Int_t k;

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

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


     TTree *TR=rl->TreeTR();
     TBranch *branch=TR->GetBranch("TrackReferences");
     if (branch==0) {
        ::Error("GoodTracksTRD","No TRD track references !");
        delete rl;
        return 5;
     }
     TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
     branch->SetAddress(&refs);
     Int_t nr=TR->GetEntries();

     //Preselect the "good" track candidates (TRD info only)
     Int_t nt=0;
     for (Int_t r=0; r<nr; r++) {
         refs->Clear();
         TR->GetEntry(r);
         Int_t n=refs->GetEntriesFast();
 
         if (n<=1) continue;

	 AliTrackReference *ref0=(AliTrackReference *)refs->UncheckedAt(0);
         if (ref0->LocalX() > 300.) continue;

	 AliTrackReference *refn=0x0;
	 for (Int_t iref=n-1; iref>=0; --iref) {
	   refn = (AliTrackReference *)refs->UncheckedAt(iref);
	   if (refn->LocalX() > 363. &&
	       refn->DetectorId() == AliTrackReference::kTRD) break;
	   refn = 0x0;
	 }
	 if (!refn) continue;
         if (TMath::Abs(ref0->Alpha() - refn->Alpha()) > 1e-5) continue;   

         new((*trdRefs)[nt++]) AliTrackReference(*refn);
     }

     //Check if the candidates are "good" from the TPC point of view
     tpcTree->GetEvent(e);
     Int_t ntpc=tpcRefs->GetEntriesFast();
     for (Int_t t=0; t<nt; t++) {
	 AliTrackReference *ref=(AliTrackReference *)trdRefs->UncheckedAt(t);
         Int_t lab=ref->Label();
         for (k=0; k<ntpc; k++) {
             AliTrackReference 
                *tpcRef=(AliTrackReference *)tpcRefs->UncheckedAt(k);
             if (tpcRef->Label()==lab) {
	        ref->SetLength(tpcRef->GetLength());
                break;
             } 
        }
	if (k==ntpc) delete trdRefs->RemoveAt(t); 
     }
     trdRefs->Compress();
     trdTree.Fill();

     trdRefs->Clear();
     tpcRefs->Clear();

   } //*** end of the loop over generated events

   trdTree.Write();
   trdFile->Close();

   delete tpcTree;
   tpcFile->Close();

   delete rl;
   return 0;
}


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