LDMX Software
LinearTruthTracking.cxx
1#include "Tracking/Reco/LinearTruthTracking.h"
2
3#include "Eigen/Dense"
4
5namespace tracking {
6namespace reco {
7
9 framework::Process& process)
10 : TrackingGeometryUser(name, process) {}
11
13 // Output seed name
14 out_trk_collection_ = parameters.get<std::string>("out_trk_collection",
15 "LinearRecoilTruthTracks");
16
17 // Input strip hits_
19 parameters.get<std::string>("input_hits_collection", "RecoilSimHits");
21 parameters.get<std::string>("input_rec_hits_collection", "EcalRecHits");
22
23 input_pass_name_ = parameters.get<std::string>("input_pass_name");
24
25 ecal_first_layer_z_threshold_ =
26 parameters.get<double>("ecal_first_layer_z_threshold");
27} // configure
28
30 auto start = std::chrono::high_resolution_clock::now();
31 std::vector<ldmx::StraightTrack> straight_truth_tracks;
32 n_events_++;
33
34 const auto& recoil_hits = event.getCollection<ldmx::SimTrackerHit>(
35 input_hits_collection_, input_pass_name_);
36 const auto& ecal_rec_hit = event.getCollection<ldmx::EcalHit>(
37 input_rec_hits_collection_, input_pass_name_);
38
39 std::vector<std::array<double, 3>> first_layer_ecal_rec_hits;
40
41 for (const auto& x_ecal : ecal_rec_hit) {
42 if (x_ecal.getZPos() < ecal_first_layer_z_threshold_) {
43 first_layer_ecal_rec_hits.push_back(
44 {x_ecal.getZPos(), x_ecal.getXPos(), x_ecal.getYPos()});
45 } // if first layer_ of Ecal
46 } // for positions in ecalRecHit
47
48 // Check if we would fit empty seeds.
49 if (recoil_hits.size() < 2) {
50 n_missing_++;
51 n_truth_ += straight_truth_tracks.size();
52 event.add(out_trk_collection_, straight_truth_tracks);
53 return;
54 } // if not enough hits_
55
56 std::vector<ldmx::SimTrackerHit> truth_points;
57
58 // Find only recoil electron out of the hit collection
59 for (const auto& point : recoil_hits) {
60 if (point.getTrackID() == 1) {
61 truth_points.push_back(point);
62 } // if trackID == 1
63 } // for point in hit collection
64
65 // should only find 1 track (recoil e-) or 0 track (no recoil e-)
66 if (truth_points.size() > 1) {
67 straight_truth_tracks.push_back(
68 truthTracker(truth_points, first_layer_ecal_rec_hits));
69 } // if we have truth points
70 else {
71 n_empty_++;
72 n_truth_ += straight_truth_tracks.size();
73 return;
74 } // else no truth points
75
76 n_truth_ += straight_truth_tracks.size();
77 event.add(out_trk_collection_, straight_truth_tracks);
78
79 auto end = std::chrono::high_resolution_clock::now();
80 auto diff = end - start;
81 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
82
83 first_layer_ecal_rec_hits.clear();
84 straight_truth_tracks.clear();
85
86} // produce
87
88ldmx::StraightTrack LinearTruthTracking::truthTracker(
89 const std::vector<ldmx::SimTrackerHit>& points,
90 std::vector<std::array<double, 3>>& ecal_points) {
92
93 // We have pre-selected hits_ corresponding to the Recoil e-, so we only fit
94 // these (no RecHit)
95 auto [m_x, b_x, m_y, b_y] = fit3DLine(points);
96
97 trk.setSlopeX(m_x);
98 trk.setInterceptX(b_x);
99 trk.setSlopeY(m_y);
100 trk.setInterceptY(b_y);
101 trk.setTheta(std::atan2(m_y, std::sqrt(1 + m_x * m_x)));
102 trk.setPhi(std::atan2(m_x, 1.0));
103 trk.setTargetLocation(0.0, b_x, b_y);
104 trk.setNhits(points.size());
105 trk.setNdf(points.size() - 2);
106
107 trk.setTrackID(points[0].getTrackID());
108 trk.setPdgID(points[0].getPdgID());
109 trk.setTruthProb(1.0);
110
111 trk.setAllTruthSensorPoints(points);
112
113 trk.setChi2(globalChiSquare(points, m_x, m_y, b_x, b_y));
114
115 // Z position from the first point in ecalPoints
116 if (ecal_points.size() > 0) {
117 double ecal_first_layer_z = ecal_points[0][0];
118
119 // Extrapolate track to first layer_ of Ecal
120 std::array<double, 3> extrapolated_point = {ecal_first_layer_z,
121 m_x * ecal_first_layer_z + b_x,
122 m_y * ecal_first_layer_z + b_y};
123
124 const std::array<double, 3>* closest_rec_hit = nullptr;
125 double min_distance = std::numeric_limits<double>::max();
126
127 // Loop through ecalPoints to find the closest recHit to our track
128 for (const auto& ecal_rec_hit : ecal_points) {
129 double temp_distance =
130 calculateDistance(extrapolated_point, ecal_rec_hit);
131
132 if (temp_distance < min_distance) {
133 min_distance = temp_distance;
134 closest_rec_hit = &ecal_rec_hit;
135 } // if compare min distance
136 } // for recHits
137
138 if (closest_rec_hit != nullptr) {
139 trk.setFirstLayerEcalRecHit(*closest_rec_hit);
140 trk.setDistancetoRecHit(min_distance);
141 trk.setEcalLayer1Location(extrapolated_point);
142 } // if set trk info
143 } // make sure there are ecalPoints to work with
144
145 return trk;
146} // truthTracker
147
149 ldmx_log(info) << "AVG Time / Event: " << std::fixed << std::setprecision(1)
150 << processing_time_ / n_events_ << " ms";
151 ldmx_log(info) << "Number of Events without enough seed points: "
152 << n_missing_;
153 ldmx_log(info) << "Number of Events with only 1 trackID == 1 point: "
154 << n_empty_;
155 ldmx_log(info) << "Total Truth Tracks / Event: " << n_truth_ << "/"
156 << n_events_;
157} // onProcessEnd
158
159double LinearTruthTracking::calculateDistance(
160 const std::array<double, 3>& point1, const std::array<double, 3>& point2) {
161 return sqrt(pow(point1[1] - point2[1], 2) + pow(point1[2] - point2[2], 2));
162} // calculateDistance
163
164std::tuple<double, double, double, double> LinearTruthTracking::fit3DLine(
165 const std::vector<ldmx::SimTrackerHit>& points) {
166 std::vector<double> z_vals, x_vals, y_vals;
167
168 for (const auto& point : points) {
169 const auto& position = point.getPosition();
170 x_vals.push_back(position[0]);
171 y_vals.push_back(position[1]);
172 z_vals.push_back(position[2]);
173 }
174
175 size_t num_points = points.size();
176 Eigen::MatrixXd a_mat(num_points, 2);
177 Eigen::VectorXd x_vec(num_points), y_vec(num_points);
178
179 for (size_t i = 0; i < num_points; ++i) {
180 a_mat(i, 0) = z_vals[i];
181 a_mat(i, 1) = 1;
182 x_vec(i) = x_vals[i];
183 y_vec(i) = y_vals[i];
184 }
185
186 Eigen::Matrix2d at_a = a_mat.transpose() * a_mat;
187 Eigen::Vector2d at_x = a_mat.transpose() * x_vec;
188 Eigen::Vector2d at_y = a_mat.transpose() * y_vec;
189
190 Eigen::Vector2d x_params = at_a.ldlt().solve(at_x);
191 Eigen::Vector2d y_params = at_a.ldlt().solve(at_y);
192
193 // Return (mx, bx, my, by)
194 return {x_params(0), x_params(1), y_params(0), y_params(1)};
195}
196
197double LinearTruthTracking::globalChiSquare(
198 const std::vector<ldmx::SimTrackerHit>& points, double m_x, double m_y,
199 double b_x, double b_y) {
200 double chi2_x = 0, chi2_y = 0;
201
202 for (const auto& point : points) {
203 const auto& position = point.getPosition();
204 double x_pos = position[0];
205 double y_pos = position[1];
206 double z_pos = position[2];
207
208 double x_fit = m_x * z_pos + b_x;
209 double y_fit = m_y * z_pos + b_y;
210
211 chi2_x += std::pow((x_pos - x_fit), 2);
212 chi2_y += std::pow((y_pos - y_fit), 2);
213 }
214
215 return chi2_x + chi2_y;
216} // globalChiSquare
217
218} // namespace reco
219} // namespace tracking
220
#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.
Definition Event.h:42
Class which represents the process under execution.
Definition Process.h:37
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
Stores reconstructed hit information from the ECAL.
Definition EcalHit.h:19
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 &parameters) 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...