CMS 3D CMS Logo

Classes | Public Member Functions | Private Member Functions | Private Attributes

HLTMuonValidator Class Reference

#include <HLTMuonValidator.h>

Inheritance diagram for HLTMuonValidator:
edm::EDAnalyzer

List of all members.

Classes

struct  matchesByDescendingPt
struct  MatchStruct

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void beginRun (const edm::Run &, const edm::EventSetup &)
 HLTMuonValidator (const edm::ParameterSet &)

Private Member Functions

void analyzePath (const edm::Event &, const std::string &, const std::string &, const std::vector< MatchStruct >, edm::Handle< trigger::TriggerEventWithRefs >)
void bookHist (std::string, std::string, std::string, std::string)
void findMatches (std::vector< MatchStruct > &, std::vector< l1extra::L1MuonParticleRef >, std::vector< std::vector< const reco::RecoChargedCandidate * > >)
void initializeHists ()

Private Attributes

double cutMaxEta_
unsigned int cutMotherId_
std::vector< double > cutsDr_
std::map< std::string, double > cutsMinPt_
DQMStoredbe_
std::map< std::string,
MonitorElement * > 
elements_
std::map< std::string,
std::vector< std::string > > 
filterLabels_
std::string genMuonCut_
StringCutObjectSelector
< reco::GenParticle > * 
genMuonSelector_
std::string genParticleLabel_
HLTConfigProvider hltConfig_
std::set< std::string > hltPaths_
std::vector< std::string > hltPathsToCheck_
std::string hltProcessName_
std::string l1CandLabel_
L1MuonMatcherAlgo l1Matcher_
std::string l2CandLabel_
std::string l3CandLabel_
std::vector< double > parametersEta_
std::vector< double > parametersPhi_
std::vector< double > parametersTurnOn_
std::string recMuonCut_
std::string recMuonLabel_
StringCutObjectSelector
< reco::Muon > * 
recMuonSelector_
std::map< std::string,
std::vector< std::string > > 
stepLabels_

Detailed Description

Generate histograms for muon trigger efficiencies Documentation available on the CMS TWiki: https://twiki.cern.ch/twiki/bin/view/CMS/MuonHLTOfflinePerformance

Date:
2011/04/29 21:41:56
Revision:
1.12
Author:
J. Klukas, M. Vander Donckt, J. Alcaraz

Definition at line 60 of file HLTMuonValidator.h.


Constructor & Destructor Documentation

HLTMuonValidator::HLTMuonValidator ( const edm::ParameterSet pset)

Definition at line 37 of file HLTMuonValidator.cc.

References cutsDr_, dbe_, genMuonCut_, genMuonSelector_, genParticleLabel_, edm::ParameterSet::getParameter(), hltPathsToCheck_, hltProcessName_, l1CandLabel_, l2CandLabel_, l3CandLabel_, cmsCodeRules::cppFunctionSkipper::operator, parametersEta_, parametersPhi_, parametersTurnOn_, recMuonCut_, recMuonLabel_, recMuonSelector_, and DQMStore::setVerbose().

                                                            :
  l1Matcher_(pset)
{

  hltProcessName_  = pset.getParameter< string         >("hltProcessName");
  hltPathsToCheck_ = pset.getParameter< vector<string> >("hltPathsToCheck");

  genParticleLabel_ = pset.getParameter<string>("genParticleLabel");
      recMuonLabel_ = pset.getParameter<string>(    "recMuonLabel");
       l1CandLabel_ = pset.getParameter<string>(     "l1CandLabel");
       l2CandLabel_ = pset.getParameter<string>(     "l2CandLabel");
       l3CandLabel_ = pset.getParameter<string>(     "l3CandLabel");

  cutsDr_      = pset.getParameter< vector<double> >("cutsDr"     );

  parametersEta_    = pset.getParameter< vector<double> >("parametersEta");
  parametersPhi_    = pset.getParameter< vector<double> >("parametersPhi");
  parametersTurnOn_ = pset.getParameter< vector<double> >("parametersTurnOn");

  genMuonCut_ = pset.getParameter<string>("genMuonCut");
  recMuonCut_ = pset.getParameter<string>("recMuonCut");

  genMuonSelector_ = 0;
  recMuonSelector_ = 0;

  dbe_ = Service<DQMStore>().operator->();
  dbe_->setVerbose(0);

}

Member Function Documentation

void HLTMuonValidator::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 212 of file HLTMuonValidator.cc.

