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 constBField = 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 debugTransform = false;
60 auto transformPos = [this, debugTransform](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 (debugTransform) {
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 transformBField = [rotation, scale, debugTransform](
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 (debugTransform) {
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 transformPos, transformBField));
121
122 auto acts_loggingLevel = Acts::Logging::FATAL;
123 if (debug_acts_) acts_loggingLevel = Acts::Logging::VERBOSE;
124
125 // Setup the steppers
126 const auto stepper = Acts::EigenStepper<>{map};
127 const auto const_stepper = Acts::EigenStepper<>{constBField};
128 const auto multi_stepper = Acts::MultiEigenStepperLoop{map};
129
130 // Setup the navigator
131 Acts::Navigator::Config navCfg{geometry().getTG()};
132 navCfg.resolveMaterial = true;
133 navCfg.resolvePassive = true;
134 navCfg.resolveSensitive = true;
135 const Acts::Navigator navigator(navCfg);
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_loggingLevel));
144
145 // Setup the finder / fitters
146 ckf_ = std::make_unique<std::decay_t<decltype(*ckf_)>>(
147 *propagator_, Acts::getDefaultLogger("CKF", acts_loggingLevel));
148 trk_extrap_ = std::make_shared<std::decay_t<decltype(*trk_extrap_)>>(
149 *propagator_, geometry_context(), magnetic_field_context());
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 loggingLevel = Acts::Logging::DEBUG;
166 ACTS_LOCAL_LOGGER(
167 Acts::getDefaultLogger("LDMX Tracking Geometry Maker", loggingLevel));
168
169 // Move this at the start of the producer
170 Acts::PropagatorOptions<Acts::StepperPlainOptions,
171 Acts::NavigatorPlainOptions, ActionList, AbortList>
172 propagator_options(geometry_context(), magnetic_field_context());
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& mInteractor =
181 propagator_options.actionList.get<Acts::MaterialInteractor>();
182 mInteractor.multipleScattering = true;
183 mInteractor.energyLoss = true;
184 mInteractor.recordInteractions = false;
185
186 // The logger can be switched to sterile, e.g. for timing logging
187 auto& sLogger =
188 propagator_options.actionList.get<Acts::detail::SteppingLogger>();
189 sLogger.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_maxSteps_;
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> truthMatchingTool = nullptr;
213 std::map<int, ldmx::SimParticle> particleMap;
214
215 if (event.exists("SimParticles")) {
216 ldmx_log(debug) << "Setting up track truth matching tool";
217 particleMap = event.getMap<int, ldmx::SimParticle>("SimParticles");
218 truthMatchingTool = std::make_shared<tracking::sim::TruthMatchingTool>(
219 particleMap, measurements);
220 }
221
222 // The mapping between the geometry identifier
223 // and the IndexsourceLink that points to the hit
224 const auto geoId_sl_map = makeGeoIdSourceLinkMap(tg, measurements);
225
226 auto hits = std::chrono::high_resolution_clock::now();
227 profiling_map_["hits"] +=
228 std::chrono::duration<double, std::milli>(hits - setup).count();
229
230 // ============ Setup the CKF ============
231
232 // Retrieve the seeds
233 const std::vector<ldmx::Track> seed_tracks =
234 event.getCollection<ldmx::Track>(seed_coll_name_, input_pass_name_);
235
236 ldmx_log(info) << "Number of " << seed_coll_name_
237 << " seed tracks = " << seed_tracks.size();
238
239 if (seed_tracks.empty()) {
240 std::vector<ldmx::Track> empty;
241 ldmx_log(info) << "No seed tracks, returning...";
242 event.add(out_trk_collection_, empty);
243 return;
244 }
245
246 // Run the CKF on each seed and produce a track candidate
247 std::vector<Acts::BoundTrackParameters> startParameters;
248
249 ldmx_log(debug) << "Transform the seed track to bound parameters";
250 int seed_track_index{0};
251 for (auto& seed : seed_tracks) {
252 // Transform the seed track to bound parameters
253 std::shared_ptr<Acts::PerigeeSurface> perigeeSurface =
254 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
255 seed.getPerigeeX(), seed.getPerigeeY(), seed.getPerigeeZ()));
256
257 Acts::BoundVector paramVec;
258 paramVec << seed.getD0(), seed.getZ0(), seed.getPhi(), seed.getTheta(),
259 seed.getQoP(), seed.getT();
260
261 Acts::BoundSquareMatrix covMat =
262 tracking::sim::utils::unpackCov(seed.getPerigeeCov());
263
264 ldmx_log(debug) << " For seed index = " << seed_track_index
265 << ": Perigee X / Y / Z = " << seed.getPerigeeX() << " / "
266 << seed.getPerigeeY() << " / " << seed.getPerigeeZ()
267 << ", D0 = " << paramVec[0] << ", Z0 = " << paramVec[1]
268 << ", Phi = " << paramVec[2] << ", Theta = " << paramVec[3]
269 << ", QoP = " << paramVec[4] << ", Time = " << paramVec[5];
270
271 ldmx_log(debug) << " Cov matrix diagonal (" << covMat(0, 0) << ", "
272 << covMat(1, 1) << ", " << covMat(2, 2) << ")";
273
274 // need to set particle hypothesis...set to electron for now...
275 auto partHypo{Acts::SinglyChargedParticleHypothesis::electron()};
276 startParameters.push_back(
277 Acts::BoundTrackParameters(perigeeSurface, paramVec, covMat, partHypo));
278
279 // This is a global variable for performance checks
280 nseeds_++;
281 // This is just to index the seed we are looking at
282 seed_track_index++;
283 } // loop on seeds
284
285 auto seeds = std::chrono::high_resolution_clock::now();
286 profiling_map_["seeds"] +=
287 std::chrono::duration<double, std::milli>(seeds - hits).count();
288
289 Acts::GainMatrixUpdater kfUpdater;
290
291 // configuration for the measurement selector. Empty geometry identifier means
292 // applicable to all the detector elements
293
294 Acts::MeasurementSelector::Config measurementSelectorCfg = {
295 // global default: no chi2 cut, only one measurement per surface
296 {Acts::GeometryIdentifier(), {{}, {outlier_pval_}, {1u}}},
297 };
298
299 Acts::MeasurementSelector measSel{measurementSelectorCfg};
300
301 tracking::sim::LdmxMeasurementCalibrator calibrator{measurements};
302
303 Acts::CombinatorialKalmanFilterExtensions<TrackContainer> ckf_extensions;
304
305 if (use1Dmeasurements_) {
306 ckf_extensions.calibrator
308 Acts::VectorMultiTrajectory>>(&calibrator);
309 } else {
310 ckf_extensions.calibrator
312 Acts::VectorMultiTrajectory>>(&calibrator);
313 }
314
315 ckf_extensions.updater.connect<
316 &Acts::GainMatrixUpdater::operator()<Acts::VectorMultiTrajectory>>(
317 &kfUpdater);
318
319 ckf_extensions.measurementSelector
320 .connect<&Acts::MeasurementSelector::select<Acts::VectorMultiTrajectory>>(
321 &measSel);
322
323 ldmx_log(debug) << "SourceLinkAccessor...";
324
325 // Create source link accessor and connect delegate
326 struct SourceLinkAccIt {
327 using BaseIt = decltype(geoId_sl_map.begin());
328 BaseIt it;
329
330#pragma GCC diagnostic push
331#pragma GCC diagnostic ignored "-Wunused-local-typedefs"
332
333 using difference_type = typename BaseIt::difference_type;
334 using iterator_category = typename BaseIt::iterator_category;
335 // using value_type = typename BaseIt::value_type::second_type;
336 using value_type = Acts::SourceLink;
337 using pointer = typename BaseIt::pointer;
338 using reference = value_type&;
339#pragma GCC diagnostic pop
340
341 SourceLinkAccIt& operator++() {
342 ++it;
343 return *this;
344 }
345 bool operator==(const SourceLinkAccIt& other) const {
346 return it == other.it;
347 }
348 bool operator!=(const SourceLinkAccIt& other) const {
349 return !(*this == other);
350 }
351 // const value_type& operator*() const { return it->second; }
352
353 // by value
354 value_type operator*() const { return value_type{it->second}; }
355 };
356
357 auto sourceLinkAccessor = [&](const Acts::Surface& surface)
358 -> std::pair<SourceLinkAccIt, SourceLinkAccIt> {
359 auto [begin, end] = geoId_sl_map.equal_range(surface.geometryId());
360 return {SourceLinkAccIt{begin}, SourceLinkAccIt{end}};
361 };
362
363 Acts::SourceLinkAccessorDelegate<SourceLinkAccIt> sourceLinkAccessorDelegate;
364 sourceLinkAccessorDelegate.connect<&decltype(sourceLinkAccessor)::operator(),
365 decltype(sourceLinkAccessor)>(
366 &sourceLinkAccessor);
367
368 ldmx_log(debug) << "Setting up surfaces...";
369
370 std::shared_ptr<const Acts::PerigeeSurface> origin_surface =
371 Acts::Surface::makeShared<Acts::PerigeeSurface>(
372 Acts::Vector3(0., 0., 0.));
373
374 ldmx_log(debug) << "About to run CKF...";
375
376 // run the CKF for all initial track states
377 auto ckf_setup = std::chrono::high_resolution_clock::now();
378 profiling_map_["ckf_setup"] +=
379 std::chrono::duration<double, std::milli>(ckf_setup - seeds).count();
380
381 Acts::VectorTrackContainer vtc;
382 Acts::VectorMultiTrajectory mtj;
383 Acts::TrackContainer tc{vtc, mtj};
384
385 // The number of track candidates (i.e. startParameters.size()) is always
386 // the same as the number of seed tracks
387 ldmx_log(debug) << "Loop on the track candidates";
388 for (size_t trackId = 0u; trackId < startParameters.size(); ++trackId) {
389 ldmx_log(debug) << "---------------------------";
390 ldmx_log(debug) << "Candidate Track ID = " << trackId;
391 // Define the CKF options here:
392 const Acts::CombinatorialKalmanFilterOptions<SourceLinkAccIt,
393 TrackContainer>
394 ckfOptions(TrackingGeometryUser::geometry_context(),
395 TrackingGeometryUser::magnetic_field_context(),
396 TrackingGeometryUser::calibration_context(),
397 sourceLinkAccessorDelegate, ckf_extensions,
398 propagator_options, true /* multiple scattering */,
399 false /* energy loss */);
400
401 ldmx_log(debug) << " Checking options: multiple scattering = "
402 << ckfOptions.multipleScattering
403 << " energy loss = " << ckfOptions.energyLoss;
404 auto results =
405 ckf_->findTracks(startParameters.at(trackId), ckfOptions, tc);
406
407 auto start_params = startParameters.at(trackId).parameters().transpose();
408 ldmx_log(debug)
409 << " Checking CKF success for track candidate with params: "
410 << " D0 = " << start_params[0] << " Z0 = " << start_params[1]
411 << ", Phi = " << start_params[2] << " Theta = " << start_params[3]
412 << ", QoP = " << start_params[4] << " Time = " << start_params[5];
413 if (not results.ok()) {
414 ldmx_log(debug) << " CKF failed!";
415 continue;
416 } else {
417 ldmx_log(debug) << " CKF succeded!";
418 }
419
420 auto& tracksFromSeed = results.value();
421 if (tracksFromSeed.size() != 1) {
422 ldmx_log(info) << " tracksFromSeed.size = " << tracksFromSeed.size();
423 }
424 // For now it seems this loop is only looping on a single element
425 for (auto& track : tracksFromSeed) {
426 // do the track smoothing...this is not done in the CKF code anymore
427 Acts::smoothTrack(geometry_context(), track); // from TrackHelpers
428 // make the empty ldmx::Track() and track state at target
429 ldmx::Track trk = ldmx::Track();
430 ldmx::Track::TrackState tsAtTarget;
431
432 // Extrapolate to the target
433 bool success = trk_extrap_->TrackStateAtSurface(
434 track, target_surface, tsAtTarget, ldmx::TrackStateType::AtTarget);
435 // Check if the extrapolation to target succeeded
436 if (success) {
437 ldmx_log(debug) << " Successfully obtained TrackState at target";
438
439 ldmx_log(debug) << " Parameters At Target: Loc0 = "
440 << tsAtTarget.params[0] << ", Loc1 "
441 << tsAtTarget.params[1]
442 << ", phi = " << tsAtTarget.params[2]
443 << ", theta = " << tsAtTarget.params[3] << ", QoP "
444 << tsAtTarget.params[4];
445
446 trk.addTrackState(tsAtTarget);
447 } else {
448 ldmx_log(info) << " Could not extrapolate to target! nhits = "
449 << track.nMeasurements() << " Printing track states: ";
450 for (const auto ts : track.trackStatesReversed()) {
451 if (ts.hasSmoothed()) {
452 ldmx_log(info) << " Track momentum for smoothed track state = "
453 << 1 / ts.smoothed()[Acts::eBoundQOverP];
454 ldmx_log(info) << " Parameters: " << ts.smoothed().transpose();
455 } else {
456 ldmx_log(info) << " Track state not smoothed!";
457 }
458 }
459 ldmx_log(info) << " ...skipping this track candidate...";
460 continue;
461 }
462
463 // get the BoundTrackParameters at the target
464 // ...use to fill in the Acts::TrackProxy object
465 // This isn't really necessary, since we can take
466 // most everything for making the ldmx::track
467 // from tsAtTarget...maybe useful for something?
468 // -->one thing this does is allow Acts to
469 // calculate the momentum 3-vector for you
470 Acts::BoundTrackParameters boundStateAtTarget =
471 tracking::sim::utils::btp(tsAtTarget, target_surface, 11);
472 track.setReferenceSurface(target_surface);
473 track.parameters() = boundStateAtTarget.parameters();
474
475 // ldmx_log(debug) << " typeid(track).name() = " << typeid(track).name();
476 // These are the parameters at the target surface
477 const Acts::BoundVector& track_pars = track.parameters();
478 // const Acts::BoundMatrix& trk_cov = track.covariance();
479
480 ldmx_log(debug) << " Track parameters (at target): Loc0 = "
481 << track_pars[Acts::eBoundLoc0]
482 << ", Loc1 = " << track_pars[Acts::eBoundLoc1]
483 << ", Theta = " << track_pars[Acts::eBoundTheta]
484 << ", Phi = " << track_pars[Acts::eBoundPhi]
485 << ", chi2 = " << track.chi2()
486 << ", nHits = " << track.nMeasurements();
487
488 // the target...it's not really perigee anymore.
489 trk.setPerigeeLocation(0, 0, 0);
490 trk.setPerigeeParameters(tsAtTarget.params);
491 trk.setPerigeeCov(tsAtTarget.cov);
492
493 trk.setChi2(track.chi2());
494 trk.setNhits(track.nMeasurements());
495 // trk.setNdf(track.nDoF());
496 // TODO Switch back to nDoF when Acts is fixed.
497 trk.setNdf(track.nMeasurements() - 5);
498 trk.setNsharedHits(track.nSharedHits());
499
500 ldmx_log(debug) << " Track momentum: px = " << track.momentum()[0]
501 << " py = " << track.momentum()[1]
502 << " pz = " << track.momentum()[2];
503 trk.setMomentum(track.momentum()[0], track.momentum()[1],
504 track.momentum()[2]);
505
506 // At least min_hits_ hits and p > 50 MeV
507 if ((trk.getNhits() <= min_hits_) || (abs(1. / trk.getQoP()) <= 0.05)) {
508 ldmx_log(debug)
509 << " > Track candidate did NOT meet the requirements: Nhits = "
510 << trk.getNhits() << " and p = " << abs(1. / trk.getQoP())
511 << " GeV";
512 // No point of doing the extrapolations if the track candidate anyway
513 // wont be kept
514 continue;
515 }
516
517 // Add measurements to the final track
518 ldmx_log(debug) << " Add measurements to the final track from "
519 << track.nTrackStates() << " TrackStates with "
520 << track.nMeasurements() << " measurements";
521
522 int trk_state_index{0};
523 for (const auto ts : track.trackStatesReversed()) {
524 // Check TrackStates Quality
525 ldmx_log(debug) << " Checking Track State index = "
526 << trk_state_index << " at location "
527 << ts.referenceSurface()
528 .transform(geometry_context())
529 .translation()
530 .transpose();
531
532 if (ts.hasSmoothed()) {
533 ldmx_log(debug) << " Smoothed track parameters: "
534 << ts.smoothed().transpose();
535 // ldmx_log(debug) << " Smoothed covariance mtx:\n" <<
536 // ts.smoothedCovariance();
537 }
538
539 // Check if the track state is a measurement
540 auto typeFlags = ts.typeFlags();
541
542 if (typeFlags.test(Acts::TrackStateFlag::MeasurementFlag) &&
543 ts.hasUncalibratedSourceLink()) {
545 ts.getUncalibratedSourceLink()
546 .template get<ActsExamples::IndexSourceLink>();
547
548 ldmx::Measurement ldmx_meas = measurements.at(sl.index());
549 ldmx_log(debug) << " Adding measurement to ldmx::track with "
550 "source link index = "
551 << sl.index();
552 ldmx_log(trace) << " Measurement:\n" << ldmx_meas;
553 trk.addMeasurementIndex(sl.index());
554 } else {
555 ldmx_log(debug) << " This TrackState is not a measurement";
556 }
557 trk_state_index++;
558 }
559
560 ldmx_log(debug) << " Starting extrapolations";
561 // Extrapolations
562 // To ECAL
563 const double ECAL_SCORING_PLANE = 240.5;
564 Acts::Vector3 pos(ECAL_SCORING_PLANE, 0., 0.);
565 Acts::Translation3 surf_translation(pos);
566 Acts::Transform3 surf_transform(surf_translation * surf_rotation);
567 const std::shared_ptr<Acts::PlaneSurface> ecal_surface =
568 Acts::Surface::makeShared<Acts::PlaneSurface>(surf_transform);
569
570 // To HCAL
571 // Let's extrapolate to 540 mm which is the mid of the side HCAL
572 Acts::Vector3 pos_hcal(540.0, 0., 0.);
573 Acts::Translation3 surf_translation_hcal(pos_hcal);
574 Acts::Transform3 surf_transform_hcal(surf_translation_hcal *
575 surf_rotation);
576 const std::shared_ptr<Acts::PlaneSurface> hcal_surface =
577 Acts::Surface::makeShared<Acts::PlaneSurface>(surf_transform_hcal);
578
579 // Beam Origin unbounded surface
580 const std::shared_ptr<Acts::Surface> beamOrigin_surface =
581 tracking::sim::utils::unboundSurface(-700);
582
583 if (taggerTracking_) {
584 ldmx_log(debug) << " Beam Origin Extrapolation";
585 ldmx::Track::TrackState tsAtBeamOrigin;
586 success = trk_extrap_->TrackStateAtSurface(
587 track, beamOrigin_surface, tsAtBeamOrigin,
588 ldmx::TrackStateType::AtBeamOrigin);
589
590 if (success) {
591 trk.addTrackState(tsAtBeamOrigin);
592 ldmx_log(debug)
593 << " Successfully obtained TrackState at beam origin";
594 }
595 }
596
597 // Recoil Extrapolation to ECAL only
598 if (!taggerTracking_) {
599 ldmx_log(debug) << " Ecal Extrapolation";
601 success = trk_extrap_->TrackStateAtSurface(
602 track, ecal_surface, tsAtEcal, ldmx::TrackStateType::AtECAL);
603
604 if (success) {
605 trk.addTrackState(tsAtEcal);
606 ldmx_log(debug) << " Successfully obtained TrackState at Ecal";
607 ldmx_log(debug) << " Parameters At Ecal: Loc0 = "
608 << tsAtEcal.params[0]
609 << ", Loc1 = " << tsAtEcal.params[1]
610 << ", phi = " << tsAtEcal.params[2]
611 << ", theta = " << tsAtEcal.params[3]
612 << ", QoP = " << tsAtEcal.params[4];
613 } else {
614 ldmx_log(info) << " Could not extrapolate to ECAL!! Please check "
615 "the track states";
616 }
617 }
618
619 // Recoil Extrapolation to HCAL only
620 if (!taggerTracking_) {
621 ldmx_log(debug) << " Hcal Extrapolation";
623 success = trk_extrap_->TrackStateAtSurface(
624 track, hcal_surface, tsAtHcal, ldmx::TrackStateType::AtHCAL);
625
626 if (success) {
627 trk.addTrackState(tsAtHcal);
628 ldmx_log(debug) << " Successfully obtained TrackState at Hcal";
629 ldmx_log(debug) << " Parameters At Hcal: Loc0 = "
630 << tsAtHcal.params[0]
631 << ", Loc1 = " << tsAtHcal.params[1]
632 << ", phi = " << tsAtHcal.params[2]
633 << ", theta = " << tsAtHcal.params[3]
634 << ", QoP = " << tsAtHcal.params[4];
635 } else {
636 ldmx_log(info) << " Could not extrapolate to HCAL!! Please check "
637 "the track states";
638 }
639 }
640
641 // Truth matching
642 if (truthMatchingTool) {
643 auto truthInfo = truthMatchingTool->TruthMatch(trk);
644 trk.setTrackID(truthInfo.trackID);
645 trk.setPdgID(truthInfo.pdgID);
646 trk.setTruthProb(truthInfo.truthProb);
647 }
648
649 // Adding the track candidate to the track collection
650 ldmx_log(debug)
651 << " > Adding the track candidate to the track collection";
652 tracks.push_back(trk);
653 ntracks_++;
654 } // // loop on tracksFromSeed (which usually has 1 element)
655 } // loop seed track parameters (i.e. track candidates)
656
657 ldmx_log(info) << "Number of CKF tracks " << tracks.size();
658
659 auto ckf_run = std::chrono::high_resolution_clock::now();
660 profiling_map_["ckf_run"] +=
661 std::chrono::duration<double, std::milli>(ckf_run - ckf_setup).count();
662
663 // Calculating Shared Hits
664 auto sharedHits = computeSharedHits(tracks, measurements, tg,
665 tracking::sim::utils::sourceLinkHash,
666 tracking::sim::utils::sourceLinkEquality);
667 for (std::size_t iTrack = 0; iTrack < sharedHits.size(); ++iTrack) {
668 tracks[iTrack].setNsharedHits(sharedHits[iTrack].size());
669 for (auto idx : sharedHits[iTrack]) {
670 tracks[iTrack].addSharedIndex(idx);
671 }
672 }
673
674 auto result_loop = std::chrono::high_resolution_clock::now();
675 profiling_map_["result_loop"] +=
676 std::chrono::duration<double, std::milli>(result_loop - ckf_run).count();
677
678 // Add the tracks to the event
679 event.add(out_trk_collection_, tracks);
680
681 auto end = std::chrono::high_resolution_clock::now();
682 // long long microseconds =
683 // std::chrono::duration_cast<std::chrono::microseconds>(end-start).count();
684 auto diff = end - start;
685 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
686}
687
689 if (use1Dmeasurements_)
690 ldmx_log(debug) << "Use1Dmeasurements = " << std::boolalpha
691 << use1Dmeasurements_;
692 if (remove_stereo_)
693 ldmx_log(debug) << "Remove_stereo = " << std::boolalpha << remove_stereo_;
694}
695
697 ldmx_log(info) << "found " << ntracks_ << " tracks / " << nseeds_
698 << " nseeds";
699 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(1)
700 << processing_time_ / nevents_ << " ms";
701 ldmx_log(info) << "Breakdown::";
702 ldmx_log(info) << "setup Avg Time/Event = " << std::fixed
703 << std::setprecision(3) << profiling_map_["setup"] / nevents_
704 << " ms";
705 ldmx_log(info) << "hits Avg Time/Event = " << std::fixed
706 << std::setprecision(2) << profiling_map_["hits"] / nevents_
707 << " ms";
708 ldmx_log(info) << "seeds Avg Time/Event = " << std::fixed
709 << std::setprecision(3) << profiling_map_["seeds"] / nevents_
710 << " ms";
711 ldmx_log(info) << "cf_setup Avg Time/Event = " << std::fixed
712 << std::setprecision(3)
713 << profiling_map_["ckf_setup"] / nevents_ << " ms";
714 ldmx_log(info) << "ckf_run Avg Time/Event = " << std::fixed
715 << std::setprecision(3) << profiling_map_["ckf_run"] / nevents_
716 << " ms";
717 ldmx_log(info) << "result_loop Avg Time/Event = " << std::fixed
718 << std::setprecision(1)
719 << profiling_map_["result_loop"] / nevents_ << " ms";
720}
721
723 dumpobj_ = parameters.getParameter<bool>("dumpobj", 0);
724 pionstates_ = parameters.getParameter<int>("pionstates", 0);
725
726 bfield_ = parameters.getParameter<double>("bfield", -1.5);
727 const_b_field_ = parameters.getParameter<bool>("const_b_field", false);
728 field_map_ = parameters.getParameter<std::string>("field_map");
729 propagator_step_size_ =
730 parameters.getParameter<double>("propagator_step_size", 200.);
731 propagator_maxSteps_ =
732 parameters.getParameter<int>("propagator_maxSteps", 10000);
733 measurement_collection_ = parameters.getParameter<std::string>(
734 "measurement_collection", "TaggerMeasurements");
735 outlier_pval_ = parameters.getParameter<double>("outlier_pval_", 3.84);
736
737 debug_acts_ = parameters.getParameter<bool>("debug_acts", false);
738
739 remove_stereo_ = parameters.getParameter<bool>("remove_stereo", false);
740 use1Dmeasurements_ = parameters.getParameter<bool>("use1Dmeasurements", true);
741 min_hits_ = parameters.getParameter<int>("min_hits", 7);
742
743 // Ckf specific options
744 use_extrapolate_location_ =
745 parameters.getParameter<bool>("use_extrapolate_location", true);
746 extrapolate_location_ = parameters.getParameter<std::vector<double>>(
747 "extrapolate_location", {0., 0., 0.});
748 use_seed_perigee_ = parameters.getParameter<bool>("use_seed_perigee", false);
749
750 // seeds from the event
751 seed_coll_name_ =
752 parameters.getParameter<std::string>("seed_coll_name", "seedTracks");
753
754 // output track collection
755 out_trk_collection_ =
756 parameters.getParameter<std::string>("out_trk_collection", "Tracks");
757
758 // keep track on which system tracking is running
759 taggerTracking_ = parameters.getParameter<bool>("taggerTracking", true);
760
761 // BField Systematics
762 map_offset_ =
763 parameters.getParameter<std::vector<double>>("map_offset_", {0., 0., 0.});
764
765 input_pass_name_ =
766 parameters.getParameter<std::string>("input_pass_name", "");
767}
768
769auto CKFProcessor::makeGeoIdSourceLinkMap(
771 const std::vector<ldmx::Measurement>& measurements)
772 -> std::unordered_multimap<Acts::GeometryIdentifier,
774 std::unordered_multimap<Acts::GeometryIdentifier,
776 geoId_sl_map;
777
778 ldmx_log(debug) << "The makeGeoIdSourceLinkMap has " << measurements.size()
779 << " measurements";
780
781 // Check the hits associated to the surfaces
782 for (unsigned int i_meas = 0; i_meas < measurements.size(); i_meas++) {
783 ldmx::Measurement meas = measurements.at(i_meas);
784 unsigned int layerid = meas.getLayerID();
785
786 const Acts::Surface* hit_surface = tg.getSurface(layerid);
787
788 if (hit_surface) {
789 // Transform the ldmx space point from global to local and store the
790 // information
791
792 ActsExamples::IndexSourceLink idx_sl(hit_surface->geometryId(), i_meas);
793 // mg aug 2024 ... these don't print statements
794 // don't compile using v36 in Acts...figure out later
795 /*
796 ldmx_log(debug)
797 << "Insert measurement on surface located at::"
798 << hit_surface->transform(geometry_context()).translation();
799 ldmx_log(debug) << "and geoId::" << hit_surface->geometryId();
800
801 ldmx_log(debug) << "Surface info::"
802 << std::tie(*hit_surface, geometry_context());
803 */
804 geoId_sl_map.insert(std::make_pair(hit_surface->geometryId(), idx_sl));
805
806 } else
807 ldmx_log(debug) << getName() << "::HIT " << i_meas << " at layer"
808 << (measurements.at(i_meas)).getLayerID()
809 << " is not associated to any surface?!";
810 }
811
812 return geoId_sl_map;
813}
814
815template <typename geometry_t, typename source_link_hash_t,
816 typename source_link_equality_t>
817std::vector<std::vector<std::size_t>> CKFProcessor::computeSharedHits(
818 std::vector<ldmx::Track> tracks, std::vector<ldmx::Measurement> meas_coll,
819 geometry_t& tg, source_link_hash_t&& sourceLinkHash,
820 source_link_equality_t&& sourceLinkEquality) const {
821 auto measurementIndexMap =
822 std::unordered_map<Acts::SourceLink, std::size_t, source_link_hash_t,
823 source_link_equality_t>(0, sourceLinkHash,
824 sourceLinkEquality);
825
826 std::vector<std::vector<std::size_t>> measurementsPerTrack;
827 boost::container::flat_map<std::size_t,
828 boost::container::flat_set<std::size_t>>
829 tracksPerMeasurement;
830 std::vector<std::size_t> sharedMeasurementsPerTrack;
831 auto numberOfTracks = 0;
832
833 // Iterate through all input tracks, collect their properties like measurement
834 // count and chi2 and fill the measurement map in order to relate tracks to
835 // each other if they have shared hits.
836 for (const auto& track : tracks) {
837 // Kick out tracks that do not fulfill our initial requirements
838 // if (track.getNhits() < nMeasurementsMin_) {
839 // continue;
840 // }
841
842 std::vector<std::size_t> measurements;
843 for (auto imeas : track.getMeasurementsIdxs()) {
844 auto meas = meas_coll.at(imeas);
845 const Acts::Surface* hit_surface = tg.getSurface(meas.getLayerID());
846 // Store the index source link
847 ActsExamples::IndexSourceLink idx_sl(hit_surface->geometryId(), imeas);
848 Acts::SourceLink sourceLink = Acts::SourceLink(idx_sl);
849
850 auto emplace = measurementIndexMap.try_emplace(
851 sourceLink, measurementIndexMap.size());
852 measurements.push_back(emplace.first->second);
853 }
854
855 measurementsPerTrack.push_back(std::move(measurements));
856
857 ++numberOfTracks;
858 }
859
860 // Now we relate measurements to tracks
861 for (std::size_t iTrack = 0; iTrack < numberOfTracks; ++iTrack) {
862 for (auto iMeasurement : measurementsPerTrack[iTrack]) {
863 tracksPerMeasurement[iMeasurement].insert(iTrack);
864 }
865 }
866
867 // Finally, we can accumulate the number of shared measurements per track
868 sharedMeasurementsPerTrack = std::vector<std::size_t>(numberOfTracks, 0);
869
870 std::vector<std::vector<std::size_t>> sharedMeasurementIdxsPerTrack;
871 for (std::size_t iTrack = 0; iTrack < numberOfTracks; ++iTrack) {
872 std::vector<std::size_t> sharedMeasurementIdxs;
873 for (auto iMeasurement : measurementsPerTrack[iTrack]) {
874 if (tracksPerMeasurement[iMeasurement].size() > 1) {
875 ++sharedMeasurementsPerTrack[iTrack];
876 sharedMeasurementIdxs.push_back(iMeasurement);
877 }
878 }
879 sharedMeasurementIdxsPerTrack.push_back(sharedMeasurementIdxs);
880 }
881 return sharedMeasurementIdxsPerTrack;
882}
883
884} // namespace reco
885} // namespace tracking
886
887DECLARE_PRODUCER_NS(tracking::reco, CKFProcessor)
#define DECLARE_PRODUCER_NS(NS, 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
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 calibrate_1d(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...