CMS 3D CMS Logo

HLTElectronMissingHitsFilter.cc
Go to the documentation of this file.
1 
14 
19 
21  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
22  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
23  electronTag_ = iConfig.getParameter< edm::InputTag >("electronProducer");
24  electronToken_ = consumes<reco::ElectronCollection>(electronTag_);
25  barrelcut_ = iConfig.getParameter<int> ("barrelcut");
26  endcapcut_ = iConfig.getParameter<int> ("endcapcut");
27  ncandcut_ = iConfig.getParameter<int> ("ncandcut");
28 }
29 
31 
35  desc.add<edm::InputTag>("candTag", edm::InputTag("hltSingleElectronOneOEMinusOneOPFilter"));
36  desc.add<edm::InputTag>("electronProducer", edm::InputTag("hltL1NonIsoHLTNonIsoSingleElectronEt15LTIPixelMatchFilte"));
37  desc.add<int>("barrelcut", 0);
38  desc.add<int>("endcapcut", 0);
39  desc.add<int>("ncandcut", 1);
40  descriptions.add("hltElectronMissingHitsFilter", desc);
41 }
42 
44 
45  using namespace trigger;
46 
47  if (saveTags())
48  filterproduct.addCollectionTag(electronTag_);
49 
51  iEvent.getByToken(candToken_,PrevFilterOutput);
52 
53  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
54  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
55  if(recoecalcands.empty())
56  PrevFilterOutput->getObjects(TriggerPhoton,recoecalcands);
57 
59  iEvent.getByToken(electronToken_,electronHandle);
60 
61  int n(0);
62 
64  for (auto & recoecalcand : recoecalcands) {
65 
66  reco::SuperClusterRef scCand = recoecalcand->superCluster();
67  for(auto iElectron = electronHandle->begin(); iElectron != electronHandle->end(); iElectron++) {
68  reco::ElectronRef electronref(reco::ElectronRef(electronHandle, iElectron - electronHandle->begin()));
69  const reco::SuperClusterRef scEle = electronref->superCluster();
70  if(scCand == scEle) {
71 
72  int missinghits = 0;
73  if (electronref->gsfTrack().isNonnull()){
74  missinghits = electronref->gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
75  } else if (electronref->track().isNonnull()){
76  missinghits = electronref->track()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
77  }else{
78  std::cerr << "Electron without track..." << std::endl;
79  }
80 
81  if(fabs(electronref->eta()) < 1.479) {
82  if (missinghits < barrelcut_) {
83  n++;
84  filterproduct.addObject(TriggerElectron, electronref);
85  }
86  }
87 
88  if(fabs(electronref->eta()) > 1.479) {
89  if (missinghits < endcapcut_) {
90  n++;
91  filterproduct.addObject(TriggerElectron, electronref);
92  }
93  }
94  }
95  }
96  }
97 
98  bool accept(n >= ncandcut_);
99 
100  return accept;
101 }
102 
103 // declare this class as a framework plugin
T getParameter(std::string const &) const
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
HLTElectronMissingHitsFilter(const edm::ParameterSet &)
edm::EDGetTokenT< reco::ElectronCollection > electronToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
~HLTElectronMissingHitsFilter() override
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)