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