CMS 3D CMS Logo

HLTMhtProducer.cc
Go to the documentation of this file.
1 
10 
15 
18 
19 // Constructor
21  : usePt_(iConfig.getParameter<bool>("usePt")),
22  excludePFMuons_(iConfig.getParameter<bool>("excludePFMuons")),
23  minNJet_(iConfig.getParameter<int>("minNJet")),
24  minPtJet_(iConfig.getParameter<double>("minPtJet")),
25  maxEtaJet_(iConfig.getParameter<double>("maxEtaJet")),
26  jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
27  pfCandidatesLabel_(iConfig.getParameter<edm::InputTag>("pfCandidatesLabel")) {
28  m_theJetToken = consumes<reco::CandidateView>(jetsLabel_);
29  if (pfCandidatesLabel_.label().empty())
30  excludePFMuons_ = false;
31  if (excludePFMuons_)
32  m_thePFCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
33 
34  // Register the products
35  produces<reco::METCollection>();
36 }
37 
38 // Destructor
40 
41 // Fill descriptions
43  // Current default is for hltPFMET
45  desc.add<bool>("usePt", true);
46  desc.add<bool>("excludePFMuons", false);
47  desc.add<int>("minNJet", 0);
48  desc.add<double>("minPtJet", 0.);
49  desc.add<double>("maxEtaJet", 999.);
50  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4PFJets"));
51  desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
52  descriptions.add("hltMhtProducer", desc);
53 }
54 
55 // Produce the products
57  // Create a pointer to the products
58  std::unique_ptr<reco::METCollection> result(new reco::METCollection());
59 
61  iEvent.getByToken(m_theJetToken, jets);
62 
64  if (excludePFMuons_)
66 
67  int nj = 0;
68  double sumet = 0., mhx = 0., mhy = 0.;
69 
70  for (auto const& aJet : *jets) {
71  double const pt = usePt_ ? aJet.pt() : aJet.et();
72  double const eta = aJet.eta();
73  double const phi = aJet.phi();
74  double const px = usePt_ ? aJet.px() : aJet.et() * cos(phi);
75  double const py = usePt_ ? aJet.py() : aJet.et() * sin(phi);
76 
77  if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
78  mhx -= px;
79  mhy -= py;
80  sumet += pt;
81  ++nj;
82  }
83  }
84 
85  if (excludePFMuons_) {
86  for (auto const& aCand : *pfCandidates) {
87  if (std::abs(aCand.pdgId()) == 13) {
88  mhx += aCand.px();
89  mhy += aCand.py();
90  }
91  }
92  }
93 
94  if (nj < minNJet_) {
95  sumet = 0;
96  mhx = 0;
97  mhy = 0;
98  }
99 
100  reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx * mhx + mhy * mhy));
101  reco::MET::Point vtx(0, 0, 0);
102  reco::MET mht(sumet, p4, vtx);
103  result->push_back(mht);
104 
105  // Put the products into the Event
106  iEvent.put(std::move(result));
107 }
int minNJet_
Minimum number of jets passing pt and eta requirements.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::PFCandidateCollection > m_thePFCandidateToken
double maxEtaJet_
Maximum (abs) eta requirement for jets.
edm::EDGetTokenT< reco::CandidateView > m_theJetToken
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
edm::InputTag jetsLabel_
Input jet, PFCandidate collections.
std::string const & label() const
Definition: InputTag.h:36
edm::InputTag pfCandidatesLabel_
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:22
int iEvent
Definition: GenABIO.cc:224
Definition: MET.h:41
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double minPtJet_
Minimum pt requirement for jets.
bool usePt_
Use pt; otherwise, use et.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
HLT enums.
HLTMhtProducer(const edm::ParameterSet &iConfig)
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
def move(src, dest)
Definition: eostools.py:511
~HLTMhtProducer() override