LDMX Software
dqm::EcalSPTrackCompare Class Reference

Compare the primary electron at the ECal scoring plane with ECAL tracks reconstructed by ACTS CKF. More...

#include <EcalSPTrackCompare.h>

Public Member Functions

 EcalSPTrackCompare (const std::string &name, framework::Process &process)
 
void configure (framework::config::Parameters &ps) override
 Callback for the EventProcessor to configure itself from the given set of parameters.
 
void analyze (const framework::Event &event) override
 Process the event and make histograms or summaries.
 
- 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.
 
virtual void onProcessEnd ()
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
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 ecal_sp_coll_name_ {"EcalScoringPlaneHits"}
 
std::string ecal_sp_pass_name_ {""}
 
std::string track_collection_ {"EcalTracks"}
 
std::string track_pass_name_ {""}
 

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

Compare the primary electron at the ECal scoring plane with ECAL tracks reconstructed by ACTS CKF.

Finds the primary beam electron (track ID 1, PDG 11) at scoring plane 31 (ECal front face), then compares its position and direction with each reconstructed EcalTrack. Fills residual histograms for the closest matching track.

Definition at line 21 of file EcalSPTrackCompare.h.

Constructor & Destructor Documentation

◆ EcalSPTrackCompare()

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

Definition at line 23 of file EcalSPTrackCompare.h.

24 : Analyzer(name, process) {}
virtual void process(Event &event) final
Processing an event for an Analyzer is calling analyze.
Analyzer(const std::string &name, Process &process)
Class constructor.

Member Function Documentation

◆ analyze()

void dqm::EcalSPTrackCompare::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 18 of file EcalSPTrackCompare.cxx.

18 {
19 // --- Find the primary electron at ECal scoring plane (plane 31) ---
20 if (!event.exists(ecal_sp_coll_name_, ecal_sp_pass_name_)) {
21 ldmx_log(debug) << "SP collection " << ecal_sp_coll_name_ << " not found";
22 return;
23 }
24
25 const auto& sp_hits = event.getCollection<ldmx::SimTrackerHit>(
26 ecal_sp_coll_name_, ecal_sp_pass_name_);
27
28 ldmx_log(debug) << "SP hits: " << sp_hits.size();
29
30 // Find the highest-momentum forward-going electron at plane 31.
31 // In dark brem events the recoil electron has trackID != 1.
32 const ldmx::SimTrackerHit* sp_electron = nullptr;
33 double best_pz = 0;
34 for (const auto& hit : sp_hits) {
35 ldmx::SimSpecialID hit_id(hit.getID());
36 if (hit_id.plane() != 31) continue;
37 if (hit.getPdgID() != 11) continue;
38 const auto p = hit.getMomentum();
39 if (p[2] <= 0) continue;
40 if (p[2] > best_pz) {
41 best_pz = p[2];
42 sp_electron = &hit;
43 }
44 }
45
46 if (!sp_electron) {
47 ldmx_log(debug) << "No forward electron found at ECal SP";
48 return;
49 }
50
51 // Extract SP electron kinematics
52 const auto sp_pos = sp_electron->getPosition();
53 const auto sp_mom = sp_electron->getMomentum();
54 double sp_x = sp_pos[0];
55 double sp_y = sp_pos[1];
56
57 ldmx_log(debug) << "SP electron: trackID=" << sp_electron->getTrackID()
58 << " pos=(" << sp_x << ", " << sp_y << ") pz=" << best_pz;
59 double sp_px = sp_mom[0];
60 double sp_py = sp_mom[1];
61 double sp_pz = sp_mom[2];
62 double sp_p = std::sqrt(sp_px * sp_px + sp_py * sp_py + sp_pz * sp_pz);
63 double sp_pt = std::sqrt(sp_px * sp_px + sp_py * sp_py);
64 double sp_theta = std::atan2(sp_pt, sp_pz);
65 double sp_phi = std::atan2(sp_py, sp_px);
66 if (sp_phi < 0) sp_phi += 2.0 * M_PI; // Shift to [0, 2π]
67
68 histograms_.fill("sp_x", sp_x);
69 histograms_.fill("sp_y", sp_y);
70 histograms_.fill("sp_p", sp_p);
71 histograms_.fill("sp_pt", sp_pt);
72 histograms_.fill("sp_theta", sp_theta);
73 histograms_.fill("sp_phi", sp_phi);
74
75 // --- Get reconstructed tracks ---
76 if (!event.exists(track_collection_, track_pass_name_)) {
77 histograms_.fill("n_tracks", 0);
78 histograms_.fill("has_match", 0);
79 return;
80 }
81
82 const auto& tracks =
83 event.getCollection<ldmx::Track>(track_collection_, track_pass_name_);
84
85 int n_tracks = tracks.size();
86 histograms_.fill("n_tracks", n_tracks);
87
88 if (n_tracks == 0) {
89 histograms_.fill("has_match", 0);
90 return;
91 }
92
93 // --- Find closest track to SP electron ---
94 int best_idx = -1;
95 double best_dr = 1e9;
96
97 for (int i = 0; i < n_tracks; ++i) {
98 const auto& track = tracks[i];
99 double trk_theta = track.getTheta();
100 double trk_phi = track.getPhi();
101
102 // Angular distance
103 double dtheta = trk_theta - sp_theta;
104 double dphi = trk_phi - sp_phi;
105 // Wrap dphi to [-pi, pi]
106 while (dphi > M_PI) dphi -= 2.0 * M_PI;
107 while (dphi < -M_PI) dphi += 2.0 * M_PI;
108
109 double dr = std::sqrt(dtheta * dtheta + dphi * dphi);
110 if (dr < best_dr) {
111 best_dr = dr;
112 best_idx = i;
113 }
114 }
115
116 histograms_.fill("has_match", 1);
117
118 // --- Fill residuals for the closest track ---
119 const auto& best = tracks[best_idx];
120
121 // Get track parameters at perigee
122 double trk_d0 = best.getD0();
123 double trk_phi = best.getPhi();
124 double trk_theta = best.getTheta();
125 double trk_qop = best.getQoP();
126 double trk_p = (std::abs(trk_qop) > 1e-9) ? 1.0 / std::abs(trk_qop) : 0.0;
127
128 // Track position at perigee (reference point is at ECAL front face z)
129 auto perigee = best.getPerigeeLocation();
130 double trk_x = perigee[0] + (-trk_d0 * std::sin(trk_phi));
131 double trk_y = perigee[1] + (trk_d0 * std::cos(trk_phi));
132
133 ldmx_log(debug) << "Track perigee: (" << perigee[0] << ", " << perigee[1]
134 << ", " << perigee[2] << ") d0=" << trk_d0;
135 ldmx_log(debug) << "Track pos: (" << trk_x << ", " << trk_y
136 << ") p=" << trk_p << " theta=" << trk_theta
137 << " phi=" << trk_phi;
138
139 // Derived track quantities
140 double trk_pt = trk_p * std::sin(trk_theta);
141
142 // Residuals
143 double dx = trk_x - sp_x;
144 double dy = trk_y - sp_y;
145 double dtheta = trk_theta - sp_theta;
146 double dphi = trk_phi - sp_phi;
147 while (dphi > M_PI) dphi -= 2.0 * M_PI;
148 while (dphi < -M_PI) dphi += 2.0 * M_PI;
149 double dp = trk_p - sp_p;
150 double dpt = trk_pt - sp_pt;
151 double dp_frac = (sp_p > 1.0) ? dp / sp_p : 0.0;
152
153 // Track properties
154 histograms_.fill("trk_x", trk_x);
155 histograms_.fill("trk_y", trk_y);
156 histograms_.fill("trk_p", trk_p);
157 histograms_.fill("trk_pt", trk_pt);
158 histograms_.fill("trk_theta", trk_theta);
159 histograms_.fill("trk_phi", trk_phi);
160
161 // Residuals
162 histograms_.fill("delta_x", dx);
163 histograms_.fill("delta_y", dy);
164 histograms_.fill("delta_theta", dtheta);
165 histograms_.fill("delta_phi", dphi);
166 histograms_.fill("delta_p", dp);
167 histograms_.fill("delta_pt", dpt);
168 histograms_.fill("delta_p_frac", dp_frac);
169 histograms_.fill("delta_r", std::sqrt(dx * dx + dy * dy));
170 histograms_.fill("delta_angle", best_dr);
171
172 // 2D correlations
173 histograms_.fill("sp_x_vs_trk_x", sp_x, trk_x);
174 histograms_.fill("sp_y_vs_trk_y", sp_y, trk_y);
175 histograms_.fill("sp_p_vs_trk_p", sp_p, trk_p);
176 histograms_.fill("sp_theta_vs_trk_theta", sp_theta, trk_theta);
177}
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:105
void fill(const std::string &name, const T &val)
Fill a 1D histogram.
Implements detector ids for special simulation-derived hits like scoring planes.
Represents a simulated tracker hit in the simulation.
std::vector< float > getPosition() const
Get the XYZ position of the hit [mm].
std::vector< double > getMomentum() const
Get the XYZ momentum of the particle at the position at which the hit took place [MeV].
int getTrackID() const
Get the Sim particle track ID of the hit.
Implementation of a track object.
Definition Track.h:53

