LDMX Software
tracking::dqm::TrackingRecoDQM Class Reference

Public Member Functions

 TrackingRecoDQM (const std::string &name, framework::Process &process)
 
 ~TrackingRecoDQM ()=default
 Destructor.
 
void analyze (const framework::Event &event) override
 Process the event and make histograms or summaries.
 
void trackMonitoring (const std::vector< ldmx::Track > &tracks, const std::vector< ldmx::Measurement > &measurements, const std::string title, const bool &doDetail, const bool &doTruth)
 
void efficiencyPlots (const std::vector< ldmx::Track > &tracks, const std::vector< ldmx::Measurement > &measurements, const std::string &title)
 
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.
 
void configure (framework::config::Parameters &parameters) override
 Configure the analyzer using the given user specified parameters.
 
void onProcessEnd () override
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
void sortTracks (const std::vector< ldmx::Track > &tracks, std::vector< ldmx::Track > &uniqueTracks, std::vector< ldmx::Track > &duplicateTracks, std::vector< ldmx::Track > &fakeTracks)
 
- Public Member Functions inherited from framework::Analyzer
 Analyzer (const std::string &name, Process &process)
 Class constructor.
 
virtual void process (Event &event) final
 Processing an event for an Analyzer is calling analyze.
 
virtual void beforeNewRun (ldmx::RunHeader &run_header) final
 Don't allow Analyzers to add parameters to the run header.
 
- Public Member Functions inherited from framework::EventProcessor
 DECLARE_FACTORY (EventProcessor, EventProcessor *, const std::string &, Process &)
 declare that we have a factory for this class
 
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()=default
 Class destructor.
 
virtual void onNewRun (const ldmx::RunHeader &run_header)
 Callback for the EventProcessor to take any necessary action when the run being processed changes.
 
virtual void onFileOpen (EventFile &event_file)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &event_file)
 Callback for the EventProcessor to take any necessary action when a event input ROOT file is closed.
 
virtual void onProcessStart ()
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
template<class T >
const T & getCondition (const std::string &condition_name)
 Access a conditions object for the current event.
 
TDirectory * getHistoDirectory ()
 Access/create a directory in the histogram file for this event processor to create histograms and analysis tuples.
 
void setStorageHint (framework::StorageControl::Hint hint)
 Mark the current event as having the given storage control hint from this module_.
 
void setStorageHint (framework::StorageControl::Hint hint, const std::string &purposeString)
 Mark the current event as having the given storage control hint from this module and the given purpose string.
 
int getLogFrequency () const
 Get the current logging frequency from the process.
 
int getRunNumber () const
 Get the run number from the process.
 
std::string getName () const
 Get the processor name.
 
void createHistograms (const std::vector< framework::config::Parameters > &histos)
 Internal function which is used to create histograms passed from the python configuration @parma histos vector of Parameters that configure histograms to create.
 

Private Attributes

std::string track_collection_
 
std::string truth_collection_
 
std::string measurement_collection_
 
std::string measurement_passname_
 
std::string ecal_sp_events_passname_
 
std::string ecal_sp_passname_
 
std::string target_sp_events_passname_
 
std::string target_sp_passname_
 
std::string track_collection_events_passname_
 
std::string track_passname_
 
std::string truth_events_passname_
 
std::string truth_passname_
 
std::string title_ {"tagger_trk_"}
 
double track_prob_cut_ {0.5}
 
std::string subdetector_ {"Tagger"}
 
bool do_truth_comparison_ {false}
 
std::vector< std::string > track_states_
 
std::shared_ptr< ldmx::Tracks > truth_track_collection_ {nullptr}
 
std::shared_ptr< std::vector< ldmx::SimTrackerHit > > ecal_scoring_hits_ {nullptr}
 
std::shared_ptr< std::vector< ldmx::SimTrackerHit > > target_scoring_hits_
 
std::vector< ldmx::Trackunique_tracks_
 
std::vector< ldmx::Trackduplicate_tracks_
 
std::vector< ldmx::Trackfake_tracks_
 
std::map< int, int > pidmap_
 

Additional Inherited Members

- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 
- Protected Attributes inherited from framework::EventProcessor
HistogramPool histograms_
 helper object for making and filling histograms
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger the_log_
 The logger for this EventProcessor.
 

Detailed Description

Definition at line 28 of file TrackingRecoDQM.h.

Constructor & Destructor Documentation

◆ TrackingRecoDQM()

tracking::dqm::TrackingRecoDQM::TrackingRecoDQM ( const std::string & name,
framework::Process & process )
inline

Definition at line 30 of file TrackingRecoDQM.h.

Base class for a module which does not produce a data product.
virtual void process(Event &event) final
Processing an event for an Analyzer is calling analyze.

Member Function Documentation

◆ analyze()

void tracking::dqm::TrackingRecoDQM::analyze ( const framework::Event & event)
overridevirtual

Process the event and make histograms or summaries.

Parameters
eventThe Event to analyze

Implements framework::Analyzer.

Definition at line 39 of file TrackingRecoDQM.cxx.

