CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BTagEntry.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <sstream>
3 
6 
9  std::string measurement_type,
10  std::string sys_type,
11  JetFlavor jf,
12  float eta_min,
13  float eta_max,
14  float pt_min,
15  float pt_max,
16  float discr_min,
17  float discr_max
18 ):
19  operatingPoint(op),
20  measurementType(measurement_type),
21  sysType(sys_type),
22  jetFlavor(jf),
23  etaMin(eta_min),
24  etaMax(eta_max),
25  ptMin(pt_min),
26  ptMax(pt_max),
27  discrMin(discr_min),
28  discrMax(discr_max)
29 {
31  measurementType.begin(), ::tolower);
32  std::transform(sysType.begin(), sysType.end(),
33  sysType.begin(), ::tolower);
34 }
35 
37 {
38  // make tokens
39  std::stringstream buff(csvLine);
40  std::vector<std::string> vec;
41  std::string token;
42  while (std::getline(buff, token, ","[0])) {
43  token = BTagEntry::trimStr(token);
44  if (token.empty()) {
45  continue;
46  }
47  vec.push_back(token);
48  }
49  if (vec.size() != 11) {
50  throw cms::Exception("BTagCalibration")
51  << "Invalid csv line; num tokens != 11: "
52  << csvLine;
53  }
54 
55  // clean string values
56  char chars[] = " \"\n";
57  for (unsigned int i = 0; i < strlen(chars); ++i) {
58  vec[1].erase(remove(vec[1].begin(),vec[1].end(),chars[i]),vec[1].end());
59  vec[2].erase(remove(vec[2].begin(),vec[2].end(),chars[i]),vec[2].end());
60  vec[10].erase(remove(vec[10].begin(),vec[10].end(),chars[i]),vec[10].end());
61  }
62 
63  // make formula
64  formula = vec[10];
65  TF1 f1("", formula.c_str()); // compile formula to check validity
66  if (f1.IsZombie()) {
67  throw cms::Exception("BTagCalibration")
68  << "Invalid csv line; formula does not compile: "
69  << csvLine;
70  }
71 
72  // make parameters
73  unsigned op = stoi(vec[0]);
74  if (op > 3) {
75  throw cms::Exception("BTagCalibration")
76  << "Invalid csv line; OperatingPoint > 3: "
77  << csvLine;
78  }
79  unsigned jf = stoi(vec[3]);
80  if (jf > 2) {
81  throw cms::Exception("BTagCalibration")
82  << "Invalid csv line; JetFlavor > 2: "
83  << csvLine;
84  }
87  vec[1],
88  vec[2],
90  stof(vec[4]),
91  stof(vec[5]),
92  stof(vec[6]),
93  stof(vec[7]),
94  stof(vec[8]),
95  stof(vec[9])
96  );
97 }
98 
100  formula(func),
101  params(p)
102 {
103  TF1 f1("", formula.c_str()); // compile formula to check validity
104  if (f1.IsZombie()) {
105  throw cms::Exception("BTagCalibration")
106  << "Invalid func string; formula does not compile: "
107  << func;
108  }
109 }
110 
112  formula(std::string(func->GetExpFormula("p").Data())),
113  params(p)
114 {
115  if (func->IsZombie()) {
116  throw cms::Exception("BTagCalibration")
117  << "Invalid TF1 function; function is zombie: "
118  << func->GetName();
119  }
120 }
121 
122 // Creates chained step functions like this:
123 // "<prevous_bin> : x<bin_high_bound ? bin_value : <next_bin>"
124 // e.g. "x<0 ? 1 : x<1 ? 2 : x<2 ? 3 : 4"
126  int nbins = hist->GetNbinsX();
127  TAxis const* axis = hist->GetXaxis();
128  std::stringstream buff;
129  buff << "x<" << axis->GetBinLowEdge(1) << " ? 0. : "; // default value
130  for (int i=1; i<nbins+1; ++i) {
131  char tmp_buff[50];
132  sprintf(tmp_buff,
133  "x<%g ? %g : ", // %g is the smaller one of %e or %f
134  axis->GetBinUpEdge(i),
135  hist->GetBinContent(i));
136  buff << tmp_buff;
137  }
138  buff << 0.; // default value
139  return buff.str();
140 }
141 
142 // Creates step functions making a binary search tree:
143 // "x<mid_bin_bound ? (<left side tree>) : (<right side tree>)"
144 // e.g. "x<2 ? (x<1 ? (x<0 ? 0:0.1) : (1)) : (x<4 ? (x<3 ? 2:3) : (0))"
145 std::string th1ToFormulaBinTree(const TH1* hist, int start=0, int end=-1) {
146  if (end == -1) { // initialize
147  start = 0.;
148  end = hist->GetNbinsX()+1;
149  TH1* h2 = (TH1*) hist->Clone();
150  h2->SetBinContent(start, 0); // kill underflow
151  h2->SetBinContent(end, 0); // kill overflow
153  delete h2;
154  return res;
155  }
156  if (start == end) { // leave is reached
157  char tmp_buff[20];
158  sprintf(tmp_buff, "%g", hist->GetBinContent(start));
159  return std::string(tmp_buff);
160  }
161  if (start == end - 1) { // no parenthesis for neighbors
162  char tmp_buff[70];
163  sprintf(tmp_buff,
164  "x<%g ? %g:%g",
165  hist->GetXaxis()->GetBinUpEdge(start),
166  hist->GetBinContent(start),
167  hist->GetBinContent(end));
168  return std::string(tmp_buff);
169  }
170 
171  // top-down recursion
172  std::stringstream buff;
173  int mid = (end-start)/2 + start;
174  char tmp_buff[25];
175  sprintf(tmp_buff,
176  "x<%g ? (",
177  hist->GetXaxis()->GetBinUpEdge(mid));
178  buff << tmp_buff
179  << th1ToFormulaBinTree(hist, start, mid)
180  << ") : ("
181  << th1ToFormulaBinTree(hist, mid+1, end)
182  << ")";
183  return buff.str();
184 }
185 
187  params(p)
188 {
189  int nbins = hist->GetNbinsX();
190  TAxis const* axis = hist->GetXaxis();
191 
192  // overwrite bounds with histo values
194  params.discrMin = axis->GetBinLowEdge(1);
195  params.discrMax = axis->GetBinUpEdge(nbins);
196  } else {
197  params.ptMin = axis->GetBinLowEdge(1);
198  params.ptMax = axis->GetBinUpEdge(nbins);
199  }
200 
201  // balanced full binary tree height = ceil(log(2*n_leaves)/log(2))
202  // breakes even around 10, but lower values are more propable in pt-spectrum
203  if (nbins < 15) {
204  formula = th1ToFormulaLin(hist);
205  } else {
207  }
208 
209  // compile formula to check validity
210  TF1 f1("", formula.c_str());
211  if (f1.IsZombie()) {
212  throw cms::Exception("BTagCalibration")
213  << "Invalid histogram; formula does not compile (>150 bins?): "
214  << hist->GetName();
215  }
216 }
217 
219 {
220  return "OperatingPoint, "
221  "measurementType, "
222  "sysType, "
223  "jetFlavor, "
224  "etaMin, "
225  "etaMax, "
226  "ptMin, "
227  "ptMax, "
228  "discrMin, "
229  "discrMax, "
230  "formula \n";
231 }
232 
234 {
235  std::stringstream buff;
236  buff << params.operatingPoint
237  << ", " << params.measurementType
238  << ", " << params.sysType
239  << ", " << params.jetFlavor
240  << ", " << params.etaMin
241  << ", " << params.etaMax
242  << ", " << params.ptMin
243  << ", " << params.ptMax
244  << ", " << params.discrMin
245  << ", " << params.discrMax
246  << ", \"" << formula
247  << "\" \n";
248  return buff.str();
249 }
250 
252  size_t s = str.find_first_not_of(" \n\r\t");
253  size_t e = str.find_last_not_of (" \n\r\t");
254 
255  if((std::string::npos == s) || (std::string::npos == e))
256  return "";
257  else
258  return str.substr(s, e-s+1);
259 }
std::string sysType
Definition: BTagEntry.h:41
int i
Definition: DBlmapReader.cc:9
Parameters(OperatingPoint op=OP_TIGHT, std::string measurement_type="comb", std::string sys_type="central", JetFlavor jf=FLAV_B, float eta_min=-99999., float eta_max=99999., float pt_min=0., float pt_max=99999., float discr_min=0., float discr_max=99999.)
Definition: BTagEntry.cc:7
std::string th1ToFormulaBinTree(const TH1 *hist, int start=0, int end=-1)
Definition: BTagEntry.cc:145
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
static std::string trimStr(std::string str)
Definition: BTagEntry.cc:251
JetFlavor jetFlavor
Definition: BTagEntry.h:42
static std::string makeCSVHeader()
Definition: BTagEntry.cc:218
std::string th1ToFormulaLin(const TH1 *hist)
Definition: BTagEntry.cc:125
std::string makeCSVLine() const
Definition: BTagEntry.cc:233
vector< ParameterSet > Parameters
OperatingPoint operatingPoint
Definition: BTagEntry.h:39
std::string formula
Definition: BTagEntry.h:78
#define end
Definition: vmac.h:37
OperatingPoint
Definition: BTagEntry.h:27
Parameters params
Definition: BTagEntry.h:79
BTagEntry()
Definition: BTagEntry.h:67
#define begin
Definition: vmac.h:30
std::string measurementType
Definition: BTagEntry.h:40