CMS 3D CMS Logo

FastTrackerRecHitMaskProducer.cc
Go to the documentation of this file.
1 // system includes
2 #include <memory>
3 
4 // framework stuff
11 
12 // data formats
18 
19 // other
21 
22 
24 {
25  public:
26 
28 
30 
31  virtual void produce(edm::Event&, const edm::EventSetup&) override;
32 
33 
34  private:
35 
36  // an alias
37  using QualityMaskCollection = std::vector<unsigned char>;
38 
39  // tokens
42 
45 
47 
49 
50 };
51 
53  : minNumberOfLayersWithMeasBeforeFiltering_(iConfig.getParameter<int>("minNumberOfLayersWithMeasBeforeFiltering"))
54  , trackQuality_(reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality")))
55  , trajectories_(iConfig.getParameter<edm::InputTag>("trajectories"),consumesCollector())
56  , recHits_(consumes<FastTrackerRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHits")))
57 {
58 
59  produces<std::vector<bool> >();
60 
61  auto const & classifier = iConfig.getParameter<edm::InputTag>("trackClassifier");
62  if ( !classifier.label().empty())
63  srcQuals = consumes<QualityMaskCollection>(classifier);
64 
65  auto const & oldHitRemovalInfo = iConfig.getParameter<edm::InputTag>("oldHitRemovalInfo");
66  if (!oldHitRemovalInfo.label().empty()) {
67  oldHitMaskToken_ = consumes<std::vector<bool> >(oldHitRemovalInfo);
68  }
69 }
70 
72 {
73  // the product
74 
75  std::unique_ptr<std::vector<bool> > collectedHits(new std::vector<bool>());
76 
77  // input
78 
80  iEvent.getByToken(recHits_,recHits);
81 
83  edm::Handle<std::vector<bool> > oldHitMasks;
84  iEvent.getByToken(oldHitMaskToken_,oldHitMasks);
85  collectedHits->insert(collectedHits->begin(),oldHitMasks->begin(),oldHitMasks->end());
86  }
87  collectedHits->resize(recHits->size(),false);
88 
89  auto const & tracks = trajectories_.tracks(iEvent);
90 
91  QualityMaskCollection const * pquals=nullptr;
92  if (!srcQuals.isUninitialized()) {
94  iEvent.getByToken(srcQuals, hqual);
95  pquals = hqual.product();
96  }
97 
98  // required quality
99  unsigned char qualMask = ~0;
101 
102  // loop over tracks and mask hits of selected tracks
103  for (auto i=0U; i<tracks.size(); ++i){
104  const reco::Track & track = tracks[i];
105  bool goodTk = (pquals) ? (*pquals)[i] & qualMask : track.quality(trackQuality_);
106  if ( !goodTk) continue;
108  for (auto hitIt = track.recHitsBegin() ; hitIt != track.recHitsEnd(); ++hitIt) {
109  if(!(*hitIt)->isValid())
110  continue;
111  const FastTrackerRecHit & hit = static_cast<const FastTrackerRecHit &>(*(*hitIt));
112  // note: for matched hits nIds() returns 2, otherwise 1
113  for(unsigned id_index = 0;id_index < hit.nIds();id_index++){
114  (*collectedHits)[unsigned(hit.id(id_index))] = true;
115  }
116  }
117 
118  }
119 
120  iEvent.put(std::move(collectedHits));
121 
122 }
123 
T getParameter(std::string const &) const
std::vector< unsigned char > QualityMaskCollection
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
virtual size_t nIds() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
TrackQuality
track quality
Definition: TrackBase.h:151
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
size_type size() const
Definition: OwnVector.h:264
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:518
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:230
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
FastTrackerRecHitMaskProducer(const edm::ParameterSet &)
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
reco::TrackCollection const & tracks(edm::Event &evt) const
const reco::TrackBase::TrackQuality trackQuality_
virtual void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: Handle.h:81
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
edm::EDGetTokenT< QualityMaskCollection > srcQuals
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:505
fixed size matrix
HLT enums.
const TrackCollectionTokens trajectories_
edm::EDGetTokenT< FastTrackerRecHitCollection > recHits_
virtual int32_t id(size_t i=0) const
bool isUninitialized() const
Definition: EDGetToken.h:73
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< std::vector< bool > > oldHitMaskToken_
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109