CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTPFJetIDProducer.cc
Go to the documentation of this file.
1 
10 
15 
16 
17 // Constructor
19  minPt_ (iConfig.getParameter<double>("minPt")),
20  maxEta_ (iConfig.getParameter<double>("maxEta")),
21  CHF_ (iConfig.getParameter<double>("CHF")),
22  NHF_ (iConfig.getParameter<double>("NHF")),
23  CEF_ (iConfig.getParameter<double>("CEF")),
24  NEF_ (iConfig.getParameter<double>("NEF")),
25  NCH_ (iConfig.getParameter<int>("NCH")),
26  NTOT_ (iConfig.getParameter<int>("NTOT")),
27  inputTag_ (iConfig.getParameter<edm::InputTag>("jetsInput")) {
28  m_thePFJetToken = consumes<reco::PFJetCollection>(inputTag_);
29 
30  // Register the products
31  produces<reco::PFJetCollection>();
32 }
33 
34 // Destructor
36 
37 // Fill descriptions
40  desc.add<double>("minPt", 20.);
41  desc.add<double>("maxEta", 1e99);
42  desc.add<double>("CHF", -99.);
43  desc.add<double>("NHF", 99.);
44  desc.add<double>("CEF", 99.);
45  desc.add<double>("NEF", 99.);
46  desc.add<int>("NCH", 0);
47  desc.add<int>("NTOT", 0);
48  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT4PFJets"));
49  descriptions.add("hltPFJetIDProducer", desc);
50 }
51 
52 // Produce the products
54 
55  // Create a pointer to the products
56  std::auto_ptr<reco::PFJetCollection> result (new reco::PFJetCollection());
57 
59  iEvent.getByToken(m_thePFJetToken, pfjets);
60 
61  for (reco::PFJetCollection::const_iterator j = pfjets->begin(); j != pfjets->end(); ++j) {
62  bool pass = false;
63  double pt = j->pt();
64  double eta = j->eta();
65 
66  if (!(pt > 0.)) continue; // skip jets with zero or negative pt
67 
68  if (pt < minPt_) {
69  pass = true;
70 
71  } else if (std::abs(eta) >= maxEta_) {
72  pass = true;
73 
74  } else {
75  double chf = j->chargedHadronEnergyFraction();
76  //double nhf = j->neutralHadronEnergyFraction() + j->HFHadronEnergyFraction();
77  double nhf = j->neutralHadronEnergyFraction();
78  double cef = j->chargedEmEnergyFraction();
79  double nef = j->neutralEmEnergyFraction();
80  int nch = j->chargedMultiplicity();
81  int ntot = j->numberOfDaughters();
82 
83  pass = true;
84  pass = pass && (ntot > NTOT_);
85  pass = pass && (nef < NEF_);
86  pass = pass && (nhf < NHF_); //to be revisited if PF-reconstruction at HLT changes to be in line with offline code
87  pass = pass && (cef < CEF_ || std::abs(eta) >= 2.4);
88  pass = pass && (chf > CHF_ || std::abs(eta) >= 2.4);
89  pass = pass && (nch > NCH_ || std::abs(eta) >= 2.4);
90  }
91 
92  if (pass) result->push_back(*j);
93  }
94 
95  // Put the products into the Event
96  iEvent.put(result);
97 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
HLTPFJetIDProducer(const edm::ParameterSet &iConfig)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::PFJetCollection > m_thePFJetToken
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
tuple result
Definition: query.py:137
double NHF_
neutral hadron fraction
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double NEF_
neutral EM fraction
int NTOT_
number of constituents
double CHF_
charged hadron fraction
int NCH_
number of charged constituents
edm::InputTag inputTag_
input PFJet collection
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)