CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
AcceptJet Class Reference

#include <AcceptJet.h>

Public Member Functions

 AcceptJet (const double &etaMin_, const double &etaMax_, const double &ptMin_, const double &ptMax_, const double &pMin_, const double &pMax_, const double &ratioMin_, const double &ratioMax_, const bool &doJetID_)
 
bool operator() (const reco::Jet &jet, const int &jetFlavour, const edm::Handle< reco::SoftLeptonTagInfoCollection > &infos, const double jec) const
 Returns true if jet and associated parton satisfy kinematic cuts. More...
 
void setDoJetID (bool b)
 
void setEtaMax (double d)
 
void setEtaMin (double d)
 Set cut parameters. More...
 
void setPRecJetMax (double d)
 
void setPRecJetMin (double d)
 
void setPtRecJetMax (double d)
 
void setPtRecJetMin (double d)
 
void setRatioMax (double d)
 
void setRatioMin (double d)
 

Protected Member Functions

double ratio (const reco::Jet &jet, const edm::Handle< reco::SoftLeptonTagInfoCollection > &infos) const
 Finds the ratio of the momentum of any leptons in the jet to jet energy. More...
 

Protected Attributes

bool doJetID
 
double etaMax
 
double etaMin
 
double pRecJetMax
 
double pRecJetMin
 
double ptRecJetMax
 
double ptRecJetMin
 
double ratioMax
 
double ratioMin
 

Detailed Description

Decide if jet and associated parton satisfy desired kinematic cuts.

Definition at line 14 of file AcceptJet.h.

Constructor & Destructor Documentation

AcceptJet::AcceptJet ( const double &  etaMin_,
const double &  etaMax_,
const double &  ptMin_,
const double &  ptMax_,
const double &  pMin_,
const double &  pMax_,
const double &  ratioMin_,
const double &  ratioMax_,
const bool &  doJetID_ 
)

Definition at line 10 of file AcceptJet.cc.

11  :
12  etaMin(etaMin_), etaMax(etaMax_), ptRecJetMin(ptMin_), ptRecJetMax(ptMax_), pRecJetMin(pMin_),
13  pRecJetMax(pMax_), ratioMin(ratioMin_), ratioMax(ratioMax_), doJetID(doJetID_) {}
double pRecJetMax
Definition: AcceptJet.h:61
double etaMax
Definition: AcceptJet.h:46
double ptRecJetMin
Definition: AcceptJet.h:57
double ratioMin
Definition: AcceptJet.h:63
double ptRecJetMax
Definition: AcceptJet.h:58
bool doJetID
Definition: AcceptJet.h:67
double etaMin
Definition: AcceptJet.h:45
double pRecJetMin
Definition: AcceptJet.h:60
double ratioMax
Definition: AcceptJet.h:64

Member Function Documentation

bool AcceptJet::operator() ( const reco::Jet jet,
const int &  jetFlavour,
const edm::Handle< reco::SoftLeptonTagInfoCollection > &  infos,
const double  jec 
) const

Returns true if jet and associated parton satisfy kinematic cuts.

Definition at line 16 of file AcceptJet.cc.

References reco::PFJet::chargedEmEnergy(), reco::PFJet::chargedHadronEnergy(), reco::btau::chargedHadronEnergyFraction, reco::PFJet::chargedMultiplicity(), doJetID, reco::LeafCandidate::energy(), reco::LeafCandidate::eta(), etaMax, etaMin, reco::PFJet::getPFConstituents(), edm::HandleBase::isValid(), patTestJEC_cfi::jec, metsig::jet, LogDebug, jets_cff::nConstituents, reco::PFJet::neutralEmEnergy(), reco::PFJet::neutralHadronEnergy(), reco::btau::neutralHadronEnergyFraction, reco::LeafCandidate::p(), pRecJetMax, pRecJetMin, reco::LeafCandidate::pt(), ptRecJetMax, ptRecJetMin, ratio(), ratioMax, and ratioMin.

