LDMX Software
LinearSeedFinder.cxx
1#include "Tracking/Reco/LinearSeedFinder.h"
2
3#include "Eigen/Dense"
4
5namespace tracking {
6namespace reco {
7
8LinearSeedFinder::LinearSeedFinder(const std::string& name,
9 framework::Process& process)
10 : TrackingGeometryUser(name, process) {}
11
13 truth_matching_tool_ = std::make_shared<tracking::sim::TruthMatchingTool>();
14}
15
17 // Output seed name
18 out_seed_collection_ = parameters.get<std::string>(
19 "out_seed_collection", getName() + "LinearRecoilSeedTracks");
20
21 // Input strip hits_
23 parameters.get<std::string>("input_hits_collection", "DigiRecoilSimHits");
25 parameters.get<std::string>("input_rec_hits_collection", "EcalRecHits");
26
27 input_pass_name_ = parameters.get<std::string>("input_pass_name", "");
28
29 sim_particles_passname_ =
30 parameters.get<std::string>("sim_particles_passname");
31
32 sim_particles_events_passname_ =
33 parameters.get<std::string>("sim_particles_events_passname");
34
35 // the uncertainty is sigma_x = 6 microns and sigma_y = 20./sqrt(12)
36 recoil_uncertainty_ =
37 parameters.get<std::vector<double>>("recoil_uncertainty", {0.006, 0.085});
38 ecal_uncertainty_ = parameters.get<double>("ecal_uncertainty", {3.87});
39 ecal_distance_threshold_ = parameters.get<double>("ecal_distance_threshold");
40 ecal_first_layer_z_threshold_ =
41 parameters.get<double>("ecal_first_layer_z_threshold");
42
43 layer12_midpoint_ = parameters.get<double>("layer12_midpoint");
44 layer23_midpoint_ = parameters.get<double>("layer23_midpoint");
45 layer34_midpoint_ = parameters.get<double>("layer34_midpoint");
46}
47
49 auto start = std::chrono::high_resolution_clock::now();
50 std::vector<ldmx::StraightTrack> straight_seed_tracks;
51 n_events_++;
52 auto tg{geometry()};
53
54 const auto& recoil_hits = event.getCollection<ldmx::Measurement>(
55 input_hits_collection_, input_pass_name_);
56 const auto& ecal_rec_hit = event.getCollection<ldmx::EcalHit>(
57 input_rec_hits_collection_, input_pass_name_);
58
59 std::vector<std::array<double, 3>> first_layer_ecal_rec_hits;
60
61 // Find RecHits at first layer_ of ECal
62 for (const auto& x_ecal : ecal_rec_hit) {
63 if (x_ecal.getZPos() < ecal_first_layer_z_threshold_) {
64 first_layer_ecal_rec_hits.push_back(
65 {x_ecal.getZPos(), x_ecal.getXPos(), x_ecal.getYPos()});
66 } // if first layer_ of Ecal
67 } // for positions in ecalRecHit
68
69 // Check if we would fit empty seeds, if so: end tracking
70 if ((recoil_hits.size() < 2) || (first_layer_ecal_rec_hits.empty()) ||
71 (uniqueLayersHit(recoil_hits) < 2)) {
72 n_missing_++;
73 n_seeds_ += straight_seed_tracks.size();
74 event.add(out_seed_collection_, straight_seed_tracks);
75 return;
76 }
77
78 // Setup truth map
79 std::map<int, ldmx::SimParticle> particle_map;
80 if (event.exists("SimParticles", sim_particles_events_passname_)) {
81 particle_map = event.getMap<int, ldmx::SimParticle>(
82 "SimParticles", sim_particles_passname_);
83 truth_matching_tool_->setup(particle_map, recoil_hits);
84 }
85
86 std::vector<
87 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
88 first_two_layers;
89 std::vector<
90 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
91 second_two_layers;
92
93 const auto& recoil_sim_hits = event.getCollection<ldmx::SimTrackerHit>(
94 "RecoilSimHits", input_pass_name_);
95 const auto& scoring_hits = event.getCollection<ldmx::SimTrackerHit>(
96 "TargetScoringPlaneHits", input_pass_name_);
97
98 // Index all sim hits_ by track ID
99 std::unordered_map<int, std::vector<const ldmx::SimTrackerHit*>>
100 sim_hits_by_track_id;
101 for (const auto& hit : recoil_sim_hits) {
102 sim_hits_by_track_id[hit.getTrackID()].push_back(&hit);
103 } // for sim hits_
104
105 // Index scoring hits_ by track ID (one positive scoring plane hit / ID)
106 std::unordered_map<int, const ldmx::SimTrackerHit*> scoring_hit_map;
107 for (const auto& sp_hit : scoring_hits) {
108 if (sp_hit.getPosition()[2] > 0)
109 scoring_hit_map[sp_hit.getTrackID()] = &sp_hit;
110 } // for sp hits_
111
112 for (const auto& point : recoil_hits) {
113 // x is in tracking coordinates, z is in ldmx coordinates
114 float x = point.getGlobalPosition()[0];
115 int track_id = point.getTrackIds()[0];
116
117 // get the key value = track_id
118 auto sim_range_it = sim_hits_by_track_id.find(track_id);
119 if (sim_range_it == sim_hits_by_track_id.end()) continue;
120
121 // access map value at track_id
122 const auto& sim_hits = sim_range_it->second;
123
124 for (const auto* sim_hit : sim_hits) {
125 float z = sim_hit->getPosition()[2];
126
127 if (x < layer12_midpoint_) {
128 if (z < layer12_midpoint_) {
129 auto sp_it = scoring_hit_map.find(track_id);
130 if (sp_it != scoring_hit_map.end()) {
131 first_two_layers.emplace_back(point, *sim_hit, *sp_it->second);
132 break;
133 } // add the associated scoring plane hit (will be needed for 3D
134 // reconstruction)
135 } // associate 1st layer_ sim hit
136 } // check if recoil hit is 1st layer_
137 else if (x < layer23_midpoint_) {
138 if (z > layer12_midpoint_ && z < layer23_midpoint_) {
139 first_two_layers.emplace_back(point, *sim_hit, ldmx::SimTrackerHit());
140 break;
141 } // associate 2nd layer_ sim hit
142 } // check if recoil hit is 2nd layer_
143 else if (x < layer34_midpoint_) {
144 if (z > layer23_midpoint_ && z < layer34_midpoint_) {
145 auto sp_it = scoring_hit_map.find(track_id);
146 if (sp_it != scoring_hit_map.end()) {
147 second_two_layers.emplace_back(point, *sim_hit, *sp_it->second);
148 break;
149 } // add the associated scoring plane hit (will be needed for 3D
150 // reconstruction)
151 } // associate 3rd layer_ sim hits_
152 } // check if recoil hit is 3rd layer_
153 else {
154 if (z > layer34_midpoint_) {
155 second_two_layers.emplace_back(point, *sim_hit,
157 break;
158 } // associate 4th layer_ sim hits_
159 } // check if recoil hits_ is 4th layer_
160
161 } // loop through simhits
162 } // loop through recoil hits_
163
164 // Reconstruct 3D sensor points on which to do fitting
165 auto first_sensor_combos = processMeasurements(first_two_layers, tg);
166 auto second_sensor_combos = processMeasurements(second_two_layers, tg);
167
168 for (const auto& [first_combo_3d_point, first_layer_one, first_layer_two] :
169 first_sensor_combos) {
170 std::tuple<std::array<double, 3>, ldmx::Measurement,
171 std::optional<ldmx::Measurement>>
172 first_sensor_point;
173
174 if (first_layer_two.has_value()) {
175 first_sensor_point = {first_combo_3d_point,
176 std::get<ldmx::Measurement>(first_layer_one),
177 std::get<ldmx::Measurement>(*first_layer_two)};
178 } // check whether we did reconstruction or...
179 else {
180 first_sensor_point = {first_combo_3d_point,
181 std::get<ldmx::Measurement>(first_layer_one),
182 std::nullopt};
183 } //...we are taking only one layer_ as the measurement (axial or stereo)
184
185 for (const auto& [second_combo_3d_point, second_layer_one,
186 second_layer_two] : second_sensor_combos) {
187 std::tuple<std::array<double, 3>, ldmx::Measurement,
188 std::optional<ldmx::Measurement>>
189 second_sensor_point;
190 if (second_layer_two.has_value()) {
191 second_sensor_point = {second_combo_3d_point,
192 std::get<ldmx::Measurement>(second_layer_one),
193 std::get<ldmx::Measurement>(*second_layer_two)};
194 } // check whether we did reconstruction or...
195 else {
196 second_sensor_point = {second_combo_3d_point,
197 std::get<ldmx::Measurement>(second_layer_one),
198 std::nullopt};
199 } //...we are taking only one layer_ as the measurement (axial or stereo)
200
201 for (const auto& rec_hit : first_layer_ecal_rec_hits) {
202 // Do fitting on 2 sensor + 1 recHit combinations = 1 degree of freedom
203 // for linear fit
204 ldmx::StraightTrack seed_track =
205 seedTracker(first_sensor_point, second_sensor_point, rec_hit);
206
207 // Seed passed RecHit distance check, add it
208 if (seed_track.getChi2() > 0.0) {
209 straight_seed_tracks.push_back(seed_track);
210 } // if chi2 > 0
211 } // for rec_hits
212 } // for second recoil tracker
213 } // for first recoil tracker
214
215 n_seeds_ += straight_seed_tracks.size();
216 event.add(out_seed_collection_, straight_seed_tracks);
217
218 auto end = std::chrono::high_resolution_clock::now();
219
220 auto diff = end - start;
221 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
222
223 first_layer_ecal_rec_hits.clear();
224 straight_seed_tracks.clear();
225
226} // produce
227
228ldmx::StraightTrack LinearSeedFinder::seedTracker(
229 const std::tuple<std::array<double, 3>, ldmx::Measurement,
230 std::optional<ldmx::Measurement>>
231 recoil_one,
232 const std::tuple<std::array<double, 3>, ldmx::Measurement,
233 std::optional<ldmx::Measurement>>
234 recoil_two,
235 const std::array<double, 3> ecal_one) {
236 auto [sensor1, layer1, layer2] = recoil_one;
237 auto [sensor2, layer3, layer4] = recoil_two;
238 std::vector<ldmx::Measurement> all_points;
239
240 // TODO: in the case where we don't have all 4 hits_, we will be fitting a
241 // sensor (weighted average of two layers) + single layer_
242 // TODO: or fitting two single layers. Currently, the single layer_ point has
243 // the uncertainty of a sensor assigned to it,
244 // TODO: but this is not a realistic uncertainty for a single layer_...
245 // IF all layers are well-defined, this sequence will add layer1, 2, 3, 4 to
246 // the allPoints vector
247 all_points.push_back(layer1);
248
249 // if layer2 doesn't exist (has_value == False), then the layer1 we added
250 // is either layer1 or 2, depending on which one has_value
251 if (layer2.has_value()) {
252 all_points.push_back(*layer2);
253 }
254
255 all_points.push_back(layer3);
256
257 // if layer4 doesn't exist (has_value == False), then the layer3 we added
258 // is either layer3 or 4, depending on which one has_value
259 if (layer4.has_value()) {
260 all_points.push_back(*layer4);
261 }
262
263 // Fit the 3 points to a 3D straight line, find track location at first layer_
264 // of Ecal, check distance to recHit used in fitting
265 // m = slope ; b = intercept
266 auto [m_x, b_x, m_y, b_y, seed_cov] = fit3DLine(sensor1, sensor2, ecal_one);
267 std::array<double, 3> temp_extrapolated_point = {
268 ecal_one[0], m_x * ecal_one[0] + b_x, m_y * ecal_one[0] + b_y};
269 double temp_distance = calculateDistance(temp_extrapolated_point, ecal_one);
270
272
273 if (temp_distance < ecal_distance_threshold_) {
274 trk.setSlopeX(m_x);
275 trk.setInterceptX(b_x);
276 trk.setSlopeY(m_y);
277 trk.setInterceptY(b_y);
278 trk.setTheta(std::atan2(m_y, std::sqrt(1 + m_x * m_x)));
279 trk.setPhi(std::atan2(m_x, 1.0));
280
281 trk.setAllSensorPoints(all_points);
282 trk.setFirstSensorPosition(sensor1);
283 trk.setSecondSensorPosition(sensor2);
284 trk.setFirstLayerEcalRecHit(ecal_one);
285 trk.setDistancetoRecHit(temp_distance);
286
287 trk.setTargetLocation(0.0, b_x, b_y);
288 trk.setEcalLayer1Location(temp_extrapolated_point);
289 trk.setChi2(
290 globalChiSquare(sensor1, sensor2, ecal_one, m_x, m_y, b_x, b_y));
291 trk.setNhits(3);
292 trk.setNdf(1);
293
294 trk.setCov(seed_cov);
295
296 // truth matching
297 if (truth_matching_tool_->configured()) {
298 auto truth_info = truth_matching_tool_->truthMatch(all_points);
299 trk.setTrackID(truth_info.track_id_);
300 trk.setPdgID(truth_info.pdg_id_);
301 trk.setTruthProb(truth_info.truth_prob_);
302 }
303
304 return trk;
305
306 } // if (track is close enough to EcalRecHit)
307 else {
308 trk.setChi2(-1);
309 return trk;
310 } // else (does not pass the threshold)
311} // SeedTracker
312
314 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(1)
315 << processing_time_ / n_events_ << " ms";
316 ldmx_log(info) << "Total Seeds/Events: " << n_seeds_ << "/" << n_events_;
317 ldmx_log(info) << "not enough seed points " << n_missing_;
318
319} // onProcessEnd
320
321std::array<double, 3> LinearSeedFinder::getPointAtZ(
322 std::array<double, 3> target, std::array<double, 3> measurement,
323 double z_target) {
324 double slope_x = (measurement[1] - target[0]) / (measurement[0] - target[2]);
325 double slope_y = (measurement[2] - target[1]) / (measurement[0] - target[2]);
326
327 double intercept_x = target[0] - slope_x * target[2];
328 double intercept_y = target[1] - slope_y * target[2];
329
330 double x_target = slope_x * z_target + intercept_x;
331 double y_target = slope_y * z_target + intercept_y;
332
333 return {z_target, x_target, y_target};
334}
335
336std::vector<std::tuple<
337 std::array<double, 3>,
338 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>,
339 std::optional<std::tuple<ldmx::Measurement, ldmx::SimTrackerHit,
341LinearSeedFinder::processMeasurements(
342 const std::vector<std::tuple<ldmx::Measurement, ldmx::SimTrackerHit,
343 ldmx::SimTrackerHit>>& measurements,
344 const geo::TrackersTrackingGeometry& tg) {
345 std::vector<
346 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
347 axial_measurements;
348 std::vector<
349 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
350 stereo_measurements;
351 std::vector<std::tuple<
352 std::array<double, 3>,
353 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>,
354 std::optional<std::tuple<ldmx::Measurement, ldmx::SimTrackerHit,
356 points_with_measurement;
357 Acts::Vector3 dummy{0., 0., 0.};
358
359 // Separate measurements into axial and stereo based on layerID
360 for (const auto& [measurement, sim_hit, scoring_hit] : measurements) {
361 if (measurement.getLayerID() % 2 == 0) {
362 axial_measurements.emplace_back(measurement, sim_hit, scoring_hit);
363 } else {
364 stereo_measurements.emplace_back(measurement, sim_hit, scoring_hit);
365 }
366 }
367
368 // If there are measurements in both axial and stereo layer_:
369 // Iterate over axial and stereo measurements to compute 3D points
370 if (!axial_measurements.empty() && !stereo_measurements.empty()) {
371 for (const auto& axial : axial_measurements) {
372 const auto& [axial_meas, axial_hit, axial_sp] = axial;
373
374 for (const auto& stereo : stereo_measurements) {
375 const auto& [stereo_meas, stereo_hit, stereo_sp] = stereo;
376
377 const Acts::Surface* axial_surface =
378 tg.getSurface(axial_meas.getLayerID());
379 const Acts::Surface* stereo_surface =
380 tg.getSurface(stereo_meas.getLayerID());
381
382 if (!axial_surface || !stereo_surface) continue;
383
384 std::vector<ldmx::SimTrackerHit> sim_hits = {axial_hit, stereo_hit};
385 Acts::Vector3 space_point =
386 simple3DHitV2(axial_meas, *axial_surface, stereo_meas,
387 *stereo_surface, axial_sp, sim_hits);
388
389 points_with_measurement.push_back(
390 {convertToLdmxStdArray(space_point), axial, stereo});
391 }
392 }
393 } else if (!axial_measurements.empty()) {
394 // if there are only axial measurements, take them to be our sensor points
395 for (const auto& axial : axial_measurements) {
396 const auto& [axial_meas, axial_hit, axial_sp] = axial;
397 const Acts::Surface* axial_surface =
398 tg.getSurface(axial_meas.getLayerID());
399
400 Acts::Vector3 axial_meas_hit = axial_surface->localToGlobal(
401 geometryContext(),
402 Acts::Vector2(axial_meas.getLocalPosition()[0], 0.0), dummy);
403
404 points_with_measurement.push_back(
405 {convertToLdmxStdArray(axial_meas_hit), axial, std::nullopt});
406 }
407 } else if (!stereo_measurements.empty()) {
408 for (const auto& stereo : stereo_measurements) {
409 const auto& [stereo_meas, stereo_hit, stereo_sp] = stereo;
410 const Acts::Surface* stereo_surface =
411 tg.getSurface(stereo_meas.getLayerID());
412
413 Acts::Vector3 stereo_meas_hit = stereo_surface->localToGlobal(
414 geometryContext(),
415 Acts::Vector2(stereo_meas.getLocalPosition()[0], 0.0), dummy);
416
417 points_with_measurement.push_back(
418 {convertToLdmxStdArray(stereo_meas_hit), stereo, std::nullopt});
419 }
420 }
421 return points_with_measurement;
422}
423
424// ACTS saves its arrays like (x, y, z)
425std::array<double, 3> LinearSeedFinder::convertToLdmxStdArray(
426 const Acts::Vector3& vec) {
427 return {vec.x(), vec.y(), vec.z()};
428}
429
430// Helper function to calculate unit vector by taking advantage of the
431// localToGlobal transformation
432std::tuple<Acts::Vector3, Acts::Vector3, Acts::Vector3>
433LinearSeedFinder::getSurfaceVectors(const Acts::Surface& surface) {
434 Acts::Vector3 dummy{0., 0., 0.};
435 Acts::Vector3 u =
436 surface.localToGlobal(geometryContext(), Acts::Vector2(1, 0), dummy) -
437 surface.center(geometryContext());
438 Acts::Vector3 v =
439 surface.localToGlobal(geometryContext(), Acts::Vector2(0, 1), dummy) -
440 surface.center(geometryContext());
441 Acts::Vector3 w = u.cross(v).normalized();
442 return {u.normalized(), v.normalized(), w};
443}
444
445// estimate the 3d position of the particle as it passes through a stereo/axial
446// pair currently this uses the sim hits_ from the target and the axial sensor
447// to calculate the angle of the particle, which is needed to project the two
448// sensors to the same z. For real data, we could use the position at the
449// target for the tagger and the axial u position of the measured hit (we only
450// project in x, which, in MC, is identically x) assumptions: axial
451// u-direction is identically global-x (tracking-global y)
452// sensors' w-directions are aligned with beam (global-z,
453// tracking-global x)
454Acts::Vector3 LinearSeedFinder::simple3DHitV2(
455 const ldmx::Measurement& axial, const Acts::Surface& axial_surface,
456 const ldmx::Measurement& stereo, const Acts::Surface& stereo_surface,
457 const ldmx::SimTrackerHit& target_sp,
458 std::vector<ldmx::SimTrackerHit> pair_sim_hits) {
459 Acts::Vector3 dummy{0., 0., 0.};
460 Acts::Vector3 hit_on_target{target_sp.getPosition()[0],
461 target_sp.getPosition()[1],
462 target_sp.getPosition()[2]}; // x,y,z
463
464 Acts::Vector3 axial_true_global{pair_sim_hits[0].getPosition()[0],
465 pair_sim_hits[0].getPosition()[1],
466 pair_sim_hits[0].getPosition()[2]};
467 // stereo_true_global is unused
468 Acts::Vector3 stereo_true_global{pair_sim_hits[1].getPosition()[0],
469 pair_sim_hits[1].getPosition()[1],
470 pair_sim_hits[1].getPosition()[2]};
471
472 Acts::Vector3 simpart_path = axial_true_global - hit_on_target;
473 Acts::Vector3 simpart_unit = simpart_path.normalized();
474
475 // Get global positions for strip origins .... actually these are in
476 // tracking coordinates!
477 Acts::Vector3 axial_origin = axial_surface.center(geometryContext());
478 Acts::Vector3 stereo_origin = stereo_surface.center(geometryContext());
479
480 // the tracking-global vector difference between stereo and axial sensor
481 // centers
482 Acts::Vector3 delta_sensors = stereo_origin - axial_origin;
483
484 // calculate the displacement in tracking global x (need to generalize) by
485 // going from tracking x=axial to x=stereo
486 double dx_proj = (simpart_unit[0] / simpart_unit[2]) *
487 delta_sensors[0]; // this looks weird because simpart_unit
488 // is in global-global and delta sensors
489 // is in tracking-global
490
491 // Compute unit vectors for both hits_
492 auto [axial_u, axial_v, axial_w] = getSurfaceVectors(axial_surface);
493 auto [stereo_u, stereo_v, stereo_w] = getSurfaceVectors(stereo_surface);
494 double salpha = dotProduct(axial_v, stereo_u);
495 double cosalpha = dotProduct(axial_u, stereo_u);
496
497 // Get sensor local measured coordinates
498 // Get local position components
499 auto [axial_u_value, axial_v_value] = axial.getLocalPosition();
500 auto [stereo_u_value, stereo_v_value] = stereo.getLocalPosition();
501
502 // Manual correction, since v should always be 0 (insensitive direction)
503 axial_v_value = 0.0;
504 stereo_v_value = 0.0;
505
506 // use the dx_proj as the displacement in u of the axial measurement
507 // it's axial_u_value - dx_proj because u is in the -x direction
508 // this calculation is in the axial frame
509 double v_intercept_useproj =
510 (stereo_u_value - (axial_u_value - dx_proj) * cosalpha) / salpha;
511 double u_intercept_useproj = axial_u_value - dx_proj;
512
513 // convert to tracking global
514 Acts::Vector3 axst_global_useproj = axial_surface.localToGlobal(
515 geometryContext(),
516 Acts::Vector2(u_intercept_useproj, v_intercept_useproj), dummy);
517 Acts::Vector3 dummy_stereo_proj = stereo_surface.localToGlobal(
518 geometryContext(),
519 Acts::Vector2(u_intercept_useproj, v_intercept_useproj), dummy);
520
521 // we want the reconstructed hit to be at the z of the stereo layer
522 Acts::Vector3 reconstructed_hit{dummy_stereo_proj[0], axst_global_useproj[1],
523 axst_global_useproj[2]};
524
525 ldmx_log(debug) << "The particle projected axst measured position is "
526 "(compare with stereo sim position): "
527 << reconstructed_hit[0] << ", " << reconstructed_hit[1]
528 << ", " << reconstructed_hit[2] << "\n";
529
530 return reconstructed_hit;
531}
532
533double LinearSeedFinder::dotProduct(const Acts::Vector3& v1,
534 const Acts::Vector3& v2) {
535 return v1.dot(v2);
536}
537
538std::tuple<double, double, double, double, std::vector<double>>
539LinearSeedFinder::fit3DLine(const std::array<double, 3>& first_recoil,
540 const std::array<double, 3>& second_recoil,
541 const std::array<double, 3>& ecal) {
542 double z_pos1 = first_recoil[0], x_pos1 = first_recoil[1],
543 y_pos1 = first_recoil[2];
544 double z_pos2 = second_recoil[0], x_pos2 = second_recoil[1],
545 y_pos2 = second_recoil[2];
546 double z_pos3 = ecal[0], x_pos3 = ecal[1], y_pos3 = ecal[2];
547
548 std::array<double, 6> weights = {
549 1 / pow(recoil_uncertainty_[0], 2), 1 / pow(recoil_uncertainty_[1], 2),
550 1 / pow(recoil_uncertainty_[0], 2), 1 / pow(recoil_uncertainty_[1], 2),
551 1 / pow(ecal_uncertainty_, 2), 1 / pow(ecal_uncertainty_, 2)};
552
553 Eigen::Matrix<double, 6, 4> a_mat;
554 Eigen::Matrix<double, 6, 1> d_vec, w_vec;
555
556 // Fill the A matrix (z, 1, 0, 0) for x and (0, 0, z, 1) for y
557 a_mat << z_pos1, 1, 0, 0, 0, 0, z_pos1, 1, z_pos2, 1, 0, 0, 0, 0, z_pos2, 1,
558 z_pos3, 1, 0, 0, 0, 0, z_pos3, 1;
559
560 // Fill the d vector with x and y values
561 d_vec << x_pos1, y_pos1, x_pos2, y_pos2, x_pos3, y_pos3;
562
563 // Fill the weights vector
564 w_vec = Eigen::Matrix<double, 6, 1>(weights.data());
565
566 // Solve the weighted least squares system
567 Eigen::MatrixXd at_w_a = a_mat.transpose() * w_vec.asDiagonal() * a_mat;
568 Eigen::MatrixXd at_w_d = a_mat.transpose() * w_vec.asDiagonal() * d_vec;
569 Eigen::VectorXd param_vec = at_w_a.ldlt().solve(at_w_d);
570
571 Eigen::Matrix4d covariance_matrix = at_w_a.inverse();
572
573 // Store only the upper triangular part of the covariance matrix since it is
574 // symmetric
575 std::vector<double> covariance_vector = {
576 covariance_matrix(0, 0), covariance_matrix(0, 1), covariance_matrix(0, 2),
577 covariance_matrix(0, 3), covariance_matrix(1, 1), covariance_matrix(1, 2),
578 covariance_matrix(1, 3), covariance_matrix(2, 2), covariance_matrix(2, 3),
579 covariance_matrix(3, 3)};
580
581 // return {slope_x, intercept_x, slope_y, intercept_y, covariance}
582 return {param_vec(0), param_vec(1), param_vec(2), param_vec(3),
583 covariance_vector};
584} // fit3DLine
585
586double LinearSeedFinder::calculateDistance(
587 const std::array<double, 3>& point1, const std::array<double, 3>& point2) {
588 return sqrt(pow(point1[1] - point2[1], 2) + pow(point1[2] - point2[2], 2));
589} // calculateDistance in xy
590
591double LinearSeedFinder::globalChiSquare(
592 const std::array<double, 3>& first_sensor,
593 const std::array<double, 3>& second_sensor,
594 const std::array<double, 3>& ecal_hit, double m_x, double m_y, double b_x,
595 double b_y) {
596 double chi2_x = 0, chi2_y = 0;
597 chi2_x += pow(
598 (m_x * first_sensor[0] + b_x - first_sensor[1]) / recoil_uncertainty_[0],
599 2);
600 chi2_y += pow(
601 (m_y * first_sensor[0] + b_y - first_sensor[2]) / recoil_uncertainty_[1],
602 2);
603
604 chi2_x += pow((m_x * second_sensor[0] + b_x - second_sensor[1]) /
605 recoil_uncertainty_[0],
606 2);
607 chi2_y += pow((m_y * second_sensor[0] + b_y - second_sensor[2]) /
608 recoil_uncertainty_[1],
609 2);
610
611 chi2_x += pow((m_x * ecal_hit[0] + b_x - ecal_hit[1]) / ecal_uncertainty_, 2);
612 chi2_y += pow((m_y * ecal_hit[0] + b_y - ecal_hit[2]) / ecal_uncertainty_, 2);
613
614 return chi2_x + chi2_y;
615} // globalChiSquare
616
617int LinearSeedFinder::uniqueLayersHit(
618 const std::vector<ldmx::Measurement>& digi_points) {
619 std::vector<ldmx::Measurement> sorted_points = digi_points;
620
621 // Sort by z position in the Recoil
622 std::sort(sorted_points.begin(), sorted_points.end(),
623 [](const ldmx::Measurement& meas1, const ldmx::Measurement& meas2) {
624 return meas1.getGlobalPosition()[0] <
625 meas2.getGlobalPosition()[0];
626 });
627
628 // Remove duplicates to ensure we only keep unique z positions
629 auto last = std::unique(
630 sorted_points.begin(), sorted_points.end(),
631 [](const ldmx::Measurement& meas1, const ldmx::Measurement& meas2) {
632 return meas1.getGlobalPosition()[0] == meas2.getGlobalPosition()[0];
633 });
634
635 // return the number of unique layer_ hits_
636 return std::distance(sorted_points.begin(), last);
637} // uniqueLayersHit
638
639} // namespace reco
640} // namespace tracking
641
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
std::string getName() const
Get the processor name.
Implements an event buffer system for storing event data.
Definition Event.h:42
bool exists(const std::string &name, const std::string &passName, bool unique=true) const
Check for the existence of an object or collection with the given name and pass name in the event.
Definition Event.cxx:105
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
std::array< float, 2 > getLocalPosition() const
Definition Measurement.h:66
Class representing a simulated particle.
Definition SimParticle.h:24
Represents a simulated tracker hit in the simulation.
std::vector< float > getPosition() const
Get the XYZ position of the hit [mm].
std::string out_seed_collection_
The name of the output collection of seeds to be stored.
LinearSeedFinder(const std::string &name, framework::Process &process)
Constructor.
void configure(framework::config::Parameters &parameters) override
Configure the processor using the given user specified parameters.
void produce(framework::Event &event) override
Run the processor and create a collection of results which indicate if a charge particle can be found...
void onProcessEnd() override
Output event statistics.
std::string input_hits_collection_
The name of the input hits collection to use in finding seeds..
void onProcessStart() override
Setup the truth matching.
std::string input_rec_hits_collection_
The name of the tagger Tracks (only for Recoil Seeding)
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...