CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTTrackMETProducer.cc
Go to the documentation of this file.
1 
10 
15 
16 // Constructor
18  : usePt_(iConfig.getParameter<bool>("usePt")),
19  useTracks_(iConfig.getParameter<bool>("useTracks")),
20  usePFRecTracks_(iConfig.getParameter<bool>("usePFRecTracks")),
21  usePFCandidatesCharged_(iConfig.getParameter<bool>("usePFCandidatesCharged")),
22  usePFCandidates_(iConfig.getParameter<bool>("usePFCandidates")),
23  excludePFMuons_(iConfig.getParameter<bool>("excludePFMuons")),
24  minNJet_(iConfig.getParameter<int>("minNJet")),
25  minPtJet_(iConfig.getParameter<double>("minPtJet")),
26  maxEtaJet_(iConfig.getParameter<double>("maxEtaJet")),
27  jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
28  tracksLabel_(iConfig.getParameter<edm::InputTag>("tracksLabel")),
29  pfRecTracksLabel_(iConfig.getParameter<edm::InputTag>("pfRecTracksLabel")),
30  pfCandidatesLabel_(iConfig.getParameter<edm::InputTag>("pfCandidatesLabel")) {
31  m_theJetToken = consumes<edm::View<reco::Jet>>(jetsLabel_);
32  m_theTrackToken = consumes<reco::TrackCollection>(tracksLabel_);
33  m_theRecTrackToken = consumes<reco::PFRecTrackCollection>(pfRecTracksLabel_);
34  m_thePFCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
35 
36  // Register the products
37  produces<reco::METCollection>();
38 }
39 
40 // Destructor
42 
43 // Fill descriptions
45  // Current default is for hltPFMET
47  desc.add<bool>("usePt", true);
48  desc.add<bool>("useTracks", false);
49  desc.add<bool>("usePFRecTracks", false);
50  desc.add<bool>("usePFCandidatesCharged", true);
51  desc.add<bool>("usePFCandidates", false);
52  desc.add<bool>("excludePFMuons", false);
53  desc.add<int>("minNJet", 0);
54  desc.add<double>("minPtJet", 0.);
55  desc.add<double>("maxEtaJet", 999.);
56  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4PFJets"));
57  desc.add<edm::InputTag>("tracksLabel", edm::InputTag("hltL3Muons"));
58  desc.add<edm::InputTag>("pfRecTracksLabel", edm::InputTag("hltLightPFTracks"));
59  desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
60  descriptions.add("hltTrackMETProducer", desc);
61 }
62 
63 // Produce the products
65  // Create a pointer to the products
66  std::unique_ptr<reco::METCollection> result(new reco::METCollection());
67 
68  if (pfCandidatesLabel_.label().empty())
69  excludePFMuons_ = false;
70 
72  if (!useJets) {
73  minNJet_ = 0;
74  }
75 
77  if (useJets)
78  iEvent.getByToken(m_theJetToken, jets);
79 
81  if (useTracks_)
82  iEvent.getByToken(m_theTrackToken, tracks);
83 
85  if (usePFRecTracks_)
86  iEvent.getByToken(m_theRecTrackToken, pfRecTracks);
87 
90  iEvent.getByToken(m_thePFCandidateToken, pfCandidates);
91 
92  int nj = 0;
93  double sumet = 0., mhx = 0., mhy = 0.;
94 
95  if (useJets && !jets->empty()) {
96  for (reco::JetView::const_iterator j = jets->begin(); j != jets->end(); ++j) {
97  double pt = usePt_ ? j->pt() : j->et();
98  double eta = j->eta();
99  double phi = j->phi();
100  double px = usePt_ ? j->px() : j->et() * cos(phi);
101  double py = usePt_ ? j->py() : j->et() * sin(phi);
102 
103  if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
104  mhx -= px;
105  mhy -= py;
106  sumet += pt;
107  ++nj;
108  }
109  }
110 
111  } else if (useTracks_ && !tracks->empty()) {
112  for (auto const& j : *tracks) {
113  double pt = j.pt();
114  double px = j.px();
115  double py = j.py();
116  double eta = j.eta();
117 
118  if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
119  mhx -= px;
120  mhy -= py;
121  sumet += pt;
122  ++nj;
123  }
124  }
125 
126  } else if (usePFRecTracks_ && !pfRecTracks->empty()) {
127  for (auto const& j : *pfRecTracks) {
128  double pt = j.trackRef()->pt();
129  double px = j.trackRef()->px();
130  double py = j.trackRef()->py();
131  double eta = j.trackRef()->eta();
132 
133  if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
134  mhx -= px;
135  mhy -= py;
136  sumet += pt;
137  ++nj;
138  }
139  }
140 
141  } else if ((usePFCandidatesCharged_ || usePFCandidates_) && !pfCandidates->empty()) {
142  for (auto const& j : *pfCandidates) {
143  if (usePFCandidatesCharged_ && j.charge() == 0)
144  continue;
145  double pt = j.pt();
146  double px = j.px();
147  double py = j.py();
148  double eta = j.eta();
149 
150  if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
151  mhx -= px;
152  mhy -= py;
153  sumet += pt;
154  ++nj;
155  }
156  }
157  }
158 
159  if (excludePFMuons_) {
160  for (auto const& j : *pfCandidates) {
161  if (std::abs(j.pdgId()) == 13) {
162  mhx += j.px();
163  mhy += j.py();
164  }
165  }
166  }
167 
168  if (nj < minNJet_) {
169  sumet = 0;
170  mhx = 0;
171  mhy = 0;
172  }
173 
174  reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx * mhx + mhy * mhy));
175  reco::MET::Point vtx(0, 0, 0);
176  reco::MET mht(sumet, p4, vtx);
177  result->push_back(mht);
178 
179  // Put the products into the Event
180  iEvent.put(std::move(result));
181 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::EDGetTokenT< reco::TrackCollection > m_theTrackToken
bool usePFCandidates_
Use PF candidates as input instead of jets.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::InputTag pfRecTracksLabel_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
edm::InputTag pfCandidatesLabel_
auto const & tracks
cannot be loose
std::vector< reco::MET > METCollection
collection of MET objects
Definition: METCollection.h:22
edm::InputTag tracksLabel_
tuple result
Definition: mps_fire.py:311
int iEvent
Definition: GenABIO.cc:224
Definition: MET.h:41
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< reco::PFCandidateCollection > m_thePFCandidateToken
vector< PseudoJet > jets
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::InputTag jetsLabel_
Input jet, track, PFRecTrack, PFCandidate collections.
double minPtJet_
Minimum pt requirement for jets (or objects used as input)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::PFRecTrackCollection > m_theRecTrackToken
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
~HLTTrackMETProducer() override
HLTTrackMETProducer(const edm::ParameterSet &iConfig)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
std::string const & label() const
Definition: InputTag.h:36
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< reco::JetView > m_theJetToken
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
double maxEtaJet_
Maximum (abs) eta requirement for jets (or objects used as input)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)