CMS 3D CMS Logo

PixelTrackReconstruction.cc
Go to the documentation of this file.
1 #include <vector>
2 
19 
20 using namespace pixeltrackfitting;
21 using edm::ParameterSet;
22 
25  : theHitSetsToken(iC.consumes<RegionsSeedingHitSets>(cfg.getParameter<edm::InputTag>("SeedingHitSets"))),
26  theFitterToken(iC.consumes<PixelFitter>(cfg.getParameter<edm::InputTag>("Fitter"))),
27  theCleanerName(cfg.getParameter<std::string>("Cleaner"))
28 {
29  edm::InputTag filterTag = cfg.getParameter<edm::InputTag>("Filter");
30  if (not filterTag.label().empty()) {
31  theFilterToken = iC.consumes<PixelTrackFilter>(filterTag);
32  }
33 }
34 
36 {
37 }
38 
40  desc.add<edm::InputTag>("SeedingHitSets", edm::InputTag("pixelTracksHitTriplets"));
41  desc.add<edm::InputTag>("Fitter", edm::InputTag("pixelFitterByHelixProjections"));
42  desc.add<edm::InputTag>("Filter", edm::InputTag("pixelTrackFilterByKinematics"));
43  desc.add<std::string>("Cleaner", "pixelTrackCleanerBySharedHits");
44 }
45 
47 {
49  ev.getByToken(theHitSetsToken, hhitSets);
50  const auto& hitSets = *hhitSets;
51 
53  ev.getByToken(theFitterToken, hfitter);
54  const auto& fitter = *hfitter;
55 
56  const PixelTrackFilter *filter = nullptr;
59  ev.getByToken(theFilterToken, hfilter);
60  filter = hfilter.product();
61  }
62 
63  std::vector<const TrackingRecHit *> hits;hits.reserve(4);
64  for(const auto& regionHitSets: hitSets) {
65  const TrackingRegion& region = regionHitSets.region();
66 
67  for(const SeedingHitSet& tuplet: regionHitSets) {
69  auto nHits = tuplet.size(); hits.resize(nHits);
70  for (unsigned int iHit = 0; iHit < nHits; ++iHit) hits[iHit] = tuplet[iHit];
71 
72  // fitting
73  std::unique_ptr<reco::Track> track = fitter.run(hits, region);
74  if (!track) continue;
75 
76  if (filter) {
77  if (!(*filter)(track.get(), hits)) {
78  continue;
79  }
80  }
81 
82  // add tracks
83  tracks.emplace_back(track.release(), tuplet);
84  }
85  }
86 
87  // skip ovelrapped tracks
88  if(!theCleanerName.empty()) {
90  es.get<PixelTrackCleaner::Record>().get(theCleanerName, hcleaner);
91  const auto& cleaner = *hcleaner;
92  if(cleaner.fast())
93  cleaner.cleanTracks(tracks);
94  else
95  tracks = PixelTrackCleanerWrapper(&cleaner).clean(tracks);
96  }
97 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
pixeltrackfitting::TracksWithTTRHs clean(const pixeltrackfitting::TracksWithTTRHs &initialT_TTRHs) const
edm::EDGetTokenT< PixelTrackFilter > theFilterToken
bool ev
std::vector< TrackWithTTRHs > TracksWithTTRHs
PixelTrackReconstruction(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void run(pixeltrackfitting::TracksWithTTRHs &tah, edm::Event &ev, const edm::EventSetup &es)
T const * product() const
Definition: Handle.h:81
std::string const & label() const
Definition: InputTag.h:36
HLT enums.
T get() const
Definition: EventSetup.h:62
bool isUninitialized() const
Definition: EDGetToken.h:73
static void fillDescriptions(edm::ParameterSetDescription &desc)
edm::EDGetTokenT< RegionsSeedingHitSets > theHitSetsToken
edm::EDGetTokenT< PixelFitter > theFitterToken