39 {
40 ldmx_log(trace) << "DQM Reading in:" << track_collection_;
41
42 if (!event.exists(track_collection_, track_collection_events_passname_)) {
43 ldmx_log(error) << "TrackCollection " << track_collection_
44 << " with pass = " << track_collection_events_passname_
45 << " not in event";
46 return;
47 }
48
49 auto tracks{
50 event.getCollection<ldmx::Track>(track_collection_, track_passname_)};
51
52 if (!event.exists(measurement_collection_, measurement_passname_)) {
53 ldmx_log(error) << "Measurement collection " << measurement_collection_
54 << " with pass = " << measurement_passname_
55 << " not in event";
56 return;
57 }
58
59 auto measurements{event.getCollection<ldmx::Measurement>(
60 measurement_collection_, measurement_passname_)};
61
62 // The truth track collection
63 if (event.exists(truth_collection_, truth_events_passname_)) {
64 truth_track_collection_ = std::make_shared<ldmx::Tracks>(
65 event.getCollection<ldmx::Track>(truth_collection_, truth_passname_));
66 do_truth_comparison_ = true;
67 }
68
69 // The scoring plane hits_
70 if (event.exists("EcalScoringPlaneHits", ecal_sp_events_passname_)) {
71 ecal_scoring_hits_ = std::make_shared<std::vector<ldmx::SimTrackerHit>>(
72 event.getCollection<ldmx::SimTrackerHit>("EcalScoringPlaneHits",
73 ecal_sp_passname_));
74 }
75
76 if (event.exists("TargetScoringPlaneHits", target_sp_events_passname_)) {
77 target_scoring_hits_ = std::make_shared<std::vector<ldmx::SimTrackerHit>>(
78 event.getCollection<ldmx::SimTrackerHit>("TargetScoringPlaneHits",
79 target_sp_passname_));
80 }
81
82 ldmx_log(debug) << "Do truth comparison::" << do_truth_comparison_;
83
84 if (do_truth_comparison_) {
85 sortTracks(tracks, unique_tracks_, duplicate_tracks_, fake_tracks_);
86 } else {
87 unique_tracks_ = tracks;
88 }
89
90 ldmx_log(debug) << "Filling histograms for " << tracks.size() << " tracks";
91
92 // General Plots
93 histograms_.fill(title_ + "N_tracks", tracks.size());
94
95 if (!unique_tracks_.empty()) {
96 ldmx_log(debug) << "Track Monitoring on " << unique_tracks_.size()
97 << " Unique Tracks";
98 trackMonitoring(unique_tracks_, measurements, title_, true,
99 do_truth_comparison_);
100 }
101
102 // Fakes and duplicates
103 if (!duplicate_tracks_.empty()) {
104 ldmx_log(debug) << "Track Monitoring on " << duplicate_tracks_.size()
105 << " duplicates";
106 trackMonitoring(duplicate_tracks_, measurements, title_ + "dup_", false,
107 false);
108 }
109 if (!fake_tracks_.empty()) {
110 ldmx_log(debug) << "Track Monitoring on " << fake_tracks_.size()
111 << " fakes";
112 trackMonitoring(fake_tracks_, measurements, title_ + "fake_", false, false);
113 }
114
115 // Track Extrapolation to Ecal Monitoring
116 ldmx_log(trace) << "Track Extrapolation to Ecal Monitoring";
117 if (std::find(track_states_.begin(), track_states_.end(), "target") !=
118 track_states_.end()) {
119 trackStateMonitoring(tracks, ldmx::TrackStateType::AtTarget, "target");
120 }
121
122 if (std::find(track_states_.begin(), track_states_.end(), "ecal") !=
123 track_states_.end()) {
124 trackStateMonitoring(tracks, ldmx::TrackStateType::AtECAL, "ecal");
125 }
126
127 if (std::find(track_states_.begin(), track_states_.end(), "beamOrigin") !=
128 track_states_.end()) {
129 trackStateMonitoring(tracks, ldmx::TrackStateType::AtBeamOrigin,
130 "beamOrigin");
131 }
132
133 // Technical Efficiency plots
134 if (do_truth_comparison_) {
135 ldmx_log(trace) << "Technical Efficiency plots";
136 efficiencyPlots(tracks, measurements, title_);
137 }
138
139 // Tagger Recoil Matching
140
141 // Clear the vectors
142 ldmx_log(trace) << "Clear the vectors";
143 unique_tracks_.clear();
144 duplicate_tracks_.clear();
145 fake_tracks_.clear();
146}
HistogramPool histograms_
helper object for making and filling histograms
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.
Represents a simulated tracker hit in the simulation.
Implementation of a track object.
Definition Track.h:53
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.

References framework::Event::exists(), framework::HistogramPool::fill(), framework::Event::getCollection(), framework::EventProcessor::histograms_, and trackStateMonitoring().

◆ configure()

void tracking::dqm::TrackingRecoDQM::configure ( framework::config::Parameters & parameters)
overridevirtual

Configure the analyzer using the given user specified parameters.

Parameters
parametersSet of parameters used to configure this analyzer.

Reimplemented from framework::EventProcessor.

Definition at line 5 of file TrackingRecoDQM.cxx.

5 {
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_events_passname_ =
13 parameters.get<std::string>("ecal_sp_events_passname");
14 ecal_sp_passname_ = parameters.get<std::string>("ecal_sp_passname");
15 target_sp_events_passname_ =
16 parameters.get<std::string>("target_sp_events_passname");
17 target_sp_passname_ = parameters.get<std::string>("target_sp_passname");
18 truth_passname_ = parameters.get<std::string>("truth_passname");
19 truth_events_passname_ = parameters.get<std::string>("truth_events_passname");
20 track_passname_ = parameters.get<std::string>("track_passname");
21 track_collection_events_passname_ =
22 parameters.get<std::string>("track_collection_events_passname");
23
24 title_ = parameters.get<std::string>("title", "tagger_trk_");
25 track_prob_cut_ = parameters.get<double>("trackProb_cut", 0.5);
26 subdetector_ = parameters.get<std::string>("subdetector", "Tagger");
27 track_states_ = parameters.get<std::vector<std::string>>("trackStates", {});
28
29 pidmap_[-321] = PIDBins::kminus;
30 pidmap_[321] = PIDBins::kplus;
31 pidmap_[-211] = PIDBins::piminus;
32 pidmap_[211] = PIDBins::piplus;
33 pidmap_[11] = PIDBins::electron;
34 pidmap_[-11] = PIDBins::positron;
35 pidmap_[2212] = PIDBins::proton;
36 pidmap_[-2212] = PIDBins::antiproton;
37}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78

