LDMX Software
EcalWABRecProcessor.cxx
1
9
10namespace ecal {
11
12// 68% Electron Radii of Containment for various theta ranges
13const std::vector<float> radius_68_theta_0_to_10 = {
14 10.12233413, 9.921772, 11.38255086, 11.67991867, 13.14337347,
15 13.17120624, 16.80994665, 17.83787244, 22.44684374, 23.74239886,
16 28.60564083, 30.27889678, 34.86404888, 36.39009394, 41.29309474,
17 43.34682279, 48.55982854, 50.80565589, 55.29496257, 57.92737879,
18 60.64828824, 65.51760517, 68.26709803, 76.32877518, 84.61219467,
19 98.21320491, 110.9880892, 120.6762931, 140.6174478, 136.4979268,
20 145.579465, 154.9803228, 164.7005, 174.7399968};
21const std::vector<float> radius_68_theta_10_to_15 = {
22 10.82307758, 11.17850518, 16.2185281, 18.62488713, 22.63408229,
23 24.71769042, 30.11217538, 32.69939046, 37.99753196, 40.81619543,
24 45.89054775, 49.03066318, 54.00440948, 59.31733555, 63.40789682,
25 64.77580021, 73.00113678, 73.25561396, 78.8914776, 86.73962133,
26 97.05926327, 96.6932739, 111.6226151, 106.5960265, 109.477541,
27 128.0220711, 145.4137195, 210.3582819, 199.6355662, 184.2513208,
28 195.076552, 206.2322029, 217.7182737, 229.5347642};
29
30const std::vector<float> radius_68_theta_15_to_20 = {
31 12.79450901, 13.02698578, 21.27450933, 25.66008312, 31.78592103,
32 35.99689874, 44.37101115, 48.82709363, 55.05972458, 59.68948687,
33 65.39866214, 70.59280337, 76.06007787, 82.22695257, 87.50371819,
34 90.60099831, 96.34848268, 101.4928478, 106.7157092, 105.0540604,
35 110.0653355, 148.3428736, 133.1449443, 146.997265, 173.3954389,
36 175.4329544, 184.3003543, 259.8415751, 215.0165, 231.1288534,
37 243.4440277, 256.085458, 269.0531444, 282.3470868};
38
39const std::vector<float> radius_68_theta_20_to_30 = {
40 14.16989595, 15.4488322, 28.31044668, 37.54285657, 48.57288885,
41 57.04243339, 68.99836079, 75.33388728, 85.00572867, 91.52574074,
42 102.5044698, 106.5315986, 116.2341378, 127.1121442, 133.8866375,
43 144.5121759, 162.1726963, 160.2986579, 171.386638, 182.5653112,
44 205.5853241, 196.3113071, 200.5907513, 228.7275694, 234.0298491,
45 251.7701385, 293.9351568, 310.521898, 344.1455457, 293.3518953,
46 303.2401036, 313.128312, 323.0165203, 332.9047287};
47
48const std::vector<float> radius_68_theta_30_to_90 = {
49 22.50983127, 26.44537503, 58.24642887, 90.59076279, 130.0592014,
50 157.4611392, 184.2187293, 202.6994588, 225.3488816, 243.3454167,
51 269.2456428, 280.6119298, 303.8591523, 322.0522722, 335.1780181,
52 350.3398234, 353.7763544, 373.9942362, 382.9453608, 401.9703438,
53 441.6281859, 432.5241826, 455.2878243, 492.2888656, 502.6653722,
54 480.2334627, 566.5438302, 505.7032783, 556.1650321, 596.9112032,
55 616.1614593, 635.4117154, 654.6619715, 673.9122276};
56
58 // Set the collection name as defined in the configuration
59 sp_pass_name_ = parameters.getParameter<std::string>("sp_pass_name");
60 collection_name_ = parameters.getParameter<std::string>("collection_name");
61 rec_pass_name_ = parameters.getParameter<std::string>("rec_pass_name");
62 rec_coll_name_ = parameters.getParameter<std::string>("rec_coll_name");
63 track_pass_name_ = parameters.getParameter<std::string>("track_pass_name");
64 track_coll_name_ = parameters.getParameter<std::string>("track_coll_name");
65}
66
68 // Define start time for processing
69 auto start = std::chrono::high_resolution_clock::now();
70 nevents_++;
71
72 // Keep track of event progress where:
73 // 0: No tracks found
74 // 1: Track found but not enough info to reconstruct either electron or photon
75 // 2: Track found and enough info to reconstruct electron
76 // 3: Track found and enough info to reconstruct electron and photon
77 int progress_num = 0;
78
79 // Get the Ecal Geometry
81 ldmx::EcalGeometry::CONDITIONS_OBJECT_NAME);
82
83 // Get the collection of ecal_rec_hits and tracks
84 const std::vector<ldmx::EcalHit> ecal_rec_hits =
85 event.getCollection<ldmx::EcalHit>(rec_coll_name_, rec_pass_name_);
86 const std::vector<ldmx::StraightTrack> linear_tracks =
87 event.getCollection<ldmx::StraightTrack>(track_coll_name_,
88 track_pass_name_);
89
90 // Define variables to save recoil electron/photon information (SP)
91 std::vector<double> recoil_e_p, recoil_y_p;
92 std::vector<float> recoil_e_pos, recoil_y_pos;
93
94 // Result object that stores kinematic variables
96
97 // Define kinematic variables
98 const std::vector<float> z_hat = {0, 0, 1};
99 float true_theta_electron = -9.;
100 float true_theta_photon = -9.;
101 float true_phi_electron = -9.;
102 float true_phi_photon = -9.;
103 float rec_theta_electron = -9.;
104 float rec_theta_photon = -9.;
105 float rec_phi_electron = -9.;
106 float rec_phi_photon = -9.;
107 float true_theta_diff_electron_photon = -9.;
108 float true_phi_diff_electron_photon = -9.;
109 float rec_theta_diff_electron_photon = -9.;
110 float rec_phi_diff_electron_photon = -9.;
111 float true_rec_theta_diff_electron = -9.;
112 float true_rec_phi_diff_electron = -9.;
113 float true_rec_theta_diff_photon = -9.;
114 float true_rec_phi_diff_photon = -9.;
115 float true_electron_shower_energy = -999.;
116 float true_photon_shower_energy = -999.;
117 float rec_electron_shower_energy = -999.;
118 float rec_photon_shower_energy = -999.;
119
120 // Create lists for rec hits and electron/photon shower hits
121 std::vector<std::array<float, 6>> rec_hit_list;
122 std::vector<std::array<float, 6>> ele_hit_list;
123 std::vector<std::array<float, 6>> phot_hit_list;
124
125 // Save rec hit info to rec_hit_list
126 for (const ldmx::EcalHit& hit : ecal_rec_hits) {
127 ldmx::EcalID id(hit.getID());
128 auto pos = geometry_->getPosition(id);
129 auto [x, y, z] = std::apply(
130 [](double a, double b, double c) {
131 return std::make_tuple(static_cast<float>(a), static_cast<float>(b),
132 static_cast<float>(c));
133 },
134 pos);
135 float energy = hit.getEnergy();
136 float layer_num = id.layer();
137 rec_hit_list.push_back({x, y, z, layer_num, 0, energy});
138 }
139
140 if (event.exists("TargetScoringPlaneHits", sp_pass_name_)) {
141 //
142 // Loop through all of the sim particles and find the recoil
143 // photon/electron.
144 //
145
146 // Find Target SP hit for recoil photon/electron
147 const std::vector<ldmx::SimTrackerHit> target_sp_hits =
148 event.getCollection<ldmx::SimTrackerHit>("TargetScoringPlaneHits",
149 sp_pass_name_);
150 float photon_p_zmax = 0, electron_p_zmax = 0;
151 for (const ldmx::SimTrackerHit& sp_hit : target_sp_hits) {
152 ldmx::SimSpecialID hit_id(sp_hit.getID());
153 if (hit_id.plane() != 1 || sp_hit.getMomentum()[2] <= 0) continue;
154
155 if (sp_hit.getPdgID() == 11) {
156 if (sp_hit.getMomentum()[2] > electron_p_zmax) {
157 recoil_e_p = sp_hit.getMomentum();
158 true_electron_shower_energy = sp_hit.getEnergy();
159 recoil_e_pos = sp_hit.getPosition();
160 electron_p_zmax = recoil_e_p[2];
161 }
162 }
163 if (sp_hit.getPdgID() == 22) {
164 if (sp_hit.getMomentum()[2] > photon_p_zmax) {
165 recoil_y_p = sp_hit.getMomentum();
166 true_photon_shower_energy = sp_hit.getEnergy();
167 recoil_y_pos = sp_hit.getPosition();
168 photon_p_zmax = recoil_y_p[2];
169 }
170 }
171 }
172
173 // Calculating true theta values using SP hit parameters
174 if (recoil_y_p.size() == 3 &&
175 std::sqrt((recoil_y_p[0]) * (recoil_y_p[0]) +
176 (recoil_y_p[1]) * (recoil_y_p[1]) +
177 (recoil_y_p[2]) * (recoil_y_p[2])) != 0 &&
178 recoil_e_p.size() == 3 &&
179 std::sqrt((recoil_e_p[0]) * (recoil_e_p[0]) +
180 (recoil_e_p[1]) * (recoil_e_p[1]) +
181 (recoil_e_p[2]) * (recoil_e_p[2])) != 0) {
182 true_theta_electron =
183 (180 / std::numbers::pi) *
184 std::acos(std::inner_product(recoil_e_p.begin(), recoil_e_p.end(),
185 z_hat.begin(), 0.0) /
186 std::sqrt((recoil_e_p[0]) * (recoil_e_p[0]) +
187 (recoil_e_p[1]) * (recoil_e_p[1]) +
188 (recoil_e_p[2]) * (recoil_e_p[2])));
189 true_theta_photon =
190 (180 / std::numbers::pi) *
191 std::acos(std::inner_product(recoil_y_p.begin(), recoil_y_p.end(),
192 z_hat.begin(), 0.0) /
193 std::sqrt((recoil_y_p[0]) * (recoil_y_p[0]) +
194 (recoil_y_p[1]) * (recoil_y_p[1]) +
195 (recoil_y_p[2]) * (recoil_y_p[2])));
196 }
197
198 // Calculating true phi values using SP hit parameters
199 if (recoil_y_p.size() == 3 && recoil_y_p[2] != 0 &&
200 recoil_e_p.size() == 3 && recoil_e_p[2] != 0) {
201 true_phi_electron =
202 (180 / std::numbers::pi) * std::atan(recoil_e_p[1] / recoil_e_p[0]);
203 if (recoil_e_p[1] < 0) {
204 true_phi_electron += 180;
205 }
206 if (recoil_e_p[0] < 0 && recoil_e_p[1] > 0) {
207 true_phi_electron += 360;
208 }
209 true_phi_photon =
210 (180 / std::numbers::pi) * std::atan(recoil_y_p[1] / recoil_y_p[0]);
211 if (recoil_y_p[1] < 0) {
212 true_phi_photon += 180;
213 }
214 if (recoil_y_p[0] < 0 && recoil_y_p[1] > 0) {
215 true_phi_photon += 360;
216 }
217 }
218
219 // Calculating true delta_phi/delta_theta using SP hit parameters
220 if (recoil_y_p.size() == 3 && recoil_e_p.size() == 3) {
221 std::array<double, 2> phi_diff_electron_arr = {recoil_e_p[0],
222 recoil_e_p[1]};
223 std::array<double, 2> phi_diff_photon_arr = {recoil_y_p[0],
224 recoil_y_p[1]};
225 std::array<double, 2> theta_diff_electron_arr = {recoil_e_p[2],
226 recoil_e_p[0]};
227 std::array<double, 2> theta_diff_photon_arr = {recoil_y_p[2],
228 recoil_y_p[0]};
229
230 true_theta_diff_electron_photon =
231 (180 / std::numbers::pi) *
232 std::acos(
233 std::inner_product(theta_diff_electron_arr.begin(),
234 theta_diff_electron_arr.end(),
235 theta_diff_photon_arr.begin(), 0.0) /
236 (std::sqrt((theta_diff_electron_arr[0]) *
237 (theta_diff_electron_arr[0]) +
238 (theta_diff_electron_arr[1]) *
239 (theta_diff_electron_arr[1])) *
240 std::sqrt(
241 (theta_diff_photon_arr[0]) * (theta_diff_photon_arr[0]) +
242 (theta_diff_photon_arr[1]) * (theta_diff_photon_arr[1]))));
243 true_phi_diff_electron_photon =
244 (180 / std::numbers::pi) *
245 std::acos(
246 std::inner_product(phi_diff_electron_arr.begin(),
247 phi_diff_electron_arr.end(),
248 phi_diff_photon_arr.begin(), 0.0) /
249 (std::sqrt(
250 (phi_diff_electron_arr[0]) * (phi_diff_electron_arr[0]) +
251 (phi_diff_electron_arr[1]) * (phi_diff_electron_arr[1])) *
252 std::sqrt((phi_diff_photon_arr[0]) * (phi_diff_photon_arr[0]) +
253 (phi_diff_photon_arr[1]) * (phi_diff_photon_arr[1]))));
254 }
255 }
256
257 // Defining variables to save best fit results
258 std::pair<Eigen::VectorXd, Eigen::VectorXd> linear_fit_coeffs;
259 std::tuple<Eigen::VectorXd, float, int, Eigen::MatrixXd, int> best_x_result;
260 std::tuple<Eigen::VectorXd, float, int, Eigen::MatrixXd, int> best_y_result;
261 std::get<1>(best_x_result) = 10e99;
262 std::get<1>(best_y_result) = 10e99;
263
264 // Looping over tracks to find best fit to hits
265 std::vector<float> ele_roc = radius_68_theta_30_to_90;
266 for (const ldmx::StraightTrack& track : linear_tracks) {
267 progress_num = 1;
268 // Determining the RoC value to use based on recoil electron (track) theta
269 std::vector<double> track_vec = {track.getSlopeX(), track.getSlopeY(), 1};
270 float track_theta =
271 (180 / std::numbers::pi) *
272 std::acos(std::inner_product(track_vec.begin(), track_vec.end(),
273 z_hat.begin(), 0.0) /
274 std::sqrt((track_vec[0]) * (track_vec[0]) +
275 (track_vec[1]) * (track_vec[1]) + 1));
276 if (track_theta <= 10) {
277 ele_roc = radius_68_theta_0_to_10;
278 } else if (track_theta > 10 && track_theta <= 15) {
279 ele_roc = radius_68_theta_10_to_15;
280 } else if (track_theta > 15 && track_theta <= 20) {
281 ele_roc = radius_68_theta_15_to_20;
282 } else if (track_theta > 20 && track_theta <= 30) {
283 ele_roc = radius_68_theta_20_to_30;
284 }
285
286 // Labeling hits as electron (1) or photon (0)
287 for (std::array<float, 6>& hit : rec_hit_list) {
288 if (std::sqrt(
289 (hit[0] - (track.getSlopeX() * hit[2] + track.getInterceptX())) *
290 (hit[0] -
291 (track.getSlopeX() * hit[2] + track.getInterceptX())) +
292 (hit[1] - (track.getSlopeY() * hit[2] + track.getInterceptY())) *
293 (hit[1] - (track.getSlopeY() * hit[2] +
294 track.getInterceptY()))) < ele_roc[hit[3]]) {
295 hit[4] = 1;
296 }
297 }
298 // Create vectors to hold electron/photon hits specifically
299 std::vector<float> ele_hit_list_x, ele_hit_list_y, ele_hit_list_z;
300 std::vector<float> phot_hit_list_x, phot_hit_list_y, phot_hit_list_z;
301
302 // Use labels to sort hits as electron/photon and calculate shower energies
303 for (const auto& hit : rec_hit_list) {
304 if (hit[4] == 1) {
305 ele_hit_list.push_back(hit);
306 ele_hit_list_x.push_back(hit[0]);
307 ele_hit_list_y.push_back(hit[1]);
308 ele_hit_list_z.push_back(hit[2]);
309 } else if (hit[4] == 0) {
310 phot_hit_list.push_back(hit);
311 phot_hit_list_x.push_back(hit[0]);
312 phot_hit_list_y.push_back(hit[1]);
313 phot_hit_list_z.push_back(hit[2]);
314 }
315 }
316
317 // Fit both photon/electron or just electron hits based on # of viable
318 // showers
319 if (phot_hit_list.size() >= 3 && ele_hit_list.size() >= 3) {
320 progress_num =
321 3; // Set progress_num to 3 to halt electron-only reconstruction
322
323 // Generate guesses and error vectors for vertex constrained fit
324 std::vector<double> x_guess = {
325 track.getSlopeX(),
326 (phot_hit_list.back()[0] - phot_hit_list[0][0]) /
327 (phot_hit_list.back()[2] - phot_hit_list[0][2]),
328 track.getInterceptX()};
329 std::vector<double> y_guess = {
330 track.getSlopeY(),
331 (phot_hit_list.back()[1] - phot_hit_list[0][1]) /
332 (phot_hit_list.back()[2] - phot_hit_list[0][2]),
333 track.getInterceptY()};
334 std::vector<float> phot_hit_error(phot_hit_list.size(),
335 0.456435464588 * 4.816);
336 std::vector<float> ele_hit_error(ele_hit_list.size(),
337 0.456435464588 * 4.816);
338
339 int max_iter = 200;
340 // Carry out fit
341 std::tuple<Eigen::VectorXd, float, int, Eigen::MatrixXd, int> x_result =
342 fit2DTracksConstrained(ele_hit_list_z, ele_hit_list_x, ele_hit_error,
343 phot_hit_list_z, phot_hit_list_x,
344 phot_hit_error, x_guess, max_iter, 0, 0.001,
345 10.0);
346
347 std::tuple<Eigen::VectorXd, float, int, Eigen::MatrixXd, int> y_result =
348 fit2DTracksConstrained(ele_hit_list_z, ele_hit_list_y, ele_hit_error,
349 phot_hit_list_z, phot_hit_list_y,
350 phot_hit_error, y_guess, max_iter, 0, 0.001,
351 40.0);
352
353 // Update best fit variables if current fit is an improvement
354 if ((std::get<1>(x_result) + std::get<1>(y_result)) / 2 <
355 (std::get<1>(best_x_result) + std::get<1>(best_y_result)) / 2) {
356 best_x_result = x_result;
357 best_y_result = y_result;
358 }
359 }
360
361 // If there isn't enough info for both electron and photon reconstruction,
362 // reconstruct electron if possible
363 else if (ele_hit_list.size() >= 3) {
364 if (progress_num != 3) {
365 progress_num = 2; // Set progress_num to 3 to indicate electron-only
366 // reconstruction
367 linear_fit_coeffs =
368 polyfitXYvsZ(ele_hit_list_x, ele_hit_list_y, ele_hit_list_z, 1);
369 }
370 }
371 }
372
373 // Calculate kinematic variables for electron and/or photon
374 // based on # of viable showers (with reconstruction information)
375 if (std::get<0>(best_x_result).size() != 0) {
376 rec_electron_shower_energy = 0;
377 rec_photon_shower_energy = 0;
378 for (const auto& hit : ele_hit_list) {
379 rec_electron_shower_energy += hit[5];
380 }
381 for (const auto& hit : phot_hit_list) {
382 rec_photon_shower_energy += hit[5];
383 }
384
385 std::vector<double> ele_params = {std::get<0>(best_x_result)(0),
386 std::get<0>(best_y_result)(0)};
387 std::vector<double> phot_params = {std::get<0>(best_x_result)(1),
388 std::get<0>(best_y_result)(1)};
389 std::vector<double> ele_params_x = {std::get<0>(best_x_result)(0)};
390 std::vector<double> phot_params_x = {std::get<0>(best_x_result)(1)};
391
392 rec_theta_electron =
393 (180 / std::numbers::pi) *
394 std::acos(1 / std::sqrt((ele_params[0]) * (ele_params[0]) +
395 (ele_params[1]) * (ele_params[1]) + 1));
396 rec_theta_photon =
397 (180 / std::numbers::pi) *
398 std::acos(1 / std::sqrt((phot_params[0]) * (phot_params[0]) +
399 (phot_params[1]) * (phot_params[1]) + 1));
400
401 rec_phi_electron =
402 (180 / std::numbers::pi) * std::atan(ele_params[1] / ele_params[0]);
403 if (ele_params[1] < 0) {
404 rec_phi_electron += 180;
405 }
406 if (ele_params[0] < 0 && ele_params[1] > 0) {
407 rec_phi_electron += 360;
408 }
409 rec_phi_photon =
410 (180 / std::numbers::pi) * std::atan(phot_params[1] / phot_params[0]);
411 if (phot_params[1] < 0) {
412 rec_phi_photon += 180;
413 }
414 if (phot_params[0] < 0 && phot_params[1] > 0) {
415 rec_phi_photon += 360;
416 }
417
418 rec_theta_diff_electron_photon =
419 (180 / std::numbers::pi) *
420 std::acos(std::inner_product(ele_params_x.begin(), ele_params_x.end(),
421 phot_params_x.begin(), 1.0) /
422 (std::sqrt((ele_params_x[0]) * (ele_params_x[0]) + (1)) *
423 std::sqrt((phot_params_x[0]) * (phot_params_x[0]) + 1)));
424 rec_phi_diff_electron_photon =
425 (180 / std::numbers::pi) *
426 std::acos(std::inner_product(ele_params.begin(), ele_params.end(),
427 phot_params.begin(), 0.0) /
428 (std::sqrt((ele_params[0]) * (ele_params[0]) +
429 (ele_params[1]) * (ele_params[1])) *
430 std::sqrt((phot_params[0]) * (phot_params[0]) +
431 (phot_params[1]) * (phot_params[1]))));
432
433 if (recoil_y_p.size() == 3 && recoil_y_p[2] != 0 &&
434 recoil_e_p.size() == 3 && recoil_e_p[2] != 0) {
435 true_rec_theta_diff_electron =
436 (180 / std::numbers::pi) *
437 std::acos(std::inner_product(ele_params_x.begin(), ele_params_x.end(),
438 recoil_e_p.begin(), recoil_e_p[2]) /
439 (std::sqrt((ele_params_x[0]) * (ele_params_x[0]) + 1) *
440 std::sqrt((recoil_e_p[0]) * (recoil_e_p[0]) +
441 (recoil_e_p[2]) * (recoil_e_p[2]))));
442 true_rec_theta_diff_photon =
443 (180 / std::numbers::pi) *
444 std::acos(std::inner_product(phot_params_x.begin(),
445 phot_params_x.end(), recoil_y_p.begin(),
446 recoil_y_p[2]) /
447 (std::sqrt((phot_params_x[0]) * (phot_params_x[0]) + 1) *
448 std::sqrt((recoil_y_p[0]) * (recoil_y_p[0]) +
449 (recoil_y_p[2]) * (recoil_y_p[2]))));
450 true_rec_phi_diff_electron =
451 (180 / std::numbers::pi) *
452 std::acos(std::inner_product(ele_params.begin(), ele_params.end(),
453 recoil_e_p.begin(), 0) /
454 (std::sqrt((ele_params[0]) * (ele_params[0]) +
455 (ele_params[1]) * (ele_params[1])) *
456 std::sqrt((recoil_e_p[0]) * (recoil_e_p[0]) +
457 (recoil_e_p[1]) * (recoil_e_p[1]))));
458 true_rec_phi_diff_photon =
459 (180 / std::numbers::pi) *
460 std::acos(std::inner_product(phot_params.begin(), phot_params.end(),
461 recoil_y_p.begin(), 0) /
462 (std::sqrt((phot_params[0]) * (phot_params[0]) +
463 (phot_params[1]) * (phot_params[1])) *
464 std::sqrt((recoil_y_p[0]) * (recoil_y_p[0]) +
465 (recoil_y_p[1]) * (recoil_y_p[1]))));
466 }
467 } else if (progress_num == 2) {
468 rec_electron_shower_energy = 0;
469 for (const auto& hit : ele_hit_list) {
470 rec_electron_shower_energy += hit[5];
471 }
472
473 std::vector<double> ele_params = {linear_fit_coeffs.first(1),
474 linear_fit_coeffs.second(1)};
475 std::vector<double> ele_params_x = {linear_fit_coeffs.first(1)};
476
477 rec_theta_electron =
478 (180 / std::numbers::pi) *
479 std::acos(1 / std::sqrt((ele_params[0]) * (ele_params[0]) +
480 (ele_params[1]) * (ele_params[1]) + 1));
481 rec_phi_electron =
482 (180 / std::numbers::pi) * std::atan(ele_params[1] / ele_params[0]);
483 if (ele_params[1] < 0) {
484 rec_phi_electron += 180;
485 }
486 if (ele_params[0] < 0 && ele_params[1] > 0) {
487 rec_phi_electron += 360;
488 }
489
490 if (recoil_y_p.size() == 3 && recoil_y_p[2] != 0 &&
491 recoil_e_p.size() == 3 && recoil_e_p[2] != 0) {
492 true_rec_theta_diff_electron =
493 (180 / std::numbers::pi) *
494 std::acos(std::inner_product(ele_params_x.begin(), ele_params_x.end(),
495 recoil_e_p.begin(), recoil_e_p[2]) /
496 (std::sqrt((ele_params_x[0]) * (ele_params_x[0]) + 1) *
497 std::sqrt((recoil_e_p[0]) * (recoil_e_p[0]) +
498 (recoil_e_p[2]) * (recoil_e_p[2]))));
499 true_rec_phi_diff_electron =
500 (180 / std::numbers::pi) *
501 std::acos(std::inner_product(ele_params.begin(), ele_params.end(),
502 recoil_e_p.begin(), 0) /
503 (std::sqrt((ele_params[0]) * (ele_params[0]) +
504 (ele_params[1]) * (ele_params[1])) *
505 std::sqrt((recoil_e_p[0]) * (recoil_e_p[0]) +
506 (recoil_e_p[1]) * (recoil_e_p[1]))));
507 }
508
509 // Set photon variables to non-physical value
510 // that corresponds to electron-only reco case
511 rec_theta_photon = -5.;
512 rec_phi_photon = -5.;
513 rec_theta_diff_electron_photon = -5.;
514 rec_phi_diff_electron_photon = -5.;
515 true_rec_theta_diff_photon = -5.;
516 true_rec_phi_diff_photon = -5.;
517 }
518
519 // Setting output object equal to calculated variables
520 result.setVariables(
521 true_theta_electron, true_theta_photon, true_phi_electron,
522 true_phi_photon, rec_theta_electron, rec_theta_photon, rec_phi_electron,
523 rec_phi_photon, true_theta_diff_electron_photon,
524 true_phi_diff_electron_photon, rec_theta_diff_electron_photon,
525 rec_phi_diff_electron_photon, true_rec_theta_diff_electron,
526 true_rec_phi_diff_electron, true_rec_theta_diff_photon,
527 true_rec_phi_diff_photon, true_electron_shower_energy,
528 true_photon_shower_energy, rec_electron_shower_energy,
529 rec_photon_shower_energy, progress_num);
530
531 event.add(collection_name_, result);
532
533 // Calculate processing time for event
534 auto end = std::chrono::high_resolution_clock::now();
535 auto diff = end - start;
536 processing_time_ += std::chrono::duration<float, std::milli>(diff).count();
537}
538
540 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(2)
541 << processing_time_ / nevents_ << " ms";
542}
543
544std::tuple<Eigen::VectorXd, float, int, Eigen::MatrixXd, int>
545EcalWABRecProcessor::fit2DTracksConstrained(
546 const std::vector<float>& x1, const std::vector<float>& y1,
547 const std::vector<float>& s1, const std::vector<float>& x2,
548 const std::vector<float>& y2, const std::vector<float>& s2,
549 const std::vector<double>& guess, int max_iter, int verbose, float d_chisq,
550 float abs_lim) {
551 /*
552 Function that fits two 2D tracks with a vertex constraint (same intercepts).
553 The fitted model is initially defined as:
554 y1 = par[0] * x1 + abs_lim * tanh(par[2] / abs_lim)
555 y2 = par[1] * x2 + abs_lim * tanh(par[2] / abs_lim)
556 After updating parameters the fitted values are computed as:
557 y1 = par[0] * (x1 - par[2])
558 y2 = par[1] * (x2 - par[2])
559
560 Inputs:
561 x1, y1, s1 : measured coordinates and errors for track 1
562 x2, y2, s2 : measured coordinates and errors for track 2
563 guess : initial guess for the parameter vector (size 3)
564 max_iter : maximum number of iterations (default 20)
565 verbose : level of verbose (default 0)
566 d_chisq : stopping criterion for chi-squared improvement (default
567 0.001) abs_lim : a parameter used in the fit (default 10) Returns: A tuple
568 containing: par : fitted parameters (Eigen::VectorXd of size 3) chisq :
569 chi-squared at minimum (float) ndof : number of degrees of freedom (int)
570 cov : covariance matrix (Eigen::MatrixXd 3x3)
571 niter : number of iterations used (int)
572 */
573 // Copy the initial guess into a 3-element parameter vector
574 Eigen::VectorXd par(3);
575 par(0) = guess[0];
576 par(1) = guess[1];
577 par(2) = guess[2];
578
579 // Determine number of points in each track and total
580 int n1 = x1.size();
581 int n2 = x2.size();
582 int n = n1 + n2;
583
584 // Concatenate x, y, and s into Eigen vectors of size n.
585 Eigen::VectorXd x(n), y(n), s(n);
586 for (int i = 0; i < n1; ++i) {
587 x(i) = x1[i];
588 y(i) = y1[i];
589 s(i) = s1[i];
590 }
591 for (int i = 0; i < n2; ++i) {
592 x(n1 + i) = x2[i];
593 y(n1 + i) = y2[i];
594 s(n1 + i) = s2[i];
595 }
596
597 // Build the weight matrix W = diag(1/s_i^2)
598 Eigen::MatrixXd W = Eigen::MatrixXd::Zero(n, n);
599 for (int i = 0; i < n; ++i) {
600 W(i, i) = 1.0 / ((s(i)) * (s(i)));
601 }
602
603 float chi_sq = 0.0;
604 float old_chi_sq = 1e12; // a large initial value
605 int n_iter = 0;
606 Eigen::MatrixXd cov(3, 3); // covariance matrix
607
608 // Iterative fitting loop
609 for (int iter = 0; iter < max_iter; ++iter) {
610 n_iter = iter + 1;
611
612 // Compute fitted y coordinates for each track using the current parameters.
613 // For track 1: y1_fit = par[0] * x1 + abs_lim * tanh(par[2]/abs_lim)
614 // For track 2: y2_fit = par[1] * x2 + abs_lim * tanh(par[2]/abs_lim)
615 Eigen::VectorXd y1_fit(n1), y2_fit(n2);
616 float tanh_term = std::tanh(par(2) / abs_lim);
617 for (int i = 0; i < n1; ++i) {
618 y1_fit(i) = par(0) * x1[i] + abs_lim * tanh_term;
619 }
620 for (int i = 0; i < n2; ++i) {
621 y2_fit(i) = par(1) * x2[i] + abs_lim * tanh_term;
622 }
623
624 // Concatenate the fitted values
625 Eigen::VectorXd y_fit(n);
626 for (int i = 0; i < n1; ++i) {
627 y_fit(i) = y1_fit(i);
628 }
629 for (int i = 0; i < n2; ++i) {
630 y_fit(n1 + i) = y2_fit(i);
631 }
632
633 // Compute chi-squared: sum_i [ (y_fit[i]-y[i])^2 / s[i]^2 ]
634 chi_sq = 0.0;
635 for (int i = 0; i < n; ++i) {
636 float diff = y_fit(i) - y(i);
637 chi_sq += (diff * diff) / ((s(i)) * (s(i)));
638 }
639
640 if (verbose > 0) {
641 ldmx_log(debug) << "Before iteration " << iter << ", chi_sq = " << chi_sq;
642 ldmx_log(debug) << "Track 1 residuals: ";
643 for (int i = 0; i < n1; ++i) {
644 ldmx_log(debug) << (y1_fit(i) - y1[i]);
645 }
646 ldmx_log(debug) << "Track 2 residuals: ";
647 for (int i = 0; i < n2; ++i) {
648 ldmx_log(debug) << (y2_fit(i) - y2[i]);
649 }
650 }
651
652 // Compute the derivatives (Jacobian components)
653 // For track 1:
654 // dy1/dpar0 = x1, dy1/dpar1 = 0, dy1/dpar2 =
655 // (1/cosh(par[2]/abs_lim))^2 (constant for all points)
656 // For track 2:
657 // dy2/dpar0 = 0, dy2/dpar1 = x2, dy2/dpar2 =
658 // (1/cosh(par[2]/abs_lim))^2
659 Eigen::VectorXd dy1_dpar_0(n1), dy1_dpar_1 = Eigen::VectorXd::Zero(n1),
660 dy1_dpar_2(n1);
661 Eigen::VectorXd dy2_dpar_0 = Eigen::VectorXd::Zero(n2), dy2_dpar_1(n2),
662 dy2_dpar_2(n2);
663
664 float d_term = 1.0 / std::cosh(par(2) / abs_lim);
665 d_term = d_term * d_term; // square it
666 for (int i = 0; i < n1; ++i) {
667 dy1_dpar_0(i) = x1[i];
668 dy1_dpar_2(i) = d_term;
669 }
670 for (int i = 0; i < n2; ++i) {
671 dy2_dpar_1(i) = x2[i];
672 dy2_dpar_2(i) = d_term;
673 }
674
675 // Concatenate the derivatives for both tracks into full vectors of length
676 // n.
677 Eigen::VectorXd dy_dpar_0(n), dy_dpar_1(n), dy_dpar_2(n);
678 for (int i = 0; i < n1; ++i) {
679 dy_dpar_0(i) = dy1_dpar_0(i);
680 dy_dpar_1(i) = dy1_dpar_1(i);
681 dy_dpar_2(i) = dy1_dpar_2(i);
682 }
683 for (int i = 0; i < n2; ++i) {
684 dy_dpar_0(n1 + i) = dy2_dpar_0(i);
685 dy_dpar_1(n1 + i) = dy2_dpar_1(i);
686 dy_dpar_2(n1 + i) = dy2_dpar_2(i);
687 }
688
689 // Build the "A" matrix (the Jacobian) in its transposed form (3 x n)
690 Eigen::MatrixXd a_trans(3, n);
691 a_trans.row(0) = dy_dpar_0.transpose();
692 a_trans.row(1) = dy_dpar_1.transpose();
693 a_trans.row(2) = dy_dpar_2.transpose();
694
695 // The Jacobian (n x 3) is the transpose of a_trans.
696 Eigen::MatrixXd a = a_trans.transpose();
697
698 // The residual vector (difference between measured and fitted y values)
699 Eigen::VectorXd dy_vec = y - y_fit;
700
701 // Compute the (3 x 3) matrix: M = a_trans * W * a
702 Eigen::MatrixXd temp = a_trans * W; // 3 x n
703 Eigen::MatrixXd temp2 = temp * a; // 3 x 3
704
705 // Add a regularization term to ensure numerical stability.
706 Eigen::MatrixXd reg =
707 1e-10 * Eigen::MatrixXd::Identity(temp2.rows(), temp2.cols());
708 Eigen::MatrixXd temp2_reg = temp2 + reg;
709
710 // Invert the matrix to obtain the covariance matrix.
711 cov = temp2_reg.inverse();
712
713 // Compute the parameter correction: dpar = cov * a_trans * W * dy_vec
714 Eigen::MatrixXd temp4 = cov * a_trans; // 3 x n
715 Eigen::MatrixXd temp5 = temp4 * W; // 3 x n
716 Eigen::VectorXd dpar = temp5 * dy_vec; // 3 x 1
717
718 // Update the parameters
719 par += dpar;
720
721 // After the update, the fitted y values are recalculated with a different
722 // formula:
723 // y1_fit = par[0]*(x1 - par[2])
724 // y2_fit = par[1]*(x2 - par[2])
725 for (int i = 0; i < n1; ++i) {
726 y1_fit(i) = par(0) * (x1[i] - par(2));
727 }
728 for (int i = 0; i < n2; ++i) {
729 y2_fit(i) = par(1) * (x2[i] - par(2));
730 }
731 for (int i = 0; i < n1; ++i) {
732 y_fit(i) = y1_fit(i);
733 }
734 for (int i = 0; i < n2; ++i) {
735 y_fit(n1 + i) = y2_fit(i);
736 }
737
738 // Recompute chi-squared with the updated fitted values.
739 float new_chi_sq = 0.0;
740 for (int i = 0; i < n; ++i) {
741 float diff = y_fit(i) - y(i);
742 new_chi_sq += (diff * diff) / ((s(i)) * (s(i)));
743 }
744 chi_sq = new_chi_sq;
745
746 // Check for convergence
747 if (iter > 0) {
748 if (std::abs(chi_sq - old_chi_sq) < d_chisq) {
749 break;
750 }
751 }
752 old_chi_sq = chi_sq;
753 } // end for loop
754
755 if (verbose > 0) {
756 ldmx_log(debug) << "At the end chi_sq = " << chi_sq;
757 ldmx_log(debug) << "Scaled residuals for track 1:";
758 for (int i = 0; i < n1; ++i) {
759 float fit_val = par(0) * (x1[i] - par(2));
760 ldmx_log(debug) << 10000 * (fit_val - y1[i]);
761 }
762 ldmx_log(debug) << "Scaled residuals for track 2:";
763 for (int i = 0; i < n2; ++i) {
764 float fit_val = par(1) * (x2[i] - par(2));
765 ldmx_log(debug) << 10000 * (fit_val - y2[i]);
766 }
767 }
768
769 int ndof = n1 + n2 - 3;
770 return std::make_tuple(par, chi_sq, ndof, cov, n_iter);
771}
772
773std::pair<Eigen::VectorXd, Eigen::VectorXd> EcalWABRecProcessor::polyfitXYvsZ(
774 const std::vector<float>& x, const std::vector<float>& y,
775 const std::vector<float>& z, int degree) {
776 /*
777 Function that fits two polynomials (x vs. z and y vs. z) to 3D hit position
778 data using a least-squares method. The fitted models are defined as: x = a₀
779 + a₁ * z + a₂ * z² + ... + aₙ * zⁿ
780 y = b₀ + b₁ * z + b₂ * z² + ... + bₙ * zⁿ
781 where n is the specified polynomial degree.
782
783 Inputs:
784 x, y, z : measured coordinates for the tracks;
785 x and y are the dependent variables, and z is the independent
786 variable (all provided as std::vector<float>) degree : degree of the
787 polynomial to be fitted (int)
788
789 Returns:
790 A pair containing:
791 first : polynomial coefficients for the x vs. z fit (Eigen::VectorXd)
792 second : polynomial coefficients for the y vs. z fit (Eigen::VectorXd)
793
794 Notes:
795 The polynomial is represented with the constant term first (i.e., [a₀, a₁,
796 ..., aₙ]), so the linear term (slope) is located at index 1.
797 */
798 const size_t n = z.size();
799 if (n == 0 || x.size() != n || y.size() != n) {
800 throw std::invalid_argument(
801 "Vectors x, y, and z must be non-empty and have the same size.");
802 }
803
804 // Construct the Vandermonde (design) matrix A (n x (degree + 1)):
805 // Each row i: [1, z[i], z[i]^2, ..., z[i]^degree]
806 Eigen::MatrixXd A(n, degree + 1);
807 for (size_t i = 0; i < n; ++i) {
808 float term = 1.0;
809 for (int j = 0; j <= degree; ++j) {
810 A(i, j) = term;
811 term *= z[i];
812 }
813 }
814
815 // Map the x and y data into Eigen vectors.
816 Eigen::VectorXd bx(n), by(n);
817 for (size_t i = 0; i < n; ++i) {
818 bx(i) = x[i];
819 by(i) = y[i];
820 }
821
822 // Solve the least-squares problems:
823 // A * coeffsX ≈ bx and A * coeffsY ≈ by
824 Eigen::VectorXd coeffsX = A.colPivHouseholderQr().solve(bx);
825 Eigen::VectorXd coeffsY = A.colPivHouseholderQr().solve(by);
826
827 return {coeffsX, coeffsY};
828}
829} // namespace ecal
830
831DECLARE_PRODUCER_NS(ecal, EcalWABRecProcessor);
Class that reconstructs important kinematic variables for WAB studies.
#define DECLARE_PRODUCER_NS(NS, CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
void produce(framework::Event &event) override
Process the event and put new data products into it.
std::string collection_name_
Name of the collection which will contain the results.
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
Callback for the EventProcessor to configure itself from the given set of parameters.
const T & getCondition(const std::string &condition_name)
Access a conditions object for the current event.
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
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
Translation between real-space positions and cell IDs within the ECal.
std::tuple< double, double, double > getPosition(EcalID id) const
Get a cell's position from its ID number.
Stores reconstructed hit information from the ECAL.
Definition EcalHit.h:19
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
Definition EcalID.h:20
Implements detector ids for special simulation-derived hits like scoring planes.
int plane() const
Get the value of the plane field from the ID, if it is a scoring plane.
Represents a simulated tracker hit in the simulation.