Skip to content

Commit

Permalink
Merge pull request #6378 from HeinAtCERN/btag-cond-format
Browse files Browse the repository at this point in the history
Btag cond format
  • Loading branch information
cmsbuild committed Nov 20, 2014
2 parents 573a490 + cf29be7 commit eccca54
Show file tree
Hide file tree
Showing 28 changed files with 1,871 additions and 15 deletions.
4 changes: 4 additions & 0 deletions CondCore/BTauPlugins/src/plugin.cc
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,9 @@
#include "CondFormats/DataRecord/interface/BTauGenericMVAJetTagComputerRcd.h"
#include "CondFormats/DataRecord/interface/TauTagMVAComputerRcd.h"

#include "CondFormats/BTauObjects/interface/BTagCalibration.h"
#include "CondFormats/DataRecord/interface/BTagCalibrationRcd.h"

using namespace PhysicsTools::Calibration;


Expand All @@ -19,3 +22,4 @@ REGISTER_PLUGIN(BTauGenericMVAJetTagComputerRcd, MVAComputerContainer);
REGISTER_PLUGIN(TauTagMVAComputerRcd, MVAComputerContainer);
REGISTER_PLUGIN(BTagTrackProbability2DRcd,TrackProbabilityCalibration);
REGISTER_PLUGIN(BTagTrackProbability3DRcd,TrackProbabilityCalibration);
REGISTER_PLUGIN(BTagCalibrationRcd, BTagCalibration);
54 changes: 54 additions & 0 deletions CondFormats/BTauObjects/interface/BTagCalibration.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,54 @@
#ifndef BTagCalibration_H
#define BTagCalibration_H

/**
* BTagCalibration
*
* The 'hierarchy' of stored information is this:
* - by tagger (BTagCalibration)
* - by operating point or reshape bin
* - by jet parton flavor
* - by type of measurement
* - by systematic
* - by eta bin
* - as 1D-function dependent of pt or discriminant
*
************************************************************/

#include <map>
#include <vector>
#include <string>
#include <istream>
#include <ostream>

#include "CondFormats/Serialization/interface/Serializable.h"
#include "CondFormats/BTauObjects/interface/BTagEntry.h"

class BTagCalibration
{
public:
BTagCalibration() {}
BTagCalibration(const std::string &tagger);
BTagCalibration(const std::string &tagger, const std::string &filename);
~BTagCalibration() {}

std::string tagger() const {return tagger_;}

void addEntry(const BTagEntry &entry);
const std::vector<BTagEntry>& getEntries(const BTagEntry::Parameters &par) const;

void readCSV(istream &s);
void readCSV(const std::string &s);
void makeCSV(ostream &s) const;
std::string makeCSV() const;

protected:
static std::string token(const BTagEntry::Parameters &par);

std::string tagger_;
std::map<std::string, std::vector<BTagEntry> > data_;

COND_SERIALIZABLE;
};

#endif // BTagCalibration_H
52 changes: 52 additions & 0 deletions CondFormats/BTauObjects/interface/BTagCalibrationReader.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,52 @@
#ifndef BTagCalibrationReader_H
#define BTagCalibrationReader_H

/**
* BTagCalibrationReader
*
* Helper class to pull out a specific set of BTagEntry's out of a
* BTagCalibration. TF1 functions are set up at initialization time.
*
************************************************************/

#include <map>
#include <string>
#include <vector>
#include <TF1.h>

#include "CondFormats/BTauObjects/interface/BTagEntry.h"
#include "CondFormats/BTauObjects/interface/BTagCalibration.h"

class BTagCalibrationReader
{
public:
BTagCalibrationReader() {}
BTagCalibrationReader(const BTagCalibration* c,
BTagEntry::OperatingPoint op,
std::string measurementType="comb",
std::string sysType="central");
~BTagCalibrationReader() {}

double eval(BTagEntry::JetFlavor jf,
float eta,
float pt,
float discr=0.) const;

protected:
struct TmpEntry {
float etaMin;
float etaMax;
float ptMin;
float ptMax;
float discrMin;
float discrMax;
TF1 func;
};
void setupTmpData(const BTagCalibration* c);

BTagEntry::Parameters params;
std::map<BTagEntry::JetFlavor, std::vector<TmpEntry> > tmpData_;
std::vector<bool> useAbsEta;
};

#endif // BTagCalibrationReader_H
84 changes: 84 additions & 0 deletions CondFormats/BTauObjects/interface/BTagEntry.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,84 @@
#ifndef BTagEntry_H
#define BTagEntry_H

/**
*
* 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.
*
************************************************************/

#include <string>
#include <TF1.h>
#include <TH1.h>

#include "CondFormats/Serialization/interface/Serializable.h"

class BTagEntry
{
public:
enum OperatingPoint {
OP_LOOSE=0,
OP_MEDIUM=1,
OP_TIGHT=2,
OP_RESHAPING=3,
};
enum JetFlavor {
FLAV_B=0,
FLAV_C=1,
FLAV_UDSG=2,
};
struct Parameters {
OperatingPoint operatingPoint;
std::string measurementType;
std::string sysType;
JetFlavor jetFlavor;
float etaMin;
float etaMax;
float ptMin;
float ptMax;
float discrMin;
float discrMax;

// default constructor
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.
);

COND_SERIALIZABLE;
};

