CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
buildtestMPlex.cc
Go to the documentation of this file.
7 
8 #include "oneapi/tbb/parallel_for.h"
9 
10 #include <memory>
11 
12 namespace mkfit {
13 
14  inline bool sortByHitsChi2(const std::pair<Track, TrackState> &cand1, const std::pair<Track, TrackState> &cand2) {
15  if (cand1.first.nFoundHits() == cand2.first.nFoundHits())
16  return cand1.first.chi2() < cand2.first.chi2();
17 
18  return cand1.first.nFoundHits() > cand2.first.nFoundHits();
19  }
20 
21  inline bool sortByPhi(const Hit &hit1, const Hit &hit2) {
22  return std::atan2(hit1.y(), hit1.x()) < std::atan2(hit2.y(), hit2.x());
23  }
24 
25  inline bool sortByEta(const Hit &hit1, const Hit &hit2) { return hit1.eta() < hit2.eta(); }
26 
27  inline bool sortTracksByEta(const Track &track1, const Track &track2) { return track1.momEta() < track2.momEta(); }
28 
29  inline bool sortTracksByPhi(const Track &track1, const Track &track2) { return track1.momPhi() < track2.momPhi(); }
30 
32  const std::vector<std::vector<Track>> &m_track_candidates;
33 
34  sortTracksByPhiStruct(std::vector<std::vector<Track>> *track_candidates) : m_track_candidates(*track_candidates) {}
35 
36  bool operator()(const std::pair<int, int> &track1, const std::pair<int, int> &track2) {
37  return m_track_candidates[track1.first][track1.second].posPhi() <
38  m_track_candidates[track2.first][track2.second].posPhi();
39  }
40  };
41 
42  // within a layer with a "reasonable" geometry, ordering by Z is the same as eta
43  inline bool sortByZ(const Hit &hit1, const Hit &hit2) { return hit1.z() < hit2.z(); }
44 
45  //==============================================================================
46  // NaN and Silly track parameter check
47  //==============================================================================
48 
49  namespace {
50 
51  int check_nan_n_silly(TrackVec &tracks, const char *prefix) {
52  int count = 0;
53  for (auto &t : tracks) {
54  if (t.hasSillyValues(Const::nan_n_silly_print_bad_cands_bkfit, false, prefix)) {
55  ++count;
56  }
57  }
58  return count;
59  }
60 
61  void check_nan_n_silly_candidates(Event &ev) {
62  // MIMI -- nan_n_silly_per_layer_count is in MkBuilder, could be in MkJob.
63  // if (Const::nan_n_silly_check_cands_every_layer)
64  // {
65  // int sc = (int) ev.nan_n_silly_per_layer_count_;
66  // if (sc > 0)
67  // printf("Nan'n'Silly: Number of silly candidates over all layers = %d\n", sc);
68  // }
70  int sc = check_nan_n_silly(ev.candidateTracks_, "Pre-bkfit silly check");
71  if (sc > 0)
72  printf("Nan'n'Silly: Number of silly pre-bkfit candidates = %d\n", sc);
73  }
74  }
75 
76  void check_nan_n_silly_bkfit(Event &ev) {
78  int sc = check_nan_n_silly(ev.fitTracks_, "Post-bkfit silly check");
79  if (sc > 0)
80  printf("Nan'n'Silly: Number of silly post-bkfit candidates = %d\n", sc);
81  }
82  }
83 
84  } // namespace
85 
86  //==============================================================================
87  // runBuildTestPlexDumbCMSSW
88  //==============================================================================
89 
90  void runBuildingTestPlexDumbCMSSW(Event &ev, const EventOfHits &eoh, MkBuilder &builder) {
91  const IterationConfig &itconf = Config::ItrInfo[0];
92 
93  MkJob job({Config::TrkInfo, itconf, eoh});
94 
95  builder.begin_event(&job, &ev, __func__);
96 
99  }
100 
101  builder.end_event();
102  }
103 
104  //==============================================================================
105  // runBuildTestPlexBestHit
106  //==============================================================================
107 
108  double runBuildingTestPlexBestHit(Event &ev, const EventOfHits &eoh, MkBuilder &builder) {
109  const IterationConfig &itconf = Config::ItrInfo[0];
110 
111  const bool validation_on = (Config::sim_val || Config::quality_val);
112 
113  if (validation_on) {
114  TrackVec seeds1;
115 
116  unsigned int algorithms[] = {4}; //only initialStep
117 
118  for (auto const &s : ev.seedTracks_) {
119  //keep seeds form the first iteration for processing
120  if (std::find(algorithms, algorithms + 1, s.algoint()) != algorithms + 1)
121  seeds1.push_back(s);
122  }
123  ev.seedTracks_.swap(seeds1); //necessary for the validation - PrepareSeeds
124  ev.relabel_bad_seedtracks(); //necessary for the validation - PrepareSeeds
125  }
126 
127  IterationMaskIfc mask_ifc;
128 
129  // To disable hit-masks, pass nullptr in place of &mask_ifc to MkJob ctor
130  // and optionally comment out ev.fill_hitmask_bool_vectors() call.
131 
133 
134  MkJob job({Config::TrkInfo, itconf, eoh, &mask_ifc});
135 
136  builder.begin_event(&job, &ev, __func__);
137 
138  // CCCC builder.PrepareSeeds();
139 
140  // EventOfCandidates event_of_cands;
142 
143 #ifdef USE_VTUNE_PAUSE
144  __SSC_MARK(0x111); // use this to resume Intel SDE at the same point
145  __itt_resume();
146 #endif
147 
148  double time = dtime();
149 
150  builder.findTracksBestHit();
151 
152  time = dtime() - time;
153 
154 #ifdef USE_VTUNE_PAUSE
155  __itt_pause();
156  __SSC_MARK(0x222); // use this to pause Intel SDE at the same point
157 #endif
158 
159  // Hack, get the tracks out.
160  ev.candidateTracks_ = builder.ref_tracks();
161 
162  // For best hit, the candidateTracks_ vector is the direct input to the backward fit so only need to do find_duplicates once
164  //Mark tracks as duplicates; if within CMSSW, remove duplicate tracks before backward fit
165  if (Config::removeDuplicates) {
167  }
168  }
169 
170  // now do backwards fit... do we want to time this section?
171  if (Config::backwardFit) {
172  builder.backwardFitBH();
173  ev.fitTracks_ = builder.ref_tracks();
174  }
175 
176  if (Config::quality_val) {
177  StdSeq::Quality qval;
178  qval.quality_val(&ev);
179  } else if (Config::sim_val || Config::cmssw_val) {
180  StdSeq::root_val(&ev);
181  }
182 
183  builder.end_event();
184 
185  // ev.print_tracks(ev.candidateTracks_, true);
186 
187  return time;
188  }
189 
190  //==============================================================================
191  // runBuildTestPlex Combinatorial: Standard TBB
192  //==============================================================================
193 
194  double runBuildingTestPlexStandard(Event &ev, const EventOfHits &eoh, MkBuilder &builder) {
195  const IterationConfig &itconf = Config::ItrInfo[0];
196 
197  const bool validation_on = (Config::sim_val || Config::quality_val);
198 
199  if (validation_on) {
200  TrackVec seeds1;
201 
202  unsigned int algorithms[] = {4}; //only initialStep
203 
204  for (auto const &s : ev.seedTracks_) {
205  //keep seeds form the first iteration for processing
206  if (std::find(algorithms, algorithms + 1, s.algoint()) != algorithms + 1)
207  seeds1.push_back(s);
208  }
209  ev.seedTracks_.swap(seeds1); //necessary for the validation - PrepareSeeds
210  ev.relabel_bad_seedtracks(); //necessary for the validation - PrepareSeeds
211  }
212 
213  IterationMaskIfc mask_ifc;
214 
215  // To disable hit-masks, pass nullptr in place of &mask_ifc to MkJob ctor
216  // and optionally comment out ev.fill_hitmask_bool_vectors() call.
217 
219 
220  MkJob job({Config::TrkInfo, itconf, eoh, &mask_ifc});
221 
222  builder.begin_event(&job, &ev, __func__);
223 
224  // CCCC builder.PrepareSeeds();
225 
227 
228 #ifdef USE_VTUNE_PAUSE
229  __SSC_MARK(0x111); // use this to resume Intel SDE at the same point
230  __itt_resume();
231 #endif
232 
233  double time = dtime();
234 
235  builder.findTracksStandard();
236 
237  time = dtime() - time;
238 
239 #ifdef USE_VTUNE_PAUSE
240  __itt_pause();
241  __SSC_MARK(0x222); // use this to pause Intel SDE at the same point
242 #endif
243 
244  check_nan_n_silly_candidates(ev);
245 
246  // first store candidate tracks
248 
249  // now do backwards fit... do we want to time this section?
250  if (Config::backwardFit) {
251  // Using the TrackVec version until we home in on THE backward fit etc.
252  // builder.backwardFit();
253  builder.select_best_comb_cands();
254  builder.backwardFitBH();
255  ev.fitTracks_ = builder.ref_tracks();
256 
257  check_nan_n_silly_bkfit(ev);
258  }
259 
261 
262  if (Config::quality_val) {
263  StdSeq::Quality qval;
264  qval.quality_val(&ev);
265  } else if (Config::sim_val || Config::cmssw_val) {
266  StdSeq::root_val(&ev);
267  }
268 
269  builder.end_event();
270 
271  // ev.print_tracks(ev.candidateTracks_, true);
272 
273  return time;
274  }
275 
276  //==============================================================================
277  // runBuildTestPlex Combinatorial: CloneEngine TBB
278  //==============================================================================
279 
280  double runBuildingTestPlexCloneEngine(Event &ev, const EventOfHits &eoh, MkBuilder &builder) {
281  const IterationConfig &itconf = Config::ItrInfo[0];
282 
283  const bool validation_on = (Config::sim_val || Config::quality_val);
284 
285  if (validation_on) {
286  TrackVec seeds1;
287 
288  unsigned int algorithms[] = {4}; //only initialStep
289 
290  for (auto const &s : ev.seedTracks_) {
291  //keep seeds form the first iteration for processing
292  if (std::find(algorithms, algorithms + 1, s.algoint()) != algorithms + 1)
293  seeds1.push_back(s);
294  }
295  ev.seedTracks_.swap(seeds1); //necessary for the validation - PrepareSeeds
296  ev.relabel_bad_seedtracks(); //necessary for the validation - PrepareSeeds
297  }
298 
299  IterationMaskIfc mask_ifc;
300 
301  // To disable hit-masks, pass nullptr in place of &mask_ifc to MkJob ctor
302  // and optionally comment out ev.fill_hitmask_bool_vectors() call.
303 
305 
306  MkJob job({Config::TrkInfo, itconf, eoh, &mask_ifc});
307 
308  builder.begin_event(&job, &ev, __func__);
309 
310  // CCCC builder.PrepareSeeds();
311 
313 
314 #ifdef USE_VTUNE_PAUSE
315  __SSC_MARK(0x111); // use this to resume Intel SDE at the same point
316  __itt_resume();
317 #endif
318 
319  double time = dtime();
320 
321  builder.findTracksCloneEngine();
322 
323  time = dtime() - time;
324 
325 #ifdef USE_VTUNE_PAUSE
326  __itt_pause();
327  __SSC_MARK(0x222); // use this to pause Intel SDE at the same point
328 #endif
329 
330  check_nan_n_silly_candidates(ev);
331 
332  // first store candidate tracks - needed for BH backward fit and root_validation
334 
335  // now do backwards fit... do we want to time this section?
336  if (Config::backwardFit) {
337  // a) TrackVec version:
338  builder.select_best_comb_cands();
339  builder.backwardFitBH();
340  ev.fitTracks_ = builder.ref_tracks();
341 
342  // b) Version that runs on CombCand / TrackCand
343  // builder.backwardFit();
344  // builder.quality_store_tracks(ev.fitTracks_);
345 
346  check_nan_n_silly_bkfit(ev);
347  }
348 
350 
351  // validation section
352  if (Config::quality_val) {
353  StdSeq::Quality qval;
354  qval.quality_val(&ev);
355  } else if (Config::sim_val || Config::cmssw_val) {
356  StdSeq::root_val(&ev);
357  }
358 
359  builder.end_event();
360 
361  // ev.print_tracks(ev.candidateTracks_, true);
362 
363  return time;
364  }
365 
366  //==============================================================================
367  // runBtpCe_MultiIter
368  //
369  // Prototype for running multiple iterations, sequentially, using the same builder.
370  // For cmmsw seeds
371  //
372  // There is, in general, a mess in how tracks are processed, marked, or copied out
373  // in various validation scenarios and export flags.
374  //
375  // In particular, MkBuilder::PrepareSeeds does a lot of things to whole / complete
376  // event,seedTracks_ -- probably this would need to be split into common / and
377  // per-iteration part.
378  // - MkBuilder::prep_*** functions also mostly do not belong there (prep_sim is called from
379  // PrepareSeeds() for cmssw seeds).
380  //
381  // At this point we need to think about what should happen to Event before all the iterations and
382  // after all the iterations ... from the Validation perspective.
383  // And if we care about doing too muich work for seeds that will never get processed.
384  //==============================================================================
385 
386  std::vector<double> runBtpCe_MultiIter(Event &ev, const EventOfHits &eoh, MkBuilder &builder, int n) {
387  std::vector<double> timevec;
388  if (n <= 0)
389  return timevec;
390  timevec.resize(n + 1, 0.0);
391 
392  const bool validation_on = (Config::sim_val || Config::quality_val);
393 
394  TrackVec seeds_used;
395  TrackVec seeds1;
396 
397  unsigned int algorithms[] = {4, 22, 23, 5, 24, 7, 8, 9, 10, 6}; //9 iterations
398 
399  if (validation_on) {
400  for (auto const &s : ev.seedTracks_) {
401  //keep seeds form the first n iterations for processing
402  if (std::find(algorithms, algorithms + n, s.algoint()) != algorithms + n)
403  seeds1.push_back(s);
404  }
405  ev.seedTracks_.swap(seeds1); //necessary for the validation - PrepareSeeds
406  ev.relabel_bad_seedtracks(); //necessary for the validation - PrepareSeeds
407  }
408 
409  IterationMaskIfc mask_ifc;
410  TrackVec seeds;
411  TrackVec tmp_tvec;
412 
413  for (int it = 0; it <= n - 1; ++it) {
414  const IterationConfig &itconf = Config::ItrInfo[it];
415 
416  // To disable hit-masks, pass nullptr in place of &mask_ifc to MkJob ctor
417  // and optionally comment out ev.fill_hitmask_bool_vectors() call.
418 
420 
421  MkJob job({Config::TrkInfo, itconf, eoh, &mask_ifc});
422 
423  builder.begin_event(&job, &ev, __func__);
424 
425  { // We could partition seeds once, store beg, end for each iteration in a map or vector.
426  seeds.clear();
427  int nc = 0;
428  for (auto &s : ev.seedTracks_) {
429  if (s.algoint() == itconf.m_track_algorithm) {
430  if (itconf.m_requires_seed_hit_sorting) {
431  s.sortHitsByLayer();
432  }
433  seeds.push_back(s);
434  ++nc;
435  } else if (nc > 0)
436  break;
437  }
438  }
439 
440  if (itconf.m_requires_dupclean_tight)
441  StdSeq::clean_cms_seedtracks_iter(&seeds, itconf, eoh.refBeamSpot());
442 
443  builder.seed_post_cleaning(seeds);
444 
445  // Add protection in case no seeds are found for iteration
446  if (seeds.size() <= 0)
447  continue;
448 
449  builder.find_tracks_load_seeds(seeds);
450 
451  double time = dtime();
452 
453  builder.findTracksCloneEngine();
454 
455  timevec[it] = dtime() - time;
456  timevec[n] += timevec[it];
457 
458  // Print min and max size of hots vectors of CombCands.
459  // builder.find_min_max_hots_size();
460 
461  if (validation_on)
462  seeds_used.insert(seeds_used.end(), seeds.begin(), seeds.end()); //cleaned seeds need to be stored somehow
463 
466  if (Algo(itconf.m_track_algorithm) == Algo::pixelPairStep) {
467  builder.filter_comb_cands([&](const TrackCand &t) { return StdSeq::qfilter_n_hits_pixseed(t, 3); });
468  } else if (Algo(itconf.m_track_algorithm) == Algo::pixelLessStep) {
469  builder.filter_comb_cands(
470  [&](const TrackCand &t) { return StdSeq::qfilter_pixelLessFwd(t, eoh.refBeamSpot(), Config::TrkInfo); });
471  } else {
472  builder.filter_comb_cands(
473  [&](const TrackCand &t) { return StdSeq::qfilter_n_hits(t, itconf.m_params.minHitsQF); });
474  }
475  }
476 
477  builder.select_best_comb_cands();
478 
479  {
480  builder.export_tracks(tmp_tvec);
481  StdSeq::find_and_remove_duplicates(tmp_tvec, itconf);
482  ev.candidateTracks_.reserve(ev.candidateTracks_.size() + tmp_tvec.size());
483  for (auto &&t : tmp_tvec)
484  ev.candidateTracks_.emplace_back(std::move(t));
485  tmp_tvec.clear();
486  }
487 
488  // now do backwards fit... do we want to time this section?
489  if (Config::backwardFit) {
490  // a) TrackVec version:
491  // builder.backwardFitBH();
492 
493  // b) Version that runs on CombCand / TrackCand
494  const bool do_backward_search = Config::backwardSearch && itconf.m_backward_search;
495 
496  // We copy seed-hits into Candidates ... now we have to remove them so backward fit stops
497  // before reaching seeding region. Ideally, we wouldn't add them in the first place but
498  // if we want to export full tracks above we need to hold on to them (alternatively, we could
499  // have a pointer to seed track in CombCandidate and copy them from there).
500  if (do_backward_search) {
502  }
503 
504  builder.backwardFit();
505 
506  if (do_backward_search) {
507  builder.beginBkwSearch();
509  builder.endBkwSearch();
510  }
511 
515  builder.filter_comb_cands(
516  [&](const TrackCand &t) { return StdSeq::qfilter_n_layers(t, eoh.refBeamSpot(), Config::TrkInfo); });
517  } else if (Algo(itconf.m_track_algorithm) == Algo::pixelLessStep) {
518  builder.filter_comb_cands([&](const TrackCand &t) {
520  });
521  }
522  }
523 
524  builder.filter_comb_cands([&](const TrackCand &t) { return StdSeq::qfilter_nan_n_silly(t); });
525 
526  builder.select_best_comb_cands(true); // true -> clear m_tracks as they were already filled once above
527 
529  builder.export_tracks(ev.fitTracks_);
530  }
531 
532  builder.end_event();
533  }
534 
535  // MIMI - Fake back event pointer for final processing (that should be done elsewhere)
536  MkJob job({Config::TrkInfo, Config::ItrInfo[0], eoh});
537  builder.begin_event(&job, &ev, __func__);
538 
539  if (validation_on) {
541  //swap for the cleaned seeds
542  ev.seedTracks_.swap(seeds_used);
543  }
544 
545  check_nan_n_silly_candidates(ev);
546 
548  check_nan_n_silly_bkfit(ev);
549 
550  // validation section
551  if (Config::quality_val) {
552  StdSeq::Quality qval;
553  qval.quality_val(&ev);
554  } else if (Config::sim_val || Config::cmssw_val) {
555  StdSeq::root_val(&ev);
556  }
557 
558  // ev.print_tracks(ev.candidateTracks_, true);
559 
560  // MIMI Unfake.
561  builder.end_event();
562 
563  // In CMSSW runOneIter we now release memory for comb-cands:
564  builder.release_memory();
565 
566  return timevec;
567  }
568 
569 } // end namespace mkfit
void find_and_remove_duplicates(TrackVec &tracks, const IterationConfig &itconf)
Definition: MkStdSeqs.cc:591
bool qfilter_n_hits_pixseed(const TRACK &t, int nMinHits)
Definition: MkStdSeqs.h:50
void find_tracks_load_seeds_BH(const TrackVec &in_seeds)
Definition: MkBuilder.cc:431
bool qfilter_pixelLessFwd(const TRACK &t, const BeamSpot &bspot, const TrackerInfo &tk_info)
quality filter tuned for pixelLess iteration during forward search
Definition: MkStdSeqs.h:78
IterationParams m_params
void export_tracks(TrackVec &out_vec)
Definition: MkBuilder.cc:369
void findTracksBestHit(SteeringParams::IterationType_e iteration_dir=SteeringParams::IT_FwdSearch)
Definition: MkBuilder.cc:447
void find_duplicates(TrackVec &tracks)
Definition: MkStdSeqs.cc:335
const BeamSpot & refBeamSpot() const
bool qfilter_n_hits(const TRACK &t, int nMinHits)
Definition: MkStdSeqs.h:42
void relabel_bad_seedtracks()
Definition: Event.cc:778
constexpr bool nan_n_silly_check_cands_post_bkfit
Definition: Config.h:63
bool qfilter_n_layers(const TRACK &t, const BeamSpot &bspot, const TrackerInfo &trk_inf)
Definition: MkStdSeqs.h:56
constexpr bool nan_n_silly_print_bad_cands_bkfit
Definition: Config.h:64
int filter_comb_cands(std::function< filter_track_cand_foo > filter)
Definition: MkBuilder.cc:277
double runBuildingTestPlexStandard(Event &ev, const EventOfHits &eoh, MkBuilder &builder)
constexpr bool nan_n_silly_check_cands_pre_bkfit
Definition: Config.h:62
void compactifyHitStorageForBestCand(bool remove_seed_hits, int backward_fit_min_hits)
Definition: MkBuilder.h:106
bool sortTracksByPhi(const Track &track1, const Track &track2)
void root_val(Event *event)
TrackVec seedTracks_
Definition: Event.h:61
float y() const
Definition: Hit.h:162
void findTracksCloneEngine(SteeringParams::IterationType_e iteration_dir=SteeringParams::IT_FwdSearch)
Definition: MkBuilder.cc:900
auto const & tracks
cannot be loose
void begin_event(MkJob *job, Event *ev, const char *build_type)
Definition: MkBuilder.cc:177
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
bool sortByHitsChi2(const std::pair< Track, TrackState > &cand1, const std::pair< Track, TrackState > &cand2)
void backwardFitBH()
Definition: MkBuilder.cc:1132
const TrackVec & ref_tracks() const
Definition: MkBuilder.h:114
void prep_simtracks(Event *event)
float z() const
Definition: Hit.h:163
void findTracksStandard(SteeringParams::IterationType_e iteration_dir=SteeringParams::IT_FwdSearch)
Definition: MkBuilder.cc:721
const std::vector< std::vector< Track > > & m_track_candidates
sortTracksByPhiStruct(std::vector< std::vector< Track >> *track_candidates)
void release_memory()
Definition: MkBuilder.cc:205
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
tuple Algo
Definition: ntupleEnum.py:25
TrackVec & ref_tracks_nc()
Definition: MkBuilder.h:115
def move
Definition: eostools.py:511
bool qfilter_pixelLessBkwd(const TRACK &t, const BeamSpot &bspot, const TrackerInfo &tk_info)
quality filter tuned for pixelLess iteration during backward search
Definition: MkStdSeqs.h:104
TrackerInfo TrkInfo
void seed_post_cleaning(TrackVec &tv)
Definition: MkBuilder.cc:380
void fill_hitmask_bool_vectors(int track_algo, std::vector< std::vector< bool >> &layer_masks)
Definition: Event.cc:805
void select_best_comb_cands(bool clear_m_tracks=false, bool remove_missing_hits=false)
Definition: MkBuilder.cc:348
void runBuildingTestPlexDumbCMSSW(Event &ev, const EventOfHits &eoh, MkBuilder &builder)
IterationsInfo ItrInfo
float momEta() const
Definition: Track.h:173
bool sortByZ(const Hit &hit1, const Hit &hit2)
std::vector< Track > TrackVec
double runBuildingTestPlexCloneEngine(Event &ev, const EventOfHits &eoh, MkBuilder &builder)
float x() const
Definition: Hit.h:161
void find_tracks_load_seeds(const TrackVec &in_seeds)
Definition: MkBuilder.cc:596
void beginBkwSearch()
Definition: MkBuilder.h:110
float eta() const
Definition: Hit.h:168
double dtime()
bool qfilter_nan_n_silly(const TRACK &t)
Definition: MkStdSeqs.h:128
std::vector< double > runBtpCe_MultiIter(Event &ev, const EventOfHits &eoh, MkBuilder &builder, int n)
void export_best_comb_cands(TrackVec &out_vec, bool remove_missing_hits=false)
Definition: MkBuilder.cc:354
bool operator()(const std::pair< int, int > &track1, const std::pair< int, int > &track2)
TrackVec fitTracks_
Definition: Event.h:61
void handle_duplicates(Event *event)
TrackVec candidateTracks_
Definition: Event.h:61
void endBkwSearch()
Definition: MkBuilder.h:111
TrackAlgorithm
track algorithm; copy from TrackBase.h to keep in standalone builds
Definition: Track.h:275
void root_val_dumb_cmssw(Event *event)
void quality_val(Event *event)
bool sortTracksByEta(const Track &track1, const Track &track2)
std::vector< std::vector< bool > > m_mask_vector
bool sortByPhi(const Hit &hit1, const Hit &hit2)
float momPhi() const
Definition: Track.h:172
Definition: fakeMenu.h:6
double runBuildingTestPlexBestHit(Event &ev, const EventOfHits &eoh, MkBuilder &builder)
bool sortByEta(const Hit &hit1, const Hit &hit2)
int clean_cms_seedtracks_iter(TrackVec *seed_ptr, const IterationConfig &itrcfg, const BeamSpot &bspot)
Definition: MkStdSeqs.cc:84