Skip to content
Snippets Groups Projects

Resolve "Add layer energy-deposit profile to energy resolution benchmark"

Compare and
4 files
+ 144
40
Compare changes
  • Side-by-side
  • Inline
Files
4
@@ -18,6 +18,10 @@
#include "TH1D.h"
#include "TGraphErrors.h"
#include "DD4hep/Detector.h"
#include "DDG4/Geant4Data.h"
#include "DDRec/CellIDPositionConverter.h"
using ROOT::RDataFrame;
using namespace ROOT::VecOps;
@@ -25,7 +29,6 @@ using namespace ROOT::VecOps;
void save_canvas(TCanvas* c, std::string label)
{
c->SaveAs(fmt::format("results/energy_scan/{}.png",label).c_str());
c->SaveAs(fmt::format("results/energy_scan/{}.pdf",label).c_str());
}
void save_canvas(TCanvas* c, std::string label, double E)
@@ -44,7 +47,7 @@ void save_canvas(TCanvas* c, std::string var_label, std::string E_label, std::st
save_canvas(c, label_with_E);
}
std::tuple <double, double, double, double> extract_sampling_fraction_parameters(std::string particle_label, std::string E_label)
std::tuple <double, double, double, double, std::vector<double>, std::vector<double>> extract_sampling_fraction_parameters(std::string particle_label, std::string E_label, dd4hep::Detector& detector)
{
std::string input_fname = fmt::format("sim_output/energy_scan/{}/sim_emcal_barrel_{}.root", E_label, particle_label);
ROOT::EnableImplicitMT();
@@ -73,6 +76,25 @@ std::tuple <double, double, double, double> extract_sampling_fraction_parameters
return sampled / thrown;
};
// Energy deposited in layers
//dd4hep::rec::CellIDPositionConverter cellid_converter(detector);
auto decoder = detector.readout("EcalBarrelHits").idSpec().decoder();
fmt::print("{}\n", decoder->fieldDescription());
auto layer_index = decoder->index("layer");
fmt::print(" layer index is {}.\n", layer_index);
auto Layer = [=](const std::vector<dd4pod::CalorimeterHitData>& evt) {
std::vector<double> vec_layer_ids;
for (const auto& i: evt) {
auto layer_id = decoder->get(i.cellID, layer_index);
vec_layer_ids.push_back(layer_id);
}
RVec<double> layer_ids(vec_layer_ids.data(), vec_layer_ids.size());
return layer_ids;
};
// Define variables
auto d1 = d0.Define("Ethr", Ethr, {"mcparticles"})
.Define("nhits", nhits, {"EcalBarrelHits"})
@@ -91,8 +113,61 @@ std::tuple <double, double, double, double> extract_sampling_fraction_parameters
{"hEsim", "Energy Deposit; Energy Deposit [GeV]; Events", 500, 0.0, 0.5},
"Esim");
auto hfsam = d1.Histo1D(
{"hfsam", "Sampling Fraction; Sampling Fraction; Events", 200, 0.0, 0.05},
{"hfsam", "Sampling Fraction; Sampling Fraction; Events", 800, 0.0, 0.2},
"fsam");
TH1D* hno_edep = new TH1D("hno_edep", "% of events with no dE", 20, 1, 21);
std::vector<double> mean_layer;
std::vector<double> rms_layer;
// Draw energy per layer
for(int layer=1; layer<21; layer++) {
auto Esim_layer = [=](const std::vector<dd4pod::CalorimeterHitData>& evt) {
auto layer_edep = 0.0;
for (const auto& i: evt) {
if (decoder->get(i.cellID, layer_index) == layer) {
layer_edep += i.energyDeposit;
}
}
//std::cout << "layer: " << layer << " layer_edep: " << layer_edep << std::endl;
return layer_edep;
};
auto d2 = d0.Define(fmt::format("Esim_layer_{}",layer).c_str(), Esim_layer, {"EcalBarrelHits"});
std::cout << "Layer to process: " << layer << std::endl;
auto hEsim_layer = d2.Histo1D({fmt::format("hEsim_layer_{}",layer).c_str(),fmt::format("Energy Deposit in layer {}; Energy Deposit [Gev]; Events",layer).c_str(), 500, 0.0, 0.01}, fmt::format("Esim_layer_{}",layer));
TCanvas* clayer = new TCanvas(fmt::format("c1_layer{}",layer).c_str(), fmt::format("c1_layer{}",layer).c_str(), 700, 500);
clayer->cd();
clayer->SetLogy(1);
auto h = hEsim_layer->DrawCopy();
double no_edep = (h->GetBinContent(1)/h->GetEntries())*100;
std::cout<<h->GetBinContent(1)<<" "<<h->GetEntries()<<std::endl;
double up_fit = h->GetMean() + 3*h->GetStdDev();
double down_fit = h->GetMean() - 3*h->GetStdDev();
h->SetLineWidth(2);
h->SetLineColor(kBlue);
h->GetXaxis()->SetRangeUser(0.,up_fit);
h->GetXaxis()->SetRangeUser(h->GetXaxis()->GetBinUpEdge(1), h->GetXaxis()->GetBinUpEdge(h->GetXaxis()->GetNbins()));
mean_layer.push_back(h->GetMean());
rms_layer.push_back(h->GetStdDev());
h->SetRangeUser();
save_canvas(clayer, fmt::format("Esim_layer{}",layer), E_label, particle_label);
hno_edep->Fill(layer,no_edep);
}
TCanvas* clayer = new TCanvas("c_noedep", "c_noedep", 700, 500);
hno_edep->GetXaxis()->SetTitle("layer");
hno_edep->GetYaxis()->SetTitle("% of events with no dE");
hno_edep->SetLineWidth(2);
hno_edep->SetLineColor(kBlue);
hno_edep->Draw("hist");
save_canvas(clayer, "Layer_nodep", E_label, particle_label);
// Event Counts
auto nevents_thrown = d1.Count();
@@ -151,7 +226,7 @@ std::tuple <double, double, double, double> extract_sampling_fraction_parameters
double mean_err = gaus->GetParError(1);
double sigma_err = gaus->GetParError(2);
save_canvas(c4, "fsam", E_label, particle_label);
return std::make_tuple(mean, sigma, mean_err, sigma_err);
return std::make_tuple(mean, sigma, mean_err, sigma_err, mean_layer, rms_layer);
}
}
@@ -190,38 +265,64 @@ void emcal_barrel_energy_scan_analysis(std::string particle_label = "electron")
gStyle->SetPadLeftMargin(0.14);
gStyle->SetPadRightMargin(0.14);
auto scanned_energies = read_scanned_energies(fmt::format("sim_output/emcal_barrel_energy_scan_points_{}.txt", particle_label));
auto scanned_energies = read_scanned_energies(fmt::format("sim_output/emcal_barrel_energy_scan_points_{}.txt", particle_label));
//Take detector layers
std::string detector_path = "";
std::string detector_name = "athena";
if(std::getenv("DETECTOR_PATH")) {
detector_path = std::getenv("DETECTOR_PATH");
}
if(std::getenv("JUGGLER_DETECTOR")) {
detector_name = std::getenv("JUGGLER_DETECTOR");
}
dd4hep::Detector& detector = dd4hep::Detector::getInstance();
detector.fromCompact(fmt::format("{}/{}.xml", detector_path, detector_name));
TGraphErrors gr_fsam(scanned_energies.size()-1);
TGraphErrors gr_fsam_res(scanned_energies.size()-1);
for (const auto& E_label : scanned_energies) {
auto [fsam, fsam_res, fsam_err, fsam_res_err] = extract_sampling_fraction_parameters(particle_label, E_label);
auto E = std::stod(E_label);
gr_fsam.SetPoint(gr_fsam.GetN(),E,100*fsam);
gr_fsam.SetPointError(gr_fsam.GetN()-1,0., 100*fsam_err);
gr_fsam_res.SetPoint(gr_fsam_res.GetN(),E,100.0*(fsam_res/fsam));
auto fsam_res_rel_err = 100.0*(sqrt(pow((fsam_res_err/fsam),2)+pow((fsam_err*fsam_res)/(fsam*fsam),2)));
gr_fsam_res.SetPointError(gr_fsam_res.GetN()-1,0.,fsam_res_rel_err);
}
TGraphErrors gr_fsam(scanned_energies.size()-1);
TGraphErrors gr_fsam_res(scanned_energies.size()-1);
TGraphErrors gr_edep_layer1(scanned_energies.size()-1);
for (const auto& E_label : scanned_energies) {
auto [fsam, fsam_res, fsam_err, fsam_res_err, mean_layer,rms_layer] = extract_sampling_fraction_parameters(particle_label, E_label, detector);
auto E = std::stod(E_label);
gr_fsam.SetPoint(gr_fsam.GetN(),E,100*fsam);
gr_fsam.SetPointError(gr_fsam.GetN()-1,0., 100*fsam_err);
gr_fsam_res.SetPoint(gr_fsam_res.GetN(),E,100.0*(fsam_res/fsam));
auto fsam_res_rel_err = 100.0*(sqrt(pow((fsam_res_err/fsam),2)+pow((fsam_err*fsam_res)/(fsam*fsam),2)));
gr_fsam_res.SetPointError(gr_fsam_res.GetN()-1,0.,fsam_res_rel_err);
gr_edep_layer1.SetPoint(gr_edep_layer1.GetN(), E, mean_layer[0]);
gr_edep_layer1.SetPointError(gr_edep_layer1.GetN()-1, 0., rms_layer[0]);
std::cout<<"mean layer 0: "<<mean_layer[0]<<" 1: "<<mean_layer[1]<<std::endl;
std::cout<<"rms layer 0: "<<rms_layer[0]<<" 1: "<<rms_layer[1]<<std::endl;
}
TCanvas* c5 = new TCanvas("c5", "c5", 700, 500);
c5->cd();
gr_fsam.SetTitle("Sampling Fraction Scan;True Energy [GeV];Sampling Fraction [%]");
gr_fsam.SetMarkerStyle(20);
gr_fsam.Fit("pol0", "", "", 2., 20.);
gr_fsam.Draw("APE");
save_canvas(c5, fmt::format("emcal_barrel_{}_fsam_scan", particle_label));
TCanvas* c6 = new TCanvas("c6", "c6", 700, 500);
c6->cd();
TF1* func_res = new TF1("func_res", "[0]/sqrt(x) + [1]", 0.25, 20.);
func_res->SetLineWidth(2);
func_res->SetLineColor(kRed);
gr_fsam_res.SetTitle("Energy Resolution;True Energy [GeV];#Delta E/E [%]");
gr_fsam_res.SetMarkerStyle(20);
gr_fsam_res.Fit(func_res,"R");
gr_fsam_res.Draw("APE");
save_canvas(c6,fmt::format("emcal_barrel_{}_fsam_scan_res", particle_label));
TCanvas* c5 = new TCanvas("c5", "c5", 700, 500);
c5->cd();
gr_fsam.SetTitle("Sampling Fraction Scan;True Energy [GeV];Sampling Fraction [%]");
gr_fsam.SetMarkerStyle(20);
gr_fsam.Fit("pol0", "", "", 2., 20.);
gr_fsam.Draw("APE");
save_canvas(c5, fmt::format("emcal_barrel_{}_fsam_scan", particle_label));
TCanvas* c6 = new TCanvas("c6", "c6", 700, 500);
c6->cd();
TF1* func_res = new TF1("func_res", "[0]/sqrt(x) + [1]", 0.25, 20.);
func_res->SetLineWidth(2);
func_res->SetLineColor(kRed);
gr_fsam_res.SetTitle("Energy Resolution;True Energy [GeV];#Delta E/E [%]");
gr_fsam_res.SetMarkerStyle(20);
gr_fsam_res.Fit(func_res,"R");
gr_fsam_res.Draw("APE");
save_canvas(c6,fmt::format("emcal_barrel_{}_fsam_scan_res", particle_label));
TCanvas* c7 = new TCanvas("c7", "c7", 700, 500);
c7->cd();
gr_edep_layer1.SetTitle("Edep 1st layer;True Energy [GeV];Edep Mean");
gr_edep_layer1.SetMarkerStyle(20);
gr_edep_layer1.Draw("APE");
save_canvas(c7, fmt::format("emcal_barrel_{}_edep_layer1_scan", particle_label));
}
Loading