ROOT logo
/****************************************************************************
 *           Very important, delicate and rather obscure macro.             *
 *                                                                          *
 *               Creates list of "trackable" tracks,                        *
 *             calculates efficiency, resolutions etc.                      *
 *  The ESD tracks must be in an appropriate state (kITSin or kITSrefit)    *
 *                                                                          *
 *           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 "AliStack.h"
  #include "AliHeader.h"
  #include "AliTrackReference.h"
  #include "AliRunLoader.h"
  #include "AliRun.h"
  #include "AliESDEvent.h"
  #include "AliESDtrack.h"

  #include "AliITSRecPoint.h"
  #include "AliITSLoader.h"
#endif

Int_t GoodTracksITS(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 AliITSComparisonV2
(Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
   gBenchmark->Start("AliITSComparisonV2");

   ::Info("AliITSComparisonV2.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","Transverse impact parameter",30,-777,777); 
   hmpt->SetFillColor(6);

   TH1F *hz=(TH1F*)gROOT->FindObject("hz");
   if (!hz) hz=new TH1F("hz","Longitudinal impact parameter",30,-777,777); 



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

   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/GoodTracksITS.root",dir);

   TFile *refFile=TFile::Open(fname,"old");
   if (!refFile || !refFile->IsOpen()) {
   ::Info("AliITSComparisonV2.C","Marking good tracks (will take a while)...");
     if (GoodTracksITS(dir)) {
        ::Error("AliITSComparisonV2.C","Can't generate the reference file !");
        return 1;
     }
   }
   refFile=TFile::Open(fname,"old");
   if (!refFile || !refFile->IsOpen()) {
     ::Error("AliITSComparisonV2.C","Can't open the reference file !");
     return 1;
   }   
  
   TTree *itsTree=(TTree*)refFile->Get("itsTree");
   if (!itsTree) {
     ::Error("AliITSComparisonV2.C","Can't get the reference tree !");
     return 2;
   }
   TBranch *branch=itsTree->GetBranch("ITS");
   if (!branch) {
     ::Error("AliITSComparisonV2.C","Can't get the ITS 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/AliESDits.root",dir);
      ef=TFile::Open(fname);
      if ((!ef)||(!ef->IsOpen())) {
         ::Error("AliITSComparisonV2.C","Can't open AliESDits.root !");
         return 4;
      }
   }
   AliESDEvent* event = new AliESDEvent();
   TTree* esdTree = (TTree*) ef->Get("esdTree");
   if (!esdTree) {
      ::Error("AliITSComparison.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 (itsTree->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=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());

        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::kITSrefit)==0) continue;

           Int_t lbl=t->GetLabel();
           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); 
        }

        Double_t alpha=esd->GetAlpha(),xv,par[5]; 
        esd->GetExternalParameters(xv,par);
        Float_t phi=TMath::ASin(par[2]) + alpha;
        if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
        if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
        Float_t lam=TMath::ATan(par[3]); 
        Float_t pt_1=TMath::Abs(par[4]);

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

        Float_t d,z; esd->GetImpactParameters(d,z);
        hmpt->Fill(10000*d);
        hz->Fill(10000*z);

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

        Float_t mom=1./(pt_1*TMath::Cos(lam));
        Float_t dedx=esd->GetITSsignal();
        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 itsTree;
   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("(micron)");
   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("AliITSComparisonV2.root","RECREATE");
   c1->Write();
   c2->Write();
   fc.Close();

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

   return 0;
}



Int_t GoodTracksITS(const Char_t *dir) {
   if (gAlice) { 
       delete AliRunLoader::Instance();
       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("GoodTracksITS","Can't start session !");
      return 1;
   }

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

   AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
   if (itsl == 0x0) {
       ::Error("GoodTracksITS","Can not find the ITSLoader");
       delete rl;
       return 4;
   }
   itsl->LoadRecPoints();
  

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

   sprintf(fname,"%s/GoodTracksTPC.root",dir);
   TFile *tpcFile=TFile::Open(fname);
   if ((!tpcFile)||(!tpcFile->IsOpen())) {
       ::Error("GoodTracksITS","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("GoodTracksITS","Can't get the TPC reference tree !");
       delete rl;
       return 6;
   }
   TBranch *tpcBranch=tpcTree->GetBranch("TPC");
   if (!tpcBranch) {
      ::Error("GoodTracksITS","Can't get the TPC reference branch !");
      delete rl;
      return 7;
   }
   tpcBranch->SetAddress(&tpcRefs);

   sprintf(fname,"%s/GoodTracksITS.root",dir);
   TFile *itsFile=TFile::Open(fname,"recreate");
   TClonesArray dummy2("AliTrackReference",1000), *itsRefs=&dummy2;
   TTree itsTree("itsTree","Tree with info about the reconstructable ITS tracks");
   itsTree.Branch("ITS",&itsRefs);

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

     rl->GetEvent(e);  itsFile->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 (k=0; k<np; k++) good[k]=0;

     TTree *cTree=itsl->TreeR();
     if (!cTree) {
        ::Error("GoodTracksITS","Can't get the cluster tree !"); 
        delete rl;
        return 8;
     }
     TBranch *branch=cTree->GetBranch("ITSRecPoints");
     if (!branch) {
        ::Error("GoodTracksITS","Can't get the clusters branch !"); 
        delete rl;
        return 9;
     }
     TClonesArray dummy("AliITSRecPoint",10000), *clusters=&dummy;
     branch->SetAddress(&clusters);

     Int_t entr=(Int_t)cTree->GetEntries();
     for (k=0; k<entr; k++) {
         cTree->GetEvent(k);
         Int_t ncl=clusters->GetEntriesFast(); if (ncl==0) continue;
         while (ncl--) {
            AliITSRecPoint *pnt=(AliITSRecPoint*)clusters->UncheckedAt(ncl);

            Int_t lay=pnt->GetLayer();
            if (lay<0 || lay>5) {
               ::Error("GoodTracksITS","Wrong layer !");
               delete rl;
               return 10;
            }

            Int_t l0=pnt->GetLabel(0);
	       if (l0>=np) {
// 		 cerr<<"Wrong label: "<<l0<<endl;
		 continue;
	       }
            Int_t l1=pnt->GetLabel(1);
	       if (l1>=np) {
// 		 cerr<<"Wrong label: "<<l1<<endl;
		 continue;
	       }
            Int_t l2=pnt->GetLabel(2);
	       if (l2>=np) {
// 		 cerr<<"Wrong label: "<<l2<<endl;
		 continue;
	       }
            Int_t mask=1<<lay;
            if (l0>=0) good[l0]|=mask; 
            if (l1>=0) good[l1]|=mask; 
            if (l2>=0) good[l2]|=mask;
         }
         clusters->Clear();
     }
   

     //****** select tracks which are "good" enough
     AliStack* stack = rl->Stack();

     tpcTree->GetEvent(e);
     Int_t nk=tpcRefs->GetEntriesFast();
     Int_t nt=0;
     for (k=0; k<nk; k++) {
        AliTrackReference *tpcRef=(AliTrackReference *)tpcRefs->UncheckedAt(k);
        Int_t lab=tpcRef->Label();
        if (good[lab] != 0x3F) continue;
        TParticle *p = (TParticle*)stack->Particle(lab);
        if (p == 0x0) {
           cerr<<"Can not get particle "<<lab<<endl;
           continue;
        }

	AliTrackReference *ref=new((*itsRefs)[nt]) AliTrackReference(*tpcRef);
	ref->SetMomentum(p->Px(),p->Py(),p->Pz());
	ref->SetPosition(p->Vx(),p->Vy(),p->Vz());
        nt++;
     }
     tpcRefs->Clear();

     itsTree.Fill();
     itsRefs->Clear();

     delete[] good;

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

   itsTree.Write();
   itsFile->Close();

   delete tpcTree;
   tpcFile->Close();

   delete rl;
   return 0;
}


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