CMS 3D CMS Logo

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  maxCF_ (iConfig.getParameter<double>("maxCF")),
26  NCH_ (iConfig.getParameter<int>("NCH")),
27  NTOT_ (iConfig.getParameter<int>("NTOT")),
28  inputTag_ (iConfig.getParameter<edm::InputTag>("jetsInput")) {
29  m_thePFJetToken = consumes<reco::PFJetCollection>(inputTag_);
30 
31  // Register the products
32  produces<reco::PFJetCollection>();
33 }
34 
35 // Destructor
37 
38 // Fill descriptions
41  desc.add<double>("minPt", 20.);
42  desc.add<double>("maxEta", 1e99);
43  desc.add<double>("CHF", -99.);
44  desc.add<double>("NHF", 99.);
45  desc.add<double>("CEF", 99.);
46  desc.add<double>("NEF", 99.);
47  desc.add<double>("maxCF", 99.);
48  desc.add<int>("NCH", 0);
49  desc.add<int>("NTOT", 0);
50  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT4PFJets"));
51  descriptions.add("hltPFJetIDProducer", desc);
52 }
53 
54 // Produce the products
56 
57  // Create a pointer to the products
58  std::unique_ptr<reco::PFJetCollection> result (new reco::PFJetCollection());
59 
61  iEvent.getByToken(m_thePFJetToken, pfjets);
62 
63  for (auto const & j : *pfjets) {
64  bool pass = false;
65  double pt = j.pt();
66  double eta = j.eta();
67  double abseta = std::abs(eta);
68 
69  if (!(pt > 0.)) continue; // skip jets with zero or negative pt
70 
71  if (pt < minPt_) {
72  pass = true;
73 
74  } else if (abseta >= maxEta_) {
75  pass = true;
76 
77  } else {
78  double chf = j.chargedHadronEnergyFraction();
79  //double nhf = j->neutralHadronEnergyFraction() + j->HFHadronEnergyFraction();
80  double nhf = j.neutralHadronEnergyFraction();
81  double cef = j.chargedEmEnergyFraction();
82  double nef = j.neutralEmEnergyFraction();
83  double cftot= chf + cef + j.chargedMuEnergyFraction();
84  int nch = j.chargedMultiplicity();
85  int ntot = j.numberOfDaughters();
86 
87  pass = true;
88  pass = pass && (ntot > NTOT_);
89  pass = pass && (nef < NEF_);
90  pass = pass && (nhf < NHF_ || abseta >= 2.4); //NHF-cut does not work in HF anymore with recent PF
91  pass = pass && (cef < CEF_ || abseta >= 2.4);
92  pass = pass && (chf > CHF_ || abseta >= 2.4);
93  pass = pass && (nch > NCH_ || abseta >= 2.4);
94  pass = pass && (cftot < maxCF_ || abseta >= 2.4);
95  }
96 
97  if (pass) result->push_back(j);
98  }
99 
100  // Put the products into the Event
101  iEvent.put(std::move(result));
102 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HLTPFJetIDProducer(const edm::ParameterSet &iConfig)
~HLTPFJetIDProducer() override
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::PFJetCollection > m_thePFJetToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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)
HLT enums.
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
def move(src, dest)
Definition: eostools.py:511