LDMX Software
TrackingRecoDQM.cxx
1#include "Tracking/dqm/TrackingRecoDQM.h"
2
3namespace tracking::dqm {
4
6 track_collection_ = parameters.get<std::string>("track_collection");
7 truth_collection_ = parameters.get<std::string>("truth_collection");
8 measurement_collection_ =
9 parameters.get<std::string>("measurement_collection");
10 measurement_passname_ = parameters.get<std::string>("measurement_passname");
11
12 ecal_sp_coll_name_ = parameters.get<std::string>("ecal_sp_coll_name");
13 ecal_sp_passname_ = parameters.get<std::string>("ecal_sp_passname");
14 target_sp_coll_name_ = parameters.get<std::string>("target_sp_coll_name");
15 target_sp_passname_ = parameters.get<std::string>("target_sp_passname");
16 truth_passname_ = parameters.get<std::string>("truth_passname");
17 truth_events_passname_ = parameters.get<std::string>("truth_events_passname");
18 track_passname_ = parameters.get<std::string>("track_passname");
19 track_collection_events_passname_ =
20 parameters.get<std::string>("track_collection_events_passname");
21
22 title_ = parameters.get<std::string>("title", "tagger_trk_");
23 track_prob_cut_ = parameters.get<double>("trackProb_cut", 0.5);
24 subdetector_ = parameters.get<std::string>("subdetector", "Tagger");
25 track_states_ = parameters.get<std::vector<std::string>>("trackStates", {});
26
27 pidmap_[-321] = PIDBins::kminus;
28 pidmap_[321] = PIDBins::kplus;
29 pidmap_[-211] = PIDBins::piminus;
30 pidmap_[211] = PIDBins::piplus;
31 pidmap_[11] = PIDBins::electron;
32 pidmap_[-11] = PIDBins::positron;
33 pidmap_[2212] = PIDBins::proton;
34 pidmap_[-2212] = PIDBins::antiproton;
35}
36
38 ldmx_log(trace) << "DQM Reading in:" << track_collection_;
39
40 if (!event.exists(track_collection_, track_collection_events_passname_)) {
41 ldmx_log(error) << "TrackCollection " << track_collection_
42 << " with pass = " << track_collection_events_passname_
43 << " not in event";
44 return;
45 }
46
47 auto tracks{
48 event.getCollection<ldmx::Track>(track_collection_, track_passname_)};
49
50 if (!event.exists(measurement_collection_, measurement_passname_)) {
51 ldmx_log(error) << "Measurement collection " << measurement_collection_
52 << " with pass = " << measurement_passname_
53 << " not in event";
54 return;
55 }
56
57 auto measurements{event.getCollection<ldmx::Measurement>(
58 measurement_collection_, measurement_passname_)};
59
60 // The truth track collection
61 if (event.exists(truth_collection_, truth_events_passname_)) {
62 truth_track_collection_ = std::make_shared<ldmx::Tracks>(
63 event.getCollection<ldmx::Track>(truth_collection_, truth_passname_));
64 do_truth_comparison_ = true;
65 }
66
67 // The scoring plane hits_
68 if (event.exists(ecal_sp_coll_name_, ecal_sp_passname_)) {
69 ecal_scoring_hits_ = std::make_shared<std::vector<ldmx::SimTrackerHit>>(
70 event.getCollection<ldmx::SimTrackerHit>(ecal_sp_coll_name_,
71 ecal_sp_passname_));
72 }
73
74 if (event.exists(target_sp_coll_name_, target_sp_passname_)) {
75 target_scoring_hits_ = std::make_shared<std::vector<ldmx::SimTrackerHit>>(
76 event.getCollection<ldmx::SimTrackerHit>(target_sp_coll_name_,
77 target_sp_passname_));
78 }
79
80 ldmx_log(debug) << "Do truth comparison::" << do_truth_comparison_;
81
82 if (do_truth_comparison_) {
83 sortTracks(tracks, unique_tracks_, duplicate_tracks_, fake_tracks_);
84 } else {
85 unique_tracks_ = tracks;
86 }
87
88 ldmx_log(debug) << "Filling histograms for " << tracks.size() << " tracks";
89
90 // General Plots
91 histograms_.fill(title_ + "N_tracks", tracks.size());
92
93 if (!unique_tracks_.empty()) {
94 ldmx_log(debug) << "Track Monitoring on " << unique_tracks_.size()
95 << " Unique Tracks";
96 trackMonitoring(unique_tracks_, measurements, title_, true,
97 do_truth_comparison_);
98 }
99
100 // Fakes and duplicates
101 if (!duplicate_tracks_.empty()) {
102 ldmx_log(debug) << "Track Monitoring on " << duplicate_tracks_.size()
103 << " duplicates";
104 trackMonitoring(duplicate_tracks_, measurements, title_ + "dup_", false,
105 false);
106 }
107 if (!fake_tracks_.empty()) {
108 ldmx_log(debug) << "Track Monitoring on " << fake_tracks_.size()
109 << " fakes";
110 trackMonitoring(fake_tracks_, measurements, title_ + "fake_", false, false);
111 }
112
113 // Track Extrapolation to Ecal Monitoring
114 ldmx_log(trace) << "Track Extrapolation to Ecal Monitoring";
115 if (std::find(track_states_.begin(), track_states_.end(), "target") !=
116 track_states_.end()) {
117 trackStateMonitoring(tracks, ldmx::TrackStateType::AtTarget, "target");
118 }
119
120 if (std::find(track_states_.begin(), track_states_.end(), "ecal") !=
121 track_states_.end()) {
122 trackStateMonitoring(tracks, ldmx::TrackStateType::AtECAL, "ecal");
123 }
124
125 if (std::find(track_states_.begin(), track_states_.end(), "beamOrigin") !=
126 track_states_.end()) {
127 trackStateMonitoring(tracks, ldmx::TrackStateType::AtBeamOrigin,
128 "beamOrigin");
129 }
130
131 // Technical Efficiency plots
132 if (do_truth_comparison_) {
133 ldmx_log(trace) << "Technical Efficiency plots";
134 efficiencyPlots(tracks, measurements, title_);
135 }
136
137 // Tagger Recoil Matching
138
139 // Clear the vectors
140 ldmx_log(trace) << "Clear the vectors";
141 unique_tracks_.clear();
142 duplicate_tracks_.clear();
143 fake_tracks_.clear();
144}
145
147 // Produce the efficiency plots. (TODO::Switch to TEfficiency instead)
148}
149
150void TrackingRecoDQM::efficiencyPlots(
151 const std::vector<ldmx::Track>& tracks,
152 const std::vector<ldmx::Measurement>& measurements,
153 const std::string& title) {
154 // Do all truth track plots - denominator
155
156 histograms_.fill(title + "truth_N_tracks", truth_track_collection_->size());
157 for (auto& truth_trk : *(truth_track_collection_)) {
158 auto truth_phi = truth_trk.getPhi();
159 auto truth_d0 = truth_trk.getD0();
160 auto truth_z0 = truth_trk.getZ0();
161 auto truth_theta = truth_trk.getTheta();
162 auto truth_qop = truth_trk.getQoP();
163 auto truth_p = 1. / abs(truth_trk.getQoP());
164 auto truth_n_hits = truth_trk.getNhits();
165
166 std::vector<double> truth_mom = truth_trk.getMomentum();
167
168 // Polar angle
169 // The momentum in the plane transverse wrt the beam axis
170 auto truth_pt_beam =
171 std::sqrt(truth_mom[1] * truth_mom[1] + truth_mom[2] * truth_mom[2]);
172
173 auto truth_beam_angle = std::atan2(truth_pt_beam, truth_mom[0]);
174
175 histograms_.fill(title + "truth_nHits", truth_n_hits);
176 histograms_.fill(title + "truth_d0", truth_d0);
177 histograms_.fill(title + "truth_z0", truth_z0);
178 histograms_.fill(title + "truth_phi", truth_phi);
179 histograms_.fill(title + "truth_theta", truth_theta);
180 histograms_.fill(title + "truth_qop", truth_qop);
181 histograms_.fill(title + "truth_p", truth_p);
182 histograms_.fill(title + "truth_beam_angle", truth_beam_angle);
183
184 if (pidmap_.count(truth_trk.getPdgID()) != 0) {
185 histograms_.fill(title + "truth_PID", pidmap_[truth_trk.getPdgID()]);
186
187 // TODO do this properly.
188
189 if (pidmap_[truth_trk.getPdgID()] == PIDBins::kminus) {
190 histograms_.fill(title + "truth_kminus_p", truth_p);
191 }
192
193 if (pidmap_[truth_trk.getPdgID()] == PIDBins::kplus) {
194 histograms_.fill(title + "truth_kplus_p", truth_p);
195 }
196
197 if (pidmap_[truth_trk.getPdgID()] == PIDBins::piminus) {
198 histograms_.fill(title + "truth_piminus_p", truth_p);
199 }
200
201 if (pidmap_[truth_trk.getPdgID()] == PIDBins::piplus) {
202 histograms_.fill(title + "truth_piplus_p", truth_p);
203 }
204
205 if (pidmap_[truth_trk.getPdgID()] == PIDBins::electron) {
206 histograms_.fill(title + "truth_electron_p", truth_p);
207 }
208
209 if (pidmap_[truth_trk.getPdgID()] == PIDBins::positron) {
210 histograms_.fill(title + "truth_positron_p", truth_p);
211 }
212
213 if (pidmap_[truth_trk.getPdgID()] == PIDBins::proton) {
214 histograms_.fill(title + "truth_proton_p", truth_p);
215 }
216 }
217
218 } // loop on truth tracks
219
220 for (auto& track : tracks) {
221 // Match the tracks to truth
222 ldmx::Track* truth_trk = nullptr;
223
224 auto it = std::find_if(truth_track_collection_->begin(),
225 truth_track_collection_->end(),
226 [&](const ldmx::Track& tt) {
227 return tt.getTrackID() == track.getTrackID();
228 });
229
230 double track_truth_prob = track.getTruthProb();
231
232 if (it != truth_track_collection_->end() &&
233 track_truth_prob >= track_prob_cut_)
234 truth_trk = &(*it);
235
236 // Match not found
237 if (!truth_trk) return;
238
239 auto truth_phi = truth_trk->getPhi();
240 auto truth_d0 = truth_trk->getD0();
241 auto truth_z0 = truth_trk->getZ0();
242 auto truth_theta = truth_trk->getTheta();
243 auto truth_qop = truth_trk->getQoP();
244 auto truth_p = 1. / abs(truth_trk->getQoP());
245 std::vector<double> truth_mom = truth_trk->getMomentum();
246
247 // Polar angle
248 // The momentum in the plane transverse wrt the beam axis
249 auto truth_pt_beam =
250 std::sqrt(truth_mom[1] * truth_mom[1] + truth_mom[2] * truth_mom[2]);
251
252 auto truth_beam_angle = std::atan2(truth_pt_beam, truth_mom[0]);
253
254 // Fill reco plots for efficiencies - numerator. The quantities are truth
255 histograms_.fill(title + "match_prob", track_truth_prob);
256 histograms_.fill(title + "match_d0", truth_d0);
257 histograms_.fill(title + "match_z0", truth_z0);
258 histograms_.fill(title + "match_phi", truth_phi);
259 histograms_.fill(title + "match_theta", truth_theta);
260 histograms_.fill(title + "match_p", truth_p);
261 histograms_.fill(title + "match_qop", truth_qop);
262 histograms_.fill(title + "match_beam_angle", truth_beam_angle);
263 histograms_.fill(title + "match_nHits", measurements.size());
264 auto dedx_measurements = track.getDedxMeasurements();
265 auto measurement_idxs = track.getMeasurementsIdxs();
266 for (size_t i = 0; i < measurement_idxs.size(); ++i) {
267 histograms_.fill(title + "match_layers_hit",
268 measurements.at(measurement_idxs[i]).getLayer());
269 // Add histogram of the measurement dE/dx
270 if (i < dedx_measurements.size()) {
271 histograms_.fill(title + "match_measurement_dedx",
272 dedx_measurements[i]);
273 }
274 }
275
276 // For some particles
277
278 if (pidmap_.count(truth_trk->getPdgID()) != 0) {
279 histograms_.fill(title + "match_PID", pidmap_[truth_trk->getPdgID()]);
280
281 // TODO do this properly.
282
283 if (pidmap_[truth_trk->getPdgID()] == PIDBins::kminus) {
284 histograms_.fill(title + "match_kminus_p", truth_p);
285 }
286
287 if (pidmap_[truth_trk->getPdgID()] == PIDBins::kplus) {
288 histograms_.fill(title + "match_kplus_p", truth_p);
289 }
290
291 if (pidmap_[truth_trk->getPdgID()] == PIDBins::piminus) {
292 histograms_.fill(title + "match_piminus_p", truth_p);
293 }
294
295 if (pidmap_[truth_trk->getPdgID()] == PIDBins::piplus) {
296 histograms_.fill(title + "match_piplus_p", truth_p);
297 }
298
299 if (pidmap_[truth_trk->getPdgID()] == PIDBins::electron) {
300 histograms_.fill(title + "match_electron_p", truth_p);
301 }
302
303 if (pidmap_[truth_trk->getPdgID()] == PIDBins::positron) {
304 histograms_.fill(title + "match_positron_p", truth_p);
305 }
306
307 if (pidmap_[truth_trk->getPdgID()] == PIDBins::proton) {
308 histograms_.fill(title + "match_proton_p", truth_p);
309 }
310 }
311 } // Loop on tracks
312
313} // Efficiency plots
314
315void TrackingRecoDQM::trackMonitoring(
316 const std::vector<ldmx::Track>& tracks,
317 const std::vector<ldmx::Measurement>& measurements, const std::string title,
318 const bool& doDetail, const bool& doTruth) {
319 for (auto& track : tracks) {
320 // Perigee track parameters
321 auto trk_d0 = track.getD0();
322 auto trk_z0 = track.getZ0();
323 auto trk_qop = track.getQoP();
324 auto trk_theta = track.getTheta();
325 auto trk_phi = track.getPhi();
326 auto trk_p = 1. / abs(trk_qop);
327 auto dedx_measurements = track.getDedxMeasurements();
328 auto measurement_idxs = track.getMeasurementsIdxs();
329 for (size_t i = 0; i < measurement_idxs.size(); ++i) {
330 histograms_.fill(title + "layers_hit",
331 measurements.at(measurement_idxs[i]).getLayer());
332 // Add histogram of the measurement dE/dx
333 if (i < dedx_measurements.size()) {
334 histograms_.fill(title + "measurement_dedx", dedx_measurements[i]);
335 }
336 }
337
338 std::vector<double> trk_mom = track.getMomentum();
339
340 // The transverse momentum in the bending plane
341 double pt_bending =
342 std::sqrt(trk_mom[0] * trk_mom[0] + trk_mom[1] * trk_mom[1]);
343
344 // The momentum in the plane transverse wrt the beam axis
345 double pt_beam =
346 std::sqrt(trk_mom[1] * trk_mom[1] + trk_mom[2] * trk_mom[2]);
347
348 // Covariance matrix
349 Acts::BoundSquareMatrix cov =
350 tracking::sim::utils::unpackCov(track.getPerigeeCov());
351
352 double sigmad0 = sqrt(
353 cov(Acts::BoundIndices::eBoundLoc0, Acts::BoundIndices::eBoundLoc0));
354 double sigmaz0 = sqrt(
355 cov(Acts::BoundIndices::eBoundLoc1, Acts::BoundIndices::eBoundLoc1));
356 double sigmaphi =
357 sqrt(cov(Acts::BoundIndices::eBoundPhi, Acts::BoundIndices::eBoundPhi));
358 double sigmatheta = sqrt(
359 cov(Acts::BoundIndices::eBoundTheta, Acts::BoundIndices::eBoundTheta));
360 double sigmaqop = sqrt(cov(Acts::BoundIndices::eBoundQOverP,
361 Acts::BoundIndices::eBoundQOverP));
362 double sigmap = (1. / trk_qop) * (1. / trk_qop) * sigmaqop;
363
364 histograms_.fill(title + "d0", trk_d0);
365 histograms_.fill(title + "z0", trk_z0);
366 histograms_.fill(title + "qop", trk_qop);
367 histograms_.fill(title + "phi", trk_phi);
368 histograms_.fill(title + "theta", trk_theta);
369 histograms_.fill(title + "p", std::abs(1. / trk_qop));
370
371 if (doDetail) {
372 histograms_.fill(title + "px", trk_mom[0]);
373 histograms_.fill(title + "py", trk_mom[1]);
374 histograms_.fill(title + "pz", trk_mom[2]);
375
376 histograms_.fill(title + "pt_bending", pt_bending);
377 histograms_.fill(title + "pt_beam", pt_beam);
378
379 histograms_.fill(title + "nHits", track.getNhits());
380 histograms_.fill(title + "Chi2", track.getChi2());
381 histograms_.fill(title + "ndf", track.getNdf());
382 histograms_.fill(title + "Chi2_per_ndf",
383 track.getChi2() / track.getNdf());
384 histograms_.fill(title + "nShared", track.getNsharedHits());
385
386 histograms_.fill(title + "d0_err", sigmad0);
387 histograms_.fill(title + "z0_err", sigmaz0);
388 histograms_.fill(title + "phi_err", sigmaphi);
389 histograms_.fill(title + "theta_err", sigmatheta);
390 histograms_.fill(title + "qop_err", sigmaqop);
391 histograms_.fill(title + "p_err", sigmap);
392
393 // 2D Error plots
394 double p = std::abs(1. / trk_qop);
395 histograms_.fill(title + "d0_err_vs_p", p, sigmad0);
396 histograms_.fill(title + "z0_err_vs_p", std::abs(1. / trk_qop), sigmaz0);
397 histograms_.fill(title + "p_err_vs_p", std::abs(1. / trk_qop), sigmap);
398
399 if (track.getNhits() == 8)
400 histograms_.fill(title + "p_err_vs_p_8hits", p, sigmap);
401 else if (track.getNhits() == 9)
402 histograms_.fill(title + "p_err_vs_p_9hits", p, sigmap);
403 else if (track.getNhits() == 10)
404 histograms_.fill(title + "p_err_vs_p_10hits", p, sigmap);
405 }
406
407 if (doTruth) {
408 // Match to the truth track
409 ldmx::Track* truth_trk = nullptr;
410
411 auto it = std::find_if(truth_track_collection_->begin(),
412 truth_track_collection_->end(),
413 [&](const ldmx::Track& tt) {
414 return tt.getTrackID() == track.getTrackID();
415 });
416
417 double track_truth_prob = track.getTruthProb();
418
419 if (it != truth_track_collection_->end() &&
420 track_truth_prob >= track_prob_cut_)
421 truth_trk = &(*it);
422
423 // Found matched track
424 if (truth_trk) {
425 auto truth_d0 = truth_trk->getD0();
426 auto truth_z0 = truth_trk->getZ0();
427 auto truth_phi = truth_trk->getPhi();
428 auto truth_theta = truth_trk->getTheta();
429 auto truth_qop = truth_trk->getQoP();
430 auto truth_p = 1. / abs(truth_trk->getQoP());
431 std::vector<double> truth_mom = truth_trk->getMomentum();
432 // Polar angle
433 // The momentum in the plane transverse wrt the beam axis
434 double truth_pt_beam = std::sqrt(truth_mom[1] * truth_mom[1] +
435 truth_mom[2] * truth_mom[2]);
436
437 // histograms_.fill(title+"truth_d0", truth_d0);
438 // histograms_.fill(title+"truth_z0", truth_z0);
439 // histograms_.fill(title+"truth_phi", truth_phi);
440 // histograms_.fill(title+"truth_theta",truth_theta);
441 // histograms_.fill(title+"truth_qop", truth_qop);
442 // histograms_.fill(title+"truth_p", truth_p);
443
444 double res_d0 = trk_d0 - truth_d0;
445 double res_z0 = trk_z0 - truth_z0;
446 double res_phi = trk_phi - truth_phi;
447 double res_theta = trk_theta - truth_theta;
448 double res_qop = trk_qop - truth_qop;
449 double res_p = trk_p - truth_p;
450 double res_pt_beam = pt_beam - truth_pt_beam;
451
452 histograms_.fill(title + "res_d0", res_d0);
453 histograms_.fill(title + "res_z0", res_z0);
454 histograms_.fill(title + "res_phi", res_phi);
455 histograms_.fill(title + "res_theta", res_theta);
456 histograms_.fill(title + "res_qop", res_qop);
457 histograms_.fill(title + "res_p", res_p);
458 histograms_.fill(title + "res_pt_beam", res_pt_beam);
459
460 double pull_d0 = res_d0 / sigmad0;
461 double pull_z0 = res_z0 / sigmaz0;
462 double pull_phi = res_phi / sigmaphi;
463 double pull_theta = res_theta / sigmatheta;
464 double pull_qop = res_qop / sigmaqop;
465 double pull_p = res_p / sigmap;
466
467 histograms_.fill(title + "pull_d0", pull_d0);
468 histograms_.fill(title + "pull_z0", pull_z0);
469 histograms_.fill(title + "pull_phi", pull_phi);
470 histograms_.fill(title + "pull_theta", pull_theta);
471 histograms_.fill(title + "pull_qop", pull_qop);
472 histograms_.fill(title + "pull_p", pull_p);
473
474 // Error plots from residuals
475
476 histograms_.fill(title + "res_p_vs_p", truth_p, res_p);
477
478 histograms_.fill(title + "res_qop_vs_p", truth_p, res_qop);
479 histograms_.fill(title + "res_d0_vs_p", truth_p, res_d0);
480 histograms_.fill(title + "res_z0_vs_p", truth_p, res_z0);
481 histograms_.fill(title + "res_phi_vs_p", truth_p, res_phi);
482 histograms_.fill(title + "res_theta_vs_p", truth_p, res_theta);
483
484 histograms_.fill(title + "pull_qop_vs_p", truth_p, pull_qop);
485 histograms_.fill(title + "pull_d0_vs_p", truth_p, pull_d0);
486 histograms_.fill(title + "pull_z0_vs_p", truth_p, pull_z0);
487 histograms_.fill(title + "pull_phi_vs_p", truth_p, pull_phi);
488 histograms_.fill(title + "pull_theta_vs_p", truth_p, pull_theta);
489
490 if (track.getNhits() == 8)
491 histograms_.fill(title + "res_p_vs_p_8hits", truth_p, res_p);
492 else if (track.getNhits() == 9)
493 histograms_.fill(title + "res_p_vs_p_9hits", truth_p, res_p);
494 else if (track.getNhits() == 10)
495 histograms_.fill(title + "res_p_vs_p_10hits", truth_p, res_p);
496
497 histograms_.fill(title + "res_pt_beam_vs_p", truth_pt_beam,
498 res_pt_beam);
499
500 } // found matched track
501 } // do TruthComparison
502 } // loop on tracks
503
504} // Track Monitoring
505
506void TrackingRecoDQM::trackStateMonitoring(const ldmx::Tracks& tracks,
507 ldmx::TrackStateType ts_type,
508 const std::string& ts_title) {
509 for (auto& track : tracks) {
510 // Match the tracks to truth
511 ldmx::Track* truth_trk = nullptr;
512
513 auto it = std::find_if(truth_track_collection_->begin(),
514 truth_track_collection_->end(),
515 [&](const ldmx::Track& tt) {
516 return tt.getTrackID() == track.getTrackID();
517 });
518
519 double track_truth_prob = track.getTruthProb();
520
521 if (it != truth_track_collection_->end() &&
522 track_truth_prob >= track_prob_cut_)
523 truth_trk = &(*it);
524
525 // Match not found, skip track
526 if (!truth_trk) continue;
527
528 // TruthTrack doesn't have the right amount of states
529
530 auto trk_ts = track.getTrackState(ts_type);
531 auto truth_ts = truth_trk->getTrackState(ts_type);
532
533 if (!trk_ts.has_value()) continue;
534
535 if (!truth_ts.has_value()) continue;
536
537 ldmx::Track::TrackState& truth_target_state = truth_ts.value();
538 ldmx::Track::TrackState& target_state = trk_ts.value();
539
540 ldmx_log(debug) << "Unpacking covariance matrix";
541 Acts::BoundSquareMatrix cov =
542 tracking::sim::utils::unpackCov(target_state.cov_);
543
544 [[maybe_unused]] double sigmaloc0 = sqrt(
545 cov(Acts::BoundIndices::eBoundLoc0, Acts::BoundIndices::eBoundLoc0));
546 [[maybe_unused]] double sigmaloc1 = sqrt(
547 cov(Acts::BoundIndices::eBoundLoc1, Acts::BoundIndices::eBoundLoc1));
548 [[maybe_unused]] double sigmaphi =
549 sqrt(cov(Acts::BoundIndices::eBoundPhi, Acts::BoundIndices::eBoundPhi));
550 [[maybe_unused]] double sigmatheta = sqrt(
551 cov(Acts::BoundIndices::eBoundTheta, Acts::BoundIndices::eBoundTheta));
552 [[maybe_unused]] double sigmaqop = sqrt(cov(
553 Acts::BoundIndices::eBoundQOverP, Acts::BoundIndices::eBoundQOverP));
554
555 double trk_qop = track.getQoP();
556 double trk_p = 1. / abs(trk_qop);
557
558 double track_state_loc0 = target_state.params_[0];
559 double track_state_loc1 = target_state.params_[1];
560 [[maybe_unused]] double track_state_phi = target_state.params_[2];
561 [[maybe_unused]] double track_state_theta = target_state.params_[3];
562 [[maybe_unused]] double track_state_p = target_state.params_[4];
563
564 double truth_state_loc0 = truth_target_state.params_[0];
565 double truth_state_loc1 = truth_target_state.params_[1];
566 [[maybe_unused]] double truth_state_phi = truth_target_state.params_[2];
567 [[maybe_unused]] double truth_state_theta = truth_target_state.params_[3];
568 [[maybe_unused]] double truth_state_p = truth_target_state.params_[4];
569
570 // Check that the track state is filled
571 if (target_state.params_.size() < 5) continue;
572
573 histograms_.fill(title_ + "trk_" + ts_title + "_loc0", track_state_loc0);
574 histograms_.fill(title_ + "trk_" + ts_title + "_loc1", track_state_loc1);
575 histograms_.fill(title_ + ts_title + "_sp_hit_X", truth_state_loc0);
576 histograms_.fill(title_ + ts_title + "_sp_hit_Y", truth_state_loc1);
577
578 // TH1F The difference(residual) between end_loc0 and sp_hit_X
579 histograms_.fill(title_ + "trk_" + ts_title + "_loc0-sp_hit_X",
580 track_state_loc0 - truth_state_loc0);
581 histograms_.fill(title_ + "trk_" + ts_title + "_loc1-sp_hit_Y",
582 track_state_loc1 - truth_state_loc1);
583
584 // TH1F The pulls of loc0 and loc1
585 histograms_.fill(title_ + ts_title + "_Pulls_of_loc0",
586 (track_state_loc0 - truth_state_loc0) / sigmaloc0);
587 histograms_.fill(title_ + ts_title + "_Pulls_of_loc1",
588 (track_state_loc1 - truth_state_loc1) / sigmaloc1);
589
590 // TODO:: TH1F The pulls of phi, theta, qop
591
592 // TH2F residual vs Nhits
593 histograms_.fill(title_ + ts_title + "_res_loc0-vs-N_hits",
594 track.getNhits(), track_state_loc0 - truth_state_loc0);
595 histograms_.fill(title_ + ts_title + "_res_loc1-vs-N_hits",
596 track.getNhits(), track_state_loc1 - truth_state_loc1);
597
598 // TH2F pulls vs Nhits
599 histograms_.fill(title_ + ts_title + "_pulls_loc0-vs-N_hits",
600 track.getNhits(),
601 (track_state_loc0 - truth_state_loc0) / sigmaloc0);
602 histograms_.fill(title_ + ts_title + "_pulls_loc1-vs-N_hits",
603 track.getNhits(),
604 (track_state_loc1 - truth_state_loc1) / sigmaloc1);
605
606 // TH2F residual vs trk_p
607 histograms_.fill(title_ + ts_title + "_res_loc0-vs-trk_p", trk_p,
608 track_state_loc0 - truth_state_loc0);
609 histograms_.fill(title_ + ts_title + "_res_loc1-vs-trk_p", trk_p,
610 track_state_loc1 - truth_state_loc1);
611
612 // TH2F pulls vs trk_p
613 histograms_.fill(title_ + ts_title + "_pulls_loc0-vs-trk_p", trk_p,
614 (track_state_loc0 - truth_state_loc0) / sigmaloc0);
615 histograms_.fill(title_ + ts_title + "_pulls_loc1-vs-trk_p", trk_p,
616 (track_state_loc1 - truth_state_loc1) / sigmaloc1);
617
618 } // loop on tracks
619}
620
621void TrackingRecoDQM::sortTracks(const std::vector<ldmx::Track>& tracks,
622 std::vector<ldmx::Track>& uniqueTracks,
623 std::vector<ldmx::Track>& duplicateTracks,
624 std::vector<ldmx::Track>& fakeTracks) {
625 // Create a copy of the const vector so we can sort it
626 std::vector<ldmx::Track> sorted_tracks = tracks;
627
628 // Sort the vector of Track objects based on their trackID member
629 std::sort(sorted_tracks.begin(), sorted_tracks.end(),
630 [](ldmx::Track& t1, ldmx::Track& t2) {
631 return t1.getTrackID() < t2.getTrackID();
632 });
633
634 // Loop over the sorted vector of Track objects
635 for (size_t i = 0; i < sorted_tracks.size(); i++) {
636 if (sorted_tracks[i].getTruthProb() < track_prob_cut_)
637 fakeTracks.push_back(sorted_tracks[i]);
638 else { // not a fake track
639 // If this is the first Track object with this trackID, add it to the
640 // uniqueTracks vector directly
641 if (uniqueTracks.size() == 0 ||
642 sorted_tracks[i].getTrackID() != sorted_tracks[i - 1].getTrackID()) {
643 uniqueTracks.push_back(sorted_tracks[i]);
644 }
645 // Otherwise, add it to the duplicateTracks vector if its truthProb is
646 // lower than the existing Track object Otherwise, if the truthProbability
647 // is higher than the track stored in uniqueTracks, put it in uniqueTracks
648 // and move the uniqueTracks.back to duplicateTracks.
649 else if (sorted_tracks[i].getTruthProb() >
650 uniqueTracks.back().getTruthProb()) {
651 duplicateTracks.push_back(uniqueTracks.back());
652 uniqueTracks.back() = sorted_tracks[i];
653 }
654 // Otherwise, add it to the duplicateTracks vector
655 else {
656 duplicateTracks.push_back(sorted_tracks[i]);
657 }
658 } // a real track
659 } // loop on sorted tracks
660 // The total number of elements in the uniqueTracks and duplicateTracks
661 // vectors should be equal to the number of elements in the original tracks
662 // vector
663 if (uniqueTracks.size() + duplicateTracks.size() + fakeTracks.size() !=
664 tracks.size()) {
665 std::cerr << "Error: unique and duplicate tracks vectors do not add up to "
666 "original tracks vector";
667 return;
668 }
669
670 // Iterate through the uniqueTracks vector and duplicateTracks vector
671 ldmx_log(trace) << "Unique tracks:";
672 for (const ldmx::Track& track : uniqueTracks) {
673 ldmx_log(trace) << "\tTrack ID: " << track.getTrackID()
674 << ", Truth Prob: " << track.getTruthProb();
675 }
676 ldmx_log(trace) << "Duplicate tracks:";
677 for (const ldmx::Track& track : duplicateTracks) {
678 ldmx_log(trace) << "\tTrack ID: " << track.getTrackID()
679 << ", Truth Prob: " << track.getTruthProb();
680 }
681 ldmx_log(trace) << "Fake tracks:";
682 for (const ldmx::Track& track : fakeTracks) {
683 ldmx_log(trace) << "\tTrack ID: " << track.getTrackID()
684 << ", Truth Prob: " << track.getTruthProb();
685 }
686}
687} // namespace tracking::dqm
688
#define DECLARE_ANALYZER(CLASS)
Macro which allows the framework to construct an analyzer given its name during configuration.
HistogramPool histograms_
helper object for making and filling histograms
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:92
const std::vector< ContentType > & getCollection(const std::string &collectionName, const std::string &passName) const
Get a collection (std::vector) of objects from the event bus.
Definition Event.h:400
void fill(const std::string &name, const T &val)
Fill a 1D histogram.
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
Represents a simulated tracker hit in the simulation.
Implementation of a track object.
Definition Track.h:53
void onProcessEnd() override
Callback for the EventProcessor to take any necessary action when the processing of events finishes,...
void configure(framework::config::Parameters &parameters) override
Configure the analyzer using the given user specified parameters.
void analyze(const framework::Event &event) override
Process the event and make histograms or summaries.
void trackStateMonitoring(const ldmx::Tracks &tracks, ldmx::TrackStateType ts_type, const std::string &ts_title)
Monitoring plots for tracks extrapolated to the ECAL Scoring plane.