References framework::Event::exists(), framework::HistogramPool::fill(), ldmx::SimTrackerHit::getMomentum(), ldmx::SimTrackerHit::getPosition(), ldmx::SimTrackerHit::getTrackID(), framework::EventProcessor::histograms_, and ldmx::SimSpecialID::plane().

◆ configure()

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

Callback for the EventProcessor to configure itself from the given set of parameters.

The parameters a processor has access to are the member variables of the python class in the sequence that has class_name equal to the EventProcessor class name.

For an example, look at MyProcessor.

Parameters
parametersParameters for configuration.

Reimplemented from framework::EventProcessor.

Definition at line 11 of file EcalSPTrackCompare.cxx.

11 {
12 ecal_sp_coll_name_ = ps.get<std::string>("ecal_sp_coll_name");
13 ecal_sp_pass_name_ = ps.get<std::string>("ecal_sp_pass_name");
14 track_collection_ = ps.get<std::string>("track_collection");
15 track_pass_name_ = ps.get<std::string>("track_pass_name");
16}

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

Member Data Documentation

◆ ecal_sp_coll_name_

std::string dqm::EcalSPTrackCompare::ecal_sp_coll_name_ {"EcalScoringPlaneHits"}
private

Definition at line 31 of file EcalSPTrackCompare.h.

31{"EcalScoringPlaneHits"};

◆ ecal_sp_pass_name_

std::string dqm::EcalSPTrackCompare::ecal_sp_pass_name_ {""}
private

Definition at line 32 of file EcalSPTrackCompare.h.

32{""};

◆ track_collection_

std::string dqm::EcalSPTrackCompare::track_collection_ {"EcalTracks"}
private

Definition at line 33 of file EcalSPTrackCompare.h.

33{"EcalTracks"};

◆ track_pass_name_

std::string dqm::EcalSPTrackCompare::track_pass_name_ {""}
private

Definition at line 34 of file EcalSPTrackCompare.h.

34{""};

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