CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTCaloJetIDProducer.cc
Go to the documentation of this file.
1 
10 
15 
16 // Constructor
18  : min_N90_(iConfig.getParameter<int>("min_N90")),
19  min_N90hits_(iConfig.getParameter<int>("min_N90hits")),
20  min_EMF_(iConfig.getParameter<double>("min_EMF")),
21  max_EMF_(iConfig.getParameter<double>("max_EMF")),
22  inputTag_(iConfig.getParameter<edm::InputTag>("jetsInput")),
23  jetIDParams_(iConfig.getParameter<edm::ParameterSet>("JetIDParams")),
24  jetIDHelper_(jetIDParams_, consumesCollector()) {
25  m_theCaloJetToken = consumes<reco::CaloJetCollection>(inputTag_);
26 
27  // Register the products
28  produces<reco::CaloJetCollection>();
29 }
30 
31 // Destructor
33 
34 // Fill descriptions
37  desc.add<int>("min_N90", -2);
38  desc.add<int>("min_N90hits", 2);
39  desc.add<double>("min_EMF", 1e-6);
40  desc.add<double>("max_EMF", 999.);
41  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT4CaloJets"));
42 
44  descNested.add<bool>("useRecHits", true);
45  descNested.add<edm::InputTag>("hbheRecHitsColl", edm::InputTag("hltHbhereco"));
46  descNested.add<edm::InputTag>("hoRecHitsColl", edm::InputTag("hltHoreco"));
47  descNested.add<edm::InputTag>("hfRecHitsColl", edm::InputTag("hltHfreco"));
48  descNested.add<edm::InputTag>("ebRecHitsColl", edm::InputTag("hltEcalRecHit", "EcalRecHitsEB"));
49  descNested.add<edm::InputTag>("eeRecHitsColl", edm::InputTag("hltEcalRecHit", "EcalRecHitsEE"));
50  desc.add<edm::ParameterSetDescription>("JetIDParams", descNested);
51 
52  descriptions.add("hltCaloJetIDProducer", desc);
53 }
54 
55 // Produce the products
57  // Create a pointer to the products
58  std::unique_ptr<reco::CaloJetCollection> result(new reco::CaloJetCollection());
59 
61  iEvent.getByToken(m_theCaloJetToken, calojets);
62 
63  for (auto const& j : *calojets) {
64  bool pass = false;
65 
66  if (!(j.energy() > 0.))
67  continue; // skip jets with zero or negative energy
68 
69  if (std::abs(j.eta()) >= 2.6) {
70  pass = true;
71 
72  } else {
73  if (min_N90hits_ > 0)
74  jetIDHelper_.calculate(iEvent, iSetup, j);
75  if ((j.emEnergyFraction() >= min_EMF_) && (j.emEnergyFraction() <= max_EMF_) && (j.n90() >= min_N90_) &&
76  ((min_N90hits_ <= 0) || (jetIDHelper_.n90Hits() >= min_N90hits_))) {
77  pass = true;
78  }
79  }
80 
81  if (pass)
82  result->push_back(j);
83  }
84 
85  // Put the products into the Event
86  iEvent.put(std::move(result));
87 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::InputTag inputTag_
input CaloJet collection
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
double min_EMF_
minimum EMF
double max_EMF_
maximum EMF
void calculate(const edm::Event &event, const edm::EventSetup &setup, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:93
~HLTCaloJetIDProducer() override
tuple result
Definition: mps_fire.py:311
int iEvent
Definition: GenABIO.cc:224
HLTCaloJetIDProducer(const edm::ParameterSet &iConfig)
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::helper::JetIDHelper jetIDHelper_
A helper to calculates calo jet ID variables.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int min_N90hits_
mininum N90hits
edm::EDGetTokenT< reco::CaloJetCollection > m_theCaloJetToken
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects