CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Private Member Functions | Friends
BTagEntry Class Reference

#include <BTagEntry.h>

Classes

struct  Parameters
 

Public Types

enum  JetFlavor { FLAV_B =0, FLAV_C =1, FLAV_UDSG =2 }
 
enum  OperatingPoint { OP_LOOSE =0, OP_MEDIUM =1, OP_TIGHT =2, OP_RESHAPING =3 }
 

Public Member Functions

 BTagEntry ()
 
 BTagEntry (const std::string &csvLine, bool validate)
 
 BTagEntry (const std::string &func, Parameters p)
 
 BTagEntry (const TF1 *func, Parameters p)
 
 BTagEntry (const TH1 *histo, Parameters p)
 
std::string makeCSVLine () const
 
 ~BTagEntry ()
 

Static Public Member Functions

static std::string makeCSVHeader ()
 
static std::string trimStr (std::string str)
 

Public Attributes

std::string formula
 
Parameters params
 

Private Member Functions

template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Friends

class boost::serialization::access
 
template<typename CondSerializationT , typename Enabled >
struct cond::serialization::access
 

Detailed Description

BTagEntry

Represents one pt- or discriminator-dependent calibration function.

measurement_type: e.g. comb, ttbar, di-mu, boosted, ... sys_type: e.g. central, plus, minus, plus_JEC, plus_JER, ...

Everything is converted into a function, as it is easiest to store it in a txt or json file.

Definition at line 24 of file BTagEntry.h.

Member Enumeration Documentation

Enumerator
FLAV_B 
FLAV_C 
FLAV_UDSG 

Definition at line 33 of file BTagEntry.h.

33  {
34  FLAV_B=0,
35  FLAV_C=1,
36  FLAV_UDSG=2,
37  };
Enumerator
OP_LOOSE 
OP_MEDIUM 
OP_TIGHT 
OP_RESHAPING 

Definition at line 27 of file BTagEntry.h.

Constructor & Destructor Documentation

BTagEntry::BTagEntry ( )
inline
BTagEntry::BTagEntry ( const std::string &  csvLine,
bool  validate 
)

Definition at line 36 of file BTagEntry.cc.

References begin, end, Exception, connectstrParser::f1, formula, mps_fire::i, params, AlCaHLTBitMon_QueryRunRegistry::string, and trimStr().

36  {
37  // make tokens
38  std::stringstream buff(csvLine);
39  std::vector<std::string> vec;
40  std::string token;
41  while (std::getline(buff, token, ","[0])) {
42  token = BTagEntry::trimStr(token);
43  if (token.empty()) {
44  continue;
45  }
46  vec.push_back(token);
47  }
48  if (vec.size() != 11) {
49  throw cms::Exception("BTagCalibration")
50  << "Invalid csv line; num tokens != 11: "
51  << csvLine;
52  }
53 
54  // clean string values
55  char chars[] = " \"\n";
56  for (unsigned int i = 0; i < strlen(chars); ++i) {
57  vec[1].erase(remove(vec[1].begin(),vec[1].end(),chars[i]),vec[1].end());
58  vec[2].erase(remove(vec[2].begin(),vec[2].end(),chars[i]),vec[2].end());
59  vec[10].erase(remove(vec[10].begin(),vec[10].end(),chars[i]),vec[10].end());
60  }
61 
62  // make formula
63  formula = vec[10];
64  if (validate) {
65  TF1 f1("", formula.c_str()); // compile formula to check validity
66  if (f1.IsZombie()) {
67  throw cms::Exception("BTagCalibration") << "Invalid csv line; formula does not compile: " << csvLine;
68  }
69  }
70 
71  // make parameters
72  unsigned op = stoi(vec[0]);
73  if (op > 3) {
74  throw cms::Exception("BTagCalibration")
75  << "Invalid csv line; OperatingPoint > 3: "
76  << csvLine;
77  }
78  unsigned jf = stoi(vec[3]);
79  if (jf > 2) {
80  throw cms::Exception("BTagCalibration")
81  << "Invalid csv line; JetFlavor > 2: "
82  << csvLine;
83  }
86  vec[1],
87  vec[2],
89  stof(vec[4]),
90  stof(vec[5]),
91  stof(vec[6]),
92  stof(vec[7]),
93  stof(vec[8]),
94  stof(vec[9])
95  );
96 }
static std::string trimStr(std::string str)
Definition: BTagEntry.cc:250
vector< ParameterSet > Parameters
std::string formula
Definition: BTagEntry.h:78
#define end
Definition: vmac.h:39
OperatingPoint
Definition: BTagEntry.h:27
Parameters params
Definition: BTagEntry.h:79
#define begin
Definition: vmac.h:32
BTagEntry::BTagEntry ( const std::string &  func,
BTagEntry::Parameters  p 
)

Definition at line 98 of file BTagEntry.cc.

References Exception, connectstrParser::f1, formula, and jets_cff::func.

98  :
99  formula(func),
100  params(p)
101 {
102  TF1 f1("", formula.c_str()); // compile formula to check validity
103  if (f1.IsZombie()) {
104  throw cms::Exception("BTagCalibration")
105  << "Invalid func string; formula does not compile: "
106  << func;
107  }
108 }
std::string formula
Definition: BTagEntry.h:78
Parameters params
Definition: BTagEntry.h:79
BTagEntry::BTagEntry ( const TF1 *  func,
BTagEntry::Parameters  p 
)

Definition at line 110 of file BTagEntry.cc.

References Exception.

