1#include "Tracking/Reco/LinearSeedFinder.h"
13 truth_matching_tool_ = std::make_shared<tracking::sim::TruthMatchingTool>();
19 "out_seed_collection",
getName() +
"LinearRecoilSeedTracks");
23 parameters.
get<std::string>(
"input_hits_collection",
"DigiRecoilSimHits");
25 parameters.
get<std::string>(
"input_rec_hits_collection",
"EcalRecHits");
27 input_pass_name_ = parameters.
get<std::string>(
"input_pass_name",
"");
29 sim_particles_passname_ =
30 parameters.
get<std::string>(
"sim_particles_passname");
32 sim_particles_events_passname_ =
33 parameters.
get<std::string>(
"sim_particles_events_passname");
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");
43 layer12_midpoint_ = parameters.
get<
double>(
"layer12_midpoint");
44 layer23_midpoint_ = parameters.
get<
double>(
"layer23_midpoint");
45 layer34_midpoint_ = parameters.
get<
double>(
"layer34_midpoint");
49 auto start = std::chrono::high_resolution_clock::now();
50 std::vector<ldmx::StraightTrack> straight_seed_tracks;
56 const auto& ecal_rec_hit =
event.getCollection<
ldmx::EcalHit>(
59 std::vector<std::array<double, 3>> first_layer_ecal_rec_hits;
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()});
70 if ((recoil_hits.size() < 2) || (first_layer_ecal_rec_hits.empty()) ||
71 (uniqueLayersHit(recoil_hits) < 2)) {
73 n_seeds_ += straight_seed_tracks.size();
79 std::map<int, ldmx::SimParticle> particle_map;
80 if (event.
exists(
"SimParticles", sim_particles_events_passname_)) {
82 "SimParticles", sim_particles_passname_);
83 truth_matching_tool_->setup(particle_map, recoil_hits);
87 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
90 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
94 "RecoilSimHits", input_pass_name_);
96 "TargetScoringPlaneHits", input_pass_name_);
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);
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;
112 for (
const auto& point : recoil_hits) {
114 float x = point.getGlobalPosition()[0];
115 int track_id = point.getTrackIds()[0];
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;
122 const auto& sim_hits = sim_range_it->second;
124 for (
const auto* sim_hit : sim_hits) {
125 float z = sim_hit->getPosition()[2];
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);
137 else if (x < layer23_midpoint_) {
138 if (z > layer12_midpoint_ && z < layer23_midpoint_) {
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);
154 if (z > layer34_midpoint_) {
155 second_two_layers.emplace_back(point, *sim_hit,
165 auto first_sensor_combos = processMeasurements(first_two_layers, tg);
166 auto second_sensor_combos = processMeasurements(second_two_layers, tg);
168 for (
const auto& [first_combo_3d_point, first_layer_one, first_layer_two] :
169 first_sensor_combos) {
171 std::optional<ldmx::Measurement>>
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)};
180 first_sensor_point = {first_combo_3d_point,
181 std::get<ldmx::Measurement>(first_layer_one),
185 for (
const auto& [second_combo_3d_point, second_layer_one,
186 second_layer_two] : second_sensor_combos) {
188 std::optional<ldmx::Measurement>>
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)};
196 second_sensor_point = {second_combo_3d_point,
197 std::get<ldmx::Measurement>(second_layer_one),
201 for (
const auto& rec_hit : first_layer_ecal_rec_hits) {
205 seedTracker(first_sensor_point, second_sensor_point, rec_hit);
208 if (seed_track.getChi2() > 0.0) {
209 straight_seed_tracks.push_back(seed_track);
215 n_seeds_ += straight_seed_tracks.size();
218 auto end = std::chrono::high_resolution_clock::now();
220 auto diff = end - start;
221 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
223 first_layer_ecal_rec_hits.clear();
224 straight_seed_tracks.clear();
230 std::optional<ldmx::Measurement>>
233 std::optional<ldmx::Measurement>>
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;
247 all_points.push_back(layer1);
251 if (layer2.has_value()) {
252 all_points.push_back(*layer2);
255 all_points.push_back(layer3);
259 if (layer4.has_value()) {
260 all_points.push_back(*layer4);
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);
273 if (temp_distance < ecal_distance_threshold_) {
275 trk.setInterceptX(b_x);
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));
281 trk.setAllSensorPoints(all_points);
282 trk.setFirstSensorPosition(sensor1);
283 trk.setSecondSensorPosition(sensor2);
284 trk.setFirstLayerEcalRecHit(ecal_one);
285 trk.setDistancetoRecHit(temp_distance);
287 trk.setTargetLocation(0.0, b_x, b_y);
288 trk.setEcalLayer1Location(temp_extrapolated_point);
290 globalChiSquare(sensor1, sensor2, ecal_one, m_x, m_y, b_x, b_y));
294 trk.setCov(seed_cov);
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_);
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_;
321std::array<double, 3> LinearSeedFinder::getPointAtZ(
322 std::array<double, 3> target, std::array<double, 3> measurement,
324 double slope_x = (measurement[1] - target[0]) / (measurement[0] - target[2]);
325 double slope_y = (measurement[2] - target[1]) / (measurement[0] - target[2]);
327 double intercept_x = target[0] - slope_x * target[2];
328 double intercept_y = target[1] - slope_y * target[2];
330 double x_target = slope_x * z_target + intercept_x;
331 double y_target = slope_y * z_target + intercept_y;
333 return {z_target, x_target, y_target};
336std::vector<std::tuple<
337 std::array<double, 3>,
338 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>,
341LinearSeedFinder::processMeasurements(
344 const geo::TrackersTrackingGeometry& tg) {
346 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
349 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>>
351 std::vector<std::tuple<
352 std::array<double, 3>,
353 std::tuple<ldmx::Measurement, ldmx::SimTrackerHit, ldmx::SimTrackerHit>,
356 points_with_measurement;
357 Acts::Vector3 dummy{0., 0., 0.};
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);
364 stereo_measurements.emplace_back(measurement, sim_hit, scoring_hit);
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;
374 for (
const auto& stereo : stereo_measurements) {
375 const auto& [stereo_meas, stereo_hit, stereo_sp] = stereo;
377 const Acts::Surface* axial_surface =
378 tg.getSurface(axial_meas.getLayerID());
379 const Acts::Surface* stereo_surface =
380 tg.getSurface(stereo_meas.getLayerID());
382 if (!axial_surface || !stereo_surface)
continue;
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);
389 points_with_measurement.push_back(
390 {convertToLdmxStdArray(space_point), axial, stereo});
393 }
else if (!axial_measurements.empty()) {
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());
400 Acts::Vector3 axial_meas_hit = axial_surface->localToGlobal(
402 Acts::Vector2(axial_meas.getLocalPosition()[0], 0.0), dummy);
404 points_with_measurement.push_back(
405 {convertToLdmxStdArray(axial_meas_hit), axial, std::nullopt});
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());
413 Acts::Vector3 stereo_meas_hit = stereo_surface->localToGlobal(
415 Acts::Vector2(stereo_meas.getLocalPosition()[0], 0.0), dummy);
417 points_with_measurement.push_back(
418 {convertToLdmxStdArray(stereo_meas_hit), stereo, std::nullopt});
421 return points_with_measurement;
425std::array<double, 3> LinearSeedFinder::convertToLdmxStdArray(
426 const Acts::Vector3& vec) {
427 return {vec.x(), vec.y(), vec.z()};
432std::tuple<Acts::Vector3, Acts::Vector3, Acts::Vector3>
433LinearSeedFinder::getSurfaceVectors(
const Acts::Surface& surface) {
434 Acts::Vector3 dummy{0., 0., 0.};
436 surface.localToGlobal(geometryContext(), Acts::Vector2(1, 0), dummy) -
437 surface.center(geometryContext());
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};
454Acts::Vector3 LinearSeedFinder::simple3DHitV2(
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],
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]};
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]};
472 Acts::Vector3 simpart_path = axial_true_global - hit_on_target;
473 Acts::Vector3 simpart_unit = simpart_path.normalized();
477 Acts::Vector3 axial_origin = axial_surface.center(geometryContext());
478 Acts::Vector3 stereo_origin = stereo_surface.center(geometryContext());
482 Acts::Vector3 delta_sensors = stereo_origin - axial_origin;
486 double dx_proj = (simpart_unit[0] / simpart_unit[2]) *
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);
504 stereo_v_value = 0.0;
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;
514 Acts::Vector3 axst_global_useproj = axial_surface.localToGlobal(
516 Acts::Vector2(u_intercept_useproj, v_intercept_useproj), dummy);
517 Acts::Vector3 dummy_stereo_proj = stereo_surface.localToGlobal(
519 Acts::Vector2(u_intercept_useproj, v_intercept_useproj), dummy);
522 Acts::Vector3 reconstructed_hit{dummy_stereo_proj[0], axst_global_useproj[1],
523 axst_global_useproj[2]};
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";
530 return reconstructed_hit;
533double LinearSeedFinder::dotProduct(
const Acts::Vector3& v1,
534 const Acts::Vector3& v2) {
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];
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)};
553 Eigen::Matrix<double, 6, 4> a_mat;
554 Eigen::Matrix<double, 6, 1> d_vec, w_vec;
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;
561 d_vec << x_pos1, y_pos1, x_pos2, y_pos2, x_pos3, y_pos3;
564 w_vec = Eigen::Matrix<double, 6, 1>(weights.data());
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);
571 Eigen::Matrix4d covariance_matrix = at_w_a.inverse();
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)};
582 return {param_vec(0), param_vec(1), param_vec(2), param_vec(3),
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));
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,
596 double chi2_x = 0, chi2_y = 0;
598 (m_x * first_sensor[0] + b_x - first_sensor[1]) / recoil_uncertainty_[0],
601 (m_y * first_sensor[0] + b_y - first_sensor[2]) / recoil_uncertainty_[1],
604 chi2_x += pow((m_x * second_sensor[0] + b_x - second_sensor[1]) /
605 recoil_uncertainty_[0],
607 chi2_y += pow((m_y * second_sensor[0] + b_y - second_sensor[2]) /
608 recoil_uncertainty_[1],
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);
614 return chi2_x + chi2_y;
617int LinearSeedFinder::uniqueLayersHit(
618 const std::vector<ldmx::Measurement>& digi_points) {
619 std::vector<ldmx::Measurement> sorted_points = digi_points;
622 std::sort(sorted_points.begin(), sorted_points.end(),
624 return meas1.getGlobalPosition()[0] <
625 meas2.getGlobalPosition()[0];
629 auto last = std::unique(
630 sorted_points.begin(), sorted_points.end(),
632 return meas1.getGlobalPosition()[0] == meas2.getGlobalPosition()[0];
636 return std::distance(sorted_points.begin(), last);
#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.
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.
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.
std::array< float, 2 > getLocalPosition() const
Class representing a simulated particle.
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 ¶meters) 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...