Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-17 22:59:01

0001 #include "RecoTracker/MkFitCMS/standalone/Shell.h"
0002 
0003 #include "RecoTracker/MkFitCore/src/Debug.h"
0004 
0005 // #include "RecoTracker/MkFitCore/src/Matriplex/MatriplexCommon.h"
0006 
0007 #include "RecoTracker/MkFitCMS/interface/runFunctions.h"
0008 
0009 #include "RecoTracker/MkFitCore/interface/HitStructures.h"
0010 #include "RecoTracker/MkFitCore/interface/MkBuilder.h"
0011 #include "RecoTracker/MkFitCore/src/MkFitter.h"
0012 #include "RecoTracker/MkFitCMS/interface/MkStdSeqs.h"
0013 #include "RecoTracker/MkFitCMS/standalone/MkStandaloneSeqs.h"
0014 
0015 #include "RecoTracker/MkFitCore/interface/Config.h"
0016 #include "RecoTracker/MkFitCore/standalone/ConfigStandalone.h"
0017 
0018 #include "RecoTracker/MkFitCore/standalone/Event.h"
0019 
0020 #include "RecoTracker/MkFitCore/interface/TrackerInfo.h"
0021 
0022 #include "TROOT.h"
0023 #include "TRint.h"
0024 
0025 #ifdef WITH_REVE
0026 #include "TRandom.h"
0027 #include "ROOT/REveJetCone.hxx"
0028 #include "ROOT/REveManager.hxx"
0029 #include "ROOT/REveScene.hxx"
0030 #include "ROOT/REveBoxSet.hxx"
0031 #endif
0032 
0033 #include "oneapi/tbb/task_arena.h"
0034 
0035 #include <vector>
0036 
0037 // clang-format off
0038 
0039 namespace {
0040   constexpr int algos[] = {4, 22, 23, 5, 24, 7, 8, 9, 10, 6};  // 10 iterations
0041   constexpr int n_algos = sizeof(algos) / sizeof(int);
0042 
0043   const char* b2a(bool b) { return b ? "true" : "false"; }
0044 }
0045 
0046 namespace mkfit {
0047 
0048   Shell::Shell(std::vector<DeadVec> &dv, const std::string &in_file, int start_ev)
0049     : m_deadvectors(dv)
0050   {
0051     m_eoh = new EventOfHits(Config::TrkInfo);
0052     m_builder = new MkBuilder(Config::silent);
0053 
0054     m_backward_fit = Config::backwardFit;
0055 
0056     m_data_file = new DataFile;
0057     m_event = new Event(0, Config::TrkInfo.n_layers());
0058 
0059     if ( ! in_file.empty() && Config::nEvents > 0) {
0060       m_evs_in_file = m_data_file->openRead(in_file, Config::TrkInfo.n_layers());
0061       GoToEvent(start_ev);
0062     } else {
0063       printf("Shell initialized but the %s, running on an empty Event.\n",
0064       in_file.empty() ? "input-file not specified" : "requested number of events to process is 0");
0065     }
0066   }
0067 
0068   Shell::~Shell() {
0069     delete m_event;
0070     delete m_data_file;
0071     delete m_builder;
0072     delete m_eoh;
0073     delete gApplication;
0074   }
0075 
0076   void Shell::Run() {
0077     std::vector<const char *> argv = { "mkFit", "-l" };
0078     int argc = argv.size();
0079     gApplication = new TRint("mkFit-shell", &argc, const_cast<char**>(argv.data()));
0080 
0081     char buf[256];
0082     sprintf(buf, "mkfit::Shell &s = * (mkfit::Shell*) %p;", this);
0083     gROOT->ProcessLine(buf);
0084     printf("Shell &s variable is set: ");
0085     gROOT->ProcessLine("s");
0086 
0087     gApplication->Run(true);
0088     printf("Shell::Run finished\n");
0089   }
0090 
0091   void Shell::Status() {
0092     printf("On event %d, selected iteration index %d, algo %d - %s\n"
0093           "  debug = %s, use_dead_modules = %s\n"
0094            "  clean_seeds = %s, backward_fit = %s, remove_duplicates = %s\n",
0095            m_event->evtID(), m_it_index, algos[m_it_index], TrackBase::algoint_to_cstr(algos[m_it_index]),
0096            b2a(g_debug), b2a(Config::useDeadModules),
0097            b2a(m_clean_seeds), b2a(m_backward_fit), b2a(m_remove_duplicates));
0098   }
0099 
0100   TrackerInfo* Shell::tracker_info() { return &Config::TrkInfo; }
0101 
0102   //===========================================================================
0103   // Event navigation / processing
0104   //===========================================================================
0105 
0106   void Shell::GoToEvent(int eid) {
0107     if (eid < 1) {
0108       fprintf(stderr, "Requested event %d is less than 1 -- 1 is the first event, %d is total number of events in file\n",
0109              eid, m_evs_in_file);
0110       throw std::runtime_error("event out of range");
0111     }
0112     if (eid > m_evs_in_file) {
0113       fprintf(stderr, "Requested event %d is grater than total number of events in file %d\n",
0114              eid, m_evs_in_file);
0115       throw std::runtime_error("event out of range");
0116     }
0117 
0118     int pos = m_event->evtID();
0119     if (eid > pos) {
0120       m_data_file->skipNEvents(eid - pos - 1);
0121     } else {
0122       m_data_file->rewind();
0123       m_data_file->skipNEvents(eid - 1);
0124     }
0125     m_event->resetCurrentSeedTracks(); // left after ProcessEvent() for debugging etc
0126     m_event->reset(eid);
0127     m_event->read_in(*m_data_file);
0128     StdSeq::loadHitsAndBeamSpot(*m_event, *m_eoh);
0129     if (Config::useDeadModules) {
0130       StdSeq::loadDeads(*m_eoh, m_deadvectors);
0131     }
0132 
0133     printf("At event %d\n", eid);
0134   }
0135 
0136   void Shell::NextEvent(int skip) {
0137     GoToEvent(m_event->evtID() + skip);
0138   }
0139 
0140   void Shell::ProcessEvent(SeedSelect_e seed_select, int selected_seed, int count) {
0141     // count is only used for SS_IndexPreCleaning and SS_IndexPostCleaning.
0142     //       There are no checks for upper bounds, ie, if requested seeds beyond the first one exist.
0143 
0144     const IterationConfig &itconf = Config::ItrInfo[m_it_index];
0145     IterationMaskIfc mask_ifc;
0146     m_event->fill_hitmask_bool_vectors(itconf.m_track_algorithm, mask_ifc.m_mask_vector);
0147 
0148     m_seeds.clear();
0149     m_tracks.clear();
0150 
0151     {
0152       int n_algo = 0; // seeds are grouped by algo
0153       for (auto &s : m_event->seedTracks_) {
0154         if (s.algoint() == itconf.m_track_algorithm) {
0155           if (seed_select == SS_UseAll || seed_select == SS_IndexPostCleaning) {
0156             m_seeds.push_back(s);
0157           } else if (seed_select == SS_Label && s.label() == selected_seed) {
0158             m_seeds.push_back(s);
0159             if (--count <= 0)
0160               break;
0161           } else if (seed_select == SS_IndexPreCleaning && n_algo >= selected_seed) {
0162             m_seeds.push_back(s);
0163             if (--count <= 0)
0164               break;
0165           }
0166           ++n_algo;
0167         } else if (n_algo > 0)
0168           break;
0169       }
0170     }
0171 
0172     printf("Shell::ProcessEvent running over %d seeds\n", (int) m_seeds.size());
0173 
0174     // Equivalent to run_OneIteration(...) without MkBuilder::release_memory().
0175     // If seed_select == SS_IndexPostCleaning the given seed is picked after cleaning.
0176     {
0177       const TrackerInfo &trackerInfo = Config::TrkInfo;
0178       const EventOfHits &eoh = *m_eoh;
0179       const IterationMaskIfcBase &it_mask_ifc = mask_ifc;
0180       MkBuilder &builder = *m_builder;
0181       TrackVec &seeds = m_seeds;
0182       TrackVec &out_tracks = m_tracks;
0183       bool do_seed_clean = m_clean_seeds;
0184       bool do_backward_fit = m_backward_fit;
0185       bool do_remove_duplicates = m_remove_duplicates;
0186 
0187       MkJob job({trackerInfo, itconf, eoh, eoh.refBeamSpot(), &it_mask_ifc});
0188 
0189       builder.begin_event(&job, m_event, __func__);
0190 
0191       // Seed cleaning not done on all iterations.
0192       do_seed_clean = m_clean_seeds && itconf.m_seed_cleaner;
0193 
0194       if (do_seed_clean) {
0195         itconf.m_seed_cleaner(seeds, itconf, eoh.refBeamSpot());
0196         printf("Shell::ProcessEvent post seed-cleaning: %d seeds\n", (int) m_seeds.size());
0197       } else {
0198         printf("Shell::ProcessEvent no seed-cleaning\n");
0199       }
0200 
0201       // Check nans in seeds -- this should not be needed when Slava fixes
0202       // the track parameter coordinate transformation.
0203       builder.seed_post_cleaning(seeds);
0204 
0205       if (seed_select == SS_IndexPostCleaning) {
0206         int seed_size = (int) seeds.size();
0207         if (selected_seed >= 0 && selected_seed < seed_size) {
0208           if (selected_seed + count >= seed_size) {
0209             count = seed_size - selected_seed;
0210             printf("  -- selected seed_index + count > seed vector size after cleaning -- trimming count to %d\n",
0211                    count);
0212           }
0213           for (int i = 0; i < count; ++i)
0214             seeds[i] = seeds[selected_seed + i];
0215           seeds.resize(count);
0216         } else {
0217           seeds.clear();
0218         }
0219       }
0220 
0221       if (seeds.empty()) {
0222         if (seed_select != SS_UseAll)
0223           printf("Shell::ProcessEvent requested seed not found among seeds of the selected iteration.\n");
0224         else
0225           printf("Shell::ProcessEvent no seeds found.\n");
0226         return;
0227       }
0228 
0229       if (itconf.m_requires_seed_hit_sorting) {
0230         for (auto &s : seeds)
0231           s.sortHitsByLayer();  // sort seed hits for the matched hits (I hope it works here)
0232       }
0233 
0234       m_event->setCurrentSeedTracks(seeds);
0235 
0236       builder.find_tracks_load_seeds(seeds, do_seed_clean);
0237 
0238       builder.findTracksCloneEngine();
0239 
0240       printf("Shell::ProcessEvent post fwd search: %d comb-cands\n", builder.ref_eocc().size());
0241 
0242       // Pre backward-fit filtering.
0243       filter_candidates_func pre_filter;
0244       if (do_backward_fit && itconf.m_pre_bkfit_filter)
0245         pre_filter = [&](const TrackCand &tc, const MkJob &jb) -> bool {
0246           return itconf.m_pre_bkfit_filter(tc, jb) && StdSeq::qfilter_nan_n_silly<TrackCand>(tc, jb);
0247         };
0248       else if (itconf.m_pre_bkfit_filter)
0249         pre_filter = itconf.m_pre_bkfit_filter;
0250       else if (do_backward_fit)
0251         pre_filter = StdSeq::qfilter_nan_n_silly<TrackCand>;
0252       // pre_filter can be null if we are not doing backward fit as nan_n_silly will be run below.
0253       if (pre_filter)
0254         builder.filter_comb_cands(pre_filter, true);
0255 
0256       printf("Shell::ProcessEvent post pre-bkf-filter (%s) and nan-filter (%s) filter: %d comb-cands\n",
0257              b2a(bool(itconf.m_pre_bkfit_filter)), b2a(do_backward_fit), builder.ref_eocc().size());
0258 
0259       job.switch_to_backward();
0260 
0261       if (do_backward_fit) {
0262         if (itconf.m_backward_search) {
0263           builder.compactifyHitStorageForBestCand(itconf.m_backward_drop_seed_hits, itconf.m_backward_fit_min_hits);
0264         }
0265 
0266         builder.backwardFit();
0267 
0268         if (itconf.m_backward_search) {
0269           builder.beginBkwSearch();
0270           builder.findTracksCloneEngine(SteeringParams::IT_BkwSearch);
0271         }
0272 
0273         printf("Shell::ProcessEvent post backward fit / search: %d comb-cands\n", builder.ref_eocc().size());
0274       }
0275 
0276       // Post backward-fit filtering.
0277       filter_candidates_func post_filter;
0278       if (do_backward_fit && itconf.m_post_bkfit_filter)
0279         post_filter = [&](const TrackCand &tc, const MkJob &jb) -> bool {
0280           return itconf.m_post_bkfit_filter(tc, jb) && StdSeq::qfilter_nan_n_silly<TrackCand>(tc, jb);
0281         };
0282       else
0283         post_filter = StdSeq::qfilter_nan_n_silly<TrackCand>;
0284       // post_filter is always at least doing nan_n_silly filter.
0285       builder.filter_comb_cands(post_filter, true);
0286 
0287       printf("Shell::ProcessEvent post post-bkf-filter (%s) and nan-filter (true): %d comb-cands\n",
0288              b2a(do_backward_fit && itconf.m_post_bkfit_filter), builder.ref_eocc().size());
0289 
0290       if (do_backward_fit && itconf.m_backward_search)
0291         builder.endBkwSearch();
0292 
0293       builder.export_best_comb_cands(out_tracks, true);
0294 
0295       if (do_remove_duplicates && itconf.m_duplicate_cleaner) {
0296         itconf.m_duplicate_cleaner(out_tracks, itconf);
0297       }
0298 
0299       printf("Shell::ProcessEvent post remove-duplicates: %d comb-cands\n", (int) out_tracks.size());
0300 
0301       // Do not clear ... useful for debugging / printouts!
0302       // m_event->resetCurrentSeedTracks();
0303 
0304       builder.end_event();
0305     }
0306 
0307     printf("Shell::ProcessEvent found %d tracks, number of seeds at end %d\n",
0308            (int) m_tracks.size(), (int) m_seeds.size());
0309   }
0310 
0311   //===========================================================================
0312   // Iteration selection
0313   //===========================================================================
0314 
0315   void Shell::SelectIterationIndex(int itidx) {
0316     if (itidx < 0 || itidx >= n_algos) {
0317       fprintf(stderr, "Requested iteration index out of range [%d, %d)", 0, n_algos);
0318       throw std::runtime_error("iteration index out of range");
0319     }
0320     m_it_index = itidx;
0321   }
0322 
0323   void Shell::SelectIterationAlgo(int algo) {
0324     for (int i = 0; i < n_algos; ++i) {
0325       if (algos[i] == algo) {
0326         m_it_index = i;
0327         return;
0328       }
0329     }
0330     fprintf(stderr, "Requested algo %d not found", algo);
0331     throw std::runtime_error("algo not found");
0332   }
0333 
0334   void Shell::PrintIterations() {
0335     printf("Shell::PrintIterations selected index = %d, %d iterations hardcoded as\n",
0336             m_it_index, n_algos);
0337     for (int i = 0; i < n_algos; ++i)
0338       printf("%d %2d %s\n", i, algos[i], TrackBase::algoint_to_cstr(algos[i]));
0339   }
0340 
0341   //===========================================================================
0342   // Flags / status setters
0343   //===========================================================================
0344 
0345   void Shell::SetDebug(bool b) { g_debug = b; }
0346   void Shell::SetCleanSeeds(bool b) { m_clean_seeds = b; }
0347   void Shell::SetBackwardFit(bool b) { m_backward_fit = b; }
0348   void Shell::SetRemoveDuplicates(bool b) { m_remove_duplicates = b; }
0349   void Shell::SetUseDeadModules(bool b) { Config::useDeadModules = b; }
0350 
0351   //===========================================================================
0352   // Analysis helpers
0353   //===========================================================================
0354 
0355   /*
0356     sim tracks are written to .bin files with a label equal to its own index.
0357     reco tracks labels are seed indices.
0358     seed labels are sim track indices
0359     --
0360     mkfit labels are seed indices in given iteration after cleaning (at seed load-time).
0361           This is no longer true -- was done like that in branch where this code originated from.
0362           It seems the label is the same as seed label.
0363   */
0364 
0365   int Shell::LabelFromHits(Track &t, bool replace, float good_frac) {
0366     std::map<int, int> lab_cnt;
0367     for (int hi = 0; hi < t.nTotalHits(); ++hi) {
0368       auto hot = t.getHitOnTrack(hi);
0369       if (hot.index < 0)
0370         continue;
0371       const Hit &h = m_event->layerHits_[hot.layer][hot.index];
0372       int hl = m_event->simHitsInfo_[h.mcHitID()].mcTrackID_;
0373       if (hl >= 0)
0374         ++lab_cnt[hl];
0375     }
0376     int max_c = -1, max_l = -1;
0377     for (auto& x : lab_cnt) {
0378       if (x.second > max_c) {
0379         max_l = x.first;
0380         max_c = x.second;
0381       }
0382     }
0383     bool success = max_c >= good_frac * t.nFoundHits();
0384     int relabel = success ? max_l : -1;
0385     // printf("found_hits=%d, best_lab %d (%d hits), existing label=%d (replace flag=%s)\n",
0386     //        t.nFoundHits(), max_l, max_c, t.label(), b2a(replace));
0387     if (replace)
0388         t.setLabel(relabel);
0389     return relabel;
0390   }
0391 
0392   void Shell::FillByLabelMaps_CkfBase() {
0393     Event &ev = *m_event;
0394     const int track_algo = Config::ItrInfo[m_it_index].m_track_algorithm;
0395 
0396     m_ckf_map.clear();
0397     m_sim_map.clear();
0398     m_seed_map.clear();
0399     m_mkf_map.clear();
0400 
0401     // Pick ckf tracks with right algo and a good label.
0402     int rec_algo_match = 0;
0403     for (auto &t : ev.cmsswTracks_) {
0404       if (t.algoint() != track_algo)
0405         continue;
0406       ++rec_algo_match;
0407       int label = LabelFromHits(t, false, 0.5);
0408       if (label >= 0) {
0409         m_ckf_map.insert(std::make_pair(label, &t));
0410       }
0411     }
0412 
0413     // Pick sim tracks with labels found by ckf.
0414     for (auto &t : ev.simTracks_) {
0415       if (t.label() >= 0 && m_ckf_map.find(t.label()) != m_ckf_map.end()) {
0416         m_sim_map.insert(std::make_pair(t.label(), &t));
0417       }
0418     }
0419 
0420     // Pick seeds with right algo and a label found by ckf.
0421     for (auto &t : ev.seedTracks_) {
0422       if (t.algoint() == track_algo && t.label() >= 0 && m_ckf_map.find(t.label()) != m_ckf_map.end()) {
0423         m_seed_map.insert(std::make_pair(t.label(), &t));
0424       }
0425     }
0426     // Some seeds seem to be labeled -1, try fixing when not otherwise found.
0427     for (auto &t : ev.seedTracks_) {
0428       if (t.algoint() == track_algo && t.label() == -1) {
0429         int lab = LabelFromHits(t, true, 0.5);
0430         if (lab >= 0 && m_seed_map.find(lab) == m_seed_map.end()) {
0431           if (m_ckf_map.find(lab) != m_ckf_map.end()) {
0432             m_seed_map.insert(std::make_pair(t.label(), &t));
0433             printf("Saved seed with label -1 -> %d\n", lab);
0434           }
0435         }
0436       }
0437     }
0438 
0439     // Pick mkfit tracks, label by 
0440     for (auto &t : m_tracks) {
0441       int label = LabelFromHits(t, false, 0.5);
0442       if (label >= 0) {
0443         m_mkf_map.insert(std::make_pair(label, &t));
0444       }
0445     }
0446 
0447     printf("Shell::FillByLabelMaps reporting tracks with label >= 0, algo=%d (%s): "
0448            "ckf: %d of %d (same algo=%d)), sim: %d of %d, seed: %d of %d, mkfit: %d w/label of %d\n",
0449            track_algo, TrackBase::algoint_to_cstr(track_algo),
0450            (int) m_ckf_map.size(), (int) ev.cmsswTracks_.size(), rec_algo_match,
0451            (int) m_sim_map.size(), (int) ev.simTracks_.size(),
0452            (int) m_seed_map.size(), (int) m_seeds.size(),
0453            (int) m_mkf_map.size(), (int) m_tracks.size()
0454     );
0455   }
0456 
0457   bool Shell::CheckMkFitLayerPlanVsReferenceHits(const Track &mkft, const Track &reft, const std::string &name) {
0458     // Check if all hit-layers of a reference track reft are in the mkfit layer plan.
0459     // Returns true if all layers are in the plan.
0460     // String name is printed in front of label, expected to be SIMK or CKF.
0461     const IterationConfig &itconf = Config::ItrInfo[m_it_index];
0462     auto lp = itconf.m_steering_params[ mkft.getEtaRegion() ].m_layer_plan;
0463     bool ret = true;
0464     for (int hi = 0; hi < reft.nTotalHits(); ++hi) {
0465       auto hot = reft.getHitOnTrack(hi);
0466       if (std::find_if(lp.begin(), lp.end(), [=](auto &x){ return x.m_layer == hot.layer; }) == lp.end())
0467       {
0468         printf("CheckMkfLayerPlanVsCkfHits: layer %d not in layer plan for region %d, %s label=%d\n",
0469                 hot.layer, mkft.getEtaRegion(), name.c_str(), reft.label());
0470         ret = false;
0471       }
0472     }
0473     return ret;
0474   }
0475 
0476   //===========================================================================
0477   // Analysis drivers / main functions / Comparators
0478   //===========================================================================
0479 
0480   void Shell::Compare() {
0481     Event &ev = *m_event;
0482     const IterationConfig &itconf = Config::ItrInfo[m_it_index];
0483 
0484     FillByLabelMaps_CkfBase();
0485 
0486     printf("------------------------------------------------------\n");
0487 
0488     const bool print_all_def = false;
0489     int mkf_cnt=0, less_hits=0, more_hits=0;
0490 
0491     // TOBTEC: look for rec-seeds with hits in tob1 and 2 only.
0492     int tot_cnt = 0, no_mkf_cnt = 0;
0493 
0494     for (auto& [l, simt_ptr]: m_sim_map)
0495     {
0496       auto &simt = * simt_ptr;
0497       auto &ckft = * m_ckf_map[l];
0498       auto mi = m_mkf_map.find(l);
0499 
0500       bool print_all = print_all_def;
0501 
0502       // TOBTEC: look for rec-seeds with hits in tob1 and 2 only.
0503       bool select = true;
0504       {
0505         auto &ckf_seed = ev.seedTracks_[ckft.label()];
0506         for (int hi = 0; hi < ckf_seed.nTotalHits(); ++hi) {
0507           const HitOnTrack hot = ckf_seed.getHitOnTrack(hi);
0508           if (hot.index >= 0 && (hot.layer < 10 || hot.layer > 13)) {
0509             select = false;
0510             break;
0511           }
0512         }
0513       }
0514       if ( ! select) continue;
0515 
0516       ++tot_cnt;
0517       //print_all = true;
0518 
0519       if (mi != m_mkf_map.end())
0520       {
0521         auto &mkft = * mi->second;
0522         mkf_cnt++;
0523         if (mkft.nFoundHits() < ckft.nFoundHits()) ++less_hits;
0524         if (mkft.nFoundHits() > ckft.nFoundHits()) ++more_hits;
0525 
0526         CheckMkFitLayerPlanVsReferenceHits(mkft, ckft, "CKF");
0527         // CheckMkFitLayerPlanVsReferenceHits(mkft, simt, "SIM");
0528 
0529         (void) print_all;
0530         if (/* itconf.m_track_algorithm == 10 ||*/ print_all) {
0531           // ckf label is wrong when validation is on (even quality val) for mixedTriplet, pixelless and tobtec
0532           // as seed tracks get removed for non-mkfit iterations and indices from rec-tracks are no longer valid.
0533           auto &ckf_seed = ev.seedTracks_[ckft.label()];
0534           auto &mkf_seed = m_seeds[mkft.label()];
0535           print("ckf  ", 0, ckft, ev);
0536           print("mkfit", 0, mkft, ev);
0537           print("sim  ", 0, simt, ev);
0538           print("ckf seed", 0, ckf_seed, ev);
0539           print("mkf seed", 0, mkf_seed, ev);
0540           printf("------------------------------------------------------\n");
0541 
0542           TrackVec ssss;
0543           ssss.push_back(mkf_seed);
0544 
0545           IterationSeedPartition pppp(1);
0546           IterationConfig::get_seed_partitioner("phase1:1:debug")(Config::TrkInfo, ssss, *m_eoh, pppp);
0547 
0548           printf("------------------------------------------------------\n");
0549           printf("\n");
0550         }
0551       }
0552       else
0553       {
0554         printf("\n!!!!! No mkfit track with this label.\n\n");
0555         ++no_mkf_cnt;
0556 
0557         auto &ckf_seed = ev.seedTracks_[ckft.label()];
0558         print("ckf ", 0, ckft, ev);
0559         print("sim ", 0, simt, ev);
0560         print("ckf seed", 0, ckf_seed, ev);
0561         auto smi = m_seed_map.find(l);
0562         if (smi != m_seed_map.end())
0563           print("seed with matching label", 0, *smi->second, ev);
0564         printf("------------------------------------------------------\n");
0565       }
0566     }
0567 
0568     printf("mkFit found %d, matching_label=%d, less_hits=%d, more_hits=%d  [algo=%d (%s)]\n",
0569            (int) ev.fitTracks_.size(), mkf_cnt, less_hits, more_hits,
0570            itconf.m_track_algorithm, TrackBase::algoint_to_cstr(itconf.m_track_algorithm));
0571 
0572     if (tot_cnt > 0) {
0573       printf("\ntobtec tob1/2 tot=%d no_mkf=%d (%f%%)\n",
0574             tot_cnt, no_mkf_cnt, 100.0 * no_mkf_cnt / tot_cnt);
0575     } else {
0576       printf("\nNo CKF tracks with seed hits in TOB1/2 found (need iteration idx 8, TobTec?)\n");
0577     }
0578 
0579     printf("-------------------------------------------------------------------------------------------\n");
0580     printf("-------------------------------------------------------------------------------------------\n");
0581     printf("\n");
0582   }
0583 
0584   //===========================================================================
0585   // Visualization helpers
0586   //===========================================================================
0587 
0588 #ifdef WITH_REVE
0589 
0590   void Shell::ShowTracker() {
0591     namespace REX = ROOT::Experimental;
0592     auto eveMng = REX::REveManager::Create();
0593     eveMng->AllowMultipleRemoteConnections(false, false);
0594 
0595     {
0596       REX::REveElement *holder = new REX::REveElement("Jets");
0597 
0598       int N_Jets = 4;
0599       TRandom &r = *gRandom;
0600 
0601       //const Double_t kR_min = 240;
0602       const Double_t kR_max = 250;
0603       const Double_t kZ_d   = 300;
0604       for (int i = 0; i < N_Jets; i++)
0605       {
0606           auto jet = new REX::REveJetCone(Form("Jet_%d",i ));
0607           jet->SetCylinder(2*kR_max, 2*kZ_d);
0608           jet->AddEllipticCone(r.Uniform(-0.5, 0.5), r.Uniform(0, TMath::TwoPi()),
0609                               0.1, 0.2);
0610           jet->SetFillColor(kRed + 4);
0611           jet->SetLineColor(kBlack);
0612           jet->SetMainTransparency(90);
0613 
0614           holder->AddElement(jet);
0615       }
0616       eveMng->GetEventScene()->AddElement(holder);
0617     }
0618 
0619     auto &ti = Config::TrkInfo;
0620     for (int l = 0; l < ti.n_layers(); ++l) {
0621       auto &li = ti[l];
0622       auto* bs = new REX::REveBoxSet(Form("Layer %d", l));
0623       bs->Reset(REX::REveBoxSet::kBT_InstancedScaledRotated, true, li.n_modules());
0624       bs->SetMainColorPtr(new Color_t);
0625       bs->UseSingleColor();
0626       if (li.is_pixel())
0627         bs->SetMainColor(li.is_barrel() ? kBlue - 3 : kCyan - 3);
0628       else
0629         bs->SetMainColor(li.is_barrel() ? kMagenta - 3 : kGreen - 3);
0630 
0631       float t[16];
0632       t[3] = t[7] = t[11] = 0;
0633       t[15] = 1;
0634       for (int m = 0; m < li.n_modules(); ++m) {
0635         auto &mi = li.module_info(m);
0636         auto &si = li.module_shape(mi.shapeid);
0637 
0638         auto &x = mi.xdir;
0639         t[0] = x[0] * si.dx1;
0640         t[1] = x[1] * si.dx1;
0641         t[2] = x[2] * si.dx1;
0642         auto y = mi.calc_ydir();
0643         t[4] = y[0] * si.dy;
0644         t[5] = y[1] * si.dy;
0645         t[6] = y[2] * si.dy;
0646         auto &z = mi.zdir;
0647         t[8] = z[0] * si.dz;
0648         t[9] = z[1] * si.dz;
0649         t[10] = z[2] * si.dz;
0650         auto &p = mi.pos;
0651         t[12] = p[0];
0652         t[13] = p[1];
0653         t[14] = p[2];
0654 
0655         bs->AddInstanceMat4(t);
0656       }
0657       bs->SetMainTransparency(60);
0658       bs->RefitPlex();
0659 
0660       eveMng->GetEventScene()->AddElement(bs);
0661     }
0662     eveMng->Show();
0663   }
0664 
0665 #endif // WITH_REVE
0666 
0667 }