References framework::config::Parameters::get().

◆ efficiencyPlots()

void tracking::dqm::TrackingRecoDQM::efficiencyPlots ( const std::vector< ldmx::Track > & tracks,
const std::vector< ldmx::Measurement > & measurements,
const std::string & title )

Definition at line 152 of file TrackingRecoDQM.cxx.

155 {
156 // Do all truth track plots - denominator
157
158 histograms_.fill(title + "truth_N_tracks", truth_track_collection_->size());
159 for (auto& truth_trk : *(truth_track_collection_)) {
160 auto truth_phi = truth_trk.getPhi();
161 auto truth_d0 = truth_trk.getD0();
162 auto truth_z0 = truth_trk.getZ0();
163 auto truth_theta = truth_trk.getTheta();
164 auto truth_qop = truth_trk.getQoP();
165 auto truth_p = 1. / abs(truth_trk.getQoP());
166 auto truth_n_hits = truth_trk.getNhits();
167
168 std::vector<double> truth_mom = truth_trk.getMomentum();
169
170 // Polar angle
171 // The momentum in the plane transverse wrt the beam axis
172 auto truth_pt_beam =
173 std::sqrt(truth_mom[1] * truth_mom[1] + truth_mom[2] * truth_mom[2]);
174
175 auto truth_beam_angle = std::atan2(truth_pt_beam, truth_mom[0]);
176
177 histograms_.fill(title + "truth_nHits", truth_n_hits);
178 histograms_.fill(title + "truth_d0", truth_d0);
179 histograms_.fill(title + "truth_z0", truth_z0);
180 histograms_.fill(title + "truth_phi", truth_phi);
181 histograms_.fill(title + "truth_theta", truth_theta);
182 histograms_.fill(title + "truth_qop", truth_qop);
183 histograms_.fill(title + "truth_p", truth_p);
184 histograms_.fill(title + "truth_beam_angle", truth_beam_angle);
185
186 if (pidmap_.count(truth_trk.getPdgID()) != 0) {
187 histograms_.fill(title + "truth_PID", pidmap_[truth_trk.getPdgID()]);
188
189 // TODO do this properly.
190
191 if (pidmap_[truth_trk.getPdgID()] == PIDBins::kminus) {
192 histograms_.fill(title + "truth_kminus_p", truth_p);
193 }
194
195 if (pidmap_[truth_trk.getPdgID()] == PIDBins::kplus) {
196 histograms_.fill(title + "truth_kplus_p", truth_p);
197 }
198
199 if (pidmap_[truth_trk.getPdgID()] == PIDBins::piminus) {
200 histograms_.fill(title + "truth_piminus_p", truth_p);
201 }
202
203 if (pidmap_[truth_trk.getPdgID()] == PIDBins::piplus) {
204 histograms_.fill(title + "truth_piplus_p", truth_p);
205 }
206
207 if (pidmap_[truth_trk.getPdgID()] == PIDBins::electron) {
208 histograms_.fill(title + "truth_electron_p", truth_p);
209 }
210
211 if (pidmap_[truth_trk.getPdgID()] == PIDBins::positron) {
212 histograms_.fill(title + "truth_positron_p", truth_p);
213 }
214
215 if (pidmap_[truth_trk.getPdgID()] == PIDBins::proton) {
216 histograms_.fill(title + "truth_proton_p", truth_p);
217 }
218 }
219
220 } // loop on truth tracks
221
222 for (auto& track : tracks) {
223 // Match the tracks to truth
224 ldmx::Track* truth_trk = nullptr;
225
226 auto it = std::find_if(truth_track_collection_->begin(),
227 truth_track_collection_->end(),
228 [&](const ldmx::Track& tt) {
229 return tt.getTrackID() == track.getTrackID();
230 });
231
232 double track_truth_prob = track.getTruthProb();
233
234 if (it != truth_track_collection_->end() &&
235 track_truth_prob >= track_prob_cut_)
236 truth_trk = &(*it);
237
238 // Match not found
239 if (!truth_trk) return;
240
241 auto truth_phi = truth_trk->getPhi();
242 auto truth_d0 = truth_trk->getD0();
243 auto truth_z0 = truth_trk->getZ0();
244 auto truth_theta = truth_trk->getTheta();
245 auto truth_qop = truth_trk->getQoP();
246 auto truth_p = 1. / abs(truth_trk->getQoP());
247 std::vector<double> truth_mom = truth_trk->getMomentum();
248
249 // Polar angle
250 // The momentum in the plane transverse wrt the beam axis
251 auto truth_pt_beam =
252 std::sqrt(truth_mom[1] * truth_mom[1] + truth_mom[2] * truth_mom[2]);
253
254 auto truth_beam_angle = std::atan2(truth_pt_beam, truth_mom[0]);
255
256 // Fill reco plots for efficiencies - numerator. The quantities are truth
257 histograms_.fill(title + "match_prob", track_truth_prob);
258 histograms_.fill(title + "match_d0", truth_d0);
259 histograms_.fill(title + "match_z0", truth_z0);
260 histograms_.fill(title + "match_phi", truth_phi);
261 histograms_.fill(title + "match_theta", truth_theta);
262 histograms_.fill(title + "match_p", truth_p);
263 histograms_.fill(title + "match_qop", truth_qop);
264 histograms_.fill(title + "match_beam_angle", truth_beam_angle);
265 histograms_.fill(title + "match_nHits", measurements.size());
266 auto dedx_measurements = track.getDedxMeasurements();
267 auto measurement_idxs = track.getMeasurementsIdxs();
268 for (size_t i = 0; i < measurement_idxs.size(); ++i) {
269 histograms_.fill(title + "match_layers_hit",
270 measurements.at(measurement_idxs[i]).getLayer());
271 // Add histogram of the measurement dE/dx
272 if (i < dedx_measurements.size()) {
273 histograms_.fill(title + "match_measurement_dedx",
274 dedx_measurements[i]);
275 }
276 }
277
278 // For some particles
279
280 if (pidmap_.count(truth_trk->getPdgID()) != 0) {
281 histograms_.fill(title + "match_PID", pidmap_[truth_trk->getPdgID()]);
282
283 // TODO do this properly.
284
285 if (pidmap_[truth_trk->getPdgID()] == PIDBins::kminus) {
286 histograms_.fill(title + "match_kminus_p", truth_p);
287 }
288
289 if (pidmap_[truth_trk->getPdgID()] == PIDBins::kplus) {
290 histograms_.fill(title + "match_kplus_p", truth_p);
291 }
292
293 if (pidmap_[truth_trk->getPdgID()] == PIDBins::piminus) {
294 histograms_.fill(title + "match_piminus_p", truth_p);
295 }
296
297 if (pidmap_[truth_trk->getPdgID()] == PIDBins::piplus) {
298 histograms_.fill(title + "match_piplus_p", truth_p);
299 }
300
301 if (pidmap_[truth_trk->getPdgID()] == PIDBins::electron) {
302 histograms_.fill(title + "match_electron_p", truth_p);
303 }
304
305 if (pidmap_[truth_trk->getPdgID()] == PIDBins::positron) {
306 histograms_.fill(title + "match_positron_p", truth_p);
307 }
308
309 if (pidmap_[truth_trk->getPdgID()] == PIDBins::proton) {
310 histograms_.fill(title + "match_proton_p", truth_p);
311 }
312 }
313 } // Loop on tracks
314
315} // Efficiency plots

