Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:32

0001 /** \class HLTCaloJetIDProducer
0002  *
0003  * See header file for documentation
0004  *
0005  *  \author a Jet/MET person
0006  *
0007  */
0008 
0009 #include "HLTrigger/JetMET/interface/HLTCaloJetIDProducer.h"
0010 
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "DataFormats/Common/interface/Handle.h"
0015 
0016 // Constructor
0017 HLTCaloJetIDProducer::HLTCaloJetIDProducer(const edm::ParameterSet& iConfig)
0018     : min_N90_(iConfig.getParameter<int>("min_N90")),
0019       min_N90hits_(iConfig.getParameter<int>("min_N90hits")),
0020       min_EMF_(iConfig.getParameter<double>("min_EMF")),
0021       max_EMF_(iConfig.getParameter<double>("max_EMF")),
0022       inputTag_(iConfig.getParameter<edm::InputTag>("jetsInput")),
0023       jetIDParams_(iConfig.getParameter<edm::ParameterSet>("JetIDParams")),
0024       jetIDHelper_(jetIDParams_, consumesCollector()) {
0025   m_theCaloJetToken = consumes<reco::CaloJetCollection>(inputTag_);
0026 
0027   // Register the products
0028   produces<reco::CaloJetCollection>();
0029 }
0030 
0031 // Destructor
0032 HLTCaloJetIDProducer::~HLTCaloJetIDProducer() = default;
0033 
0034 // Fill descriptions
0035 void HLTCaloJetIDProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0036   edm::ParameterSetDescription desc;
0037   desc.add<int>("min_N90", -2);
0038   desc.add<int>("min_N90hits", 2);
0039   desc.add<double>("min_EMF", 1e-6);
0040   desc.add<double>("max_EMF", 999.);
0041   desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT4CaloJets"));
0042 
0043   edm::ParameterSetDescription descNested;
0044   descNested.add<bool>("useRecHits", true);
0045   descNested.add<edm::InputTag>("hbheRecHitsColl", edm::InputTag("hltHbhereco"));
0046   descNested.add<edm::InputTag>("hoRecHitsColl", edm::InputTag("hltHoreco"));
0047   descNested.add<edm::InputTag>("hfRecHitsColl", edm::InputTag("hltHfreco"));
0048   descNested.add<edm::InputTag>("ebRecHitsColl", edm::InputTag("hltEcalRecHit", "EcalRecHitsEB"));
0049   descNested.add<edm::InputTag>("eeRecHitsColl", edm::InputTag("hltEcalRecHit", "EcalRecHitsEE"));
0050   desc.add<edm::ParameterSetDescription>("JetIDParams", descNested);
0051 
0052   descriptions.add("hltCaloJetIDProducer", desc);
0053 }
0054 
0055 // Produce the products
0056 void HLTCaloJetIDProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0057   // Create a pointer to the products
0058   std::unique_ptr<reco::CaloJetCollection> result(new reco::CaloJetCollection());
0059 
0060   edm::Handle<reco::CaloJetCollection> calojets;
0061   iEvent.getByToken(m_theCaloJetToken, calojets);
0062 
0063   for (auto const& j : *calojets) {
0064     bool pass = false;
0065 
0066     if (!(j.energy() > 0.))
0067       continue;  // skip jets with zero or negative energy
0068 
0069     if (std::abs(j.eta()) >= 2.6) {
0070       pass = true;
0071 
0072     } else {
0073       if (min_N90hits_ > 0)
0074         jetIDHelper_.calculate(iEvent, iSetup, j);
0075       if ((j.emEnergyFraction() >= min_EMF_) && (j.emEnergyFraction() <= max_EMF_) && (j.n90() >= min_N90_) &&
0076           ((min_N90hits_ <= 0) || (jetIDHelper_.n90Hits() >= min_N90hits_))) {
0077         pass = true;
0078       }
0079     }
0080 
0081     if (pass)
0082       result->push_back(j);
0083   }
0084 
0085   // Put the products into the Event
0086   iEvent.put(std::move(result));
0087 }