Skip to content
Snippets Groups Projects
Commit 53d458b1 authored by Finn King's avatar Finn King
Browse files

AnalysisEfficiency: Changing the way a method for the fake rate estimate is...

AnalysisEfficiency: Changing the way a method for the fake rate estimate is selected and configured.
parent b56a0ac7
No related branches found
No related tags found
No related merge requests found
...@@ -27,8 +27,8 @@ AnalysisEfficiency::AnalysisEfficiency(Configuration& config, std::shared_ptr<De ...@@ -27,8 +27,8 @@ AnalysisEfficiency::AnalysisEfficiency(Configuration& config, std::shared_ptr<De
config_.setDefault<XYVector>("inpixel_cut_edge", {Units::get(5.0, "um"), Units::get(5.0, "um")}); config_.setDefault<XYVector>("inpixel_cut_edge", {Units::get(5.0, "um"), Units::get(5.0, "um")});
config_.setDefault<double>("masked_pixel_distance_cut", 1.); config_.setDefault<double>("masked_pixel_distance_cut", 1.);
config_.setDefault<double>("spatial_cut_sensoredge", 1.); config_.setDefault<double>("spatial_cut_sensoredge", 1.);
config_.setDefault<double>("fake_rate_radius", -1.); config_.setDefault<FakeRateMethod>("fake_rate_method", FakeRateMethod::RADIUS);
config_.setDefault<double>("fake_rate_sensoredge", -1); config_.setDefault<double>("fake_rate_distance", 2.);
config_.setDefault<int>("n_charge_bins", 1000); config_.setDefault<int>("n_charge_bins", 1000);
config_.setDefault<double>("charge_histo_range", 1000.0); config_.setDefault<double>("charge_histo_range", 1000.0);
...@@ -39,8 +39,8 @@ AnalysisEfficiency::AnalysisEfficiency(Configuration& config, std::shared_ptr<De ...@@ -39,8 +39,8 @@ AnalysisEfficiency::AnalysisEfficiency(Configuration& config, std::shared_ptr<De
m_inpixelEdgeCut = config_.get<XYVector>("inpixel_cut_edge"); m_inpixelEdgeCut = config_.get<XYVector>("inpixel_cut_edge");
m_maskedPixelDistanceCut = config_.get<int>("masked_pixel_distance_cut"); m_maskedPixelDistanceCut = config_.get<int>("masked_pixel_distance_cut");
spatial_cut_sensoredge = config_.get<double>("spatial_cut_sensoredge"); spatial_cut_sensoredge = config_.get<double>("spatial_cut_sensoredge");
m_fake_rate_radius = config_.get<double>("fake_rate_radius"); m_fake_rate_method = config_.get<FakeRateMethod>("fake_rate_method");
m_fake_rate_sensoredge = config_.get<double>("fake_rate_sensoredge"); m_fake_rate_distance = config_.get<double>("fake_rate_distance");
m_n_charge_bins = config_.get<int>("n_charge_bins"); m_n_charge_bins = config_.get<int>("n_charge_bins");
m_charge_histo_range = config_.get<double>("charge_histo_range"); m_charge_histo_range = config_.get<double>("charge_histo_range");
} }
...@@ -567,8 +567,8 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard) ...@@ -567,8 +567,8 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard)
} }
// fake rate analysis // fake rate analysis
if(m_fake_rate_radius > 0) { if(m_fake_rate_method == FakeRateMethod::RADIUS) {
LOG_ONCE(STATUS) << "Estimating fake rate based on radial cut around pixels."; LOG_ONCE(STATUS) << "Estimating fake rate based on radial cut around pixels (RADIUS method).";
int fake_hits = 0; int fake_hits = 0;
int fake_clusters = 0; int fake_clusters = 0;
...@@ -583,7 +583,7 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard) ...@@ -583,7 +583,7 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard)
// discard tracks without intercept, using a tolerance defined // discard tracks without intercept, using a tolerance defined
// by the radial cut which we will use later. // by the radial cut which we will use later.
if(m_detector->hasIntercept(track.get(), -m_fake_rate_radius)) { if(m_detector->hasIntercept(track.get(), -m_fake_rate_distance)) {
continue; continue;
} }
...@@ -593,7 +593,7 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard) ...@@ -593,7 +593,7 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard)
double norm_xdistance = (interceptLocal.X() - cluster->local().x()) / m_detector->getPitch().X(); double norm_xdistance = (interceptLocal.X() - cluster->local().x()) / m_detector->getPitch().X();
double norm_ydistance = (interceptLocal.Y() - cluster->local().y()) / m_detector->getPitch().Y(); double norm_ydistance = (interceptLocal.Y() - cluster->local().y()) / m_detector->getPitch().Y();
double norm_rdistance = sqrt(norm_xdistance * norm_xdistance + norm_ydistance * norm_ydistance); double norm_rdistance = sqrt(norm_xdistance * norm_xdistance + norm_ydistance * norm_ydistance);
if(norm_rdistance < m_fake_rate_radius) { if(norm_rdistance < m_fake_rate_distance) {
track_too_close = true; track_too_close = true;
break; break;
} }
...@@ -619,9 +619,9 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard) ...@@ -619,9 +619,9 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard)
fakePixelPerEventVsTime->Fill(static_cast<double>(Units::convert(event->start(), "s")), fake_hits); fakePixelPerEventVsTime->Fill(static_cast<double>(Units::convert(event->start(), "s")), fake_hits);
fakePixelPerEventVsTimeLong->Fill(static_cast<double>(Units::convert(event->start(), "s")), fake_hits); fakePixelPerEventVsTimeLong->Fill(static_cast<double>(Units::convert(event->start(), "s")), fake_hits);
hFakeClusterPerEvent->Fill(fake_clusters); hFakeClusterPerEvent->Fill(fake_clusters);
}
} else if(m_fake_rate_sensoredge > 0) { if(m_fake_rate_method == FakeRateMethod::EDGE) {
LOG_ONCE(STATUS) << "Estimating fake rate based on events without DUT intercepting tracks."; LOG_ONCE(STATUS) << "Estimating fake rate based on events without DUT intercepting tracks (EDGE method).";
bool track_in_active = false; bool track_in_active = false;
int fake_hits = 0; int fake_hits = 0;
...@@ -632,7 +632,7 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard) ...@@ -632,7 +632,7 @@ StatusCode AnalysisEfficiency::run(const std::shared_ptr<Clipboard>& clipboard)
// check if one of the tracks intercepts the dut // check if one of the tracks intercepts the dut
// A positive value of m_fake_rate_radius means that we want to // A positive value of m_fake_rate_radius means that we want to
// check an area that is larger than the sensor matrix. // check an area that is larger than the sensor matrix.
if(m_detector->hasIntercept(track.get(), -m_fake_rate_sensoredge)) { if(m_detector->hasIntercept(track.get(), -m_fake_rate_distance)) {
track_in_active = true; track_in_active = true;
break; break;
} }
......
...@@ -96,8 +96,13 @@ namespace corryvreckan { ...@@ -96,8 +96,13 @@ namespace corryvreckan {
TH1D* hFakeClusterCharge; TH1D* hFakeClusterCharge;
TH1D* hFakeClusterSize; TH1D* hFakeClusterSize;
double m_chi2ndofCut, m_timeCutFrameEdge, m_inpixelBinSize, spatial_cut_sensoredge, m_fake_rate_radius, enum class FakeRateMethod {
m_fake_rate_sensoredge, m_charge_histo_range; RADIUS,
EDGE,
} m_fake_rate_method;
double m_chi2ndofCut, m_timeCutFrameEdge, m_inpixelBinSize, spatial_cut_sensoredge, m_fake_rate_distance,
m_charge_histo_range;
int m_n_charge_bins; int m_n_charge_bins;
XYVector m_inpixelEdgeCut; XYVector m_inpixelEdgeCut;
int m_maskedPixelDistanceCut = 1; int m_maskedPixelDistanceCut = 1;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment