CMS 3D CMS Logo

SoftLeptonTagPlotter.cc
Go to the documentation of this file.
4 
5 #include <sstream>
6 #include <string>
7 
8 using namespace std ;
9 using namespace RecoBTag;
10 
11 static const string ordinal[9] = { "1st", "2nd", "3rd", "4th", "5th", "6th", "7th", "8th", "9th" };
12 
14  const EtaPtBin & etaPtBin, const edm::ParameterSet& pSet,
15  unsigned int mc, bool wf, DQMStore::IBooker & ibook) :
16  BaseTagInfoPlotter(tagName, etaPtBin), mcPlots_(mc), willFinalize_(wf)
17 {
18  const std::string softLepDir(theExtensionString.substr(1));
19 
20  if(willFinalize_) return;
21 
22  for (int i = 0; i < s_leptons; i++) {
23  std::string s;
24  s += ordinal[i] + " lepton ";
25  m_leptonId.push_back(std::make_unique<FlavourHistograms<double>>(
26  s + "id",
27  "Lepton identification discriminaint",
28  60, -0.1, 1.1, false, false, true, "b", softLepDir, mcPlots_, ibook));
29  m_leptonPt.push_back(std::make_unique<FlavourHistograms<double>>(
30  s + "pT",
31  "Lepton transverse moementum",
32  100, 0.0, 20.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
33  m_sip2dsig.push_back(std::make_unique<FlavourHistograms<double>>(
34  s + "sip2dsig",
35  "Lepton signed 2D impact parameter significance",
36  100, -20.0, 30.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
37  m_sip3dsig.push_back(std::make_unique<FlavourHistograms<double>>(
38  s + "sip3dsig",
39  "Lepton signed 3D impact parameter significance",
40  100, -20.0, 30.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
41  m_sip2d.push_back(std::make_unique<FlavourHistograms<double>>(
42  s + "sip2d",
43  "Lepton signed 2D impact parameter",
44  100, -20.0, 30.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
45  m_sip3d.push_back(std::make_unique<FlavourHistograms<double>>(
46  s + "sip3d",
47  "Lepton signed 3D impact parameter",
48  100, -20.0, 30.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
49  m_ptRel.push_back(std::make_unique<FlavourHistograms<double>>(
50  s + "pT rel",
51  "Lepton transverse momentum relative to jet axis",
52  100, 0.0, 10.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
53  m_p0Par.push_back(std::make_unique<FlavourHistograms<double>>(
54  s + "p0 par",
55  "Lepton moementum along jet axis in the B rest frame",
56  100, 0.0, 10.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
57  m_etaRel.push_back(std::make_unique<FlavourHistograms<double>>(
58  s + "eta rel",
59  "Lepton pseudorapidity relative to jet axis",
60  100, -5.0, 25.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
61  m_deltaR.push_back(std::make_unique<FlavourHistograms<double>>(
62  s + "delta R",
63  "Lepton pseudoangular distance from jet axis",
64  100, 0.0, 0.6, false, false, true, "b", softLepDir, mcPlots_, ibook));
65  m_ratio.push_back(std::make_unique<FlavourHistograms<double>>(
66  s + "energy ratio",
67  "Ratio of lepton momentum to jet energy",
68  100, 0.0, 2.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
69  m_ratioRel.push_back(std::make_unique<FlavourHistograms<double>>(
70  s + "parallel energy ratio",
71  "Ratio of lepton momentum along the jet axis to jet energy",
72  100, 0.0, 2.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
73  }
74 }
75 
77 
78 void SoftLeptonTagPlotter::analyzeTag( const reco::BaseTagInfo * baseTagInfo, double jec, int jetFlavour, float w/*=1*/)
79 {
80  const reco::SoftLeptonTagInfo * tagInfo =
81  dynamic_cast<const reco::SoftLeptonTagInfo *>(baseTagInfo);
82 
83  if (!tagInfo) {
84  throw cms::Exception("Configuration")
85  << "BTagPerformanceAnalyzer: Extended TagInfo not of type SoftLeptonTagInfo. " << endl;
86  }
87 
88  int n_leptons = tagInfo->leptons();
89 
90  for (int i = 0; i != n_leptons && i != s_leptons; ++i) {
91  const reco::SoftLeptonProperties& properties = tagInfo->properties(i);
92  m_leptonPt[i]->fill(jetFlavour, tagInfo->lepton(i)->pt(), w);
93  m_leptonId[i]->fill(jetFlavour, properties.quality(), w);
94  m_sip2dsig[i]->fill(jetFlavour, properties.sip2dsig, w);
95  m_sip3dsig[i]->fill(jetFlavour, properties.sip3dsig, w);
96  m_sip2d[i]->fill(jetFlavour, properties.sip2d, w);
97  m_sip3d[i]->fill(jetFlavour, properties.sip3d, w);
98  m_ptRel[i]->fill(jetFlavour, properties.ptRel, w);
99  m_p0Par[i]->fill(jetFlavour, properties.p0Par, w);
100  m_etaRel[i]->fill(jetFlavour, properties.etaRel, w);
101  m_deltaR[i]->fill(jetFlavour, properties.deltaR, w);
102  m_ratio[i]->fill(jetFlavour, properties.ratio, w);
103  m_ratioRel[i]->fill(jetFlavour, properties.ratioRel, w);
104  }
105 }
106 
108 {
109  if(willFinalize_) return;
110 
111  const std::string cName("SoftLeptonPlots" + theExtensionString);
112  setTDRStyle()->cd();
113  TCanvas canvas(cName.c_str(), cName.c_str(), 600, 900);
114  canvas.UseCurrentStyle();
115  canvas.Divide(2,3);
116  canvas.Print((name + cName + ".ps[").c_str());
117  for (int i = 0; i < s_leptons; i++) {
118  canvas.cd(1)->Clear();
119  m_leptonId[i]->plot();
120  canvas.cd(2)->Clear();
121  m_leptonPt[i]->plot();
122  canvas.cd(3)->Clear();
123  m_sip2d[i]->plot();
124  canvas.cd(4)->Clear();
125  m_sip3d[i]->plot();
126  canvas.cd(5)->Clear();
127  m_sip2dsig[i]->plot();
128  canvas.cd(6)->Clear();
129  m_sip3dsig[i]->plot();
130  canvas.Print((name + cName + ".ps").c_str());
131 
132  canvas.cd(1)->Clear();
133  m_etaRel[i]->plot();
134  canvas.cd(2)->Clear();
135  m_deltaR[i]->plot();
136  canvas.cd(3)->Clear();
137  m_ratio[i]->plot();
138  canvas.cd(4)->Clear();
139  m_ratioRel[i]->plot();
140  canvas.cd(5)->Clear();
141  m_ptRel[i]->plot();
142  canvas.cd(6)->Clear();
143  m_p0Par[i]->plot();
144  canvas.Print((name + cName + ".ps").c_str());
145  }
146  canvas.Print((name + cName + ".ps]").c_str());
147 }
148 
149 
151 {
152  if(willFinalize_) return;
153 
154  for (int i=0; i != s_leptons; ++i) {
155  m_leptonId[i]->epsPlot( name );
156  m_leptonPt[i]->epsPlot( name );
157  m_sip2d[i]->epsPlot( name );
158  m_sip3d[i]->epsPlot( name );
159  m_sip2dsig[i]->epsPlot( name );
160  m_sip3dsig[i]->epsPlot( name );
161  m_ptRel[i]->epsPlot( name );
162  m_p0Par[i]->epsPlot( name );
163  m_etaRel[i]->epsPlot( name );
164  m_deltaR[i]->epsPlot( name );
165  m_ratio[i]->epsPlot( name );
166  m_ratioRel[i]->epsPlot( name );
167  }
168 }
169 
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_deltaR
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_sip2d
void psPlot(const std::string &name) override
float quality(Quality::Generic qual, bool throwIfUndefined=true) const
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_ptRel
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_ratio
const double w
Definition: UKUtility.cc:23
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_sip3d
~SoftLeptonTagPlotter(void) override
const std::string theExtensionString
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_ratioRel
void analyzeTag(const reco::BaseTagInfo *baseTagInfo, double jec, int jetFlavour, float w) override
const SoftLeptonProperties & properties(size_t i) const
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_p0Par
SoftLeptonTagPlotter(const std::string &tagName, const EtaPtBin &etaPtBin, const edm::ParameterSet &pSet, unsigned int mc, bool willFinalize, DQMStore::IBooker &ibook)
static const string ordinal[9]
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_etaRel
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_sip2dsig
def setTDRStyle()
Definition: plotscripts.py:89
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_leptonId
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_leptonPt
Definition: Tools.h:23
static const int s_leptons
std::vector< std::unique_ptr< FlavourHistograms< double > > > m_sip3dsig
def canvas(sub, attr)
Definition: svgfig.py:482
void epsPlot(const std::string &name) override