![]() |
![]() |
00001 #include "DQMOffline/RecoB/interface/SoftLeptonTagPlotter.h" 00002 #include "DataFormats/TrackReco/interface/Track.h" 00003 #include "DQMOffline/RecoB/interface/Tools.h" 00004 00005 #include <sstream> 00006 #include <string> 00007 00008 using namespace std ; 00009 using namespace RecoBTag; 00010 00011 static const string ordinal[9] = { "1st", "2nd", "3rd", "4th", "5th", "6th", "7th", "8th", "9th" }; 00012 00013 SoftLeptonTagPlotter::SoftLeptonTagPlotter(const std::string & tagName, 00014 const EtaPtBin & etaPtBin, const edm::ParameterSet& pSet, const bool& mc, const bool& update) : 00015 BaseTagInfoPlotter(tagName, etaPtBin), mcPlots_(mc) 00016 { 00017 const std::string softLepDir(theExtensionString.substr(1)); 00018 00019 for (int i = 0; i < s_leptons; i++) { 00020 std::ostringstream s(""); 00021 s << ordinal[i] << " lepton "; 00022 m_leptonId[i] = new FlavourHistograms<double> ( 00023 s.str() + "id", 00024 "Lepton identification discriminaint", 00025 60, -0.1, 1.1, false, false, true, "b", update,softLepDir,mcPlots_ ); 00026 m_leptonPt[i] = new FlavourHistograms<double> ( 00027 s.str() + "pT", 00028 "Lepton transverse moementum", 00029 100, 0.0, 20.0, false, false, true, "b", update,softLepDir,mcPlots_); 00030 m_sip2d[i] = new FlavourHistograms<double> ( 00031 s.str() + "sip2d", 00032 "Lepton signed 2D impact parameter significance", 00033 100, -20.0, 30.0, false, false, true, "b", update,softLepDir,mcPlots_); 00034 m_sip3d[i] = new FlavourHistograms<double> ( 00035 s.str() + "sip3d", 00036 "Lepton signed 3D impact parameter significance", 00037 100, -20.0, 30.0, false, false, true, "b", update,softLepDir,mcPlots_); 00038 m_ptRel[i] = new FlavourHistograms<double> ( 00039 s.str() + "pT rel", 00040 "Lepton transverse moementum relative to jet axis", 00041 100, 0.0, 10.0, false, false, true, "b", update,softLepDir,mcPlots_); 00042 m_p0Par[i] = new FlavourHistograms<double> ( 00043 s.str() + "p0 par", 00044 "Lepton moementum along jet axis in the B rest frame", 00045 100, 0.0, 10.0, false, false, true, "b", update,softLepDir,mcPlots_); 00046 m_etaRel[i] = new FlavourHistograms<double> ( 00047 s.str() + "eta rel", 00048 "Lepton pseudorapidity relative to jet axis", 00049 100, -5.0, 25.0, false, false, true, "b", update,softLepDir,mcPlots_); 00050 m_deltaR[i] = new FlavourHistograms<double> ( 00051 s.str() + "delta R", 00052 "Lepton pseudoangular distance from jet axis", 00053 100, 0.0, 0.6, false, false, true, "b", update,softLepDir,mcPlots_); 00054 m_ratio[i] = new FlavourHistograms<double> ( 00055 s.str() + "energy ratio", 00056 "Ratio of lepton momentum to jet energy", 00057 100, 0.0, 2.0, false, false, true, "b", update,softLepDir,mcPlots_); 00058 m_ratioRel[i] = new FlavourHistograms<double> ( 00059 s.str() + "parallel energy ratio", 00060 "Ratio of lepton momentum along the jet axis to jet energy", 00061 100, 0.0, 2.0, false, false, true, "b", update,softLepDir,mcPlots_); 00062 } 00063 } 00064 00065 SoftLeptonTagPlotter::~SoftLeptonTagPlotter () 00066 { 00067 for (int i = 0; i != s_leptons; ++i) { 00068 delete m_leptonId[i]; 00069 delete m_leptonPt[i]; 00070 delete m_sip2d[i]; 00071 delete m_sip3d[i]; 00072 delete m_ptRel[i]; 00073 delete m_p0Par[i]; 00074 delete m_etaRel[i]; 00075 delete m_deltaR[i]; 00076 delete m_ratio[i]; 00077 delete m_ratioRel[i]; 00078 } 00079 } 00080 00081 void SoftLeptonTagPlotter::analyzeTag( const reco::BaseTagInfo * baseTagInfo, 00082 const int & jetFlavour ) 00083 { 00084 00085 const reco::SoftLeptonTagInfo * tagInfo = 00086 dynamic_cast<const reco::SoftLeptonTagInfo *>(baseTagInfo); 00087 00088 if (!tagInfo) { 00089 throw cms::Exception("Configuration") 00090 << "BTagPerformanceAnalyzer: Extended TagInfo not of type SoftLeptonTagInfo. " << endl; 00091 } 00092 00093 int n_leptons = tagInfo->leptons(); 00094 00095 for (int i = 0; i != n_leptons && i != s_leptons; ++i) { 00096 const reco::SoftLeptonProperties& properties = tagInfo->properties(i); 00097 m_leptonPt[i]->fill( jetFlavour, tagInfo->lepton(i)->pt() ); 00098 m_leptonId[i]->fill( jetFlavour, properties.quality() ); 00099 m_sip2d[i]->fill( jetFlavour, properties.sip2d ); 00100 m_sip3d[i]->fill( jetFlavour, properties.sip3d ); 00101 m_ptRel[i]->fill( jetFlavour, properties.ptRel ); 00102 m_p0Par[i]->fill( jetFlavour, properties.p0Par ); 00103 m_etaRel[i]->fill( jetFlavour, properties.etaRel ); 00104 m_deltaR[i]->fill( jetFlavour, properties.deltaR ); 00105 m_ratio[i]->fill( jetFlavour, properties.ratio ); 00106 m_ratioRel[i]->fill( jetFlavour, properties.ratioRel ); 00107 } 00108 } 00109 00110 void SoftLeptonTagPlotter::psPlot(const std::string & name) 00111 { 00112 const std::string cName("SoftLeptonPlots" + theExtensionString); 00113 setTDRStyle()->cd(); 00114 TCanvas canvas(cName.c_str(), cName.c_str(), 600, 900); 00115 canvas.UseCurrentStyle(); 00116 canvas.Divide(2,3); 00117 canvas.Print((name + cName + ".ps[").c_str()); 00118 for (int i = 0; i < s_leptons; i++) { 00119 canvas.cd(1)->Clear(); 00120 m_leptonId[i]->plot(); 00121 canvas.cd(2)->Clear(); 00122 m_leptonPt[i]->plot(); 00123 canvas.cd(3)->Clear(); 00124 m_sip2d[i]->plot(); 00125 canvas.cd(4)->Clear(); 00126 m_sip3d[i]->plot(); 00127 canvas.cd(5)->Clear(); 00128 m_ptRel[i]->plot(); 00129 canvas.cd(6)->Clear(); 00130 m_p0Par[i]->plot(); 00131 canvas.Print((name + cName + ".ps").c_str()); 00132 00133 canvas.cd(1)->Clear(); 00134 m_etaRel[i]->plot(); 00135 canvas.cd(2)->Clear(); 00136 m_deltaR[i]->plot(); 00137 canvas.cd(3)->Clear(); 00138 m_ratio[i]->plot(); 00139 canvas.cd(4)->Clear(); 00140 m_ratioRel[i]->plot(); 00141 canvas.cd(5)->Clear(); 00142 canvas.cd(6)->Clear(); 00143 canvas.Print((name + cName + ".ps").c_str()); 00144 } 00145 canvas.Print((name + cName + ".ps]").c_str()); 00146 } 00147 00148 00149 void SoftLeptonTagPlotter::epsPlot(const std::string & name) 00150 { 00151 for (int i=0; i != s_leptons; ++i) { 00152 m_leptonId[i]->epsPlot( name ); 00153 m_leptonPt[i]->epsPlot( name ); 00154 m_sip2d[i]->epsPlot( name ); 00155 m_sip3d[i]->epsPlot( name ); 00156 m_ptRel[i]->epsPlot( name ); 00157 m_p0Par[i]->epsPlot( name ); 00158 m_etaRel[i]->epsPlot( name ); 00159 m_deltaR[i]->epsPlot( name ); 00160 m_ratio[i]->epsPlot( name ); 00161 m_ratioRel[i]->epsPlot( name ); 00162 } 00163 } 00164