LDMX Software
SeedFinderProcessor.cxx
1#include "Tracking/Reco/SeedFinderProcessor.h"
2
3#include "Acts/Definitions/TrackParametrization.hpp"
4#include "Acts/Seeding/EstimateTrackParamsFromSeed.hpp"
5#include "Eigen/Dense"
6#include "Tracking/Sim/TrackingUtils.h"
7
8/* This processor takes in input a set of 3D space points and builds seedTracks
9 * using the ACTS algorithm which is based on the ATLAS 3-space point conformal
10 * fit.
11 *
12 */
13
14using Eigen::MatrixXd;
15using Eigen::VectorXd;
16
17namespace tracking {
18namespace reco {
19
21 framework::Process& process)
22 : TrackingGeometryUser(name, process) {
23 // TODO REMOVE FROM DEFAULT
24 /*
25 output_file_ = new TFile("seeder.root", "RECREATE");
26 output_tree_ = new TTree("seeder", "seeder");
27
28 output_tree_->Branch("nevents", &nevents_);
29 output_tree_->Branch("xhit", &xhit_);
30 output_tree_->Branch("yhit", &yhit_);
31 output_tree_->Branch("zhit", &zhit_);
32
33 output_tree_->Branch("b0", &b0_);
34 output_tree_->Branch("b1", &b1_);
35 output_tree_->Branch("b2", &b2_);
36 output_tree_->Branch("b3", &b3_);
37 output_tree_->Branch("b4", &b4_);
38 */
39}
40
42 truth_matching_tool_ = std::make_shared<tracking::sim::TruthMatchingTool>();
43}
44
46 // Output seed name
47 out_seed_collection_ = parameters.get<std::string>("out_seed_collection",
48 getName() + "SeedTracks");
49
50 // Input strip hits
52 parameters.get<std::string>("input_hits_collection", "TaggerSimHits");
53
54 // Tagger tracks - only for Recoil Seed finding
56 parameters.get<std::string>("tagger_trks_collection", "TaggerTracks");
57
59 parameters.get<std::vector<double>>("perigee_location", {-700, 0., 0.});
60 pmin_ = parameters.get<double>("pmin", 0.05 * Acts::UnitConstants::GeV);
61 pmax_ = parameters.get<double>("pmax", 8 * Acts::UnitConstants::GeV);
62 d0max_ = parameters.get<double>("d0max", -15. * Acts::UnitConstants::mm);
63 d0min_ = parameters.get<double>("d0min", -45. * Acts::UnitConstants::mm);
64 z0max_ = parameters.get<double>("z0max", 60. * Acts::UnitConstants::mm);
65 phicut_ = parameters.get<double>("phicut", 0.1);
66 thetacut_ = parameters.get<double>("thetacut", 0.2);
67 loc0cut_ = parameters.get<double>("loc0cut", 0.1);
68 loc1cut_ = parameters.get<double>("loc1cut", 0.3);
70 parameters.get<std::vector<std::string>>("strategies", {"0,1,2,3,4"});
71 inflate_factors_ = parameters.get<std::vector<double>>(
72 "inflate_factors", {10., 10., 10., 10., 10., 10.});
73 bfield_ = parameters.get<double>("bfield", 1.5);
74 input_pass_name_ = parameters.get<std::string>("input_pass_name");
75 sim_particles_coll_name_ =
76 parameters.get<std::string>("sim_particles_coll_name");
77 sim_particles_passname_ =
78 parameters.get<std::string>("sim_particles_passname");
79 tagger_trks_event_collection_passname_ =
80 parameters.get<std::string>("tagger_trks_event_collection_passname");
81 sim_particles_event_passname_ =
82 parameters.get<std::string>("sim_particles_event_passname");
83 u_error_ = parameters.get<double>("u_error");
84 v_error_ = parameters.get<double>("v_error");
85}
86
88 // tg is unused, should it be? FIXME
89 // const auto& tg{geometry()};
90 auto start = std::chrono::high_resolution_clock::now();
91 ldmx::Tracks seed_tracks;
92
93 nevents_++;
94
95 // check if SimParticleMap is available for truth matching
96 std::map<int, ldmx::SimParticle> particle_map;
97
98 const std::vector<ldmx::Measurement> measurements =
99 event.getCollection<ldmx::Measurement>(input_hits_collection_,
100 input_pass_name_);
101
102 std::vector<ldmx::Track> tagger_tracks;
104 tagger_trks_event_collection_passname_)) {
105 tagger_tracks = event.getCollection<ldmx::Track>(tagger_trks_collection_,
106 input_pass_name_);
107 }
108
109 // Create an unbound surface at the target
110 std::shared_ptr<Acts::Surface> tgt_surf =
111 tracking::sim::utils::unboundSurface(0.);
112
113 // Create the pseudomeasurements at the target
114
115 ldmx::Measurements target_pseudo_meas;
116
117 for (auto tagtrk : tagger_tracks) {
118 auto ts = tagtrk.getTrackState(ldmx::TrackStateType::AtTarget);
119
120 // The covariance matrix passed to the pseudo measurement is considered as
121 // uncorrelated. This is an approx that considers that loc-u and loc-v from
122 // the track have small correlation.
123
124 if (ts.has_value()) {
125 auto track_state = ts.value();
126
127 Acts::BoundSquareMatrix cov =
128 tracking::sim::utils::unpackCov(track_state.cov_);
129 double locu = track_state.params_[0];
130 double locv = track_state.params_[1];
131 double covuu =
132 cov(Acts::BoundIndices::eBoundLoc0, Acts::BoundIndices::eBoundLoc0);
133 double covvv =
134 cov(Acts::BoundIndices::eBoundLoc1, Acts::BoundIndices::eBoundLoc1);
135
136 ldmx::Measurement pseudo_meas;
137 pseudo_meas.setLocalPosition(locu, locv);
138 Acts::Vector3 dummy{0., 0., 0.};
139 Acts::Vector2 local_pos{locu, locv};
140 Acts::Vector3 global_pos =
141 tgt_surf->localToGlobal(geometryContext(), local_pos, dummy);
142
143 pseudo_meas.setGlobalPosition(global_pos(0), global_pos(1),
144 global_pos(2));
145 pseudo_meas.setTime(0.);
146 pseudo_meas.setLocalCovariance(covuu, covvv);
147
148 target_pseudo_meas.push_back(pseudo_meas);
149 }
150 }
151
152 if (event.exists(sim_particles_coll_name_, sim_particles_event_passname_)) {
153 particle_map = event.getMap<int, ldmx::SimParticle>(
154 sim_particles_coll_name_, sim_particles_passname_);
155 truth_matching_tool_->setup(particle_map, measurements);
156 }
157
158 ldmx_log(debug) << "Preparing the strategies";
159
160 groups_map_.clear();
161 // set the seeding strategy
162 // strategy is a list of layers from which to make the seed
163 // this must include 5 layers; layer_ numbering starts at 0.
164 // std::vector<int> strategy = {9,10,11,12,13};
165 std::vector<int> strategy = {0, 1, 2, 3, 4};
166 bool success = groupStrips(measurements, strategy);
167 if (success) findSeedsFromMap(seed_tracks, target_pseudo_meas);
168
169 // currently, we only use a single strategy but eventually
170 // we will use more. Below is an example of how to add them
171 /*
172 groups_map.clear();
173 strategy = {9,10,11,12,13};
174 success = GroupStrips(measurements,strategy);
175 if (success)
176 FindSeedsFromMap(seed_tracks, target_pseudo_meas);
177 */
178
179 groups_map_.clear();
180 // output_tree_->Fill();
181 ntracks_ += seed_tracks.size();
182 event.add(out_seed_collection_, seed_tracks);
183
184 auto end = std::chrono::high_resolution_clock::now();
185
186 // long long microseconds =
187 // std::chrono::duration_cast<std::chrono::microseconds>(end-start).count();
188
189 auto diff = end - start;
190 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
191
192 // Seed finding using 2D Hits
193 // - The hits should keep track if they are already associated to a track or
194 // not. This can be used for subsequent passes of seed-finding
195
196 // This should go into a digitization producer, which takes care of producing
197 // measurements from:
198 // - raw hits in data
199 // - sim hits in MC
200 // Step 0: Get the sim hits and project them on the surfaces to mimic 2d
201 // hits Step 1: Smear the hits and associate an uncertainty to those
202 // measurements.
203
204 xhit_.clear();
205 yhit_.clear();
206 zhit_.clear();
207
208 b0_.clear();
209 b1_.clear();
210 b2_.clear();
211 b3_.clear();
212 b4_.clear();
213
214} // produce
215
216// Seed finder from Robert's in HPS
217// https://github.com/JeffersonLab/hps-java/blob/47712878302eb0c0374d077a208a6f8f0e2c3dc6/tracking/src/main/java/org/hps/recon/tracking/kalman/SeedTrack.java
218// Adapted to possible 3D hit points.
219
220// yOrigin is the location along the beam about which we fit the seed helix
221// perigee_location is where the track parameters will be extracted
222
223// while this takes in a target measurement (from tagger, this is pmeas_tgt)
224// this code doesn't do anything with it yet.
225
226ldmx::Track SeedFinderProcessor::seedTracker(
227 const ldmx::Measurements& vmeas, double xOrigin,
228 const Acts::Vector3& perigee_location,
229 const ldmx::Measurements& pmeas_tgt) {
230 // Fit a straight line in the non-bending plane and a parabola in the bending
231 // plane
232
233 // Each measurement is treated as a 3D point, where the v direction is in the
234 // center of the strip with sigma equal to the length of the strip / sqrt(12).
235 // In this way it's easier to incorporate the tagger track extrapolation to
236 // the fit
237
238 Acts::ActsMatrix<5, 5> a = Acts::ActsMatrix<5, 5>::Zero();
239 Acts::ActsVector<5> y = Acts::ActsVector<5>::Zero();
240
241 for (auto meas : vmeas) {
242 double xmeas = meas.getGlobalPosition()[0] - xOrigin;
243
244 // Get the surface
245 const Acts::Surface* hit_surface = geometry().getSurface(meas.getLayerID());
246
247 // Get the global to local transformation
248 auto rot = hit_surface->transform(geometryContext()).rotation();
249 auto tr = hit_surface->transform(geometryContext()).translation();
250
251 auto rotl2g = rot.transpose();
252
253 // Only for saving purposes
254 Acts::Vector2 loc{meas.getLocalPosition()[0], 0.};
255
256 xhit_.push_back(xmeas);
257 yhit_.push_back(meas.getGlobalPosition()[1]);
258 zhit_.push_back(meas.getGlobalPosition()[2]);
259
260 Acts::ActsMatrix<2, 5> a_i;
261
262 a_i(0, 0) = rotl2g(0, 1);
263 a_i(0, 1) = rotl2g(0, 1) * xmeas;
264 a_i(0, 2) = rotl2g(0, 1) * xmeas * xmeas;
265 a_i(0, 3) = rotl2g(0, 2);
266 a_i(0, 4) = rotl2g(0, 2) * xmeas;
267
268 a_i(1, 0) = rotl2g(1, 1);
269 a_i(1, 1) = rotl2g(1, 1) * xmeas;
270 a_i(1, 2) = rotl2g(1, 1) * xmeas * xmeas;
271 a_i(1, 3) = rotl2g(1, 2);
272 a_i(1, 4) = rotl2g(1, 2) * xmeas;
273
274 // Fill the yprime vector
275 Acts::Vector2 offset = (rot.transpose() * tr).topRows<2>();
276 Acts::Vector2 xoffset = {rotl2g(0, 0) * xmeas, rotl2g(1, 0) * xmeas};
277
278 loc(0) = meas.getLocalPosition()[0];
279 loc(1) = 0.;
280 // weight matrix
281 Acts::ActsMatrix<2, 2> w_i = Acts::ActsMatrix<2, 2>::Zero();
282
283 w_i(0, 0) = 1. / (u_error_ * u_error_);
284 w_i(1, 1) = 1. / (v_error_ * v_error_);
285
286 Acts::Vector2 yprime_i = loc + offset - xoffset;
287 y += (a_i.transpose()) * w_i * yprime_i;
288
289 Acts::ActsMatrix<2, 5> wa_i = (w_i * a_i);
290 a += a_i.transpose() * wa_i;
291 }
292
293 Acts::ActsVector<5> b;
294 b = a.inverse() * y;
295
296 b0_.push_back(b(0));
297 b1_.push_back(b(1));
298 b2_.push_back(b(2));
299 b3_.push_back(b(3));
300 b4_.push_back(b(4));
301
302 // Acts::ActsVector<5> hlx = Acts::ActsVector<5>::Zero();
303 Acts::ActsVector<3> ref{0., 0., 0.};
304
305 double relative_perigee_x = perigee_location(0) - xOrigin;
306
307 std::shared_ptr<const Acts::PerigeeSurface> seed_perigee =
308 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
309 relative_perigee_x, perigee_location(1), perigee_location(2)));
310
311 // in mm
312 Acts::Vector3 seed_pos{relative_perigee_x,
313 b(0) + b(1) * relative_perigee_x +
314 b(2) * relative_perigee_x * relative_perigee_x,
315 b(3) + b(4) * relative_perigee_x};
316 Acts::Vector3 dir{1, b(1) + 2 * b(2) * relative_perigee_x, b(4)};
317 dir /= dir.norm();
318
319 // Momentum at xmeas
320 // R in meters, p in GeV
321 double p = 0.3 * bfield_ * (1. / (2. * abs(b(2)))) * 0.001;
322 // std::cout<<"Momentum "<< p*dir << std::endl;
323
324 // Convert it to MeV since that's what TrackUtils assumes
325 Acts::Vector3 seed_mom = p * dir / Acts::UnitConstants::MeV;
326 Acts::ActsScalar q =
327 b(2) < 0 ? -1 * Acts::UnitConstants::e : +1 * Acts::UnitConstants::e;
328
329 // Linear intersection with the perigee line. TODO:: Use propagator instead
330 // Project the position on the surface.
331 // This is mainly necessary for the perigee surface, where
332 // the mean might not fulfill the perigee condition.
333
334 // mg Aug 2024 .. interect has changed, but just remove boundary check
335 // and change intersection to intersections
336 // auto intersection =
337 // (*seed_perigee).intersect(geometry_context(), seed_pos, dir, false);
338
339 // Acts::FreeVector seed_free = tracking::sim::utils::toFreeParameters(
340 // intersection.intersection.position, seed_mom, q);
341
342 auto intersection =
343 (*seed_perigee).intersect(geometryContext(), seed_pos, dir);
344
345 Acts::FreeVector seed_free = tracking::sim::utils::toFreeParameters(
346 intersection.intersections()[0].position(), seed_mom, q);
347
348 auto bound_params = Acts::transformFreeToBoundParameters(
349 seed_free, *seed_perigee, geometryContext())
350 .value();
351
352 ldmx_log(trace) << "bound parameters at perigee location" << bound_params;
353
354 Acts::BoundVector stddev;
355 // sigma set to 75% of momentum
356 double sigma_p = 0.75 * p * Acts::UnitConstants::GeV;
357 stddev[Acts::eBoundLoc0] =
358 inflate_factors_[Acts::eBoundLoc0] * 2 * Acts::UnitConstants::mm;
359 stddev[Acts::eBoundLoc1] =
360 inflate_factors_[Acts::eBoundLoc1] * 5 * Acts::UnitConstants::mm;
361 stddev[Acts::eBoundPhi] =
362 inflate_factors_[Acts::eBoundPhi] * 5 * Acts::UnitConstants::degree;
363 stddev[Acts::eBoundTheta] =
364 inflate_factors_[Acts::eBoundTheta] * 5 * Acts::UnitConstants::degree;
365 stddev[Acts::eBoundQOverP] =
366 inflate_factors_[Acts::eBoundQOverP] * (1. / p) * (1. / p) * sigma_p;
367 stddev[Acts::eBoundTime] =
368 inflate_factors_[Acts::eBoundTime] * 1000 * Acts::UnitConstants::ns;
369
370 ldmx_log(debug)
371 << "Making covariance matrix as diagonal matrix with inflated terms";
372 Acts::BoundSquareMatrix bound_cov = stddev.cwiseProduct(stddev).asDiagonal();
373
374 ldmx_log(debug) << "...now putting together the seed track ...";
375
376 ldmx::Track trk = ldmx::Track();
377 trk.setPerigeeLocation(perigee_location(0), perigee_location(1),
378 perigee_location(2));
379 trk.setChi2(0.);
380 trk.setNhits(5);
381 trk.setNdf(0);
382 trk.setNsharedHits(0);
383 std::vector<double> v_seed_params(
384 (bound_params).data(),
385 bound_params.data() + bound_params.rows() * bound_params.cols());
386 std::vector<double> v_seed_cov;
387 tracking::sim::utils::flatCov(bound_cov, v_seed_cov);
388 trk.setPerigeeParameters(v_seed_params);
389 trk.setPerigeeCov(v_seed_cov);
390
391 // Store the global position and momentum at the perigee
392 // TODO:: The eFreePos0 is wrong due to the linear intersection.
393 // YZ are ~ correct
394 trk.setPosition(seed_free[Acts::eFreePos0], seed_free[Acts::eFreePos1],
395 seed_free[Acts::eFreePos2]);
396 trk.setMomentum(seed_free[Acts::eFreeDir0], seed_free[Acts::eFreeDir1],
397 seed_free[Acts::eFreeDir2]);
398
399 ldmx_log(debug)
400 << "...making the ParticleHypothesis ...assume electron for now";
401 auto part_hypo{Acts::SinglyChargedParticleHypothesis::electron()};
402
403 ldmx_log(debug) << "Making BoundTrackParameters seedParameters";
404 Acts::BoundTrackParameters seed_parameters(
405 seed_perigee, std::move(bound_params), bound_cov, part_hypo);
406
407 ldmx_log(debug) << "Returning seed track";
408 return trk;
409}
410
412 // output_file_->cd();
413 // output_tree_->Write();
414 // output_file_->Close();
415 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(1)
416 << processing_time_ / nevents_ << " ms";
417 ldmx_log(info) << "Total Seeds/Events: " << ntracks_ << "/" << nevents_;
418 ldmx_log(info) << "Seeds discarded due to multiple hits on layers "
419 << ndoubles_;
420 ldmx_log(info) << "not enough seed points " << nmissing_;
421 ldmx_log(info) << " nfailpmin=" << nfailpmin_;
422 ldmx_log(info) << " nfailpmax=" << nfailpmax_;
423 ldmx_log(info) << " nfaild0max=" << nfaild0max_;
424 ldmx_log(info) << " nfaild0min=" << nfaild0min_;
425 ldmx_log(info) << " nfailphicut=" << nfailphi_;
426 ldmx_log(info) << " nfailthetacut=" << nfailtheta_;
427 ldmx_log(info) << " nfailz0max=" << nfailz0max_;
428}
429
430// Given a strategy, group the hits according to some options
431// Not a good algorithm. The best would be to organize all the hits in sensors
432// *first* then only select the hits that we are interested into. TODO!
433
434bool SeedFinderProcessor::groupStrips(
435 const std::vector<ldmx::Measurement>& measurements,
436 const std::vector<int> strategy) {
437 // std::cout<<"Using stratedy"<<std::endl;
438 // for (auto& e : strategy) {
439 // std::cout<<e<<" ";
440 //}
441 // std::cout<<std::endl;
442
443 for (auto& meas : measurements) {
444 ldmx_log(trace) << meas;
445
446 if (std::find(strategy.begin(), strategy.end(), meas.getLayer()) !=
447 strategy.end()) {
448 ldmx_log(debug) << "Adding measurement from layer_ = " << meas.getLayer();
449 groups_map_[meas.getLayer()].push_back(&meas);
450 }
451
452 } // loop meas
453
454 if (groups_map_.size() < 5)
455 return false;
456 else
457 return true;
458}
459
460// For each strategy, form all the possible combinatorics and form a seedTrack
461// for each of those This will reshuffle all points. (issue?) Will sort the
462// meas_for_seed vector
463
464void SeedFinderProcessor::findSeedsFromMap(ldmx::Tracks& seeds,
465 const ldmx::Measurements& pmeas) {
466 std::map<int, std::vector<const ldmx::Measurement*>>::iterator groups_iter =
467 groups_map_.begin();
468 // Vector of iterators
469 constexpr size_t k = 5;
470 std::vector<std::vector<const ldmx::Measurement*>::iterator> it;
471 it.reserve(k);
472
473 unsigned int ikey = 0;
474 for (auto& key : groups_map_) {
475 it[ikey] = key.second.begin();
476 ikey++;
477 }
478
479 // K vectors in an array v[0],v[1].... v[K-1]
480
481 // Loop over all combinations
482 while (it[0] != groups_iter->second.end()) {
483 // process the pointed-to elements
484
485 /*
486 for (int j=0; j<K; j++) {
487 const ldmx::Measurement* meas = (*(it[j]));
488 std::cout<<meas->getGlobalPosition()[0]<<","
489 <<meas->getGlobalPosition()[1]<<","
490 <<meas->getGlobalPosition()[2]<<","<<std::endl;
491 }
492 */
493
494 std::vector<ldmx::Measurement> meas_for_seeds;
495 meas_for_seeds.reserve(5);
496
497 ldmx_log(debug) << " Grouping ";
498
499 for (int j = 0; j < k; j++) {
500 const ldmx::Measurement* meas = (*(it[j]));
501 meas_for_seeds.push_back(*meas);
502 }
503
504 std::sort(meas_for_seeds.begin(), meas_for_seeds.end(),
505 [](const ldmx::Measurement& m1, const ldmx::Measurement& m2) {
506 return m1.getGlobalPosition()[0] < m2.getGlobalPosition()[0];
507 });
508
509 if (meas_for_seeds.size() < 5) {
510 nmissing_++;
511 return;
512 }
513
514 ldmx_log(debug) << "making seedTrack";
515
516 Acts::Vector3 perigee{perigee_location_[0], perigee_location_[1],
518
519 ldmx::Track seed_track =
520 seedTracker(meas_for_seeds, meas_for_seeds.at(2).getGlobalPosition()[0],
521 perigee, pmeas);
522
523 bool fail = false;
524
525 // Remove failed fits
526 if (1. / abs(seed_track.getQoP()) < pmin_) {
527 nfailpmin_++;
528 fail = true;
529 } else if (1. / abs(seed_track.getQoP()) > pmax_) {
530 nfailpmax_++;
531 fail = true;
532 }
533
534 // Remove large part of fake tracks and duplicates with the following cuts
535 // for various compatibility checks.
536
537 else if (abs(seed_track.getZ0()) > z0max_) {
538 nfailz0max_++;
539 fail = true;
540 } else if (seed_track.getD0() < d0min_) {
541 nfaild0min_++;
542 fail = true;
543 } else if (seed_track.getD0() > d0max_) {
544 nfaild0max_++;
545 fail = true;
546 } else if (abs(seed_track.getPhi()) > phicut_) {
547 fail = true;
548 nfailphi_++;
549 } else if (abs(seed_track.getTheta() - piover2_) > thetacut_) {
550 fail = true;
551 nfailtheta_++;
552 }
553
554 // If I didn't use the target pseudo measurements in the track finding
555 // I can use them for compatibility with the tagger track
556
557 // TODO this should protect against running this check on tagger seeder.
558 // This is true only if this seeder is not run twice on the tagger after
559 // already having tagger tracks available.
560 if (pmeas.size() > 0) {
561 // I can have multiple target pseudo measurements
562 // A seed is rejected if it is found incompatible with all the target
563 // extrapolations
564
565 // This is set but unused, eventually we will use tagger track position at
566 // target to inform recoil tracking bool tgt_compatible = false;
567 for (auto tgt_pseudomeas : pmeas) {
568 // The d0/z0 are in a frame with the same orientation of the target
569 // surface
570 double delta_loc0 =
571 seed_track.getD0() - tgt_pseudomeas.getLocalPosition()[0];
572 double delta_loc1 =
573 seed_track.getZ0() - tgt_pseudomeas.getLocalPosition()[1];
574
575 if (abs(delta_loc0) < loc0cut_ && abs(delta_loc1) < loc1cut_) {
576 // found at least 1 compatible target location
577 // tgt_compatible = true;
578 break;
579 }
580 }
581 } // pmeas > 0
582
583 if (!fail) {
584 if (truth_matching_tool_->configured()) {
585 auto truth_info = truth_matching_tool_->truthMatch(meas_for_seeds);
586 seed_track.setTrackID(truth_info.track_id_);
587 seed_track.setPdgID(truth_info.pdg_id_);
588 seed_track.setTruthProb(truth_info.truth_prob_);
589 }
590
591 seeds.push_back(seed_track);
592 }
593
594 else {
595 b0_.pop_back();
596 b1_.pop_back();
597 b2_.pop_back();
598 b3_.pop_back();
599 b4_.pop_back();
600 }
601
602 // Go to next combination
603 ldmx_log(debug) << "Go to the next combination";
604
605 ++it[k - 1];
606 for (int i = k - 1;
607 (i > 0) && (it[i] == (std::next(groups_iter, i))->second.end()); --i) {
608 it[i] = std::next(groups_iter, i)->second.begin();
609 ++it[i - 1];
610 }
611 }
612} // find seeds
613
614} // namespace reco
615} // namespace tracking
616
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
std::string getName() const
Get the processor name.
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:37
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
void setLocalPosition(const float &meas_u, const float &meas_v)
Set the local position i.e.
Definition Measurement.h:59
void setGlobalPosition(const float &meas_x, const float &meas_y, const float &meas_z)
Set the global position i.e.
Definition Measurement.h:40
void setLocalCovariance(const float &cov_uu, const float &cov_vv)
Set cov(U,U) and cov(V, V).
Definition Measurement.h:75
void setTime(const float &meas_t)
Set the measurement time in ns.
Definition Measurement.h:91
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:161
SeedFinderProcessor(const std::string &name, framework::Process &process)
Constructor.
std::string out_seed_collection_
The name of the output collection of seeds to be stored.
double pmax_
Maximum cut on the momentum of the seeds.
std::vector< std::string > strategies_
List of stragies for seed finding.
std::string input_hits_collection_
The name of the input hits collection to use in finding seeds..
void onProcessEnd() override
Callback for the EventProcessor to take any necessary action when the processing of events finishes,...
void produce(framework::Event &event) override
Run the processor and create a collection of results which indicate if a charge particle can be found...
void onProcessStart() override
Callback for the EventProcessor to take any necessary action when the processing of events starts,...
double pmin_
Minimum cut on the momentum of the seeds.
std::string tagger_trks_collection_
The name of the tagger Tracks (only for Recoil Seeding)
std::vector< double > perigee_location_
Location of the perigee for the helix track parameters.
double d0max_
Max d0 allowed for the seeds.
void configure(framework::config::Parameters &parameters) override
Configure the processor using the given user specified parameters.
double d0min_
Min d0 allowed for the seeds.
double z0max_
Max z0 allowed for the seeds.
a helper base class providing some methods to shorten access to common conditions used within the tra...
The measurement calibrator can be a function or a class/struct able to retrieve the sim hits containe...