Go to the documentation of this file.00001
00002 #include "HLTrigger/JetMET/interface/HLTHtMhtProducer.h"
00003
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "FWCore/Framework/interface/EventSetup.h"
00007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00009 #include "FWCore/Utilities/interface/InputTag.h"
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "DataFormats/Common/interface/View.h"
00012
00013 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00014 #include "DataFormats/METReco/interface/METCollection.h"
00015 #include "DataFormats/METReco/interface/MET.h"
00016 #include "DataFormats/TrackReco/interface/Track.h"
00017 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00018
00019
00020 HLTHtMhtProducer::HLTHtMhtProducer(const edm::ParameterSet & iConfig) :
00021 usePt_ ( iConfig.getParameter<bool>("usePt") ),
00022 useTracks_ ( iConfig.getParameter<bool>("useTracks") ),
00023 excludePFMuons_ ( iConfig.getParameter<bool>("excludePFMuons") ),
00024 minNJetHt_ ( iConfig.getParameter<int>("minNJetHt") ),
00025 minNJetMht_ ( iConfig.getParameter<int>("minNJetMht") ),
00026 minPtJetHt_ ( iConfig.getParameter<double>("minPtJetHt") ),
00027 minPtJetMht_ ( iConfig.getParameter<double>("minPtJetMht") ),
00028 maxEtaJetHt_ ( iConfig.getParameter<double>("maxEtaJetHt") ),
00029 maxEtaJetMht_ ( iConfig.getParameter<double>("maxEtaJetMht") ),
00030 jetsLabel_ ( iConfig.getParameter<edm::InputTag>("jetsLabel") ),
00031 tracksLabel_ ( iConfig.getParameter<edm::InputTag>("tracksLabel") ),
00032 pfCandidatesLabel_ ( iConfig.getParameter<edm::InputTag>("pfCandidatesLabel") )
00033 {
00034 produces<reco::METCollection>();
00035 }
00036
00037
00038 HLTHtMhtProducer::~HLTHtMhtProducer() {
00039 }
00040
00041
00042 void HLTHtMhtProducer::fillDescriptions(edm::ConfigurationDescriptions & descriptions) {
00043 edm::ParameterSetDescription desc;
00044 desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltCaloJetCorrected"));
00045 desc.add<bool>("usePt", true);
00046 desc.add<int>("minNJetHt", 0);
00047 desc.add<int>("minNJetMht", 0);
00048 desc.add<double>("minPtJetHt", 40);
00049 desc.add<double>("minPtJetMht", 30);
00050 desc.add<double>("maxEtaJetHt", 3);
00051 desc.add<double>("maxEtaJetMht", 999);
00052 desc.add<bool>("useTracks", false);
00053 desc.add<edm::InputTag>("tracksLabel", edm::InputTag("hltL3Muons"));
00054 desc.add<bool>("excludePFMuons", false);
00055 desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
00056 descriptions.add("hltHtMhtProducer", desc);
00057 }
00058
00059
00060 void HLTHtMhtProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00061
00062 std::auto_ptr<reco::METCollection> metobject(new reco::METCollection());
00063
00064
00065 edm::Handle<edm::View<reco::Jet> > jets;
00066 iEvent.getByLabel(jetsLabel_, jets);
00067 edm::Handle<reco::TrackCollection> tracks;
00068 if (useTracks_) iEvent.getByLabel(tracksLabel_, tracks);
00069 edm::Handle<reco::PFCandidateCollection> pfCandidates;
00070 if (excludePFMuons_) iEvent.getByLabel(pfCandidatesLabel_, pfCandidates);
00071
00072 int nj_ht = 0, nj_mht = 0;
00073 double ht=0.;
00074 double mhtx=0., mhty=0.;
00075
00076
00077 for(edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); jet++ ) {
00078 double mom = (usePt_ ? jet->pt() : jet->et());
00079 if (mom > minPtJetHt_ and std::abs(jet->eta()) < maxEtaJetHt_) {
00080 ht += mom;
00081 ++nj_ht;
00082 }
00083 if (mom > minPtJetMht_ and std::abs(jet->eta()) < maxEtaJetMht_) {
00084 mhtx -= mom*cos(jet->phi());
00085 mhty -= mom*sin(jet->phi());
00086 ++nj_mht;
00087 }
00088 }
00089 if (useTracks_) {
00090 for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); track++) {
00091 if (track->pt() > minPtJetHt_ and std::abs(track->eta()) < maxEtaJetHt_) {
00092 ht += track->pt();
00093 }
00094 if (track->pt() > minPtJetMht_ and std::abs(track->eta()) < maxEtaJetMht_) {
00095 mhtx -= track->px();
00096 mhty -= track->py();
00097 }
00098 }
00099 }
00100 if (excludePFMuons_) {
00101 reco::PFCandidateCollection::const_iterator i (pfCandidates->begin());
00102 for (; i != pfCandidates->end(); ++i) {
00103 if(abs(i->pdgId())==13){
00104 mhtx += (i->pt())*cos(i->phi());
00105 mhty += (i->pt())*sin(i->phi());
00106 }
00107 }
00108 }
00109
00110 if (nj_ht < minNJetHt_ ) { ht = 0; }
00111 if (nj_mht < minNJetMht_) { mhtx = 0; mhty = 0; }
00112
00113 metobject->push_back(
00114 reco::MET(
00115 ht,
00116 reco::MET::LorentzVector(mhtx, mhty, 0, sqrt(mhtx*mhtx + mhty*mhty)),
00117 reco::MET::Point()
00118 )
00119 );
00120
00121 iEvent.put(metobject);
00122 }