ROOT logo
//*************************************************************************
//                Example of using the PID classes
//  which calculates the efficiency and contamination of the combined PID.
//
//      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
//*************************************************************************

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

  #include <AliStack.h>
  #include <AliRunLoader.h>
  #include <AliRun.h>
  #include <AliTPCpidESD.h>
  #include <AliESDEvent.h>
#endif

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

static Int_t allpisel=0;
static Int_t allkasel=0;
static Int_t allprsel=0;
static Int_t allnsel=0;


Double_t tpcBethe(Double_t *x, Double_t *par) {
  //
  // This is the TPC Bethe-Bloch function normalised to 1 at the minimum
  //
  Double_t p=x[0];   //momentum
  Double_t m=par[0]; //mass
  Double_t bg=p/m;
  return AliTPCpidESD::Bethe(bg);  
}


Int_t AliPIDComparison(const Char_t *dir=".") { 
   gBenchmark->Start("AliPIDComparison");

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

   Double_t pi=0.2,pa=3;

   TH2F *tpcHist=(TH2F*)gROOT->FindObject("tpcHist");
   if (!tpcHist)
     tpcHist=new TH2F("tpcHist","TPC dE/dX vs momentum",100,pi,pa,100,0.,6.);
   tpcHist->SetXTitle("p (GeV/c)"); tpcHist->SetYTitle("dE/dx (Arb. Units)");
   tpcHist->SetMarkerStyle(8); 
   tpcHist->SetMarkerSize(0.3);
 
   const Char_t *hname[]={
    "piG","piR","piF","piGood","piFake",
    "kaG","kaR","kaF","kaGood","kaFake",
    "prG","prR","prF","prGood","prFake"
   };
   Int_t nh=sizeof(hname)/sizeof(const Char_t *);
   TH1F **hprt=new TH1F*[nh]; 

   for (Int_t i=0; i<nh; i++) {
     hprt[i]=(TH1F*)gROOT->FindObject(hname[i]);
     if (hprt[i]==0) {hprt[i]=new TH1F(hname[i],"",20,pi,pa);hprt[i]->Sumw2();}
   }
   TH1F *piG=hprt[0];
   TH1F *piR=hprt[1];
   TH1F *piF=hprt[2];
   TH1F *piGood=hprt[3]; 
        piGood->SetTitle("Combined PID for pions");
        piGood->SetLineColor(4); piGood->SetXTitle("p (GeV/c)");
   TH1F *piFake=hprt[4]; 
        piFake->SetLineColor(2);
   TH1F *kaG=hprt[5];
   TH1F *kaR=hprt[6];
   TH1F *kaF=hprt[7];
   TH1F *kaGood=hprt[8];
        kaGood->SetTitle("Combined PID for kaons"); 
        kaGood->SetLineColor(4); kaGood->SetXTitle("p (GeV/c)");
   TH1F *kaFake=hprt[9]; 
        kaFake->SetLineColor(2);
   TH1F *prG=hprt[10];
   TH1F *prR=hprt[11];
   TH1F *prF=hprt[12];
   TH1F *prGood=hprt[13];
        prGood->SetTitle("Combined PID for protons"); 
        prGood->SetLineColor(4); prGood->SetXTitle("p (GeV/c)");
   TH1F *prFake=hprt[14]; 
        prFake->SetLineColor(2);

   delete[] hprt;

   Char_t fname[100];

   if (gAlice) {
      delete AliRunLoader::Instance();
      delete gAlice;
      gAlice=0;
   }
   sprintf(fname,"%s/galice.root",dir);
   AliRunLoader *rl = AliRunLoader::Open(fname);
   if (rl == 0x0) {
      ::Error("AliPIDComparison.C","Can not open session !");
      return 1;
   }
   if (rl->LoadgAlice()) {
      ::Error("AliPIDComparison.C","LoadgAlice returned error !");
      delete rl;
      return 1;
   }
   if (rl->LoadHeader()) {
      ::Error("AliPIDComparison.C","LoadHeader returned error !");
      delete rl;
      return 1;
   }
   rl->LoadKinematics();

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


   //****** Tentative particle type "concentrations"
   Double_t c[5]={0.01, 0.01, 0.85, 0.10, 0.05};
   AliPID::SetPriors(c);


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

      rl->GetEvent(e);
 
      e++;

      Int_t ntrk=event->GetNumberOfTracks();
      cerr<<"Number of ESD tracks : "<<ntrk<<endl; 

      Int_t pisel=0,kasel=0,prsel=0,nsel=0;

      AliStack *stack = rl->Stack();

      while (ntrk--) {
        AliESDtrack *t=event->GetTrack(ntrk);

	//*** Some track quality cuts ****
        if (t->GetITSclusters(0) < 6 ) continue;
        if (t->GetTPCclusters(0) < 60) continue;
        //if (t->GetTRDclusters(0) < 60) continue;

        if (!t->IsOn(AliESDtrack::kITSpid)) continue;
        if (!t->IsOn(AliESDtrack::kTPCpid)) continue;
        //if (!t->IsOn(AliESDtrack::kTRDpid)) continue;
        if (!t->IsOn(AliESDtrack::kTOFpid)) continue;
        {
           nsel++;

           Double_t p=t->GetInnerParam()->GetP();
	   Double_t dedx=t->GetTPCsignal()/47.;
	   tpcHist->Fill(p,dedx,1);

           Int_t lab=TMath::Abs(t->GetLabel());
           TParticle *part=stack->Particle(lab);
           Int_t code=part->GetPdgCode();

           Double_t r[10]; t->GetESDpid(r);

           AliPID pid(r);

           Double_t w[10];
           w[0]=pid.GetProbability(AliPID::kElectron);
           w[1]=pid.GetProbability(AliPID::kMuon);
           w[2]=pid.GetProbability(AliPID::kPion);
           w[3]=pid.GetProbability(AliPID::kKaon);
           w[4]=pid.GetProbability(AliPID::kProton);

           if (TMath::Abs(code)==2212) prR->Fill(p);
           if (w[4]>w[3] && w[4]>w[2] && w[4]>w[1] && w[4]>w[0]) {//proton
	      prsel++;
	      prG->Fill(p);
              if (TMath::Abs(code)!=2212) prF->Fill(p);
           }

           if (TMath::Abs(code)==321) kaR->Fill(p);
           if (w[3]>w[4] && w[3]>w[2] && w[3]>w[1] && w[3]>w[0]) {//kaon
	      kasel++;
	      kaG->Fill(p);
              if (TMath::Abs(code)!=321) kaF->Fill(p);
           }

	   if (TMath::Abs(code)==211) piR->Fill(p);
	   if (w[2]>w[4] && w[2]>w[3] && w[2]>w[0] && w[2]>w[1]) {//pion
	      pisel++;
	      piG->Fill(p);
	      if (TMath::Abs(code)!=211) piF->Fill(p);
           }
	   /*
           if (TMath::Abs(code)==11) piR->Fill(p);
           if (w[0]>w[4] && w[0]>w[3] && w[0]>w[3] && w[0]>w[1]) {//electron
	      pisel++;
	      piG->Fill(p);
              if (TMath::Abs(code)!=11) piF->Fill(p);
           }
	   */
	}
      }
      cout<<"Number of selected ESD tracks : "<<nsel<<endl;
      cout<<"Number of selected pion ESD tracks : "<<pisel<<endl;
      cout<<"Number of selected kaon ESD tracks : "<<kasel<<endl;
      cout<<"Number of selected proton ESD tracks : "<<prsel<<endl;

      allnsel+=nsel; allpisel+=pisel; allkasel+=kasel; allprsel+=prsel;

   } // ***** End of the loop over events

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

   TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
   if (!c1) {
      c1=new TCanvas("c1","",0,0,600,1000);
      c1->Divide(1,4);
   }

   c1->cd(1);
   tpcHist->Draw();

   TF1 *fpi=(TF1*)gROOT->FindObject("piBethe");
   if (!fpi) {
      fpi=new TF1("piBethe",tpcBethe,pi,pa,1);
      fpi->SetLineColor(2);
      fpi->SetParameter(0,AliPID::ParticleMass(AliPID::kPion));
   }
   fpi->Draw("same");

   TF1 *fka=(TF1*)gROOT->FindObject("kaBethe");
   if (!fka) {
      fka=new TF1("kaBethe",tpcBethe,pi,pa,1);
      fka->SetLineColor(3);
      fka->SetParameter(0,AliPID::ParticleMass(AliPID::kKaon));
   }
   fka->Draw("same");

   TF1 *fpr=(TF1*)gROOT->FindObject("prBethe");
   if (!fpr) {
      fpr=new TF1("prBethe",tpcBethe,pi,pa,1);
      fpr->SetLineColor(4);
      fpr->SetParameter(0,AliPID::ParticleMass(AliPID::kProton));
   }
   fpr->Draw("same");


   c1->cd(2);
   //piG->Sumw2(); piF->Sumw2(); piR->Sumw2();
   piGood->Add(piG,piF,1,-1);
   piGood->Divide(piGood,piR,1,1,"b");
   piFake->Divide(piF,piG,1,1,"b");
   piGood->SetMinimum(0);
   piGood->Draw("hist");
   piFake->Draw("same");

   c1->cd(3);
   //kaG->Sumw2(); kaF->Sumw2(); kaR->Sumw2();
   kaGood->Add(kaG,kaF,1,-1);
   kaGood->Divide(kaGood,kaR,1,1,"b");
   kaFake->Divide(kaF,kaG,1,1,"b");
   kaGood->SetMinimum(0);
   kaGood->Draw("hist");
   kaFake->Draw("same");

   c1->cd(4);
   //prG->Sumw2(); prF->Sumw2(); prR->Sumw2();
   prGood->Add(prG,prF,1,-1);
   prGood->Divide(prGood,prR,1,1,"b");
   prFake->Divide(prF,prG,1,1,"b");
   prGood->SetMinimum(0);
   prGood->Draw("hist");
   prFake->Draw("same");

   c1->Update();

   cout<<endl;
   cout<<endl;
   e=(Int_t)piR->GetEntries();
   Int_t o=(Int_t)piG->GetEntries();
   if (e*o) cout<<"Efficiency (contamination) for pions : "<<
      piGood->GetEntries()/e<<'('<<piF->GetEntries()/o<<')'<<endl; 
   e=(Int_t)kaR->GetEntries();
   o=(Int_t)kaG->GetEntries();
   if (e*o) cout<<"Efficiency (contamination) for kaons : "<<
      kaGood->GetEntries()/e<<'('<<kaF->GetEntries()/o<<')'<<endl; 
   e=(Int_t)prR->GetEntries();
   o=(Int_t)prG->GetEntries();
   if (e*o) cout<<"Efficiency (contamination) for protons : "<<
      prGood->GetEntries()/e<<'('<<prF->GetEntries()/o<<')'<<endl;
   cout<<endl; 

   cout<<"Total number of selected ESD tracks : "<<allnsel<<endl;
   cout<<"Total number of selected pion ESD tracks : "<<allpisel<<endl;
   cout<<"Total number of selected kaon ESD tracks : "<<allkasel<<endl;
   cout<<"Total number of selected proton ESD tracks : "<<allprsel<<endl;

   ef->Close();
   TFile fc("AliPIDComparison.root","RECREATE");
   c1->Write();
   fc.Close();

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

   delete rl;
   return 0;
}
 AliPIDComparison.C:1
 AliPIDComparison.C:2
 AliPIDComparison.C:3
 AliPIDComparison.C:4
 AliPIDComparison.C:5
 AliPIDComparison.C:6
 AliPIDComparison.C:7
 AliPIDComparison.C:8
 AliPIDComparison.C:9
 AliPIDComparison.C:10
 AliPIDComparison.C:11
 AliPIDComparison.C:12
 AliPIDComparison.C:13
 AliPIDComparison.C:14
 AliPIDComparison.C:15
 AliPIDComparison.C:16
 AliPIDComparison.C:17
 AliPIDComparison.C:18
 AliPIDComparison.C:19
 AliPIDComparison.C:20
 AliPIDComparison.C:21
 AliPIDComparison.C:22
 AliPIDComparison.C:23
 AliPIDComparison.C:24
 AliPIDComparison.C:25
 AliPIDComparison.C:26
 AliPIDComparison.C:27
 AliPIDComparison.C:28
 AliPIDComparison.C:29
 AliPIDComparison.C:30
 AliPIDComparison.C:31
 AliPIDComparison.C:32
 AliPIDComparison.C:33
 AliPIDComparison.C:34
 AliPIDComparison.C:35
 AliPIDComparison.C:36
 AliPIDComparison.C:37
 AliPIDComparison.C:38
 AliPIDComparison.C:39
 AliPIDComparison.C:40
 AliPIDComparison.C:41
 AliPIDComparison.C:42
 AliPIDComparison.C:43
 AliPIDComparison.C:44
 AliPIDComparison.C:45
 AliPIDComparison.C:46
 AliPIDComparison.C:47
 AliPIDComparison.C:48
 AliPIDComparison.C:49
 AliPIDComparison.C:50
 AliPIDComparison.C:51
 AliPIDComparison.C:52
 AliPIDComparison.C:53
 AliPIDComparison.C:54
 AliPIDComparison.C:55
 AliPIDComparison.C:56
 AliPIDComparison.C:57
 AliPIDComparison.C:58
 AliPIDComparison.C:59
 AliPIDComparison.C:60
 AliPIDComparison.C:61
 AliPIDComparison.C:62
 AliPIDComparison.C:63
 AliPIDComparison.C:64
 AliPIDComparison.C:65
 AliPIDComparison.C:66
 AliPIDComparison.C:67
 AliPIDComparison.C:68
 AliPIDComparison.C:69
 AliPIDComparison.C:70
 AliPIDComparison.C:71
 AliPIDComparison.C:72
 AliPIDComparison.C:73
 AliPIDComparison.C:74
 AliPIDComparison.C:75
 AliPIDComparison.C:76
 AliPIDComparison.C:77
 AliPIDComparison.C:78
 AliPIDComparison.C:79
 AliPIDComparison.C:80
 AliPIDComparison.C:81
 AliPIDComparison.C:82
 AliPIDComparison.C:83
 AliPIDComparison.C:84
 AliPIDComparison.C:85
 AliPIDComparison.C:86
 AliPIDComparison.C:87
 AliPIDComparison.C:88
 AliPIDComparison.C:89
 AliPIDComparison.C:90
 AliPIDComparison.C:91
 AliPIDComparison.C:92
 AliPIDComparison.C:93
 AliPIDComparison.C:94
 AliPIDComparison.C:95
 AliPIDComparison.C:96
 AliPIDComparison.C:97
 AliPIDComparison.C:98
 AliPIDComparison.C:99
 AliPIDComparison.C:100
 AliPIDComparison.C:101
 AliPIDComparison.C:102
 AliPIDComparison.C:103
 AliPIDComparison.C:104
 AliPIDComparison.C:105
 AliPIDComparison.C:106
 AliPIDComparison.C:107
 AliPIDComparison.C:108
 AliPIDComparison.C:109
 AliPIDComparison.C:110
 AliPIDComparison.C:111
 AliPIDComparison.C:112
 AliPIDComparison.C:113
 AliPIDComparison.C:114
 AliPIDComparison.C:115
 AliPIDComparison.C:116
 AliPIDComparison.C:117
 AliPIDComparison.C:118
 AliPIDComparison.C:119
 AliPIDComparison.C:120
 AliPIDComparison.C:121
 AliPIDComparison.C:122
 AliPIDComparison.C:123
 AliPIDComparison.C:124
 AliPIDComparison.C:125
 AliPIDComparison.C:126
 AliPIDComparison.C:127
 AliPIDComparison.C:128
 AliPIDComparison.C:129
 AliPIDComparison.C:130
 AliPIDComparison.C:131
 AliPIDComparison.C:132
 AliPIDComparison.C:133
 AliPIDComparison.C:134
 AliPIDComparison.C:135
 AliPIDComparison.C:136
 AliPIDComparison.C:137
 AliPIDComparison.C:138
 AliPIDComparison.C:139
 AliPIDComparison.C:140
 AliPIDComparison.C:141
 AliPIDComparison.C:142
 AliPIDComparison.C:143
 AliPIDComparison.C:144
 AliPIDComparison.C:145
 AliPIDComparison.C:146
 AliPIDComparison.C:147
 AliPIDComparison.C:148
 AliPIDComparison.C:149
 AliPIDComparison.C:150
 AliPIDComparison.C:151
 AliPIDComparison.C:152
 AliPIDComparison.C:153
 AliPIDComparison.C:154
 AliPIDComparison.C:155
 AliPIDComparison.C:156
 AliPIDComparison.C:157
 AliPIDComparison.C:158
 AliPIDComparison.C:159
 AliPIDComparison.C:160
 AliPIDComparison.C:161
 AliPIDComparison.C:162
 AliPIDComparison.C:163
 AliPIDComparison.C:164
 AliPIDComparison.C:165
 AliPIDComparison.C:166
 AliPIDComparison.C:167
 AliPIDComparison.C:168
 AliPIDComparison.C:169
 AliPIDComparison.C:170
 AliPIDComparison.C:171
 AliPIDComparison.C:172
 AliPIDComparison.C:173
 AliPIDComparison.C:174
 AliPIDComparison.C:175
 AliPIDComparison.C:176
 AliPIDComparison.C:177
 AliPIDComparison.C:178
 AliPIDComparison.C:179
 AliPIDComparison.C:180
 AliPIDComparison.C:181
 AliPIDComparison.C:182
 AliPIDComparison.C:183
 AliPIDComparison.C:184
 AliPIDComparison.C:185
 AliPIDComparison.C:186
 AliPIDComparison.C:187
 AliPIDComparison.C:188
 AliPIDComparison.C:189
 AliPIDComparison.C:190
 AliPIDComparison.C:191
 AliPIDComparison.C:192
 AliPIDComparison.C:193
 AliPIDComparison.C:194
 AliPIDComparison.C:195
 AliPIDComparison.C:196
 AliPIDComparison.C:197
 AliPIDComparison.C:198
 AliPIDComparison.C:199
 AliPIDComparison.C:200
 AliPIDComparison.C:201
 AliPIDComparison.C:202
 AliPIDComparison.C:203
 AliPIDComparison.C:204
 AliPIDComparison.C:205
 AliPIDComparison.C:206
 AliPIDComparison.C:207
 AliPIDComparison.C:208
 AliPIDComparison.C:209
 AliPIDComparison.C:210
 AliPIDComparison.C:211
 AliPIDComparison.C:212
 AliPIDComparison.C:213
 AliPIDComparison.C:214
 AliPIDComparison.C:215
 AliPIDComparison.C:216
 AliPIDComparison.C:217
 AliPIDComparison.C:218
 AliPIDComparison.C:219
 AliPIDComparison.C:220
 AliPIDComparison.C:221
 AliPIDComparison.C:222
 AliPIDComparison.C:223
 AliPIDComparison.C:224
 AliPIDComparison.C:225
 AliPIDComparison.C:226
 AliPIDComparison.C:227
 AliPIDComparison.C:228
 AliPIDComparison.C:229
 AliPIDComparison.C:230
 AliPIDComparison.C:231
 AliPIDComparison.C:232
 AliPIDComparison.C:233
 AliPIDComparison.C:234
 AliPIDComparison.C:235
 AliPIDComparison.C:236
 AliPIDComparison.C:237
 AliPIDComparison.C:238
 AliPIDComparison.C:239
 AliPIDComparison.C:240
 AliPIDComparison.C:241
 AliPIDComparison.C:242
 AliPIDComparison.C:243
 AliPIDComparison.C:244
 AliPIDComparison.C:245
 AliPIDComparison.C:246
 AliPIDComparison.C:247
 AliPIDComparison.C:248
 AliPIDComparison.C:249
 AliPIDComparison.C:250
 AliPIDComparison.C:251
 AliPIDComparison.C:252
 AliPIDComparison.C:253
 AliPIDComparison.C:254
 AliPIDComparison.C:255
 AliPIDComparison.C:256
 AliPIDComparison.C:257
 AliPIDComparison.C:258
 AliPIDComparison.C:259
 AliPIDComparison.C:260
 AliPIDComparison.C:261
 AliPIDComparison.C:262
 AliPIDComparison.C:263
 AliPIDComparison.C:264
 AliPIDComparison.C:265
 AliPIDComparison.C:266
 AliPIDComparison.C:267
 AliPIDComparison.C:268
 AliPIDComparison.C:269
 AliPIDComparison.C:270
 AliPIDComparison.C:271
 AliPIDComparison.C:272
 AliPIDComparison.C:273
 AliPIDComparison.C:274
 AliPIDComparison.C:275
 AliPIDComparison.C:276
 AliPIDComparison.C:277
 AliPIDComparison.C:278
 AliPIDComparison.C:279
 AliPIDComparison.C:280
 AliPIDComparison.C:281
 AliPIDComparison.C:282
 AliPIDComparison.C:283
 AliPIDComparison.C:284
 AliPIDComparison.C:285
 AliPIDComparison.C:286
 AliPIDComparison.C:287
 AliPIDComparison.C:288
 AliPIDComparison.C:289
 AliPIDComparison.C:290
 AliPIDComparison.C:291
 AliPIDComparison.C:292
 AliPIDComparison.C:293
 AliPIDComparison.C:294
 AliPIDComparison.C:295
 AliPIDComparison.C:296
 AliPIDComparison.C:297
 AliPIDComparison.C:298
 AliPIDComparison.C:299
 AliPIDComparison.C:300
 AliPIDComparison.C:301
 AliPIDComparison.C:302
 AliPIDComparison.C:303
 AliPIDComparison.C:304
 AliPIDComparison.C:305
 AliPIDComparison.C:306
 AliPIDComparison.C:307
 AliPIDComparison.C:308
 AliPIDComparison.C:309
 AliPIDComparison.C:310
 AliPIDComparison.C:311
 AliPIDComparison.C:312
 AliPIDComparison.C:313
 AliPIDComparison.C:314
 AliPIDComparison.C:315
 AliPIDComparison.C:316
 AliPIDComparison.C:317
 AliPIDComparison.C:318
 AliPIDComparison.C:319
 AliPIDComparison.C:320
 AliPIDComparison.C:321
 AliPIDComparison.C:322
 AliPIDComparison.C:323
 AliPIDComparison.C:324
 AliPIDComparison.C:325
 AliPIDComparison.C:326
 AliPIDComparison.C:327
 AliPIDComparison.C:328
 AliPIDComparison.C:329
 AliPIDComparison.C:330
 AliPIDComparison.C:331
 AliPIDComparison.C:332
 AliPIDComparison.C:333
 AliPIDComparison.C:334
 AliPIDComparison.C:335
 AliPIDComparison.C:336
 AliPIDComparison.C:337
 AliPIDComparison.C:338