Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTHtMhtProducer
0002  *
0003  * See header file for documentation
0004  *
0005  *  \author Steven Lowette
0006  *
0007  */
0008 
0009 #include "HLTrigger/JetMET/interface/HLTHtMhtProducer.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 #include "DataFormats/METReco/interface/MET.h"
0017 #include "DataFormats/METReco/interface/METFwd.h"
0018 
0019 // Constructor
0020 HLTHtMhtProducer::HLTHtMhtProducer(const edm::ParameterSet& iConfig)
0021     : usePt_(iConfig.getParameter<bool>("usePt")),
0022       excludePFMuons_(iConfig.getParameter<bool>("excludePFMuons")),
0023       minNJetHt_(iConfig.getParameter<int>("minNJetHt")),
0024       minNJetMht_(iConfig.getParameter<int>("minNJetMht")),
0025       minPtJetHt_(iConfig.getParameter<double>("minPtJetHt")),
0026       minPtJetMht_(iConfig.getParameter<double>("minPtJetMht")),
0027       maxEtaJetHt_(iConfig.getParameter<double>("maxEtaJetHt")),
0028       maxEtaJetMht_(iConfig.getParameter<double>("maxEtaJetMht")),
0029       jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
0030       pfCandidatesLabel_(iConfig.getParameter<edm::InputTag>("pfCandidatesLabel")) {
0031   m_theJetToken = consumes<reco::CandidateView>(jetsLabel_);
0032   m_thePFCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
0033 
0034   // Register the products
0035   produces<reco::METCollection>();
0036 }
0037 
0038 // Destructor
0039 HLTHtMhtProducer::~HLTHtMhtProducer() = default;
0040 
0041 // Fill descriptions
0042 void HLTHtMhtProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0043   // Current default is for hltHtMht
0044   edm::ParameterSetDescription desc;
0045   desc.add<bool>("usePt", false);
0046   desc.add<bool>("excludePFMuons", false);
0047   desc.add<int>("minNJetHt", 0);
0048   desc.add<int>("minNJetMht", 0);
0049   desc.add<double>("minPtJetHt", 40.);
0050   desc.add<double>("minPtJetMht", 30.);
0051   desc.add<double>("maxEtaJetHt", 3.);
0052   desc.add<double>("maxEtaJetMht", 5.);
0053   desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltCaloJetL1FastJetCorrected"));
0054   desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
0055   descriptions.add("hltHtMhtProducer", desc);
0056 }
0057 
0058 // Produce the products
0059 void HLTHtMhtProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0060   // Create a pointer to the products
0061   std::unique_ptr<reco::METCollection> result(new reco::METCollection());
0062 
0063   if (pfCandidatesLabel_.label().empty())
0064     excludePFMuons_ = false;
0065 
0066   edm::Handle<reco::CandidateView> jets;
0067   iEvent.getByToken(m_theJetToken, jets);
0068 
0069   edm::Handle<reco::PFCandidateCollection> pfCandidates;
0070   if (excludePFMuons_)
0071     iEvent.getByToken(m_thePFCandidateToken, pfCandidates);
0072 
0073   int nj_ht = 0, nj_mht = 0;
0074   double ht = 0., mhx = 0., mhy = 0.;
0075 
0076   for (auto const& aJet : *jets) {
0077     double const pt = usePt_ ? aJet.pt() : aJet.et();
0078     double const eta = aJet.eta();
0079     double const phi = aJet.phi();
0080     double const px = usePt_ ? aJet.px() : aJet.et() * cos(phi);
0081     double const py = usePt_ ? aJet.py() : aJet.et() * sin(phi);
0082 
0083     if (pt > minPtJetHt_ && std::abs(eta) < maxEtaJetHt_) {
0084       ht += pt;
0085       ++nj_ht;
0086     }
0087 
0088     if (pt > minPtJetMht_ && std::abs(eta) < maxEtaJetMht_) {
0089       mhx -= px;
0090       mhy -= py;
0091       ++nj_mht;
0092     }
0093   }
0094 
0095   if (excludePFMuons_) {
0096     for (auto const& aCand : *pfCandidates) {
0097       if (std::abs(aCand.pdgId()) == 13) {
0098         mhx += aCand.px();
0099         mhy += aCand.py();
0100       }
0101     }
0102   }
0103 
0104   if (nj_ht < minNJetHt_) {
0105     ht = 0;
0106   }
0107   if (nj_mht < minNJetMht_) {
0108     mhx = 0;
0109     mhy = 0;
0110   }
0111 
0112   reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx * mhx + mhy * mhy));
0113   reco::MET::Point vtx(0, 0, 0);
0114   reco::MET htmht(ht, p4, vtx);
0115   result->push_back(htmht);
0116 
0117   // Put the products into the Event
0118   iEvent.put(std::move(result));
0119 }