Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTTrackMETProducer
0002  *
0003  * See header file for documentation
0004  *
0005  *  \author Steven Lowette
0006  *
0007  */
0008 
0009 #include "HLTrigger/JetMET/interface/HLTTrackMETProducer.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 HLTTrackMETProducer::HLTTrackMETProducer(const edm::ParameterSet& iConfig)
0018     : usePt_(iConfig.getParameter<bool>("usePt")),
0019       useTracks_(iConfig.getParameter<bool>("useTracks")),
0020       usePFRecTracks_(iConfig.getParameter<bool>("usePFRecTracks")),
0021       usePFCandidatesCharged_(iConfig.getParameter<bool>("usePFCandidatesCharged")),
0022       usePFCandidates_(iConfig.getParameter<bool>("usePFCandidates")),
0023       excludePFMuons_(iConfig.getParameter<bool>("excludePFMuons")),
0024       minNJet_(iConfig.getParameter<int>("minNJet")),
0025       minPtJet_(iConfig.getParameter<double>("minPtJet")),
0026       maxEtaJet_(iConfig.getParameter<double>("maxEtaJet")),
0027       jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
0028       tracksLabel_(iConfig.getParameter<edm::InputTag>("tracksLabel")),
0029       pfRecTracksLabel_(iConfig.getParameter<edm::InputTag>("pfRecTracksLabel")),
0030       pfCandidatesLabel_(iConfig.getParameter<edm::InputTag>("pfCandidatesLabel")) {
0031   m_theJetToken = consumes<edm::View<reco::Jet>>(jetsLabel_);
0032   m_theTrackToken = consumes<reco::TrackCollection>(tracksLabel_);
0033   m_theRecTrackToken = consumes<reco::PFRecTrackCollection>(pfRecTracksLabel_);
0034   m_thePFCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
0035 
0036   // Register the products
0037   produces<reco::METCollection>();
0038 }
0039 
0040 // Destructor
0041 HLTTrackMETProducer::~HLTTrackMETProducer() = default;
0042 
0043 // Fill descriptions
0044 void HLTTrackMETProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0045   // Current default is for hltPFMET
0046   edm::ParameterSetDescription desc;
0047   desc.add<bool>("usePt", true);
0048   desc.add<bool>("useTracks", false);
0049   desc.add<bool>("usePFRecTracks", false);
0050   desc.add<bool>("usePFCandidatesCharged", true);
0051   desc.add<bool>("usePFCandidates", false);
0052   desc.add<bool>("excludePFMuons", false);
0053   desc.add<int>("minNJet", 0);
0054   desc.add<double>("minPtJet", 0.);
0055   desc.add<double>("maxEtaJet", 999.);
0056   desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4PFJets"));
0057   desc.add<edm::InputTag>("tracksLabel", edm::InputTag("hltL3Muons"));
0058   desc.add<edm::InputTag>("pfRecTracksLabel", edm::InputTag("hltLightPFTracks"));
0059   desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
0060   descriptions.add("hltTrackMETProducer", desc);
0061 }
0062 
0063 // Produce the products
0064 void HLTTrackMETProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0065   // Create a pointer to the products
0066   std::unique_ptr<reco::METCollection> result(new reco::METCollection());
0067 
0068   if (pfCandidatesLabel_.label().empty())
0069     excludePFMuons_ = false;
0070 
0071   bool useJets = !useTracks_ && !usePFRecTracks_ && !usePFCandidatesCharged_ && !usePFCandidates_;
0072   if (!useJets) {
0073     minNJet_ = 0;
0074   }
0075 
0076   edm::Handle<reco::JetView> jets;
0077   if (useJets)
0078     iEvent.getByToken(m_theJetToken, jets);
0079 
0080   edm::Handle<reco::TrackCollection> tracks;
0081   if (useTracks_)
0082     iEvent.getByToken(m_theTrackToken, tracks);
0083 
0084   edm::Handle<reco::PFRecTrackCollection> pfRecTracks;
0085   if (usePFRecTracks_)
0086     iEvent.getByToken(m_theRecTrackToken, pfRecTracks);
0087 
0088   edm::Handle<reco::PFCandidateCollection> pfCandidates;
0089   if (excludePFMuons_ || usePFCandidatesCharged_ || usePFCandidates_)
0090     iEvent.getByToken(m_thePFCandidateToken, pfCandidates);
0091 
0092   int nj = 0;
0093   double sumet = 0., mhx = 0., mhy = 0.;
0094 
0095   if (useJets && !jets->empty()) {
0096     for (reco::JetView::const_iterator j = jets->begin(); j != jets->end(); ++j) {
0097       double pt = usePt_ ? j->pt() : j->et();
0098       double eta = j->eta();
0099       double phi = j->phi();
0100       double px = usePt_ ? j->px() : j->et() * cos(phi);
0101       double py = usePt_ ? j->py() : j->et() * sin(phi);
0102 
0103       if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
0104         mhx -= px;
0105         mhy -= py;
0106         sumet += pt;
0107         ++nj;
0108       }
0109     }
0110 
0111   } else if (useTracks_ && !tracks->empty()) {
0112     for (auto const& j : *tracks) {
0113       double pt = j.pt();
0114       double px = j.px();
0115       double py = j.py();
0116       double eta = j.eta();
0117 
0118       if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
0119         mhx -= px;
0120         mhy -= py;
0121         sumet += pt;
0122         ++nj;
0123       }
0124     }
0125 
0126   } else if (usePFRecTracks_ && !pfRecTracks->empty()) {
0127     for (auto const& j : *pfRecTracks) {
0128       double pt = j.trackRef()->pt();
0129       double px = j.trackRef()->px();
0130       double py = j.trackRef()->py();
0131       double eta = j.trackRef()->eta();
0132 
0133       if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
0134         mhx -= px;
0135         mhy -= py;
0136         sumet += pt;
0137         ++nj;
0138       }
0139     }
0140 
0141   } else if ((usePFCandidatesCharged_ || usePFCandidates_) && !pfCandidates->empty()) {
0142     for (auto const& j : *pfCandidates) {
0143       if (usePFCandidatesCharged_ && j.charge() == 0)
0144         continue;
0145       double pt = j.pt();
0146       double px = j.px();
0147       double py = j.py();
0148       double eta = j.eta();
0149 
0150       if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
0151         mhx -= px;
0152         mhy -= py;
0153         sumet += pt;
0154         ++nj;
0155       }
0156     }
0157   }
0158 
0159   if (excludePFMuons_) {
0160     for (auto const& j : *pfCandidates) {
0161       if (std::abs(j.pdgId()) == 13) {
0162         mhx += j.px();
0163         mhy += j.py();
0164       }
0165     }
0166   }
0167 
0168   if (nj < minNJet_) {
0169     sumet = 0;
0170     mhx = 0;
0171     mhy = 0;
0172   }
0173 
0174   reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx * mhx + mhy * mhy));
0175   reco::MET::Point vtx(0, 0, 0);
0176   reco::MET mht(sumet, p4, vtx);
0177   result->push_back(mht);
0178 
0179   // Put the products into the Event
0180   iEvent.put(std::move(result));
0181 }