![]() |
![]() |
#include <BTagPerformanceAnalyzerMC.h>
Top level steering routine for b tag performance analysis.
Definition at line 45 of file BTagPerformanceAnalyzerMC.h.
typedef std::map<edm::RefToBase<reco::Jet>, unsigned int, JetRefCompare> BTagPerformanceAnalyzerMC::FlavourMap [private] |
Definition at line 68 of file BTagPerformanceAnalyzerMC.h.
typedef std::pair<reco::Jet, reco::JetFlavour> BTagPerformanceAnalyzerMC::JetWithFlavour [private] |
Definition at line 67 of file BTagPerformanceAnalyzerMC.h.
typedef std::map<edm::RefToBase<reco::Jet>, reco::JetFlavour::Leptons, JetRefCompare> BTagPerformanceAnalyzerMC::LeptonMap [private] |
Definition at line 69 of file BTagPerformanceAnalyzerMC.h.
BTagPerformanceAnalyzerMC::BTagPerformanceAnalyzerMC | ( | const edm::ParameterSet & | pSet | ) | [explicit] |
Definition at line 20 of file BTagPerformanceAnalyzerMC.cc.
References bookHistos(), electronPlots, edm::ParameterSet::getParameter(), jetMatcher, muonPlots, MatchJet::setThreshold(), and tauPlots.
: partonKinematics(pSet.getParameter< bool >("partonKinematics")), ptPartonMin(pSet.getParameter<double>("ptPartonMin")), ptPartonMax(pSet.getParameter<double>("ptPartonMax")), jetSelector( pSet.getParameter<double>("etaMin"), pSet.getParameter<double>("etaMax"), pSet.getParameter<double>("ptRecJetMin"), pSet.getParameter<double>("ptRecJetMax"), 0.0, 99999.0, pSet.getParameter<double>("ratioMin"), pSet.getParameter<double>("ratioMax") ), etaRanges(pSet.getParameter< vector<double> >("etaRanges")), ptRanges(pSet.getParameter< vector<double> >("ptRanges")), produceEps(pSet.getParameter< bool >("produceEps")), producePs(pSet.getParameter< bool >("producePs")), inputFile(pSet.getParameter<std::string>( "inputfile" )), update(pSet.getParameter<bool>( "update" )), allHisto(pSet.getParameter<bool>( "allHistograms" )), finalize(pSet.getParameter< bool >("finalizePlots")), finalizeOnly(pSet.getParameter< bool >("finalizeOnly")), jetMCSrc(pSet.getParameter<edm::InputTag>("jetMCSrc")), slInfoTag(pSet.getParameter<edm::InputTag>("softLeptonInfo")), moduleConfig(pSet.getParameter< vector<edm::ParameterSet> >("tagConfig")), mcPlots_(pSet.getParameter< bool >("mcPlots")), makeDiffPlots_(pSet.getParameter< bool >("differentialPlots")), jetCorrector(pSet.getParameter<std::string>("jetCorrection")), jetMatcher(pSet.getParameter<edm::ParameterSet>("recJetMatching")) { double ptRecJetMin = pSet.getParameter<double>("ptRecJetMin"); jetMatcher.setThreshold(0.25 * ptRecJetMin); switch(pSet.getParameter<unsigned int>("leptonPlots")) { case 11: electronPlots = true; muonPlots = false; tauPlots = false; break; case 13: muonPlots = true; electronPlots = false; tauPlots = false; break; case 15: tauPlots = true; electronPlots = false; tauPlots = false; break; default: electronPlots = false; muonPlots = false; tauPlots = false; } bookHistos(pSet); }
BTagPerformanceAnalyzerMC::~BTagPerformanceAnalyzerMC | ( | ) |
Definition at line 235 of file BTagPerformanceAnalyzerMC.cc.
References binJetTagPlotters, binTagCorrelationPlotters, binTagInfoPlotters, differentialPlots, finalize, makeDiffPlots_, and mcPlots_.
{ for (vector<vector<JetTagPlotter*> >::iterator iJetLabel = binJetTagPlotters.begin(); iJetLabel != binJetTagPlotters.end(); ++iJetLabel) for (vector<JetTagPlotter*>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) delete *iPlotter; if (finalize && mcPlots_ && makeDiffPlots_) { for(vector<vector<BTagDifferentialPlot*> >::iterator iJetLabel = differentialPlots.begin(); iJetLabel != differentialPlots.end(); ++iJetLabel) for (vector<BTagDifferentialPlot *>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++ iPlotter) delete *iPlotter; } for (vector<vector<TagCorrelationPlotter*> >::iterator iJetLabel = binTagCorrelationPlotters.begin(); iJetLabel != binTagCorrelationPlotters.end(); ++iJetLabel) for (vector<TagCorrelationPlotter* >::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) delete *iPlotter; for (vector<vector<BaseTagInfoPlotter*> >::iterator iJetLabel = binTagInfoPlotters.begin(); iJetLabel != binTagInfoPlotters.end(); ++iJetLabel) for (vector<BaseTagInfoPlotter*>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) delete *iPlotter; }
void BTagPerformanceAnalyzerMC::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
needed for lepton specific plots
needed for lepton specific plots
Implements edm::EDAnalyzer.
Definition at line 263 of file BTagPerformanceAnalyzerMC.cc.
References abs, edm::AssociationVector< KeyRefProd, CVal, KeyRef, SizeType, KeyReferenceHelper >::begin(), binJetTagPlotters, binTagCorrelationPlotters, binTagInfoPlotters, binTagInfoPlottersToModuleConfig, metsig::electron, electronPlots, edm::AssociationVector< KeyRefProd, CVal, KeyRef, SizeType, KeyReferenceHelper >::end(), eventInitialized, Exception, finalizeOnly, edm::Event::getByLabel(), getJetWithFlavour(), edm::ProductID::id(), reco::BaseTagInfo::jet(), jetMCSrc, jetSelector, jetTagInputTags, L1TDQM_cfg::labels, EgammaValidation_Wenu_cff::leptons, LogDebug, moduleConfig, metsig::muon, muonPlots, partonKinematics, edm::Handle< T >::product(), findQualityFiles::size, edm::AssociationVector< KeyRefProd, CVal, KeyRef, SizeType, KeyReferenceHelper >::size(), slInfoTag, tagCorrelationInputTags, tagInfoInputTags, BaseTagInfoPlotter::tagInfoRequirements(), metsig::tau, tauPlots, and tiDataFormatType.
{ eventInitialized = false; if (finalizeOnly) return; edm::Handle<JetFlavourMatchingCollection> jetMC; FlavourMap flavours; LeptonMap leptons; iEvent.getByLabel(jetMCSrc, jetMC); for (JetFlavourMatchingCollection::const_iterator iter = jetMC->begin(); iter != jetMC->end(); ++iter) { unsigned int fl = std::abs(iter->second.getFlavour()); flavours.insert(std::make_pair(iter->first, fl)); const reco::JetFlavour::Leptons &lep = iter->second.getLeptons(); leptons.insert(std::make_pair(iter->first, lep)); } edm::Handle<reco::SoftLeptonTagInfoCollection> infoHandle; iEvent.getByLabel(slInfoTag, infoHandle); // Look first at the jetTags for (unsigned int iJetLabel = 0; iJetLabel != jetTagInputTags.size(); ++iJetLabel) { edm::Handle<reco::JetTagCollection> tagHandle; iEvent.getByLabel(jetTagInputTags[iJetLabel], tagHandle); const reco::JetTagCollection & tagColl = *(tagHandle.product()); LogDebug("Info") << "Found " << tagColl.size() << " B candidates in collection " << jetTagInputTags[iJetLabel]; int plotterSize = binJetTagPlotters[iJetLabel].size(); for (JetTagCollection::const_iterator tagI = tagColl.begin(); tagI != tagColl.end(); ++tagI) { // Identify parton associated to jet. if (flavours[tagI->first] == 5 && ((electronPlots && !leptons[tagI->first].electron) || (muonPlots && !leptons[tagI->first].muon) || (tauPlots && !leptons[tagI->first].tau))) continue; JetWithFlavour jetWithFlavour; if (!getJetWithFlavour(tagI->first, flavours, jetWithFlavour, iSetup)) continue; if (!jetSelector(jetWithFlavour.first, std::abs(jetWithFlavour.second.getFlavour()), infoHandle)) continue; for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) { bool inBin = false; if (partonKinematics) inBin = binJetTagPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.second.getLorentzVector().Eta(), jetWithFlavour.second.getLorentzVector().Pt()); else inBin = binJetTagPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.first); // Fill histograms if in desired pt/rapidity bin. if (inBin) binJetTagPlotters[iJetLabel][iPlotter]->analyzeTag(jetWithFlavour.first, tagI->second, std::abs(jetWithFlavour.second.getFlavour())); } } } // Now look at Tag Correlations for (unsigned int iJetLabel = 0; iJetLabel != tagCorrelationInputTags.size(); ++iJetLabel) { const std::pair<edm::InputTag, edm::InputTag>& inputTags = tagCorrelationInputTags[iJetLabel]; edm::Handle<reco::JetTagCollection> tagHandle1; iEvent.getByLabel(inputTags.first, tagHandle1); const reco::JetTagCollection& tagColl1 = *(tagHandle1.product()); edm::Handle<reco::JetTagCollection> tagHandle2; iEvent.getByLabel(inputTags.second, tagHandle2); const reco::JetTagCollection& tagColl2 = *(tagHandle2.product()); int plotterSize = binTagCorrelationPlotters[iJetLabel].size(); for (JetTagCollection::const_iterator tagI = tagColl1.begin(); tagI != tagColl1.end(); ++tagI) { if (flavours[tagI->first] == 5 && ((electronPlots && !leptons[tagI->first].electron) || (muonPlots && !leptons[tagI->first].muon) || (tauPlots && !leptons[tagI->first].tau))) continue; JetWithFlavour jetWithFlavour; if (!getJetWithFlavour(tagI->first, flavours, jetWithFlavour, iSetup)) continue; if (!jetSelector(jetWithFlavour.first, std::abs(jetWithFlavour.second.getFlavour()), infoHandle)) continue; for(int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) { bool inBin = false; if (partonKinematics) inBin = binTagCorrelationPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.second.getLorentzVector().Eta(), jetWithFlavour.second.getLorentzVector().Pt()); else inBin = binTagCorrelationPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.first); if(inBin) { double discr2 = tagColl2[tagI->first]; binTagCorrelationPlotters[iJetLabel][iPlotter]->analyzeTags(tagI->second, discr2, std::abs(jetWithFlavour.second.getFlavour())); } } } } // Now look at the TagInfos for (unsigned int iJetLabel = 0; iJetLabel != tiDataFormatType.size(); ++iJetLabel) { int plotterSize = binTagInfoPlotters[iJetLabel].size(); for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) binTagInfoPlotters[iJetLabel][iPlotter]->setEventSetup(iSetup); vector<edm::InputTag> & inputTags = tagInfoInputTags[iJetLabel]; if (inputTags.empty()) { // deferred retrieval of input tags BaseTagInfoPlotter *firstPlotter = binTagInfoPlotters[iJetLabel][0]; int iModule = binTagInfoPlottersToModuleConfig[firstPlotter]; vector<string> labels = firstPlotter->tagInfoRequirements(); if (labels.empty()) labels.push_back("label"); for (vector<string>::const_iterator iLabels = labels.begin(); iLabels != labels.end(); ++iLabels) { edm::InputTag inputTag = moduleConfig[iModule].getParameter<InputTag>(*iLabels); inputTags.push_back(inputTag); } } unsigned int nInputTags = inputTags.size(); vector< edm::Handle< View<BaseTagInfo> > > tagInfoHandles(nInputTags); edm::ProductID jetProductID; unsigned int nTagInfos = 0; for (unsigned int iInputTags = 0; iInputTags < inputTags.size(); ++iInputTags) { edm::Handle< View<BaseTagInfo> > & tagInfoHandle = tagInfoHandles[iInputTags]; iEvent.getByLabel(inputTags[iInputTags], tagInfoHandle); unsigned int size = tagInfoHandle->size(); LogDebug("Info") << "Found " << size << " B candidates in collection " << inputTags[iInputTags]; edm::ProductID thisProductID = (size > 0) ? (*tagInfoHandle)[0].jet().id() : edm::ProductID(); if (iInputTags == 0) { jetProductID = thisProductID; nTagInfos = size; } else if (jetProductID != thisProductID) throw cms::Exception("Configuration") << "TagInfos are referencing a different jet collection." << endl; else if (nTagInfos != size) throw cms::Exception("Configuration") << "TagInfo collections are having a different size." << endl; } for (unsigned int iTagInfos = 0; iTagInfos < nTagInfos; ++iTagInfos) { vector<const BaseTagInfo*> baseTagInfos(nInputTags); edm::RefToBase<Jet> jetRef; for (unsigned int iTagInfo = 0; iTagInfo < nInputTags; iTagInfo++) { const BaseTagInfo &baseTagInfo = (*tagInfoHandles[iTagInfo])[iTagInfos]; if (iTagInfo == 0) jetRef = baseTagInfo.jet(); else if (baseTagInfo.jet() != jetRef) throw cms::Exception("Configuration") << "TagInfos pointing to different jets." << endl; baseTagInfos[iTagInfo] = &baseTagInfo; } // Identify parton associated to jet. if (flavours[jetRef] == 5 && ((electronPlots && !leptons[jetRef].electron) || (muonPlots && !leptons[jetRef].muon) || (tauPlots && !leptons[jetRef].tau))) continue; JetWithFlavour jetWithFlavour; if (!getJetWithFlavour(jetRef, flavours, jetWithFlavour, iSetup)) continue; if (!jetSelector(jetWithFlavour.first, std::abs(jetWithFlavour.second.getFlavour()), infoHandle)) continue; for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) { bool inBin = false; if (partonKinematics) inBin = binTagInfoPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.second.getLorentzVector().Eta(), jetWithFlavour.second.getLorentzVector().Pt()); else inBin = binTagInfoPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(*jetRef); // Fill histograms if in desired pt/rapidity bin. if (inBin) binTagInfoPlotters[iJetLabel][iPlotter]->analyzeTag(baseTagInfos, std::abs(jetWithFlavour.second.getFlavour())); } } } }
void BTagPerformanceAnalyzerMC::bookHistos | ( | const edm::ParameterSet & | pSet | ) | [private] |
Definition at line 233 of file DTDataIntegrityTest.cc.
References DTDataIntegrityTest::dbe, and DQMStore::setCurrentFolder().
Referenced by BTagPerformanceAnalyzerMC().
{ stringstream dduId_s; dduId_s << dduId; dbe->setCurrentFolder("DT/00-DataIntegrity/FED" + dduId_s.str()); string histoName; }
void BTagPerformanceAnalyzerMC::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 496 of file BTagPerformanceAnalyzerMC.cc.
References begin, binJetTagPlotters, binTagInfoPlotters, differentialPlots, epsBaseName, finalize, makeDiffPlots_, produceEps, producePs, psBaseName, and plotscripts::setTDRStyle().
{ if (!finalize) return; setTDRStyle(); for (unsigned int iJetLabel = 0; iJetLabel != binJetTagPlotters.size(); ++iJetLabel) { int plotterSize = binJetTagPlotters[iJetLabel].size(); for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) { binJetTagPlotters[iJetLabel][iPlotter]->finalize(); // binJetTagPlotters[iJetLabel][iPlotter]->write(allHisto); if (producePs) (*binJetTagPlotters[iJetLabel][iPlotter]).psPlot(psBaseName); if (produceEps) (*binJetTagPlotters[iJetLabel][iPlotter]).epsPlot(epsBaseName); } if(makeDiffPlots_) { for (vector<BTagDifferentialPlot *>::iterator iPlotter = differentialPlots[iJetLabel].begin(); iPlotter != differentialPlots[iJetLabel].end(); ++ iPlotter) { (*iPlotter)->process(); if (producePs) (*iPlotter)->psPlot(psBaseName); if (produceEps) (*iPlotter)->epsPlot(epsBaseName); // (**iPlotter).write(allHisto); } } } for (vector<vector<BaseTagInfoPlotter*> >::iterator iJetLabel = binTagInfoPlotters.begin(); iJetLabel != binTagInfoPlotters.end(); ++iJetLabel) { for (vector<BaseTagInfoPlotter*>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) { (*iPlotter)->finalize(); // binTagInfoPlotters[iJetLabel][iPlotter]->write(allHisto); if (producePs) (*iPlotter)->psPlot(psBaseName); if (produceEps) (*iPlotter)->epsPlot(epsBaseName); } } }
EtaPtBin BTagPerformanceAnalyzerMC::getEtaPtBin | ( | const int & | iEta, |
const int & | iPt | ||
) | [private] |
Definition at line 204 of file BTagPerformanceAnalyzerMC.cc.
References etaRanges, funct::false, ptRanges, and funct::true.
{ // DEFINE BTagBin: bool etaActive_ , ptActive_; double etaMin_, etaMax_, ptMin_, ptMax_ ; if ( iEta != -1 ) { etaActive_ = true ; etaMin_ = etaRanges[iEta] ; etaMax_ = etaRanges[iEta+1] ; } else { etaActive_ = false ; etaMin_ = etaRanges[0] ; etaMax_ = etaRanges[etaRanges.size() - 1] ; } if ( iPt != -1 ) { ptActive_ = true ; ptMin_ = ptRanges[iPt] ; ptMax_ = ptRanges[iPt+1] ; } else { ptActive_ = false ; ptMin_ = ptRanges[0] ; ptMax_ = ptRanges[ptRanges.size() - 1] ; } return EtaPtBin(etaActive_ , etaMin_ , etaMax_ , ptActive_ , ptMin_ , ptMax_ ); }
bool BTagPerformanceAnalyzerMC::getJetWithFlavour | ( | edm::RefToBase< reco::Jet > | caloRef, |
FlavourMap | flavours, | ||
JetWithFlavour & | jetWithFlavour, | ||
const edm::EventSetup & | es | ||
) | [private] |
Definition at line 454 of file BTagPerformanceAnalyzerMC.cc.
References eventInitialized, i, edm::RefToBase< T >::id(), edm::ProductID::id(), edm::RefToBase< T >::isNull(), jetCorrector, jetMatcher, LogTrace, MatchJet::matchCollections(), edm::RefToBaseVector< T >::push_back(), and CorrectJet::setEventSetup().
Referenced by analyze().
{ edm::ProductID recProdId = jetRef.id(); edm::ProductID refProdId = (flavours.begin() == flavours.end()) ? recProdId : flavours.begin()->first.id(); if (!eventInitialized) { jetCorrector.setEventSetup(es); if (recProdId != refProdId) { edm::RefToBaseVector<Jet> refJets; for(FlavourMap::const_iterator iter = flavours.begin(); iter != flavours.end(); ++iter) refJets.push_back(iter->first); const edm::RefToBaseProd<Jet> recJetsProd(jetRef); edm::RefToBaseVector<Jet> recJets; for(unsigned int i = 0; i < recJetsProd->size(); i++) recJets.push_back(edm::RefToBase<Jet>(recJetsProd, i)); jetMatcher.matchCollections(refJets, recJets, es); } eventInitialized = true; } if (recProdId != refProdId) { jetRef = jetMatcher(jetRef); if (jetRef.isNull()) return false; } jetWithFlavour.first = jetCorrector(*jetRef); jetWithFlavour.second = reco::JetFlavour(jetWithFlavour.first.p4(), math::XYZPoint (0,0,0), flavours[jetRef]); LogTrace("Info") << "Found jet with flavour "<<jetWithFlavour.second.getFlavour()<<endl; LogTrace("Info") << jetWithFlavour.first.p()<<" , "<< jetWithFlavour.first.pt()<<" - " << jetWithFlavour.second.getLorentzVector().P()<<" , "<< jetWithFlavour.second.getLorentzVector().Pt()<<endl; return true; }
bool BTagPerformanceAnalyzerMC::allHisto [private] |
Definition at line 83 of file BTagPerformanceAnalyzerMC.h.
std::vector< std::vector<JetTagPlotter*> > BTagPerformanceAnalyzerMC::binJetTagPlotters [private] |
Definition at line 89 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), endJob(), and ~BTagPerformanceAnalyzerMC().
std::vector< std::vector<TagCorrelationPlotter*> > BTagPerformanceAnalyzerMC::binTagCorrelationPlotters [private] |
Definition at line 90 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), and ~BTagPerformanceAnalyzerMC().
std::vector< std::vector<BaseTagInfoPlotter*> > BTagPerformanceAnalyzerMC::binTagInfoPlotters [private] |
Definition at line 91 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), endJob(), and ~BTagPerformanceAnalyzerMC().
std::map<BaseTagInfoPlotter*, size_t> BTagPerformanceAnalyzerMC::binTagInfoPlottersToModuleConfig [private] |
Definition at line 99 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
std::vector< std::vector<BTagDifferentialPlot*> > BTagPerformanceAnalyzerMC::differentialPlots [private] |
Definition at line 96 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob(), and ~BTagPerformanceAnalyzerMC().
bool BTagPerformanceAnalyzerMC::electronPlots [private] |
Definition at line 107 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), and BTagPerformanceAnalyzerMC().
std::string BTagPerformanceAnalyzerMC::epsBaseName [private] |
Definition at line 82 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob().
std::vector<double> BTagPerformanceAnalyzerMC::etaRanges [private] |
Definition at line 80 of file BTagPerformanceAnalyzerMC.h.
Referenced by getEtaPtBin().
bool BTagPerformanceAnalyzerMC::eventInitialized [private] |
Definition at line 106 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), and getJetWithFlavour().
bool BTagPerformanceAnalyzerMC::finalize [private] |
Definition at line 84 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob(), and ~BTagPerformanceAnalyzerMC().
bool BTagPerformanceAnalyzerMC::finalizeOnly [private] |
Definition at line 85 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
std::string BTagPerformanceAnalyzerMC::inputFile [private] |
Definition at line 82 of file BTagPerformanceAnalyzerMC.h.
Definition at line 103 of file BTagPerformanceAnalyzerMC.h.
Referenced by getJetWithFlavour().
Definition at line 104 of file BTagPerformanceAnalyzerMC.h.
Referenced by BTagPerformanceAnalyzerMC(), and getJetWithFlavour().
Definition at line 86 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
Definition at line 79 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
std::vector<edm::InputTag> BTagPerformanceAnalyzerMC::jetTagInputTags [private] |
Definition at line 92 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
bool BTagPerformanceAnalyzerMC::makeDiffPlots_ [private] |
Definition at line 101 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob(), and ~BTagPerformanceAnalyzerMC().
bool BTagPerformanceAnalyzerMC::mcPlots_ [private] |
Definition at line 101 of file BTagPerformanceAnalyzerMC.h.
Referenced by ~BTagPerformanceAnalyzerMC().
std::vector<edm::ParameterSet> BTagPerformanceAnalyzerMC::moduleConfig [private] |
Definition at line 98 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
bool BTagPerformanceAnalyzerMC::muonPlots [private] |
Definition at line 107 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), and BTagPerformanceAnalyzerMC().
bool BTagPerformanceAnalyzerMC::partonKinematics [private] |
Definition at line 77 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
bool BTagPerformanceAnalyzerMC::produceEps [private] |
Definition at line 81 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob().
bool BTagPerformanceAnalyzerMC::producePs [private] |
Definition at line 81 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob().
std::string BTagPerformanceAnalyzerMC::psBaseName [private] |
Definition at line 82 of file BTagPerformanceAnalyzerMC.h.
Referenced by endJob().
double BTagPerformanceAnalyzerMC::ptPartonMax [private] |
Definition at line 78 of file BTagPerformanceAnalyzerMC.h.
double BTagPerformanceAnalyzerMC::ptPartonMin [private] |
Definition at line 78 of file BTagPerformanceAnalyzerMC.h.
std::vector<double> BTagPerformanceAnalyzerMC::ptRanges [private] |
Definition at line 80 of file BTagPerformanceAnalyzerMC.h.
Referenced by getEtaPtBin().
Definition at line 87 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
std::vector< std::pair<edm::InputTag, edm::InputTag> > BTagPerformanceAnalyzerMC::tagCorrelationInputTags [private] |
Definition at line 93 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
std::vector< std::vector<edm::InputTag> > BTagPerformanceAnalyzerMC::tagInfoInputTags [private] |
Definition at line 94 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
bool BTagPerformanceAnalyzerMC::tauPlots [private] |
Definition at line 107 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze(), and BTagPerformanceAnalyzerMC().
std::vector<std::string> BTagPerformanceAnalyzerMC::tiDataFormatType [private] |
Definition at line 76 of file BTagPerformanceAnalyzerMC.h.
Referenced by analyze().
bool BTagPerformanceAnalyzerMC::update [private] |
Definition at line 83 of file BTagPerformanceAnalyzerMC.h.