Skip to content
Snippets Groups Projects

Resolve "pi0 resolution for ECal barrel"

Merged Marshall Scott requested to merge 19-pi0-resolution-for-ecal-barrel into master
Compare and Show latest version
1 file
+ 11
11
Compare changes
  • Side-by-side
  • Inline
@@ -121,17 +121,8 @@ void emcal_barrel_pions_analysis(const char* input_fname = "sim_output/sim_emcal
.Define("nhits", nhits, {"EcalBarrelHits"})
.Define("Esim", Esim, {"EcalBarrelHits"})
.Define("fsam", fsam, {"Esim","Ethr"})
.Define("dE", eResol, {"Esim","Ethr"})
.Define("dE_rel", eResol_rel, {"Esim","Ethr"})
.Define("pid", getpid, {"mcparticles"})
;
// Set sampling Fraction, ideally this will be taken from a json file
auto svec = d1.Mean("fsam");
for (auto &ele : svec){cout << ele << endl;}
auto d2 = d1.Define("dE", eResol, {"Esim","Ethr"})
.Define("dE_rel", eResol_rel, {"Esim","Ethr"})
;
// Define Histograms
auto hEthr = d1.Histo1D({"hEthr", "Thrown Energy; Thrown Energy [GeV]; Events", 100, 0.0, 7.5}, "Ethr");
@@ -140,6 +131,15 @@ void emcal_barrel_pions_analysis(const char* input_fname = "sim_output/sim_emcal
auto hfsam = d1.Histo1D({"hfsam", "Sampling Fraction; Sampling Fraction; Events", 100, 0.0, 0.1}, "fsam");
auto hpid = d1.Histo1D({"hpid", "PID; PID; Count", 100, -220, 220}, "pid");
// Set sampling Fraction, ideally this will be taken from a json file
samp_frac = hfsam -> GetMean();
auto d2 = d1.Define("dE", eResol, {"Esim","Ethr"})
.Define("dE_rel", eResol_rel, {"Esim","Ethr"})
;
// Event Counts
auto nevents_thrown = d1.Count();
std::cout << "Number of Thrown Events: " << (*nevents_thrown) << "\n";
@@ -194,8 +194,8 @@ void emcal_barrel_pions_analysis(const char* input_fname = "sim_output/sim_emcal
c5->SaveAs("results/emcal_barrel_pions_pid.pdf");
//Energy Resolution Work
auto hdE = d1.Histo1D({"hdE", "dE; dE[GeV]; Events", 20, -7.5, 7.5}, "dE");//changed from 100
auto hdE_rel = d1.Histo1D({"hdE_rel", "dE Relative; dE Relative; Events", 20, -2, 2}, "dE_rel");//changed from 100
auto hdE = d2.Histo1D({"hdE", "dE; dE[GeV]; Events", 20, -7.5, 7.5}, "dE");//changed from 100
auto hdE_rel = d2.Histo1D({"hdE_rel", "dE Relative; dE Relative; Events", 20, -2, 2}, "dE_rel");//changed from 100
hdE_rel->Fit("gaus", "", "", -2, 2);
double* res = hdE_rel->GetFunction("gaus")->GetParameters();
Loading