File indexing completed on 2026-05-19 07:37:07
0001
0002
0003
0004
0005
0006 #include "ROOT/RDataFrame.hxx"
0007 #include <iostream>
0008 #include <algorithm>
0009 #include <string>
0010
0011 #include "edm4hep/MCParticleCollection.h"
0012 #include "edm4hep/SimCalorimeterHitCollection.h"
0013
0014 #include "common_bench/benchmark.h"
0015
0016 #include <boost/range/combine.hpp>
0017
0018 #include "DD4hep/IDDescriptor.h"
0019 #include "DD4hep/Readout.h"
0020 #include "DD4hep/Segmentations.h"
0021
0022 #include "TCanvas.h"
0023 #include "TStyle.h"
0024 #include "TMath.h"
0025 #include "TH1.h"
0026 #include "TF1.h"
0027 #include "TH1D.h"
0028 #include "TFitResult.h"
0029 #include "TLegend.h"
0030 #include "TString.h"
0031 #include "TGraph.h"
0032 #include "TGraph2D.h"
0033 #include "TGraphErrors.h"
0034 #include "TLine.h"
0035 #include "TError.h"
0036
0037 R__LOAD_LIBRARY(libfmt.so)
0038 #include "fmt/core.h"
0039 #include "DD4hep/Detector.h"
0040 #include "DDG4/Geant4Data.h"
0041 #include "DDRec/CellIDPositionConverter.h"
0042 #include "emcal_barrel_common_functions.h"
0043 #include <nlohmann/json.hpp>
0044
0045 using ROOT::RDataFrame;
0046 using namespace ROOT::VecOps;
0047
0048 void emcal_barrel_pion_rejection_analysis(
0049 const char* input_fname1 = "sim_output/sim_emcal_barrel_piRej_electron.edm4hep.root",
0050 const char* input_fname2 = "sim_output/sim_emcal_barrel_piRej_piminus.edm4hep.root"
0051 )
0052 {
0053
0054 gErrorIgnoreLevel = kFatal;
0055
0056
0057 gROOT->SetStyle("Plain");
0058 gStyle->SetOptFit(1);
0059 gStyle->SetLineWidth(2);
0060 gStyle->SetPadTickX(1);
0061 gStyle->SetPadTickY(1);
0062 gStyle->SetPadGridX(1);
0063 gStyle->SetPadGridY(1);
0064 gStyle->SetPadLeftMargin(0.14);
0065 gStyle->SetPadRightMargin(0.14);
0066
0067 ROOT::EnableImplicitMT();
0068 ROOT::RDataFrame d0("events", {input_fname1, input_fname2});
0069
0070
0071 if (! d0.HasColumn("EcalBarrelScFiHits")) {
0072 std::cout << "EcalBarrelScFiHits is required" << std::endl;
0073 return;
0074 }
0075
0076
0077 std::string detector_path = "";
0078 std::string detector_name = "athena";
0079 if(std::getenv("DETECTOR_PATH")) {
0080 detector_path = std::getenv("DETECTOR_PATH");
0081 }
0082 if(std::getenv("DETECTOR_CONFIG")) {
0083 detector_name = std::getenv("DETECTOR_CONFIG");
0084 }
0085
0086
0087
0088
0089
0090
0091
0092
0093
0094
0095
0096
0097 int layerNum;
0098 int dep_min = 1;
0099 int dep_max = 6;
0100
0101
0102 dd4hep::Detector& detector = dd4hep::Detector::getInstance();
0103 detector.fromCompact(fmt::format("{}/{}.xml", detector_path, detector_name));
0104
0105 auto decoder = detector.readout("EcalBarrelImagingHits").idSpec().decoder();
0106 auto decoderScFi = detector.readout("EcalBarrelScFiHits").idSpec().decoder();
0107 auto layer_index = decoder->index("layer");
0108 auto layer_indexScFi = decoderScFi->index("layer");
0109
0110
0111 auto Ethr = [](std::vector<edm4hep::MCParticleData> const& input) {
0112 return TMath::Sqrt(input[2].momentum.x*input[2].momentum.x + input[2].momentum.y*input[2].momentum.y + input[2].momentum.z*input[2].momentum.z + input[2].mass*input[2].mass);
0113 };
0114
0115
0116 auto Pthr = [](std::vector<edm4hep::MCParticleData> const& input) {
0117 return TMath::Sqrt(input[2].momentum.x*input[2].momentum.x + input[2].momentum.y*input[2].momentum.y + input[2].momentum.z*input[2].momentum.z);
0118 };
0119
0120
0121 auto Eta = [](std::vector<edm4hep::MCParticleData> const& input) {
0122 double E = TMath::Sqrt(input[2].momentum.x*input[2].momentum.x + input[2].momentum.y*input[2].momentum.y + input[2].momentum.z*input[2].momentum.z + input[2].mass*input[2].mass);
0123 return 0.5*TMath::Log((E + input[2].momentum.z) / (E - input[2].momentum.z));
0124 };
0125
0126
0127 auto pT = [](std::vector<edm4hep::MCParticleData> const& input) {
0128 return TMath::Sqrt(input[2].momentum.x*input[2].momentum.x + input[2].momentum.y*input[2].momentum.y);
0129 };
0130
0131
0132 auto nhits = [] (const std::vector<edm4hep::SimCalorimeterHitData>& evt) {return (int) evt.size(); };
0133
0134
0135 auto Esim = [](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0136 double total_edep = 0.0;
0137 for (const auto& i: evt){
0138 total_edep += i.energy;
0139 }
0140 return total_edep;
0141 };
0142
0143
0144 auto Esim_dep2 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0145 auto total_edep = 0.0;
0146 for (const auto& i: evt) {
0147 if( decoder->get(i.cellID, layer_index) < 3 ){
0148 total_edep += i.energy;
0149 }
0150 }
0151 return total_edep;
0152 };
0153
0154
0155
0156 auto Esim_dep3 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0157 auto total_edep = 0.0;
0158 for (const auto& i: evt) {
0159 if( decoder->get(i.cellID, layer_index) < 4 ){
0160 total_edep += i.energy;
0161 }
0162 }
0163 return total_edep;
0164 };
0165
0166
0167 auto Esim_dep4 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0168 auto total_edep = 0.0;
0169 for (const auto& i: evt) {
0170 if( decoder->get(i.cellID, layer_index) < 5 ){
0171 total_edep += i.energy;
0172 }
0173 }
0174 return total_edep;
0175 };
0176
0177
0178 auto Esim_dep5 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0179 auto total_edep = 0.0;
0180 for (const auto& i: evt) {
0181 if( decoder->get(i.cellID, layer_index) < 6 ){
0182 total_edep += i.energy;
0183 }
0184 }
0185 return total_edep;
0186 };
0187
0188
0189 auto Esim_dep6 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0190 auto total_edep = 0.0;
0191 for (const auto& i: evt) {
0192 if( decoder->get(i.cellID, layer_index) < 7 ){
0193 total_edep += i.energy;
0194 }
0195 }
0196 return total_edep;
0197 };
0198
0199
0200 auto Esim_dep6_ScFi = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0201 auto total_edep = 0.0;
0202 for (const auto& i: evt) {
0203 if( decoderScFi->get(i.cellID, layer_indexScFi) < 7 ){
0204 total_edep += i.energy;
0205 }
0206 }
0207 return total_edep;
0208 };
0209
0210
0211 auto Esim_dep7 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0212 auto total_edep = 0.0;
0213 for (const auto& i: evt) {
0214 if( decoder->get(i.cellID, layer_index) < 8 ){
0215 total_edep += i.energy;
0216 }
0217 }
0218 return total_edep;
0219 };
0220
0221 auto Esim_dep8 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0222 auto total_edep = 0.0;
0223 for (const auto& i: evt) {
0224 if( decoder->get(i.cellID, layer_index) < 9 ){
0225 total_edep += i.energy;
0226 }
0227 }
0228 return total_edep;
0229 };
0230
0231
0232 auto Esim_dep9 = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0233 auto total_edep = 0.0;
0234 for (const auto& i: evt) {
0235 if( decoder->get(i.cellID, layer_index) < 10 ){
0236 total_edep += i.energy;
0237 }
0238 }
0239 return total_edep;
0240 };
0241
0242
0243
0244 auto Esim_dep = [=](const std::vector<edm4hep::SimCalorimeterHitData>& evt) {
0245 std::vector<double>res(20);
0246 for (const auto& i: evt) {
0247 res[decoder->get(i.cellID, layer_index)] += i.energy;
0248 }
0249 return res;
0250 };
0251
0252
0253 auto Esim_dep_sum = [&dep_min, &dep_max](const std::vector<double>& dep) {
0254 double res = 0;
0255 for (int i = dep_min; i < dep_max + 1; i++) {
0256 if (i >= dep.size()) continue;
0257 res += dep[i];
0258 }
0259 return res;
0260 };
0261
0262
0263 auto Esim_depN = [&layerNum](const std::vector<double>& dep) {
0264 return (layerNum < dep.size() ? dep[layerNum] : 0.0);
0265 };
0266
0267
0268 auto fsam = [](const double& sampled, const double& thrown) {
0269 return sampled / thrown;
0270 };
0271
0272
0273 auto fEp = [](const double& E_front, const double& mom) {
0274 return E_front / mom;
0275 };
0276
0277
0278 auto getpid = [](std::vector<edm4hep::MCParticleData> const& input) {
0279 return input[2].PDG;
0280 };
0281
0282
0283 auto getdau = [](std::vector<edm4hep::MCParticleData> const& input){
0284 return input[2].daughters_begin;
0285 };
0286
0287
0288 auto is_electron = [](std::vector<edm4hep::MCParticleData> const& input){
0289 return (input[2].PDG == 11 ? true : false);
0290 };
0291
0292
0293 auto is_piMinus = [](std::vector<edm4hep::MCParticleData> const& input){
0294 return (input[2].PDG == -211 ? true : false);
0295 };
0296
0297
0298 auto EDep_bool = [](std::vector<double> const& input){
0299 return (input[0] > input[1]);
0300 };
0301
0302 auto Diff = [](const double &esim, const double &edep){
0303 return esim-edep;
0304 };
0305
0306
0307
0308 auto d1 = d0.Define("Ethr", Ethr, {"MCParticles"})
0309 .Define("Pthr", Pthr, {"MCParticles"})
0310 .Define("nhits", nhits, {"EcalBarrelImagingHits"})
0311 .Define("Esim", Esim, {"EcalBarrelImagingHits"})
0312 .Define("EsimScFi", Esim, {"EcalBarrelScFiHits"})
0313 .Define("EsimOverP", fEp, {"Esim", "Pthr"})
0314 .Define("EsimScFiOverP", fEp, {"EsimScFi", "Pthr"})
0315 .Define("EsimTot", "EsimScFi+Esim")
0316 .Define("EsimTotOverP", fEp, {"EsimTot", "Pthr"})
0317 .Define("fsam", fsam, {"Esim","Ethr"})
0318 .Define("pid", getpid, {"MCParticles"})
0319 .Define("EDep", Esim_dep, {"EcalBarrelImagingHits"})
0320 .Define("EDepSum", Esim_dep_sum, {"EDep"})
0321 .Define("EDepN", Esim_depN, {"EDep"})
0322 .Define("EDep2", Esim_dep2, {"EcalBarrelImagingHits"})
0323 .Define("EDep3", Esim_dep3, {"EcalBarrelImagingHits"})
0324 .Define("EDep4", Esim_dep4, {"EcalBarrelImagingHits"})
0325 .Define("EDep5", Esim_dep5, {"EcalBarrelImagingHits"})
0326 .Define("EDep6", Esim_dep6, {"EcalBarrelImagingHits"})
0327 .Define("EDep6OverP", fEp, {"EDep6", "Pthr"})
0328 .Define("EOverP", fEp, {"EDep3", "Pthr"})
0329 .Define("Eta", Eta, {"MCParticles"})
0330 .Define("pT", pT, {"MCParticles"})
0331 .Define("EDepOverP", fEp, {"EDepN", "Pthr"})
0332 .Define("EDepOverPT", fEp, {"EDepN", "pT"})
0333 .Define("EDepSumOverP", fEp, {"EDepSum", "Pthr"})
0334 .Define("EDepSumOverPT", fEp, {"EDepSum", "pT"})
0335 .Define("EDepFrac", fEp, {"EDepSum", "Esim"})
0336 ;
0337
0338
0339 dep_min = 1;
0340 dep_max = 6;
0341 auto d_ele = d1.Filter(is_electron, {"MCParticles"});
0342 auto d_pim = d1.Filter(is_piMinus, {"MCParticles"});
0343
0344
0345 std::string currentCut = "(EDep6OverP>2.5e-3)&&(EDep6>5e-3)";
0346
0347
0348
0349 std::vector<std::string> var = {"Esim [GeV];", "EsimTot [GeV];", "EDep6 [GeV];", "EDep6/p;", "pT [GeV];", "#eta;", "EsimScFi [GeV]", "EsimScFi/p"};
0350 std::vector<std::string> var_save = {"Esim", "EsimTot", "EDep6", "EDep6OverP", "pT", "eta", "EsimScFi", "EsimScFiOverP"};
0351 std::vector<std::string> col = {"Esim", "EsimTot", "EDep6", "EDep6OverP", "pT", "Eta", "EsimScFi", "EsimScFiOverP"};
0352 std::vector<std::vector<double>> h1Ranges = {{0,0.2}, {0, 0.2}, {0,0.25}, {0, 0.02}, {0, 18}, {-1, 1}, {0,0.2}, {0,0.2}};
0353 for (int i = 0; i < var.size(); i++){
0354 std::string title = "#pi^{-}, e^{-};" + var[i] + " Events";
0355 auto he = d_ele.Histo1D({"he", title.c_str(), 100, h1Ranges[i][0], h1Ranges[i][1]}, col[i]);
0356 auto hp = d_pim.Histo1D({"hp", title.c_str(), 100, h1Ranges[i][0], h1Ranges[i][1]}, col[i]);
0357
0358 hp->GetYaxis()->SetTitleOffset(1.4);
0359 he->SetLineWidth(2);
0360 he->SetLineColor(kRed);
0361 hp->SetLineWidth(2);
0362 hp->SetLineColor(kBlue);
0363 auto c = new TCanvas("c", "c", 700, 500);
0364 auto leng = new TLegend(0.7, 0.7, 0.9, 0.9);
0365 if (var[i] != "EsimScFi/p"){
0366 hp->DrawClone();
0367 he->DrawClone("same");
0368 }
0369 else {
0370 he->DrawClone();
0371 hp->DrawClone("same");
0372 }
0373 c->Update();
0374
0375 leng->AddEntry(he.GetPtr(),"e^{-}","l");
0376 leng->AddEntry(hp.GetPtr(),"#pi^{-}","l");
0377 leng->Draw();
0378 c->SaveAs(("results/emcal_barrel_pion_rej_uncut_comb_" + var_save[i] + ".png").c_str());
0379 }
0380
0381
0382
0383
0384
0385 std::string cutEEta;
0386 std::vector<std::vector<double>> EBins = {{0,2}, {2, 4}, {4, 6}, {6, 9}, {9, 12}, {12, 19}};
0387 for (int i = 0; i < EBins.size(); i++){
0388 std::string minCut = "Pthr>="+std::to_string(EBins[i][0]);
0389 std::string maxCut = "Pthr<"+std::to_string(EBins[i][1]);
0390 cutEEta += "(" + minCut + "&&" + maxCut + "&&";
0391
0392 for (int j = -1; j < 1; j++){
0393 std::string title = "#pi^{-}, e^{-}";
0394 title += fmt::format(" : {} < E < {}", EBins[i][0], EBins[i][1]);
0395 title += fmt::format(" & {} < #eta < {};EDep6/p; Events", j, j+1);
0396 std::string etaCutMin = fmt::format("Eta>={}", j);
0397 std::string etaCutMax = fmt::format("Eta<{}",j+1);
0398 cutEEta += "(" + etaCutMin + "&&" + etaCutMax;
0399 auto he = d_ele.Filter(minCut).Filter(maxCut).Filter(etaCutMin).Filter(etaCutMax).Histo1D({"he", title.c_str(), 50, 0, 0.02}, "EDep6OverP");
0400 auto hp = d_pim.Filter(minCut).Filter(maxCut).Filter(etaCutMin).Filter(etaCutMax).Histo1D({"hp", title.c_str(), 50, 0, 0.02}, "EDep6OverP");
0401 auto hecopy = he->DrawCopy();
0402 hecopy->Fit("gaus", "", "", 0, hecopy->GetMaximum());
0403 TF1* gaus = hecopy->GetFunction("gaus");
0404 if (gaus != nullptr) {
0405 double* res = gaus->GetParameters();
0406 cutEEta += fmt::format("&&EDep6OverP>={}", res[1] - 2.0*res[2]);
0407 cutEEta += fmt::format("&&EDep6OverP<{})||",res[1] + 3.0*res[2]);
0408 } else {
0409 cutEEta += ")||";
0410 std::cerr << "Warning: Gaussian fit failed for E bin " << i << ", eta bin " << j << std::endl;
0411 }
0412
0413 hp->GetYaxis()->SetTitleOffset(1.4);
0414 he->SetLineWidth(2);
0415 he->SetLineColor(kRed);
0416 hp->SetLineWidth(2);
0417 hp->SetLineColor(kBlue);
0418 auto c = new TCanvas("c", "c", 700, 500);
0419 auto leng = new TLegend(0.7, 0.7, 0.9, 0.9);
0420 hp->DrawClone();
0421 he->DrawClone("same");
0422 c->Update();
0423
0424 leng->AddEntry(he.GetPtr(),"e^{-}","l");
0425 leng->AddEntry(hp.GetPtr(),"#pi^{-}","l");
0426 leng->Draw();
0427 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_uncut_comb_E{}Eta{}.png", i, j+1)).c_str());
0428 }
0429 cutEEta.pop_back();
0430 cutEEta.pop_back();
0431 cutEEta += ")||";
0432 }
0433 cutEEta.pop_back();
0434 cutEEta.pop_back();
0435 currentCut = cutEEta;
0436
0437
0438 auto d_pim_cut = d_pim.Filter(currentCut.c_str());
0439 auto d_ele_cut = d_ele.Filter(currentCut.c_str());
0440
0441
0442 std::vector<double> E = {5, 10, 18};
0443 std::vector<double> ledges5 = {2.8, 0.4, 0.3, 0.5};
0444 std::vector<double> ledges10 = {1.4, 0.5, 0.6, 1.0};
0445 std::vector<double> ledges18 = {0.9, 0.9, 1.0, 1.8};
0446 std::vector<double> maxRate5 = {0.1, 100, 500, 1000};
0447 std::vector<double> maxRate10 = {10, 400, 800, 1000};
0448 std::vector<double> maxRate18 = {200, 800, 1000, 100};
0449 std::vector<vector<double>> maxRate = {maxRate5, maxRate10, maxRate18};
0450 std::vector<vector<double>> lowEdges = {ledges5, ledges10, ledges18};
0451 double suppression = 1e-4;
0452 std::vector<vector<double>> rejRatios = lowEdges;
0453 std::vector<vector<double>> effEle = lowEdges;
0454 std::vector<vector<double>> effPim = lowEdges;
0455 std::vector<std::string> etaBin = {"Eta >= -3.5 && Eta < -2.0", "Eta >= -2.0 && Eta < -1.0", "Eta >= -1.0 && Eta < 0", "Eta >= 0 && Eta < 1.0"};
0456 std::vector<std::string> etaTitle = {"-3.5 < #eta < -2.0", "-2.0 < #eta < -1.0", "-1.0 < #eta < 0", "0 < #eta < 1.0"};
0457
0458
0459 std::vector<double> pBins = {0.1,0.2,0.3,0.4,0.5,1,2,3,4,5,10,12,14,16,18};
0460 EBins = {{0,2}, {2, 4}, {4, 6}, {6, 9}, {9, 12}, {12, 19}};
0461 auto tg = new TGraphErrors();
0462
0463 for (int i = 0; i < EBins.size(); i++){
0464 std::string filter = currentCut;
0465 filter += "&&(Pthr>=" + std::to_string(EBins[i][0]) + "&&Pthr<" + std::to_string(EBins[i][1]) + ")";
0466 double numer = (double)*d_ele_cut.Filter(filter.c_str()).Count();
0467 double denom = (double)*d_pim_cut.Filter(filter.c_str()).Count();
0468 double error = std::sqrt(std::pow(numer / denom, 2.0)*(1.0/numer + 1.0/denom));
0469 double ratio = numer / denom;
0470 if (denom == 0){ratio = 1; error = 1;}
0471 tg->SetPoint(i, 0.5*(EBins[i][0] + EBins[i][1]), ratio);
0472 tg->SetPointError(i, 0.5*(EBins[i][1] - EBins[i][0]), error);
0473 }
0474 double e_eff = (double)*d_ele_cut.Count() / (double)*d_ele.Count();
0475 tg->SetTitle(("#pi Rejection with #varepsilon_{e} = "+ std::to_string(e_eff)).c_str());
0476 tg->GetXaxis()->SetTitle("p [GeV]");
0477 tg->GetYaxis()->SetTitle("R_{e/#pi}");
0478 tg->SetMarkerColor(kBlue);
0479 tg->SetMarkerStyle(20);
0480
0481 auto cp = new TCanvas("cp", "cp");
0482 cp->SetLogy();
0483 cp->SetLogx();
0484 tg->DrawClone("ap");
0485 cp->SaveAs("results/emcal_barrel_pion_rej_RatioRej.png");
0486
0487
0488
0489
0490 dep_min = 1;
0491 dep_max = 6;
0492 for (int i = 0; i < 3; i++){
0493 for (int j = 2; j < 4; j++){
0494
0495
0496 std::string pCut = "Pthr>=" + std::to_string(lowEdges[i][j]) + "&&Pthr<" + std::to_string(E[i]);
0497 auto e_eta = d_ele.Filter(etaBin[j]).Filter(pCut);
0498 auto p_eta = d_pim.Filter(etaBin[j]).Filter(pCut);
0499
0500
0501 std::string title = "e^{-} (E = " + std::to_string((int)E[i]) + " GeV) : " + etaTitle[j] + "; p [GeV]; Events";
0502 auto he = e_eta.Histo1D({"he", title.c_str(), 100, lowEdges[i][j], E[i]}, "Pthr");
0503 he->SetLineColor(kBlue);
0504 auto he_cut = e_eta.Filter(currentCut).Histo1D({"he_cut", title.c_str(), 100, lowEdges[i][j], E[i]}, "Pthr");
0505 he_cut->GetYaxis()->SetTitleOffset(1.4);
0506 he_cut->SetLineWidth(2);
0507 he_cut->SetLineColor(kRed);
0508
0509 title = "#pi^{-} (E = " + std::to_string((int)E[i]) + " GeV) : " + etaTitle[j] + "; p [GeV]; Events";
0510 auto hp = p_eta.Histo1D({"hp", title.c_str(), 100, lowEdges[i][j], E[i]}, "Pthr");
0511 hp->SetLineColor(kBlue);
0512 auto hp_cut = p_eta.Filter(currentCut).Histo1D({"hp", title.c_str(), 100, lowEdges[i][j], E[i]}, "Pthr");
0513 hp_cut->GetYaxis()->SetTitleOffset(1.4);
0514 hp_cut->SetLineWidth(2);
0515 hp_cut->SetLineColor(kRed);
0516
0517 auto c = new TCanvas("c", "c", 700, 500);
0518 c->SetLogy(1);
0519 he->DrawClone();
0520 he_cut->DrawClone("same");
0521 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_mom_ele_E{}_eta{}.png", (int)E[i], j)).c_str());
0522 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_mom_ele_E{}_eta{}.pdf", (int)E[i], j)).c_str());
0523 c->Clear();
0524
0525 hp->DrawClone();
0526 hp_cut->DrawClone("same");
0527 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_mom_pim_E{}_eta{}.png", (int)E[i], j)).c_str());
0528 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_mom_pim_E{}_eta{}.pdf", (int)E[i], j)).c_str());
0529 c->Clear();
0530
0531
0532
0533 rejRatios[i][j] = (double)hp_cut->Integral() / (double)he_cut->Integral();
0534 effPim[i][j] = (double)hp_cut->Integral() / (double)hp->Integral();
0535 effEle[i][j] = (double)he_cut->Integral() / (double)he->Integral();
0536
0537 hp_cut->Divide(he.GetPtr());
0538 title = "#pi^{-}/e^{-} (E = " + std::to_string((int)E[i]) + " GeV) : " + etaTitle[j];
0539 hp_cut->SetTitle(title.c_str());
0540 hp_cut->SetLineColor(kBlack);
0541 hp_cut->DrawClone();
0542 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_ratio_pim_E{}_eta{}.png", (int)E[i], j)).c_str());
0543 c->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_ratio_pim_E{}_eta{}.pdf", (int)E[i], j)).c_str());
0544 c->Clear();
0545
0546
0547 std::vector<std::string> endStr = {";pT [GeV]; Events", ";EDep6/p; Events"};
0548 std::vector<std::string> var_save_loc = {"pT", "EDep6OverP"};
0549 std::vector<std::string> col_loc = {"pT", "EDep6OverP"};
0550 std::vector<std::vector<double>> h1Ranges = {{0, E[i]}, {0, 0.02}};
0551 for (int k = 0; k < 2; k++){
0552 auto cl = new TCanvas("cl", "cl", 700, 500);
0553 title = "e^{-} (E = " + std::to_string((int)E[i]) + " GeV) : " + etaTitle[j] + endStr[k];
0554 auto he1 = e_eta.Histo1D({"he", title.c_str(), 100, h1Ranges[k][0], h1Ranges[k][1]}, col_loc[k]);
0555 he1->SetLineColor(kBlue);
0556 auto he1_cut = e_eta.Filter(currentCut).Histo1D({"he_cut", title.c_str(), 100, h1Ranges[k][0], h1Ranges[k][1]}, col_loc[k]);
0557 he1_cut->GetYaxis()->SetTitleOffset(1.4);
0558 he1_cut->SetLineWidth(2);
0559 he1_cut->SetLineColor(kRed);
0560
0561 title = "#pi^{-} (E = " + std::to_string((int)E[i]) + " GeV) : " + etaTitle[j] + endStr[k];
0562 auto hp1 = p_eta.Histo1D({"hp", title.c_str(), 100, h1Ranges[k][0], h1Ranges[k][1]}, col_loc[k]);
0563 hp1->SetLineColor(kBlue);
0564 auto hp1_cut = p_eta.Filter(currentCut).Histo1D({"hp_cut", title.c_str(), 100, h1Ranges[k][0], h1Ranges[k][1]}, col_loc[k]);
0565 hp1_cut->GetYaxis()->SetTitleOffset(1.4);
0566 hp1_cut->SetLineWidth(2);
0567 hp1_cut->SetLineColor(kRed);
0568
0569 cl->SetLogy(1);
0570 he1->DrawClone();
0571 he1_cut->DrawClone("same");
0572 cl->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_{}_ele_E{}_eta{}.png", col_loc[k], (int)E[i], j)).c_str());
0573 cl->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_{}_ele_E{}_eta{}.pdf", col_loc[k], (int)E[i], j)).c_str());
0574 cl->Clear();
0575
0576 hp1->DrawClone();
0577 hp1_cut->DrawClone("same");
0578 cl->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_{}_pim_E{}_eta{}.png", col_loc[k], (int)E[i], j)).c_str());
0579 cl->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_{}_pim_E{}_eta{}.pdf", col_loc[k], (int)E[i], j)).c_str());
0580 cl->Clear();
0581
0582
0583 title = "#pi^{-}, e^{-} (E = " + std::to_string((int)E[i]) + " GeV) : " + etaTitle[j] + endStr[k];
0584 hp1_cut->SetLineColor(kBlue);
0585 he1_cut->SetLineColor(kRed);
0586 he1_cut->SetTitle(title.c_str());
0587
0588 auto leng = new TLegend(0.7, 0.7, 0.9, 0.9);
0589 he1_cut->DrawClone();
0590 hp1_cut->DrawClone("same");
0591 cl->Update();
0592
0593 leng->AddEntry(he1_cut.GetPtr(),"e^{-}","l");
0594 leng->AddEntry(hp1_cut.GetPtr(),"#pi^{-}","l");
0595 leng->Draw();
0596 cl->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_{}_comb_E{}_eta{}.png", col_loc[k], (int)E[i], j)).c_str());
0597 cl->SaveAs((fmt::format("results/emcal_barrel_pion_rej_cut_{}_comb_E{}_eta{}.pdf", col_loc[k], (int)E[i], j)).c_str());
0598
0599 }
0600 }
0601 }
0602
0603
0604
0605 std::string test_tag = "Barrel_emcal_pion_rejection";
0606
0607 std:string detectorEle = "Barrel_emcal";
0608
0609 for (int i = 0; i < etaTitle.size(); i++){
0610 etaTitle[i].erase(std::remove(etaTitle[i].begin(), etaTitle[i].end(), '#'), etaTitle[i].end());
0611 std::replace(etaTitle[i].begin(), etaTitle[i].end(), 'e', 'E');
0612 }
0613
0614
0615 common_bench::Test pion_rejection_E18_Eta2{
0616 {{"name", fmt::format("{}_E{}_EtaBin{}", test_tag, (int)E[0], 2)},
0617 {"title", "Pion Rejection1"},
0618 {"description", fmt::format("Pion rejection with E = {}, and {}", (int)E[0], etaTitle[2])},
0619 {"quantity", "pi-/e-"},
0620 {"cut", currentCut},
0621 {"e- efficiency", std::to_string(effEle[0][2])},
0622 {"pi- efficiency", std::to_string(effPim[0][2])},
0623 {"target", std::to_string(suppression * maxRate[0][2])}
0624 }
0625 };
0626 suppression * maxRate[0][2] >= rejRatios[0][2] ? pion_rejection_E18_Eta2.pass(rejRatios[0][2]) : pion_rejection_E18_Eta2.fail(rejRatios[0][2]);
0627
0628
0629 common_bench::Test pion_rejection_E18_Eta3{
0630 {{"name", fmt::format("{}_E{}_EtaBin{}", test_tag, (int)E[0], 3)},
0631 {"title", "Pion Rejection"},
0632 {"description", fmt::format("Pion rejection with E = {}, and {}", (int)E[0], etaTitle[3])},
0633 {"quantity", "pi-/e-"},
0634 {"cut", currentCut},
0635 {"e- efficiency", std::to_string(effEle[0][3])},
0636 {"pi- efficiency", std::to_string(effPim[0][3])},
0637 {"target", std::to_string(suppression * maxRate[0][3])}
0638 }
0639 };
0640 suppression * maxRate[0][3] >= rejRatios[0][3] ? pion_rejection_E18_Eta3.pass(rejRatios[0][3]) : pion_rejection_E18_Eta3.fail(rejRatios[0][3]);
0641
0642
0643 common_bench::Test pion_rejection_E10_Eta2{
0644 {{"name", fmt::format("{}_E{}_EtaBin{}", test_tag, (int)E[1], 2)},
0645 {"title", "Pion Rejection"},
0646 {"description", fmt::format("Pion rejection with E = {}, and {}", (int)E[1], etaTitle[2])},
0647 {"quantity", "pi-/e-"},
0648 {"cut", currentCut},
0649 {"e- efficiency", std::to_string(effEle[1][2])},
0650 {"pi- efficiency", std::to_string(effPim[1][2])},
0651 {"target", std::to_string(suppression * maxRate[1][2])}
0652 }
0653 };
0654 suppression * maxRate[1][2] >= rejRatios[1][2] ? pion_rejection_E10_Eta2.pass(rejRatios[1][2]) : pion_rejection_E10_Eta2.fail(rejRatios[1][2]);
0655
0656
0657 common_bench::Test pion_rejection_E10_Eta3{
0658 {{"name", fmt::format("{}_E{}_EtaBin{}", test_tag, (int)E[1], 3)},
0659 {"title", "Pion Rejection"},
0660 {"description", fmt::format("Pion rejection with E = {}, and {}", (int)E[1], etaTitle[3])},
0661 {"quantity", "pi-/e-"},
0662 {"e- efficiency", std::to_string(effEle[1][3])},
0663 {"pi- efficiency", std::to_string(effPim[1][3])},
0664 {"target", std::to_string(suppression * maxRate[1][3])}
0665 }
0666 };
0667 suppression * maxRate[1][3] >= rejRatios[1][3] ? pion_rejection_E10_Eta3.pass(rejRatios[1][3]) : pion_rejection_E10_Eta3.fail(rejRatios[1][3]);
0668
0669
0670 common_bench::Test pion_rejection_E5_Eta2{
0671 {{"name", fmt::format("{}_E{}_EtaBin{}", test_tag, (int)E[2], 2)},
0672 {"title", "Pion Rejection"},
0673 {"description", fmt::format("Pion rejection with E = {}, and {}", (int)E[2], etaTitle[2])},
0674 {"quantity", "pi-/e-"},
0675 {"cut", currentCut},
0676 {"e- efficiency", std::to_string(effEle[2][2])},
0677 {"pi- efficiency", std::to_string(effPim[2][2])},
0678 {"target", std::to_string(suppression * maxRate[2][2])}
0679 }
0680 };
0681 suppression * maxRate[2][2] >= rejRatios[2][2] ? pion_rejection_E5_Eta2.pass(rejRatios[2][2]) : pion_rejection_E5_Eta2.fail(rejRatios[2][2]);
0682
0683
0684 common_bench::Test pion_rejection_E5_Eta3{
0685 {{"name", fmt::format("{}_E{}_EtaBin{}", test_tag, (int)E[2], 3)},
0686 {"title", "Pion Rejection"},
0687 {"description", fmt::format("Pion rejection with E = {}, and {}", (int)E[2], etaTitle[3])},
0688 {"quantity", "pi-/e-"},
0689 {"cut", currentCut},
0690 {"e- efficiency", std::to_string(effEle[2][3])},
0691 {"pi- efficiency", std::to_string(effPim[2][3])},
0692 {"target", std::to_string(suppression * maxRate[2][3])}
0693 }
0694 };
0695 suppression * maxRate[2][3] >= rejRatios[2][3] ? pion_rejection_E5_Eta3.pass(rejRatios[2][3]) : pion_rejection_E5_Eta3.fail(rejRatios[2][3]);
0696
0697
0698 common_bench::write_test({pion_rejection_E18_Eta2,
0699 pion_rejection_E18_Eta3,
0700 pion_rejection_E10_Eta2,
0701 pion_rejection_E10_Eta3,
0702 pion_rejection_E5_Eta2,
0703 pion_rejection_E5_Eta3},
0704 fmt::format("results/{}_pion_rej.json", detectorEle));
0705 }