17 {
18 
19  // temporary fudge to correct for double loop error
20  // jetPartonMomentum /= 2.0;
21 
22  if ( fabs(jet.eta()) < etaMin ||
23  fabs(jet.eta()) > etaMax ) return false;
24 
25 // if ( jetFlavour.underlyingParton4Vec().P() < pPartonMin ||
26 // jetFlavour.underlyingParton4Vec().P() > pPartonMax ) accept = false;
27 //
28 // if ( jetFlavour.underlyingParton4Vec().Pt() < ptPartonMin ||
29 // jetFlavour.underlyingParton4Vec().Pt() > ptPartonMax ) accept = false;
30 
31  if ( jet.pt()*jec < ptRecJetMin ||
32  jet.pt()*jec > ptRecJetMax ) return false;
33 
34  if ( jet.p()*jec < pRecJetMin ||
35  jet.p()*jec > pRecJetMax ) return false;
36 
37  if ( !infos.isValid() ) {
38  LogDebug("infos not valid") << "A valid SoftLeptonTagInfoCollection was not found!"
39  << " Skipping ratio check.";
40  }
41  else {
42  double pToEratio = ratio( jet, infos );
43  if ( pToEratio/jec < ratioMin ||
44  pToEratio/jec > ratioMax ) return false;
45  }
46 
47  if(doJetID){
48  const reco::PFJet * pfjet = dynamic_cast<const reco::PFJet *>(&jet);
49  if(pfjet) {
51  double neutralEmEnergyFraction = pfjet->neutralEmEnergy()/(jet.energy()*jec);
52  int nConstituents = pfjet->getPFConstituents().size();
54  int chargedMultiplicity = pfjet->chargedMultiplicity();
55  double chargedEmEnergyFraction = pfjet->chargedEmEnergy()/(jet.energy()*jec);
56  if(!(neutralHadronEnergyFraction < 0.99
57  && neutralEmEnergyFraction < 0.99
58  && nConstituents > 1
59  && chargedHadronEnergyFraction > 0.0
60  && chargedMultiplicity > 0.0
61  && chargedEmEnergyFraction < 0.99))
62  return false; //2012 values
63  }
64  //Looks to not work -> to be fixed
65  //else std::cout<<"Jets are not PF jets, put 'doJetID' to False."<<std::endl;
66  }
67 
68  return true;
69 }
#define LogDebug(id)
double pRecJetMax
Definition: AcceptJet.h:61
double etaMax
Definition: AcceptJet.h:46
double eta() const final
momentum pseudorapidity
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:142
double ptRecJetMin
Definition: AcceptJet.h:57
double ratioMin
Definition: AcceptJet.h:63
double ptRecJetMax
Definition: AcceptJet.h:58
bool doJetID
Definition: AcceptJet.h:67
double etaMin
Definition: AcceptJet.h:45
double pRecJetMin
Definition: AcceptJet.h:60
double pt() const final
transverse momentum
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:155
double ratioMax
Definition: AcceptJet.h:64
Jets made from PFObjects.
Definition: PFJet.h:21
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:150
double energy() const final
energy
double ratio(const reco::Jet &jet, const edm::Handle< reco::SoftLeptonTagInfoCollection > &infos) const
Finds the ratio of the momentum of any leptons in the jet to jet energy.
Definition: AcceptJet.cc:71
bool isValid() const
Definition: HandleBase.h:74
nConstituents
Definition: jets_cff.py:165
double p() const final
magnitude of momentum vector
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:102
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:98
double AcceptJet::ratio ( const reco::Jet jet,
const edm::Handle< reco::SoftLeptonTagInfoCollection > &  infos 
) const
protected

Finds the ratio of the momentum of any leptons in the jet to jet energy.

Definition at line 71 of file AcceptJet.cc.

References reco::deltaR(), MillePedeFileConverter_cfg::e, reco::LeafCandidate::eta(), gen::k, and reco::LeafCandidate::phi().

Referenced by operator()(), and setDoJetID().