◆ onProcessEnd()

void tracking::dqm::TrackingRecoDQM::onProcessEnd ( )
overridevirtual

Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.

Reimplemented from framework::EventProcessor.

Definition at line 148 of file TrackingRecoDQM.cxx.

148 {
149 // Produce the efficiency plots. (TODO::Switch to TEfficiency instead)
150}

◆ sortTracks()

void tracking::dqm::TrackingRecoDQM::sortTracks ( const std::vector< ldmx::Track > & tracks,
std::vector< ldmx::Track > & uniqueTracks,
std::vector< ldmx::Track > & duplicateTracks,
std::vector< ldmx::Track > & fakeTracks )

Definition at line 623 of file TrackingRecoDQM.cxx.

626 {
627 // Create a copy of the const vector so we can sort it
628 std::vector<ldmx::Track> sorted_tracks = tracks;
629
630 // Sort the vector of Track objects based on their trackID member
631 std::sort(sorted_tracks.begin(), sorted_tracks.end(),
632 [](ldmx::Track& t1, ldmx::Track& t2) {
633 return t1.getTrackID() < t2.getTrackID();
634 });
635
636 // Loop over the sorted vector of Track objects
637 for (size_t i = 0; i < sorted_tracks.size(); i++) {
638 if (sorted_tracks[i].getTruthProb() < track_prob_cut_)
639 fakeTracks.push_back(sorted_tracks[i]);
640 else { // not a fake track
641 // If this is the first Track object with this trackID, add it to the
642 // uniqueTracks vector directly
643 if (uniqueTracks.size() == 0 ||
644 sorted_tracks[i].getTrackID() != sorted_tracks[i - 1].getTrackID()) {
645 uniqueTracks.push_back(sorted_tracks[i]);
646 }
647 // Otherwise, add it to the duplicateTracks vector if its truthProb is
648 // lower than the existing Track object Otherwise, if the truthProbability
649 // is higher than the track stored in uniqueTracks, put it in uniqueTracks
650 // and move the uniqueTracks.back to duplicateTracks.
651 else if (sorted_tracks[i].getTruthProb() >
652 uniqueTracks.back().getTruthProb()) {
653 duplicateTracks.push_back(uniqueTracks.back());
654 uniqueTracks.back() = sorted_tracks[i];
655 }
656 // Otherwise, add it to the duplicateTracks vector
657 else {
658 duplicateTracks.push_back(sorted_tracks[i]);
659 }
660 } // a real track
661 } // loop on sorted tracks
662 // The total number of elements in the uniqueTracks and duplicateTracks
663 // vectors should be equal to the number of elements in the original tracks
664 // vector
665 if (uniqueTracks.size() + duplicateTracks.size() + fakeTracks.size() !=
666 tracks.size()) {
667 std::cerr << "Error: unique and duplicate tracks vectors do not add up to "
668 "original tracks vector";
669 return;
670 }
671
672 // Iterate through the uniqueTracks vector and duplicateTracks vector
673 ldmx_log(trace) << "Unique tracks:";
674 for (const ldmx::Track& track : uniqueTracks) {
675 ldmx_log(trace) << "\tTrack ID: " << track.getTrackID()
676 << ", Truth Prob: " << track.getTruthProb();
677 }
678 ldmx_log(trace) << "Duplicate tracks:";
679 for (const ldmx::Track& track : duplicateTracks) {
680 ldmx_log(trace) << "\tTrack ID: " << track.getTrackID()
681 << ", Truth Prob: " << track.getTruthProb();
682 }
683 ldmx_log(trace) << "Fake tracks:";
684 for (const ldmx::Track& track : fakeTracks) {
685 ldmx_log(trace) << "\tTrack ID: " << track.getTrackID()
686 << ", Truth Prob: " << track.getTruthProb();
687 }
688}

