Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "HLTriggerOffline/Btag/interface/HLTBTagPerformanceAnalyzer.h"
0002 #include <algorithm>
0003 #include <set>
0004 
0005 using namespace edm;
0006 using namespace reco;
0007 
0008 namespace {
0009   // find the index of the object key of an association vector closest to a given
0010   // jet, within a given distance
0011   template <typename T, typename V>
0012   int closestJet(const RefToBase<reco::Jet> jet, const edm::AssociationVector<T, V> &association, double distance) {
0013     int closest = -1;
0014     for (unsigned int i = 0; i < association.size(); ++i) {
0015       double d = ROOT::Math::VectorUtil::DeltaR(jet->momentum(), association[i].first->momentum());
0016       if (d < distance) {
0017         distance = d;
0018         closest = i;
0019       }
0020     }
0021     return closest;
0022   }
0023 
0024   std::set<std::string> const keepSetJet{"jetNSecondaryVertices",
0025                                          "jetNSelectedTracks",
0026                                          "jetNTracks",
0027                                          "Jet_JP",
0028                                          "chargedHadronEnergyFraction",
0029                                          "neutralHadronEnergyFraction",
0030                                          "photonEnergyFraction",
0031                                          "electronEnergyFraction",
0032                                          "muonEnergyFraction",
0033                                          "chargedHadronMultiplicity",
0034                                          "neutralHadronMultiplicity",
0035                                          "photonMultiplicity",
0036                                          "electronMultiplicity",
0037                                          "muonMultiplicity",
0038                                          "hadronMultiplicity",
0039                                          "hadronPhotonMultiplicity",
0040                                          "totalMultiplicity"};
0041 
0042   std::set<std::string> const keepSetTrack{
0043       "trackChi2",       "trackNTotalHits",    "trackNPixelHits",   "trackSip3dVal",    "trackSip3dSig",
0044       "trackSip2dVal",   "trackSip2dSig",      "trackPtRel",        "trackDeltaR",      "trackPtRatio",
0045       "trackSip3dSig_0", "trackSip3dSig_1",    "trackSip3dSig_2",   "trackSip3dSig_3",  "trackMomentum",
0046       "trackEta",        "trackPhi",           "trackDecayLenVal",  "trackDecayLenSig", "trackJetDistVal",
0047       "trackJetDistSig", "trackSumJetEtRatio", "trackSumJetDeltaR", "trackEtaRel"
0048 
0049   };
0050   std::set<std::string> const keepSetVtx{"vertexMass",
0051                                          "vertexNTracks"
0052                                          "vertexFitProb",
0053                                          "vertexCategory",
0054                                          "vertexEnergyRatio",
0055                                          "vertexJetDeltaR",
0056                                          "vertexBoostOverSqrtJetPt",
0057                                          "flightDistance1dVal",
0058                                          "flightDistance1dSig",
0059                                          "flightDistance2dVal",
0060                                          "flightDistance2dSig",
0061                                          "flightDistance3dVal",
0062                                          "flightDistance3dSig"};
0063 
0064   auto initializeKeepSet() {
0065     std::set<std::string> ret;
0066     // Make a combined set of inputs avaiable
0067     std::set_union(std::begin(keepSetJet),
0068                    std::end(keepSetJet),
0069                    std::begin(keepSetTrack),
0070                    std::end(keepSetTrack),
0071                    std::inserter(ret, std::begin(ret)));
0072     std::set_union(std::begin(ret),
0073                    std::end(ret),
0074                    std::begin(keepSetVtx),
0075                    std::end(keepSetVtx),
0076                    std::inserter(ret, std::begin(ret)));
0077     return ret;
0078   }
0079   std::set<std::string> const keepSet = initializeKeepSet();
0080 }  // namespace
0081 
0082 // constructors and destructor
0083 HLTBTagPerformanceAnalyzer::HLTBTagPerformanceAnalyzer(const edm::ParameterSet &iConfig) {
0084   mainFolder_ = iConfig.getParameter<std::string>("mainFolder");
0085   hlTriggerResults_ = consumes<edm::TriggerResults>(iConfig.getParameter<InputTag>("TriggerResults"));
0086   JetTagCollection_ =
0087       edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag>>("JetTag"),
0088                             [this](edm::InputTag const &tag) { return mayConsume<reco::JetTagCollection>(tag); });
0089   //        shallowTagInfosTokenCalo_ =
0090   //        consumes<std::vector<reco::ShallowTagInfo> >
0091   //        (edm::InputTag("hltDeepCombinedSecondaryVertexBJetTagsInfosCalo"));
0092   shallowTagInfosTokenPf_ =
0093       consumes<std::vector<reco::ShallowTagInfo>>(edm::InputTag("hltDeepCombinedSecondaryVertexBJetTagsInfos"));
0094   m_mcPartons = consumes<JetFlavourMatchingCollection>(iConfig.getParameter<InputTag>("mcPartons"));
0095   hltPathNames_ = iConfig.getParameter<std::vector<std::string>>("HLTPathNames");
0096   edm::ParameterSet mc = iConfig.getParameter<edm::ParameterSet>("mcFlavours");
0097   m_mcLabels = mc.getParameterNamesForType<std::vector<unsigned int>>();
0098 
0099   EDConsumerBase::labelsForToken(m_mcPartons, label);
0100   m_mcPartons_Label = label.module;
0101 
0102   for (unsigned int i = 0; i < JetTagCollection_.size(); i++) {
0103     EDConsumerBase::labelsForToken(JetTagCollection_[i], label);
0104     JetTagCollection_Label.push_back(label.module);
0105   }
0106 
0107   EDConsumerBase::labelsForToken(hlTriggerResults_, label);
0108   hlTriggerResults_Label = label.module;
0109 
0110   for (unsigned int i = 0; i < m_mcLabels.size(); ++i)
0111     m_mcFlavours.push_back(mc.getParameter<std::vector<unsigned int>>(m_mcLabels[i]));
0112   m_mcMatching = m_mcPartons_Label != "none";
0113 
0114   m_mcRadius = 0.3;
0115 
0116   HCALSpecialsNames[HEP17] = "HEP17";
0117   HCALSpecialsNames[HEP18] = "HEP18";
0118   HCALSpecialsNames[HEM17] = "HEM17";
0119 }
0120 
0121 HLTBTagPerformanceAnalyzer::~HLTBTagPerformanceAnalyzer() {
0122   // do anything here that needs to be done at desctruction time
0123   // (e.g. close files, deallocate resources etc.)
0124 }
0125 
0126 void HLTBTagPerformanceAnalyzer::dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) {
0127   triggerConfChanged_ = true;
0128   EDConsumerBase::labelsForToken(hlTriggerResults_, label);
0129 
0130   hltConfigProvider_.init(iRun, iSetup, label.process, triggerConfChanged_);
0131   const std::vector<std::string> &allHltPathNames = hltConfigProvider_.triggerNames();
0132 
0133   // fill hltPathIndexs_ with the trigger number of each hltPathNames_
0134   for (size_t trgs = 0; trgs < hltPathNames_.size(); trgs++) {
0135     unsigned int found = 1;
0136     int it_mem = -1;
0137     for (size_t it = 0; it < allHltPathNames.size(); ++it) {
0138       found = allHltPathNames.at(it).find(hltPathNames_[trgs]);
0139       if (found == 0) {
0140         it_mem = (int)it;
0141       }
0142     }  // for allallHltPathNames
0143     hltPathIndexs_.push_back(it_mem);
0144   }  // for hltPathNames_
0145 
0146   // fill _isfoundHLTs for each hltPathNames_
0147   for (size_t trgs = 0; trgs < hltPathNames_.size(); trgs++) {
0148     if (hltPathIndexs_[trgs] < 0) {
0149       _isfoundHLTs.push_back(false);
0150     } else {
0151       _isfoundHLTs.push_back(true);
0152     }
0153   }
0154 }
0155 
0156 void HLTBTagPerformanceAnalyzer::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0157   bool trigRes = false;
0158   bool MCOK = false;
0159   using namespace edm;
0160 
0161   // get triggerResults
0162   Handle<TriggerResults> TriggerResulsHandler;
0163   Handle<reco::JetFlavourMatchingCollection> h_mcPartons;
0164   if (hlTriggerResults_Label.empty() || hlTriggerResults_Label == "NULL") {
0165     edm::LogInfo("NoTriggerResults") << "TriggerResults ==> Empty";
0166     return;
0167   }
0168   iEvent.getByToken(hlTriggerResults_, TriggerResulsHandler);
0169   if (TriggerResulsHandler.isValid())
0170     trigRes = true;
0171   if (!trigRes) {
0172     edm::LogInfo("NoTriggerResults") << "TriggerResults ==> not readable";
0173     return;
0174   }
0175   const TriggerResults &triggerResults = *(TriggerResulsHandler.product());
0176 
0177   // get partons
0178   if (m_mcMatching && !m_mcPartons_Label.empty() && m_mcPartons_Label != "NULL") {
0179     iEvent.getByToken(m_mcPartons, h_mcPartons);
0180     if (h_mcPartons.isValid())
0181       MCOK = true;
0182   }
0183 
0184   // fill the 1D and 2D DQM plot
0185   Handle<reco::JetTagCollection> JetTagHandler;
0186   for (unsigned int ind = 0; ind < hltPathNames_.size(); ind++) {
0187     bool BtagOK = false;
0188     JetTagMap JetTag;
0189     if (!_isfoundHLTs[ind])
0190       continue;  // if the hltPath is not in the event, skip the event
0191     if (!triggerResults.accept(hltPathIndexs_[ind]))
0192       continue;  // if the hltPath was not accepted skip the event
0193 
0194     // get JetTagCollection
0195     if (!JetTagCollection_Label[ind].empty() && JetTagCollection_Label[ind] != "NULL") {
0196       iEvent.getByToken(JetTagCollection_[ind], JetTagHandler);
0197       iEvent.getByToken(shallowTagInfosTokenPf_, shallowTagInfosPf);
0198       //                        iEvent.getByToken(shallowTagInfosTokenCalo_,
0199       //                        shallowTagInfosCalo);
0200       if (JetTagHandler.isValid())
0201         BtagOK = true;
0202     }
0203 
0204     // fill JetTag map
0205     if (BtagOK)
0206       for (auto iter = JetTagHandler->begin(); iter != JetTagHandler->end(); iter++) {
0207         JetTag.insert(JetTagMap::value_type(iter->first, iter->second));
0208       }
0209     else {
0210       edm::LogInfo("NoCollection") << "Collection " << JetTagCollection_Label[ind] << " ==> not found";
0211       return;
0212     }
0213     // fill Inputs for All
0214     if (shallowTagInfosPf.isValid()) {
0215       for (auto &info : *(shallowTagInfosPf)) {
0216         TaggingVariableList vars = info.taggingVariables();
0217         for (auto entry = vars.begin(); entry != vars.end(); ++entry) {
0218           if (keepSet.find(TaggingVariableTokens[entry->first]) !=
0219               keepSet.end()) {  // if Input name in defined list to keep
0220             try {
0221               H1_.at(ind)[TaggingVariableTokens[entry->first]]->Fill(std::fmax(0.0, entry->second));
0222             } catch (const std::exception &e) {
0223               continue;
0224             }
0225           } else
0226             continue;
0227         }
0228       }
0229     } else {
0230       edm::LogInfo("NoCollection") << "No shallowTagInfosPf collection";
0231     }
0232     // fill tagging
0233     for (auto &BtagJT : JetTag) {
0234       std::map<HCALSpecials, bool> inmodule;
0235       inmodule[HEP17] = (BtagJT.first->phi() >= -0.87) && (BtagJT.first->phi() < -0.52) && (BtagJT.first->eta() > 1.3);
0236       inmodule[HEP18] = (BtagJT.first->phi() >= -0.52) && (BtagJT.first->phi() < -0.17) && (BtagJT.first->eta() > 1.3);
0237       inmodule[HEM17] = (BtagJT.first->phi() >= -0.87) && (BtagJT.first->phi() < -0.52) && (BtagJT.first->eta() < -1.3);
0238 
0239       // fill 1D btag plot for 'all'
0240       H1_.at(ind)[JetTagCollection_Label[ind]]->Fill(std::fmax(0.0, BtagJT.second));
0241       for (const auto &i : HCALSpecialsNames) {
0242         if (inmodule[i.first])
0243           H1mod_.at(ind)[JetTagCollection_Label[ind]][i.first]->Fill(std::fmax(0.0, BtagJT.second));
0244       }
0245       if (MCOK) {
0246         int m = closestJet(BtagJT.first, *h_mcPartons, m_mcRadius);
0247         unsigned int flavour = (m != -1) ? abs((*h_mcPartons)[m].second.getFlavour()) : 0;
0248         for (unsigned int i = 0; i < m_mcLabels.size(); ++i) {
0249           std::string flavour_str = m_mcLabels[i];
0250           flavours_t flav_collection = m_mcFlavours[i];
0251           auto it = std::find(flav_collection.begin(), flav_collection.end(), flavour);
0252           if (it == flav_collection.end())
0253             continue;
0254           std::string label = JetTagCollection_Label[ind] + "__";
0255           label += flavour_str;
0256           H1_.at(ind)[label]->Fill(std::fmax(0.0, BtagJT.second));  // fill 1D btag plot for 'b,c,uds'
0257           for (const auto &j : HCALSpecialsNames) {
0258             if (inmodule[j.first])
0259               H1mod_.at(ind)[label][j.first]->Fill(
0260                   std::fmax(0.0, BtagJT.second));  // fill 1D btag plot for 'b,c,uds' in
0261                                                    // modules (HEP17 etc.)
0262           }
0263           label = JetTagCollection_Label[ind] + "___";
0264           label += flavour_str;
0265           std::string labelEta = label;
0266           std::string labelPhi = label;
0267           std::string labelEtaPhi = label;
0268           std::string labelEtaPhi_threshold = label;
0269           label += "_disc_pT";
0270           H2_.at(ind)[label]->Fill(std::fmax(0.0, BtagJT.second),
0271                                    BtagJT.first->pt());  // fill 2D btag, jetPt plot for 'b,c,uds'
0272           for (const auto &j : HCALSpecialsNames) {
0273             if (inmodule[j.first])
0274               H2mod_.at(ind)[label][j.first]->Fill(std::fmax(0.0, BtagJT.second), BtagJT.first->pt());
0275           }
0276           labelEta += "_disc_eta";
0277           H2Eta_.at(ind)[labelEta]->Fill(std::fmax(0.0, BtagJT.second),
0278                                          BtagJT.first->eta());  // fill 2D btag, jetEta plot for 'b,c,uds'
0279           labelPhi += "_disc_phi";
0280           H2Phi_.at(ind)[labelPhi]->Fill(std::fmax(0.0, BtagJT.second),
0281                                          BtagJT.first->phi());  // fill 2D btag, jetPhi plot for 'b,c,uds'
0282           labelEtaPhi += "_eta_phi";
0283           H2EtaPhi_.at(ind)[labelEtaPhi]->Fill(BtagJT.first->eta(),
0284                                                BtagJT.first->phi());  // fill 2D btag, jetPhi plot for 'b,c,uds'
0285           labelEtaPhi_threshold += "_eta_phi_disc05";
0286           if (BtagJT.second > 0.5) {
0287             H2EtaPhi_threshold_.at(ind)[labelEtaPhi_threshold]->Fill(
0288                 BtagJT.first->eta(),
0289                 BtagJT.first->phi());  // fill 2D btag, jetPhi plot for 'b,c,uds'
0290           }
0291         }  /// for flavour
0292       }    /// if MCOK
0293     }      /// for BtagJT
0294   }        // for triggers
0295 }
0296 
0297 //// ------------ method called once each job just before starting event loop
0298 ///------------
0299 void HLTBTagPerformanceAnalyzer::bookHistograms(DQMStore::IBooker &ibooker,
0300                                                 edm::Run const &iRun,
0301                                                 edm::EventSetup const &iSetup) {
0302   // book the DQM plots for each path and for each flavour
0303   using namespace std;
0304   assert(hltPathNames_.size() == JetTagCollection_.size());
0305   std::string dqmFolder;
0306   for (unsigned int ind = 0; ind < hltPathNames_.size(); ind++) {
0307     float btagL = 0.;
0308     float btagU = 1.;
0309     int btagBins = 100;
0310     dqmFolder = Form("%s/Discriminator/%s", mainFolder_.c_str(), hltPathNames_[ind].c_str());
0311     H1_.push_back(std::map<std::string, MonitorElement *>());
0312     H2_.push_back(std::map<std::string, MonitorElement *>());
0313     H1mod_.push_back(std::map<std::string, std::map<HCALSpecials, MonitorElement *>>());
0314     H2mod_.push_back(std::map<std::string, std::map<HCALSpecials, MonitorElement *>>());
0315     H2Eta_.push_back(std::map<std::string, MonitorElement *>());
0316     H2Phi_.push_back(std::map<std::string, MonitorElement *>());
0317     H2EtaPhi_.push_back(std::map<std::string, MonitorElement *>());
0318     H2EtaPhi_threshold_.push_back(std::map<std::string, MonitorElement *>());
0319     ibooker.setCurrentFolder(dqmFolder);
0320 
0321     // book 1D btag plot for 'all'
0322     if (!JetTagCollection_Label[ind].empty() && JetTagCollection_Label[ind] != "NULL") {
0323       H1_.back()[JetTagCollection_Label[ind]] = ibooker.book1D(
0324           JetTagCollection_Label[ind] + "_all", JetTagCollection_Label[ind] + "_all", btagBins, btagL, btagU);
0325       H1_.back()[JetTagCollection_Label[ind]]->setAxisTitle(JetTagCollection_Label[ind] + "discriminant", 1);
0326       // Input storing
0327       ibooker.setCurrentFolder(dqmFolder + "/inputs");
0328       ibooker.setCurrentFolder(dqmFolder + "/inputs/Jet");
0329       for (int i = 0; i < 100; i++) {
0330         if (keepSetJet.find(TaggingVariableTokens[i]) != keepSetJet.end()) {  // if input name in defined set
0331           std::string inpt = TaggingVariableTokens[i];
0332           H1_.back()[inpt] = ibooker.book1D(inpt, inpt, 105, -5, 100.);
0333           H1_.back()[inpt]->setAxisTitle(inpt, 1);
0334         } else
0335           continue;
0336       }
0337       ibooker.setCurrentFolder(dqmFolder + "/inputs/Track");
0338       for (int i = 0; i < 100; i++) {
0339         if (keepSetTrack.find(TaggingVariableTokens[i]) != keepSetTrack.end()) {  // if input name in defined set
0340           std::string inpt = TaggingVariableTokens[i];
0341           H1_.back()[inpt] = ibooker.book1D(inpt, inpt, 105, -5, 100.);
0342           H1_.back()[inpt]->setAxisTitle(inpt, 1);
0343         } else
0344           continue;
0345       }
0346       ibooker.setCurrentFolder(dqmFolder + "/inputs/Vertex");
0347       for (int i = 0; i < 100; i++) {
0348         if (keepSetVtx.find(TaggingVariableTokens[i]) != keepSetVtx.end()) {  // if input name in defined set
0349           std::string inpt = TaggingVariableTokens[i];
0350           H1_.back()[inpt] = ibooker.book1D(inpt, inpt, 105, -5, 100.);
0351           H1_.back()[inpt]->setAxisTitle(inpt, 1);
0352         } else
0353           continue;
0354       }
0355 
0356       for (const auto &i : HCALSpecialsNames) {
0357         ibooker.setCurrentFolder(dqmFolder + "/" + i.second);
0358         H1mod_.back()[JetTagCollection_Label[ind]][i.first] = ibooker.book1D(
0359             JetTagCollection_Label[ind] + "_all", JetTagCollection_Label[ind] + "_all", btagBins, btagL, btagU);
0360         H1mod_.back()[JetTagCollection_Label[ind]][i.first]->setAxisTitle(JetTagCollection_Label[ind] + "discriminant",
0361                                                                           1);
0362       }
0363       ibooker.setCurrentFolder(dqmFolder);
0364     }
0365     int nBinsPt = 60;
0366     double pTmin = 30;
0367     double pTMax = 330;
0368     int nBinsPhi = 54;
0369     double phimin = -M_PI;
0370     double phiMax = M_PI;
0371     int nBinsEta = 40;
0372     double etamin = -2.4;
0373     double etaMax = 2.4;
0374 
0375     for (unsigned int i = 0; i < m_mcLabels.size(); ++i) {
0376       std::string flavour = m_mcLabels[i];
0377       std::string label;
0378       std::string labelEta;
0379       std::string labelPhi;
0380       std::string labelEtaPhi;
0381       std::string labelEtaPhi_threshold;
0382       if (!JetTagCollection_Label[ind].empty() && JetTagCollection_Label[ind] != "NULL") {
0383         label = JetTagCollection_Label[ind] + "__";
0384         label += flavour;
0385 
0386         // book 1D btag plot for 'b,c,light,g'
0387         H1_.back()[label] = ibooker.book1D(
0388             label, Form("%s %s", JetTagCollection_Label[ind].c_str(), flavour.c_str()), btagBins, btagL, btagU);
0389         H1_.back()[label]->setAxisTitle("disc", 1);
0390         for (const auto &j : HCALSpecialsNames) {
0391           ibooker.setCurrentFolder(dqmFolder + "/" + j.second);
0392           H1mod_.back()[label][j.first] = ibooker.book1D(
0393               label, Form("%s %s", JetTagCollection_Label[ind].c_str(), flavour.c_str()), btagBins, btagL, btagU);
0394           H1mod_.back()[label][j.first]->setAxisTitle("disc", 1);
0395         }
0396         ibooker.setCurrentFolder(dqmFolder);
0397         label = JetTagCollection_Label[ind] + "___";
0398         labelEta = label;
0399         labelPhi = label;
0400         labelEtaPhi = label;
0401         labelEtaPhi_threshold = label;
0402         label += flavour + "_disc_pT";
0403         labelEta += flavour + "_disc_eta";
0404         labelPhi += flavour + "_disc_phi";
0405         labelEtaPhi += flavour + "_eta_phi";
0406         labelEtaPhi_threshold += flavour + "_eta_phi_disc05";
0407 
0408         // book 2D btag plot for 'b,c,light,g'
0409         H2_.back()[label] = ibooker.book2D(label, label, btagBins, btagL, btagU, nBinsPt, pTmin, pTMax);
0410         H2_.back()[label]->setAxisTitle("pT", 2);
0411         H2_.back()[label]->setAxisTitle("disc", 1);
0412         for (const auto &j : HCALSpecialsNames) {
0413           ibooker.setCurrentFolder(dqmFolder + "/" + j.second);
0414           H2mod_.back()[label][j.first] = ibooker.book2D(label, label, btagBins, btagL, btagU, nBinsPt, pTmin, pTMax);
0415           H2mod_.back()[label][j.first]->setAxisTitle("pT", 2);
0416           H2mod_.back()[label][j.first]->setAxisTitle("disc", 1);
0417         }
0418         ibooker.setCurrentFolder(dqmFolder);
0419         H2Eta_.back()[labelEta] = ibooker.book2D(labelEta, labelEta, btagBins, btagL, btagU, nBinsEta, etamin, etaMax);
0420         H2Eta_.back()[labelEta]->setAxisTitle("eta", 2);
0421         H2Eta_.back()[labelEta]->setAxisTitle("disc", 1);
0422         H2Phi_.back()[labelPhi] = ibooker.book2D(labelPhi, labelPhi, btagBins, btagL, btagU, nBinsPhi, phimin, phiMax);
0423         H2Phi_.back()[labelPhi]->setAxisTitle("phi", 2);
0424         H2Phi_.back()[labelPhi]->setAxisTitle("disc", 1);
0425         H2EtaPhi_.back()[labelEtaPhi] =
0426             ibooker.book2D(labelEtaPhi, labelEtaPhi, nBinsEta, etamin, etaMax, nBinsPhi, phimin, phiMax);
0427         H2EtaPhi_.back()[labelEtaPhi]->setAxisTitle("phi", 2);
0428         H2EtaPhi_.back()[labelEtaPhi]->setAxisTitle("eta", 1);
0429         H2EtaPhi_threshold_.back()[labelEtaPhi_threshold] = ibooker.book2D(
0430             labelEtaPhi_threshold, labelEtaPhi_threshold, nBinsEta, etamin, etaMax, nBinsPhi, phimin, phiMax);
0431         H2EtaPhi_threshold_.back()[labelEtaPhi_threshold]->setAxisTitle("phi", 2);
0432         H2EtaPhi_threshold_.back()[labelEtaPhi_threshold]->setAxisTitle("eta", 1);
0433       }
0434     }  /// for mc.size()
0435   }    /// for hltPathNames_.size()
0436 }
0437 
0438 // define this as a plug-in
0439 DEFINE_FWK_MODULE(HLTBTagPerformanceAnalyzer);