References analyzePath(), edm::HandleBase::failedToGet(), genMuonCut_, genMuonSelector_, genParticleLabel_, genParticleCandidates2GenParticles_cfi::genParticles, edm::Event::getByLabel(), hltPaths_, i, edm::HandleBase::isValid(), LogTrace, recMuonCut_, recMuonLabel_, recMuonSelector_, python::multivaluedict::sort(), and LaserTracksInput_cfi::source.

{

  static int eventNumber = 0;
  eventNumber++;
  LogTrace("HLTMuonVal") << "In HLTMuonValidator::analyze,  " 
                         << "Event: " << eventNumber;

  Handle<          TriggerEventWithRefs> rawTriggerEvent;
  Handle<                MuonCollection> recMuons;
  Handle<         GenParticleCollection> genParticles;

  iEvent.getByLabel("hltTriggerSummaryRAW", rawTriggerEvent);
  if (rawTriggerEvent.failedToGet())
    {LogError("HLTMuonVal") << "No trigger summary found"; return;}
  iEvent.getByLabel(    recMuonLabel_, recMuons     );
  iEvent.getByLabel(genParticleLabel_, genParticles );

  vector<string> sources;
  if (genParticles.isValid()) sources.push_back("gen");
  if (    recMuons.isValid()) sources.push_back("rec");

  for (size_t sourceNo = 0; sourceNo < sources.size(); sourceNo++) {

    string source = sources[sourceNo];
    
    // If this is the first event, initialize selectors
    if (!genMuonSelector_) genMuonSelector_ =
      new StringCutObjectSelector<reco::GenParticle>(genMuonCut_);
    if (!recMuonSelector_) recMuonSelector_ =
      new StringCutObjectSelector<reco::Muon       >(recMuonCut_);

    // Make each good gen/rec muon into the base cand for a MatchStruct
    vector<MatchStruct> matches;
    if (source == "gen" && genParticles.isValid())
      for (size_t i = 0; i < genParticles->size(); i++)
        if ((*genMuonSelector_)(genParticles->at(i)))
          matches.push_back(MatchStruct(& genParticles->at(i)));
    if (source == "rec" && recMuons.isValid())
      for (size_t i = 0; i < recMuons->size(); i++)
        if ((*recMuonSelector_)(recMuons->at(i)))
          matches.push_back(MatchStruct(& recMuons->at(i)));
    
    // Sort the MatchStructs by pT for later filling of turn-on curve
    sort(matches.begin(), matches.end(), matchesByDescendingPt());

    set<string>::iterator iPath;
    for (iPath = hltPaths_.begin(); iPath != hltPaths_.end(); iPath++)
      analyzePath(iEvent, * iPath, source, matches, rawTriggerEvent);

  } // End loop over sources

}
void HLTMuonValidator::analyzePath ( const edm::Event ,
const std::string &  ,
const std::string &  ,
const std::vector< MatchStruct ,
edm::Handle< trigger::TriggerEventWithRefs  
) [private]

Referenced by analyze().

void HLTMuonValidator::beginJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 70 of file HLTMuonValidator.cc.

{
}
void HLTMuonValidator::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 77 of file HLTMuonValidator.cc.

References filterLabels_, hltConfig_, hltPaths_, hltPathsToCheck_, hltProcessName_, i, L1MuonMatcherAlgo::init(), HLTConfigProvider::init(), initializeHists(), j, l1Matcher_, listBenchmarks::pattern, inputsource_file_cfi::runNumber, and HLTConfigProvider::triggerNames().

{

  static int runNumber = 0;
  runNumber++;

  l1Matcher_.init(iSetup);

  bool changedConfig;
  if (!hltConfig_.init(iRun, iSetup, hltProcessName_, changedConfig)) {
    LogError("HLTMuonVal") << "Initialization of HLTConfigProvider failed!!"; 
    return;
  }

  hltPaths_.clear();
  filterLabels_.clear();
  vector<string> validTriggerNames = hltConfig_.triggerNames();
  for (size_t i = 0; i < hltPathsToCheck_.size(); i++) {
    TPRegexp pattern(hltPathsToCheck_[i]);
    for (size_t j = 0; j < validTriggerNames.size(); j++)
      if (TString(validTriggerNames[j]).Contains(pattern))
        hltPaths_.insert(validTriggerNames[j]);
    // Add fake HLT path where we bypass all filters
    if (TString("NoFilters").Contains(pattern))
      hltPaths_.insert("NoFilters");
  }

  initializeHists();

}
void HLTMuonValidator::bookHist ( std::string  ,
std::string  ,
std::string  ,
std::string   
) [private]

Definition at line 476 of file HLTMuonValidator.cc.

References DQMStore::book1D(), dbe_, prof2calltree::edges, elements_, h, i, label, max(), min, mergeVDriftHistosByStation::name, parametersEta_, parametersPhi_, parametersTurnOn_, LaserTracksInput_cfi::source, and indexGen::title.

Referenced by initializeHists().

{

  string sourceUpper = source; 
  sourceUpper[0] = toupper(sourceUpper[0]);
  string name  = source + "Pass" + type + "_" + label;
  string rootName = path + "_" + name;
  TH1F * h;

  if (type.find("MaxPt") != string::npos) {
    string desc = (type == "MaxPt1") ? "Leading" : "Next-to-Leading";
    string title = "pT of " + desc + " " + sourceUpper + " Muon "+
                   "matched to " + label;
    const size_t nBins = parametersTurnOn_.size() - 1;
    float * edges = new float[nBins + 1];
    for (size_t i = 0; i < nBins + 1; i++) edges[i] = parametersTurnOn_[i];
    h = new TH1F(rootName.c_str(), title.c_str(), nBins, edges);
  }

  else {
    string symbol = (type == "Eta") ? "#eta" : "#phi";
    string title  = symbol + " of " + sourceUpper + " Muons " +
                    "matched to " + label;
    vector<double> params = (type == "Eta") ? parametersEta_ : parametersPhi_; 
    int    nBins = (int)params[0];
    double min   = params[1];
    double max   = params[2];
    h = new TH1F(rootName.c_str(), title.c_str(), nBins, min, max);
  }

  h->Sumw2();
  elements_[rootName] = dbe_->book1D(name, h);
  delete h;

}
void HLTMuonValidator::findMatches ( std::vector< MatchStruct > &  ,
std::vector< l1extra::L1MuonParticleRef ,
std::vector< std::vector< const reco::RecoChargedCandidate * > >   
) [private]

Definition at line 391 of file HLTMuonValidator.cc.

References begin, cutsDr_, deltaR(), reco::Candidate::eta(), i, edm::eventsetup::heterocontainer::insert(), j, kNull, and reco::Candidate::phi().

{

  set<size_t>::iterator it;

  set<size_t> indicesL1;
  for (size_t i = 0; i < candsL1.size(); i++) 
    indicesL1.insert(i);

  vector< set<size_t> > indicesHlt(candsHlt.size());
  for (size_t i = 0; i < candsHlt.size(); i++)
    for (size_t j = 0; j < candsHlt[i].size(); j++)
      indicesHlt[i].insert(j);

  for (size_t i = 0; i < matches.size(); i++) {

    const Candidate * cand = matches[i].candBase;

    double bestDeltaR = cutsDr_[0];
    size_t bestMatch = kNull;
    for (it = indicesL1.begin(); it != indicesL1.end(); it++) {
     if (candsL1[*it].isAvailable()) {
      double dR = deltaR(cand->eta(), cand->phi(),
                         candsL1[*it]->eta(), candsL1[*it]->phi());
      if (dR < bestDeltaR) {
        bestMatch = *it;
        bestDeltaR = dR;
      }
      // TrajectoryStateOnSurface propagated;
      // float dR = 999., dPhi = 999.;
      // bool isValid = l1Matcher_.match(* cand, * candsL1[*it], 
      //                                 dR, dPhi, propagated);
      // if (isValid && dR < bestDeltaR) {
      //   bestMatch = *it;
      //   bestDeltaR = dR;
      // }
     } else {
       LogWarning("HLTMuonValidator")
         << "Ref candsL1[*it]: product not available "
         << *it;
     }
    }
    if (bestMatch != kNull)
      matches[i].candL1 = & * candsL1[bestMatch];
    indicesL1.erase(bestMatch);

    matches[i].candHlt.assign(candsHlt.size(), 0);
    for (size_t j = 0; j < candsHlt.size(); j++) {
      size_t level = (candsHlt.size() == 4) ? (j < 2) ? 2 : 3 :
                     (candsHlt.size() == 2) ? (j < 1) ? 2 : 3 :
                     2;
      bestDeltaR = cutsDr_[level - 2];
      bestMatch = kNull;
      for (it = indicesHlt[j].begin(); it != indicesHlt[j].end(); it++) {
        double dR = deltaR(cand->eta(), cand->phi(),
                           candsHlt[j][*it]->eta(), candsHlt[j][*it]->phi());
        if (dR < bestDeltaR) {
          bestMatch = *it;
          bestDeltaR = dR;
        }
      }
      if (bestMatch != kNull)
        matches[i].candHlt[j] = candsHlt[j][bestMatch];
      indicesHlt[j].erase(bestMatch);
    }

//     cout << "    Muon: " << cand->eta() << ", ";
//     if (matches[i].candL1) cout << matches[i].candL1->eta() << ", ";
//     else cout << "none, ";
//     for (size_t j = 0; j < candsHlt.size(); j++) 
//       if (matches[i].candHlt[j]) cout << matches[i].candHlt[j]->eta() << ", ";
//       else cout << "none, ";
//     cout << endl;

  }

}
void HLTMuonValidator::initializeHists ( ) [private]

Definition at line 111 of file HLTMuonValidator.cc.

References cmsMakeMELists::baseDir, DQMStore::bookFloat(), bookHist(), cutsMinPt_, dbe_, elements_, filterLabels_, spr::find(), DQMStore::get(), hltConfig_, hltPaths_, i, j, HLTConfigProvider::moduleLabels(), path(), DQMStore::setCurrentFolder(), LaserTracksInput_cfi::source, stepLabels_, and dtDQMClient_cfg::threshold.

Referenced by beginRun().

{

  vector<string> sources(2);
  sources[0] = "gen";
  sources[1] = "rec";

  set<string>::iterator iPath;

  for (iPath = hltPaths_.begin(); iPath != hltPaths_.end(); iPath++) {
 
    string path = * iPath;

    if (path == "NoFilters") {
      filterLabels_[path].push_back("hltL1sL1SingleMuOpenL1SingleMu0");
      // The L2 and L3 collections will be built manually later on
      filterLabels_[path].push_back("");
      filterLabels_[path].push_back("");
    }

    vector<string> moduleLabels;
    if (path != "NoFilters") moduleLabels = hltConfig_.moduleLabels(path);

    for (size_t i = 0; i < moduleLabels.size(); i++)
      if (moduleLabels[i].find("Filtered") != string::npos)
        filterLabels_[path].push_back(moduleLabels[i]);

    double cutMaxEta = 2.4;


    // Choose a pT cut for gen/rec muons based on the pT cut in the path
    unsigned int threshold = 3;
    TPRegexp ptRegexp("Mu([0-9]*)");
    TObjArray * regexArray = ptRegexp.MatchS(path);
    if (regexArray->GetEntriesFast() == 2) {
      threshold = atoi(((TObjString *)regexArray->At(1))->GetString());
    }
    delete regexArray;
    // We select a whole number min pT cut slightly above the path's final 
    // pt threshold, then subtract a bit to let through particle gun muons with
    // exact integer pT:
    double cutMinPt = ceil(threshold * 1.1) - 0.01;
    if (cutMinPt < 0. || path == "NoFilters") cutMinPt = 0.;
    cutsMinPt_[path] = cutMinPt;

    string baseDir = "HLT/Muon/Distributions/";
    string pathSansSuffix = path;
    if (path.rfind("_v") < path.length())
      pathSansSuffix = path.substr(0, path.rfind("_v"));
    dbe_->setCurrentFolder(baseDir + pathSansSuffix);

    if (dbe_->get(baseDir + path + "/CutMinPt") == 0) {

      elements_[path + "_" + "CutMinPt" ] = dbe_->bookFloat("CutMinPt" );
      elements_[path + "_" + "CutMaxEta"] = dbe_->bookFloat("CutMaxEta");
      elements_[path + "_" + "CutMinPt" ]->Fill(cutMinPt);
      elements_[path + "_" + "CutMaxEta"]->Fill(cutMaxEta);

      // Standardize the names that will be applied to each step
      const int nFilters = filterLabels_[path].size();
      stepLabels_[path].push_back("All");
      stepLabels_[path].push_back("L1");
      if (nFilters == 2) {
        stepLabels_[path].push_back("L2");
      }
      if (nFilters == 3) {
        stepLabels_[path].push_back("L2");
        stepLabels_[path].push_back("L3");
      }
      if (nFilters == 5) {
        stepLabels_[path].push_back("L2");
        stepLabels_[path].push_back("L2Iso");
        stepLabels_[path].push_back("L3");
        stepLabels_[path].push_back("L3Iso");
      }

      //     string l1Name = path + "_L1Quality";
      //     elements_[l1Name.c_str()] = 
      //       dbe_->book1D("L1Quality", "Quality of L1 Muons", 8, 0, 8);
      //     for (size_t i = 0; i < 8; i++)
      //       elements_[l1Name.c_str()]->setBinLabel(i + 1, Form("%i", i));

      for (size_t i = 0; i < sources.size(); i++) {
        string source = sources[i];
        for (size_t j = 0; j < stepLabels_[path].size(); j++) {
          bookHist(path, stepLabels_[path][j], source, "Eta");
          bookHist(path, stepLabels_[path][j], source, "Phi");
          bookHist(path, stepLabels_[path][j], source, "MaxPt1");
          bookHist(path, stepLabels_[path][j], source, "MaxPt2");
        }
      }

    }

  }

}

Member Data Documentation

double HLTMuonValidator::cutMaxEta_ [private]

Definition at line 125 of file HLTMuonValidator.h.

unsigned int HLTMuonValidator::cutMotherId_ [private]

Definition at line 126 of file HLTMuonValidator.h.

std::vector<double> HLTMuonValidator::cutsDr_ [private]

Definition at line 127 of file HLTMuonValidator.h.

Referenced by findMatches(), and HLTMuonValidator().

std::map<std::string, double> HLTMuonValidator::cutsMinPt_ [private]

Definition at line 123 of file HLTMuonValidator.h.

Referenced by initializeHists().

Definition at line 138 of file HLTMuonValidator.h.

Referenced by bookHist(), HLTMuonValidator(), and initializeHists().

std::map<std::string, MonitorElement *> HLTMuonValidator::elements_ [private]

Definition at line 139 of file HLTMuonValidator.h.

Referenced by bookHist(), and initializeHists().

std::map<std::string, std::vector<std::string> > HLTMuonValidator::filterLabels_ [private]

Definition at line 111 of file HLTMuonValidator.h.

Referenced by beginRun(), and initializeHists().

std::string HLTMuonValidator::genMuonCut_ [private]

Definition at line 128 of file HLTMuonValidator.h.

Referenced by analyze(), and HLTMuonValidator().

Definition at line 131 of file HLTMuonValidator.h.

Referenced by analyze(), and HLTMuonValidator().

std::string HLTMuonValidator::genParticleLabel_ [private]

Definition at line 113 of file HLTMuonValidator.h.

Referenced by analyze(), and HLTMuonValidator().

Definition at line 134 of file HLTMuonValidator.h.

Referenced by beginRun(), and initializeHists().

std::set<std::string> HLTMuonValidator::hltPaths_ [private]

Definition at line 110 of file HLTMuonValidator.h.

Referenced by analyze(), beginRun(), and initializeHists().

std::vector<std::string> HLTMuonValidator::hltPathsToCheck_ [private]

Definition at line 109 of file HLTMuonValidator.h.

Referenced by beginRun(), and HLTMuonValidator().

std::string HLTMuonValidator::hltProcessName_ [private]

Definition at line 107 of file HLTMuonValidator.h.

Referenced by beginRun(), and HLTMuonValidator().

std::string HLTMuonValidator::l1CandLabel_ [private]

Definition at line 115 of file HLTMuonValidator.h.

Referenced by HLTMuonValidator().

Definition at line 136 of file HLTMuonValidator.h.

Referenced by beginRun().

std::string HLTMuonValidator::l2CandLabel_ [private]

Definition at line 116 of file HLTMuonValidator.h.

Referenced by HLTMuonValidator().

std::string HLTMuonValidator::l3CandLabel_ [private]

Definition at line 117 of file HLTMuonValidator.h.

Referenced by HLTMuonValidator().

std::vector<double> HLTMuonValidator::parametersEta_ [private]

Definition at line 119 of file HLTMuonValidator.h.

Referenced by bookHist(), and HLTMuonValidator().

std::vector<double> HLTMuonValidator::parametersPhi_ [private]

Definition at line 120 of file HLTMuonValidator.h.

Referenced by bookHist(), and HLTMuonValidator().

std::vector<double> HLTMuonValidator::parametersTurnOn_ [private]

Definition at line 121 of file HLTMuonValidator.h.

Referenced by bookHist(), and HLTMuonValidator().

std::string HLTMuonValidator::recMuonCut_ [private]

Definition at line 129 of file HLTMuonValidator.h.

Referenced by analyze(), and HLTMuonValidator().

std::string HLTMuonValidator::recMuonLabel_ [private]

Definition at line 114 of file HLTMuonValidator.h.

Referenced by analyze(), and HLTMuonValidator().

Definition at line 132 of file HLTMuonValidator.h.

Referenced by analyze(), and HLTMuonValidator().

std::map<std::string, std::vector<std::string> > HLTMuonValidator::stepLabels_ [private]

Definition at line 140 of file HLTMuonValidator.h.

Referenced by initializeHists().