◆ trackMonitoring()

void tracking::dqm::TrackingRecoDQM::trackMonitoring ( const std::vector< ldmx::Track > & tracks,
const std::vector< ldmx::Measurement > & measurements,
const std::string title,
const bool & doDetail,
const bool & doTruth )

Definition at line 317 of file TrackingRecoDQM.cxx.

320 {
321 for (auto& track : tracks) {
322 // Perigee track parameters
323 auto trk_d0 = track.getD0();
324 auto trk_z0 = track.getZ0();
325 auto trk_qop = track.getQoP();
326 auto trk_theta = track.getTheta();
327 auto trk_phi = track.getPhi();
328 auto trk_p = 1. / abs(trk_qop);
329 auto dedx_measurements = track.getDedxMeasurements();
330 auto measurement_idxs = track.getMeasurementsIdxs();
331 for (size_t i = 0; i < measurement_idxs.size(); ++i) {
332 histograms_.fill(title + "layers_hit",
333 measurements.at(measurement_idxs[i]).getLayer());
334 // Add histogram of the measurement dE/dx
335 if (i < dedx_measurements.size()) {
336 histograms_.fill(title + "measurement_dedx", dedx_measurements[i]);
337 }
338 }
339
340 std::vector<double> trk_mom = track.getMomentum();
341
342 // The transverse momentum in the bending plane
343 double pt_bending =
344 std::sqrt(trk_mom[0] * trk_mom[0] + trk_mom[1] * trk_mom[1]);
345
346 // The momentum in the plane transverse wrt the beam axis
347 double pt_beam =
348 std::sqrt(trk_mom[1] * trk_mom[1] + trk_mom[2] * trk_mom[2]);
349
350 // Covariance matrix
351 Acts::BoundSquareMatrix cov =
352 tracking::sim::utils::unpackCov(track.getPerigeeCov());
353
354 double sigmad0 = sqrt(
355 cov(Acts::BoundIndices::eBoundLoc0, Acts::BoundIndices::eBoundLoc0));
356 double sigmaz0 = sqrt(
357 cov(Acts::BoundIndices::eBoundLoc1, Acts::BoundIndices::eBoundLoc1));
358 double sigmaphi =
359 sqrt(cov(Acts::BoundIndices::eBoundPhi, Acts::BoundIndices::eBoundPhi));
360 double sigmatheta = sqrt(
361 cov(Acts::BoundIndices::eBoundTheta, Acts::BoundIndices::eBoundTheta));
362 double sigmaqop = sqrt(cov(Acts::BoundIndices::eBoundQOverP,
363 Acts::BoundIndices::eBoundQOverP));
364 double sigmap = (1. / trk_qop) * (1. / trk_qop) * sigmaqop;
365
366 histograms_.fill(title + "d0", trk_d0);
367 histograms_.fill(title + "z0", trk_z0);
368 histograms_.fill(title + "qop", trk_qop);
369 histograms_.fill(title + "phi", trk_phi);
370 histograms_.fill(title + "theta", trk_theta);
371 histograms_.fill(title + "p", std::abs(1. / trk_qop));
372
373 if (doDetail) {
374 histograms_.fill(title + "px", trk_mom[0]);
375 histograms_.fill(title + "py", trk_mom[1]);
376 histograms_.fill(title + "pz", trk_mom[2]);
377
378 histograms_.fill(title + "pt_bending", pt_bending);
379 histograms_.fill(title + "pt_beam", pt_beam);
380
381 histograms_.fill(title + "nHits", track.getNhits());
382 histograms_.fill(title + "Chi2", track.getChi2());
383 histograms_.fill(title + "ndf", track.getNdf());
384 histograms_.fill(title + "Chi2_per_ndf",
385 track.getChi2() / track.getNdf());
386 histograms_.fill(title + "nShared", track.getNsharedHits());
387
388 histograms_.fill(title + "d0_err", sigmad0);
389 histograms_.fill(title + "z0_err", sigmaz0);
390 histograms_.fill(title + "phi_err", sigmaphi);
391 histograms_.fill(title + "theta_err", sigmatheta);
392 histograms_.fill(title + "qop_err", sigmaqop);
393 histograms_.fill(title + "p_err", sigmap);
394
395 // 2D Error plots
396 double p = std::abs(1. / trk_qop);
397 histograms_.fill(title + "d0_err_vs_p", p, sigmad0);
398 histograms_.fill(title + "z0_err_vs_p", std::abs(1. / trk_qop), sigmaz0);
399 histograms_.fill(title + "p_err_vs_p", std::abs(1. / trk_qop), sigmap);
400
401 if (track.getNhits() == 8)
402 histograms_.fill(title + "p_err_vs_p_8hits", p, sigmap);
403 else if (track.getNhits() == 9)
404 histograms_.fill(title + "p_err_vs_p_9hits", p, sigmap);
405 else if (track.getNhits() == 10)
406 histograms_.fill(title + "p_err_vs_p_10hits", p, sigmap);
407 }
408
409 if (doTruth) {
410 // Match to the truth track
411 ldmx::Track* truth_trk = nullptr;
412
413 auto it = std::find_if(truth_track_collection_->begin(),
414 truth_track_collection_->end(),
415 [&](const ldmx::Track& tt) {
416 return tt.getTrackID() == track.getTrackID();
417 });
418
419 double track_truth_prob = track.getTruthProb();
420
421 if (it != truth_track_collection_->end() &&
422 track_truth_prob >= track_prob_cut_)
423 truth_trk = &(*it);
424
425 // Found matched track
426 if (truth_trk) {
427 auto truth_d0 = truth_trk->getD0();
428 auto truth_z0 = truth_trk->getZ0();
429 auto truth_phi = truth_trk->getPhi();
430 auto truth_theta = truth_trk->getTheta();
431 auto truth_qop = truth_trk->getQoP();
432 auto truth_p = 1. / abs(truth_trk->getQoP());
433 std::vector<double> truth_mom = truth_trk->getMomentum();
434 // Polar angle
435 // The momentum in the plane transverse wrt the beam axis
436 double truth_pt_beam = std::sqrt(truth_mom[1] * truth_mom[1] +
437 truth_mom[2] * truth_mom[2]);
438
439 // histograms_.fill(title+"truth_d0", truth_d0);
440 // histograms_.fill(title+"truth_z0", truth_z0);
441 // histograms_.fill(title+"truth_phi", truth_phi);
442 // histograms_.fill(title+"truth_theta",truth_theta);
443 // histograms_.fill(title+"truth_qop", truth_qop);
444 // histograms_.fill(title+"truth_p", truth_p);
445
446 double res_d0 = trk_d0 - truth_d0;
447 double res_z0 = trk_z0 - truth_z0;
448 double res_phi = trk_phi - truth_phi;
449 double res_theta = trk_theta - truth_theta;
450 double res_qop = trk_qop - truth_qop;
451 double res_p = trk_p - truth_p;
452 double res_pt_beam = pt_beam - truth_pt_beam;
453
454 histograms_.fill(title + "res_d0", res_d0);
455 histograms_.fill(title + "res_z0", res_z0);
456 histograms_.fill(title + "res_phi", res_phi);
457 histograms_.fill(title + "res_theta", res_theta);
458 histograms_.fill(title + "res_qop", res_qop);
459 histograms_.fill(title + "res_p", res_p);
460 histograms_.fill(title + "res_pt_beam", res_pt_beam);
461
462 double pull_d0 = res_d0 / sigmad0;
463 double pull_z0 = res_z0 / sigmaz0;
464 double pull_phi = res_phi / sigmaphi;
465 double pull_theta = res_theta / sigmatheta;
466 double pull_qop = res_qop / sigmaqop;
467 double pull_p = res_p / sigmap;
468
469 histograms_.fill(title + "pull_d0", pull_d0);
470 histograms_.fill(title + "pull_z0", pull_z0);
471 histograms_.fill(title + "pull_phi", pull_phi);
472 histograms_.fill(title + "pull_theta", pull_theta);
473 histograms_.fill(title + "pull_qop", pull_qop);
474 histograms_.fill(title + "pull_p", pull_p);
475
476 // Error plots from residuals
477
478 histograms_.fill(title + "res_p_vs_p", truth_p, res_p);
479
480 histograms_.fill(title + "res_qop_vs_p", truth_p, res_qop);
481 histograms_.fill(title + "res_d0_vs_p", truth_p, res_d0);
482 histograms_.fill(title + "res_z0_vs_p", truth_p, res_z0);
483 histograms_.fill(title + "res_phi_vs_p", truth_p, res_phi);
484 histograms_.fill(title + "res_theta_vs_p", truth_p, res_theta);
485
486 histograms_.fill(title + "pull_qop_vs_p", truth_p, pull_qop);
487 histograms_.fill(title + "pull_d0_vs_p", truth_p, pull_d0);
488 histograms_.fill(title + "pull_z0_vs_p", truth_p, pull_z0);
489 histograms_.fill(title + "pull_phi_vs_p", truth_p, pull_phi);
490 histograms_.fill(title + "pull_theta_vs_p", truth_p, pull_theta);
491
492 if (track.getNhits() == 8)
493 histograms_.fill(title + "res_p_vs_p_8hits", truth_p, res_p);
494 else if (track.getNhits() == 9)
495 histograms_.fill(title + "res_p_vs_p_9hits", truth_p, res_p);
496 else if (track.getNhits() == 10)
497 histograms_.fill(title + "res_p_vs_p_10hits", truth_p, res_p);
498
499 histograms_.fill(title + "res_pt_beam_vs_p", truth_pt_beam,
500 res_pt_beam);
501
502 } // found matched track
503 } // do TruthComparison
504 } // loop on tracks
505
506} // Track Monitoring

