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.
4 
5 #include<iostream>
6 
7 using namespace std;
8 
9 AcceptJet::AcceptJet(const double& etaMin_, const double& etaMax_, const double& ptMin_, const double& ptMax_,
10  const double& pMin_, const double& pMax_, const double& ratioMin_, const double& ratioMax_) :
11  etaMin(etaMin_), etaMax(etaMax_), ptRecJetMin(ptMin_), ptRecJetMax(ptMax_), pRecJetMin(pMin_),
12  pRecJetMax(pMax_), ratioMin(ratioMin_), ratioMax(ratioMax_) {}
13 
14 
15 bool AcceptJet::operator() (const reco::Jet & jet, const int & jetFlavour, const edm::Handle<reco::SoftLeptonTagInfoCollection> & infos) const
16 {
17 
18  // temporary fudge to correct for double loop error
19  // jetPartonMomentum /= 2.0;
20 
21  if ( fabs(jet.eta()) < etaMin ||
22  fabs(jet.eta()) > etaMax ) return false;
23 
24 // if ( jetFlavour.underlyingParton4Vec().P() < pPartonMin ||
25 // jetFlavour.underlyingParton4Vec().P() > pPartonMax ) accept = false;
26 //
27 // if ( jetFlavour.underlyingParton4Vec().Pt() < ptPartonMin ||
28 // jetFlavour.underlyingParton4Vec().Pt() > ptPartonMax ) accept = false;
29 
30  if ( jet.pt() < ptRecJetMin ||
31  jet.pt() > ptRecJetMax ) return false;
32 
33  if ( jet.p() < pRecJetMin ||
34  jet.p() > pRecJetMax ) return false;
35 
36  if ( !infos.isValid() ) {
37  edm::LogWarning("infos not valid") << "A valid SoftLeptonTagInfoCollection was not found!"
38  << " Skipping ratio check.";
39  }
40  else {
41  double pToEratio = ratio( jet, infos );
42  if ( pToEratio < ratioMin ||
43  pToEratio > ratioMax ) return false;
44  }
45 
46  return true;
47 }
48 
50 {
51  double jetRatio = 0.0;
52  reco::SoftLeptonTagInfoCollection::const_iterator infoiter = infos->begin();
53  for( ; infoiter != infos->end(); ++infoiter)
54  {
55  if( reco::deltaR(jet.eta(), jet.phi(), infoiter->jet()->eta(), infoiter->jet()->phi()) > 1e-4 )
56  continue;
57 
58  if( infoiter->leptons() == 0 )
59  break;
60 
61  for( unsigned int k = 0; k != infoiter->leptons(); ++k )
62  {
63  double tempRatio = infoiter->properties(k).ratio;
64  if( tempRatio > jetRatio )
65  jetRatio = tempRatio;
66  }
67  }
68 
69  return jetRatio;
70 }
double pRecJetMax
Definition: AcceptJet.h:59
double etaMax
Definition: AcceptJet.h:44
double ptRecJetMin
Definition: AcceptJet.h:55
double ratioMin
Definition: AcceptJet.h:61
virtual double p() const GCC11_FINAL
magnitude of momentum vector
double ptRecJetMax
Definition: AcceptJet.h:56
double etaMin
Definition: AcceptJet.h:43
Base class for all types of Jets.
Definition: Jet.h:21
double pRecJetMin
Definition: AcceptJet.h:58
double ratioMax
Definition: AcceptJet.h:62
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
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_)
Definition: AcceptJet.cc:9
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
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:49
bool isValid() const
Definition: HandleBase.h:76
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
int k[5][pyjets_maxn]
virtual float pt() const GCC11_FINAL
transverse momentum
bool operator()(const reco::Jet &jet, const int &jetFlavour, const edm::Handle< reco::SoftLeptonTagInfoCollection > &infos) const
Returns true if jet and associated parton satisfy kinematic cuts.
Definition: AcceptJet.cc:15