72 {
73  double jetRatio = 0.0;
74  reco::SoftLeptonTagInfoCollection::const_iterator infoiter = infos->begin();
75  for( ; infoiter != infos->end(); ++infoiter)
76  {
77  if( reco::deltaR(jet.eta(), jet.phi(), infoiter->jet()->eta(), infoiter->jet()->phi()) > 1e-4 )
78  continue;
79 
80  if( infoiter->leptons() == 0 )
81  break;
82 
83  for( unsigned int k = 0; k != infoiter->leptons(); ++k )
84  {
85  double tempRatio = infoiter->properties(k).ratio;
86  if( tempRatio > jetRatio )
87  jetRatio = tempRatio;
88  }
89  }
90 
91  return jetRatio;
92 }
double eta() const final
momentum pseudorapidity
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
int k[5][pyjets_maxn]
double phi() const final
momentum azimuthal angle
void AcceptJet::setDoJetID ( bool  b)
inline

Definition at line 34 of file AcceptJet.h.

References b, doJetID, and ratio().

34 { doJetID = b ; }
bool doJetID
Definition: AcceptJet.h:67
double b
Definition: hdecay.h:120
void AcceptJet::setEtaMax ( double  d)
inline

Definition at line 25 of file AcceptJet.h.

References edmIntegrityCheck::d, and etaMax.

25 { etaMax = d ; }
double etaMax
Definition: AcceptJet.h:46
void AcceptJet::setEtaMin ( double  d)
inline

Set cut parameters.

Definition at line 24 of file AcceptJet.h.

References edmIntegrityCheck::d, and etaMin.

24 { etaMin = d ; }
double etaMin
Definition: AcceptJet.h:45
void AcceptJet::setPRecJetMax ( double  d)
inline

Definition at line 31 of file AcceptJet.h.

References edmIntegrityCheck::d, and pRecJetMax.

void AcceptJet::setPRecJetMin ( double  d)
inline

Definition at line 30 of file AcceptJet.h.

References edmIntegrityCheck::d, and pRecJetMin.

30 { pRecJetMin = d ; }
double pRecJetMin
Definition: AcceptJet.h:60
void AcceptJet::setPtRecJetMax ( double  d)
inline

Definition at line 29 of file AcceptJet.h.

References edmIntegrityCheck::d, and ptRecJetMax.

29 { ptRecJetMax = d ; }
double ptRecJetMax
Definition: AcceptJet.h:58
void AcceptJet::setPtRecJetMin ( double  d)
inline

Definition at line 28 of file AcceptJet.h.

References edmIntegrityCheck::d, and ptRecJetMin.

28 { ptRecJetMin = d ; }
double ptRecJetMin
Definition: AcceptJet.h:57
void AcceptJet::setRatioMax ( double  d)
inline

Definition at line 33 of file AcceptJet.h.

References edmIntegrityCheck::d, and ratioMax.

33 { ratioMax = d ; }
double ratioMax
Definition: AcceptJet.h:64
void AcceptJet::setRatioMin ( double  d)
inline

Definition at line 32 of file AcceptJet.h.

References edmIntegrityCheck::d, and ratioMin.

32 { ratioMin = d ; }
double ratioMin
Definition: AcceptJet.h:63

Member Data Documentation

bool AcceptJet::doJetID
protected

Definition at line 67 of file AcceptJet.h.

Referenced by operator()(), and setDoJetID().

double AcceptJet::etaMax
protected

Definition at line 46 of file AcceptJet.h.

Referenced by operator()(), and setEtaMax().

double AcceptJet::etaMin
protected

Definition at line 45 of file AcceptJet.h.

Referenced by operator()(), and setEtaMin().

double AcceptJet::pRecJetMax
protected

Definition at line 61 of file AcceptJet.h.

Referenced by operator()(), and setPRecJetMax().

double AcceptJet::pRecJetMin
protected

Definition at line 60 of file AcceptJet.h.

Referenced by operator()(), and setPRecJetMin().

double AcceptJet::ptRecJetMax
protected

Definition at line 58 of file AcceptJet.h.

Referenced by operator()(), and setPtRecJetMax().

double AcceptJet::ptRecJetMin
protected

Definition at line 57 of file AcceptJet.h.

Referenced by operator()(), and setPtRecJetMin().

double AcceptJet::ratioMax
protected

Definition at line 64 of file AcceptJet.h.

Referenced by operator()(), and setRatioMax().

double AcceptJet::ratioMin
protected

Definition at line 63 of file AcceptJet.h.

Referenced by operator()(), and setRatioMin().