◆ trackStateMonitoring()

void tracking::dqm::TrackingRecoDQM::trackStateMonitoring ( const ldmx::Tracks & tracks,
ldmx::TrackStateType ts_type,
const std::string & ts_title )

Monitoring plots for tracks extrapolated to the ECAL Scoring plane.

This aims Tracks will be truth matched first to get the trackID. The hit with the trackID

Definition at line 508 of file TrackingRecoDQM.cxx.

510 {
511 for (auto& track : tracks) {
512 // Match the tracks to truth
513 ldmx::Track* truth_trk = nullptr;
514
515 auto it = std::find_if(truth_track_collection_->begin(),
516 truth_track_collection_->end(),
517 [&](const ldmx::Track& tt) {
518 return tt.getTrackID() == track.getTrackID();
519 });
520
521 double track_truth_prob = track.getTruthProb();
522
523 if (it != truth_track_collection_->end() &&
524 track_truth_prob >= track_prob_cut_)
525 truth_trk = &(*it);
526
527 // Match not found, skip track
528 if (!truth_trk) continue;
529
530 // TruthTrack doesn't have the right amount of states
531
532 auto trk_ts = track.getTrackState(ts_type);
533 auto truth_ts = truth_trk->getTrackState(ts_type);
534
535 if (!trk_ts.has_value()) continue;
536
537 if (!truth_ts.has_value()) continue;
538
539 ldmx::Track::TrackState& truth_target_state = truth_ts.value();
540 ldmx::Track::TrackState& target_state = trk_ts.value();
541
542 ldmx_log(debug) << "Unpacking covariance matrix";
543 Acts::BoundSquareMatrix cov =
544 tracking::sim::utils::unpackCov(target_state.cov_);
545
546 [[maybe_unused]] double sigmaloc0 = sqrt(
547 cov(Acts::BoundIndices::eBoundLoc0, Acts::BoundIndices::eBoundLoc0));
548 [[maybe_unused]] double sigmaloc1 = sqrt(
549 cov(Acts::BoundIndices::eBoundLoc1, Acts::BoundIndices::eBoundLoc1));
550 [[maybe_unused]] double sigmaphi =
551 sqrt(cov(Acts::BoundIndices::eBoundPhi, Acts::BoundIndices::eBoundPhi));
552 [[maybe_unused]] double sigmatheta = sqrt(
553 cov(Acts::BoundIndices::eBoundTheta, Acts::BoundIndices::eBoundTheta));
554 [[maybe_unused]] double sigmaqop = sqrt(cov(
555 Acts::BoundIndices::eBoundQOverP, Acts::BoundIndices::eBoundQOverP));
556
557 double trk_qop = track.getQoP();
558 double trk_p = 1. / abs(trk_qop);
559
560 double track_state_loc0 = target_state.params_[0];
561 double track_state_loc1 = target_state.params_[1];
562 [[maybe_unused]] double track_state_phi = target_state.params_[2];
563 [[maybe_unused]] double track_state_theta = target_state.params_[3];
564 [[maybe_unused]] double track_state_p = target_state.params_[4];
565
566 double truth_state_loc0 = truth_target_state.params_[0];
567 double truth_state_loc1 = truth_target_state.params_[1];
568 [[maybe_unused]] double truth_state_phi = truth_target_state.params_[2];
569 [[maybe_unused]] double truth_state_theta = truth_target_state.params_[3];
570 [[maybe_unused]] double truth_state_p = truth_target_state.params_[4];
571
572 // Check that the track state is filled
573 if (target_state.params_.size() < 5) continue;
574
575 histograms_.fill(title_ + "trk_" + ts_title + "_loc0", track_state_loc0);
576 histograms_.fill(title_ + "trk_" + ts_title + "_loc1", track_state_loc1);
577 histograms_.fill(title_ + ts_title + "_sp_hit_X", truth_state_loc0);
578 histograms_.fill(title_ + ts_title + "_sp_hit_Y", truth_state_loc1);
579
580 // TH1F The difference(residual) between end_loc0 and sp_hit_X
581 histograms_.fill(title_ + "trk_" + ts_title + "_loc0-sp_hit_X",
582 track_state_loc0 - truth_state_loc0);
583 histograms_.fill(title_ + "trk_" + ts_title + "_loc1-sp_hit_Y",
584 track_state_loc1 - truth_state_loc1);
585
586 // TH1F The pulls of loc0 and loc1
587 histograms_.fill(title_ + ts_title + "_Pulls_of_loc0",
588 (track_state_loc0 - truth_state_loc0) / sigmaloc0);
589 histograms_.fill(title_ + ts_title + "_Pulls_of_loc1",
590 (track_state_loc1 - truth_state_loc1) / sigmaloc1);
591
592 // TODO:: TH1F The pulls of phi, theta, qop
593
594 // TH2F residual vs Nhits
595 histograms_.fill(title_ + ts_title + "_res_loc0-vs-N_hits",
596 track.getNhits(), track_state_loc0 - truth_state_loc0);
597 histograms_.fill(title_ + ts_title + "_res_loc1-vs-N_hits",
598 track.getNhits(), track_state_loc1 - truth_state_loc1);
599
600 // TH2F pulls vs Nhits
601 histograms_.fill(title_ + ts_title + "_pulls_loc0-vs-N_hits",
602 track.getNhits(),
603 (track_state_loc0 - truth_state_loc0) / sigmaloc0);
604 histograms_.fill(title_ + ts_title + "_pulls_loc1-vs-N_hits",
605 track.getNhits(),
606 (track_state_loc1 - truth_state_loc1) / sigmaloc1);
607
608 // TH2F residual vs trk_p
609 histograms_.fill(title_ + ts_title + "_res_loc0-vs-trk_p", trk_p,
610 track_state_loc0 - truth_state_loc0);
611 histograms_.fill(title_ + ts_title + "_res_loc1-vs-trk_p", trk_p,
612 track_state_loc1 - truth_state_loc1);
613
614 // TH2F pulls vs trk_p
615 histograms_.fill(title_ + ts_title + "_pulls_loc0-vs-trk_p", trk_p,
616 (track_state_loc0 - truth_state_loc0) / sigmaloc0);
617 histograms_.fill(title_ + ts_title + "_pulls_loc1-vs-trk_p", trk_p,
618 (track_state_loc1 - truth_state_loc1) / sigmaloc1);
619
620 } // loop on tracks
621}

