CMS 3D CMS Logo

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