CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTEgammaEtFilter.cc
Go to the documentation of this file.
1 
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  relaxed_ = iConfig.getUntrackedParameter<bool> ("relaxed",true) ;
28  L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand");
29  L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand");
30  inputToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(inputTag_);
31 }
32 
33 void
37  desc.add<edm::InputTag>("inputTag",edm::InputTag("HLTEgammaL1MatchFilter"));
38  desc.add<edm::InputTag>("L1IsoCand",edm::InputTag("hltL1IsoRecoEcalCandidate"));
39  desc.add<edm::InputTag>("L1NonIsoCand",edm::InputTag("hltL1NonIsoRecoEcalCandidate"));
40  desc.addUntracked<bool>("relaxed",true);
41  desc.add<double>("etcutEB", 1.0);
42  desc.add<double>("etcutEE", 1.0);
43  desc.add<int>("ncandcut", 1);
44  descriptions.add("hltEgammaEtFilter",desc);
45 }
46 
48 
49 
50 // ------------ method called to produce the data ------------
51 bool
53 {
54  using namespace trigger;
55 
56  // The filter object
57  if (saveTags()) {
58  filterproduct.addCollectionTag(L1IsoCollTag_);
59  if (relaxed_) filterproduct.addCollectionTag(L1NonIsoCollTag_);
60  }
61 
62  // Ref to Candidate object to be recorded in filter object
64 
65  // get hold of filtered candidates
66  //edm::Handle<reco::HLTFilterObjectWithRefs> recoecalcands;
68  iEvent.getByToken (inputToken_,PrevFilterOutput);
69 
70  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands; // vref with your specific C++ collection type
71  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
72  if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton,recoecalcands); //we dont know if its type trigger cluster or trigger photon
73 
74  // look at all candidates, check cuts and add to filter object
75  int n(0);
76 
77  for (unsigned int i=0; i<recoecalcands.size(); i++) {
78 
79  ref = recoecalcands[i] ;
80 
81  if( ( fabs(ref->eta()) < 1.479 && ref->et() >= etcutEB_ ) || ( fabs(ref->eta()) >= 1.479 && ref->et() >= etcutEE_ ) ){
82  n++;
83  // std::cout << "Passed eta: " << ref->eta() << std::endl;
84  filterproduct.addObject(TriggerCluster, ref);
85  }
86  }
87 
88 
89  // filter decision
90  bool accept(n>=ncandcut_);
91 
92  return accept;
93 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken_
int i
Definition: DBlmapReader.cc:9
edm::InputTag L1NonIsoCollTag_
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::InputTag inputTag_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
int iEvent
Definition: GenABIO.cc:230
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::InputTag L1IsoCollTag_
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)