BTagEntry() {}
BTagEntry(const std::string &csvLine);
BTagEntry(const std::string &func, Parameters p);
BTagEntry(const TF1* func, Parameters p);
BTagEntry(const TH1* histo, Parameters p);
~BTagEntry() {}
static std::string makeCSVHeader();
std::string makeCSVLine() const;
static std::string trimStr(std::string str);

// public, no getters needed
std::string formula;
Parameters params;

COND_SERIALIZABLE;
};

#endif // BTagEntry_H
88 changes: 88 additions & 0 deletions CondFormats/BTauObjects/src/BTagCalibration.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,88 @@
#include <fstream>
#include <sstream>

#include "CondFormats/BTauObjects/interface/BTagCalibration.h"
#include "FWCore/Utilities/interface/Exception.h"


BTagCalibration::BTagCalibration(const std::string &taggr):
tagger_(taggr)
{}

BTagCalibration::BTagCalibration(const std::string &taggr,
const std::string &filename):
tagger_(taggr)
{
std::ifstream ifs(filename);
readCSV(ifs);
ifs.close();
}

void BTagCalibration::addEntry(const BTagEntry &entry)
{
data_[token(entry.params)].push_back(entry);
}

const std::vector<BTagEntry>& BTagCalibration::getEntries(
const BTagEntry::Parameters &par) const
{
auto tok = token(par);
if (!data_.count(tok)) {
throw cms::Exception("BTagCalibration")
<< "(OperatingPoint, measurementType, sysType) not available: "
<< tok;
}
return data_.at(tok);
}

void BTagCalibration::readCSV(const std::string &s)
{
std::stringstream buff(s);
readCSV(buff);
}

void BTagCalibration::readCSV(istream &s)
{
std::string line;

// firstline might be the header
getline(s,line);
if (line.find("OperatingPoint") == std::string::npos) {
addEntry(BTagEntry(line));
}

while (getline(s,line)) {
line = BTagEntry::trimStr(line);
if (line.empty()) { // skip empty lines
continue;
}
addEntry(BTagEntry(line));
}
}

void BTagCalibration::makeCSV(ostream &s) const
{
s << BTagEntry::makeCSVHeader();
for (auto i = data_.cbegin(); i != data_.cend(); ++i) {
auto vec = i->second;
for (auto j = vec.cbegin(); j != vec.cend(); ++j) {
s << j->makeCSVLine();
}
}
}

std::string BTagCalibration::makeCSV() const
{
std::stringstream buff;
makeCSV(buff);
return buff.str();
}

std::string BTagCalibration::token(const BTagEntry::Parameters &par)
{
std::stringstream buff;
buff << par.operatingPoint << ", "
<< par.measurementType << ", "
<< par.sysType;
return buff.str();
}
72 changes: 72 additions & 0 deletions CondFormats/BTauObjects/src/BTagCalibrationReader.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,72 @@
#include "CondFormats/BTauObjects/interface/BTagCalibrationReader.h"

BTagCalibrationReader::BTagCalibrationReader(const BTagCalibration* c,
BTagEntry::OperatingPoint op,
std::string measurementType,
std::string sysType):
params(BTagEntry::Parameters(op, measurementType, sysType)),
useAbsEta(true)
{
setupTmpData(c);
}

double BTagCalibrationReader::eval(BTagEntry::JetFlavor jf,
float eta,
float pt,
float discr) const
{
bool use_discr = (params.operatingPoint == BTagEntry::OP_RESHAPING);
if (useAbsEta[jf] && eta < 0) {
eta = -eta;
}

// search linearly through eta, pt and discr ranges and eval
// future: find some clever data structure based on intervals
const auto &entries = tmpData_.at(jf);
for (unsigned i=0; i<entries.size(); ++i) {
const BTagCalibrationReader::TmpEntry &e = entries.at(i);
if (
e.etaMin <= eta && eta < e.etaMax // find eta
&& e.ptMin <= pt && pt < e.ptMax // check pt
){
if (use_discr) { // discr. reshaping?
if (e.discrMin <= discr && discr < e.discrMax) { // check discr
return e.func.Eval(discr);
}
} else {
return e.func.Eval(pt);
}
}
}

return 0.; // default value
}

void BTagCalibrationReader::setupTmpData(const BTagCalibration* c)
{
useAbsEta = std::vector<bool>(4, true);
const auto &entries = c->getEntries(params);
for (unsigned i=0; i<entries.size(); ++i) {
const BTagEntry &be = entries[i];
BTagCalibrationReader::TmpEntry te;
te.etaMin = be.params.etaMin;
te.etaMax = be.params.etaMax;
te.ptMin = be.params.ptMin;
te.ptMax = be.params.ptMax;
te.discrMin = be.params.discrMin;
te.discrMax = be.params.discrMax;

if (params.operatingPoint == BTagEntry::OP_RESHAPING) {
te.func = TF1("", be.formula.c_str(),
be.params.discrMin, be.params.discrMax);
} else {
te.func = TF1("", be.formula.c_str(),
be.params.ptMin, be.params.ptMax);
}

tmpData_[be.params.jetFlavor].push_back(te);
if (te.etaMin < 0) {
useAbsEta[be.params.jetFlavor] = false;
}
}
}

0 comments on commit eccca54

Please sign in to comment.