CMS 3D CMS Logo

PixelTrackReconstruction.cc
Go to the documentation of this file.
2 
9 
11 
13 
15 
18 
22 
24 
25 #include <vector>
26 
27 using namespace pixeltrackfitting;
28 using namespace ctfseeding;
29 using edm::ParameterSet;
30 
33  : theHitSetsToken(iC.consumes<RegionsSeedingHitSets>(cfg.getParameter<edm::InputTag>("SeedingHitSets"))),
34  theFitterToken(iC.consumes<PixelFitter>(cfg.getParameter<edm::InputTag>("Fitter"))),
35  theCleanerName(cfg.getParameter<std::string>("Cleaner"))
36 {
37  edm::InputTag filterTag = cfg.getParameter<edm::InputTag>("Filter");
38  if(filterTag.label() != "") {
39  theFilterToken = iC.consumes<PixelTrackFilter>(filterTag);
40  }
41 }
42 
44 {
45 }
46 
48  desc.add<edm::InputTag>("SeedingHitSets", edm::InputTag("pixelTracksHitTriplets"));
49  desc.add<edm::InputTag>("Fitter", edm::InputTag("pixelFitterByHelixProjections"));
50  desc.add<edm::InputTag>("Filter", edm::InputTag("pixelTrackFilterByKinematics"));
51  desc.add<std::string>("Cleaner", "pixelTrackCleanerBySharedHits");
52 }
53 
55 {
57  ev.getByToken(theHitSetsToken, hhitSets);
58  const auto& hitSets = *hhitSets;
59 
61  ev.getByToken(theFitterToken, hfitter);
62  const auto& fitter = *hfitter;
63 
64  const PixelTrackFilter *filter = nullptr;
67  ev.getByToken(theFilterToken, hfilter);
68  filter = hfilter.product();
69  }
70 
71  std::vector<const TrackingRecHit *> hits;hits.reserve(4);
72  for(const auto& regionHitSets: hitSets) {
73  const TrackingRegion& region = regionHitSets.region();
74 
75  for(const SeedingHitSet& tuplet: regionHitSets) {
77  auto nHits = tuplet.size(); hits.resize(nHits);
78  for (unsigned int iHit = 0; iHit < nHits; ++iHit) hits[iHit] = tuplet[iHit];
79 
80  // fitting
81  std::unique_ptr<reco::Track> track = fitter.run(hits, region);
82  if (!track) continue;
83 
84  if (filter) {
85  if (!(*filter)(track.get(), hits)) {
86  continue;
87  }
88  }
89 
90  // add tracks
91  tracks.emplace_back(track.release(), tuplet);
92  }
93  }
94 
95  // skip ovelrapped tracks
96  if(!theCleanerName.empty()) {
98  es.get<PixelTrackCleaner::Record>().get(theCleanerName, hcleaner);
99  const auto& cleaner = *hcleaner;
100  if(cleaner.fast())
101  cleaner.cleanTracks(tracks);
102  else
103  tracks = PixelTrackCleanerWrapper(&cleaner).clean(tracks);
104  }
105 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
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
const T & get() const
Definition: EventSetup.h:55
std::string const & label() const
Definition: InputTag.h:36
HLT enums.
bool isUninitialized() const
Definition: EDGetToken.h:73
static void fillDescriptions(edm::ParameterSetDescription &desc)
edm::EDGetTokenT< RegionsSeedingHitSets > theHitSetsToken
edm::EDGetTokenT< PixelFitter > theFitterToken