CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTMETCleanerUsingJetID.cc
Go to the documentation of this file.
1 
11 
12 
13 // Constructor
15  : minPt_ (iConfig.getParameter<double>("minPt")),
16  maxEta_ (iConfig.getParameter<double>("maxEta")),
17  metLabel_ (iConfig.getParameter<edm::InputTag>("metLabel")),
18  jetsLabel_ (iConfig.getParameter<edm::InputTag>("jetsLabel")),
19  goodJetsLabel_ (iConfig.getParameter<edm::InputTag>("goodJetsLabel")) {
20  m_theMETToken = consumes<reco::CaloMETCollection>(metLabel_);
21  m_theJetToken = consumes<reco::CaloJetCollection>(jetsLabel_);
22  m_theGoodJetToken = consumes<reco::CaloJetCollection>(goodJetsLabel_);
23 
24  // Register the products
25  produces<reco::CaloMETCollection>();
26 }
27 
28 // Destructor
30 
31 // Fill descriptions
34  desc.add<double>("minPt", 20.);
35  desc.add<double>("maxEta", 5.);
36  desc.add<edm::InputTag>("metLabel", edm::InputTag("hltMet"));
37  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4CaloJets"));
38  desc.add<edm::InputTag>("goodJetsLabel", edm::InputTag("hltCaloJetIDPassed"));
39  descriptions.add("hltMETCleanerUsingJetID",desc);
40 }
41 
42 // Produce the products
44 
45  // Create a pointer to the products
46  std::auto_ptr<reco::CaloMETCollection> result(new reco::CaloMETCollection);
47 
51 
52  iEvent.getByToken(m_theMETToken, met);
53  iEvent.getByToken(m_theJetToken, jets);
54  iEvent.getByToken(m_theGoodJetToken, goodJets);
55 
56  double mex_jets = 0.;
57  double mey_jets = 0.;
58  double sumet_jets = 0.;
59  if (jets->size() > 0 ) {
60  for(reco::CaloJetCollection::const_iterator j = jets->begin(); j != jets->end(); ++j) {
61  double pt = j->pt() ;
62  double eta = j->eta();
63  double px = j->px() ;
64  double py = j->py() ;
65 
66  if (pt > minPt_ && std::abs(eta) < maxEta_) {
67  mex_jets -= px;
68  mey_jets -= py;
69  sumet_jets += pt;
70  }
71  }
72  }
73 
74  double mex_goodJets = 0.;
75  double mey_goodJets = 0.;
76  double sumet_goodJets = 0.;
77  if (goodJets->size() > 0) {
78  for(reco::CaloJetCollection::const_iterator j = goodJets->begin(); j != goodJets->end(); ++j) {
79  double pt = j->pt() ;
80  double eta = j->eta();
81  double px = j->px() ;
82  double py = j->py() ;
83 
84  if (pt > minPt_ && std::abs(eta) < maxEta_) {
85  mex_goodJets -= px;
86  mey_goodJets -= py;
87  sumet_goodJets += pt;
88  }
89  }
90  }
91 
92  if (met->size() > 0) {
93  double mex_diff = mex_goodJets - mex_jets;
94  double mey_diff = mey_goodJets - mey_jets;
95  //double sumet_diff = sumet_goodJets - sumet_jets; // cannot set sumet...
96  reco::Candidate::LorentzVector p4_clean(met->front().px()+mex_diff, mey_diff+met->front().py(), 0, sqrt((met->front().px()+mex_diff)*(met->front().px() +mex_diff)+(met->front().py()+mey_diff)*(met->front().py() +mey_diff)));
97 
98  reco::CaloMET cleanmet = met->front() ;
99  cleanmet.setP4(p4_clean);
100  result->push_back(cleanmet);
101  }
102 
103  iEvent.put( result );
104 }
edm::InputTag jetsLabel_
Input tag for the &#39;all jets&#39; collection.
edm::InputTag goodJetsLabel_
Input tag for the &#39;good jets&#39; collection.
HLTMETCleanerUsingJetID(const edm::ParameterSet &iConfig)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
double maxEta_
Maximum (abs) eta requirement for jets.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::CaloMETCollection > m_theMETToken
tuple result
Definition: mps_fire.py:95
edm::InputTag metLabel_
Input tag for the MET collection.
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::CaloJetCollection > m_theGoodJetToken
double minPt_
Minimum pt requirement for jets.
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< reco::CaloMET > CaloMETCollection
collection of CaloMET objects
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
edm::EDGetTokenT< reco::CaloJetCollection > m_theJetToken