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.
5 
7  jetsInput_ (iConfig.getParameter<edm::InputTag>("jetsInput")),
8  min_NHEF_ (iConfig.getParameter<double>("min_NHEF")),
9  max_NHEF_ (iConfig.getParameter<double>("max_NHEF")),
10  min_NEMF_ (iConfig.getParameter<double>("min_NEMF")),
11  max_NEMF_ (iConfig.getParameter<double>("max_NEMF")),
12  min_CEMF_ (iConfig.getParameter<double>("min_CEMF")),
13  max_CEMF_ (iConfig.getParameter<double>("max_CEMF")),
14  min_CHEF_ (iConfig.getParameter<double>("min_CHEF")),
15  max_CHEF_ (iConfig.getParameter<double>("max_CHEF")),
16  min_pt_ (iConfig.getParameter<double>("min_pt"))
17 {
18  m_thePFJetToken = consumes<reco::PFJetCollection>(jetsInput_);
19  produces< reco::PFJetCollection > ();
20 }
21 
24  desc.add<edm::InputTag>("jetsInput",edm::InputTag("hltAntiKT5PFJets"));
25  desc.add<double>("min_NHEF",-999.0);
26  desc.add<double>("max_NHEF",999.0);
27  desc.add<double>("min_NEMF",-999.0);
28  desc.add<double>("max_NEMF",999.0);
29  desc.add<double>("min_CEMF",-999.0);
30  desc.add<double>("max_CEMF",999.0);
31  desc.add<double>("min_CHEF",-999.0);
32  desc.add<double>("max_CHEF",999.0);
33  desc.add<double>("min_pt",30.0);
34  descriptions.add("hltPFJetIDProducer", desc);
35 }
36 
38 {
39 
40 }
41 
43 {
44 
45 }
46 
48 {
49 
51  iEvent.getByToken(m_thePFJetToken,pfjets);
52 
53  std::auto_ptr<reco::PFJetCollection> result (new reco::PFJetCollection);
54 
55  for (reco::PFJetCollection::const_iterator pfjetc = pfjets->begin();
56  pfjetc != pfjets->end(); ++pfjetc) {
57 
58  if (std::abs(pfjetc->eta()) >= 2.4) {
59  result->push_back(*pfjetc);
60  } else {
61  if (pfjetc->pt()<min_pt_) result->push_back(*pfjetc);
62  else if ((pfjetc->neutralHadronEnergyFraction() >= min_NHEF_) && (pfjetc->neutralHadronEnergyFraction() <= max_NHEF_) &&
63  (pfjetc->neutralEmEnergyFraction() >= min_NEMF_) && (pfjetc->neutralEmEnergyFraction() <= max_NEMF_) &&
64  (pfjetc->chargedEmEnergyFraction() >= min_CEMF_) && (pfjetc->chargedEmEnergyFraction() <= max_CEMF_) &&
65  (pfjetc->chargedHadronEnergyFraction() >= min_CHEF_) && (pfjetc->chargedHadronEnergyFraction() <= max_CHEF_)) {
66  result->push_back(*pfjetc);
67  }
68  }
69  } // pfjetc
70 
71  iEvent.put( result);
72 
73 }
edm::InputTag jetsInput_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void produce(edm::Event &, const edm::EventSetup &)
virtual void beginJob()
int iEvent
Definition: GenABIO.cc:243
edm::EDGetTokenT< reco::PFJetCollection > m_thePFJetToken
HLTPFJetIDProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)