test
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  : usePt_ (iConfig.getParameter<bool>("usePt")),
16  minPt_ (iConfig.getParameter<double>("minPt")),
17  maxEta_ (iConfig.getParameter<double>("maxEta")),
18  metLabel_ (iConfig.getParameter<edm::InputTag>("metLabel")),
19  jetsLabel_ (iConfig.getParameter<edm::InputTag>("jetsLabel")),
20  goodJetsLabel_ (iConfig.getParameter<edm::InputTag>("goodJetsLabel")) {
21  m_theMETToken = consumes<reco::CaloMETCollection>(metLabel_);
22  m_theJetToken = consumes<reco::CaloJetCollection>(jetsLabel_);
23  m_theGoodJetToken = consumes<reco::CaloJetCollection>(goodJetsLabel_);
24 
25  // Register the products
26  produces<reco::CaloMETCollection>();
27 }
28 
29 // Destructor
31 
32 // Fill descriptions
35  desc.add<bool>("usePt", false);
36  desc.add<double>("minPt", 20.);
37  desc.add<double>("maxEta", 5.);
38  desc.add<edm::InputTag>("metLabel", edm::InputTag("hltMet"));
39  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4CaloJets"));
40  desc.add<edm::InputTag>("goodJetsLabel", edm::InputTag("hltCaloJetIDPassed"));
41  descriptions.add("hltMETCleanerUsingJetID",desc);
42 }
43 
44 // Produce the products
46 
47  // Create a pointer to the products
48  std::auto_ptr<reco::CaloMETCollection> result(new reco::CaloMETCollection);
49 
53 
54  iEvent.getByToken(m_theMETToken, met);
55  iEvent.getByToken(m_theJetToken, jets);
56  iEvent.getByToken(m_theGoodJetToken, goodJets);
57 
58  double mex_jets = 0.;
59  double mey_jets = 0.;
60  double sumet_jets = 0.;
61  if (jets->size() > 0 ) {
62  for(reco::CaloJetCollection::const_iterator j = jets->begin(); j != jets->end(); ++j) {
63  double pt = usePt_ ? j->pt() : j->et();
64  double eta = j->eta();
65  double phi = j->phi();
66  double px = usePt_ ? j->px() : j->et() * cos(phi);
67  double py = usePt_ ? j->py() : j->et() * sin(phi);
68 
69  if (pt > minPt_ && std::abs(eta) < maxEta_) {
70  mex_jets -= px;
71  mey_jets -= py;
72  sumet_jets += pt;
73  }
74  }
75  }
76 
77  double mex_goodJets = 0.;
78  double mey_goodJets = 0.;
79  double sumet_goodJets = 0.;
80  if (goodJets->size() > 0) {
81  for(reco::CaloJetCollection::const_iterator j = goodJets->begin(); j != goodJets->end(); ++j) {
82  double pt = usePt_ ? j->pt() : j->et();
83  double eta = j->eta();
84  double phi = j->phi();
85  double px = usePt_ ? j->px() : j->pt() * cos(phi);
86  double py = usePt_ ? j->py() : j->pt() * sin(phi);
87 
88  if (pt > minPt_ && std::abs(eta) < maxEta_) {
89  mex_goodJets -= px;
90  mey_goodJets -= py;
91  sumet_goodJets += pt;
92  }
93  }
94  }
95 
96  if (met->size() > 0) {
97  double mex_diff = mex_goodJets - mex_jets;
98  double mey_diff = mey_goodJets - mey_jets;
99  //double sumet_diff = sumet_goodJets - sumet_jets; // cannot set sumet...
100  reco::Candidate::LorentzVector p4_diff(mex_diff, mey_diff, 0, sqrt(mex_diff*mex_diff + mey_diff*mey_diff));
101 
102  reco::CaloMET cleanmet = met->front();
103  cleanmet.setP4(cleanmet.p4() + p4_diff);
104  result->push_back(cleanmet);
105  }
106 
107  iEvent.put( result );
108 }
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)
tuple met
____________________________________________________________________________||
Definition: CaloMET_cfi.py:4
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double maxEta_
Maximum (abs) eta requirement for jets.
virtual void setP4(const LorentzVector &p4)
set 4-momentum
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T eta() const
edm::EDGetTokenT< reco::CaloMETCollection > m_theMETToken
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:116
T sqrt(T t)
Definition: SSEVec.h:48
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
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:41
edm::EDGetTokenT< reco::CaloJetCollection > m_theJetToken
bool usePt_
Use pt; otherwise, use et.
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: DDAxes.h:10