CMS 3D CMS Logo

EtaPtBin.cc
Go to the documentation of this file.
2 
3 #include <algorithm>
4 #include <sstream>
5 
6 
7 
8 
9 EtaPtBin::EtaPtBin ( const bool& etaActive_ , const double& etaMin_ , const double& etaMax_ ,
10  const bool& ptActive_ , const double& ptMin_ , const double& ptMax_ )
11  : etaActive ( etaActive_ ) , etaMin ( etaMin_ ) , etaMax ( etaMax_ ) ,
12  ptActive ( ptActive_ ) , ptMin ( ptMin_ ) , ptMax ( ptMax_ ) {
13 
15  ptActive , ptMin , ptMax );
16 }
17 
18 
19 
21  ( const bool& etaActive_ , const double& etaMin_ , const double& etaMax_ ,
22  const bool& ptActive_ , const double& ptMin_ , const double& ptMax_)
23 {
24  // create string only from the active parts
25  std::stringstream stream ( "" );
26 
27  if ( etaActive_ ) {
28  stream << "_ETA_" << etaMin_ << "-" << etaMax_;
29  }
30 
31  if ( ptActive_ ) {
32  stream << "_PT_" << ptMin_ << "-" << ptMax_;
33  }
34  if (!(etaActive_||ptActive_)) stream << "_GLOBAL";
35 
36  std::string descr(stream.str());
37  // remove blanks which are introduced when adding doubles
38  std::remove(descr.begin(), descr.end(), ' ');
39  std::replace(descr.begin(), descr.end(), '.' , 'v' );
40 
41  return descr;
42 }
43 
44 bool EtaPtBin::inBin(const reco::Jet & jet, const double jec) const
45 {
46  return inBin(jet.eta(), jet.pt()*jec);
47 }
48 
49 // bool EtaPtBin::inBin(const BTagMCTools::JetFlavour & jetFlavour) const
50 // {
51 // return inBin(jetFlavour.underlyingParton4Vec().Eta(),
52 // jetFlavour.underlyingParton4Vec().Pt());
53 // }
54 
55 
56 bool EtaPtBin::inBin (const double & eta , const double & pt ) const {
57  if ( etaActive ) {
58  if ( fabs(eta) < etaMin ) return false;
59  if ( fabs(eta) > etaMax ) return false;
60  }
61 
62  if ( ptActive ) {
63  if ( pt < ptMin ) return false;
64  if ( pt > ptMax ) return false;
65  }
66 
67  return true;
68 }
double eta() const final
momentum pseudorapidity
bool inBin(const reco::Jet &jet, const double jec) const
Definition: EtaPtBin.cc:44
Base class for all types of Jets.
Definition: Jet.h:20
double ptMin
Definition: EtaPtBin.h:58
std::string descriptionString
Definition: EtaPtBin.h:63
def replace(string, replacements)
double pt() const final
transverse momentum
static std::string buildDescriptionString(const bool &etaActive_, const double &etaMin_, const double &etaMax_, const bool &ptActive_, const double &ptMin_, const double &ptMax_)
Definition: EtaPtBin.cc:21
bool inBin(const double &eta, const double &pt) const
Check if jet/parton are within rapidity/pt cuts.
Definition: EtaPtBin.cc:56
double etaMax
Definition: EtaPtBin.h:55
double etaMin
Definition: EtaPtBin.h:54
EtaPtBin(const bool &etaActive_, const double &etaMin_, const double &etaMax_, const bool &ptActive_, const double &ptMin_, const double &ptMax_)
Definition: EtaPtBin.cc:9
bool ptActive
Definition: EtaPtBin.h:57
bool etaActive
Definition: EtaPtBin.h:53
double ptMax
Definition: EtaPtBin.h:59
def remove(d, key, TELL=False)
Definition: MatrixUtil.py:209