1#include "Tracking/Reco/LinearTruthTracking.h"
15 "LinearRecoilTruthTracks");
19 parameters.
get<std::string>(
"input_hits_collection",
"RecoilSimHits");
21 parameters.
get<std::string>(
"input_recHits_collection",
"EcalRecHits");
23 input_pass_name_ = parameters.
get<std::string>(
"input_pass_name");
25 ecal_first_layer_z_threshold_ =
26 parameters.
get<
double>(
"ecal_first_layer_z_threshold");
30 auto start = std::chrono::high_resolution_clock::now();
31 std::vector<ldmx::StraightTrack> straight_truth_tracks;
34 const std::vector<ldmx::SimTrackerHit> recoil_hits =
37 const std::vector<ldmx::EcalHit> ecal_rec_hit =
41 std::vector<std::array<double, 3>> first_layer_ecal_rec_hits;
43 for (
const auto& x_ecal : ecal_rec_hit) {
44 if (x_ecal.getZPos() < ecal_first_layer_z_threshold_) {
45 first_layer_ecal_rec_hits.push_back(
46 {x_ecal.getZPos(), x_ecal.getXPos(), x_ecal.getYPos()});
51 if (recoil_hits.size() < 2) {
53 n_truth_ += straight_truth_tracks.size();
58 std::vector<ldmx::SimTrackerHit> truth_points;
61 for (
const auto& point : recoil_hits) {
62 if (point.getTrackID() == 1) {
63 truth_points.push_back(point);
68 if (truth_points.size() > 1) {
69 straight_truth_tracks.push_back(
70 truthTracker(truth_points, first_layer_ecal_rec_hits));
74 n_truth_ += straight_truth_tracks.size();
78 n_truth_ += straight_truth_tracks.size();
81 auto end = std::chrono::high_resolution_clock::now();
82 auto diff = end - start;
83 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
85 first_layer_ecal_rec_hits.clear();
86 straight_truth_tracks.clear();
91 const std::vector<ldmx::SimTrackerHit>& points,
92 std::vector<std::array<double, 3>>& ecal_points) {
97 auto [m_x, b_x, m_y, b_y] = fit3DLine(points);
100 trk.setInterceptX(b_x);
102 trk.setInterceptY(b_y);
103 trk.setTheta(std::atan2(m_y, std::sqrt(1 + m_x * m_x)));
104 trk.setPhi(std::atan2(m_x, 1.0));
105 trk.setTargetLocation(0.0, b_x, b_y);
106 trk.setNhits(points.size());
107 trk.setNdf(points.size() - 2);
109 trk.setTrackID(points[0].getTrackID());
110 trk.setPdgID(points[0].getPdgID());
111 trk.setTruthProb(1.0);
113 trk.setAllTruthSensorPoints(points);
115 trk.setChi2(globalChiSquare(points, m_x, m_y, b_x, b_y));
118 if (ecal_points.size() > 0) {
119 double ecal_first_layer_z = ecal_points[0][0];
122 std::array<double, 3> extrapolated_point = {ecal_first_layer_z,
123 m_x * ecal_first_layer_z + b_x,
124 m_y * ecal_first_layer_z + b_y};
126 const std::array<double, 3>* closest_rec_hit =
nullptr;
127 double min_distance = std::numeric_limits<double>::max();
130 for (
const auto& ecal_rec_hit : ecal_points) {
131 double temp_distance =
132 calculateDistance(extrapolated_point, ecal_rec_hit);
134 if (temp_distance < min_distance) {
135 min_distance = temp_distance;
136 closest_rec_hit = &ecal_rec_hit;
140 if (closest_rec_hit !=
nullptr) {
141 trk.setFirstLayerEcalRecHit(*closest_rec_hit);
142 trk.setDistancetoRecHit(min_distance);
143 trk.setEcalLayer1Location(extrapolated_point);
151 ldmx_log(info) <<
"AVG Time / Event: " << std::fixed << std::setprecision(1)
152 << processing_time_ / n_events_ <<
" ms";
153 ldmx_log(info) <<
"Number of Events without enough seed points: "
155 ldmx_log(info) <<
"Number of Events with only 1 trackID == 1 point: "
157 ldmx_log(info) <<
"Total Truth Tracks / Event: " << n_truth_ <<
"/"
161double LinearTruthTracking::calculateDistance(
162 const std::array<double, 3>& point1,
const std::array<double, 3>& point2) {
163 return sqrt(pow(point1[1] - point2[1], 2) + pow(point1[2] - point2[2], 2));
166std::tuple<double, double, double, double> LinearTruthTracking::fit3DLine(
167 const std::vector<ldmx::SimTrackerHit>& points) {
168 std::vector<double> z_vals, x_vals, y_vals;
170 for (
const auto& point : points) {
171 const auto& position = point.getPosition();
172 x_vals.push_back(position[0]);
173 y_vals.push_back(position[1]);
174 z_vals.push_back(position[2]);
177 size_t num_points = points.size();
178 Eigen::MatrixXd a_mat(num_points, 2);
179 Eigen::VectorXd x_vec(num_points), y_vec(num_points);
181 for (
size_t i = 0; i < num_points; ++i) {
182 a_mat(i, 0) = z_vals[i];
184 x_vec(i) = x_vals[i];
185 y_vec(i) = y_vals[i];
188 Eigen::Matrix2d at_a = a_mat.transpose() * a_mat;
189 Eigen::Vector2d at_x = a_mat.transpose() * x_vec;
190 Eigen::Vector2d at_y = a_mat.transpose() * y_vec;
192 Eigen::Vector2d x_params = at_a.ldlt().solve(at_x);
193 Eigen::Vector2d y_params = at_a.ldlt().solve(at_y);
196 return {x_params(0), x_params(1), y_params(0), y_params(1)};
199double LinearTruthTracking::globalChiSquare(
200 const std::vector<ldmx::SimTrackerHit>& points,
double m_x,
double m_y,
201 double b_x,
double b_y) {
202 double chi2_x = 0, chi2_y = 0;
204 for (
const auto& point : points) {
205 const auto& position = point.getPosition();
206 double x_pos = position[0];
207 double y_pos = position[1];
208 double z_pos = position[2];
210 double x_fit = m_x * z_pos + b_x;
211 double y_fit = m_y * z_pos + b_y;
213 chi2_x += std::pow((x_pos - x_fit), 2);
214 chi2_y += std::pow((y_pos - y_fit), 2);
217 return chi2_x + chi2_y;
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Implements an event buffer system for storing event data.
Class which represents the process under execution.
Class encapsulating parameters for configuring a processor.
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Stores reconstructed hit information from the ECAL.
Represents a simulated tracker hit in the simulation.
std::string input_rec_hits_collection_
The name of the tagger Tracks (only for Recoil Seeding)
void produce(framework::Event &event) override
Run the processor and create a collection of results which indicate if a charge particle can be found...
std::string input_hits_collection_
The name of the input hits collection to use in finding seeds..
LinearTruthTracking(const std::string &name, framework::Process &process)
Constructor.
std::string out_trk_collection_
The name of the output collection of seeds to be stored.
void onProcessEnd() override
Output event statistics.
void configure(framework::config::Parameters ¶meters) override
Configure the processor using the given user specified parameters.
a helper base class providing some methods to shorten access to common conditions used within the tra...
The measurement calibrator can be a function or a class/struct able to retrieve the sim hits containe...