CMS 3D CMS Logo

HLTMETCleanerUsingJetID.cc
Go to the documentation of this file.
1 
11 
12 // Constructor
14  : minPt_(iConfig.getParameter<double>("minPt")),
15  maxEta_(iConfig.getParameter<double>("maxEta")),
16  metLabel_(iConfig.getParameter<edm::InputTag>("metLabel")),
17  jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
18  goodJetsLabel_(iConfig.getParameter<edm::InputTag>("goodJetsLabel")) {
19  m_theMETToken = consumes<reco::CaloMETCollection>(metLabel_);
20  m_theJetToken = consumes<reco::CaloJetCollection>(jetsLabel_);
21  m_theGoodJetToken = consumes<reco::CaloJetCollection>(goodJetsLabel_);
22 
23  // Register the products
24  produces<reco::CaloMETCollection>();
25 }
26 
27 // Destructor
29 
30 // Fill descriptions
33  desc.add<double>("minPt", 20.);
34  desc.add<double>("maxEta", 5.);
35  desc.add<edm::InputTag>("metLabel", edm::InputTag("hltMet"));
36  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4CaloJets"));
37  desc.add<edm::InputTag>("goodJetsLabel", edm::InputTag("hltCaloJetIDPassed"));
38  descriptions.add("hltMETCleanerUsingJetID", desc);
39 }
40 
41 // Produce the products
43  // Create a pointer to the products
44  std::unique_ptr<reco::CaloMETCollection> result(new reco::CaloMETCollection);
45 
49 
50  iEvent.getByToken(m_theMETToken, met);
51  iEvent.getByToken(m_theJetToken, jets);
52  iEvent.getByToken(m_theGoodJetToken, goodJets);
53 
54  double mex_jets = 0.;
55  double mey_jets = 0.;
56  //double sumet_jets = 0.;
57  if (!jets->empty()) {
58  for (auto const& j : *jets) {
59  double pt = j.pt();
60  double eta = j.eta();
61  double px = j.px();
62  double py = j.py();
63 
64  if (pt > minPt_ && std::abs(eta) < maxEta_) {
65  mex_jets -= px;
66  mey_jets -= py;
67  //sumet_jets += pt;
68  }
69  }
70  }
71 
72  double mex_goodJets = 0.;
73  double mey_goodJets = 0.;
74  //double sumet_goodJets = 0.;
75  if (!goodJets->empty()) {
76  for (auto const& j : *goodJets) {
77  double pt = j.pt();
78  double eta = j.eta();
79  double px = j.px();
80  double py = j.py();
81 
82  if (pt > minPt_ && std::abs(eta) < maxEta_) {
83  mex_goodJets -= px;
84  mey_goodJets -= py;
85  //sumet_goodJets += pt;
86  }
87  }
88  }
89 
90  if (!met->empty()) {
91  double mex_diff = mex_goodJets - mex_jets;
92  double mey_diff = mey_goodJets - mey_jets;
93  //double sumet_diff = sumet_goodJets - sumet_jets; // cannot set sumet...
94  reco::Candidate::LorentzVector p4_clean(met->front().px() + mex_diff,
95  mey_diff + met->front().py(),
96  0,
97  sqrt((met->front().px() + mex_diff) * (met->front().px() + mex_diff) +
98  (met->front().py() + mey_diff) * (met->front().py() + mey_diff)));
99 
100  reco::CaloMET cleanmet = met->front();
101  cleanmet.setP4(p4_clean);
102  result->push_back(cleanmet);
103  }
104 
105  iEvent.put(std::move(result));
106 }
~HLTMETCleanerUsingJetID() override
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)
double maxEta_
Maximum (abs) eta requirement for jets.
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::CaloMETCollection > m_theMETToken
edm::InputTag metLabel_
Input tag for the MET collection.
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::CaloJetCollection > m_theGoodJetToken
double minPt_
Minimum pt requirement for jets.
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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:36
HLT enums.
edm::EDGetTokenT< reco::CaloJetCollection > m_theJetToken
def move(src, dest)
Definition: eostools.py:511