CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTMhtProducer.cc
Go to the documentation of this file.
1 
10 
15 
16 
17 // Constructor
19  usePt_ ( iConfig.getParameter<bool>("usePt") ),
20  excludePFMuons_ ( iConfig.getParameter<bool>("excludePFMuons") ),
21  minNJet_ ( iConfig.getParameter<int>("minNJet") ),
22  minPtJet_ ( iConfig.getParameter<double>("minPtJet") ),
23  maxEtaJet_ ( iConfig.getParameter<double>("maxEtaJet") ),
24  jetsLabel_ ( iConfig.getParameter<edm::InputTag>("jetsLabel") ),
25  pfCandidatesLabel_ ( iConfig.getParameter<edm::InputTag>("pfCandidatesLabel") ) {
26  m_theJetToken = consumes<edm::View<reco::Jet>>(jetsLabel_);
27  if (pfCandidatesLabel_.label() == "") excludePFMuons_ = false;
28  if (excludePFMuons_) m_thePFCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
29 
30  // Register the products
31  produces<reco::METCollection>();
32 }
33 
34 // Destructor
36 
37 // Fill descriptions
39  // Current default is for hltPFMET
41  desc.add<bool>("usePt", true);
42  desc.add<bool>("excludePFMuons", false);
43  desc.add<int>("minNJet",0);
44  desc.add<double>("minPtJet", 0.);
45  desc.add<double>("maxEtaJet", 999.);
46  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4PFJets"));
47  desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
48  descriptions.add("hltMhtProducer", desc);
49 }
50 
51 // Produce the products
53 
54  // Create a pointer to the products
55  std::auto_ptr<reco::METCollection> result(new reco::METCollection());
56 
58  iEvent.getByToken(m_theJetToken, jets);
59 
61  if (excludePFMuons_)
62  iEvent.getByToken(m_thePFCandidateToken, pfCandidates);
63 
64  int nj = 0;
65  double sumet = 0., mhx = 0., mhy = 0.;
66 
67  if (jets->size() > 0) {
68  for(reco::JetView::const_iterator j = jets->begin(); j != jets->end(); ++j) {
69  double pt = usePt_ ? j->pt() : j->et();
70  double eta = j->eta();
71  double phi = j->phi();
72  double px = usePt_ ? j->px() : j->et() * cos(phi);
73  double py = usePt_ ? j->py() : j->et() * sin(phi);
74 
75  if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
76  mhx -= px;
77  mhy -= py;
78  sumet += pt;
79  ++nj;
80  }
81  }
82  }
83 
84  if (excludePFMuons_) {
85  for (reco::PFCandidateCollection::const_iterator j = pfCandidates->begin(); j != pfCandidates->end(); ++j) {
86  if (std::abs(j->pdgId()) == 13) {
87  mhx += j->px();
88  mhy += j->py();
89  }
90  }
91  }
92 
93  if (nj < minNJet_) { sumet = 0; mhx = 0; mhy = 0; }
94 
95  reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx*mhx + mhy*mhy));
96  reco::MET::Point vtx(0, 0, 0);
97  reco::MET mht(sumet, p4, vtx);
98  result->push_back(mht);
99 
100  // Put the products into the Event
101  iEvent.put(result);
102 }
int minNJet_
Minimum number of jets passing pt and eta requirements.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::PFCandidateCollection > m_thePFCandidateToken
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
edm::EDGetTokenT< reco::JetView > m_theJetToken
double maxEtaJet_
Maximum (abs) eta requirement for jets.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
edm::InputTag jetsLabel_
Input jet, PFCandidate collections.
edm::InputTag pfCandidatesLabel_
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:23
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
Definition: MET.h:42
T sqrt(T t)
Definition: SSEVec.h:48
double p4[4]
Definition: TauolaWrapper.h:92
vector< PseudoJet > jets
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
double minPtJet_
Minimum pt requirement for jets.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool usePt_
Use pt; otherwise, use et.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
std::string const & label() const
Definition: InputTag.h:43
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
HLTMhtProducer(const edm::ParameterSet &iConfig)
math::XYZPoint Point
point in the space
Definition: Candidate.h:41