References framework::HistogramPool::fill(), and framework::EventProcessor::histograms_.

Referenced by analyze().

Member Data Documentation

◆ do_truth_comparison_

bool tracking::dqm::TrackingRecoDQM::do_truth_comparison_ {false}
private

Definition at line 90 of file TrackingRecoDQM.h.

90{false};

◆ duplicate_tracks_

std::vector<ldmx::Track> tracking::dqm::TrackingRecoDQM::duplicate_tracks_
private

Definition at line 108 of file TrackingRecoDQM.h.

◆ ecal_scoring_hits_

std::shared_ptr<std::vector<ldmx::SimTrackerHit> > tracking::dqm::TrackingRecoDQM::ecal_scoring_hits_ {nullptr}
private

Definition at line 97 of file TrackingRecoDQM.h.

97{nullptr};

◆ ecal_sp_events_passname_

std::string tracking::dqm::TrackingRecoDQM::ecal_sp_events_passname_
private

Definition at line 78 of file TrackingRecoDQM.h.

◆ ecal_sp_passname_

std::string tracking::dqm::TrackingRecoDQM::ecal_sp_passname_
private

Definition at line 79 of file TrackingRecoDQM.h.

◆ fake_tracks_

std::vector<ldmx::Track> tracking::dqm::TrackingRecoDQM::fake_tracks_
private

