CMS 3D CMS Logo

HLTElectronEoverpFilterRegional.cc
Go to the documentation of this file.
1 
9 
11 
15 
20 
21 //
22 // constructors and destructor
23 //
25 {
26  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
27  electronIsolatedProducer_ = iConfig.getParameter< edm::InputTag > ("electronIsolatedProducer");
28  electronNonIsolatedProducer_ = iConfig.getParameter< edm::InputTag > ("electronNonIsolatedProducer");
29  eoverpbarrelcut_ = iConfig.getParameter<double> ("eoverpbarrelcut");
30  eoverpendcapcut_ = iConfig.getParameter<double> ("eoverpendcapcut");
31  ncandcut_ = iConfig.getParameter<int> ("ncandcut");
32  doIsolated_ = iConfig.getParameter<bool> ("doIsolated");
33  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
34  electronIsolatedToken_ = consumes<reco::ElectronCollection>(electronIsolatedProducer_);
35  if(!doIsolated_) electronNonIsolatedToken_ = consumes<reco::ElectronCollection>(electronNonIsolatedProducer_);
36 }
37 
39 
40 void
44  desc.add<edm::InputTag>("candTag",edm::InputTag("hltElectronPixelMatchFilter"));
45  desc.add<edm::InputTag>("electronIsolatedProducer",edm::InputTag("pixelMatchElectronsForHLT"));
46  desc.add<edm::InputTag>("electronNonIsolatedProducer",edm::InputTag("pixelMatchElectronsForHLT"));
47  desc.add<double>("eoverpbarrelcut",1.5);
48  desc.add<double>("eoverpendcapcut",2.45);
49  desc.add<int>("ncandcut",1);
50  desc.add<bool>("doIsolated",true);
51  descriptions.add("hltElectronEoverpFilter",desc);
52 }
53 
54 // ------------ method called to produce the data ------------
55 bool
57 {
58 
59  // The filter object
60  using namespace trigger;
61  if (saveTags()) {
64  }
65 
67  iEvent.getByToken (candToken_,PrevFilterOutput);
68 
69  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
70  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
71  if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton, recoecalcands);
72  // Get the HLT electrons from EgammaHLTPixelMatchElectronProducers
73  edm::Handle<reco::ElectronCollection> electronIsolatedHandle;
74  iEvent.getByToken(electronIsolatedToken_,electronIsolatedHandle);
75 
76  edm::Handle<reco::ElectronCollection> electronNonIsolatedHandle;
77  if(!doIsolated_) {
78  iEvent.getByToken(electronNonIsolatedToken_,electronNonIsolatedHandle);
79  }
80 
81  // look at all candidates, check cuts and add to filter object
82  int n(0);
83 
84  //loop over all the RecoCandidates from the previous filter,
85  // associate them with the corresponding Electron object
86  //(the matching is done checking the super clusters)
87  // and put into the event a Ref to the Electron objects that passes the
88  // selections
89  for (auto & recoecalcand : recoecalcands) {
90  reco::SuperClusterRef recr2 = recoecalcand->superCluster();
91 
92  //loop over the electrons to find the matching one
93  for(auto iElectron = electronIsolatedHandle->begin(); iElectron != electronIsolatedHandle->end(); iElectron++){
94 
95  reco::ElectronRef electronref(reco::ElectronRef(electronIsolatedHandle,iElectron - electronIsolatedHandle->begin()));
96  const reco::SuperClusterRef theClus = electronref->superCluster();
97 
98  if(&(*recr2) == &(*theClus)) {
99 
100  float elecEoverp = 0;
101  const math::XYZVector trackMom = electronref->track()->momentum();
102  if( trackMom.R() != 0) elecEoverp =
103  electronref->superCluster()->energy()/ trackMom.R();
104 
105  if( fabs(electronref->eta()) < 1.5 ){
106  if ( elecEoverp < eoverpbarrelcut_) {
107  n++;
108  filterproduct.addObject(TriggerElectron, electronref);
109  }
110  }
111  if( fabs(electronref->eta()) > 1.5 ){
112  if ( elecEoverp < eoverpendcapcut_) {
113  n++;
114  filterproduct.addObject(TriggerElectron, electronref);
115  }
116  }
117  }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
118  }//end of loop over electrons
119 
120  if(!doIsolated_) {
121  //loop over the electrons to find the matching one
122  for(auto iElectron = electronNonIsolatedHandle->begin(); iElectron != electronNonIsolatedHandle->end(); iElectron++){
123 
124  reco::ElectronRef electronref(reco::ElectronRef(electronNonIsolatedHandle,iElectron - electronNonIsolatedHandle->begin()));
125  const reco::SuperClusterRef theClus = electronref->superCluster();
126 
127  if(&(*recr2) == &(*theClus)) {
128 
129  float elecEoverp = 0;
130  const math::XYZVector trackMom = electronref->track()->momentum();
131  if( trackMom.R() != 0) elecEoverp =
132  electronref->superCluster()->energy()/ trackMom.R();
133 
134  if( fabs(electronref->eta()) < 1.5 ){
135  if ( elecEoverp < eoverpbarrelcut_) {
136  n++;
137  filterproduct.addObject(TriggerElectron, electronref);
138  }
139  }
140  if( fabs(electronref->eta()) > 1.5 ){
141  if ( elecEoverp < eoverpendcapcut_) {
142  n++;
143  filterproduct.addObject(TriggerElectron, electronref);
144  }
145  }
146  }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
147  }//end of loop over electrons
148  }
149  }//end of loop ober candidates
150 
151  // filter decision
152  bool accept(n>=ncandcut_);
153 
154  return accept;
155 }
156 
157 // declare this class as a framework plugin
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
HLTElectronEoverpFilterRegional(const edm::ParameterSet &)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
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
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45
edm::EDGetTokenT< reco::ElectronCollection > electronNonIsolatedToken_
edm::EDGetTokenT< reco::ElectronCollection > electronIsolatedToken_