ROOT logo
void DrawSliceProton(Double_t p, Bool_t singlePiGauss){
  //TVirtualFitter::Set
  TFile *in = TFile::Open("HFEtask.root");
  TList *qa = (TList *)in->Get("HFE_QA");

  // Make Plots for TPC
  TList *pidqa = (TList *)qa->FindObject("HFEpidQA");
  AliHFEtpcPIDqa *tpcqa = (AliHFEtpcPIDqa *)pidqa->FindObject("TPCQA");
  TH2 *hTPCsig = tpcqa->MakeSpectrumNSigma(AliHFEdetPIDqa::kBeforePID);

  // define cut model
  TF1 cutmodel("cutmodel", "[0] * TMath::Exp([1]*x) + [2]", 0, 20);
  cutmodel.SetParameter(0, -2.75);
  cutmodel.SetParameter(1, -0.8757);
  cutmodel.SetParameter(2, -0.9);

  Int_t nBinsP = 2;
  Int_t pbin = hTPCsig->GetXaxis()->FindBin(p);
  TH1 *hslice = hTPCsig->ProjectionY("hslice", pbin - nBinsP, pbin + nBinsP);

  // Functions needed
  TF1 *pi1 = new TF1("pi1", "[0] * TMath::Gaus(x, [1], [2]) + [3]*TMath::Gaus(x, [4], [5]) + [6]*TMath::Gaus(x, [7], [8])", -20, 20),
      *pi2 = new TF1("pi2", "[0] * TMath::Gaus(x, [1], [2]) + [3]*TMath::Gaus(x, [4], [5]) + [6]*TMath::Gaus(x, [7], [8])", -20, 20),
      *el1 = new TF1("el1", "[0] * TMath::Gaus(x, [1], [2])", -20, 20),
      *el2 = new TF1("el1", "[0] * TMath::Gaus(x, [1], [2])", -20, 20),
      *combined = new TF1("combined",  "[0] * TMath::Gaus(x, [1], [2]) + [3]*TMath::Gaus(x, [4], [5]) +[6] * TMath::Gaus(x, [7], [8]) + [9]*TMath::Gaus(x, [10], [11])" , -20, 20);

  // Constraints
  Double_t minele = 1e0;
  Double_t maxele = 1e3;
  Double_t minpr = 1e1; 
  Double_t maxpr = 1e4;
  Double_t minpi = 1e2;
  Double_t maxpi = 1e5;
  pi1->SetParLimits(0, minpi, maxpi); pi1->SetParLimits(3, minpi, maxpi);
  pi1->SetParLimits(2, 0.4, 1.5); pi1->SetParLimits(5, 0.4, 1.5);
  pi1->SetParLimits(1, -8,-4); pi1->SetParLimits(4, -5.5,-3);
  if(singlePiGauss){
    pi1->FixParameter(3,0);
    pi1->FixParameter(4,0);
    pi1->FixParameter(5,0);
  }
  // Set Protons
  pi1->SetParLimits(6, minpr, maxpr);
  pi1->SetParLimits(7, -7., -5.8);
  pi1->SetParLimits(8, 0.1, 1.3);
  // Set Electrons
  el1->SetParLimits(0, minele, maxele); 
  el1->SetParLimits(1, -0.9., 0.);
  el1->SetParLimits(2, 0.9, 1.5); 
  el1->SetParameter(1, -0.1);
  //el1->FixParameter(1, -1.34726e-01);
  //el1->FixParameter(2, 9.77243e-01);
  hslice->Fit(pi1, "N", "", -8, -1);
  hslice->Fit(el1, "NL", "", -.8, 3);

  // Use single fits to constrain combined fit
  combined->FixParameter(1, pi1->GetParameter(1));
  combined->FixParameter(2, pi1->GetParameter(2));
  combined->FixParameter(4, pi1->GetParameter(4));
  combined->FixParameter(5, pi1->GetParameter(5));
  if(singlePiGauss) 
    combined->FixParameter(3, 0);
  else
    combined->SetParLimits(3, minpi, maxpi);
  // Protons
  combined->FixParameter(7, pi1->GetParameter(7));
  combined->FixParameter(8, pi1->GetParameter(8));
  // Electrons
  combined->FixParameter(10, el1->GetParameter(1));
  combined->FixParameter(11, el1->GetParameter(2));
  combined->SetParLimits(0, minpi, maxpi);
  combined->SetParLimits(6, minpr, maxpr);
  combined->SetParLimits(9, minele, maxele);

  hslice->Fit(combined, "N", "", -8., 3.);

  for(Int_t ipar = 0; ipar < 9; ipar++) pi2->SetParameter(ipar, combined->GetParameter(ipar));
  for(Int_t ipar = 0; ipar < 3; ipar++) el2->SetParameter(ipar, combined->GetParameter(ipar+9));

  combined->SetLineColor(kBlack);
  el2->SetLineColor(kGreen);
  pi2->SetLineColor(kBlue);
  combined->SetLineWidth(2);
  el2->SetLineWidth(2);
  pi2->SetLineWidth(2);
  el2->SetLineStyle(2);
  pi2->SetLineStyle(2);

  hslice->SetStats(kFALSE);
  hslice->SetTitle("TPC Signal");
  hslice->GetXaxis()->SetTitle("TPC dE/dx - <dE/dx>|_{electrons} [#sigma]");

  TCanvas *output = new TCanvas("output", "Slice output");
  output->SetLogy();
  output->cd();
  hslice->Draw();
  combined->Draw("same");
  el2->Draw("same");
  pi2->Draw("same");
  
  TLegend *leg = new TLegend(0.5, 0.7, 0.89, 0.89);
  leg->SetFillStyle(0);
  leg->SetBorderSize(0);
  leg->AddEntry(hslice, "Measurement", "l");
  leg->AddEntry(combined, "Electrons + Pions", "l");
  leg->AddEntry(el2, "Electrons", "l");
  leg->AddEntry(pi2, "Pions", "l");
  leg->Draw();

  TPaveText *momentum = new TPaveText(0.6, 0.45, 0.8, 0.55, "NDC");
  momentum->SetBorderSize(0);
  momentum->SetFillStyle(0);
  momentum->AddText(Form("p = %.2fGeV/c", p));
  momentum->Draw();

  // Draw also Selection Bands
  TF1 *fBandPions = new TF1(*pi2), *fBandElectrons = new TF1(*el2);
  fBandPions->SetRange(cutmodel.Eval(p), 5);
  fBandElectrons->SetRange(cutmodel.Eval(p), 5);
  
  fBandPions->SetFillColor(kBlue);
  fBandPions->SetFillStyle(3004);
  fBandElectrons->SetFillColor(kGreen);
  fBandElectrons->SetFillStyle(3005);
  fBandPions->SetLineWidth(0);
  fBandElectrons->SetLineWidth(0);
  
  fBandElectrons->Draw("same");
  fBandPions->Draw("same");

  // Now calculate the contamination
  Double_t nPions = pi2->Integral(cutmodel.Eval(p), 5);
  Double_t nCandidates = combined->Integral(cutmodel.Eval(p), 5);
  Double_t contamination = nPions/nCandidates;
  Double_t nCandidatesHisto = hslice->Integral(hslice->GetXaxis()->FindBin(cutmodel.Eval(p)), hslice->GetXaxis()->FindBin(5));
  Double_t nPionsBefore = pi1->Integral(cutmodel.Eval(p), 5);
  Double_t contaminationHisto = nPionsBefore/nCandidatesHisto;
  printf("Contamination @ %.3fGeV/c: %f\n", p, contamination);
  printf("Contamination @ %.3fGeV/c: %f (determined from double gauss and histogram)\n", p, contaminationHisto);

  TPaveText *tcont = new TPaveText(0.45, 0.35, 0.9, 0.45, "NDC");
  tcont->SetBorderSize(0);
  tcont->SetFillStyle(0);
  tcont->AddText(Form("contamination: %f", contamination));
  tcont->Draw();

  // Save Canvas
  TFile *out = new TFile(Form("slice%dpr.root", ((Int_t)(100*p))), "RECREATE");
  out->cd();
  output->Write();
  out->Close(); delete out;
}
 DrawSliceProton.C:1
 DrawSliceProton.C:2
 DrawSliceProton.C:3
 DrawSliceProton.C:4
 DrawSliceProton.C:5
 DrawSliceProton.C:6
 DrawSliceProton.C:7
 DrawSliceProton.C:8
 DrawSliceProton.C:9
 DrawSliceProton.C:10
 DrawSliceProton.C:11
 DrawSliceProton.C:12
 DrawSliceProton.C:13
 DrawSliceProton.C:14
 DrawSliceProton.C:15
 DrawSliceProton.C:16
 DrawSliceProton.C:17
 DrawSliceProton.C:18
 DrawSliceProton.C:19
 DrawSliceProton.C:20
 DrawSliceProton.C:21
 DrawSliceProton.C:22
 DrawSliceProton.C:23
 DrawSliceProton.C:24
 DrawSliceProton.C:25
 DrawSliceProton.C:26
 DrawSliceProton.C:27
 DrawSliceProton.C:28
 DrawSliceProton.C:29
 DrawSliceProton.C:30
 DrawSliceProton.C:31
 DrawSliceProton.C:32
 DrawSliceProton.C:33
 DrawSliceProton.C:34
 DrawSliceProton.C:35
 DrawSliceProton.C:36
 DrawSliceProton.C:37
 DrawSliceProton.C:38
 DrawSliceProton.C:39
 DrawSliceProton.C:40
 DrawSliceProton.C:41
 DrawSliceProton.C:42
 DrawSliceProton.C:43
 DrawSliceProton.C:44
 DrawSliceProton.C:45
 DrawSliceProton.C:46
 DrawSliceProton.C:47
 DrawSliceProton.C:48
 DrawSliceProton.C:49
 DrawSliceProton.C:50
 DrawSliceProton.C:51
 DrawSliceProton.C:52
 DrawSliceProton.C:53
 DrawSliceProton.C:54
 DrawSliceProton.C:55
 DrawSliceProton.C:56
 DrawSliceProton.C:57
 DrawSliceProton.C:58
 DrawSliceProton.C:59
 DrawSliceProton.C:60
 DrawSliceProton.C:61
 DrawSliceProton.C:62
 DrawSliceProton.C:63
 DrawSliceProton.C:64
 DrawSliceProton.C:65
 DrawSliceProton.C:66
 DrawSliceProton.C:67
 DrawSliceProton.C:68
 DrawSliceProton.C:69
 DrawSliceProton.C:70
 DrawSliceProton.C:71
 DrawSliceProton.C:72
 DrawSliceProton.C:73
 DrawSliceProton.C:74
 DrawSliceProton.C:75
 DrawSliceProton.C:76
 DrawSliceProton.C:77
 DrawSliceProton.C:78
 DrawSliceProton.C:79
 DrawSliceProton.C:80
 DrawSliceProton.C:81
 DrawSliceProton.C:82
 DrawSliceProton.C:83
 DrawSliceProton.C:84
 DrawSliceProton.C:85
 DrawSliceProton.C:86
 DrawSliceProton.C:87
 DrawSliceProton.C:88
 DrawSliceProton.C:89
 DrawSliceProton.C:90
 DrawSliceProton.C:91
 DrawSliceProton.C:92
 DrawSliceProton.C:93
 DrawSliceProton.C:94
 DrawSliceProton.C:95
 DrawSliceProton.C:96
 DrawSliceProton.C:97
 DrawSliceProton.C:98
 DrawSliceProton.C:99
 DrawSliceProton.C:100
 DrawSliceProton.C:101
 DrawSliceProton.C:102
 DrawSliceProton.C:103
 DrawSliceProton.C:104
 DrawSliceProton.C:105
 DrawSliceProton.C:106
 DrawSliceProton.C:107
 DrawSliceProton.C:108
 DrawSliceProton.C:109
 DrawSliceProton.C:110
 DrawSliceProton.C:111
 DrawSliceProton.C:112
 DrawSliceProton.C:113
 DrawSliceProton.C:114
 DrawSliceProton.C:115
 DrawSliceProton.C:116
 DrawSliceProton.C:117
 DrawSliceProton.C:118
 DrawSliceProton.C:119
 DrawSliceProton.C:120
 DrawSliceProton.C:121
 DrawSliceProton.C:122
 DrawSliceProton.C:123
 DrawSliceProton.C:124
 DrawSliceProton.C:125
 DrawSliceProton.C:126
 DrawSliceProton.C:127
 DrawSliceProton.C:128
 DrawSliceProton.C:129
 DrawSliceProton.C:130
 DrawSliceProton.C:131
 DrawSliceProton.C:132
 DrawSliceProton.C:133
 DrawSliceProton.C:134
 DrawSliceProton.C:135
 DrawSliceProton.C:136
 DrawSliceProton.C:137
 DrawSliceProton.C:138
 DrawSliceProton.C:139
 DrawSliceProton.C:140
 DrawSliceProton.C:141
 DrawSliceProton.C:142
 DrawSliceProton.C:143
 DrawSliceProton.C:144
 DrawSliceProton.C:145
 DrawSliceProton.C:146
 DrawSliceProton.C:147
 DrawSliceProton.C:148
 DrawSliceProton.C:149
 DrawSliceProton.C:150
 DrawSliceProton.C:151
 DrawSliceProton.C:152
 DrawSliceProton.C:153
 DrawSliceProton.C:154