ROOT logo
computePriors(){
  char filename[100];
  sprintf(filename,"yourfile.root");
  TFile *f = new TFile(filename);
  TTree *t = f->Get("yourtree");

  Float_t mass[7] = {0.000511,0.139570,0.493677,0.938272,1.877837,2.817402,1.408054};

  Int_t nev = t->GetEntries();

  printf("ntracks = %i\n",nev);

  const Int_t maxstep = 10;

  TH1D *hpriors[7][maxstep];
  TH1D *hratio[7];

  for(Int_t i=0;i<7;i++){
    for(Int_t j=0;j<maxstep;j++){
      char histoname[100];
      sprintf(histoname,"priors%istep%i",i,j);
      hpriors[i][j] = new TH1D(histoname,histoname,10,0,10);
      if(j==0){
	sprintf(histoname,"ratio%i",i);
	hratio[i] = new TH1D(histoname,histoname,10,0,10);
	for(Int_t k=1;k <= 100;k++){
	  hpriors[i][j]->SetBinContent(k,1);
	}
      }
    }
  }

  for(Int_t isp=0;isp < 7;isp++) hratio[isp]->Sumw2();

  // loop variable definition 
  Float_t pt, eta, weight=1.0;
  Float_t ppi, pka, ppr, pel, pmu, pde, ptr, phe, ptot;

  for(Int_t j=1; j <maxstep;j++){
    for(Int_t isp=0;isp<7;isp++) hratio[isp]->Reset();
    printf("step %i\n",j);
    for(Int_t i=0;i<nev;i++){
      t->GetEvent(i);
      Int_t sw = t->GetLeaf("kTOF")->GetValue() * t->GetLeaf("kTPC")->GetValue(); 
      if(sw){
	pt = t->GetLeaf("pt")->GetValue();
	eta = t->GetLeaf("eta")->GetValue();
	//	weight = 1.0; // if you have fill the tree with weights
	if(pt > 0 && pt < 10){
	  ppi = t->GetLeaf("tofW")->GetValue(2)*t->GetLeaf("tpcW")->GetValue(2)*hpriors[1][j-1]->Interpolate(pt);
	  pka = t->GetLeaf("tofW")->GetValue(3)*t->GetLeaf("tpcW")->GetValue(3)*hpriors[2][j-1]->Interpolate(pt);
	  ppr = t->GetLeaf("tofW")->GetValue(4)*t->GetLeaf("tpcW")->GetValue(4)*hpriors[3][j-1]->Interpolate(pt);
	  pel = t->GetLeaf("tofW")->GetValue(0)*t->GetLeaf("tpcW")->GetValue(0)*hpriors[0][j-1]->Interpolate(pt);
	  pmu = t->GetLeaf("tofW")->GetValue(1)*t->GetLeaf("tpcW")->GetValue(1)*hpriors[0][j-1]->Interpolate(pt);
	  pde = t->GetLeaf("tofW")->GetValue(5)*t->GetLeaf("tpcW")->GetValue(5)*hpriors[4][j-1]->Interpolate(pt);
	  ptr = t->GetLeaf("tofW")->GetValue(6)*t->GetLeaf("tpcW")->GetValue(6)*hpriors[5][j-1]->Interpolate(pt);
	  phe = t->GetLeaf("tofW")->GetValue(7)*t->GetLeaf("tpcW")->GetValue(7)*hpriors[6][j-1]->Interpolate(pt);
	  ptot = ppi+pka+ppr+pel+pde+ptr+phe+pmu;

	  if(ptot > 0){
	    // fill ratio for |y|<0.5
	    if(TMath::Abs(eta2y(pt,mass[1],eta)) < 0.5) hratio[1]->Fill(pt,ppi/ptot*weight);
	    if(TMath::Abs(eta2y(pt,mass[2],eta)) < 0.5) hratio[2]->Fill(pt,pka/ptot*weight);
	    if(TMath::Abs(eta2y(pt,mass[3],eta)) < 0.5) hratio[3]->Fill(pt,ppr/ptot*weight);
	    if(TMath::Abs(eta2y(pt,mass[0],eta)) < 0.5) hratio[0]->Fill(pt,pel/ptot*weight);
	    if(TMath::Abs(eta2y(pt,mass[4],eta)) < 0.5) hratio[4]->Fill(pt,pde/ptot*weight);
	    if(TMath::Abs(eta2y(pt,mass[5],eta)) < 0.5) hratio[5]->Fill(pt,ptr/ptot*weight);
	    if(TMath::Abs(eta2y(pt,mass[6],eta)) < 0.5) hratio[6]->Fill(pt,phe/ptot*weight);

	    // fill priors assuming #eta cut during the fill of the tree
	    hpriors[1][j]->Fill(pt,ppi/ptot*weight);
	    hpriors[2][j]->Fill(pt,pka/ptot*weight);
	    hpriors[3][j]->Fill(pt,ppr/ptot*weight);
	    hpriors[0][j]->Fill(pt,pel/ptot*weight);
	    hpriors[4][j]->Fill(pt,pde/ptot*weight);
	    hpriors[5][j]->Fill(pt,ptr/ptot*weight);
	    hpriors[6][j]->Fill(pt,phe/ptot*weight);
	  }
	}
      }
    }

    for(Int_t isp=0;isp < 7;isp++) hpriors[isp][j]->Sumw2();

    // Normalize abundances to the pion one
    hpriors[0][j]->Divide(hpriors[1][j]);
    hpriors[2][j]->Divide(hpriors[1][j]);
    hpriors[3][j]->Divide(hpriors[1][j]);
    hpriors[4][j]->Divide(hpriors[1][j]);
    hpriors[5][j]->Divide(hpriors[1][j]);
    hpriors[6][j]->Divide(hpriors[1][j]);
    hpriors[1][j]->Divide(hpriors[1][j]);

    // make ratios w.r.t. pions (|y|<0.5)
    hratio[0]->Divide(hratio[1]);
    hratio[2]->Divide(hratio[1]);
    hratio[3]->Divide(hratio[1]);
    hratio[4]->Divide(hratio[1]);
    hratio[5]->Divide(hratio[1]);
    hratio[6]->Divide(hratio[1]);
    hratio[1]->Divide(hratio[1]);
  }

  // Write output
  sprintf(filename,"priors.root");
  TFile *fout = new TFile(filename,"RECREATE");
  for(Int_t j=maxstep-2; j <maxstep;j++){ // last 2 steps
    hpriors[0][j]->Write();
    hpriors[1][j]->Write();
    hpriors[2][j]->Write();
    hpriors[3][j]->Write();
    hpriors[4][j]->Write();
    hpriors[5][j]->Write();
    hpriors[6][j]->Write();
  }
  hratio[0]->Write();
  hratio[1]->Write();
  hratio[2]->Write();
  hratio[3]->Write();
  hratio[4]->Write();
  hratio[5]->Write();
  hratio[6]->Write();
  fout->Close();
}
Double_t y2eta(Double_t pt, Double_t m, Double_t y){
  Double_t mt=TMath::Sqrt(m*m+pt*pt);
  return TMath::ASinH(mt/pt*TMath::SinH(y));
}
Double_t eta2y(Double_t pt, Double_t m, Double_t eta){
  Double_t mt=TMath::Sqrt(m*m+pt*pt);
  return TMath::ASinH(pt/mt*TMath::SinH(eta));
}
 computePriors.C:1
 computePriors.C:2
 computePriors.C:3
 computePriors.C:4
 computePriors.C:5
 computePriors.C:6
 computePriors.C:7
 computePriors.C:8
 computePriors.C:9
 computePriors.C:10
 computePriors.C:11
 computePriors.C:12
 computePriors.C:13
 computePriors.C:14
 computePriors.C:15
 computePriors.C:16
 computePriors.C:17
 computePriors.C:18
 computePriors.C:19
 computePriors.C:20
 computePriors.C:21
 computePriors.C:22
 computePriors.C:23
 computePriors.C:24
 computePriors.C:25
 computePriors.C:26
 computePriors.C:27
 computePriors.C:28
 computePriors.C:29
 computePriors.C:30
 computePriors.C:31
 computePriors.C:32
 computePriors.C:33
 computePriors.C:34
 computePriors.C:35
 computePriors.C:36
 computePriors.C:37
 computePriors.C:38
 computePriors.C:39
 computePriors.C:40
 computePriors.C:41
 computePriors.C:42
 computePriors.C:43
 computePriors.C:44
 computePriors.C:45
 computePriors.C:46
 computePriors.C:47
 computePriors.C:48
 computePriors.C:49
 computePriors.C:50
 computePriors.C:51
 computePriors.C:52
 computePriors.C:53
 computePriors.C:54
 computePriors.C:55
 computePriors.C:56
 computePriors.C:57
 computePriors.C:58
 computePriors.C:59
 computePriors.C:60
 computePriors.C:61
 computePriors.C:62
 computePriors.C:63
 computePriors.C:64
 computePriors.C:65
 computePriors.C:66
 computePriors.C:67
 computePriors.C:68
 computePriors.C:69
 computePriors.C:70
 computePriors.C:71
 computePriors.C:72
 computePriors.C:73
 computePriors.C:74
 computePriors.C:75
 computePriors.C:76
 computePriors.C:77
 computePriors.C:78
 computePriors.C:79
 computePriors.C:80
 computePriors.C:81
 computePriors.C:82
 computePriors.C:83
 computePriors.C:84
 computePriors.C:85
 computePriors.C:86
 computePriors.C:87
 computePriors.C:88
 computePriors.C:89
 computePriors.C:90
 computePriors.C:91
 computePriors.C:92
 computePriors.C:93
 computePriors.C:94
 computePriors.C:95
 computePriors.C:96
 computePriors.C:97
 computePriors.C:98
 computePriors.C:99
 computePriors.C:100
 computePriors.C:101
 computePriors.C:102
 computePriors.C:103
 computePriors.C:104
 computePriors.C:105
 computePriors.C:106
 computePriors.C:107
 computePriors.C:108
 computePriors.C:109
 computePriors.C:110
 computePriors.C:111
 computePriors.C:112
 computePriors.C:113
 computePriors.C:114
 computePriors.C:115
 computePriors.C:116
 computePriors.C:117
 computePriors.C:118
 computePriors.C:119
 computePriors.C:120
 computePriors.C:121
 computePriors.C:122
 computePriors.C:123
 computePriors.C:124
 computePriors.C:125
 computePriors.C:126
 computePriors.C:127
 computePriors.C:128
 computePriors.C:129
 computePriors.C:130
 computePriors.C:131
 computePriors.C:132
 computePriors.C:133