CMS 3D CMS Logo

HLTEgammaEtFilter.cc
Go to the documentation of this file.
1 
8 #include "HLTEgammaEtFilter.h"
9 
11 
15 
17 
18 //
19 // constructors and destructor
20 //
22 {
23  inputTag_ = iConfig.getParameter< edm::InputTag > ("inputTag");
24  etcutEB_ = iConfig.getParameter<double> ("etcutEB");
25  etcutEE_ = iConfig.getParameter<double> ("etcutEE");
26  ncandcut_ = iConfig.getParameter<int> ("ncandcut");
27  l1EGTag_ = iConfig.getParameter< edm::InputTag > ("l1EGCand");
28  inputToken_ = consumes<trigger::TriggerFilterObjectWithRefs> (inputTag_);
29 }
30 
31 void
35  desc.add<edm::InputTag>("inputTag", edm::InputTag("HLTEgammaL1MatchFilter"));
36  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
37  desc.add<double>("etcutEB", 1.0);
38  desc.add<double>("etcutEE", 1.0);
39  desc.add<int>("ncandcut", 1);
40  descriptions.add("hltEgammaEtFilter", desc);
41 }
42 
44 
45 
46 // ------------ method called to produce the data ------------
47 bool
49 {
50  using namespace trigger;
51 
52  // The filter object
53  if (saveTags()) {
54  filterproduct.addCollectionTag(l1EGTag_);;
55  }
56 
57  // Ref to Candidate object to be recorded in filter object
59 
60  // get hold of filtered candidates
61  //edm::Handle<reco::HLTFilterObjectWithRefs> recoecalcands;
63  iEvent.getByToken (inputToken_, PrevFilterOutput);
64 
65  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands; // vref with your specific C++ collection type
66  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
67  if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton,recoecalcands); //we dont know if its type trigger cluster or trigger photon
68 
69  // look at all candidates, check cuts and add to filter object
70  int n(0);
71 
72  for (auto & recoecalcand : recoecalcands) {
73 
74  ref = recoecalcand ;
75 
76  if( ( fabs(ref->eta()) < 1.479 && ref->et() >= etcutEB_ ) || ( fabs(ref->eta()) >= 1.479 && ref->et() >= etcutEE_ ) ){
77  n++;
78  // std::cout << "Passed eta: " << ref->eta() << std::endl;
79  filterproduct.addObject(TriggerCluster, ref);
80  }
81  }
82 
83 
84  // filter decision
85  bool accept(n>=ncandcut_);
86 
87  return accept;
88 }
89 
90 // declare this class as a framework plugin
T getParameter(std::string const &) const
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken_
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::InputTag l1EGTag_
edm::InputTag inputTag_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const 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
ParameterDescriptionBase * add(U const &iLabel, T const &value)
~HLTEgammaEtFilter() override
HLTEgammaEtFilter(const edm::ParameterSet &)
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)