110  :
111  formula(std::string(func->GetExpFormula("p").Data())),
112  params(p)
113 {
114  if (func->IsZombie()) {
115  throw cms::Exception("BTagCalibration")
116  << "Invalid TF1 function; function is zombie: "
117  << func->GetName();
118  }
119 }
std::string formula
Definition: BTagEntry.h:78
Parameters params
Definition: BTagEntry.h:79
BTagEntry::BTagEntry ( const TH1 *  histo,
BTagEntry::Parameters  p 
)

Definition at line 185 of file BTagEntry.cc.

References BTagEntry::Parameters::discrMax, BTagEntry::Parameters::discrMin, Exception, connectstrParser::f1, formula, pileupCalc::nbins, OP_RESHAPING, BTagEntry::Parameters::operatingPoint, params, BTagEntry::Parameters::ptMax, BTagEntry::Parameters::ptMin, th1ToFormulaBinTree(), and th1ToFormulaLin().

185  :
186  params(p)
187 {
188  int nbins = hist->GetNbinsX();
189  TAxis const* axis = hist->GetXaxis();
190 
191  // overwrite bounds with histo values
193  params.discrMin = axis->GetBinLowEdge(1);
194  params.discrMax = axis->GetBinUpEdge(nbins);
195  } else {
196  params.ptMin = axis->GetBinLowEdge(1);
197  params.ptMax = axis->GetBinUpEdge(nbins);
198  }
199 
200  // balanced full binary tree height = ceil(log(2*n_leaves)/log(2))
201  // breakes even around 10, but lower values are more propable in pt-spectrum
202  if (nbins < 15) {
204  } else {
206  }
207 
208  // compile formula to check validity
209  TF1 f1("", formula.c_str());
210  if (f1.IsZombie()) {
211  throw cms::Exception("BTagCalibration")
212  << "Invalid histogram; formula does not compile (>150 bins?): "
213  << hist->GetName();
214  }
215 }
std::string th1ToFormulaBinTree(const TH1 *hist, int start=0, int end=-1)
Definition: BTagEntry.cc:144
std::string th1ToFormulaLin(const TH1 *hist)
Definition: BTagEntry.cc:124
OperatingPoint operatingPoint
Definition: BTagEntry.h:39
std::string formula
Definition: BTagEntry.h:78
Parameters params
Definition: BTagEntry.h:79
BTagEntry::~BTagEntry ( )
inline

Definition at line 72 of file BTagEntry.h.

References makeCSVHeader(), makeCSVLine(), str, AlCaHLTBitMon_QueryRunRegistry::string, and trimStr().

72 {}

Member Function Documentation

std::string BTagEntry::makeCSVHeader ( )
static

Definition at line 217 of file BTagEntry.cc.

Referenced by BTagCalibration::makeCSV(), and ~BTagEntry().

218 {
219  return "OperatingPoint, "
220  "measurementType, "
221  "sysType, "
222  "jetFlavor, "
223  "etaMin, "
224  "etaMax, "
225  "ptMin, "
226  "ptMax, "
227  "discrMin, "
228  "discrMax, "
229  "formula \n";
230 }
std::string BTagEntry::makeCSVLine ( ) const

Definition at line 232 of file BTagEntry.cc.

References BTagEntry::Parameters::discrMax, BTagEntry::Parameters::discrMin, BTagEntry::Parameters::etaMax, BTagEntry::Parameters::etaMin, formula, BTagEntry::Parameters::jetFlavor, BTagEntry::Parameters::measurementType, BTagEntry::Parameters::operatingPoint, params, BTagEntry::Parameters::ptMax, BTagEntry::Parameters::ptMin, and BTagEntry::Parameters::sysType.

Referenced by ~BTagEntry().

233 {
234  std::stringstream buff;
235  buff << params.operatingPoint
236  << ", " << params.measurementType
237  << ", " << params.sysType
238  << ", " << params.jetFlavor
239  << ", " << params.etaMin
240  << ", " << params.etaMax
241  << ", " << params.ptMin
242  << ", " << params.ptMax
243  << ", " << params.discrMin
244  << ", " << params.discrMax
245  << ", \"" << formula
246  << "\" \n";
247  return buff.str();
248 }
std::string sysType
Definition: BTagEntry.h:41
JetFlavor jetFlavor
Definition: BTagEntry.h:42
OperatingPoint operatingPoint
Definition: BTagEntry.h:39
std::string formula
Definition: BTagEntry.h:78
Parameters params
Definition: BTagEntry.h:79
std::string measurementType
Definition: BTagEntry.h:40
template<class Archive >
void BTagEntry::serialize ( Archive &  ar,
const unsigned int  version 
)
private
std::string BTagEntry::trimStr ( std::string  str)
static

Definition at line 250 of file BTagEntry.cc.

References MillePedeFileConverter_cfg::e, and alignCSCRings::s.

Referenced by BTagEntry(), BTagCalibration::readCSV(), and ~BTagEntry().

250  {
251  size_t s = str.find_first_not_of(" \n\r\t");
252  size_t e = str.find_last_not_of (" \n\r\t");
253 
254  if((std::string::npos == s) || (std::string::npos == e))
255  return "";
256  else
257  return str.substr(s, e-s+1);
258 }
#define str(s)

Friends And Related Function Documentation

friend class boost::serialization::access
friend

Definition at line 81 of file BTagEntry.h.

template<typename CondSerializationT , typename Enabled >
friend struct cond::serialization::access
friend

Definition at line 81 of file BTagEntry.h.

Member Data Documentation

std::string BTagEntry::formula

Definition at line 78 of file BTagEntry.h.

Referenced by BTagEntry(), and makeCSVLine().

Parameters BTagEntry::params

Definition at line 79 of file BTagEntry.h.

Referenced by BTagCalibration::addEntry(), BTagEntry(), and makeCSVLine().