CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
pat::GenJetMatcher Class Reference

#include <SmearedJetProducerT.h>

Public Member Functions

 GenJetMatcher (const edm::ParameterSet &cfg, edm::ConsumesCollector &&collector)
 
void getTokens (const edm::Event &event)
 
template<class T >
const reco::GenJetmatch (const T &jet, double resolution)
 

Static Public Member Functions

static void fillDescriptions (edm::ParameterSetDescription &desc)
 

Private Attributes

double m_dPt_max_factor
 
double m_dR_max
 
edm::Handle< reco::GenJetCollectionm_genJets
 
edm::EDGetTokenT< reco::GenJetCollectionm_genJetsToken
 

Detailed Description

Definition at line 40 of file SmearedJetProducerT.h.

Constructor & Destructor Documentation

◆ GenJetMatcher()

pat::GenJetMatcher::GenJetMatcher ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  collector 
)
inline

Definition at line 42 of file SmearedJetProducerT.h.

43  : m_genJetsToken(collector.consumes<reco::GenJetCollection>(cfg.getParameter<edm::InputTag>("genJets"))),
44  m_dR_max(cfg.getParameter<double>("dRMax")),
45  m_dPt_max_factor(cfg.getParameter<double>("dPtMaxFactor")) {
46  // Empty
47  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< GenJet > GenJetCollection
collection of GenJet objects
edm::EDGetTokenT< reco::GenJetCollection > m_genJetsToken

Member Function Documentation

◆ fillDescriptions()

static void pat::GenJetMatcher::fillDescriptions ( edm::ParameterSetDescription desc)
inlinestatic

Definition at line 49 of file SmearedJetProducerT.h.

References submitPVResolutionJobs::desc.

Referenced by SmearedJetProducerT< T >::fillDescriptions().

49  {
50  desc.add<edm::InputTag>("genJets");
51  desc.add<double>("dRMax");
52  desc.add<double>("dPtMaxFactor");
53  }

◆ getTokens()

void pat::GenJetMatcher::getTokens ( const edm::Event event)
inline

Definition at line 55 of file SmearedJetProducerT.h.

References m_genJets, and m_genJetsToken.

55 { event.getByToken(m_genJetsToken, m_genJets); }
edm::Handle< reco::GenJetCollection > m_genJets
edm::EDGetTokenT< reco::GenJetCollection > m_genJetsToken

◆ match()

template<class T >
const reco::GenJet* pat::GenJetMatcher::match ( const T jet,
double  resolution 
)
inline

Definition at line 58 of file SmearedJetProducerT.h.

References funct::abs(), PbPb_ZMuSkimMuonDPG_cff::deltaR, HGC3DClusterGenMatchSelector_cfi::dR, l1tCaloJetHTTProducer_cfi::genJets, infinity, metsig::jet, m_dPt_max_factor, m_dR_max, m_genJets, and L1TObjectsTimingClient_cff::resolution.

58  {
60 
61  // Try to find a gen jet matching
62  // dR < m_dR_max
63  // dPt < m_dPt_max_factor * resolution
64 
65  double min_dR = std::numeric_limits<double>::infinity();
66  const reco::GenJet* matched_genJet = nullptr;
67 
68  for (const auto& genJet : genJets) {
69  double dR = deltaR(genJet, jet);
70 
71  if (dR > min_dR)
72  continue;
73 
74  if (dR < m_dR_max) {
75  double dPt = std::abs(genJet.pt() - jet.pt());
76  if (dPt > m_dPt_max_factor * resolution)
77  continue;
78 
79  min_dR = dR;
80  matched_genJet = &genJet;
81  }
82  }
83 
84  return matched_genJet;
85  }
edm::Handle< reco::GenJetCollection > m_genJets
std::vector< GenJet > GenJetCollection
collection of GenJet objects
const double infinity
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23

Member Data Documentation

◆ m_dPt_max_factor

double pat::GenJetMatcher::m_dPt_max_factor
private

Definition at line 92 of file SmearedJetProducerT.h.

Referenced by match().

◆ m_dR_max

double pat::GenJetMatcher::m_dR_max
private

Definition at line 91 of file SmearedJetProducerT.h.

Referenced by match().

◆ m_genJets

edm::Handle<reco::GenJetCollection> pat::GenJetMatcher::m_genJets
private

Definition at line 89 of file SmearedJetProducerT.h.

Referenced by getTokens(), and match().

◆ m_genJetsToken

edm::EDGetTokenT<reco::GenJetCollection> pat::GenJetMatcher::m_genJetsToken
private

Definition at line 88 of file SmearedJetProducerT.h.

Referenced by getTokens().