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 outputFile_ = new TFile("seeder.root", "RECREATE");
26 outputTree_ = new TTree("seeder", "seeder");
27
28 outputTree_->Branch("nevents", &nevents_);
29 outputTree_->Branch("xhit", &xhit_);
30 outputTree_->Branch("yhit", &yhit_);
31 outputTree_->Branch("zhit", &zhit_);
32
33 outputTree_->Branch("b0", &b0_);
34 outputTree_->Branch("b1", &b1_);
35 outputTree_->Branch("b2", &b2_);
36 outputTree_->Branch("b3", &b3_);
37 outputTree_->Branch("b4", &b4_);
38 */
39}
40
42
44 truthMatchingTool_ = std::make_shared<tracking::sim::TruthMatchingTool>();
45}
46
48 // Output seed name
49 out_seed_collection_ = parameters.getParameter<std::string>(
50 "out_seed_collection", getName() + "SeedTracks");
51
52 // Input strip hits
53 input_hits_collection_ = parameters.getParameter<std::string>(
54 "input_hits_collection", "TaggerSimHits");
55
56 // Tagger tracks - only for Recoil Seed finding
57 tagger_trks_collection_ = parameters.getParameter<std::string>(
58 "tagger_trks_collection", "TaggerTracks");
59
60 perigee_location_ = parameters.getParameter<std::vector<double>>(
61 "perigee_location", {-700, 0., 0.});
62 pmin_ =
63 parameters.getParameter<double>("pmin", 0.05 * Acts::UnitConstants::GeV);
64 pmax_ = parameters.getParameter<double>("pmax", 8 * Acts::UnitConstants::GeV);
65 d0max_ =
66 parameters.getParameter<double>("d0max", -15. * Acts::UnitConstants::mm);
67 d0min_ =
68 parameters.getParameter<double>("d0min", -45. * Acts::UnitConstants::mm);
69 z0max_ =
70 parameters.getParameter<double>("z0max", 60. * Acts::UnitConstants::mm);
71
72 phicut_ = parameters.getParameter<double>("phicut", 0.1);
73 thetacut_ = parameters.getParameter<double>("thetacut", 0.2);
74
75 loc0cut_ = parameters.getParameter<double>("loc0cut", 0.1);
76 loc1cut_ = parameters.getParameter<double>("loc1cut", 0.3);
77
78 strategies_ = parameters.getParameter<std::vector<std::string>>(
79 "strategies", {"0,1,2,3,4"});
80
81 inflate_factors_ = parameters.getParameter<std::vector<double>>(
82 "inflate_factors", {10., 10., 10., 10., 10., 10.});
83
84 bfield_ = parameters.getParameter<double>("bfield", 1.5);
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> particleMap;
97
98 const std::vector<ldmx::Measurement> measurements =
99 event.getCollection<ldmx::Measurement>(input_hits_collection_);
100
101 std::vector<ldmx::Track> tagger_tracks;
102 if (event.exists(tagger_trks_collection_)) {
103 tagger_tracks = event.getCollection<ldmx::Track>(tagger_trks_collection_);
104 }
105
106 // Create an unbound surface at the target
107 std::shared_ptr<Acts::Surface> tgt_surf =
108 tracking::sim::utils::unboundSurface(0.);
109
110 // Create the pseudomeasurements at the target
111
112 ldmx::Measurements target_pseudo_meas;
113
114 for (auto tagtrk : tagger_tracks) {
115 auto ts = tagtrk.getTrackState(ldmx::TrackStateType::AtTarget);
116
117 // The covariance matrix passed to the pseudo measurement is considered as
118 // uncorrelated. This is an approx that considers that loc-u and loc-v from
119 // the track have small correlation.
120
121 if (ts.has_value()) {
122 auto trackState = ts.value();
123
124 Acts::BoundSquareMatrix cov =
125 tracking::sim::utils::unpackCov(trackState.cov);
126 double locu = trackState.params[0];
127 double locv = trackState.params[1];
128 double covuu =
129 cov(Acts::BoundIndices::eBoundLoc0, Acts::BoundIndices::eBoundLoc0);
130 double covvv =
131 cov(Acts::BoundIndices::eBoundLoc1, Acts::BoundIndices::eBoundLoc1);
132
133 ldmx::Measurement pseudo_meas;
134 pseudo_meas.setLocalPosition(locu, locv);
135 Acts::Vector3 dummy{0., 0., 0.};
136 Acts::Vector2 local_pos{locu, locv};
137 Acts::Vector3 global_pos =
138 tgt_surf->localToGlobal(geometry_context(), local_pos, dummy);
139
140 pseudo_meas.setGlobalPosition(global_pos(0), global_pos(1),
141 global_pos(2));
142 pseudo_meas.setTime(0.);
143 pseudo_meas.setLocalCovariance(covuu, covvv);
144
145 target_pseudo_meas.push_back(pseudo_meas);
146 }
147 }
148
149 if (event.exists("SimParticles")) {
150 particleMap = event.getMap<int, ldmx::SimParticle>("SimParticles");
151 truthMatchingTool_->setup(particleMap, measurements);
152 }
153
154 ldmx_log(debug) << "Preparing the strategies";
155
156 groups_map.clear();
157 // set the seeding strategy
158 // strategy is a list of layers from which to make the seed
159 // this must include 5 layers; layer numbering starts at 0.
160 // std::vector<int> strategy = {9,10,11,12,13};
161 std::vector<int> strategy = {0, 1, 2, 3, 4};
162 bool success = GroupStrips(measurements, strategy);
163 if (success) FindSeedsFromMap(seed_tracks, target_pseudo_meas);
164
165 // currently, we only use a single strategy but eventually
166 // we will use more. Below is an example of how to add them
167 /*
168 groups_map.clear();
169 strategy = {9,10,11,12,13};
170 success = GroupStrips(measurements,strategy);
171 if (success)
172 FindSeedsFromMap(seed_tracks, target_pseudo_meas);
173 */
174
175 groups_map.clear();
176 // outputTree_->Fill();
177 ntracks_ += seed_tracks.size();
178 event.add(out_seed_collection_, seed_tracks);
179
180 auto end = std::chrono::high_resolution_clock::now();
181
182 // long long microseconds =
183 // std::chrono::duration_cast<std::chrono::microseconds>(end-start).count();
184
185 auto diff = end - start;
186 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
187
188 // Seed finding using 2D Hits
189 // - The hits should keep track if they are already associated to a track or
190 // not. This can be used for subsequent passes of seed-finding
191
192 // This should go into a digitization producer, which takes care of producing
193 // measurements from:
194 // - raw hits in data
195 // - sim hits in MC
196 // Step 0: Get the sim hits and project them on the surfaces to mimic 2d hits
197 // Step 1: Smear the hits and associate an uncertainty to those measurements.
198
199 xhit_.clear();
200 yhit_.clear();
201 zhit_.clear();
202
203 b0_.clear();
204 b1_.clear();
205 b2_.clear();
206 b3_.clear();
207 b4_.clear();
208
209} // produce
210
211// Seed finder from Robert's in HPS
212// https://github.com/JeffersonLab/hps-java/blob/47712878302eb0c0374d077a208a6f8f0e2c3dc6/tracking/src/main/java/org/hps/recon/tracking/kalman/SeedTrack.java
213// Adapted to possible 3D hit points.
214
215// yOrigin is the location along the beam about which we fit the seed helix
216// perigee_location is where the track parameters will be extracted
217
218// while this takes in a target measurement (from tagger, this is pmeas_tgt)
219// this code doesn't do anything with it yet.
220
221ldmx::Track SeedFinderProcessor::SeedTracker(
222 const ldmx::Measurements& vmeas, double xOrigin,
223 const Acts::Vector3& perigee_location,
224 const ldmx::Measurements& pmeas_tgt) {
225 // Fit a straight line in the non-bending plane and a parabola in the bending
226 // plane
227
228 // Each measurement is treated as a 3D point, where the v direction is in the
229 // center of the strip with sigma equal to the length of the strip / sqrt(12).
230 // In this way it's easier to incorporate the tagger track extrapolation to
231 // the fit
232
233 Acts::ActsMatrix<5, 5> A = Acts::ActsMatrix<5, 5>::Zero();
234 Acts::ActsVector<5> Y = Acts::ActsVector<5>::Zero();
235
236 for (auto meas : vmeas) {
237 double xmeas = meas.getGlobalPosition()[0] - xOrigin;
238
239 // Get the surface
240 const Acts::Surface* hit_surface = geometry().getSurface(meas.getLayerID());
241
242 // Get the global to local transformation
243 auto rot = hit_surface->transform(geometry_context()).rotation();
244 auto tr = hit_surface->transform(geometry_context()).translation();
245
246 auto rotl2g = rot.transpose();
247
248 // Only for saving purposes
249 Acts::Vector2 loc{meas.getLocalPosition()[0], 0.};
250
251 xhit_.push_back(xmeas);
252 yhit_.push_back(meas.getGlobalPosition()[1]);
253 zhit_.push_back(meas.getGlobalPosition()[2]);
254
255 Acts::ActsMatrix<2, 5> A_i;
256
257 A_i(0, 0) = rotl2g(0, 1);
258 A_i(0, 1) = rotl2g(0, 1) * xmeas;
259 A_i(0, 2) = rotl2g(0, 1) * xmeas * xmeas;
260 A_i(0, 3) = rotl2g(0, 2);
261 A_i(0, 4) = rotl2g(0, 2) * xmeas;
262
263 A_i(1, 0) = rotl2g(1, 1);
264 A_i(1, 1) = rotl2g(1, 1) * xmeas;
265 A_i(1, 2) = rotl2g(1, 1) * xmeas * xmeas;
266 A_i(1, 3) = rotl2g(1, 2);
267 A_i(1, 4) = rotl2g(1, 2) * xmeas;
268
269 // Fill the yprime vector
270 Acts::Vector2 offset = (rot.transpose() * tr).topRows<2>();
271 Acts::Vector2 xoffset = {rotl2g(0, 0) * xmeas, rotl2g(1, 0) * xmeas};
272
273 loc(0) = meas.getLocalPosition()[0];
274 loc(1) = 0.;
275 double uError = sqrt(vmeas[0].getLocalCovariance()[0]);
276
277 // TODO Fix vError for measurements
278 double vError = 40. / sqrt(12);
279
280 Acts::ActsMatrix<2, 2> W_i =
281 Acts::ActsMatrix<2, 2>::Zero(); // weight matrix
282
283 W_i(0, 0) = 1. / (uError * uError);
284 W_i(1, 1) = 1. / (vError * vError);
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 relativePerigeeX = perigee_location(0) - xOrigin;
306
307 std::shared_ptr<const Acts::PerigeeSurface> seed_perigee =
308 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
309 relativePerigeeX, perigee_location(1), perigee_location(2)));
310
311 // in mm
312 Acts::Vector3 seed_pos{relativePerigeeX,
313 B(0) + B(1) * relativePerigeeX +
314 B(2) * relativePerigeeX * relativePerigeeX,
315 B(3) + B(4) * relativePerigeeX};
316 Acts::Vector3 dir{1, B(1) + 2 * B(2) * relativePerigeeX, B(4)};
317 dir /= dir.norm();
318
319 // Momentum at xmeas
320 double p =
321 0.3 * bfield_ * (1. / (2. * abs(B(2)))) * 0.001; // R in meters, p in GeV
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(geometry_context(), 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, geometry_context())
350 .value();
351
352 ldmx_log(debug) << "bound parameters at perigee location" << std::endl
353 << 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 trk.setPerigeeLocation(perigee_location(0), perigee_location(1),
379 perigee_location(2));
380 trk.setChi2(0.);
381 trk.setNhits(5);
382 trk.setNdf(0);
383 trk.setNsharedHits(0);
384 std::vector<double> v_seed_params(
385 (bound_params).data(),
386 bound_params.data() + bound_params.rows() * bound_params.cols());
387 std::vector<double> v_seed_cov;
388 tracking::sim::utils::flatCov(bound_cov, v_seed_cov);
389 trk.setPerigeeParameters(v_seed_params);
390 trk.setPerigeeCov(v_seed_cov);
391
392 // Store the global position and momentum at the perigee
393 // TODO:: The eFreePos0 is wrong due to the linear intersection.
394 // YZ are ~ correct
395 trk.setPosition(seed_free[Acts::eFreePos0], seed_free[Acts::eFreePos1],
396 seed_free[Acts::eFreePos2]);
397 trk.setMomentum(seed_free[Acts::eFreeDir0], seed_free[Acts::eFreeDir1],
398 seed_free[Acts::eFreeDir2]);
399
400 ldmx_log(debug)
401 << "...making the ParticleHypothesis ...assume electron for now";
402 auto partHypo{Acts::SinglyChargedParticleHypothesis::electron()};
403
404 ldmx_log(debug) << "Making BoundTrackParameters seedParameters";
405 Acts::BoundTrackParameters seedParameters(
406 seed_perigee, std::move(bound_params), bound_cov, partHypo);
407
408 ldmx_log(debug) << "Returning seed track";
409 return trk;
410}
411
413 // outputFile_->cd();
414 // outputTree_->Write();
415 // outputFile_->Close();
416 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(1)
417 << processing_time_ / nevents_ << " ms";
418 ldmx_log(info) << "Total Seeds/Events: " << ntracks_ << "/" << nevents_;
419 ldmx_log(info) << "Seeds discarded due to multiple hits on layers "
420 << ndoubles_;
421 ldmx_log(info) << "not enough seed points " << nmissing_;
422 ldmx_log(info) << " nfailpmin=" << nfailpmin_;
423 ldmx_log(info) << " nfailpmax=" << nfailpmax_;
424 ldmx_log(info) << " nfaild0max=" << nfaild0max_;
425 ldmx_log(info) << " nfaild0min=" << nfaild0min_;
426 ldmx_log(info) << " nfailphicut=" << nfailphi_;
427 ldmx_log(info) << " nfailthetacut=" << nfailtheta_;
428 ldmx_log(info) << " nfailz0max=" << nfailz0max_;
429}
430
431// Given a strategy, group the hits according to some options
432// Not a good algorithm. The best would be to organize all the hits in sensors
433// *first* then only select the hits that we are interested into. TODO!
434
435bool SeedFinderProcessor::GroupStrips(
436 const std::vector<ldmx::Measurement>& measurements,
437 const std::vector<int> strategy) {
438 // std::cout<<"Using stratedy"<<std::endl;
439 // for (auto& e : strategy) {
440 // std::cout<<e<<" ";
441 //}
442 // std::cout<<std::endl;
443
444 for (auto& meas : measurements) {
445 ldmx_log(debug) << meas;
446
447 if (std::find(strategy.begin(), strategy.end(), meas.getLayer()) !=
448 strategy.end()) {
449 ldmx_log(debug) << "Adding measurement from layer = " << meas.getLayer();
450 groups_map[meas.getLayer()].push_back(&meas);
451 }
452
453 } // loop meas
454
455 if (groups_map.size() < 5)
456 return false;
457 else
458 return true;
459}
460
461// For each strategy, form all the possible combinatorics and form a seedTrack
462// for each of those This will reshuffle all points. (issue?) Will sort the
463// meas_for_seed vector
464
465void SeedFinderProcessor::FindSeedsFromMap(ldmx::Tracks& seeds,
466 const ldmx::Measurements& pmeas) {
467 std::map<int, std::vector<const ldmx::Measurement*>>::iterator groups_iter =
468 groups_map.begin();
469 // Vector of iterators
470 constexpr size_t K = 5;
471 std::vector<std::vector<const ldmx::Measurement*>::iterator> it;
472 it.reserve(K);
473
474 unsigned int ikey = 0;
475 for (auto& key : groups_map) {
476 it[ikey] = key.second.begin();
477 ikey++;
478 }
479
480 // K vectors in an array v[0],v[1].... v[K-1]
481
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 seedTrack =
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(seedTrack.getQoP()) < pmin_) {
527 nfailpmin_++;
528 fail = true;
529 } else if (1. / abs(seedTrack.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(seedTrack.getZ0()) > z0max_) {
538 nfailz0max_++;
539 fail = true;
540 } else if (seedTrack.getD0() < d0min_) {
541 nfaild0min_++;
542 fail = true;
543 } else if (seedTrack.getD0() > d0max_) {
544 nfaild0max_++;
545 fail = true;
546 } else if (abs(seedTrack.getPhi()) > phicut_) {
547 fail = true;
548 nfailphi_++;
549 } else if (abs(seedTrack.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 seedTrack.getD0() - tgt_pseudomeas.getLocalPosition()[0];
572 double delta_loc1 =
573 seedTrack.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 (truthMatchingTool_->configured()) {
585 auto truthInfo = truthMatchingTool_->TruthMatch(meas_for_seeds);
586 seedTrack.setTrackID(truthInfo.trackID);
587 seedTrack.setPdgID(truthInfo.pdgID);
588 seedTrack.setTruthProb(truthInfo.truthProb);
589 }
590
591 seeds.push_back(seedTrack);
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
617DECLARE_PRODUCER_NS(tracking::reco, SeedFinderProcessor)
#define DECLARE_PRODUCER_NS(NS, 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:41
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:27
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
void setTime(const float &t)
Set the measurement time in ns.
Definition Measurement.h:90
void setGlobalPosition(const float &x, const float &y, const float &z)
Set the global position i.e.
Definition Measurement.h:40
void setLocalPosition(const float &u, const float &v)
Set the local position i.e.
Definition Measurement.h:58
void setLocalCovariance(const float &cov_uu, const float &cov_vv)
Set cov(U,U) and cov(V, V).
Definition Measurement.h:74
Class representing a simulated particle.
Definition SimParticle.h:23
Implementation of a track object.
Definition Track.h:52
void setPerigeeParameters(const std::vector< double > &par)
d_0 z_0 phi_0 theta q/p t
Definition Track.h:144
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...