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