CMS 3D CMS Logo

Shell.cc
Go to the documentation of this file.
2 
4 
5 // #include "RecoTracker/MkFitCore/src/Matriplex/MatriplexCommon.h"
6 
8 
14 
17 
19 
20 #ifndef NO_ROOT
21 #include "TROOT.h"
22 #include "TRint.h"
23 #endif
24 
25 #include "oneapi/tbb/task_arena.h"
26 
27 #include <vector>
28 
29 // clang-format off
30 
31 namespace {
32  constexpr int algos[] = {4, 22, 23, 5, 24, 7, 8, 9, 10, 6}; // 10 iterations
33  constexpr int n_algos = sizeof(algos) / sizeof(int);
34 
35  const char* b2a(bool b) { return b ? "true" : "false"; }
36 }
37 
38 namespace mkfit {
39 
40  Shell::Shell(std::vector<DeadVec> &dv, const std::string &in_file, int start_ev)
41  : m_deadvectors(dv)
42  {
45 
47 
48  m_data_file = new DataFile;
50 
51  m_event = new Event(0, Config::TrkInfo.n_layers());
52  GoToEvent(start_ev);
53  }
54 
55  void Shell::Run() {
56 #ifndef NO_ROOT
57  std::vector<const char *> argv = { "mkFit", "-l" };
58  int argc = argv.size();
59  TRint rint("mkFit-shell", &argc, const_cast<char**>(argv.data()));
60 
61  char buf[256];
62  sprintf(buf, "mkfit::Shell &s = * (mkfit::Shell*) %p;", this);
63  gROOT->ProcessLine(buf);
64  printf("Shell &s variable is set\n");
65 
66  rint.Run(true);
67  printf("Shell::Run finished\n");
68 #else
69  printf("Shell::Run() no root, we rot -- erroring out. Recompile with WITH_ROOT=1\n");
70  exit(1);
71 #endif
72  }
73 
74  void Shell::Status() {
75  printf("On event %d, selected iteration index %d, algo %d - %s\n"
76  " debug = %s, use_dead_modules = %s\n"
77  " clean_seeds = %s, backward_fit = %s, remove_duplicates = %s\n",
81  }
82 
83  //===========================================================================
84  // Event navigation / processing
85  //===========================================================================
86 
87  void Shell::GoToEvent(int eid) {
88  if (eid < 1) {
89  fprintf(stderr, "Requested event %d is less than 1 -- 1 is the first event, %d is total number of events in file\n",
91  throw std::runtime_error("event out of range");
92  }
93  if (eid > m_evs_in_file) {
94  fprintf(stderr, "Requested event %d is grater than total number of events in file %d\n",
96  throw std::runtime_error("event out of range");
97  }
98 
99  int pos = m_event->evtID();
100  if (eid > pos) {
101  m_data_file->skipNEvents(eid - pos - 1);
102  } else {
103  m_data_file->rewind();
105  }
106  m_event->reset(eid);
111  }
112 
113  printf("At event %d\n", eid);
114  }
115 
116  void Shell::NextEvent(int skip) {
117  GoToEvent(m_event->evtID() + skip);
118  }
119 
120  void Shell::ProcessEvent(SeedSelect_e seed_select, int selected_seed, int count) {
121  // count is only used for SS_IndexPreCleaning and SS_IndexPostCleaning.
122  // There are no checks for upper bounds, ie, if requested seeds beyond the first one exist.
123 
124  const IterationConfig &itconf = Config::ItrInfo[m_it_index];
125  IterationMaskIfc mask_ifc;
127 
128  m_seeds.clear();
129  m_tracks.clear();
130 
131  {
132  int n_algo = 0; // seeds are grouped by algo
133  for (auto &s : m_event->seedTracks_) {
134  if (s.algoint() == itconf.m_track_algorithm) {
135  if (seed_select == SS_UseAll || seed_select == SS_IndexPostCleaning) {
136  m_seeds.push_back(s);
137  } else if (seed_select == SS_Label && s.label() == selected_seed) {
138  m_seeds.push_back(s);
139  break;
140  } else if (seed_select == SS_IndexPreCleaning && n_algo >= selected_seed) {
141  m_seeds.push_back(s);
142  if (--count <= 0)
143  break;
144  }
145  ++n_algo;
146  } else if (n_algo > 0)
147  break;
148  }
149  }
150 
151  printf("Shell::ProcessEvent running over %d seeds\n", (int) m_seeds.size());
152 
153  // Equivalent to run_OneIteration(...) without MkBuilder::release_memory().
154  // If seed_select == SS_IndexPostCleaning the given seed is picked after cleaning.
155  {
156  const TrackerInfo &trackerInfo = Config::TrkInfo;
157  const EventOfHits &eoh = *m_eoh;
158  const IterationMaskIfcBase &it_mask_ifc = mask_ifc;
161  TrackVec &out_tracks = m_tracks;
162  bool do_seed_clean = m_clean_seeds;
163  bool do_backward_fit = m_backward_fit;
164  bool do_remove_duplicates = m_remove_duplicates;
165 
166  MkJob job({trackerInfo, itconf, eoh, eoh.refBeamSpot(), &it_mask_ifc});
167 
168  builder.begin_event(&job, m_event, __func__);
169 
170  // Seed cleaning not done on all iterations.
171  do_seed_clean = m_clean_seeds && itconf.m_seed_cleaner;
172 
173  if (do_seed_clean)
174  itconf.m_seed_cleaner(seeds, itconf, eoh.refBeamSpot());
175 
176  // Check nans in seeds -- this should not be needed when Slava fixes
177  // the track parameter coordinate transformation.
179 
180  if (seed_select == SS_IndexPostCleaning) {
181  if (selected_seed >= 0 && selected_seed < (int)seeds.size()) {
182  for (int i = 0; i < count; ++i)
183  seeds[i] = seeds[selected_seed + i];
184  seeds.resize(count);
185  } else {
186  seeds.clear();
187  }
188  }
189 
190  if (seeds.empty()) {
191  if (seed_select != SS_UseAll)
192  printf("Shell::ProcessEvent requested seed not found among seeds of the selected iteration.\n");
193  else
194  printf("Shell::ProcessEvent no seeds found.\n");
195  return;
196  }
197 
198  if (itconf.m_requires_seed_hit_sorting) {
199  for (auto &s : seeds)
200  s.sortHitsByLayer(); // sort seed hits for the matched hits (I hope it works here)
201  }
202 
204 
205  builder.find_tracks_load_seeds(seeds, do_seed_clean);
206 
208 
209  printf("Shell::ProcessEvent post fwd search: %d comb-cands\n", builder.ref_eocc().size());
210 
211  // Pre backward-fit filtering.
212  filter_candidates_func pre_filter;
213  if (do_backward_fit && itconf.m_pre_bkfit_filter)
214  pre_filter = [&](const TrackCand &tc, const MkJob &jb) -> bool {
215  return itconf.m_pre_bkfit_filter(tc, jb) && StdSeq::qfilter_nan_n_silly<TrackCand>(tc, jb);
216  };
217  else if (itconf.m_pre_bkfit_filter)
218  pre_filter = itconf.m_pre_bkfit_filter;
219  else if (do_backward_fit)
220  pre_filter = StdSeq::qfilter_nan_n_silly<TrackCand>;
221  // pre_filter can be null if we are not doing backward fit as nan_n_silly will be run below.
222  if (pre_filter)
223  builder.filter_comb_cands(pre_filter, true);
224 
225  printf("Shell::ProcessEvent post pre-bkf-filter (%s) and nan-filter (%s) filter: %d comb-cands\n",
226  b2a(bool(itconf.m_pre_bkfit_filter)), b2a(do_backward_fit), builder.ref_eocc().size());
227 
228  job.switch_to_backward();
229 
230  if (do_backward_fit) {
231  if (itconf.m_backward_search) {
233  }
234 
236 
237  if (itconf.m_backward_search) {
240  }
241 
242  printf("Shell::ProcessEvent post backward fit / search: %d comb-cands\n", builder.ref_eocc().size());
243  }
244 
245  // Post backward-fit filtering.
246  filter_candidates_func post_filter;
247  if (do_backward_fit && itconf.m_post_bkfit_filter)
248  post_filter = [&](const TrackCand &tc, const MkJob &jb) -> bool {
249  return itconf.m_post_bkfit_filter(tc, jb) && StdSeq::qfilter_nan_n_silly<TrackCand>(tc, jb);
250  };
251  else
252  post_filter = StdSeq::qfilter_nan_n_silly<TrackCand>;
253  // post_filter is always at least doing nan_n_silly filter.
254  builder.filter_comb_cands(post_filter, true);
255 
256  printf("Shell::ProcessEvent post post-bkf-filter (%s) and nan-filter (true): %d comb-cands\n",
257  b2a(do_backward_fit && itconf.m_post_bkfit_filter), builder.ref_eocc().size());
258 
259  if (do_backward_fit && itconf.m_backward_search)
261 
262  builder.export_best_comb_cands(out_tracks, true);
263 
264  if (do_remove_duplicates && itconf.m_duplicate_cleaner) {
265  itconf.m_duplicate_cleaner(out_tracks, itconf);
266  }
267 
268  printf("Shell::ProcessEvent post remove-duplicates: %d comb-cands\n", (int) out_tracks.size());
269 
271 
272  builder.end_event();
273  }
274 
275  printf("Shell::ProcessEvent found %d tracks, number of seeds at end %d\n",
276  (int) m_tracks.size(), (int) m_seeds.size());
277  }
278 
279  //===========================================================================
280  // Iteration selection
281  //===========================================================================
282 
283  void Shell::SelectIterationIndex(int itidx) {
284  if (itidx < 0 || itidx >= n_algos) {
285  fprintf(stderr, "Requested iteration index out of range [%d, %d)", 0, n_algos);
286  throw std::runtime_error("iteration index out of range");
287  }
288  m_it_index = itidx;
289  }
290 
292  for (int i = 0; i < n_algos; ++i) {
293  if (algos[i] == algo) {
294  m_it_index = i;
295  return;
296  }
297  }
298  fprintf(stderr, "Requested algo %d not found", algo);
299  throw std::runtime_error("algo not found");
300  }
301 
303  printf("Shell::PrintIterations selected index = %d, %d iterations hardcoded as\n",
304  m_it_index, n_algos);
305  for (int i = 0; i < n_algos; ++i)
306  printf("%d %2d %s\n", i, algos[i], TrackBase::algoint_to_cstr(algos[i]));
307  }
308 
309  //===========================================================================
310  // Flags / status setters
311  //===========================================================================
312 
313  void Shell::SetDebug(bool b) { g_debug = b; }
318 
319  //===========================================================================
320  // Analysis helpers
321  //===========================================================================
322 
323  /*
324  sim tracks are written to .bin files with a label equal to its own index.
325  reco tracks labels are seed indices.
326  seed labels are sim track indices
327  --
328  mkfit labels are seed indices in given iteration after cleaning (at seed load-time)
329  */
330 
331  int Shell::LabelFromHits(Track &t, bool replace, float good_frac) {
332  std::map<int, int> lab_cnt;
333  for (int hi = 0; hi < t.nTotalHits(); ++hi) {
334  auto hot = t.getHitOnTrack(hi);
335  if (hot.index < 0)
336  continue;
337  const Hit &h = m_event->layerHits_[hot.layer][hot.index];
338  int hl = m_event->simHitsInfo_[h.mcHitID()].mcTrackID_;
339  if (hl >= 0)
340  ++lab_cnt[hl];
341  }
342  int max_c = -1, max_l = -1;
343  for (auto& x : lab_cnt) {
344  if (x.second > max_c) {
345  max_l = x.first;
346  max_c = x.second;
347  }
348  }
349  bool success = max_c >= good_frac * t.nFoundHits();
350  int relabel = success ? max_l : -1;
351  // printf("found_hits=%d, best_lab %d (%d hits), existing label=%d (replace flag=%s)\n",
352  // t.nFoundHits(), max_l, max_c, t.label(), b2a(replace));
353  if (replace)
354  t.setLabel(relabel);
355  return relabel;
356  }
357 
359  Event &ev = *m_event;
360  const int track_algo = Config::ItrInfo[m_it_index].m_track_algorithm;
361 
362  m_ckf_map.clear();
363  m_sim_map.clear();
364  m_seed_map.clear();
365  m_mkf_map.clear();
366 
367  // Pick ckf tracks with right algo and a good label.
368  int rec_algo_match = 0;
369  for (auto &t : ev.cmsswTracks_) {
370  if (t.algoint() != track_algo)
371  continue;
372  ++rec_algo_match;
373  int label = LabelFromHits(t, false, 0.5);
374  if (label >= 0) {
375  m_ckf_map.insert(std::make_pair(label, &t));
376  }
377  }
378 
379  // Pick sim tracks with labels found by ckf.
380  for (auto &t : ev.simTracks_) {
381  if (t.label() >= 0 && m_ckf_map.find(t.label()) != m_ckf_map.end()) {
382  m_sim_map.insert(std::make_pair(t.label(), &t));
383  }
384  }
385 
386  // Pick seeds with right algo and a label found by ckf.
387  for (auto &t : ev.seedTracks_) {
388  if (t.algoint() == track_algo && t.label() >= 0 && m_ckf_map.find(t.label()) != m_ckf_map.end()) {
389  m_seed_map.insert(std::make_pair(t.label(), &t));
390  }
391  }
392  // Some seeds seem to be labeled -1, try fixing when not otherwise found.
393  for (auto &t : ev.seedTracks_) {
394  if (t.algoint() == track_algo && t.label() == -1) {
395  int lab = LabelFromHits(t, true, 0.5);
396  if (lab >= 0 && m_seed_map.find(lab) == m_seed_map.end()) {
397  if (m_ckf_map.find(lab) != m_ckf_map.end()) {
398  m_seed_map.insert(std::make_pair(t.label(), &t));
399  printf("Saved seed with label -1 -> %d\n", lab);
400  }
401  }
402  }
403  }
404 
405  // Pick mkfit tracks, label by
406  for (auto &t : m_tracks) {
407  int label = LabelFromHits(t, false, 0.5);
408  if (label >= 0) {
409  m_mkf_map.insert(std::make_pair(label, &t));
410  }
411  }
412 
413  printf("Shell::FillByLabelMaps reporting tracks with label >= 0, algo=%d (%s): "
414  "ckf: %d of %d (same algo=%d)), sim: %d of %d, seed: %d of %d, mkfit: %d w/label of %d\n",
415  track_algo, TrackBase::algoint_to_cstr(track_algo),
416  (int) m_ckf_map.size(), (int) ev.cmsswTracks_.size(), rec_algo_match,
417  (int) m_sim_map.size(), (int) ev.simTracks_.size(),
418  (int) m_seed_map.size(), (int) m_seeds.size(),
419  (int) m_mkf_map.size(), (int) m_tracks.size()
420  );
421  }
422 
423  bool Shell::CheckMkFitLayerPlanVsReferenceHits(const Track &mkft, const Track &reft, const std::string &name) {
424  // Check if all hit-layers of a reference track reft are in the mkfit layer plan.
425  // Returns true if all layers are in the plan.
426  // String name is printed in front of label, expected to be SIMK or CKF.
427  const IterationConfig &itconf = Config::ItrInfo[m_it_index];
428  auto lp = itconf.m_steering_params[ mkft.getEtaRegion() ].m_layer_plan;
429  bool ret = true;
430  for (int hi = 0; hi < reft.nTotalHits(); ++hi) {
431  auto hot = reft.getHitOnTrack(hi);
432  if (std::find_if(lp.begin(), lp.end(), [=](auto &x){ return x.m_layer == hot.layer; }) == lp.end())
433  {
434  printf("CheckMkfLayerPlanVsCkfHits: layer %d not in layer plan for region %d, %s label=%d\n",
435  hot.layer, mkft.getEtaRegion(), name.c_str(), reft.label());
436  ret = false;
437  }
438  }
439  return ret;
440  }
441 
442  //===========================================================================
443  // Analysis drivers / main functions / Comparators
444  //===========================================================================
445 
446  void Shell::Compare() {
447  Event &ev = *m_event;
448  const IterationConfig &itconf = Config::ItrInfo[m_it_index];
449 
451 
452  printf("------------------------------------------------------\n");
453 
454  const bool print_all_def = false;
455  int mkf_cnt=0, less_hits=0, more_hits=0;
456 
457  // TOBTEC: look for rec-seeds with hits in tob1 and 2 only.
458  int tot_cnt = 0, no_mkf_cnt = 0;
459 
460  for (auto& [l, simt_ptr]: m_sim_map)
461  {
462  auto &simt = * simt_ptr;
463  auto &ckft = * m_ckf_map[l];
464  auto mi = m_mkf_map.find(l);
465 
466  bool print_all = print_all_def;
467 
468  // TOBTEC: look for rec-seeds with hits in tob1 and 2 only.
469  bool select = true;
470  {
471  auto &ckf_seed = ev.seedTracks_[ckft.label()];
472  for (int hi = 0; hi < ckf_seed.nTotalHits(); ++hi) {
473  const HitOnTrack hot = ckf_seed.getHitOnTrack(hi);
474  if (hot.index >= 0 && (hot.layer < 10 || hot.layer > 13)) {
475  select = false;
476  break;
477  }
478  }
479  }
480  if ( ! select) continue;
481 
482  ++tot_cnt;
483  //print_all = true;
484 
485  if (mi != m_mkf_map.end())
486  {
487  auto &mkft = * mi->second;
488  mkf_cnt++;
489  if (mkft.nFoundHits() < ckft.nFoundHits()) ++less_hits;
490  if (mkft.nFoundHits() > ckft.nFoundHits()) ++more_hits;
491 
492  CheckMkFitLayerPlanVsReferenceHits(mkft, ckft, "CKF");
493  // CheckMkFitLayerPlanVsReferenceHits(mkft, simt, "SIM");
494 
495  (void) print_all;
496  if (/* itconf.m_track_algorithm == 10 ||*/ print_all) {
497  // ckf label is wrong when validation is on (even quality val) for mixedTriplet, pixelless and tobtec
498  // as seed tracks get removed for non-mkfit iterations and indices from rec-tracks are no longer valid.
499  auto &ckf_seed = ev.seedTracks_[ckft.label()];
500  auto &mkf_seed = m_seeds[mkft.label()];
501  print("ckf ", 0, ckft, ev);
502  print("mkfit", 0, mkft, ev);
503  print("sim ", 0, simt, ev);
504  print("ckf seed", 0, ckf_seed, ev);
505  print("mkf seed", 0, mkf_seed, ev);
506  printf("------------------------------------------------------\n");
507 
508  TrackVec ssss;
509  ssss.push_back(mkf_seed);
510 
511  IterationSeedPartition pppp(1);
512  IterationConfig::get_seed_partitioner("phase1:1:debug")(Config::TrkInfo, ssss, *m_eoh, pppp);
513 
514  printf("------------------------------------------------------\n");
515  printf("\n");
516  }
517  }
518  else
519  {
520  printf("\n!!!!! No mkfit track with this label.\n\n");
521  ++no_mkf_cnt;
522 
523  auto &ckf_seed = ev.seedTracks_[ckft.label()];
524  print("ckf ", 0, ckft, ev);
525  print("sim ", 0, simt, ev);
526  print("ckf seed", 0, ckf_seed, ev);
527  auto smi = m_seed_map.find(l);
528  if (smi != m_seed_map.end())
529  print("seed with matching label", 0, *smi->second, ev);
530  printf("------------------------------------------------------\n");
531  }
532  }
533 
534  printf("mkFit found %d, matching_label=%d, less_hits=%d, more_hits=%d [algo=%d (%s)]\n",
535  (int) ev.fitTracks_.size(), mkf_cnt, less_hits, more_hits,
537 
538  if (tot_cnt > 0) {
539  printf("\ntobtec tob1/2 tot=%d no_mkf=%d (%f%%)\n",
540  tot_cnt, no_mkf_cnt, 100.0 * no_mkf_cnt / tot_cnt);
541  } else {
542  printf("\nNo CKF tracks with seed hits in TOB1/2 found (need iteration idx 8, TobTec?)\n");
543  }
544 
545  printf("-------------------------------------------------------------------------------------------\n");
546  printf("-------------------------------------------------------------------------------------------\n");
547  printf("\n");
548  }
549 
550 }
std::map< int, Track * > m_sim_map
Definition: Shell.h:74
void ProcessEvent(SeedSelect_e seed_select=SS_UseAll, int selected_seed=-1, int count=1)
Definition: Shell.cc:120
const EventOfCombCandidates & ref_eocc() const
Definition: MkBuilder.h:73
Shell(std::vector< DeadVec > &dv, const std::string &in_file, int start_ev)
Definition: Shell.cc:40
TrackVec m_seeds
Definition: Shell.h:68
void rewind()
Definition: Event.cc:1007
int openRead(const std::string &fname, int expected_n_layers)
Definition: Event.cc:919
const BeamSpot & refBeamSpot() const
void reset(int evtID)
Definition: Event.cc:33
void PrintIterations()
Definition: Shell.cc:302
ret
prodAgent to be discontinued
bool m_remove_duplicates
Definition: Shell.h:66
DataFile * m_data_file
Definition: Shell.h:58
void compactifyHitStorageForBestCand(bool remove_seed_hits, int backward_fit_min_hits)
Definition: MkBuilder.h:62
def replace(string, replacements)
int label() const
Definition: Track.h:188
TrackVec seedTracks_
Definition: Event.h:74
void findTracksCloneEngine(SteeringParams::IterationType_e iteration_dir=SteeringParams::IT_FwdSearch)
Definition: MkBuilder.cc:932
int evtID() const
Definition: Event.h:23
std::map< int, Track * > m_mkf_map
Definition: Shell.h:74
bool CheckMkFitLayerPlanVsReferenceHits(const Track &mkft, const Track &reft, const std::string &name)
Definition: Shell.cc:423
bool m_clean_seeds
Definition: Shell.h:64
void SetDebug(bool b)
Definition: Shell.cc:313
void begin_event(MkJob *job, Event *ev, const char *build_type)
Definition: MkBuilder.cc:194
void resetCurrentSeedTracks()
Definition: Event.cc:910
void Status()
Definition: Shell.cc:74
clean_seeds_func m_seed_cleaner
void SelectIterationIndex(int itidx)
Definition: Shell.cc:283
TEMPL(T2) struct Divides void
Definition: Factorize.h:24
EventOfHits * eoh()
Definition: Shell.h:40
void Run()
Definition: Shell.cc:55
char const * label
MkBuilder * m_builder
Definition: Shell.h:61
std::vector< DeadVec > & m_deadvectors
Definition: Shell.h:57
void loadHitsAndBeamSpot(Event &ev, EventOfHits &eoh)
void FillByLabelMaps_CkfBase()
Definition: Shell.cc:358
void skipNEvents(int n_to_skip)
Definition: Event.cc:1034
std::function< filter_candidates_cf > filter_candidates_func
Definition: FunctionTypes.h:28
EventOfHits * m_eoh
Definition: Shell.h:60
bool g_debug
Definition: Debug.cc:2
void find_tracks_load_seeds(const TrackVec &in_seeds, const bool seeds_sorted)
Definition: MkBuilder.cc:609
Definition: EPCuts.h:4
void setCurrentSeedTracks(const TrackVec &seeds)
Definition: Event.cc:865
static partition_seeds_func get_seed_partitioner(const std::string &name)
void SetCleanSeeds(bool b)
Definition: Shell.cc:314
void GoToEvent(int eid)
Definition: Shell.cc:87
bool m_backward_fit
Definition: Shell.h:65
TrackerInfo TrkInfo
void loadDeads(EventOfHits &eoh, const std::vector< DeadVec > &deadvectors)
Definition: MkStdSeqs.cc:21
void seed_post_cleaning(TrackVec &tv)
Definition: MkBuilder.cc:407
void SetRemoveDuplicates(bool b)
Definition: Shell.cc:316
void fill_hitmask_bool_vectors(int track_algo, std::vector< std::vector< bool >> &layer_masks)
Definition: Event.cc:819
MCHitInfoVec simHitsInfo_
Definition: Event.h:72
IterationsInfo ItrInfo
MkBuilder * builder()
Definition: Shell.h:41
int filter_comb_cands(filter_candidates_func filter, bool attempt_all_cands)
Definition: MkBuilder.cc:301
void NextEvent(int skip=1)
Definition: Shell.cc:116
std::vector< HitVec > layerHits_
Definition: Event.h:70
static const char * algoint_to_cstr(int algo)
Definition: Track.cc:331
TrackVec m_tracks
Definition: Shell.h:69
filter_candidates_func m_post_bkfit_filter
std::vector< Track > TrackVec
int getEtaRegion() const
Definition: Track.h:277
clean_duplicates_func m_duplicate_cleaner
void beginBkwSearch()
Definition: MkBuilder.h:66
void print(std::string_view label, const MeasurementState &s)
Definition: Hit.cc:8
double b
Definition: hdecay.h:120
int nTotalHits() const
Definition: Track.h:528
void export_best_comb_cands(TrackVec &out_vec, bool remove_missing_hits=false)
Definition: MkBuilder.cc:384
Event * m_event
Definition: Shell.h:59
void SetBackwardFit(bool b)
Definition: Shell.cc:315
int m_it_index
Definition: Shell.h:63
float x
int LabelFromHits(Track &t, bool replace, float good_frac)
Definition: Shell.cc:331
void endBkwSearch()
Definition: MkBuilder.h:67
void Compare()
Definition: Shell.cc:446
void SetUseDeadModules(bool b)
Definition: Shell.cc:317
std::vector< std::vector< bool > > m_mask_vector
void SelectIterationAlgo(int algo)
Definition: Shell.cc:291
filter_candidates_func m_pre_bkfit_filter
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::map< int, Track * > m_ckf_map
Definition: Shell.h:74
std::vector< SteeringParams > m_steering_params
int m_evs_in_file
Definition: Shell.h:62
SeedSelect_e
Definition: Shell.h:18
std::map< int, Track * > m_seed_map
Definition: Shell.h:74
void read_in(DataFile &data_file, FILE *in_fp=0)
Definition: Event.cc:207
def exit(msg="")
HitOnTrack getHitOnTrack(int posHitIdx) const
Definition: Track.h:464