CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHtMhtProducer.cc
Go to the documentation of this file.
1 
3 
11 
14 
15 
17  usePt_ ( iConfig.getParameter<bool>("usePt") ),
18  useTracks_ ( iConfig.getParameter<bool>("useTracks") ),
19  excludePFMuons_ ( iConfig.getParameter<bool>("excludePFMuons") ),
20  minNJetHt_ ( iConfig.getParameter<int>("minNJetHt") ),
21  minNJetMht_ ( iConfig.getParameter<int>("minNJetMht") ),
22  minPtJetHt_ ( iConfig.getParameter<double>("minPtJetHt") ),
23  minPtJetMht_ ( iConfig.getParameter<double>("minPtJetMht") ),
24  maxEtaJetHt_ ( iConfig.getParameter<double>("maxEtaJetHt") ),
25  maxEtaJetMht_ ( iConfig.getParameter<double>("maxEtaJetMht") ),
26  jetsLabel_ ( iConfig.getParameter<edm::InputTag>("jetsLabel") ),
27  tracksLabel_ ( iConfig.getParameter<edm::InputTag>("tracksLabel") ),
28  pfCandidatesLabel_ ( iConfig.getParameter<edm::InputTag>("pfCandidatesLabel") )
29 {
30  m_theJetToken = consumes<edm::View<reco::Jet>>(jetsLabel_);
31  m_theTrackToken = consumes<reco::TrackCollection>(tracksLabel_);
32  m_thePfCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
33  produces<reco::METCollection>();
34 }
35 
36 
38 }
39 
40 
43  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltCaloJetCorrected"));
44  desc.add<bool>("usePt", true);
45  desc.add<int>("minNJetHt", 0);
46  desc.add<int>("minNJetMht", 0);
47  desc.add<double>("minPtJetHt", 40);
48  desc.add<double>("minPtJetMht", 30);
49  desc.add<double>("maxEtaJetHt", 3);
50  desc.add<double>("maxEtaJetMht", 999);
51  desc.add<bool>("useTracks", false);
52  desc.add<edm::InputTag>("tracksLabel", edm::InputTag("hltL3Muons"));
53  desc.add<bool>("excludePFMuons", false);
54  desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
55  descriptions.add("hltHtMhtProducer", desc);
56 }
57 
58 
60 
61  std::auto_ptr<reco::METCollection> metobject(new reco::METCollection());
62 
63  //edm::Handle<reco::CaloJetCollection> jets;
65  iEvent.getByToken(m_theJetToken, jets);
67  if (useTracks_) iEvent.getByToken(m_theTrackToken, tracks);
69  if (excludePFMuons_) iEvent.getByToken(m_thePfCandidateToken, pfCandidates);
70 
71  int nj_ht = 0, nj_mht = 0;
72  double ht=0.;
73  double mhtx=0., mhty=0.;
74 
75  //for (reco::CaloJetCollection::const_iterator jet = jets->begin(); jet != jets->end(); jet++) {
76  for(edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); jet++ ) {
77  double mom = (usePt_ ? jet->pt() : jet->et());
78  if (mom > minPtJetHt_ and std::abs(jet->eta()) < maxEtaJetHt_) {
79  ht += mom;
80  ++nj_ht;
81  }
82  if (mom > minPtJetMht_ and std::abs(jet->eta()) < maxEtaJetMht_) {
83  mhtx -= mom*cos(jet->phi());
84  mhty -= mom*sin(jet->phi());
85  ++nj_mht;
86  }
87  }
88  if (useTracks_) {
89  for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); track++) {
90  if (track->pt() > minPtJetHt_ and std::abs(track->eta()) < maxEtaJetHt_) {
91  ht += track->pt();
92  }
93  if (track->pt() > minPtJetMht_ and std::abs(track->eta()) < maxEtaJetMht_) {
94  mhtx -= track->px();
95  mhty -= track->py();
96  }
97  }
98  }
99  if (excludePFMuons_) {
100  reco::PFCandidateCollection::const_iterator i (pfCandidates->begin());
101  for (; i != pfCandidates->end(); ++i) {
102  if(abs(i->pdgId())==13){
103  mhtx += (i->pt())*cos(i->phi());
104  mhty += (i->pt())*sin(i->phi());
105  }
106  }
107  }
108 
109  if (nj_ht < minNJetHt_ ) { ht = 0; }
110  if (nj_mht < minNJetMht_) { mhtx = 0; mhty = 0; }
111 
112  metobject->push_back(
113  reco::MET(
114  ht,
115  reco::MET::LorentzVector(mhtx, mhty, 0, sqrt(mhtx*mhtx + mhty*mhty)),
117  )
118  );
119 
120  iEvent.put(metobject);
121 }
int i
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::InputTag pfCandidatesLabel_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
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:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
Definition: MET.h:32
edm::EDGetTokenT< reco::PFCandidateCollection > m_thePfCandidateToken
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::EDGetTokenT< reco::TrackCollection > m_theTrackToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< edm::View< reco::Jet > > m_theJetToken
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
tuple tracks
Definition: testEve_cfg.py:39
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
edm::InputTag tracksLabel_
HLTHtMhtProducer(const edm::ParameterSet &iConfig)
edm::InputTag jetsLabel_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)