#include <HLTMuonMatchAndPlot.h>
Public Member Functions | |
void | analyze (const edm::Event &, const edm::EventSetup &) |
void | beginRun (const edm::Run &, const edm::EventSetup &) |
void | endRun (const edm::Run &, const edm::EventSetup &) |
void | fillEdges (size_t &nBins, float *&edges, std::vector< double > binning) |
template<class T > | |
void | fillMapFromPSet (map< string, T > &m, ParameterSet pset, string target) |
template<class T > | |
void | fillMapFromPSet (std::map< std::string, T > &, edm::ParameterSet, std::string) |
HLTMuonMatchAndPlot (const edm::ParameterSet &, std::string, std::vector< std::string >) | |
Constructor. | |
template<class T1 , class T2 > | |
vector< size_t > | matchByDeltaR (const vector< T1 > &collection1, const vector< T2 > &collection2, const double maxDeltaR) |
template<class T1 , class T2 > | |
std::vector< size_t > | matchByDeltaR (const std::vector< T1 > &, const std::vector< T2 > &, const double maxDeltaR=NOMATCH) |
Private Member Functions | |
void | book1D (std::string, std::string, std::string) |
void | book2D (std::string, std::string, std::string, std::string) |
reco::MuonCollection | selectedMuons (const reco::MuonCollection &, const reco::BeamSpot &, bool, const StringCutObjectSelector< reco::Muon > &, double, double) |
trigger::TriggerObjectCollection | selectedTriggerObjects (const trigger::TriggerObjectCollection &, const trigger::TriggerEvent &, const edm::ParameterSet &) |
Private Attributes | |
std::map< std::string, std::vector< double > > | binParams_ |
unsigned int | cutMinPt_ |
DQMStore * | dbe_ |
std::string | destination_ |
bool | hasProbeRecoCuts |
bool | hasTargetRecoCuts |
std::map< std::string, MonitorElement * > | hists_ |
std::string | hltPath_ |
std::string | hltProcessName_ |
std::map< std::string, edm::InputTag > | inputTags_ |
std::vector< std::string > | moduleLabels_ |
std::map< std::string, double > | plotCuts_ |
double | probeD0Cut_ |
StringCutObjectSelector < reco::Muon > | probeMuonSelector_ |
edm::ParameterSet | probeParams_ |
double | probeZ0Cut_ |
std::vector< std::string > | requiredTriggers_ |
double | targetD0Cut_ |
StringCutObjectSelector < reco::Muon > | targetMuonSelector_ |
edm::ParameterSet | targetParams_ |
double | targetZ0Cut_ |
std::string | triggerLevel_ |
Match reconstructed muons to HLT objects and plot efficiencies.
Note that this is not a true EDAnalyzer; rather, the intent is that one EDAnalyzer would call a separate instantiation of HLTMuonMatchAndPlot for each HLT path under consideration.
Documentation available on the CMS TWiki: https://twiki.cern.ch/twiki/bin/view/CMS/MuonHLTOfflinePerformance
Definition at line 66 of file HLTMuonMatchAndPlot.h.
HLTMuonMatchAndPlot::HLTMuonMatchAndPlot | ( | const edm::ParameterSet & | , |
std::string | , | ||
std::vector< std::string > | |||
) |
Constructor.
Definition at line 44 of file HLTMuonMatchAndPlot.cc.
References binParams_, cutMinPt_, dbe_, fillMapFromPSet(), hltPath_, hltProcessName_, inputTags_, edm::InputTag::instance(), edm::InputTag::label(), moduleLabels_, cppFunctionSkipper::operator, plotCuts_, and triggerLevel_.
: hltProcessName_(pset.getParameter<string>("hltProcessName")), destination_(pset.getUntrackedParameter<string>("destination")), requiredTriggers_(pset.getUntrackedParameter<vstring>("requiredTriggers")), targetParams_(pset.getParameterSet("targetParams")), probeParams_(pset.getParameterSet("probeParams")), hltPath_(hltPath), moduleLabels_(moduleLabels), hasTargetRecoCuts(targetParams_.exists("recoCuts")), hasProbeRecoCuts(probeParams_.exists("recoCuts")), targetMuonSelector_(targetParams_.getUntrackedParameter<string>("recoCuts", "")), targetZ0Cut_(targetParams_.getUntrackedParameter<double>("z0Cut",0.)), targetD0Cut_(targetParams_.getUntrackedParameter<double>("d0Cut",0.)), probeMuonSelector_(probeParams_.getUntrackedParameter<string>("recoCuts", "")), probeZ0Cut_(probeParams_.getUntrackedParameter<double>("z0Cut",0.)), probeD0Cut_(probeParams_.getUntrackedParameter<double>("d0Cut",0.)) { // Create std::map<string, T> from ParameterSets. fillMapFromPSet(inputTags_, pset, "inputTags"); fillMapFromPSet(binParams_, pset, "binParams"); fillMapFromPSet(plotCuts_, pset, "plotCuts"); // Set HLT process name for TriggerResults and TriggerSummary. InputTag & resTag = inputTags_["triggerResults"]; resTag = InputTag(resTag.label(), resTag.instance(), hltProcessName_); InputTag & sumTag = inputTags_["triggerSummary"]; sumTag = InputTag(sumTag.label(), sumTag.instance(), hltProcessName_); // Prepare the DQMStore object. dbe_ = edm::Service<DQMStore>().operator->(); dbe_->setVerbose(0); // Get the trigger level. triggerLevel_ = "L3"; TPRegexp levelRegexp("L[1-3]"); size_t nModules = moduleLabels_.size(); TObjArray * levelArray = levelRegexp.MatchS(moduleLabels_[nModules - 1]); if (levelArray->GetEntriesFast() > 0) { triggerLevel_ = ((TObjString *)levelArray->At(0))->GetString(); } delete levelArray; // Get the pT cut by parsing the name of the HLT path. cutMinPt_ = 3; TPRegexp ptRegexp("Mu([0-9]*)"); TObjArray * objArray = ptRegexp.MatchS(hltPath_); if (objArray->GetEntriesFast() >= 2) { TObjString * ptCutString = (TObjString *)objArray->At(1); cutMinPt_ = atoi(ptCutString->GetString()); cutMinPt_ = ceil(cutMinPt_ * plotCuts_["minPtFactor"]); } delete objArray; }
void HLTMuonMatchAndPlot::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) |
Definition at line 165 of file HLTMuonMatchAndPlot.cc.
References allMuons_cfi::allMuons, SiPixelRawToDigiRegional_cfi::beamSpot, reco::LeafCandidate::charge(), cutMinPt_, deltaR(), reco::TrackBase::dxy(), reco::TrackBase::dz(), EFFICIENCY_SUFFIXES, eta(), reco::LeafCandidate::eta(), trigger::TriggerObject::eta(), edm::Event::getByLabel(), hasProbeRecoCuts, hasTargetRecoCuts, hists_, HLT_8E33v2_cff::hltMuons, i, reco::Muon::innerTrack(), inputTags_, reco::Muon::isStandAloneMuon(), reco::Muon::isTrackerMuon(), edm::HandleBase::isValid(), j, gen::k, matchByDeltaR(), metsig::muon, reco::Muon::outerTrack(), reco::LeafCandidate::p4(), reco::LeafCandidate::phi(), phi, trigger::TriggerObject::phi(), plotCuts_, probeD0Cut_, probeMuonSelector_, probeZ0Cut_, reco::LeafCandidate::pt(), trigger::TriggerObject::pt(), requiredTriggers_, selectedMuons(), selectedTriggerObjects(), findQualityFiles::size, createPayload::suffix, targetD0Cut_, targetMuonSelector_, targetParams_, targetZ0Cut_, triggerLevel_, and patRefSel_triggerSelection_cff::triggerResults.
{ // Get objects from the event. Handle<MuonCollection> allMuons; iEvent.getByLabel(inputTags_["recoMuon"], allMuons); Handle<BeamSpot> beamSpot; iEvent.getByLabel(inputTags_["beamSpot"], beamSpot); Handle<TriggerEvent> triggerSummary; iEvent.getByLabel(inputTags_["triggerSummary"], triggerSummary); Handle<TriggerResults> triggerResults; edm::Handle<VertexCollection> vertices; iEvent.getByLabel("offlinePrimaryVertices", vertices); //edm::Handle<GenParticleCollection> gen; //iEvent.getByLabel("genParticles", gen); //GenParticleCollection::const_iterator g_part; //std::vector<const GenParticle *> gen_leptsp; //std::vector<const GenParticle *> gen_momsp; if(!triggerSummary.isValid()) { edm::LogError("HLTMuonMatchAndPlot")<<"Missing triggerSummary with label " << inputTags_["triggerSummary"] <<std::endl; return; } iEvent.getByLabel(inputTags_["triggerResults"], triggerResults); if(!triggerResults.isValid()) { edm::LogError("HLTMuonMatchAndPlot")<<"Missing triggerResults with label " << inputTags_["triggerResults"] <<std::endl; return; } /* if(gen != 0) { for(g_part = gen->begin(); g_part != gen->end(); g_part++){ if( abs(g_part->pdgId()) == 13) { if(!( g_part->status() ==1 || (g_part->status() ==2 && abs(g_part->pdgId())==5))) continue; bool GenMomExists (true); bool GenGrMomExists(true); if( g_part->pt() < 10.0 ) continue; if( fabs(g_part->eta()) > 2.4 ) continue; int gen_id= g_part->pdgId(); const GenParticle* gen_lept = &(*g_part); // get mother of gen_lept const GenParticle* gen_mom = static_cast<const GenParticle*> (gen_lept->mother()); if(gen_mom==NULL) GenMomExists=false; int m_id=-999; if(GenMomExists) m_id = gen_mom->pdgId(); if(m_id != gen_id || !GenMomExists); else{ int id= m_id; while(id == gen_id && GenMomExists){ gen_mom = static_cast<const GenParticle*> (gen_mom->mother()); if(gen_mom==NULL){ GenMomExists=false; } if(GenMomExists) id=gen_mom->pdgId(); } } if(GenMomExists) m_id = gen_mom->pdgId(); gen_leptsp.push_back(gen_lept); gen_momsp.push_back(gen_mom); } } } std::vector<GenParticle> gen_lepts; for(size_t i = 0; i < gen_leptsp.size(); i++) { gen_lepts.push_back(*gen_leptsp[i]); } */ // Throw out this event if it doesn't pass the required triggers. for (size_t i = 0; i < requiredTriggers_.size(); i++) { unsigned int triggerIndex = triggerResults->find(requiredTriggers_[i]); if (triggerIndex < triggerResults->size() || !triggerResults->accept(triggerIndex)) return; } // Select objects based on the configuration. MuonCollection targetMuons = selectedMuons(* allMuons, * beamSpot, hasTargetRecoCuts, targetMuonSelector_, targetD0Cut_, targetZ0Cut_); MuonCollection probeMuons = selectedMuons(* allMuons, * beamSpot, hasProbeRecoCuts, probeMuonSelector_, probeD0Cut_, probeZ0Cut_); TriggerObjectCollection allTriggerObjects = triggerSummary->getObjects(); TriggerObjectCollection hltMuons = selectedTriggerObjects(allTriggerObjects, * triggerSummary, targetParams_); // Fill plots for HLT muons. for (size_t i = 0; i < hltMuons.size(); i++) { hists_["hltPt"]->Fill(hltMuons[i].pt()); hists_["hltEta"]->Fill(hltMuons[i].eta()); hists_["hltPhi"]->Fill(hltMuons[i].phi()); } // Find the best trigger object matches for the targetMuons. vector<size_t> matches = matchByDeltaR(targetMuons, hltMuons, plotCuts_[triggerLevel_ + "DeltaR"]); // Fill plots for matched muons. bool pairalreadyconsidered = false; for (size_t i = 0; i < targetMuons.size(); i++) { Muon & muon = targetMuons[i]; // Fill plots which are not efficiencies. if (matches[i] < targetMuons.size()) { TriggerObject & hltMuon = hltMuons[matches[i]]; double ptRes = (muon.pt() - hltMuon.pt()) / muon.pt(); double etaRes = muon.eta() - hltMuon.eta(); double phiRes = muon.phi() - hltMuon.phi(); hists_["resolutionEta"]->Fill(etaRes); hists_["resolutionPhi"]->Fill(phiRes); hists_["resolutionPt"]->Fill(ptRes); hists_["deltaR"]->Fill(deltaR(muon, hltMuon)); } // Fill numerators and denominators for efficiency plots. for (size_t j = 0; j < 2; j++) { string suffix = EFFICIENCY_SUFFIXES[j]; // If no match was found, then the numerator plots don't get filled. if (suffix == "numer" && matches[i] >= targetMuons.size()) continue; if (muon.pt() > cutMinPt_) { hists_["efficiencyEta_" + suffix]->Fill(muon.eta()); hists_["efficiencyPhiVsEta_" + suffix]->Fill(muon.eta(), muon.phi()); } if (fabs(muon.eta()) < plotCuts_["maxEta"]) { hists_["efficiencyTurnOn_" + suffix]->Fill(muon.pt()); } if (muon.pt() > cutMinPt_ && fabs(muon.eta()) < plotCuts_["maxEta"]) { const Track * track = 0; if (muon.isTrackerMuon()) track = & * muon.innerTrack(); else if (muon.isStandAloneMuon()) track = & * muon.outerTrack(); if (track) { double d0 = track->dxy(beamSpot->position()); double z0 = track->dz(beamSpot->position()); hists_["efficiencyVertex_" + suffix]->Fill(vertices->size()); hists_["efficiencyPhi_" + suffix]->Fill(muon.phi()); hists_["efficiencyD0_" + suffix]->Fill(d0); hists_["efficiencyZ0_" + suffix]->Fill(z0); hists_["efficiencyCharge_" + suffix]->Fill(muon.charge()); } } } // Fill plots for tag and probe // Muon cannot be a tag because doesn't match an hlt muon if(matches[i] >= targetMuons.size()) continue; for (size_t k = 0; k < targetMuons.size(); k++) { if(k == i) continue; if(muon.pt() < 20.0) continue; Muon & theProbe = targetMuons[k]; if (muon.charge() != theProbe.charge() && !pairalreadyconsidered) { double mass = (muon.p4() + theProbe.p4()).M(); if(mass > 60 && mass < 120) { hists_["massVsEtaZ_denom"]->Fill(theProbe.eta()); hists_["massVsPtZ_denom"]->Fill(theProbe.pt()); hists_["massVsVertexZ_denom"]->Fill(vertices->size()); if(matches[k] < targetMuons.size()) { hists_["massVsEtaZ_numer"]->Fill(theProbe.eta()); hists_["massVsPtZ_numer"]->Fill(theProbe.pt()); hists_["massVsVertexZ_numer"]->Fill(vertices->size()); } pairalreadyconsidered = true; } if(mass > 1 && mass < 4) { hists_["massVsEtaJpsi_denom"]->Fill(theProbe.eta()); hists_["massVsPtJpsi_denom"]->Fill(theProbe.pt()); hists_["massVsVertexJpsi_denom"]->Fill(vertices->size()); if(matches[k] < targetMuons.size()) { hists_["massVsEtaJpsi_numer"]->Fill(theProbe.eta()); hists_["massVsPtJpsi_numer"]->Fill(theProbe.pt()); hists_["massVsVertexJpsi_numer"]->Fill(vertices->size()); } pairalreadyconsidered = true; } } } // End loop over denominator and numerator. } // End loop over targetMuons. // Plot fake rates (efficiency for HLT objects to not get matched to RECO). vector<size_t> hltMatches = matchByDeltaR(hltMuons, targetMuons, plotCuts_[triggerLevel_ + "DeltaR"]); for (size_t i = 0; i < hltMuons.size(); i++) { TriggerObject & hltMuon = hltMuons[i]; bool isFake = hltMatches[i] > hltMuons.size(); for (size_t j = 0; j < 2; j++) { string suffix = EFFICIENCY_SUFFIXES[j]; // If match is found, then numerator plots should not get filled if (suffix == "numer" && ! isFake) continue; hists_["fakerateVertex_" + suffix]->Fill(vertices->size()); hists_["fakerateEta_" + suffix]->Fill(hltMuon.eta()); hists_["fakeratePhi_" + suffix]->Fill(hltMuon.phi()); hists_["fakerateTurnOn_" + suffix]->Fill(hltMuon.pt()); } // End loop over numerator and denominator. } // End loop over hltMuons. } // End analyze() method.
void HLTMuonMatchAndPlot::beginRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | iSetup | ||
) |
Definition at line 101 of file HLTMuonMatchAndPlot.cc.
References book1D(), book2D(), dbe_, destination_, EFFICIENCY_SUFFIXES, hltPath_, i, DQMStore::setCurrentFolder(), and createPayload::suffix.
{ TPRegexp suffixPtCut("Mu[0-9]+$"); string baseDir = destination_; if (baseDir[baseDir.size() - 1] != '/') baseDir += '/'; string pathSansSuffix = hltPath_; if (hltPath_.rfind("_v") < hltPath_.length()) pathSansSuffix = hltPath_.substr(0, hltPath_.rfind("_v")); dbe_->setCurrentFolder(baseDir + pathSansSuffix); // Form is book1D(name, binningType, title) where 'binningType' is used // to fetch the bin settings from binParams_. book1D("deltaR", "deltaR", ";#Deltar(reco, HLT);"); book1D("hltPt", "pt", ";p_{T} of HLT object"); book1D("hltEta", "eta", ";#eta of HLT object"); book1D("hltPhi", "phi", ";#phi of HLT object"); book1D("resolutionEta", "resolutionEta", ";#eta^{reco}-#eta^{HLT};"); book1D("resolutionPhi", "resolutionPhi", ";#phi^{reco}-#phi^{HLT};"); book1D("resolutionPt", "resolutionRel", ";(p_{T}^{reco}-p_{T}^{HLT})/|p_{T}^{reco}|;"); for (size_t i = 0; i < 2; i++) { string suffix = EFFICIENCY_SUFFIXES[i]; book1D("efficiencyEta_" + suffix, "eta", ";#eta;"); book1D("efficiencyPhi_" + suffix, "phi", ";#phi;"); book1D("efficiencyTurnOn_" + suffix, "pt", ";p_{T};"); book1D("efficiencyD0_" + suffix, "d0", ";d0;"); book1D("efficiencyZ0_" + suffix, "z0", ";z0;"); book1D("efficiencyCharge_" + suffix, "charge", ";charge;"); book1D("efficiencyVertex_" + suffix, "NVertex", ";NVertex;"); book2D("efficiencyPhiVsEta_" + suffix, "etaCoarse", "phiCoarse", ";#eta;#phi"); book1D("fakerateEta_" + suffix, "eta", ";#eta;"); book1D("fakerateVertex_" + suffix, "NVertex", ";NVertex;"); book1D("fakeratePhi_" + suffix, "phi", ";#phi;"); book1D("fakerateTurnOn_" + suffix, "pt", ";p_{T};"); book1D("massVsEtaZ_" + suffix, "etaCoarse", ";#eta"); book1D("massVsEtaJpsi_" + suffix, "etaCoarse", ";#eta"); book1D("massVsPtZ_" + suffix, "ptCoarse", ";p_{T}"); book1D("massVsPtJpsi_" + suffix, "ptCoarse", ";p_{T}"); book1D("massVsVertexZ_" + suffix, "NVertex", ";NVertex"); book1D("massVsVertexJpsi_" + suffix, "NVertex", ";NVertex"); } }
void HLTMuonMatchAndPlot::book1D | ( | std::string | , |
std::string | , | ||
std::string | |||
) | [private] |
Referenced by beginRun().
void HLTMuonMatchAndPlot::book2D | ( | std::string | , |
std::string | , | ||
std::string | , | ||
std::string | |||
) | [private] |
Referenced by beginRun().
void HLTMuonMatchAndPlot::endRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | iSetup | ||
) |
Definition at line 158 of file HLTMuonMatchAndPlot.cc.
{ }
void HLTMuonMatchAndPlot::fillEdges | ( | size_t & | nBins, |
float *& | edges, | ||
std::vector< double > | binning | ||
) |
Definition at line 385 of file HLTMuonMatchAndPlot.cc.
Referenced by HcalSubdetDigiMonitor::book1D(), and HcalSubdetDigiMonitor::book2D().
{ if (binning.size() < 3) { LogWarning("HLTMuonVal") << "Invalid binning parameters!"; return; } // Fixed-width binning. if (binning.size() == 3) { nBins = binning[0]; edges = new float[nBins + 1]; const double min = binning[1]; const double binwidth = (binning[2] - binning[1]) / nBins; for (size_t i = 0; i <= nBins; i++) edges[i] = min + (binwidth * i); } // Variable-width binning. else { nBins = binning.size() - 1; edges = new float[nBins + 1]; for (size_t i = 0; i <= nBins; i++) edges[i] = binning[i]; } }
void HLTMuonMatchAndPlot::fillMapFromPSet | ( | std::map< std::string, T > & | , |
edm::ParameterSet | , | ||
std::string | |||
) |
Referenced by HLTMuonMatchAndPlot().
void HLTMuonMatchAndPlot::fillMapFromPSet | ( | map< string, T > & | m, |
ParameterSet | pset, | ||
string | target | ||
) |
Definition at line 418 of file HLTMuonMatchAndPlot.cc.
References edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNames(), edm::ParameterSet::getParameterSet(), edm::ParameterSet::getUntrackedParameter(), edm::ParameterSet::getUntrackedParameterSet(), and cscdqm::h::names.
{ // Get the ParameterSet with name 'target' from 'pset' ParameterSet targetPset; if (pset.existsAs<ParameterSet>(target, true)) // target is tracked targetPset = pset.getParameterSet(target); else if (pset.existsAs<ParameterSet>(target, false)) // target is untracked targetPset = pset.getUntrackedParameterSet(target); // Get the parameter names from targetPset vector<string> names = targetPset.getParameterNames(); vector<string>::const_iterator iter; for (iter = names.begin(); iter != names.end(); ++iter) { if (targetPset.existsAs<T>(* iter, true)) // target is tracked m[* iter] = targetPset.getParameter<T>(* iter); else if (targetPset.existsAs<T>(* iter, false)) // target is untracked m[* iter] = targetPset.getUntrackedParameter<T>(* iter); } }
vector<size_t> HLTMuonMatchAndPlot::matchByDeltaR | ( | const vector< T1 > & | collection1, |
const vector< T2 > & | collection2, | ||
const double | maxDeltaR | ||
) |
Definition at line 446 of file HLTMuonMatchAndPlot.cc.
References deltaR(), i, j, gen::k, NOMATCH, and query::result.
{ const size_t n1 = collection1.size(); const size_t n2 = collection2.size(); vector<size_t> result(n1, -1); vector<vector<double> > deltaRMatrix(n1, vector<double>(n2, NOMATCH)); for (size_t i = 0; i < n1; i++) for (size_t j = 0; j < n2; j++) { deltaRMatrix[i][j] = deltaR(collection1[i], collection2[j]); } // Run through the matrix n1 times to make sure we've found all matches. for (size_t k = 0; k < n1; k++) { size_t i_min = -1; size_t j_min = -1; double minDeltaR = maxDeltaR; // find the smallest deltaR for (size_t i = 0; i < n1; i++) for (size_t j = 0; j < n2; j++) if (deltaRMatrix[i][j] < minDeltaR) { i_min = i; j_min = j; minDeltaR = deltaRMatrix[i][j]; } // If a match has been made, save it and make those candidates unavailable. if (minDeltaR < maxDeltaR) { result[i_min] = j_min; deltaRMatrix[i_min] = vector<double>(n2, NOMATCH); for (size_t i = 0; i < n1; i++) deltaRMatrix[i][j_min] = NOMATCH; } } return result; }
std::vector<size_t> HLTMuonMatchAndPlot::matchByDeltaR | ( | const std::vector< T1 > & | , |
const std::vector< T2 > & | , | ||
const double | maxDeltaR = NOMATCH |
||
) |
Referenced by analyze().
MuonCollection HLTMuonMatchAndPlot::selectedMuons | ( | const reco::MuonCollection & | allMuons, |
const reco::BeamSpot & | beamSpot, | ||
bool | hasRecoCuts, | ||
const StringCutObjectSelector< reco::Muon > & | selector, | ||
double | d0Cut, | ||
double | z0Cut | ||
) | [private] |
Definition at line 490 of file HLTMuonMatchAndPlot.cc.
References reco::TrackBase::dxy(), reco::TrackBase::dz(), and reco::BeamSpot::position().
Referenced by analyze().
{ // If there is no selector (recoCuts does not exists), return an empty collection. if (!hasRecoCuts) return MuonCollection(); MuonCollection reducedMuons(allMuons); MuonCollection::iterator iter = reducedMuons.begin(); while (iter != reducedMuons.end()) { const Track * track = 0; if (iter->isTrackerMuon()) track = & * iter->innerTrack(); else if (iter->isStandAloneMuon()) track = & * iter->outerTrack(); if (track && selector(* iter) && fabs(track->dxy(beamSpot.position())) < d0Cut && fabs(track->dz(beamSpot.position())) < z0Cut) ++iter; else reducedMuons.erase(iter); } return reducedMuons; }
TriggerObjectCollection HLTMuonMatchAndPlot::selectedTriggerObjects | ( | const trigger::TriggerObjectCollection & | triggerObjects, |
const trigger::TriggerEvent & | triggerSummary, | ||
const edm::ParameterSet & | pset | ||
) | [private] |
Definition at line 520 of file HLTMuonMatchAndPlot.cc.
References edm::ParameterSet::exists(), trigger::TriggerEvent::filterIndex(), trigger::TriggerEvent::filterKeys(), edm::ParameterSet::getUntrackedParameter(), hltProcessName_, j, relativeConstraints::keys, moduleLabels_, and trigger::TriggerEvent::sizeFilters().
Referenced by analyze().
{ // If pset is empty, return an empty collection. if (!pset.exists("hltCuts")) return TriggerObjectCollection(); StringCutObjectSelector<TriggerObject> selector (pset.getUntrackedParameter<string>("hltCuts")); InputTag filterTag(moduleLabels_[moduleLabels_.size() - 1], "", hltProcessName_); size_t filterIndex = triggerSummary.filterIndex(filterTag); TriggerObjectCollection selectedObjects; if (filterIndex < triggerSummary.sizeFilters()) { const Keys &keys = triggerSummary.filterKeys(filterIndex); for (size_t j = 0; j < keys.size(); j++ ){ TriggerObject foundObject = triggerObjects[keys[j]]; if (selector(foundObject)) selectedObjects.push_back(foundObject); } } return selectedObjects; }
std::map<std::string, std::vector<double> > HLTMuonMatchAndPlot::binParams_ [private] |
Definition at line 108 of file HLTMuonMatchAndPlot.h.
Referenced by HcalSubdetDigiMonitor::book1D(), HcalSubdetDigiMonitor::book2D(), and HLTMuonMatchAndPlot().
unsigned int HLTMuonMatchAndPlot::cutMinPt_ [private] |
Definition at line 115 of file HLTMuonMatchAndPlot.h.
Referenced by analyze(), and HLTMuonMatchAndPlot().
DQMStore* HLTMuonMatchAndPlot::dbe_ [private] |
Definition at line 118 of file HLTMuonMatchAndPlot.h.
Referenced by beginRun(), HcalSubdetDigiMonitor::book1D(), HcalSubdetDigiMonitor::book2D(), and HLTMuonMatchAndPlot().
std::string HLTMuonMatchAndPlot::destination_ [private] |
Definition at line 105 of file HLTMuonMatchAndPlot.h.
Referenced by beginRun().
bool HLTMuonMatchAndPlot::hasProbeRecoCuts [private] |
Definition at line 123 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
bool HLTMuonMatchAndPlot::hasTargetRecoCuts [private] |
Definition at line 122 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
std::map<std::string, MonitorElement *> HLTMuonMatchAndPlot::hists_ [private] |
Definition at line 119 of file HLTMuonMatchAndPlot.h.
Referenced by analyze(), HcalSubdetDigiMonitor::book1D(), and HcalSubdetDigiMonitor::book2D().
std::string HLTMuonMatchAndPlot::hltPath_ [private] |
Definition at line 116 of file HLTMuonMatchAndPlot.h.
Referenced by beginRun(), and HLTMuonMatchAndPlot().
std::string HLTMuonMatchAndPlot::hltProcessName_ [private] |
Definition at line 104 of file HLTMuonMatchAndPlot.h.
Referenced by HLTMuonMatchAndPlot(), and selectedTriggerObjects().
std::map<std::string, edm::InputTag> HLTMuonMatchAndPlot::inputTags_ [private] |
Definition at line 107 of file HLTMuonMatchAndPlot.h.
Referenced by analyze(), and HLTMuonMatchAndPlot().
std::vector<std::string> HLTMuonMatchAndPlot::moduleLabels_ [private] |
Definition at line 117 of file HLTMuonMatchAndPlot.h.
Referenced by HLTMuonMatchAndPlot(), and selectedTriggerObjects().
std::map<std::string, double> HLTMuonMatchAndPlot::plotCuts_ [private] |
Definition at line 109 of file HLTMuonMatchAndPlot.h.
Referenced by analyze(), and HLTMuonMatchAndPlot().
double HLTMuonMatchAndPlot::probeD0Cut_ [private] |
Definition at line 130 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
Definition at line 128 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
Definition at line 111 of file HLTMuonMatchAndPlot.h.
double HLTMuonMatchAndPlot::probeZ0Cut_ [private] |
Definition at line 129 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
std::vector<std::string> HLTMuonMatchAndPlot::requiredTriggers_ [private] |
Definition at line 106 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
double HLTMuonMatchAndPlot::targetD0Cut_ [private] |
Definition at line 127 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
Definition at line 125 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
Definition at line 110 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
double HLTMuonMatchAndPlot::targetZ0Cut_ [private] |
Definition at line 126 of file HLTMuonMatchAndPlot.h.
Referenced by analyze().
std::string HLTMuonMatchAndPlot::triggerLevel_ [private] |
Definition at line 114 of file HLTMuonMatchAndPlot.h.
Referenced by analyze(), and HLTMuonMatchAndPlot().