LDMX Software
CKFProcessor.cxx
1#include "Tracking/Reco/CKFProcessor.h"
2
3#include "Acts/EventData/TrackContainer.hpp"
4#include "Acts/Utilities/TrackHelpers.hpp"
5#include "SimCore/Event/SimParticle.h"
6#include "Tracking/Reco/TruthMatchingTool.h"
7#include "Tracking/Sim/GeometryContainers.h"
8
9//--- C++ StdLib ---//
10#include <algorithm> //std::vector reverse
11#include <iostream>
12#include <typeinfo>
13// eN files
14#include <fstream>
15
16namespace tracking {
17namespace reco {
18
19CKFProcessor::CKFProcessor(const std::string& name, framework::Process& process)
20 : TrackingGeometryUser(name, process) {}
21
23 profiling_map_["setup"] = 0.;
24 profiling_map_["hits"] = 0.;
25 profiling_map_["seeds"] = 0.;
26 profiling_map_["ckf_setup"] = 0.;
27 profiling_map_["ckf_run"] = 0.;
28 profiling_map_["result_loop"] = 0.;
29
30 // Generate a constant magnetic field
31 Acts::Vector3 b_field(0., 0., bfield_ * Acts::UnitConstants::T);
32
33 // Setup a constant magnetic field
34 const auto const_b_field = std::make_shared<Acts::ConstantBField>(b_field);
35
36 // Define the target surface - be careful:
37 // x - downstream
38 // y - left (when looking along x)
39 // z - up
40 // Passing identity here means that your target surface is oriented in the
41 // same way
42 surf_rotation_ = Acts::RotationMatrix3::Zero();
43 // u direction along +Y
44 surf_rotation_(1, 0) = 1;
45 // v direction along +Z
46 surf_rotation_(2, 1) = 1;
47 // w direction along +X
48 surf_rotation_(0, 2) = 1;
49
50 Acts::Vector3 target_pos(0., 0., 0.);
51 Acts::Translation3 target_translation(target_pos);
52 Acts::Transform3 target_transform(target_translation * surf_rotation_);
53
54 // Unbounded surface
55 target_surface_ =
56 Acts::Surface::makeShared<Acts::PlaneSurface>(target_transform);
57
58 // Custom transformation of the interpolated bfield map
59 bool debug_transform = false;
60 auto transform_pos = [this, debug_transform](const Acts::Vector3& pos_) {
61 Acts::Vector3 rot_pos;
62 rot_pos(0) = pos_(1);
63 rot_pos(1) = pos_(2);
64 rot_pos(2) = pos_(0) + DIPOLE_OFFSET;
65
66 // Systematic effect
67 rot_pos(0) += this->map_offset_[0];
68 rot_pos(1) += this->map_offset_[1];
69 rot_pos(2) += this->map_offset_[2];
70
71 // Apply A rotation around the center of the magnet. (I guess offset first
72 // and then rotation)
73
74 if (debug_transform) {
75 std::cout << "PF::DEFAULT3 TRANSFORM" << std::endl;
76 std::cout << "PF::Check:: transforming Pos" << std::endl;
77 std::cout << pos_ << std::endl;
78 std::cout << "TO" << std::endl;
79 std::cout << rot_pos << std::endl;
80 }
81
82 return rot_pos;
83 };
84
85 Acts::RotationMatrix3 rotation = Acts::RotationMatrix3::Identity();
86 double scale = 1.;
87
88 auto transform_b_field = [rotation, scale, debug_transform](
89 const Acts::Vector3& field,
90 const Acts::Vector3& /*pos_*/) {
91 // Rotate the field in tracking coordinates
92 Acts::Vector3 rot_field;
93 rot_field(0) = field(2);
94 rot_field(1) = field(0);
95 rot_field(2) = field(1);
96
97 // Scale the field
98 rot_field = scale * rot_field;
99
100 // Rotate the field
101 rot_field = rotation * rot_field;
102
103 // A distortion scaled by position.
104 if (debug_transform) {
105 std::cout << "PF::DEFAULT3 TRANSFORM" << std::endl;
106 std::cout << "PF::Check:: transforming" << std::endl;
107 std::cout << field << std::endl;
108 std::cout << "TO" << std::endl;
109 std::cout << rot_field << std::endl;
110 }
111
112 return rot_field;
113 };
114
115 // Setup a interpolated bfield map
116 const auto map = std::make_shared<InterpolatedMagneticField3>(
117 loadDefaultBField(field_map_,
118 // default_transformPos,
119 // default_transformBField));
120 transform_pos, transform_b_field));
121
122 auto acts_logging_level = Acts::Logging::FATAL;
123 if (debug_acts_) acts_logging_level = Acts::Logging::VERBOSE;
124
125 // Setup the steppers
126 const auto stepper = Acts::EigenStepper<>{map};
127 const auto const_stepper = Acts::EigenStepper<>{const_b_field};
128 const auto multi_stepper = Acts::MultiEigenStepperLoop{map};
129
130 // Setup the navigator
131 Acts::Navigator::Config nav_cfg{geometry().getTG()};
132 nav_cfg.resolveMaterial = true;
133 nav_cfg.resolvePassive = true;
134 nav_cfg.resolveSensitive = true;
135 const Acts::Navigator navigator(nav_cfg);
136
137 // Setup the propagators
138 propagator_ =
139 const_b_field_
140 ? std::make_unique<CkfPropagator>(const_stepper, navigator)
141 : std::make_unique<CkfPropagator>(
142 stepper, navigator,
143 Acts::getDefaultLogger("ACTS_PROP", acts_logging_level));
144
145 // Setup the finder / fitters
146 ckf_ = std::make_unique<std::decay_t<decltype(*ckf_)>>(
147 *propagator_, Acts::getDefaultLogger("CKF", acts_logging_level));
148 trk_extrap_ = std::make_shared<std::decay_t<decltype(*trk_extrap_)>>(
149 *propagator_, geometryContext(), magneticFieldContext());
150}
151
153 eventnr_++;
154 // get the tracking geometry from conditions
155 auto tg{geometry()};
156
157 // TODO use global variable instead and call clear;
158
159 std::vector<ldmx::Track> tracks;
160
161 auto start = std::chrono::high_resolution_clock::now();
162
163 nevents_++;
164
165 auto logging_level = Acts::Logging::DEBUG;
166 ACTS_LOCAL_LOGGER(
167 Acts::getDefaultLogger("LDMX Tracking Geometry Maker", logging_level));
168
169 // Move this at the start of the producer
170 Acts::PropagatorOptions<Acts::StepperPlainOptions,
171 Acts::NavigatorPlainOptions, ActionList, AbortList>
172 propagator_options(geometryContext(), magneticFieldContext());
173
174 propagator_options.pathLimit = std::numeric_limits<double>::max();
175 // Activate loop protection at some pt value
176 propagator_options.loopProtection = false;
177 //(startParameters.transverseMomentum() < cfg.ptLoopers);
178
179 // Switch the material interaction on/off & eventually into logging mode
180 auto& m_interactor =
181 propagator_options.actionList.get<Acts::MaterialInteractor>();
182 m_interactor.multipleScattering = true;
183 m_interactor.energyLoss = true;
184 m_interactor.recordInteractions = false;
185
186 // The logger can be switched to sterile, e.g. for timing logging
187 auto& s_logger =
188 propagator_options.actionList.get<Acts::detail::SteppingLogger>();
189 s_logger.sterile = true;
190 // Set a maximum step size
191 propagator_options.stepping.maxStepSize =
192 propagator_step_size_ * Acts::UnitConstants::mm;
193 propagator_options.maxSteps = propagator_max_steps_;
194
195 // #######################//
196 // Kalman Filter algorithm//
197 // #######################//
198
199 // Step 1 - Form the source links
200
201 // a) Loop over the sim Hits
202
203 auto setup = std::chrono::high_resolution_clock::now();
204 profiling_map_["setup"] +=
205 std::chrono::duration<double, std::milli>(setup - start).count();
206
207 const std::vector<ldmx::Measurement> measurements =
208 event.getCollection<ldmx::Measurement>(measurement_collection_,
209 input_pass_name_);
210
211 // check if SimParticleMap is available for truth matching
212 std::shared_ptr<tracking::sim::TruthMatchingTool> truth_matching_tool =
213 nullptr;
214 std::map<int, ldmx::SimParticle> particle_map;
215
216 if (event.exists("SimParticles", sim_particles_event_passname_)) {
217 ldmx_log(debug) << "Setting up track truth matching tool";
218 particle_map = event.getMap<int, ldmx::SimParticle>(
219 "SimParticles", sim_particles_event_passname_);
220 truth_matching_tool = std::make_shared<tracking::sim::TruthMatchingTool>(
221 particle_map, measurements);
222 }
223
224 // The mapping between the geometry identifier
225 // and the IndexsourceLink that points to the hit
226 const auto geo_id_sl_map = makeGeoIdSourceLinkMap(tg, measurements);
227
228 auto hits = std::chrono::high_resolution_clock::now();
229 profiling_map_["hits"] +=
230 std::chrono::duration<double, std::milli>(hits - setup).count();
231
232 // ============ Setup the CKF ============
233
234 // Retrieve the seeds
235 const std::vector<ldmx::Track> seed_tracks =
236 event.getCollection<ldmx::Track>(seed_coll_name_, input_pass_name_);
237
238 ldmx_log(info) << "Number of " << seed_coll_name_
239 << " seed tracks = " << seed_tracks.size();
240
241 if (seed_tracks.empty()) {
242 std::vector<ldmx::Track> empty;
243 ldmx_log(info) << "No seed tracks, returning...";
244 event.add(out_trk_collection_, empty);
245 return;
246 }
247
248 // Run the CKF on each seed and produce a track candidate
249 std::vector<Acts::BoundTrackParameters> start_parameters;
250
251 ldmx_log(debug) << "Transform the seed track to bound parameters";
252 int seed_track_index{0};
253 for (auto& seed : seed_tracks) {
254 // Transform the seed track to bound parameters
255 std::shared_ptr<Acts::PerigeeSurface> perigee_surface =
256 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
257 seed.getPerigeeX(), seed.getPerigeeY(), seed.getPerigeeZ()));
258
259 Acts::BoundVector param_vec;
260 param_vec << seed.getD0(), seed.getZ0(), seed.getPhi(), seed.getTheta(),
261 seed.getQoP(), seed.getT();
262
263 Acts::BoundSquareMatrix cov_mat =
264 tracking::sim::utils::unpackCov(seed.getPerigeeCov());
265
266 ldmx_log(debug) << " For seed index_ = " << seed_track_index
267 << ": Perigee X / Y / Z = " << seed.getPerigeeX() << " / "
268 << seed.getPerigeeY() << " / " << seed.getPerigeeZ()
269 << ", D0 = " << param_vec[0] << ", Z0 = " << param_vec[1]
270 << ", Phi = " << param_vec[2]
271 << ", Theta = " << param_vec[3]
272 << ", QoP = " << param_vec[4]
273 << ", Time = " << param_vec[5];
274
275 ldmx_log(debug) << " Cov matrix diagonal (" << cov_mat(0, 0) << ", "
276 << cov_mat(1, 1) << ", " << cov_mat(2, 2) << ")";
277
278 // need to set particle hypothesis...set to electron for now...
279 auto part_hypo{Acts::SinglyChargedParticleHypothesis::electron()};
280 start_parameters.push_back(Acts::BoundTrackParameters(
281 perigee_surface, param_vec, cov_mat, part_hypo));
282
283 // This is a global variable for performance checks
284 nseeds_++;
285 // This is just to index_ the seed we are looking at
286 seed_track_index++;
287 } // loop on seeds
288
289 auto seeds = std::chrono::high_resolution_clock::now();
290 profiling_map_["seeds"] +=
291 std::chrono::duration<double, std::milli>(seeds - hits).count();
292
293 Acts::GainMatrixUpdater kf_updater;
294
295 // configuration for the measurement selector. Empty geometry identifier means
296 // applicable to all the detector elements
297
298 Acts::MeasurementSelector::Config measurement_selector_cfg = {
299 // global default: no chi2 cut, only one measurement per surface
300 {Acts::GeometryIdentifier(), {{}, {outlier_pval_}, {1u}}},
301 };
302
303 Acts::MeasurementSelector meas_sel{measurement_selector_cfg};
304
305 tracking::sim::LdmxMeasurementCalibrator calibrator{measurements};
306
307 Acts::CombinatorialKalmanFilterExtensions<TrackContainer> ckf_extensions;
308
309 if (use1_dmeasurements_) {
310 ckf_extensions.calibrator
312 Acts::VectorMultiTrajectory>>(&calibrator);
313 } else {
314 ckf_extensions.calibrator
316 Acts::VectorMultiTrajectory>>(&calibrator);
317 }
318
319 ckf_extensions.updater.connect<
320 &Acts::GainMatrixUpdater::operator()<Acts::VectorMultiTrajectory>>(
321 &kf_updater);
322
323 ckf_extensions.measurementSelector
324 .connect<&Acts::MeasurementSelector::select<Acts::VectorMultiTrajectory>>(
325 &meas_sel);
326
327 ldmx_log(debug) << "SourceLinkAccessor...";
328
329 // Create source link accessor and connect delegate
330 struct SourceLinkAccIt {
331 using BaseIt = decltype(geo_id_sl_map.begin());
332 BaseIt it_;
333
334#pragma GCC diagnostic push
335#pragma GCC diagnostic ignored "-Wunused-local-typedefs"
336
337 using difference_type = typename BaseIt::difference_type;
338 using iterator_category = typename BaseIt::iterator_category;
339 // using value_type = typename BaseIt::value_type::second_type;
340 using value_type = Acts::SourceLink;
341 using pointer = typename BaseIt::pointer;
342 using reference = value_type&;
343#pragma GCC diagnostic pop
344
345 SourceLinkAccIt& operator++() {
346 ++it_;
347 return *this;
348 }
349 bool operator==(const SourceLinkAccIt& other) const {
350 return it_ == other.it_;
351 }
352 bool operator!=(const SourceLinkAccIt& other) const {
353 return !(*this == other);
354 }
355 // const value_type& operator*() const { return it->second; }
356
357 // by value
358 value_type operator*() const { return value_type{it_->second}; }
359 };
360
361 auto source_link_accessor = [&](const Acts::Surface& surface)
362 -> std::pair<SourceLinkAccIt, SourceLinkAccIt> {
363 auto [begin, end] = geo_id_sl_map.equal_range(surface.geometryId());
364 return {SourceLinkAccIt{begin}, SourceLinkAccIt{end}};
365 };
366
367 Acts::SourceLinkAccessorDelegate<SourceLinkAccIt>
368 source_link_accessor_delegate;
369 source_link_accessor_delegate
370 .connect<&decltype(source_link_accessor)::operator(),
371 decltype(source_link_accessor)>(&source_link_accessor);
372
373 ldmx_log(debug) << "Setting up surfaces...";
374
375 std::shared_ptr<const Acts::PerigeeSurface> origin_surface =
376 Acts::Surface::makeShared<Acts::PerigeeSurface>(
377 Acts::Vector3(0., 0., 0.));
378
379 ldmx_log(debug) << "About to run CKF...";
380
381 // run the CKF for all initial track states
382 auto ckf_setup = std::chrono::high_resolution_clock::now();
383 profiling_map_["ckf_setup"] +=
384 std::chrono::duration<double, std::milli>(ckf_setup - seeds).count();
385
386 Acts::VectorTrackContainer vtc;
387 Acts::VectorMultiTrajectory mtj;
388 Acts::TrackContainer tc{vtc, mtj};
389
390 // The number of track candidates (i.e. startParameters.size()) is always
391 // the same as the number of seed tracks
392 ldmx_log(debug) << "Loop on the track candidates";
393 for (size_t track_id = 0u; track_id < start_parameters.size(); ++track_id) {
394 ldmx_log(debug) << "---------------------------";
395 ldmx_log(debug) << "Candidate Track ID = " << track_id;
396 // Define the CKF options here:
397 const Acts::CombinatorialKalmanFilterOptions<SourceLinkAccIt,
398 TrackContainer>
399 ckf_options(TrackingGeometryUser::geometryContext(),
400 TrackingGeometryUser::magneticFieldContext(),
401 TrackingGeometryUser::calibrationContext(),
402 source_link_accessor_delegate, ckf_extensions,
403 propagator_options, true /* multiple scattering */,
404 false /* energy loss */);
405
406 ldmx_log(debug) << " Checking options: multiple scattering = "
407 << ckf_options.multipleScattering
408 << " energy loss = " << ckf_options.energyLoss;
409 auto results =
410 ckf_->findTracks(start_parameters.at(track_id), ckf_options, tc);
411
412 auto start_params = start_parameters.at(track_id).parameters().transpose();
413 ldmx_log(debug)
414 << " Checking CKF success for track candidate with params: "
415 << " D0 = " << start_params[0] << " Z0 = " << start_params[1]
416 << ", Phi = " << start_params[2] << " Theta = " << start_params[3]
417 << ", QoP = " << start_params[4] << " Time = " << start_params[5];
418 if (not results.ok()) {
419 ldmx_log(debug) << " CKF failed!";
420 continue;
421 } else {
422 ldmx_log(debug) << " CKF succeded!";
423 }
424
425 auto& tracks_from_seed = results.value();
426 if (tracks_from_seed.size() != 1) {
427 ldmx_log(info) << " tracksFromSeed.size = " << tracks_from_seed.size();
428 }
429 // For now it seems this loop is only looping on a single element
430 for (auto& track : tracks_from_seed) {
431 // do the track smoothing...this is not done in the CKF code anymore
432 Acts::smoothTrack(geometryContext(), track); // from TrackHelpers
433 // make the empty ldmx::Track() and track state at target
434 ldmx::Track trk = ldmx::Track();
435 ldmx::Track::TrackState ts_at_target;
436
437 // Extrapolate to the target
438 bool success = trk_extrap_->trackStateAtSurface(
439 track, target_surface_, ts_at_target, ldmx::TrackStateType::AtTarget);
440 // Check if the extrapolation to target succeeded
441 if (success) {
442 ldmx_log(debug) << " Successfully obtained TrackState at target";
443
444 ldmx_log(debug) << " Parameters At Target: Loc0 = "
445 << ts_at_target.params_[0] << ", Loc1 "
446 << ts_at_target.params_[1]
447 << ", phi = " << ts_at_target.params_[2]
448 << ", theta = " << ts_at_target.params_[3] << ", QoP "
449 << ts_at_target.params_[4];
450
451 trk.addTrackState(ts_at_target);
452 } else {
453 ldmx_log(info) << " Could not extrapolate to target! nhits = "
454 << track.nMeasurements() << " Printing track states: ";
455 for (const auto ts : track.trackStatesReversed()) {
456 if (ts.hasSmoothed()) {
457 ldmx_log(info) << " Track momentum for smoothed track state = "
458 << 1 / ts.smoothed()[Acts::eBoundQOverP];
459 ldmx_log(info) << " Parameters: " << ts.smoothed().transpose();
460 } else {
461 ldmx_log(info) << " Track state not smoothed!";
462 }
463 }
464 ldmx_log(info) << " ...skipping this track candidate...";
465 continue;
466 }
467
468 // get the BoundTrackParameters at the target
469 // ...use to fill in the Acts::TrackProxy object
470 // This isn't really necessary, since we can take
471 // most everything for making the ldmx::track
472 // from tsAtTarget...maybe useful for something?
473 // -->one thing this does is allow Acts to
474 // calculate the momentum 3-vector for you
475 Acts::BoundTrackParameters bound_state_at_target =
476 tracking::sim::utils::btp(ts_at_target, target_surface_, 11);
477 track.setReferenceSurface(target_surface_);
478 track.parameters() = bound_state_at_target.parameters();
479
480 // ldmx_log(debug) << " typeid(track).name() = " << typeid(track).name();
481 // These are the parameters at the target surface
482 const Acts::BoundVector& track_pars = track.parameters();
483 // const Acts::BoundMatrix& trk_cov = track.covariance();
484
485 ldmx_log(debug) << " Track parameters (at target): Loc0 = "
486 << track_pars[Acts::eBoundLoc0]
487 << ", Loc1 = " << track_pars[Acts::eBoundLoc1]
488 << ", Theta = " << track_pars[Acts::eBoundTheta]
489 << ", Phi = " << track_pars[Acts::eBoundPhi]
490 << ", chi2 = " << track.chi2()
491 << ", nHits = " << track.nMeasurements();
492
493 // the target...it's not really perigee anymore.
494 trk.setPerigeeLocation(0, 0, 0);
495 trk.setPerigeeParameters(ts_at_target.params_);
496 trk.setPerigeeCov(ts_at_target.cov_);
497
498 trk.setChi2(track.chi2());
499 trk.setNhits(track.nMeasurements());
500 // trk.setNdf(track.nDoF());
501 // TODO Switch back to nDoF when Acts is fixed.
502 trk.setNdf(track.nMeasurements() - 5);
503 trk.setNsharedHits(track.nSharedHits());
504
505 ldmx_log(debug) << " Track momentum: px = " << track.momentum()[0]
506 << " py = " << track.momentum()[1]
507 << " pz = " << track.momentum()[2];
508 trk.setMomentum(track.momentum()[0], track.momentum()[1],
509 track.momentum()[2]);
510
511 // At least min_hits hits and p > 50 MeV
512 if ((trk.getNhits() <= min_hits_) || (abs(1. / trk.getQoP()) <= 0.05)) {
513 ldmx_log(debug)
514 << " > Track candidate did NOT meet the requirements: Nhits = "
515 << trk.getNhits() << " and p = " << abs(1. / trk.getQoP())
516 << " GeV";
517 // No point of doing the extrapolations if the track candidate anyway
518 // wont be kept
519 continue;
520 }
521
522 // Add measurements to the final track
523 ldmx_log(debug) << " Add measurements to the final track from "
524 << track.nTrackStates() << " TrackStates with "
525 << track.nMeasurements() << " measurements";
526
527 int trk_state_index{0};
528 for (const auto ts : track.trackStatesReversed()) {
529 // Check TrackStates Quality
530 ldmx_log(debug) << " Checking Track State index_ = "
531 << trk_state_index << " at location "
532 << ts.referenceSurface()
533 .transform(geometryContext())
534 .translation()
535 .transpose();
536
537 if (ts.hasSmoothed()) {
538 ldmx_log(debug) << " Smoothed track parameters: "
539 << ts.smoothed().transpose();
540 // ldmx_log(debug) << " Smoothed covariance mtx:\n" <<
541 // ts.smoothedCovariance();
542 }
543
544 // Check if the track state is a measurement
545 auto type_flags = ts.typeFlags();
546
547 if (type_flags.test(Acts::TrackStateFlag::MeasurementFlag) &&
548 ts.hasUncalibratedSourceLink()) {
550 ts.getUncalibratedSourceLink()
551 .template get<acts_examples::IndexSourceLink>();
552
553 ldmx::Measurement ldmx_meas = measurements.at(sl.index());
554 ldmx_log(debug) << " Adding measurement to ldmx::track with "
555 "source link index_ = "
556 << sl.index();
557 ldmx_log(trace) << " Measurement:\n" << ldmx_meas;
558 trk.addMeasurementIndex(sl.index());
559 } else {
560 ldmx_log(debug) << " This TrackState is not a measurement";
561 }
562 trk_state_index++;
563 }
564
565 ldmx_log(debug) << " Starting extrapolations";
566 // Extrapolations
567 // To ECAL
568 const double ecal_scoring_plane = 240.5;
569 Acts::Vector3 pos(ecal_scoring_plane, 0., 0.);
570 Acts::Translation3 surf_translation(pos);
571 Acts::Transform3 surf_transform(surf_translation * surf_rotation_);
572 const std::shared_ptr<Acts::PlaneSurface> ecal_surface =
573 Acts::Surface::makeShared<Acts::PlaneSurface>(surf_transform);
574
575 // Beam Origin unbounded surface
576 const std::shared_ptr<Acts::Surface> beam_origin_surface =
577 tracking::sim::utils::unboundSurface(-700);
578
579 if (tagger_tracking_) {
580 ldmx_log(debug) << " Beam Origin Extrapolation";
581 ldmx::Track::TrackState ts_at_beam_origin;
582 success = trk_extrap_->trackStateAtSurface(
583 track, beam_origin_surface, ts_at_beam_origin,
584 ldmx::TrackStateType::AtBeamOrigin);
585
586 if (success) {
587 trk.addTrackState(ts_at_beam_origin);
588 ldmx_log(debug)
589 << " Successfully obtained TrackState at beam origin";
590 }
591 }
592
593 // Recoil Extrapolation to ECAL only
594 if (!tagger_tracking_) {
595 ldmx_log(debug) << " Ecal Extrapolation";
596 ldmx::Track::TrackState ts_at_ecal;
597 success = trk_extrap_->trackStateAtSurface(
598 track, ecal_surface, ts_at_ecal, ldmx::TrackStateType::AtECAL);
599
600 if (success) {
601 trk.addTrackState(ts_at_ecal);
602 ldmx_log(debug) << " Successfully obtained TrackState at Ecal";
603 ldmx_log(debug) << " Parameters At Ecal: Loc0 = "
604 << ts_at_ecal.params_[0]
605 << ", Loc1 = " << ts_at_ecal.params_[1]
606 << ", phi = " << ts_at_ecal.params_[2]
607 << ", theta = " << ts_at_ecal.params_[3]
608 << ", QoP = " << ts_at_ecal.params_[4];
609 } else {
610 ldmx_log(info) << " Could not extrapolate to ECAL!! Please check "
611 "the track states";
612 }
613 }
614
615 // Truth matching
616 if (truth_matching_tool) {
617 auto truth_info = truth_matching_tool->truthMatch(trk);
618 trk.setTrackID(truth_info.track_id_);
619 trk.setPdgID(truth_info.pdg_id_);
620 trk.setTruthProb(truth_info.truth_prob_);
621 }
622
623 // Adding the track candidate to the track collection
624 ldmx_log(debug)
625 << " > Adding the track candidate to the track collection";
626 tracks.push_back(trk);
627 ntracks_++;
628 } // // loop on tracksFromSeed (which usually has 1 element)
629 } // loop seed track parameters (i.e. track candidates)
630
631 ldmx_log(info) << "Number of CKF tracks " << tracks.size();
632
633 auto ckf_run = std::chrono::high_resolution_clock::now();
634 profiling_map_["ckf_run"] +=
635 std::chrono::duration<double, std::milli>(ckf_run - ckf_setup).count();
636
637 // Calculating Shared Hits
638 auto shared_hits = computeSharedHits(
639 tracks, measurements, tg, tracking::sim::utils::sourceLinkHash,
640 tracking::sim::utils::sourceLinkEquality);
641 for (std::size_t i_track = 0; i_track < shared_hits.size(); ++i_track) {
642 tracks[i_track].setNsharedHits(shared_hits[i_track].size());
643 for (auto idx : shared_hits[i_track]) {
644 tracks[i_track].addSharedIndex(idx);
645 }
646 }
647
648 auto result_loop = std::chrono::high_resolution_clock::now();
649 profiling_map_["result_loop"] +=
650 std::chrono::duration<double, std::milli>(result_loop - ckf_run).count();
651
652 // Add the tracks to the event
653 event.add(out_trk_collection_, tracks);
654
655 auto end = std::chrono::high_resolution_clock::now();
656 // long long microseconds =
657 // std::chrono::duration_cast<std::chrono::microseconds>(end-start).count();
658 auto diff = end - start;
659 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
660}
661
663 if (use1_dmeasurements_)
664 ldmx_log(debug) << "Use1Dmeasurements = " << std::boolalpha
665 << use1_dmeasurements_;
666 if (remove_stereo_)
667 ldmx_log(debug) << "Remove_stereo = " << std::boolalpha << remove_stereo_;
668}
669
671 ldmx_log(info) << "found " << ntracks_ << " tracks / " << nseeds_
672 << " nseeds";
673 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(1)
674 << processing_time_ / nevents_ << " ms";
675 ldmx_log(info) << "Breakdown::";
676 ldmx_log(info) << "setup Avg Time/Event = " << std::fixed
677 << std::setprecision(3) << profiling_map_["setup"] / nevents_
678 << " ms";
679 ldmx_log(info) << "hits Avg Time/Event = " << std::fixed
680 << std::setprecision(2) << profiling_map_["hits"] / nevents_
681 << " ms";
682 ldmx_log(info) << "seeds Avg Time/Event = " << std::fixed
683 << std::setprecision(3) << profiling_map_["seeds"] / nevents_
684 << " ms";
685 ldmx_log(info) << "cf_setup Avg Time/Event = " << std::fixed
686 << std::setprecision(3)
687 << profiling_map_["ckf_setup"] / nevents_ << " ms";
688 ldmx_log(info) << "ckf_run Avg Time/Event = " << std::fixed
689 << std::setprecision(3) << profiling_map_["ckf_run"] / nevents_
690 << " ms";
691 ldmx_log(info) << "result_loop Avg Time/Event = " << std::fixed
692 << std::setprecision(1)
693 << profiling_map_["result_loop"] / nevents_ << " ms";
694}
695
697 dumpobj_ = parameters.get<bool>("dumpobj", 0);
698 pionstates_ = parameters.get<int>("pionstates", 0);
699
700 bfield_ = parameters.get<double>("bfield", -1.5);
701 const_b_field_ = parameters.get<bool>("const_b_field", false);
702 field_map_ = parameters.get<std::string>("field_map");
703 propagator_step_size_ = parameters.get<double>("propagator_step_size", 200.);
704 propagator_max_steps_ = parameters.get<int>("propagator_maxSteps", 10000);
705 measurement_collection_ = parameters.get<std::string>(
706 "measurement_collection", "TaggerMeasurements");
707 outlier_pval_ = parameters.get<double>("outlier_pval_", 3.84);
708
709 debug_acts_ = parameters.get<bool>("debug_acts", false);
710
711 remove_stereo_ = parameters.get<bool>("remove_stereo", false);
712 use1_dmeasurements_ = parameters.get<bool>("use1Dmeasurements", true);
713 min_hits_ = parameters.get<int>("min_hits", 7);
714
715 // Ckf specific options
716 use_extrapolate_location_ =
717 parameters.get<bool>("use_extrapolate_location", true);
718 extrapolate_location_ =
719 parameters.get<std::vector<double>>("extrapolate_location", {0., 0., 0.});
720 use_seed_perigee_ = parameters.get<bool>("use_seed_perigee", false);
721
722 // seeds from the event
723 seed_coll_name_ = parameters.get<std::string>("seed_coll_name", "seedTracks");
724
725 sim_particles_event_passname_ =
726 parameters.get<std::string>("sim_particles_event_passname");
727
728 // output track collection
729 out_trk_collection_ =
730 parameters.get<std::string>("out_trk_collection", "Tracks");
731
732 // keep track on which system tracking is running
733 tagger_tracking_ = parameters.get<bool>("taggerTracking", true);
734
735 // BField Systematics
736 map_offset_ =
737 parameters.get<std::vector<double>>("map_offset_", {0., 0., 0.});
738
739 input_pass_name_ = parameters.get<std::string>("input_pass_name");
740}
741
742auto CKFProcessor::makeGeoIdSourceLinkMap(
744 const std::vector<ldmx::Measurement>& measurements)
745 -> std::unordered_multimap<Acts::GeometryIdentifier,
747 std::unordered_multimap<Acts::GeometryIdentifier,
749 geo_id_sl_map;
750
751 ldmx_log(debug) << "The makeGeoIdSourceLinkMap has " << measurements.size()
752 << " measurements";
753
754 // Check the hits associated to the surfaces
755 for (unsigned int i_meas = 0; i_meas < measurements.size(); i_meas++) {
756 ldmx::Measurement meas = measurements.at(i_meas);
757 unsigned int layerid = meas.getLayerID();
758
759 const Acts::Surface* hit_surface = tg.getSurface(layerid);
760
761 if (hit_surface) {
762 // Transform the ldmx space point from global to local and store the
763 // information
764
765 acts_examples::IndexSourceLink idx_sl(hit_surface->geometryId(), i_meas);
766 // mg aug 2024 ... these don't print statements
767 // don't compile using v36 in Acts...figure out later
768 /*
769 ldmx_log(debug)
770 << "Insert measurement on surface located at::"
771 << hit_surface->transform(geometry_context()).translation();
772 ldmx_log(debug) << "and geoId::" << hit_surface->geometryId();
773
774 ldmx_log(debug) << "Surface info::"
775 << std::tie(*hit_surface, geometry_context());
776 */
777 geo_id_sl_map.insert(std::make_pair(hit_surface->geometryId(), idx_sl));
778
779 } else
780 ldmx_log(debug) << getName() << "::HIT " << i_meas << " at layer_"
781 << (measurements.at(i_meas)).getLayerID()
782 << " is not associated to any surface?!";
783 }
784
785 return geo_id_sl_map;
786}
787
788template <typename geometry_t, typename source_link_hash_t,
789 typename source_link_equality_t>
790std::vector<std::vector<std::size_t>> CKFProcessor::computeSharedHits(
791 std::vector<ldmx::Track> tracks, std::vector<ldmx::Measurement> meas_coll,
792 geometry_t& tg, source_link_hash_t&& sourceLinkHash,
793 source_link_equality_t&& sourceLinkEquality) const {
794 auto measurement_index_map =
795 std::unordered_map<Acts::SourceLink, std::size_t, source_link_hash_t,
796 source_link_equality_t>(0, sourceLinkHash,
797 sourceLinkEquality);
798
799 std::vector<std::vector<std::size_t>> measurements_per_track;
800 boost::container::flat_map<std::size_t,
801 boost::container::flat_set<std::size_t>>
802 tracks_per_measurement;
803 std::vector<std::size_t> shared_measurements_per_track;
804 auto number_of_tracks = 0;
805
806 // Iterate through all input tracks, collect their properties like measurement
807 // count and chi2 and fill the measurement map in order to relate tracks to
808 // each other if they have shared hits.
809 for (const auto& track : tracks) {
810 // Kick out tracks that do not fulfill our initial requirements
811 // if (track.getNhits() < n_measurements_min_) {
812 // continue;
813 // }
814
815 std::vector<std::size_t> measurements;
816 for (auto imeas : track.getMeasurementsIdxs()) {
817 auto meas = meas_coll.at(imeas);
818 const Acts::Surface* hit_surface = tg.getSurface(meas.getLayerID());
819 // Store the index_ source link
820 acts_examples::IndexSourceLink idx_sl(hit_surface->geometryId(), imeas);
821 Acts::SourceLink source_link = Acts::SourceLink(idx_sl);
822
823 auto emplace = measurement_index_map.try_emplace(
824 source_link, measurement_index_map.size());
825 measurements.push_back(emplace.first->second);
826 }
827
828 measurements_per_track.push_back(std::move(measurements));
829
830 ++number_of_tracks;
831 }
832
833 // Now we relate measurements to tracks
834 for (std::size_t i_track = 0; i_track < number_of_tracks; ++i_track) {
835 for (auto i_measurement : measurements_per_track[i_track]) {
836 tracks_per_measurement[i_measurement].insert(i_track);
837 }
838 }
839
840 // Finally, we can accumulate the number of shared measurements per track
841 shared_measurements_per_track = std::vector<std::size_t>(number_of_tracks, 0);
842
843 std::vector<std::vector<std::size_t>> shared_measurement_idxs_per_track;
844 for (std::size_t i_track = 0; i_track < number_of_tracks; ++i_track) {
845 std::vector<std::size_t> shared_measurement_idxs;
846 for (auto i_measurement : measurements_per_track[i_track]) {
847 if (tracks_per_measurement[i_measurement].size() > 1) {
848 ++shared_measurements_per_track[i_track];
849 shared_measurement_idxs.push_back(i_measurement);
850 }
851 }
852 shared_measurement_idxs_per_track.push_back(shared_measurement_idxs);
853 }
854 return shared_measurement_idxs_per_track;
855}
856
857} // namespace reco
858} // namespace tracking
859
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
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 which represents the process under execution.
Definition Process.h:36
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
int getLayerID() const
Run-specific configuration and data stored in its own output TTree alongside the event TTree in the o...
Definition RunHeader.h:57
Class representing a simulated particle.
Definition SimParticle.h:23
Implementation of a track object.
Definition Track.h:53
void setPerigeeParameters(const std::vector< double > &par)
d_0 z_0 phi_0 theta q/p t
Definition Track.h:156
void produce(framework::Event &event) override
Run the processor.
void configure(framework::config::Parameters &parameters) override
Configure the processor using the given user specified parameters.
void onNewRun(const ldmx::RunHeader &rh) override
onNewRun is the first function called for each processor after the conditions are fully configured an...
CKFProcessor(const std::string &name, framework::Process &process)
Constructor.
int nseeds_
n seeds and n tracks
void onProcessStart() override
Callback for the EventProcessor to take any necessary action when the processing of events starts,...
void onProcessEnd() override
Callback for the EventProcessor to take any necessary action when the processing of events finishes,...
a helper base class providing some methods to shorten access to common conditions used within the tra...
void calibrate1d(const Acts::GeometryContext &, const Acts::CalibrationContext &, const Acts::SourceLink &genericSourceLink, typename traj_t::TrackStateProxy trackState) const
Find the measurement corresponding to the source link.
void calibrate(const Acts::GeometryContext &, const Acts::CalibrationContext &, const Acts::SourceLink &genericSourceLink, typename traj_t::TrackStateProxy trackState) const
Find the measurement corresponding to the source link.
The measurement calibrator can be a function or a class/struct able to retrieve the sim hits containe...