Definition at line 110 of file TrackingRecoDQM.h.

◆ measurement_collection_

std::string tracking::dqm::TrackingRecoDQM::measurement_collection_
private

Definition at line 75 of file TrackingRecoDQM.h.

◆ measurement_passname_

std::string tracking::dqm::TrackingRecoDQM::measurement_passname_
private

Definition at line 76 of file TrackingRecoDQM.h.

◆ pidmap_

std::map<int, int> tracking::dqm::TrackingRecoDQM::pidmap_
private

Definition at line 113 of file TrackingRecoDQM.h.

◆ subdetector_

std::string tracking::dqm::TrackingRecoDQM::subdetector_ {"Tagger"}
private

Definition at line 89 of file TrackingRecoDQM.h.

89{"Tagger"};

◆ target_scoring_hits_

std::shared_ptr<std::vector<ldmx::SimTrackerHit> > tracking::dqm::TrackingRecoDQM::target_scoring_hits_
private
Initial value:
{
nullptr}

Definition at line 100 of file TrackingRecoDQM.h.

100 {
101 nullptr};

◆ target_sp_events_passname_

std::string tracking::dqm::TrackingRecoDQM::target_sp_events_passname_
private

Definition at line 80 of file TrackingRecoDQM.h.

◆ target_sp_passname_

std::string tracking::dqm::TrackingRecoDQM::target_sp_passname_
private

Definition at line 81 of file TrackingRecoDQM.h.

◆ title_

std::string tracking::dqm::TrackingRecoDQM::title_ {"tagger_trk_"}
private

Definition at line 87 of file TrackingRecoDQM.h.

87{"tagger_trk_"};

◆ track_collection_

std::string tracking::dqm::TrackingRecoDQM::track_collection_
private

Definition at line 73 of file TrackingRecoDQM.h.

◆ track_collection_events_passname_

std::string tracking::dqm::TrackingRecoDQM::track_collection_events_passname_
private

Definition at line 82 of file TrackingRecoDQM.h.

◆ track_passname_

std::string tracking::dqm::TrackingRecoDQM::track_passname_
private

Definition at line 83 of file TrackingRecoDQM.h.

◆ track_prob_cut_

double tracking::dqm::TrackingRecoDQM::track_prob_cut_ {0.5}
private

Definition at line 88 of file TrackingRecoDQM.h.

88{0.5};

◆ track_states_

std::vector<std::string> tracking::dqm::TrackingRecoDQM::track_states_
private

Definition at line 91 of file TrackingRecoDQM.h.

◆ truth_collection_

std::string tracking::dqm::TrackingRecoDQM::truth_collection_
private

Definition at line 74 of file TrackingRecoDQM.h.

◆ truth_events_passname_

std::string tracking::dqm::TrackingRecoDQM::truth_events_passname_
private

Definition at line 84 of file TrackingRecoDQM.h.

◆ truth_passname_

std::string tracking::dqm::TrackingRecoDQM::truth_passname_
private

Definition at line 85 of file TrackingRecoDQM.h.

◆ truth_track_collection_

std::shared_ptr<ldmx::Tracks> tracking::dqm::TrackingRecoDQM::truth_track_collection_ {nullptr}
private

Definition at line 94 of file TrackingRecoDQM.h.

94{nullptr};

◆ unique_tracks_

std::vector<ldmx::Track> tracking::dqm::TrackingRecoDQM::unique_tracks_
private

Definition at line 106 of file TrackingRecoDQM.h.


The documentation for this class was generated from the following files: