Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 ///////////////////////////////////////////////////////////////////////////////
0002 //                    Header file for this                                    //
0003 ////////////////////////////////////////////////////////////////////////////////
0004 #include "HLTriggerOffline/Egamma/interface/EmDQM.h"
0005 
0006 using namespace ROOT::Math::VectorUtil;
0007 
0008 ////////////////////////////////////////////////////////////////////////////////
0009 //                             Constructor                                    //
0010 ////////////////////////////////////////////////////////////////////////////////
0011 EmDQM::EmDQM(const edm::ParameterSet &pset) {
0012   // are we running in automatic configuration mode with the HLTConfigProvider
0013   // or with a per path python config file
0014   autoConfMode_ = pset.getUntrackedParameter<bool>("autoConfMode", false);
0015 
0016   // set global parameters
0017   triggerObject_ = pset.getParameter<edm::InputTag>("triggerobject");
0018   verbosity_ = pset.getUntrackedParameter<unsigned int>("verbosity", 0);
0019   genEtaAcc_ = pset.getParameter<double>("genEtaAcc");
0020   genEtAcc_ = pset.getParameter<double>("genEtAcc");
0021   ptMax_ = pset.getUntrackedParameter<double>("PtMax", 1000.);
0022   ptMin_ = pset.getUntrackedParameter<double>("PtMin", 0.);
0023   etaMax_ = pset.getUntrackedParameter<double>("EtaMax", 2.7);
0024   phiMax_ = pset.getUntrackedParameter<double>("PhiMax", 3.15);
0025   nbins_ = pset.getUntrackedParameter<unsigned int>("Nbins", 40);
0026   eta2DMax_ = pset.getUntrackedParameter<double>("Eta2DMax", 2.8);
0027   phi2DMax_ = pset.getUntrackedParameter<double>("Phi2DMax", 3.2);
0028   nbins2D_ = pset.getUntrackedParameter<unsigned int>("Nbins2D", 16);
0029   minEtForEtaEffPlot_ = pset.getUntrackedParameter<unsigned int>("minEtForEtaEffPlot", 15);
0030   useHumanReadableHistTitles_ = pset.getUntrackedParameter<bool>("useHumanReadableHistTitles", false);
0031   mcMatchedOnly_ = pset.getUntrackedParameter<bool>("mcMatchedOnly", true);
0032   noPhiPlots_ = pset.getUntrackedParameter<bool>("noPhiPlots", true);
0033   noIsolationPlots_ = pset.getUntrackedParameter<bool>("noIsolationPlots", true);
0034 
0035   if (!autoConfMode_) {
0036     paramSets.push_back(pset);
0037     isData_ = false;
0038   } else {
0039     isData_ = pset.getParameter<bool>("isData");
0040   }
0041 
0042   histoFillerEle = new HistoFiller<reco::ElectronCollection>(this);
0043   histoFillerPho = new HistoFiller<reco::RecoEcalCandidateCollection>(this);
0044   histoFillerClu = new HistoFiller<reco::RecoEcalCandidateCollection>(this);
0045   histoFillerL1Iso = new HistoFiller<l1extra::L1EmParticleCollection>(this);
0046   histoFillerL1NonIso = new HistoFiller<l1extra::L1EmParticleCollection>(this);
0047 
0048   // consumes
0049   genParticles_token = consumes<edm::View<reco::Candidate>>(edm::InputTag("genParticles"));
0050   triggerObject_token = consumes<trigger::TriggerEventWithRefs>(triggerObject_);
0051   hltResults_token = consumes<edm::TriggerResults>(edm::InputTag("TriggerResults", "", triggerObject_.process()));
0052   if (autoConfMode_) {
0053     gencutColl_fidWenu_token = mayConsume<edm::View<reco::Candidate>>(edm::InputTag("fiducialWenu"));
0054     gencutColl_fidZee_token = mayConsume<edm::View<reco::Candidate>>(edm::InputTag("fiducialZee"));
0055     gencutColl_fidTripleEle_token = mayConsume<edm::View<reco::Candidate>>(edm::InputTag("fiducialTripleEle"));
0056     gencutColl_fidGammaJet_token = mayConsume<edm::View<reco::Candidate>>(edm::InputTag("fiducialGammaJet"));
0057     gencutColl_fidDiGamma_token = mayConsume<edm::View<reco::Candidate>>(edm::InputTag("fiducialDiGamma"));
0058   } else {
0059     gencutColl_manualConf_token =
0060         consumes<edm::View<reco::Candidate>>(pset.getParameter<edm::InputTag>("cutcollection"));
0061   }
0062 }
0063 
0064 void EmDQM::dqmBeginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) {
0065   bool changed(true);
0066   if (hltConfig_.init(iRun, iSetup, triggerObject_.process(), changed)) {
0067     // if init returns TRUE, initialisation has succeeded!
0068 
0069     if (autoConfMode_) {
0070       if (verbosity_ >= OUTPUT_ALL) {
0071         // Output general information on the menu
0072         edm::LogPrint("EmDQM") << "inited=" << hltConfig_.inited();
0073         edm::LogPrint("EmDQM") << "changed=" << hltConfig_.changed();
0074         edm::LogPrint("EmDQM") << "processName=" << hltConfig_.processName();
0075         edm::LogPrint("EmDQM") << "tableName=" << hltConfig_.tableName();
0076         edm::LogPrint("EmDQM") << "size=" << hltConfig_.size();
0077         edm::LogInfo("EmDQM") << "The following filter types are not analyzed: \n"
0078                               << "\tHLTGlobalSumsMET\n"
0079                               << "\tHLTHtMhtFilter\n"
0080                               << "\tHLTMhtFilter\n"
0081                               << "\tHLTJetTag\n"
0082                               << "\tHLT1CaloJet\n"
0083                               << "\tHLT1CaloMET\n"
0084                               << "\tHLT1CaloBJet\n"
0085                               << "\tHLT1Tau\n"
0086                               << "\tHLT1PFTau\n"
0087                               << "\tPFTauSelector\n"
0088                               << "\tHLT1PFJet\n"
0089                               << "\tHLTPFJetCollectionsFilter\n"
0090                               << "\tHLTPFJetCollectionsVBFFilter\n"
0091                               << "\tHLTPFJetTag\n"
0092                               << "\tEtMinCaloJetSelector\n"
0093                               << "\tEtMinPFJetSelector\n"
0094                               << "\tLargestEtCaloJetSelector\n"
0095                               << "\tLargestEtPFJetSelector\n"
0096                               << "\tHLTEgammaTriggerFilterObjectWrapper\n"
0097                               << "\tHLTEgammaDoubleLegCombFilter\n"
0098                               << "\tHLT2ElectronTau\n"
0099                               << "\tHLT2ElectronMET\n"
0100                               << "\tHLT2ElectronPFTau\n"
0101                               << "\tHLTPMMassFilter\n"
0102                               << "\tHLTHcalTowerFilter\n"
0103                               << "\tHLT1Photon\n"
0104                               << "\tHLTRFilter\n"
0105                               << "\tHLTRHemisphere\n"
0106                               << "\tHLTElectronPFMTFilter\n"
0107                               << "\tPrimaryVertexObjectFilter\n"
0108                               << "\tHLTEgammaAllCombMassFilter\n"
0109                               << "\tHLTMuon*\n";
0110       }
0111 
0112       // All electron and photon paths
0113       std::vector<std::vector<std::string>> egammaPaths = findEgammaPaths();
0114       // std::cout << "Found " << egammaPaths[TYPE_SINGLE_ELE].size() << "
0115       // single electron paths" << std::endl; std::cout << "Found " <<
0116       // egammaPaths[TYPE_DOUBLE_ELE].size() << " double electron paths" <<
0117       // std::endl; std::cout << "Found " << egammaPaths[TYPE_TRIPLE_ELE].size()
0118       // << " triple electron paths" << std::endl; std::cout << "Found " <<
0119       // egammaPaths[TYPE_SINGLE_PHOTON].size() << " single photon paths" <<
0120       // std::endl; std::cout << "Found " <<
0121       // egammaPaths[TYPE_DOUBLE_PHOTON].size() << " double photon paths" <<
0122       // std::endl;
0123 
0124       std::vector<std::string> filterModules;
0125 
0126       for (unsigned int j = 0; j < egammaPaths.size(); j++) {
0127         if (verbosity_ >= OUTPUT_ALL) {
0128           switch (j) {
0129             case TYPE_SINGLE_ELE:
0130               edm::LogPrint("EmDQM") << "////////////////////////////////////////"
0131                                         "/\nSingle electron paths: ";
0132               break;
0133             case TYPE_DOUBLE_ELE:
0134               edm::LogPrint("EmDQM") << "////////////////////////////////////////"
0135                                         "/\nDouble electron paths: ";
0136               break;
0137             case TYPE_TRIPLE_ELE:
0138               edm::LogPrint("EmDQM") << "////////////////////////////////////////"
0139                                         "/\nTriple electron paths: ";
0140               break;
0141             case TYPE_SINGLE_PHOTON:
0142               edm::LogPrint("EmDQM") << "////////////////////////////////////////"
0143                                         "/\nSingle photon paths: ";
0144               break;
0145             case TYPE_DOUBLE_PHOTON:
0146               edm::LogPrint("EmDQM") << "////////////////////////////////////////"
0147                                         "/\nDouble photon paths: ";
0148               break;
0149           }
0150         }
0151 
0152         for (unsigned int i = 0; i < egammaPaths.at(j).size(); i++) {
0153           // get pathname of this trigger
0154           const std::string pathName = egammaPaths.at(j).at(i);
0155           if (verbosity_ >= OUTPUT_ALL)
0156             edm::LogPrint("EmDQM") << "Path: " << pathName;
0157 
0158           // get filters of the current path
0159           filterModules = getFilterModules(pathName);
0160 
0161           //--------------------
0162           edm::ParameterSet paramSet;
0163 
0164           paramSet.addUntrackedParameter("pathIndex", hltConfig_.triggerIndex(pathName));
0165           paramSet.addParameter("@module_label", hltConfig_.removeVersion(pathName) + "_DQM");
0166           // paramSet.addParameter("@module_label", pathName + "_DQM");
0167 
0168           // plotting parameters (untracked because they don't affect the
0169           // physics)
0170           double genEtMin = getPrimaryEtCut(pathName);
0171           if (genEtMin >= 0) {
0172             paramSet.addUntrackedParameter("genEtMin", genEtMin);
0173           } else {
0174             if (verbosity_ >= OUTPUT_WARNINGS)
0175               edm::LogWarning("EmDQM") << "Pathname: '" << pathName
0176                                        << "':  Unable to determine a minimum Et. Will not include "
0177                                           "this path in the validation.";
0178             continue;
0179           }
0180 
0181           // set the x axis of the et plots to some reasonable value based
0182           // on the primary et cut determined from the path name
0183           double ptMax = ptMax_;
0184           double ptMin = ptMin_;
0185           if (ptMax < (1.2 * genEtMin)) {
0186             paramSet.addUntrackedParameter<double>("PtMax", (10 * ceil(0.12 * genEtMin)));
0187             paramSet.addUntrackedParameter<double>("PtMin", (10 * ceil(0.12 * genEtMin) - ptMax + ptMin));
0188           } else {
0189             paramSet.addUntrackedParameter<double>("PtMax", ptMax);
0190             paramSet.addUntrackedParameter<double>("PtMin", ptMin);
0191           }
0192 
0193           // preselction cuts
0194           switch (j) {
0195             case TYPE_SINGLE_ELE:
0196               paramSet.addParameter<unsigned>("reqNum", 1);
0197               paramSet.addParameter<int>("pdgGen", 11);
0198               paramSet.addParameter<edm::InputTag>("cutcollection", edm::InputTag("fiducialWenu"));
0199               paramSet.addParameter<int>("cutnum", 1);
0200               break;
0201             case TYPE_DOUBLE_ELE:
0202               paramSet.addParameter<unsigned>("reqNum", 2);
0203               paramSet.addParameter<int>("pdgGen", 11);
0204               paramSet.addParameter<edm::InputTag>("cutcollection", edm::InputTag("fiducialZee"));
0205               paramSet.addParameter<int>("cutnum", 2);
0206               break;
0207             case TYPE_TRIPLE_ELE:
0208               paramSet.addParameter<unsigned>("reqNum", 3);
0209               paramSet.addParameter<int>("pdgGen", 11);
0210               paramSet.addParameter<edm::InputTag>("cutcollection", edm::InputTag("fiducialTripleEle"));
0211               paramSet.addParameter<int>("cutnum", 3);
0212               break;
0213             case TYPE_SINGLE_PHOTON:
0214               paramSet.addParameter<unsigned>("reqNum", 1);
0215               paramSet.addParameter<int>("pdgGen", 22);
0216               paramSet.addParameter<edm::InputTag>("cutcollection", edm::InputTag("fiducialGammaJet"));
0217               paramSet.addParameter<int>("cutnum", 1);
0218               break;
0219             case TYPE_DOUBLE_PHOTON:
0220               paramSet.addParameter<unsigned>("reqNum", 2);
0221               paramSet.addParameter<int>("pdgGen", 22);
0222               paramSet.addParameter<edm::InputTag>("cutcollection", edm::InputTag("fiducialDiGamma"));
0223               paramSet.addParameter<int>("cutnum", 2);
0224           }
0225           //--------------------
0226 
0227           // TODO: extend this
0228           if (isData_)
0229             paramSet.addParameter<edm::InputTag>("cutcollection", edm::InputTag("gsfElectrons"));
0230 
0231           std::vector<edm::ParameterSet> filterVPSet;
0232           edm::ParameterSet filterPSet;
0233           std::string moduleLabel;
0234 
0235           // loop over filtermodules of current trigger path
0236           for (std::vector<std::string>::iterator filter = filterModules.begin(); filter != filterModules.end();
0237                ++filter) {
0238             std::string moduleType = hltConfig_.modulePSet(*filter).getParameter<std::string>("@module_type");
0239             moduleLabel = hltConfig_.modulePSet(*filter).getParameter<std::string>("@module_label");
0240 
0241             // first check if it is a filter we are not interrested in
0242             if (moduleType == "Pythia6GeneratorFilter" || moduleType == "HLTTriggerTypeFilter" ||
0243                 moduleType == "HLTLevel1Activity" || moduleType == "HLTPrescaler" || moduleType == "HLTBool")
0244               continue;
0245 
0246             // now check for the known filter types
0247             if (moduleType == "HLTLevel1GTSeed") {
0248               filterPSet = makePSetForL1SeedFilter(moduleLabel);
0249             } else if (moduleType == "HLTEgammaL1MatchFilterRegional") {
0250               filterPSet = makePSetForL1SeedToSuperClusterMatchFilter(moduleLabel);
0251             } else if (moduleType == "HLTEgammaEtFilter") {
0252               filterPSet = makePSetForEtFilter(moduleLabel);
0253             } else if (moduleType == "HLTElectronOneOEMinusOneOPFilterRegional") {
0254               filterPSet = makePSetForOneOEMinusOneOPFilter(moduleLabel);
0255             } else if (moduleType == "HLTElectronPixelMatchFilter") {
0256               filterPSet = makePSetForPixelMatchFilter(moduleLabel);
0257             } else if (moduleType == "HLTEgammaGenericFilter") {
0258               filterPSet = makePSetForEgammaGenericFilter(moduleLabel);
0259             } else if (moduleType == "HLTEgammaGenericQuadraticFilter") {
0260               filterPSet = makePSetForEgammaGenericQuadraticFilter(moduleLabel);
0261             } else if (moduleType == "HLTElectronGenericFilter") {
0262               filterPSet = makePSetForElectronGenericFilter(moduleLabel);
0263             } else if (moduleType == "HLTEgammaDoubleEtDeltaPhiFilter") {
0264               filterPSet = makePSetForEgammaDoubleEtDeltaPhiFilter(moduleLabel);
0265             } else if (moduleType == "HLTGlobalSumsMET" || moduleType == "HLTHtMhtFilter" ||
0266                        moduleType == "HLTMhtFilter" || moduleType == "HLTJetTag" || moduleType == "HLT1CaloJet" ||
0267                        moduleType == "HLT1CaloMET" || moduleType == "HLT1CaloBJet" || moduleType == "HLT1Tau" ||
0268                        moduleType == "HLT1PFTau" || moduleType == "PFTauSelector" || moduleType == "HLT1PFJet" ||
0269                        moduleType == "HLTPFJetCollectionsFilter" || moduleType == "HLTPFJetCollectionsVBFFilter" ||
0270                        moduleType == "HLTPFJetTag" || moduleType == "EtMinCaloJetSelector" ||
0271                        moduleType == "EtMinPFJetSelector" || moduleType == "LargestEtCaloJetSelector" ||
0272                        moduleType == "LargestEtPFJetSelector" ||
0273                        moduleType == "HLTEgammaTriggerFilterObjectWrapper"  // 'fake' filter
0274                        || moduleType == "HLTEgammaDoubleLegCombFilter"      // filter does not put
0275                                                                             // anything in
0276                                                                             // TriggerEventWithRefs
0277                        || moduleType == "HLT2ElectronMET" || moduleType == "HLT2ElectronTau" ||
0278                        moduleType == "HLT2ElectronPFTau" || moduleType == "HLTPMMassFilter" ||
0279                        moduleType == "HLTHcalTowerFilter" || moduleType == "HLT1Photon" || moduleType == "HLTRFilter" ||
0280                        moduleType == "HLTRHemisphere" || moduleType == "HLTElectronPFMTFilter" ||
0281                        moduleType == "PrimaryVertexObjectFilter" || moduleType == "HLTEgammaAllCombMassFilter" ||
0282                        moduleType.find("HLTMuon") != std::string::npos)
0283               continue;
0284             else {
0285               if (verbosity_ >= OUTPUT_WARNINGS)
0286                 edm::LogWarning("EmDQM") << "No parameter set for filter '" << moduleLabel << "' with filter type '"
0287                                          << moduleType << "' added. Module will not be analyzed.";
0288               continue;
0289             }
0290 
0291             // something went wrong when the parameter set is empty.
0292             if (!filterPSet.empty()) {
0293               if (!hltConfig_.modulePSet(moduleLabel).exists("saveTags")) {
0294                 // make sure that 'theHLTOutputTypes' before an unseeded filter
0295                 // in a photon path is set to trigger::TriggerPhoton this is
0296                 // coupled to the parameter 'saveTag = true'
0297                 if (moduleLabel.find("Unseeded") != std::string::npos &&
0298                     (j == TYPE_DOUBLE_PHOTON || j == TYPE_SINGLE_PHOTON)) {
0299                   filterVPSet.back().addParameter<int>("theHLTOutputTypes", trigger::TriggerPhoton);
0300                 }
0301               }
0302               // if ncandcut is -1 (when parsing for the number of particles in
0303               // the name of the L1seed filter fails), fall back to setting
0304               // ncandcut to the number of particles needed for the given path.
0305               if (filterPSet.getParameter<int>("ncandcut") < 0) {
0306                 edm::LogPrint("EmDQM") << "No number of candidates for filter " << moduleLabel << " found. Set to "
0307                                        << paramSet.getParameter<int>("cutnum") << ", determined from path name.";
0308                 filterPSet.addParameter<int>("ncandcut", paramSet.getParameter<int>("cutnum"));
0309               } else if (filterPSet.getParameter<int>("ncandcut") > paramSet.getParameter<int>("cutnum")) {
0310                 edm::LogInfo("EmDQM") << "Changed required number of candidates from "
0311                                       << paramSet.getParameter<int>("cutnum") << " to "
0312                                       << filterPSet.getParameter<int>("ncandcut") << " for filter " << moduleLabel;
0313                 paramSet.addParameter<int>("cutnum", filterPSet.getParameter<int>("ncandcut"));
0314                 paramSet.addParameter<unsigned>("reqNum", (unsigned)filterPSet.getParameter<int>("ncandcut"));
0315               }
0316 
0317               filterVPSet.push_back(filterPSet);
0318             } else
0319               break;
0320 
0321           }  // end loop over filter modules of current trigger path
0322 
0323           // do not include this path when an empty filterPSet is detected.
0324           if (!filterPSet.empty()) {
0325             std::string lastModuleName = filterPSet.getParameter<edm::InputTag>("HLTCollectionLabels").label();
0326             if (!hltConfig_.modulePSet(lastModuleName).exists("saveTags")) {
0327               // make sure that 'theHLTOutputTypes' of the last filter of a
0328               // photon path is set to trigger::TriggerPhoton this is coupled to
0329               // the parameter 'saveTag = true'
0330               if ((j == TYPE_SINGLE_PHOTON || j == TYPE_DOUBLE_PHOTON) && pathName.rfind("Ele") == std::string::npos) {
0331                 filterVPSet.back().addParameter<int>("theHLTOutputTypes", trigger::TriggerPhoton);
0332               }
0333             }
0334             paramSet.addParameter<std::vector<edm::ParameterSet>>("filters", filterVPSet);
0335           } else {
0336             if (verbosity_ >= OUTPUT_ALL)
0337               edm::LogPrint("EmDQM") << "Will not include this path in the validation due to "
0338                                         "errors while generating the parameter set.";
0339             continue;
0340           }
0341 
0342           // dump generated parameter set
0343           // std::cout << paramSet.dump() << std::endl;
0344 
0345           paramSets.push_back(paramSet);
0346         }  // loop over all paths of this analysis type
0347 
0348       }  // loop over analysis types (single ele etc.)
0349     }
0350 
0351     hltCollectionLabelsFoundPerPath.reserve(paramSets.size());
0352     hltCollectionLabelsMissedPerPath.reserve(paramSets.size());
0353 
0354     ////////////////////////////////////////////////////////////
0355     // loop over all the trigger path parameter sets
0356     ////////////////////////////////////////////////////////////
0357     for (std::vector<edm::ParameterSet>::iterator psetIt = paramSets.begin(); psetIt != paramSets.end(); ++psetIt) {
0358       hltCollectionLabelsFoundPerPath.push_back(hltCollectionLabelsFound);
0359       hltCollectionLabelsMissedPerPath.push_back(hltCollectionLabelsMissed);
0360     }
0361 
0362     if (changed) {
0363       // The HLT config has actually changed wrt the previous Run, hence rebook
0364       // your histograms or do anything else dependent on the revised HLT config
0365     }
0366   } else {
0367     // if init returns FALSE, initialisation has NOT succeeded, which indicates
0368     // a problem with the file and/or code and needs to be investigated!
0369     if (verbosity_ >= OUTPUT_ERRORS)
0370       edm::LogError("EmDQM") << " HLT config extraction failure with process name '" << triggerObject_.process()
0371                              << "'.";
0372     // In this case, all access methods will return empty values!
0373   }
0374 }
0375 
0376 void EmDQM::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &iRun, edm::EventSetup const &iSetup) {
0377   ////////////////////////////////////////////////////////////
0378   // loop over all the trigger path parameter sets
0379   ////////////////////////////////////////////////////////////
0380   for (std::vector<edm::ParameterSet>::iterator psetIt = paramSets.begin(); psetIt != paramSets.end(); ++psetIt) {
0381     SetVarsFromPSet(psetIt);
0382 
0383     iBooker.setCurrentFolder(dirname_);
0384 
0385     ////////////////////////////////////////////////////////////
0386     //  Set up Histogram of Effiency vs Step.                 //
0387     //   theHLTCollectionLabels is a vector of InputTags      //
0388     //    from the configuration file.                        //
0389     ////////////////////////////////////////////////////////////
0390     // Et & eta distributions
0391     std::vector<MonitorElement *> etahist;
0392     std::vector<MonitorElement *> phihist;
0393     std::vector<MonitorElement *> ethist;
0394     std::vector<MonitorElement *> etahistmatch;
0395     std::vector<MonitorElement *> phihistmatch;
0396     std::vector<MonitorElement *> ethistmatch;
0397     std::vector<MonitorElement *> histEtOfHltObjMatchToGen;
0398     std::vector<MonitorElement *> histEtaOfHltObjMatchToGen;
0399     std::vector<MonitorElement *> histPhiOfHltObjMatchToGen;
0400     std::vector<MonitorElement *> etaphihist;
0401     std::vector<MonitorElement *> etaphihistmatch;
0402     std::vector<MonitorElement *> histEtaPhiOfHltObjMatchToGen;
0403     // Plots of efficiency per step
0404     MonitorElement *total;
0405     MonitorElement *totalmatch;
0406     // generator histograms
0407     MonitorElement *etgen;
0408     MonitorElement *etagen;
0409     MonitorElement *phigen;
0410     MonitorElement *etaphigen;
0411 
0412     std::string histName = "total_eff";
0413     std::string histTitle = "total events passing";
0414     if (!mcMatchedOnly_) {
0415       // This plot will have bins equal to 2+(number of
0416       //        HLTCollectionLabels in the config file)
0417       total = iBooker.book1D(
0418           histName.c_str(), histTitle.c_str(), numOfHLTCollectionLabels + 2, 0, numOfHLTCollectionLabels + 2);
0419       total->setBinLabel(numOfHLTCollectionLabels + 1, "Total");
0420       total->setBinLabel(numOfHLTCollectionLabels + 2, "Gen");
0421       for (unsigned int u = 0; u < numOfHLTCollectionLabels; u++) {
0422         total->setBinLabel(u + 1, theHLTCollectionLabels[u].label());
0423       }
0424     }
0425 
0426     histName = "total_eff_MC_matched";
0427     histTitle = "total events passing (mc matched)";
0428     totalmatch = iBooker.book1D(
0429         histName.c_str(), histTitle.c_str(), numOfHLTCollectionLabels + 2, 0, numOfHLTCollectionLabels + 2);
0430     totalmatch->setBinLabel(numOfHLTCollectionLabels + 1, "Total");
0431     totalmatch->setBinLabel(numOfHLTCollectionLabels + 2, "Gen");
0432     for (unsigned int u = 0; u < numOfHLTCollectionLabels; u++) {
0433       totalmatch->setBinLabel(u + 1, theHLTCollectionLabels[u].label());
0434     }
0435 
0436     MonitorElement *tmphisto;
0437     // MonitorElement* tmpiso;
0438 
0439     ////////////////////////////////////////////////////////////
0440     // Set up generator-level histograms                      //
0441     ////////////////////////////////////////////////////////////
0442     std::string pdgIdString;
0443     switch (pdgGen) {
0444       case 11:
0445         pdgIdString = "Electron";
0446         break;
0447       case 22:
0448         pdgIdString = "Photon";
0449         break;
0450       default:
0451         pdgIdString = "Particle";
0452     }
0453 
0454     histName = "gen_et";
0455     histTitle = "E_{T} of " + pdgIdString + "s";
0456     etgen = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, plotPtMin, plotPtMax);
0457     histName = "gen_eta";
0458     histTitle = "#eta of " + pdgIdString + "s ";
0459     etagen = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -etaMax_, etaMax_);
0460     histName = "gen_phi";
0461     histTitle = "#phi of " + pdgIdString + "s ";
0462     if (!noPhiPlots_)
0463       phigen = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -phiMax_, phiMax_);
0464     histName = "gen_etaphi";
0465     histTitle = "#eta-#phi of " + pdgIdString + "s ";
0466     etaphigen = iBooker.book2D(
0467         histName.c_str(), histTitle.c_str(), nbins2D_ - 2, -eta2DMax_, eta2DMax_, nbins2D_, -phi2DMax_, phi2DMax_);
0468     ////////////////////////////////////////////////////////////
0469     //  Set up histograms of HLT objects                      //
0470     ////////////////////////////////////////////////////////////
0471     // Determine what strings to use for histogram titles
0472     std::vector<std::string> HltHistTitle;
0473     if (theHLTCollectionHumanNames.size() == numOfHLTCollectionLabels && useHumanReadableHistTitles_) {
0474       HltHistTitle = theHLTCollectionHumanNames;
0475     } else {
0476       for (unsigned int i = 0; i < numOfHLTCollectionLabels; i++) {
0477         HltHistTitle.push_back(theHLTCollectionLabels[i].label());
0478       }
0479     }
0480 
0481     for (unsigned int i = 0; i < numOfHLTCollectionLabels; i++) {
0482       if (!mcMatchedOnly_) {
0483         // Et distribution of HLT objects passing filter i
0484         histName = theHLTCollectionLabels[i].label() + "et_all";
0485         histTitle = HltHistTitle[i] + " Et (ALL)";
0486         tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, plotPtMin, plotPtMax);
0487         ethist.push_back(tmphisto);
0488 
0489         // Eta distribution of HLT objects passing filter i
0490         histName = theHLTCollectionLabels[i].label() + "eta_all";
0491         histTitle = HltHistTitle[i] + " #eta (ALL)";
0492         tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -etaMax_, etaMax_);
0493         etahist.push_back(tmphisto);
0494 
0495         if (!noPhiPlots_) {
0496           // Phi distribution of HLT objects passing filter i
0497           histName = theHLTCollectionLabels[i].label() + "phi_all";
0498           histTitle = HltHistTitle[i] + " #phi (ALL)";
0499           tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -phiMax_, phiMax_);
0500           phihist.push_back(tmphisto);
0501         }
0502 
0503         histName = theHLTCollectionLabels[i].label() + "etaphi_all";
0504         histTitle = HltHistTitle[i] + " #eta-#phi (ALL)";
0505         tmphisto = iBooker.book2D(
0506             histName.c_str(), histTitle.c_str(), nbins2D_ - 2, -eta2DMax_, eta2DMax_, nbins2D_, -phi2DMax_, phi2DMax_);
0507         etaphihist.push_back(tmphisto);
0508 
0509         // Et distribution of HLT object that is closest delta-R match to sorted
0510         // gen particle(s)
0511         histName = theHLTCollectionLabels[i].label() + "et";
0512         histTitle = HltHistTitle[i] + " Et";
0513         tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, plotPtMin, plotPtMax);
0514         histEtOfHltObjMatchToGen.push_back(tmphisto);
0515 
0516         // eta distribution of HLT object that is closest delta-R match to
0517         // sorted gen particle(s)
0518         histName = theHLTCollectionLabels[i].label() + "eta";
0519         histTitle = HltHistTitle[i] + " eta";
0520         tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -etaMax_, etaMax_);
0521         histEtaOfHltObjMatchToGen.push_back(tmphisto);
0522 
0523         if (!noPhiPlots_) {
0524           // phi distribution of HLT object that is closest delta-R match to
0525           // sorted gen particle(s)
0526           histName = theHLTCollectionLabels[i].label() + "phi";
0527           histTitle = HltHistTitle[i] + " phi";
0528           tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -phiMax_, phiMax_);
0529           histPhiOfHltObjMatchToGen.push_back(tmphisto);
0530         }
0531 
0532         histName = theHLTCollectionLabels[i].label() + "etaphi";
0533         histTitle = HltHistTitle[i] + " eta-phi";
0534         tmphisto = iBooker.book2D(
0535             histName.c_str(), histTitle.c_str(), nbins2D_ - 2, -eta2DMax_, eta2DMax_, nbins2D_, -phi2DMax_, phi2DMax_);
0536         histEtaPhiOfHltObjMatchToGen.push_back(tmphisto);
0537       }
0538 
0539       // Et distribution of gen object matching HLT object passing filter i
0540       histName = theHLTCollectionLabels[i].label() + "et_MC_matched";
0541       histTitle = HltHistTitle[i] + " Et (MC matched)";
0542       tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, plotPtMin, plotPtMax);
0543       ethistmatch.push_back(tmphisto);
0544 
0545       // Eta distribution of gen object matching HLT object passing filter i
0546       histName = theHLTCollectionLabels[i].label() + "eta_MC_matched";
0547       histTitle = HltHistTitle[i] + " #eta (MC matched)";
0548       tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -etaMax_, etaMax_);
0549       etahistmatch.push_back(tmphisto);
0550 
0551       if (!noPhiPlots_) {
0552         // Phi distribution of gen object matching HLT object passing filter i
0553         histName = theHLTCollectionLabels[i].label() + "phi_MC_matched";
0554         histTitle = HltHistTitle[i] + " #phi (MC matched)";
0555         tmphisto = iBooker.book1D(histName.c_str(), histTitle.c_str(), nbins_, -phiMax_, phiMax_);
0556         phihistmatch.push_back(tmphisto);
0557       }
0558 
0559       histName = theHLTCollectionLabels[i].label() + "etaphi_MC_matched";
0560       histTitle = HltHistTitle[i] + " #eta-#phi (MC matched)";
0561       tmphisto = iBooker.book2D(
0562           histName.c_str(), histTitle.c_str(), nbins2D_ - 2, -eta2DMax_, eta2DMax_, nbins2D_, -phi2DMax_, phi2DMax_);
0563       etaphihistmatch.push_back(tmphisto);
0564     }
0565 
0566     // Et & eta distributions
0567     etahists.push_back(etahist);
0568     phihists.push_back(phihist);
0569     ethists.push_back(ethist);
0570     etahistmatchs.push_back(etahistmatch);
0571     phihistmatchs.push_back(phihistmatch);
0572     ethistmatchs.push_back(ethistmatch);
0573     histEtOfHltObjMatchToGens.push_back(histEtOfHltObjMatchToGen);
0574     histEtaOfHltObjMatchToGens.push_back(histEtaOfHltObjMatchToGen);
0575     histPhiOfHltObjMatchToGens.push_back(histPhiOfHltObjMatchToGen);
0576     etaphihists.push_back(etaphihist);
0577     etaphihistmatchs.push_back(etaphihistmatch);
0578     histEtaPhiOfHltObjMatchToGens.push_back(histEtaPhiOfHltObjMatchToGen);
0579     // commented out because uses data not included in HTLDEBUG and uses
0580     // Isolation distributions
0581     // etahistisos.push_back(etahistiso);
0582     // phihistisos.push_back(phihistiso);
0583     // ethistisos.push_back(ethistiso);
0584     // etahistisomatchs.push_back(etahistisomatch);
0585     // phihistisomatchs.push_back(phihistisomatch);
0586     // ethistisomatchs.push_back(ethistisomatch);
0587     // histEtIsoOfHltObjMatchToGens.push_back(histEtIsoOfHltObjMatchToGen);
0588     // histEtaIsoOfHltObjMatchToGens.push_back(histEtaIsoOfHltObjMatchToGen);
0589     // histPhiIsoOfHltObjMatchToGens.push_back(histPhiIsoOfHltObjMatchToGen);
0590 
0591     totals.push_back(total);
0592     totalmatchs.push_back(totalmatch);
0593     etgens.push_back(etgen);
0594     etagens.push_back(etagen);
0595     phigens.push_back(phigen);
0596     etaphigens.push_back(etaphigen);
0597   }
0598 }
0599 
0600 ////////////////////////////////////////////////////////////////////////////////
0601 //                                Destructor                                  //
0602 ////////////////////////////////////////////////////////////////////////////////
0603 EmDQM::~EmDQM() {}
0604 ////////////////////////////////////////////////////////////////////////////////
0605 
0606 bool EmDQM::checkGeneratedParticlesRequirement(const edm::Event &event) {
0607   ////////////////////////////////////////////////////////////
0608   // Decide if this was an event of interest.               //
0609   //  Did the highest energy particles happen               //
0610   //  to have |eta| < 2.5 ?  Then continue.                 //
0611   ////////////////////////////////////////////////////////////
0612   edm::Handle<edm::View<reco::Candidate>> genParticles;
0613   event.getByToken(genParticles_token, genParticles);
0614   if (!genParticles.isValid()) {
0615     if (verbosity_ >= OUTPUT_WARNINGS)
0616       edm::LogWarning("EmDQM") << "genParticles invalid.";
0617     return false;
0618   }
0619 
0620   std::vector<reco::LeafCandidate> allSortedGenParticles;
0621 
0622   for (edm::View<reco::Candidate>::const_iterator currentGenParticle = genParticles->begin();
0623        currentGenParticle != genParticles->end();
0624        currentGenParticle++) {
0625     // TODO: do we need to check the states here again ?
0626     // in principle, there should collections produced with the python
0627     // configuration (other than 'genParticles') which fulfill these criteria
0628     if (!(abs((*currentGenParticle).pdgId()) == pdgGen && (*currentGenParticle).status() == 1 &&
0629           (*currentGenParticle).et() > 2.0))
0630       continue;
0631 
0632     reco::LeafCandidate tmpcand(*(currentGenParticle));
0633 
0634     if (tmpcand.et() < plotEtMin)
0635       continue;
0636 
0637     allSortedGenParticles.push_back(tmpcand);
0638   }
0639 
0640   std::sort(allSortedGenParticles.begin(), allSortedGenParticles.end(), pTGenComparator_);
0641 
0642   // return false if not enough particles found
0643   if (allSortedGenParticles.size() < gencut_)
0644     return false;
0645 
0646   // additional check (this might be legacy code and we need to check
0647   // whether this should not be removed ?)
0648 
0649   // We now have a sorted collection of all generated particles
0650   // with pdgId = pdgGen.
0651   // Loop over them to see if the top gen particles have eta within acceptance
0652   // bool keepEvent = true;
0653   for (unsigned int i = 0; i < gencut_; i++) {
0654     bool inECALgap = fabs(allSortedGenParticles[i].eta()) > 1.4442 && fabs(allSortedGenParticles[i].eta()) < 1.556;
0655     if ((fabs(allSortedGenParticles[i].eta()) > genEtaAcc_) || inECALgap) {
0656       // edm::LogWarning("EmDQM") << "Throwing event away. Gen particle with
0657       // pdgId="<< allSortedGenParticles[i].pdgId() <<"; et="<<
0658       // allSortedGenParticles[i].et() <<"; and eta="<<
0659       // allSortedGenParticles[i].eta() <<" beyond acceptance.";
0660       return false;
0661     }
0662   }
0663 
0664   // all tests passed
0665   return true;
0666 }
0667 ////////////////////////////////////////////////////////////////////////////////
0668 
0669 bool EmDQM::checkRecoParticlesRequirement(const edm::Event &event) {
0670   // note that this code is very similar to the one in
0671   // checkGeneratedParticlesRequirement(..) and hopefully can be merged with it
0672   // at some point in the future
0673 
0674   edm::Handle<edm::View<reco::Candidate>> referenceParticles;
0675   // get the right data according to the trigger path currently looked at
0676   if (autoConfMode_) {
0677     switch (reqNum) {
0678       case 1:
0679         if (pdgGen == 11)
0680           event.getByToken(gencutColl_fidWenu_token, referenceParticles);
0681         else
0682           event.getByToken(gencutColl_fidGammaJet_token, referenceParticles);
0683         break;
0684       case 2:
0685         if (pdgGen == 11)
0686           event.getByToken(gencutColl_fidZee_token, referenceParticles);
0687         else
0688           event.getByToken(gencutColl_fidDiGamma_token, referenceParticles);
0689         break;
0690       case 3:
0691         event.getByToken(gencutColl_fidTripleEle_token, referenceParticles);
0692         break;
0693     }
0694   } else {
0695     event.getByToken(gencutColl_manualConf_token, referenceParticles);
0696   }
0697   if (!referenceParticles.isValid()) {
0698     if (verbosity_ >= OUTPUT_WARNINGS)
0699       edm::LogWarning("EmDQM") << "referenceParticles invalid.";
0700     return false;
0701   }
0702 
0703   std::vector<const reco::Candidate *> allSortedReferenceParticles;
0704 
0705   for (edm::View<reco::Candidate>::const_iterator currentReferenceParticle = referenceParticles->begin();
0706        currentReferenceParticle != referenceParticles->end();
0707        currentReferenceParticle++) {
0708     if (currentReferenceParticle->et() <= 2.0)
0709       continue;
0710 
0711     // Note that for determining the overall efficiency,
0712     // we should only allow
0713     //
0714     // HOWEVER: for turn-on curves, we need to let
0715     //          more electrons pass
0716     if (currentReferenceParticle->et() < plotEtMin)
0717       continue;
0718 
0719     // TODO: instead of filling a new vector we could simply count here...
0720     allSortedReferenceParticles.push_back(&(*currentReferenceParticle));
0721   }
0722 
0723   // std::sort(allSortedReferenceParticles.begin(),
0724   // allSortedReferenceParticles.end(),pTComparator_);
0725 
0726   // return false if not enough particles found
0727   return allSortedReferenceParticles.size() >= gencut_;
0728 }
0729 
0730 ////////////////////////////////////////////////////////////////////////////////
0731 //                     method called to for each event                        //
0732 ////////////////////////////////////////////////////////////////////////////////
0733 void EmDQM::analyze(const edm::Event &event, const edm::EventSetup &setup) {
0734   // loop over all the trigger path parameter sets
0735   unsigned int vPos = 0;
0736   for (std::vector<edm::ParameterSet>::iterator psetIt = paramSets.begin(); psetIt != paramSets.end();
0737        ++psetIt, ++vPos) {
0738     SetVarsFromPSet(psetIt);
0739     // get the set forthe current path
0740     hltCollectionLabelsFound = hltCollectionLabelsFoundPerPath.at(vPos);
0741     hltCollectionLabelsMissed = hltCollectionLabelsMissedPerPath.at(vPos);
0742 
0743     ////////////////////////////////////////////////////////////
0744     //           Check if there's enough gen particles        //
0745     //             of interest                                //
0746     ////////////////////////////////////////////////////////////
0747     // get the right data according to the trigger path currently looked at
0748     edm::Handle<edm::View<reco::Candidate>> cutCounter;
0749     if (autoConfMode_) {
0750       switch (reqNum) {
0751         case 1:
0752           if (pdgGen == 11)
0753             event.getByToken(gencutColl_fidWenu_token, cutCounter);
0754           else
0755             event.getByToken(gencutColl_fidGammaJet_token, cutCounter);
0756           break;
0757         case 2:
0758           if (pdgGen == 11)
0759             event.getByToken(gencutColl_fidZee_token, cutCounter);
0760           else
0761             event.getByToken(gencutColl_fidDiGamma_token, cutCounter);
0762           break;
0763         case 3:
0764           event.getByToken(gencutColl_fidTripleEle_token, cutCounter);
0765           break;
0766       }
0767     } else {
0768       event.getByToken(gencutColl_manualConf_token, cutCounter);
0769     }
0770     if (cutCounter->size() < (unsigned int)gencut_) {
0771       // edm::LogWarning("EmDQM") << "Less than "<< gencut_ <<" gen particles
0772       // with pdgId=" << pdgGen;
0773       continue;
0774     }
0775 
0776     // fill L1 and HLT info
0777     // get objects possed by each filter
0778     edm::Handle<trigger::TriggerEventWithRefs> triggerObj;
0779     event.getByToken(triggerObject_token, triggerObj);
0780     if (!triggerObj.isValid()) {
0781       if (verbosity_ >= OUTPUT_WARNINGS)
0782         edm::LogWarning("EmDQM") << "parameter triggerobject (" << triggerObject_
0783                                  << ") does not corresond to a valid TriggerEventWithRefs product. "
0784                                     "Please check especially the process name (e.g. when running "
0785                                     "over reprocessed datasets)";
0786       continue;
0787     }
0788 
0789     // Were enough high energy gen particles found?
0790     if (event.isRealData()) {
0791       // running validation on data.
0792       // TODO: we should check that the entire
0793       //       run is on the same type (all data or
0794       //       all MC). Otherwise one gets
0795       //       uninterpretable results...
0796       if (!checkRecoParticlesRequirement(event))
0797         continue;
0798     } else {  // MC
0799       // if no, throw event away
0800       if (!checkGeneratedParticlesRequirement(event))
0801         continue;
0802     }
0803 
0804     // It was an event worth keeping. Continue.
0805 
0806     ////////////////////////////////////////////////////////////
0807     //  Fill the bin labeled "Total"                          //
0808     //   This will be the number of events looked at.         //
0809     ////////////////////////////////////////////////////////////
0810     if (!mcMatchedOnly_)
0811       totals.at(vPos)->Fill(numOfHLTCollectionLabels + 0.5);
0812     totalmatchs.at(vPos)->Fill(numOfHLTCollectionLabels + 0.5);
0813 
0814     ////////////////////////////////////////////////////////////
0815     //               Fill generator info                      //
0816     ////////////////////////////////////////////////////////////
0817     // the gencut_ highest Et generator objects of the preselected type are our
0818     // matches
0819 
0820     std::vector<reco::Particle> sortedGen;
0821     for (edm::View<reco::Candidate>::const_iterator genpart = cutCounter->begin(); genpart != cutCounter->end();
0822          genpart++) {
0823       reco::Particle tmpcand(genpart->charge(), genpart->p4(), genpart->vertex(), genpart->pdgId(), genpart->status());
0824       if (tmpcand.et() >= plotEtMin) {
0825         sortedGen.push_back(tmpcand);
0826       }
0827     }
0828     std::sort(sortedGen.begin(), sortedGen.end(), pTComparator_);
0829 
0830     // Now the collection of gen particles is sorted by pt.
0831     // So, remove all particles from the collection so that we
0832     // only have the top "1 thru gencut_" particles in it
0833     if (sortedGen.size() < gencut_) {
0834       continue;
0835     }
0836     sortedGen.erase(sortedGen.begin() + gencut_, sortedGen.end());
0837 
0838     for (unsigned int i = 0; i < gencut_; i++) {
0839       etgens.at(vPos)->Fill(sortedGen[i].et());  // validity has been implicitily checked by the
0840                                                  // cut on gencut_ above
0841       if (sortedGen[i].et() > minEtForEtaEffPlot_) {
0842         etagens.at(vPos)->Fill(sortedGen[i].eta());
0843         if (!noPhiPlots_)
0844           phigens.at(vPos)->Fill(sortedGen[i].phi());
0845         etaphigens.at(vPos)->Fill(sortedGen[i].eta(), sortedGen[i].phi());
0846       }
0847     }  // END of loop over Generated particles
0848     if (gencut_ >= reqNum && !mcMatchedOnly_)
0849       totals.at(vPos)->Fill(numOfHLTCollectionLabels +
0850                             1.5);  // this isn't really needed anymore keep for backward comp.
0851     if (gencut_ >= reqNum)
0852       totalmatchs.at(vPos)->Fill(numOfHLTCollectionLabels +
0853                                  1.5);  // this isn't really needed anymore keep for backward comp.
0854 
0855     bool accepted = true;  // flags that the event has been accepted by all filters before
0856     edm::Handle<edm::TriggerResults> hltResults;
0857     event.getByToken(hltResults_token, hltResults);
0858     ////////////////////////////////////////////////////////////
0859     //            Loop over filter modules                    //
0860     ////////////////////////////////////////////////////////////
0861     for (unsigned int n = 0; n < numOfHLTCollectionLabels; n++) {
0862       // check that there are not less sortedGen particles than nCandCut
0863       // requires for this filter
0864       if (sortedGen.size() < nCandCuts.at(n)) {
0865         if (verbosity_ >= OUTPUT_ERRORS)
0866           edm::LogError("EmDQM") << "There are less generated particles than the module '"
0867                                  << theHLTCollectionLabels[n].label() << "' requires.";
0868         continue;
0869       }
0870       std::vector<reco::Particle> sortedGenForFilter(sortedGen);
0871       sortedGenForFilter.erase(sortedGenForFilter.begin() + nCandCuts.at(n), sortedGenForFilter.end());
0872 
0873       // Fill only if this filter was run.
0874       if (pathIndex != 0 &&
0875           hltConfig_.moduleIndex(pathIndex, theHLTCollectionLabels[n].label()) > hltResults->index(pathIndex))
0876         break;
0877       // These numbers are from the Parameter Set, such as:
0878       //   theHLTOutputTypes = cms.uint32(100)
0879       switch (theHLTOutputTypes[n]) {
0880         case trigger::TriggerL1NoIsoEG:  // Non-isolated Level 1
0881           histoFillerL1NonIso->fillHistos(triggerObj, event, vPos, n, sortedGenForFilter, accepted);
0882           break;
0883         case trigger::TriggerL1IsoEG:  // Isolated Level 1
0884           histoFillerL1Iso->fillHistos(triggerObj, event, vPos, n, sortedGenForFilter, accepted);
0885           break;
0886         case trigger::TriggerPhoton:  // Photon
0887           histoFillerPho->fillHistos(triggerObj, event, vPos, n, sortedGenForFilter, accepted);
0888           break;
0889         case trigger::TriggerElectron:  // Electron
0890           histoFillerEle->fillHistos(triggerObj, event, vPos, n, sortedGenForFilter, accepted);
0891           break;
0892         case trigger::TriggerCluster:  // TriggerCluster
0893           histoFillerClu->fillHistos(triggerObj, event, vPos, n, sortedGenForFilter, accepted);
0894           break;
0895         default:
0896           throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]");
0897       }
0898     }  // END of loop over filter modules
0899 
0900     // earse the dummy and fill with real set
0901     hltCollectionLabelsFoundPerPath.at(vPos) = hltCollectionLabelsFound;
0902     hltCollectionLabelsMissedPerPath.at(vPos) = hltCollectionLabelsMissed;
0903   }
0904 }
0905 
0906 ////////////////////////////////////////////////////////////////////////////////
0907 // fillHistos                                                                 //
0908 //   Called by analyze method.                                                //
0909 ////////////////////////////////////////////////////////////////////////////////
0910 template <class T>
0911 void HistoFiller<T>::fillHistos(edm::Handle<trigger::TriggerEventWithRefs> &triggerObj,
0912                                 const edm::Event &iEvent,
0913                                 unsigned int vPos,
0914                                 unsigned int n,
0915                                 std::vector<reco::Particle> &sortedGen,
0916                                 bool &accepted) {
0917   std::vector<edm::Ref<T>> recoecalcands;
0918   if ((triggerObj->filterIndex(dqm->theHLTCollectionLabels[n]) >= triggerObj->size())) {  // only process if available
0919     dqm->hltCollectionLabelsMissed.insert(dqm->theHLTCollectionLabels[n].encode());
0920     accepted = false;
0921     return;
0922   }
0923 
0924   dqm->hltCollectionLabelsFound.insert(dqm->theHLTCollectionLabels[n].encode());
0925 
0926   ////////////////////////////////////////////////////////////
0927   //      Retrieve saved filter objects                     //
0928   ////////////////////////////////////////////////////////////
0929   triggerObj->getObjects(
0930       triggerObj->filterIndex(dqm->theHLTCollectionLabels[n]), dqm->theHLTOutputTypes[n], recoecalcands);
0931   // Danger: special case, L1 non-isolated
0932   // needs to be merged with L1 iso
0933   if (dqm->theHLTOutputTypes[n] == trigger::TriggerL1NoIsoEG) {
0934     std::vector<edm::Ref<T>> isocands;
0935     triggerObj->getObjects(triggerObj->filterIndex(dqm->theHLTCollectionLabels[n]), trigger::TriggerL1IsoEG, isocands);
0936     if (!isocands.empty()) {
0937       for (unsigned int i = 0; i < isocands.size(); i++)
0938         recoecalcands.push_back(isocands[i]);
0939     }
0940   }  // END of if theHLTOutputTypes == 82
0941 
0942   if (recoecalcands.empty()) {  // stop if no object passed the previous filter
0943     accepted = false;
0944     return;
0945   }
0946 
0947   // if (recoecalcands.size() >= reqNum )
0948   if (recoecalcands.size() >= dqm->nCandCuts.at(n) && !dqm->mcMatchedOnly_)
0949     dqm->totals.at(vPos)->Fill(n + 0.5);
0950 
0951   ///////////////////////////////////////////////////
0952   // check for validity                            //
0953   // prevents crash in CMSSW_3_1_0_pre6            //
0954   ///////////////////////////////////////////////////
0955   for (unsigned int j = 0; j < recoecalcands.size(); j++) {
0956     if (!(recoecalcands.at(j).isAvailable())) {
0957       if (dqm->verbosity_ >= dqm->OUTPUT_ERRORS)
0958         edm::LogError("EmDQMInvalidRefs") << "Event content inconsistent: TriggerEventWithRefs contains "
0959                                              "invalid Refs. Invalid refs for: "
0960                                           << dqm->theHLTCollectionLabels[n].label()
0961                                           << ". The collection that this module uses may has been dropped in "
0962                                              "the event.";
0963       return;
0964     }
0965   }
0966 
0967   if (!dqm->mcMatchedOnly_) {
0968     ////////////////////////////////////////////////////////////
0969     // Loop over the Generated Particles, and find the        //
0970     // closest HLT object match.                              //
0971     ////////////////////////////////////////////////////////////
0972     // for (unsigned int i=0; i < gencut_; i++) {
0973     for (unsigned int i = 0; i < dqm->nCandCuts.at(n); i++) {
0974       math::XYZVector currentGenParticleMomentum = sortedGen[i].momentum();
0975 
0976       float closestDeltaR = 0.5;
0977       int closestEcalCandIndex = -1;
0978       for (unsigned int j = 0; j < recoecalcands.size(); j++) {
0979         float deltaR = DeltaR(recoecalcands[j]->momentum(), currentGenParticleMomentum);
0980 
0981         if (deltaR < closestDeltaR) {
0982           closestDeltaR = deltaR;
0983           closestEcalCandIndex = j;
0984         }
0985       }
0986 
0987       // If an HLT object was found within some delta-R
0988       // of this gen particle, store it in a histogram
0989       if (closestEcalCandIndex >= 0) {
0990         dqm->histEtOfHltObjMatchToGens.at(vPos).at(n)->Fill(recoecalcands[closestEcalCandIndex]->et());
0991         dqm->histEtaOfHltObjMatchToGens.at(vPos).at(n)->Fill(recoecalcands[closestEcalCandIndex]->eta());
0992         if (!dqm->noPhiPlots_)
0993           dqm->histPhiOfHltObjMatchToGens.at(vPos).at(n)->Fill(recoecalcands[closestEcalCandIndex]->phi());
0994         dqm->histEtaPhiOfHltObjMatchToGens.at(vPos).at(n)->Fill(recoecalcands[closestEcalCandIndex]->eta(),
0995                                                                 recoecalcands[closestEcalCandIndex]->phi());
0996       }  // END of if closestEcalCandIndex >= 0
0997     }
0998 
0999     ////////////////////////////////////////////////////////////
1000     //  Loop over all HLT objects in this filter step, and    //
1001     //  fill histograms.                                      //
1002     ////////////////////////////////////////////////////////////
1003     //  bool foundAllMatches = false;
1004     //  unsigned int numOfHLTobjectsMatched = 0;
1005     for (unsigned int i = 0; i < recoecalcands.size(); i++) {
1006       //// See if this HLT object has a gen-level match
1007       // float closestGenParticleDr = 99.0;
1008       // for(unsigned int j =0; j < gencut_; j++) {
1009       //  math::XYZVector currentGenParticle = sortedGen[j].momentum();
1010 
1011       //  double currentDeltaR =
1012       //  DeltaR(recoecalcands[i]->momentum(),currentGenParticle); if (
1013       //  currentDeltaR < closestGenParticleDr ) {
1014       //    closestGenParticleDr = currentDeltaR;
1015       //  }
1016       //}
1017       //// If this HLT object did not have a gen particle match, go to next HLT
1018       /// object
1019       // if ( !(fabs(closestGenParticleDr < 0.3)) ) continue;
1020 
1021       // numOfHLTobjectsMatched++;
1022       // if (numOfHLTobjectsMatched >= gencut_) foundAllMatches=true;
1023 
1024       // Fill HLT object histograms
1025       dqm->ethists.at(vPos).at(n)->Fill(recoecalcands[i]->et());
1026       dqm->etahists.at(vPos).at(n)->Fill(recoecalcands[i]->eta());
1027       if (!dqm->noPhiPlots_)
1028         dqm->phihists.at(vPos).at(n)->Fill(recoecalcands[i]->phi());
1029       dqm->etaphihists.at(vPos).at(n)->Fill(recoecalcands[i]->eta(), recoecalcands[i]->phi());
1030     }
1031   }
1032 
1033   ////////////////////////////////////////////////////////////
1034   //        Fill mc matched objects into histograms         //
1035   ////////////////////////////////////////////////////////////
1036   unsigned int matchedMcParts = 0;
1037   float mindist = 0.3;
1038   if (n == 0)
1039     mindist = 0.5;  // low L1-resolution => allow wider matching
1040   for (unsigned int i = 0; i < dqm->nCandCuts.at(n); ++i) {
1041     // match generator candidate
1042     bool matchThis = false;
1043     math::XYZVector candDir = sortedGen[i].momentum();
1044     // unsigned int closest = 0;
1045     double closestDr = 1000.;
1046     for (unsigned int trigOb = 0; trigOb < recoecalcands.size(); ++trigOb) {
1047       double dr = DeltaR(recoecalcands[trigOb]->momentum(), candDir);
1048       if (dr < closestDr) {
1049         closestDr = dr;
1050         // closest = trigOb;
1051       }
1052       if (closestDr > mindist) {  // it's not really a "match" if it's that far away
1053         // closest = -1;
1054       } else {
1055         matchedMcParts++;
1056         matchThis = true;
1057       }
1058     }
1059     if (!matchThis) {
1060       accepted = false;
1061       continue;  // only plot matched candidates
1062     }
1063     // fill coordinates of mc particle matching trigger object
1064     dqm->ethistmatchs.at(vPos).at(n)->Fill(sortedGen[i].et());
1065     if (sortedGen[i].et() > dqm->minEtForEtaEffPlot_) {
1066       dqm->etahistmatchs.at(vPos).at(n)->Fill(sortedGen[i].eta());
1067       if (!dqm->noPhiPlots_)
1068         dqm->phihistmatchs.at(vPos).at(n)->Fill(sortedGen[i].phi());
1069       dqm->etaphihistmatchs.at(vPos).at(n)->Fill(sortedGen[i].eta(), sortedGen[i].phi());
1070     }
1071   }
1072   // fill total mc matched efficiency
1073 
1074   if (matchedMcParts >= dqm->nCandCuts.at(n) && accepted == true)
1075     dqm->totalmatchs.at(vPos)->Fill(n + 0.5);
1076 }
1077 
1078 void EmDQM::dqmEndRun(edm::Run const &iRun, edm::EventSetup const &iSetup) {
1079   // loop over all the trigger path parameter sets
1080   unsigned int vPos = 0;
1081   for (std::vector<edm::ParameterSet>::iterator psetIt = paramSets.begin(); psetIt != paramSets.end();
1082        ++psetIt, ++vPos) {
1083     SetVarsFromPSet(psetIt);
1084 
1085     // print information about hltCollectionLabels which were not found
1086     // (but only those which were never found)
1087 
1088     // check which ones were never found
1089     std::vector<std::string> labelsNeverFound;
1090 
1091     for (const auto &tag : hltCollectionLabelsMissedPerPath.at(vPos)) {
1092       auto atag = edm::InputTag(tag);
1093       if ((hltCollectionLabelsFoundPerPath.at(vPos)).count(atag.encode()) == 0)
1094         // never found
1095         labelsNeverFound.push_back(atag.encode());
1096 
1097     }  // loop over all tags which were missed at least once
1098 
1099     if (labelsNeverFound.empty())
1100       continue;
1101 
1102     std::sort(labelsNeverFound.begin(), labelsNeverFound.end());
1103 
1104     // there was at least one label which was never found
1105     // (note that this could also be because the corresponding
1106     // trigger path slowly fades out to zero efficiency)
1107     if (verbosity_ >= OUTPUT_WARNINGS)
1108       edm::LogWarning("EmDQM") << "There were some HLTCollectionLabels which were never found:";
1109 
1110     for (auto const &tag : labelsNeverFound) {
1111       if (verbosity_ >= OUTPUT_ALL)
1112         edm::LogPrint("EmDQM") << "  " << tag;
1113     }
1114   }
1115 }
1116 
1117 // returns count of non-overlapping occurrences of 'sub' in 'str'
1118 int EmDQM::countSubstring(const std::string &str, const std::string &sub) {
1119   if (sub.length() == 0)
1120     return 0;
1121   int count = 0;
1122   for (size_t offset = str.find(sub); offset != std::string::npos; offset = str.find(sub, offset + sub.length())) {
1123     ++count;
1124   }
1125   return count;
1126 }
1127 
1128 //----------------------------------------------------------------------
1129 // find egamma trigger paths from trigger name
1130 std::vector<std::vector<std::string>> EmDQM::findEgammaPaths() {
1131   std::vector<std::vector<std::string>> Paths(5);
1132   // Loop over all paths in the menu
1133   for (unsigned int i = 0; i < hltConfig_.size(); i++) {
1134     std::string path = hltConfig_.triggerName(i);
1135 
1136     // Find electron and photon paths and classify them
1137     if (int(path.find("HLT_")) == 0) {  // Path should start with 'HLT_'
1138 
1139       int scCount = countSubstring(path, "_SC");
1140       int eleCount = countSubstring(path, "Ele");
1141       int doubleEleCount = countSubstring(path, "DoubleEle");
1142       int doubleSCCount = countSubstring(path, "DiSC");
1143       int tripleEleCount = countSubstring(path, "TripleEle");
1144       int photonCount = countSubstring(path, "hoton");
1145       int doublePhotonCount = countSubstring(path, "DoublePhoton") + countSubstring(path, "Diphoton");
1146 
1147       int totEleCount = 2 * tripleEleCount + doubleEleCount + eleCount + scCount + 2 * doubleSCCount;
1148       int totPhotonCount = doublePhotonCount + photonCount;
1149 
1150       if (totEleCount + totPhotonCount < 1)
1151         continue;
1152       switch (totEleCount) {
1153         case 1:
1154           Paths[TYPE_SINGLE_ELE].push_back(path);
1155           // std::cout << "Electron \t" << path << std::endl;
1156           break;
1157         case 2:
1158           Paths[TYPE_DOUBLE_ELE].push_back(path);
1159           // std::cout << "DoubleElectron \t" << path << std::endl;
1160           break;
1161         case 3:
1162           Paths[TYPE_TRIPLE_ELE].push_back(path);
1163           // std::cout << "TripleElectron \t" << path << std::endl;
1164           break;
1165       }
1166 
1167       switch (totPhotonCount) {
1168         case 1:
1169           Paths[TYPE_SINGLE_PHOTON].push_back(path);
1170           // std::cout << "Photon \t\t" << path << std::endl;
1171           break;
1172         case 2:
1173           Paths[TYPE_DOUBLE_PHOTON].push_back(path);
1174           // std::cout << "DoublePhoton \t" << path << std::endl;
1175           break;
1176       }
1177     }
1178     // std::cout << i << " triggerName: " << path << " containing " <<
1179     // hltConfig_.size(i) << " modules."<< std::endl;
1180   }
1181   return Paths;
1182 }
1183 
1184 //----------------------------------------------------------------------
1185 // get the names of filters of a given path
1186 std::vector<std::string> EmDQM::getFilterModules(const std::string &path) {
1187   std::vector<std::string> filters;
1188 
1189   // std::cout << "Pathname: " << path << std::endl;
1190 
1191   // Loop over all modules in the path
1192   for (unsigned int i = 0; i < hltConfig_.size(path); i++) {
1193     std::string module = hltConfig_.moduleLabel(path, i);
1194     std::string moduleType = hltConfig_.moduleType(module);
1195     std::string moduleEDMType = hltConfig_.moduleEDMType(module);
1196 
1197     // Find filters
1198     if (moduleEDMType == "EDFilter" ||
1199         moduleType.find("Filter") != std::string::npos) {  // older samples may not have EDMType data
1200                                                            // included
1201       filters.push_back(module);
1202       // std::cout << i << "    moduleLabel: " << module << "    moduleType: "
1203       // << moduleType << "    moduleEDMType: " << moduleEDMType << std::endl;
1204     }
1205   }
1206   return filters;
1207 }
1208 
1209 //----------------------------------------------------------------------
1210 // get the primary Et cut from the trigger name
1211 double EmDQM::getPrimaryEtCut(const std::string &path) {
1212   double minEt = -1;
1213 
1214   boost::regex reg("^HLT_.*?(Ele|hoton|EG|SC)([[:digit:]]+).*");
1215 
1216   boost::smatch what;
1217   if (boost::regex_match(path, what, reg, boost::match_extra)) {
1218     minEt = std::stod(what[2]);
1219   }
1220 
1221   return minEt;
1222 }
1223 
1224 //----------------------------------------------------------------------
1225 
1226 edm::ParameterSet EmDQM::makePSetForL1SeedFilter(const std::string &moduleName) {
1227   // generates a PSet to analyze the behaviour of an L1 seed.
1228   //
1229   // moduleName is the name of the HLT module which filters
1230   // on the L1 seed.
1231   edm::ParameterSet retPSet;
1232 
1233   retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1234   // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1235   // edm::InputTag(moduleName, "", processName_));
1236   retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1237   retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections",
1238                                                    std::vector<edm::InputTag>(1, std::string("none")));
1239   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerL1NoIsoEG);
1240 
1241   // as HLTLevel1GTSeed has no parameter ncandcut we determine the value from
1242   // the name of the filter
1243 
1244   int orCount = countSubstring(moduleName, "OR") + countSubstring(moduleName, "Or");
1245   int egCount = countSubstring(moduleName, "EG");
1246   int dEgCount = countSubstring(moduleName, "DoubleEG");
1247   int tEgCount = countSubstring(moduleName, "TripleEG");
1248 
1249   int candCount = 2 * tEgCount + dEgCount + egCount;
1250   // if L1 is and OR of triggers try the assumption that all of them are similar
1251   // first and if not take the lowest number. If this is not successful, let the
1252   // path name decide
1253   if (orCount > 0 && candCount > 0) {
1254     if (egCount % (orCount + 1) == 0 && dEgCount % (orCount + 1) == 0 && tEgCount % (orCount + 1) == 0)
1255       candCount /= (orCount + 1);
1256     else if (egCount - dEgCount - tEgCount > 0)
1257       candCount = 1;  // at least one singleEG present
1258     else if (dEgCount > 0)
1259       candCount = 2;  // at least one doubleEG and no singleEG present
1260     else if (tEgCount > 0)
1261       candCount = 3;  // only tripleEG present
1262     else
1263       candCount = -1;
1264   }
1265 
1266   switch (candCount) {
1267     case 0:
1268       retPSet.addParameter<int>("ncandcut", 0);
1269       break;
1270     case 1:
1271       retPSet.addParameter<int>("ncandcut", 1);
1272       break;
1273     case 2:
1274       retPSet.addParameter<int>("ncandcut", 2);
1275       break;
1276     case 3:
1277       retPSet.addParameter<int>("ncandcut", 3);
1278       break;
1279     default:
1280       retPSet.addParameter<int>("ncandcut", -1);
1281   }
1282 
1283   return retPSet;
1284 }
1285 
1286 //----------------------------------------------------------------------
1287 
1288 edm::ParameterSet EmDQM::makePSetForL1SeedToSuperClusterMatchFilter(const std::string &moduleName) {
1289   // generates a PSet to analyze the behaviour of L1 to supercluster match
1290   // filter.
1291   //
1292   // moduleName is the name of the HLT module which requires the match
1293   // between supercluster and L1 seed.
1294   //
1295   edm::ParameterSet retPSet;
1296 
1297   retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1298   // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1299   // edm::InputTag(moduleName, "", processName_));
1300   retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1301   retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections",
1302                                                    std::vector<edm::InputTag>(1, std::string("none")));
1303   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerCluster);
1304   retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1305 
1306   return retPSet;
1307 }
1308 
1309 //----------------------------------------------------------------------
1310 
1311 edm::ParameterSet EmDQM::makePSetForEtFilter(const std::string &moduleName) {
1312   edm::ParameterSet retPSet;
1313 
1314   retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1315   // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1316   // edm::InputTag(moduleName, "", processName_));
1317   retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1318   retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections",
1319                                                    std::vector<edm::InputTag>(1, std::string("none")));
1320   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerCluster);
1321   retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1322 
1323   return retPSet;
1324 }
1325 
1326 //----------------------------------------------------------------------
1327 
1328 edm::ParameterSet EmDQM::makePSetForOneOEMinusOneOPFilter(const std::string &moduleName) {
1329   edm::ParameterSet retPSet;
1330 
1331   retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1332   // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1333   // edm::InputTag(moduleName, "", processName_));
1334   retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1335   retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections",
1336                                                    std::vector<edm::InputTag>(1, std::string("none")));
1337   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerElectron);
1338   retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1339 
1340   return retPSet;
1341 }
1342 
1343 //----------------------------------------------------------------------
1344 
1345 edm::ParameterSet EmDQM::makePSetForPixelMatchFilter(const std::string &moduleName) {
1346   edm::ParameterSet retPSet;
1347 
1348   retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1349   // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1350   // edm::InputTag(moduleName, "", processName_));
1351   retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1352   retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections",
1353                                                    std::vector<edm::InputTag>(1, std::string("none")));
1354   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerCluster);
1355   retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1356 
1357   return retPSet;
1358 }
1359 
1360 //----------------------------------------------------------------------
1361 
1362 edm::ParameterSet EmDQM::makePSetForEgammaDoubleEtDeltaPhiFilter(const std::string &moduleName) {
1363   edm::ParameterSet retPSet;
1364 
1365   retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1366   // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1367   // edm::InputTag(moduleName, "", processName_));
1368   retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1369   retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections",
1370                                                    std::vector<edm::InputTag>(1, std::string("none")));
1371   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerCluster);
1372   retPSet.addParameter<int>("ncandcut", 2);
1373 
1374   return retPSet;
1375 }
1376 
1377 //----------------------------------------------------------------------
1378 
1379 edm::ParameterSet EmDQM::makePSetForEgammaGenericFilter(const std::string &moduleName) {
1380   edm::ParameterSet retPSet;
1381 
1382   // example usages of HLTEgammaGenericFilter are:
1383   //   R9 shape filter
1384   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolR9ShapeFilter
1385   //   cluster shape filter
1386   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolClusterShapeFilter
1387   //   Ecal isolation filter
1388   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TIghterEleIdIsolEcalIsolFilter H/E
1389   //   filter hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolHEFilter
1390   //   HCAL isolation filter
1391   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolHcalIsolFilter
1392 
1393   // infer the type of filter by the type of the producer which
1394   // generates the collection used to cut on this
1395   edm::InputTag varTag = hltConfig_.modulePSet(moduleName).getParameter<edm::InputTag>("varTag");
1396   // edm::InputTag nonIsoTag =
1397   // hltConfig_.modulePSet(moduleName).getParameter<edm::InputTag>("");
1398   // std::cout
1399   // << "varTag.label " << varTag.label() << " nonIsoTag.label " <<
1400   // nonIsoTag.label() << std::endl;
1401 
1402   std::string inputType = hltConfig_.moduleType(varTag.label());
1403   // std::cout << "inputType " << inputType << " moduleName " << moduleName <<
1404   // std::endl;
1405 
1406   //--------------------
1407   // sanity check: non-isolated path should be produced by the
1408   // same type of module
1409 
1410   // first check that the non-iso tag is non-empty
1411   // if (nonIsoTag.label().empty()) {
1412   //  edm::LogError("EmDQM") << "nonIsoTag of HLTEgammaGenericFilter '" <<
1413   //  moduleName <<  "' is empty."; return retPSet;
1414   //}
1415   // if (inputType != hltConfig_.moduleType(nonIsoTag.label())) {
1416   //  edm::LogError("EmDQM") << "C++ Type of varTag '" << inputType << "' and
1417   //  nonIsoTag '" << hltConfig_.moduleType(nonIsoTag.label()) << "' are not the
1418   //  same for HLTEgammaGenericFilter '" << moduleName <<  "'."; return retPSet;
1419   //}
1420   //--------------------
1421 
1422   // parameter saveTag determines the output type
1423   if (hltConfig_.saveTags(moduleName))
1424     retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerPhoton);
1425   else
1426     retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerCluster);
1427 
1428   std::vector<edm::InputTag> isoCollections;
1429   isoCollections.push_back(varTag);
1430   // if (!nonIsoTag.label().empty())
1431   //   isoCollections.push_back(nonIsoTag);
1432 
1433   //--------------------
1434   // the following cases seem to have identical PSets ?
1435   //--------------------
1436 
1437   if (inputType == "EgammaHLTR9IDProducer" ||                      // R9 ID
1438       inputType == "EgammaHLTClusterShapeProducer" ||              // cluster shape
1439       inputType == "EgammaHLTGsfTrackVarProducer" ||               // GSF track deta and dphi filter
1440       inputType == "EgammaHLTBcHcalIsolationProducersRegional" ||  // HCAL isolation, or H for H/E (based on CaloTowers)
1441       inputType == "EgammaHLTHcalVarProducerFromRecHit" ||         // HCAL isolation, or H for H/E (based on RecHits)
1442       inputType == "EgammaHLTEcalPFClusterIsolationProducer" ||    // ECAL PF isolation filter
1443       inputType == "EgammaHLTHcalPFClusterIsolationProducer" ||    // HCAL PF isolation filter
1444       inputType == "EgammaHLTElectronTrackIsolationProducers"      // Track isolation filter
1445   ) {
1446     retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1447     // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1448     // edm::InputTag(moduleName, "", processName_));
1449     retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1450     retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections", isoCollections);
1451     retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1452 
1453     return retPSet;
1454   }
1455 
1456   if (verbosity_ >= OUTPUT_ERRORS)
1457     edm::LogError("EmDQM") << "Can't determine what the HLTEgammaGenericFilter '" << moduleName
1458                            << "' should do: uses a collection produced by a module of C++ type '" << inputType << "'.";
1459   return edm::ParameterSet();
1460 }
1461 
1462 //----------------------------------------------------------------------
1463 
1464 edm::ParameterSet EmDQM::makePSetForEgammaGenericQuadraticFilter(const std::string &moduleName) {
1465   edm::ParameterSet retPSet;
1466 
1467   // example usages of HLTEgammaGenericFilter are:
1468   //   R9 shape filter
1469   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolR9ShapeFilter
1470   //   cluster shape filter
1471   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolClusterShapeFilter
1472   //   Ecal isolation filter
1473   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TIghterEleIdIsolEcalIsolFilter H/E
1474   //   filter hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolHEFilter
1475   //   HCAL isolation filter
1476   //   hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolHcalIsolFilter
1477 
1478   // infer the type of filter by the type of the producer which
1479   // generates the collection used to cut on this
1480   edm::InputTag varTag = hltConfig_.modulePSet(moduleName).getParameter<edm::InputTag>("varTag");
1481   // edm::InputTag nonIsoTag =
1482   // hltConfig_.modulePSet(moduleName).getParameter<edm::InputTag>("nonIsoTag");
1483   // std::cout << "varTag.label " << varTag.label() << " nonIsoTag.label " <<
1484   // nonIsoTag.label() << std::endl;
1485 
1486   std::string inputType = hltConfig_.moduleType(varTag.label());
1487   // std::cout << "inputType " << inputType << " moduleName " << moduleName <<
1488   // std::endl;
1489 
1490   //--------------------
1491   // sanity check: non-isolated path should be produced by the
1492   // same type of module
1493 
1494   // first check that the non-iso tag is non-empty
1495   // if (nonIsoTag.label().empty()) {
1496   //  edm::LogError("EmDQM") << "nonIsoTag of HLTEgammaGenericFilter '" <<
1497   //  moduleName <<  "' is empty."; return retPSet;
1498   //}
1499   // if (inputType != hltConfig_.moduleType(nonIsoTag.label())) {
1500   //  edm::LogError("EmDQM") << "C++ Type of varTag '" << inputType << "' and
1501   //  nonIsoTag '" << hltConfig_.moduleType(nonIsoTag.label()) << "' are not the
1502   //  same for HLTEgammaGenericFilter '" << moduleName <<  "'."; return retPSet;
1503   //}
1504   //--------------------
1505 
1506   // parameter saveTag determines the output type
1507   if (hltConfig_.saveTags(moduleName))
1508     retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerPhoton);
1509   else
1510     retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerCluster);
1511 
1512   std::vector<edm::InputTag> isoCollections;
1513   isoCollections.push_back(varTag);
1514   // if (!nonIsoTag.label().empty())
1515   //   isoCollections.push_back(nonIsoTag);
1516 
1517   //--------------------
1518   // the following cases seem to have identical PSets ?
1519   //--------------------
1520 
1521   if (inputType == "EgammaHLTR9IDProducer" ||                      // R9 ID
1522       inputType == "EgammaHLTClusterShapeProducer" ||              // cluster shape
1523       inputType == "EgammaHLTBcHcalIsolationProducersRegional" ||  // HCAL isolation, or H for H/E (based on CaloTowers)
1524       inputType == "EgammaHLTHcalVarProducerFromRecHit" ||         // HCAL isolation, or H for H/E (based on RecHits)
1525       inputType == "EgammaHLTEcalPFClusterIsolationProducer" ||    // ECAL PF isolation filter
1526       inputType == "EgammaHLTHcalPFClusterIsolationProducer" ||    // HCAL PF isolation filter
1527       inputType == "EgammaHLTPhotonTrackIsolationProducersRegional"  // Photon track isolation
1528   ) {
1529     retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1530     // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1531     // edm::InputTag(moduleName, "", processName_));
1532     retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1533     retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections", isoCollections);
1534     retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1535 
1536     return retPSet;
1537   }
1538 
1539   if (verbosity_ >= OUTPUT_ERRORS)
1540     edm::LogError("EmDQM") << "Can't determine what the HLTEgammaGenericQuadraticFilter '" << moduleName
1541                            << "' should do: uses a collection produced by a module of C++ type '" << inputType << "'.";
1542   return edm::ParameterSet();
1543 }
1544 
1545 //----------------------------------------------------------------------
1546 
1547 edm::ParameterSet EmDQM::makePSetForElectronGenericFilter(const std::string &moduleName) {
1548   edm::ParameterSet retPSet;
1549 
1550   // example usages of HLTElectronGenericFilter are:
1551   //
1552   // deta filter
1553   // hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolDetaFilter dphi
1554   // filter hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolDphiFilter
1555   // track isolation
1556   // hltL1NonIsoHLTNonIsoSingleElectronEt17TighterEleIdIsolTrackIsolFilter
1557 
1558   // infer the type of filter by the type of the producer which
1559   // generates the collection used to cut on this
1560   edm::InputTag varTag = hltConfig_.modulePSet(moduleName).getParameter<edm::InputTag>("varTag");
1561   // edm::InputTag nonIsoTag =
1562   // hltConfig_.modulePSet(moduleName).getParameter<edm::InputTag>("nonIsoTag");
1563   // std::cout << "varTag.label " << varTag.label() << " nonIsoTag.label " <<
1564   // nonIsoTag.label() << std::endl;
1565 
1566   std::string inputType = hltConfig_.moduleType(varTag.label());
1567   // std::cout << "inputType iso " << inputType << " inputType noniso " <<
1568   // hltConfig_.moduleType(nonIsoTag.label()) << " moduleName " << moduleName <<
1569   // std::endl;
1570 
1571   //--------------------
1572   // sanity check: non-isolated path should be produced by the
1573   // same type of module
1574   // if (nonIsoTag.label().empty()) {
1575   //  edm::LogError("EmDQM") << "nonIsoTag of HLTElectronGenericFilter '" <<
1576   //  moduleName <<  "' is empty."; return retPSet;
1577   //}
1578   // if (inputType != hltConfig_.moduleType(nonIsoTag.label())) {
1579   //  edm::LogError("EmDQM") << "C++ Type of varTag '" << inputType << "' and
1580   //  nonIsoTag '" << hltConfig_.moduleType(nonIsoTag.label()) << "' are not the
1581   //  same for HLTElectronGenericFilter '" << moduleName <<  "'."; return
1582   //  retPSet;
1583   //}
1584   //--------------------
1585 
1586   // the type of object to look for seems to be the
1587   // same for all uses of HLTEgammaGenericFilter
1588   retPSet.addParameter<int>("theHLTOutputTypes", trigger::TriggerElectron);
1589 
1590   std::vector<edm::InputTag> isoCollections;
1591   isoCollections.push_back(varTag);
1592   // if (!nonIsoTag.label().empty())
1593   //   isoCollections.push_back(nonIsoTag);
1594 
1595   //--------------------
1596   // the following cases seem to have identical PSets ?
1597   //--------------------
1598 
1599   // note that whether deta or dphi is used is determined from
1600   // the product instance (not the module label)
1601   if (inputType == "EgammaHLTElectronTrackIsolationProducers"  // track isolation
1602   ) {
1603     retPSet.addParameter<std::vector<double>>("PlotBounds", std::vector<double>(2, 0.0));
1604     // retPSet.addParameter<edm::InputTag>("HLTCollectionLabels",
1605     // edm::InputTag(moduleName, "", processName_));
1606     retPSet.addParameter<edm::InputTag>("HLTCollectionLabels", edm::InputTag(moduleName, "", triggerObject_.process()));
1607     retPSet.addParameter<std::vector<edm::InputTag>>("IsoCollections", isoCollections);
1608     retPSet.addParameter<int>("ncandcut", hltConfig_.modulePSet(moduleName).getParameter<int>("ncandcut"));
1609 
1610     return retPSet;
1611   }
1612 
1613   if (verbosity_ >= OUTPUT_ERRORS)
1614     edm::LogError("EmDQM") << "Can't determine what the HLTElectronGenericFilter '" << moduleName
1615                            << "' should do: uses a collection produced by a module of C++ type '" << inputType << "'.";
1616   return edm::ParameterSet();
1617 }
1618 
1619 //----------------------------------------------------------------------
1620 
1621 void EmDQM::SetVarsFromPSet(std::vector<edm::ParameterSet>::iterator psetIt) {
1622   dirname_ = "HLT/HLTEgammaValidation/" + psetIt->getParameter<std::string>("@module_label");
1623 
1624   pathIndex = psetIt->getUntrackedParameter<unsigned int>("pathIndex", 0);
1625   // parameters for generator study
1626   reqNum = psetIt->getParameter<unsigned int>("reqNum");
1627   pdgGen = psetIt->getParameter<int>("pdgGen");
1628   // plotting parameters (untracked because they don't affect the physics)
1629   plotEtMin = psetIt->getUntrackedParameter<double>("genEtMin", 0.);
1630   plotPtMin = psetIt->getUntrackedParameter<double>("PtMin", 0.);
1631   plotPtMax = psetIt->getUntrackedParameter<double>("PtMax", 1000.);
1632 
1633   // preselction cuts
1634   gencutCollection_ = psetIt->getParameter<edm::InputTag>("cutcollection");
1635   gencut_ = psetIt->getParameter<int>("cutnum");
1636 
1637   ////////////////////////////////////////////////////////////
1638   //         Read in the Vector of Parameter Sets.          //
1639   //           Information for each filter-step             //
1640   ////////////////////////////////////////////////////////////
1641   std::vector<edm::ParameterSet> filters = psetIt->getParameter<std::vector<edm::ParameterSet>>("filters");
1642 
1643   // empty vectors of parameters from previous paths
1644   theHLTCollectionLabels.clear();
1645   theHLTOutputTypes.clear();
1646   theHLTCollectionHumanNames.clear();
1647   plotBounds.clear();
1648   isoNames.clear();
1649   plotiso.clear();
1650   nCandCuts.clear();
1651 
1652   int i = 0;
1653   for (std::vector<edm::ParameterSet>::iterator filterconf = filters.begin(); filterconf != filters.end();
1654        filterconf++) {
1655     theHLTCollectionLabels.push_back(filterconf->getParameter<edm::InputTag>("HLTCollectionLabels"));
1656     theHLTOutputTypes.push_back(filterconf->getParameter<int>("theHLTOutputTypes"));
1657     // Grab the human-readable name, if it is not specified, use the Collection
1658     // Label
1659     theHLTCollectionHumanNames.push_back(
1660         filterconf->getUntrackedParameter<std::string>("HLTCollectionHumanName", theHLTCollectionLabels[i].label()));
1661 
1662     std::vector<double> bounds = filterconf->getParameter<std::vector<double>>("PlotBounds");
1663     // If the size of plot "bounds" vector != 2, abort
1664     assert(bounds.size() == 2);
1665     plotBounds.push_back(std::pair<double, double>(bounds[0], bounds[1]));
1666     isoNames.push_back(filterconf->getParameter<std::vector<edm::InputTag>>("IsoCollections"));
1667 
1668     // for (unsigned int i=0; i<isoNames.back().size(); i++) {
1669     //  switch(theHLTOutputTypes.back())  {
1670     //  case trigger::TriggerL1NoIsoEG:
1671     //    histoFillerL1NonIso->isoNameTokens_.push_back(consumes<edm::AssociationMap<edm::OneToValue<l1extra::L1EmParticleCollection
1672     //    , float>>>(isoNames.back()[i])); break;
1673     //  case trigger::TriggerL1IsoEG: // Isolated Level 1
1674     //    histoFillerL1Iso->isoNameTokens_.push_back(consumes<edm::AssociationMap<edm::OneToValue<l1extra::L1EmParticleCollection
1675     //    , float>>>(isoNames.back()[i])); break;
1676     //  case trigger::TriggerPhoton: // Photon
1677     //    histoFillerPho->isoNameTokens_.push_back(consumes<edm::AssociationMap<edm::OneToValue<reco::RecoEcalCandidateCollection
1678     //    , float>>>(isoNames.back()[i])); break;
1679     //  case trigger::TriggerElectron: // Electron
1680     //    histoFillerEle->isoNameTokens_.push_back(consumes<edm::AssociationMap<edm::OneToValue<reco::ElectronCollection
1681     //    , float>>>(isoNames.back()[i])); break;
1682     //  case trigger::TriggerCluster: // TriggerCluster
1683     //    histoFillerClu->isoNameTokens_.push_back(consumes<edm::AssociationMap<edm::OneToValue<reco::RecoEcalCandidateCollection
1684     //    , float>>>(isoNames.back()[i])); break;
1685     //  default:
1686     //    throw(cms::Exception("Release Validation Error") << "HLT output type
1687     //    not implemented: theHLTOutputTypes[n]" );
1688     //  }
1689     //}
1690 
1691     // If the size of the isoNames vector is not greater than zero, abort
1692     assert(!isoNames.back().empty());
1693     if (isoNames.back().at(0).label() == "none") {
1694       plotiso.push_back(false);
1695     } else {
1696       if (!noIsolationPlots_)
1697         plotiso.push_back(true);
1698       else
1699         plotiso.push_back(false);
1700     }
1701     nCandCuts.push_back(filterconf->getParameter<int>("ncandcut"));
1702     i++;
1703   }  // END of loop over parameter sets
1704 
1705   // Record number of HLTCollectionLabels
1706   numOfHLTCollectionLabels = theHLTCollectionLabels.size();
1707 }
1708 
1709 //----------------------------------------------------------------------
1710 
1711 DEFINE_FWK_MODULE(EmDQM);