Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:48

0001 #include "DQMOffline/RecoB/interface/SoftLeptonTagPlotter.h"
0002 #include "DataFormats/TrackReco/interface/Track.h"
0003 #include "DQMOffline/RecoB/interface/Tools.h"
0004 
0005 #include <sstream>
0006 #include <string>
0007 
0008 using namespace std;
0009 using namespace RecoBTag;
0010 
0011 static const string ordinal[9] = {"1st", "2nd", "3rd", "4th", "5th", "6th", "7th", "8th", "9th"};
0012 
0013 SoftLeptonTagPlotter::SoftLeptonTagPlotter(const std::string& tagName,
0014                                            const EtaPtBin& etaPtBin,
0015                                            const edm::ParameterSet& pSet,
0016                                            unsigned int mc,
0017                                            bool wf,
0018                                            DQMStore::IBooker& ibook)
0019     : BaseTagInfoPlotter(tagName, etaPtBin), mcPlots_(mc), willFinalize_(wf) {
0020   const std::string softLepDir(theExtensionString.substr(1));
0021 
0022   if (willFinalize_)
0023     return;
0024 
0025   for (int i = 0; i < s_leptons; i++) {
0026     std::string s;
0027     s += ordinal[i] + " lepton ";
0028     m_leptonId.push_back(std::make_unique<FlavourHistograms<double>>(s + "id",
0029                                                                      "Lepton identification discriminaint",
0030                                                                      60,
0031                                                                      -0.1,
0032                                                                      1.1,
0033                                                                      false,
0034                                                                      false,
0035                                                                      true,
0036                                                                      "b",
0037                                                                      softLepDir,
0038                                                                      mcPlots_,
0039                                                                      ibook));
0040     m_leptonPt.push_back(std::make_unique<FlavourHistograms<double>>(
0041         s + "pT", "Lepton transverse moementum", 100, 0.0, 20.0, false, false, true, "b", softLepDir, mcPlots_, ibook));
0042     m_sip2dsig.push_back(std::make_unique<FlavourHistograms<double>>(s + "sip2dsig",
0043                                                                      "Lepton signed 2D impact parameter significance",
0044                                                                      100,
0045                                                                      -20.0,
0046                                                                      30.0,
0047                                                                      false,
0048                                                                      false,
0049                                                                      true,
0050                                                                      "b",
0051                                                                      softLepDir,
0052                                                                      mcPlots_,
0053                                                                      ibook));
0054     m_sip3dsig.push_back(std::make_unique<FlavourHistograms<double>>(s + "sip3dsig",
0055                                                                      "Lepton signed 3D impact parameter significance",
0056                                                                      100,
0057                                                                      -20.0,
0058                                                                      30.0,
0059                                                                      false,
0060                                                                      false,
0061                                                                      true,
0062                                                                      "b",
0063                                                                      softLepDir,
0064                                                                      mcPlots_,
0065                                                                      ibook));
0066     m_sip2d.push_back(std::make_unique<FlavourHistograms<double>>(s + "sip2d",
0067                                                                   "Lepton signed 2D impact parameter",
0068                                                                   100,
0069                                                                   -20.0,
0070                                                                   30.0,
0071                                                                   false,
0072                                                                   false,
0073                                                                   true,
0074                                                                   "b",
0075                                                                   softLepDir,
0076                                                                   mcPlots_,
0077                                                                   ibook));
0078     m_sip3d.push_back(std::make_unique<FlavourHistograms<double>>(s + "sip3d",
0079                                                                   "Lepton signed 3D impact parameter",
0080                                                                   100,
0081                                                                   -20.0,
0082                                                                   30.0,
0083                                                                   false,
0084                                                                   false,
0085                                                                   true,
0086                                                                   "b",
0087                                                                   softLepDir,
0088                                                                   mcPlots_,
0089                                                                   ibook));
0090     m_ptRel.push_back(std::make_unique<FlavourHistograms<double>>(s + "pT rel",
0091                                                                   "Lepton transverse momentum relative to jet axis",
0092                                                                   100,
0093                                                                   0.0,
0094                                                                   10.0,
0095                                                                   false,
0096                                                                   false,
0097                                                                   true,
0098                                                                   "b",
0099                                                                   softLepDir,
0100                                                                   mcPlots_,
0101                                                                   ibook));
0102     m_p0Par.push_back(std::make_unique<FlavourHistograms<double>>(s + "p0 par",
0103                                                                   "Lepton moementum along jet axis in the B rest frame",
0104                                                                   100,
0105                                                                   0.0,
0106                                                                   10.0,
0107                                                                   false,
0108                                                                   false,
0109                                                                   true,
0110                                                                   "b",
0111                                                                   softLepDir,
0112                                                                   mcPlots_,
0113                                                                   ibook));
0114     m_etaRel.push_back(std::make_unique<FlavourHistograms<double>>(s + "eta rel",
0115                                                                    "Lepton pseudorapidity relative to jet axis",
0116                                                                    100,
0117                                                                    -5.0,
0118                                                                    25.0,
0119                                                                    false,
0120                                                                    false,
0121                                                                    true,
0122                                                                    "b",
0123                                                                    softLepDir,
0124                                                                    mcPlots_,
0125                                                                    ibook));
0126     m_deltaR.push_back(std::make_unique<FlavourHistograms<double>>(s + "delta R",
0127                                                                    "Lepton pseudoangular distance from jet axis",
0128                                                                    100,
0129                                                                    0.0,
0130                                                                    0.6,
0131                                                                    false,
0132                                                                    false,
0133                                                                    true,
0134                                                                    "b",
0135                                                                    softLepDir,
0136                                                                    mcPlots_,
0137                                                                    ibook));
0138     m_ratio.push_back(std::make_unique<FlavourHistograms<double>>(s + "energy ratio",
0139                                                                   "Ratio of lepton momentum to jet energy",
0140                                                                   100,
0141                                                                   0.0,
0142                                                                   2.0,
0143                                                                   false,
0144                                                                   false,
0145                                                                   true,
0146                                                                   "b",
0147                                                                   softLepDir,
0148                                                                   mcPlots_,
0149                                                                   ibook));
0150     m_ratioRel.push_back(
0151         std::make_unique<FlavourHistograms<double>>(s + "parallel energy ratio",
0152                                                     "Ratio of lepton momentum along the jet axis to jet energy",
0153                                                     100,
0154                                                     0.0,
0155                                                     2.0,
0156                                                     false,
0157                                                     false,
0158                                                     true,
0159                                                     "b",
0160                                                     softLepDir,
0161                                                     mcPlots_,
0162                                                     ibook));
0163   }
0164 }
0165 
0166 SoftLeptonTagPlotter::~SoftLeptonTagPlotter() {}
0167 
0168 void SoftLeptonTagPlotter::analyzeTag(const reco::BaseTagInfo* baseTagInfo, double jec, int jetFlavour, float w /*=1*/) {
0169   const reco::SoftLeptonTagInfo* tagInfo = dynamic_cast<const reco::SoftLeptonTagInfo*>(baseTagInfo);
0170 
0171   if (!tagInfo) {
0172     throw cms::Exception("Configuration")
0173         << "BTagPerformanceAnalyzer: Extended TagInfo not of type SoftLeptonTagInfo. " << endl;
0174   }
0175 
0176   int n_leptons = tagInfo->leptons();
0177 
0178   for (int i = 0; i != n_leptons && i != s_leptons; ++i) {
0179     const reco::SoftLeptonProperties& properties = tagInfo->properties(i);
0180     m_leptonPt[i]->fill(jetFlavour, tagInfo->lepton(i)->pt(), w);
0181     m_leptonId[i]->fill(jetFlavour, properties.quality(), w);
0182     m_sip2dsig[i]->fill(jetFlavour, properties.sip2dsig, w);
0183     m_sip3dsig[i]->fill(jetFlavour, properties.sip3dsig, w);
0184     m_sip2d[i]->fill(jetFlavour, properties.sip2d, w);
0185     m_sip3d[i]->fill(jetFlavour, properties.sip3d, w);
0186     m_ptRel[i]->fill(jetFlavour, properties.ptRel, w);
0187     m_p0Par[i]->fill(jetFlavour, properties.p0Par, w);
0188     m_etaRel[i]->fill(jetFlavour, properties.etaRel, w);
0189     m_deltaR[i]->fill(jetFlavour, properties.deltaR, w);
0190     m_ratio[i]->fill(jetFlavour, properties.ratio, w);
0191     m_ratioRel[i]->fill(jetFlavour, properties.ratioRel, w);
0192   }
0193 }
0194 
0195 void SoftLeptonTagPlotter::psPlot(const std::string& name) {
0196   if (willFinalize_)
0197     return;
0198 
0199   const std::string cName("SoftLeptonPlots" + theExtensionString);
0200   setTDRStyle()->cd();
0201   TCanvas canvas(cName.c_str(), cName.c_str(), 600, 900);
0202   canvas.UseCurrentStyle();
0203   canvas.Divide(2, 3);
0204   canvas.Print((name + cName + ".ps[").c_str());
0205   for (int i = 0; i < s_leptons; i++) {
0206     canvas.cd(1)->Clear();
0207     m_leptonId[i]->plot();
0208     canvas.cd(2)->Clear();
0209     m_leptonPt[i]->plot();
0210     canvas.cd(3)->Clear();
0211     m_sip2d[i]->plot();
0212     canvas.cd(4)->Clear();
0213     m_sip3d[i]->plot();
0214     canvas.cd(5)->Clear();
0215     m_sip2dsig[i]->plot();
0216     canvas.cd(6)->Clear();
0217     m_sip3dsig[i]->plot();
0218     canvas.Print((name + cName + ".ps").c_str());
0219 
0220     canvas.cd(1)->Clear();
0221     m_etaRel[i]->plot();
0222     canvas.cd(2)->Clear();
0223     m_deltaR[i]->plot();
0224     canvas.cd(3)->Clear();
0225     m_ratio[i]->plot();
0226     canvas.cd(4)->Clear();
0227     m_ratioRel[i]->plot();
0228     canvas.cd(5)->Clear();
0229     m_ptRel[i]->plot();
0230     canvas.cd(6)->Clear();
0231     m_p0Par[i]->plot();
0232     canvas.Print((name + cName + ".ps").c_str());
0233   }
0234   canvas.Print((name + cName + ".ps]").c_str());
0235 }
0236 
0237 void SoftLeptonTagPlotter::epsPlot(const std::string& name) {
0238   if (willFinalize_)
0239     return;
0240 
0241   for (int i = 0; i != s_leptons; ++i) {
0242     m_leptonId[i]->epsPlot(name);
0243     m_leptonPt[i]->epsPlot(name);
0244     m_sip2d[i]->epsPlot(name);
0245     m_sip3d[i]->epsPlot(name);
0246     m_sip2dsig[i]->epsPlot(name);
0247     m_sip3dsig[i]->epsPlot(name);
0248     m_ptRel[i]->epsPlot(name);
0249     m_p0Par[i]->epsPlot(name);
0250     m_etaRel[i]->epsPlot(name);
0251     m_deltaR[i]->epsPlot(name);
0252     m_ratio[i]->epsPlot(name);
0253     m_ratioRel[i]->epsPlot(name);
0254   }
0255 }