CMS 3D CMS Logo

HLTCaloJetIDProducer.cc
Go to the documentation of this file.
1 
10 
15 
16 
17 // Constructor
19  min_N90_ (iConfig.getParameter<int>("min_N90")),
20  min_N90hits_(iConfig.getParameter<int>("min_N90hits")),
21  min_EMF_ (iConfig.getParameter<double>("min_EMF")),
22  max_EMF_ (iConfig.getParameter<double>("max_EMF")),
23  inputTag_ (iConfig.getParameter<edm::InputTag>("jetsInput")),
24  jetIDParams_(iConfig.getParameter<edm::ParameterSet>("JetIDParams")),
25  jetIDHelper_(jetIDParams_,consumesCollector()) {
26  m_theCaloJetToken = consumes<reco::CaloJetCollection>(inputTag_);
27 
28  // Register the products
29  produces<reco::CaloJetCollection>();
30 }
31 
32 // Destructor
34 
35 // Fill descriptions
38  desc.add<int>("min_N90", -2);
39  desc.add<int>("min_N90hits", 2);
40  desc.add<double>("min_EMF", 1e-6);
41  desc.add<double>("max_EMF", 999.);
42  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT4CaloJets"));
43 
45  descNested.add<bool>("useRecHits", true);
46  descNested.add<edm::InputTag>("hbheRecHitsColl", edm::InputTag("hltHbhereco"));
47  descNested.add<edm::InputTag>("hoRecHitsColl", edm::InputTag("hltHoreco"));
48  descNested.add<edm::InputTag>("hfRecHitsColl", edm::InputTag("hltHfreco"));
49  descNested.add<edm::InputTag>("ebRecHitsColl", edm::InputTag("hltEcalRecHit","EcalRecHitsEB"));
50  descNested.add<edm::InputTag>("eeRecHitsColl", edm::InputTag("hltEcalRecHit","EcalRecHitsEE"));
51  desc.add<edm::ParameterSetDescription>("JetIDParams", descNested);
52 
53  descriptions.add("hltCaloJetIDProducer", desc);
54 }
55 
56 // Produce the products
58 
59  // Create a pointer to the products
60  std::unique_ptr<reco::CaloJetCollection> result (new reco::CaloJetCollection());
61 
63  iEvent.getByToken(m_theCaloJetToken, calojets);
64 
65  for (auto const & j : *calojets) {
66  bool pass = false;
67 
68  if (!(j.energy() > 0.)) continue; // skip jets with zero or negative energy
69 
70  if (std::abs(j.eta()) >= 2.6) {
71  pass = true;
72 
73  } else {
74  if (min_N90hits_ > 0) jetIDHelper_.calculate(iEvent, iSetup, j);
75  if ((j.emEnergyFraction() >= min_EMF_) &&
76  (j.emEnergyFraction() <= max_EMF_) &&
77  (j.n90() >= min_N90_) &&
78  ((min_N90hits_ <= 0) || (jetIDHelper_.n90Hits() >= min_N90hits_)) ) {
79 
80  pass = true;
81  }
82  }
83 
84  if (pass) result->push_back(j);
85  }
86 
87  // Put the products into the Event
88  iEvent.put(std::move(result));
89 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
edm::InputTag inputTag_
input CaloJet collection
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
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:100
~HLTCaloJetIDProducer() override
int iEvent
Definition: GenABIO.cc:230
HLTCaloJetIDProducer(const edm::ParameterSet &iConfig)
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
HLT enums.
edm::EDGetTokenT< reco::CaloJetCollection > m_theCaloJetToken
def move(src, dest)
Definition: eostools.py:511
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects