CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AcceptJet.cc
Go to the documentation of this file.
5 
6 #include<iostream>
7 
8 using namespace std;
9 
10 AcceptJet::AcceptJet(const double& etaMin_, const double& etaMax_, const double& ptMin_, const double& ptMax_,
11  const double& pMin_, const double& pMax_, const double& ratioMin_, const double& ratioMax_, const bool& doJetID_) :
12  etaMin(etaMin_), etaMax(etaMax_), ptRecJetMin(ptMin_), ptRecJetMax(ptMax_), pRecJetMin(pMin_),
13  pRecJetMax(pMax_), ratioMin(ratioMin_), ratioMax(ratioMax_), doJetID(doJetID_) {}
14 
15 
16 bool AcceptJet::operator() (const reco::Jet & jet, const int & jetFlavour, const edm::Handle<reco::SoftLeptonTagInfoCollection> & infos, const double jec) const
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  edm::LogWarning("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 }
70 
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 pRecJetMax
Definition: AcceptJet.h:61
virtual double energy() const GCC11_FINAL
energy
double etaMax
Definition: AcceptJet.h:46
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:138
double ptRecJetMin
Definition: AcceptJet.h:57
double ratioMin
Definition: AcceptJet.h:63
virtual double p() const GCC11_FINAL
magnitude of momentum vector
double ptRecJetMax
Definition: AcceptJet.h:58
bool doJetID
Definition: AcceptJet.h:67
double etaMin
Definition: AcceptJet.h:45
Base class for all types of Jets.
Definition: Jet.h:20
double pRecJetMin
Definition: AcceptJet.h:60
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:151
double ratioMax
Definition: AcceptJet.h:64
Jets made from PFObjects.
Definition: PFJet.h:21
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:146
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
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: AcceptJet.cc:10
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:76
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
int k[5][pyjets_maxn]
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:98
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.
Definition: AcceptJet.cc:16
virtual float pt() const GCC11_FINAL
transverse momentum
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:94