CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTPFEnergyFractionsFilter.cc
Go to the documentation of this file.
1 
12 
14 
16 
18 
21 
25 
26 //
27 // constructors and destructor
28 //
30 {
31  inputPFJetTag_ = iConfig.getParameter< edm::InputTag > ("inputPFJetTag");
32  nJet_ = iConfig.getParameter<unsigned int> ("nJet");
33  min_CEEF_ = iConfig.getParameter<double> ("min_CEEF");
34  max_CEEF_ = iConfig.getParameter<double> ("max_CEEF");
35  min_NEEF_ = iConfig.getParameter<double> ("min_NEEF");
36  max_NEEF_ = iConfig.getParameter<double> ("max_NEEF");
37  min_CHEF_ = iConfig.getParameter<double> ("min_CHEF");
38  max_CHEF_ = iConfig.getParameter<double> ("max_CHEF");
39  min_NHEF_ = iConfig.getParameter<double> ("min_NHEF");
40  max_NHEF_ = iConfig.getParameter<double> ("max_NHEF");
41  triggerType_ = iConfig.getParameter<int> ("triggerType");
42  m_thePFJetToken = consumes<reco::PFJetCollection>(inputPFJetTag_);
43 }
44 
46 
47 void
51  desc.add<edm::InputTag>("inputPFJetTag",edm::InputTag("hltAntiKT5PFJets"));
52  desc.add<unsigned int>("nJet",1);
53  desc.add<double>("min_CEEF",-99.);
54  desc.add<double>("max_CEEF",99.);
55  desc.add<double>("min_NEEF",-99.);
56  desc.add<double>("max_NEEF",99.);
57  desc.add<double>("min_CHEF",-99.);
58  desc.add<double>("max_CHEF",99.);
59  desc.add<double>("min_NHEF",-99.);
60  desc.add<double>("max_NHEF",99.);
61  desc.add<int>("triggerType",trigger::TriggerJet);
62  descriptions.add("hltPFEnergyFractionsFilter",desc);
63 }
64 
65 // ------------ method called to produce the data ------------
66 bool
68 {
69  using namespace std;
70  using namespace edm;
71  using namespace reco;
72  using namespace trigger;
73 
74  // The filter object
75  if (saveTags()) filterproduct.addCollectionTag(inputPFJetTag_);
76 
77  // PFJets
79  iEvent.getByToken(m_thePFJetToken,recopfjets);
80 
81  //Checking
82  bool accept(false);
83 
84  if(recopfjets->size() >= nJet_){
85  accept = true;
86  unsigned int countJet(0);
87  //PF information
88  PFJetCollection::const_iterator i (recopfjets->begin());
89  for(; i != recopfjets->end(); ++i ){
90  if(countJet>=nJet_) break;
91  //
92  if(i->chargedEmEnergyFraction()<min_CEEF_) accept = false;
93  if(i->chargedEmEnergyFraction()>max_CEEF_) accept = false;
94  //
95  if(i->neutralEmEnergyFraction()<min_NEEF_) accept = false;
96  if(i->neutralEmEnergyFraction()>max_NEEF_) accept = false;
97  //
98  if(i->chargedHadronEnergyFraction()<min_CHEF_) accept = false;
99  if(i->chargedHadronEnergyFraction()>max_CHEF_) accept = false;
100  //
101  if(i->neutralHadronEnergyFraction()<min_NHEF_) accept = false;
102  if(i->neutralHadronEnergyFraction()>max_NHEF_) accept = false;
103  //
104  if(accept==false) break;
105  countJet++;
106  }
107 
108  //Store NJet_ jets
109  if(accept==true){
110  countJet = 0;
111  PFJetCollection::const_iterator i (recopfjets->begin());
112  for(; i != recopfjets->end(); ++i ){
113  if(countJet>=nJet_) break;
114  filterproduct.addObject(triggerType_,PFJetRef(recopfjets,distance(recopfjets->begin(),i)));
115  countJet++;
116  }
117  }
118  }// End of (recopfjets->size() >= nJet_)
119 
120  return accept;
121 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::PFJetCollection > m_thePFJetToken
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::Ref< PFJetCollection > PFJetRef
edm references
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
HLTPFEnergyFractionsFilter(const edm::ParameterSet &)
virtual 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&lt;C&gt;)
int iEvent
Definition: GenABIO.cc:230
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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