CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HLTMuonMatchAndPlot Class Reference

#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 (std::map< std::string, T > &, edm::ParameterSet, std::string)
 
template<class T >
void fillMapFromPSet (map< string, T > &m, ParameterSet pset, string target)
 
 HLTMuonMatchAndPlot (const edm::ParameterSet &, std::string, std::vector< std::string >)
 Constructor. More...
 
template<class T1 , class T2 >
std::vector< size_t > matchByDeltaR (const std::vector< T1 > &, const std::vector< T2 > &, const double maxDeltaR=NOMATCH)
 
template<class T1 , class T2 >
vector< size_t > matchByDeltaR (const vector< T1 > &collection1, const vector< T2 > &collection2, const double maxDeltaR)
 

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_
 
DQMStoredbe_
 
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_
 

Detailed Description

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

Date:
2011/05/27 13:38:18
Revision:
1.18
Author
J. Slaunwhite, Jeff Klukas

Definition at line 66 of file HLTMuonMatchAndPlot.h.

Constructor & Destructor Documentation

HLTMuonMatchAndPlot::HLTMuonMatchAndPlot ( const edm::ParameterSet ,
std::string  ,
std::vector< std::string >   
)

Constructor.

Definition at line 44 of file HLTMuonMatchAndPlot.cc.

References binParams_, gather_cfg::cout, cutMinPt_, dbe_, fillMapFromPSet(), hltPath_, hltProcessName_, inputTags_, edm::InputTag::instance(), edm::InputTag::label(), moduleLabels_, cppFunctionSkipper::operator, plotCuts_, and triggerLevel_.

46  :
47  hltProcessName_(pset.getParameter<string>("hltProcessName")),
48  destination_(pset.getUntrackedParameter<string>("destination")),
49  requiredTriggers_(pset.getUntrackedParameter<vstring>("requiredTriggers")),
50  targetParams_(pset.getParameterSet("targetParams")),
51  probeParams_(pset.getParameterSet("probeParams")),
52  hltPath_(hltPath),
53  moduleLabels_(moduleLabels),
59  probeMuonSelector_(probeParams_.getUntrackedParameter<string>("recoCuts", "")),
60  probeZ0Cut_(probeParams_.getUntrackedParameter<double>("z0Cut",0.)),
62 {
63 
64  // Create std::map<string, T> from ParameterSets.
65  fillMapFromPSet(inputTags_, pset, "inputTags");
66  fillMapFromPSet(binParams_, pset, "binParams");
67  fillMapFromPSet(plotCuts_, pset, "plotCuts");
68 
69  // Set HLT process name for TriggerResults and TriggerSummary.
70  InputTag & resTag = inputTags_["triggerResults"];
71  resTag = InputTag(resTag.label(), resTag.instance(), hltProcessName_);
72  InputTag & sumTag = inputTags_["triggerSummary"];
73  sumTag = InputTag(sumTag.label(), sumTag.instance(), hltProcessName_);
74 
75  // Prepare the DQMStore object.
77  dbe_->setVerbose(0);
78 
79  // Get the trigger level.
80  triggerLevel_ = "L3";
81  TPRegexp levelRegexp("L[1-3]");
82  size_t nModules = moduleLabels_.size();
83  TObjArray * levelArray = levelRegexp.MatchS(moduleLabels_[nModules - 1]);
84  std::cout << moduleLabels_[nModules - 1] << std::endl;
85  if (levelArray->GetEntriesFast() > 0) {
86  triggerLevel_ = ((TObjString *)levelArray->At(0))->GetString();
87  }
88  delete levelArray;
89 
90  // Get the pT cut by parsing the name of the HLT path.
91  cutMinPt_ = 3;
92  TPRegexp ptRegexp("Mu([0-9]*)");
93  TObjArray * objArray = ptRegexp.MatchS(hltPath_);
94  if (objArray->GetEntriesFast() >= 2) {
95  TObjString * ptCutString = (TObjString *)objArray->At(1);
96  cutMinPt_ = atoi(ptCutString->GetString());
97  cutMinPt_ = ceil(cutMinPt_ * plotCuts_["minPtFactor"]);
98  }
99  delete objArray;
100 }
T getUntrackedParameter(std::string const &, T const &) const
StringCutObjectSelector< reco::Muon > probeMuonSelector_
edm::ParameterSet targetParams_
std::vector< std::string > moduleLabels_
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::map< std::string, std::vector< double > > binParams_
std::vector< std::string > requiredTriggers_
void setVerbose(unsigned level)
Definition: DQMStore.cc:393
StringCutObjectSelector< reco::Muon > targetMuonSelector_
std::map< std::string, edm::InputTag > inputTags_
std::string const & label() const
Definition: InputTag.h:25
std::vector< std::string > vstring
edm::ParameterSet probeParams_
void fillMapFromPSet(std::map< std::string, T > &, edm::ParameterSet, std::string)
tuple cout
Definition: gather_cfg.py:121
std::map< std::string, double > plotCuts_
std::string const & instance() const
Definition: InputTag.h:26

Member Function Documentation

void HLTMuonMatchAndPlot::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 166 of file HLTMuonMatchAndPlot.cc.

References allMuons_cfi::allMuons, SiPixelRawToDigiRegional_cfi::beamSpot, reco::LeafCandidate::charge(), cutMinPt_, deltaR(), reco::TrackBase::dxy(), reco::TrackBase::dz(), EFFICIENCY_SUFFIXES, trigger::TriggerObject::eta(), reco::LeafCandidate::eta(), eta(), edm::Event::getByLabel(), hasProbeRecoCuts, hasTargetRecoCuts, hists_, HLT_2011_cff::hltMuons, i, reco::Muon::innerTrack(), inputTags_, reco::Muon::isStandAloneMuon(), reco::Muon::isTrackerMuon(), edm::HandleBase::isValid(), j, gen::k, scaleCards::mass, matchByDeltaR(), metsig::muon, reco::Muon::outerTrack(), reco::LeafCandidate::p4(), phi, trigger::TriggerObject::phi(), reco::LeafCandidate::phi(), plotCuts_, probeD0Cut_, probeMuonSelector_, probeZ0Cut_, trigger::TriggerObject::pt(), reco::LeafCandidate::pt(), requiredTriggers_, selectedMuons(), selectedTriggerObjects(), findQualityFiles::size, createPayload::suffix, targetD0Cut_, targetMuonSelector_, targetParams_, targetZ0Cut_, triggerLevel_, and patRefSel_triggerSelection_cff::triggerResults.

169 {
170 
171  // Get objects from the event.
173  iEvent.getByLabel(inputTags_["recoMuon"], allMuons);
174 
176  iEvent.getByLabel(inputTags_["beamSpot"], beamSpot);
177 
178  Handle<TriggerEvent> triggerSummary;
179  iEvent.getByLabel(inputTags_["triggerSummary"], triggerSummary);
181 
183  iEvent.getByLabel("offlinePrimaryVertices", vertices);
184 
185  //edm::Handle<GenParticleCollection> gen;
186  //iEvent.getByLabel("genParticles", gen);
187  //GenParticleCollection::const_iterator g_part;
188  //std::vector<const GenParticle *> gen_leptsp;
189  //std::vector<const GenParticle *> gen_momsp;
190 
191 
192  if(!triggerSummary.isValid())
193  {
194  edm::LogError("HLTMuonMatchAndPlot")<<"Missing triggerSummary with label " << inputTags_["triggerSummary"] <<std::endl;
195  return;
196  }
197  iEvent.getByLabel(inputTags_["triggerResults"], triggerResults);
198  if(!triggerResults.isValid())
199  {
200  edm::LogError("HLTMuonMatchAndPlot")<<"Missing triggerResults with label " << inputTags_["triggerResults"] <<std::endl;
201  return;
202  }
203 
204  /*
205  if(gen != 0) {
206  for(g_part = gen->begin(); g_part != gen->end(); g_part++){
207  if( abs(g_part->pdgId()) == 13) {
208  if(!( g_part->status() ==1 || (g_part->status() ==2 && abs(g_part->pdgId())==5))) continue;
209  bool GenMomExists (true);
210  bool GenGrMomExists(true);
211  if( g_part->pt() < 10.0 ) continue;
212  if( fabs(g_part->eta()) > 2.4 ) continue;
213  int gen_id= g_part->pdgId();
214  const GenParticle* gen_lept = &(*g_part);
215  // get mother of gen_lept
216  const GenParticle* gen_mom = static_cast<const GenParticle*> (gen_lept->mother());
217  if(gen_mom==NULL) GenMomExists=false;
218  int m_id=-999;
219  if(GenMomExists) m_id = gen_mom->pdgId();
220  if(m_id != gen_id || !GenMomExists);
221  else{
222  int id= m_id;
223  while(id == gen_id && GenMomExists){
224  gen_mom = static_cast<const GenParticle*> (gen_mom->mother());
225  if(gen_mom==NULL){
226  GenMomExists=false;
227  }
228  if(GenMomExists) id=gen_mom->pdgId();
229  }
230  }
231  if(GenMomExists) m_id = gen_mom->pdgId();
232  gen_leptsp.push_back(gen_lept);
233  gen_momsp.push_back(gen_mom);
234  }
235  }
236  }
237 
238 
239  std::vector<GenParticle> gen_lepts;
240  for(size_t i = 0; i < gen_leptsp.size(); i++) {
241  gen_lepts.push_back(*gen_leptsp[i]);
242  }
243 
244 
245  */
246 
247 
248 
249  // Throw out this event if it doesn't pass the required triggers.
250  for (size_t i = 0; i < requiredTriggers_.size(); i++) {
251  unsigned int triggerIndex = triggerResults->find(requiredTriggers_[i]);
252  if (triggerIndex < triggerResults->size() ||
253  !triggerResults->accept(triggerIndex))
254  return;
255  }
256 
257  // Select objects based on the configuration.
260  TriggerObjectCollection allTriggerObjects = triggerSummary->getObjects();
262  selectedTriggerObjects(allTriggerObjects, * triggerSummary, targetParams_);
263 
264  // Fill plots for HLT muons.
265  for (size_t i = 0; i < hltMuons.size(); i++) {
266  hists_["hltPt"]->Fill(hltMuons[i].pt());
267  hists_["hltEta"]->Fill(hltMuons[i].eta());
268  hists_["hltPhi"]->Fill(hltMuons[i].phi());
269  }
270 
271  // Find the best trigger object matches for the targetMuons.
272  vector<size_t> matches = matchByDeltaR(targetMuons, hltMuons,
273  plotCuts_[triggerLevel_ + "DeltaR"]);
274 
275 
276  // Fill plots for matched muons.
277  bool pairalreadyconsidered = false;
278  for (size_t i = 0; i < targetMuons.size(); i++) {
279 
280  Muon & muon = targetMuons[i];
281 
282  // Fill plots which are not efficiencies.
283  if (matches[i] < targetMuons.size()) {
284  TriggerObject & hltMuon = hltMuons[matches[i]];
285  double ptRes = (muon.pt() - hltMuon.pt()) / muon.pt();
286  double etaRes = muon.eta() - hltMuon.eta();
287  double phiRes = muon.phi() - hltMuon.phi();
288  hists_["resolutionEta"]->Fill(etaRes);
289  hists_["resolutionPhi"]->Fill(phiRes);
290  hists_["resolutionPt"]->Fill(ptRes);
291  hists_["deltaR"]->Fill(deltaR(muon, hltMuon));
292  }
293 
294  // Fill numerators and denominators for efficiency plots.
295  for (size_t j = 0; j < 2; j++) {
296 
297  string suffix = EFFICIENCY_SUFFIXES[j];
298 
299  // If no match was found, then the numerator plots don't get filled.
300  if (suffix == "numer" && matches[i] >= targetMuons.size()) continue;
301 
302  if (muon.pt() > cutMinPt_) {
303  hists_["efficiencyEta_" + suffix]->Fill(muon.eta());
304  hists_["efficiencyPhiVsEta_" + suffix]->Fill(muon.eta(), muon.phi());
305  }
306 
307  if (fabs(muon.eta()) < plotCuts_["maxEta"]) {
308  hists_["efficiencyTurnOn_" + suffix]->Fill(muon.pt());
309  }
310 
311  if (muon.pt() > cutMinPt_ && fabs(muon.eta()) < plotCuts_["maxEta"]) {
312  const Track * track = 0;
313  if (muon.isTrackerMuon()) track = & * muon.innerTrack();
314  else if (muon.isStandAloneMuon()) track = & * muon.outerTrack();
315  if (track) {
316  double d0 = track->dxy(beamSpot->position());
317  double z0 = track->dz(beamSpot->position());
318  hists_["efficiencyVertex_" + suffix]->Fill(vertices->size());
319  hists_["efficiencyPhi_" + suffix]->Fill(muon.phi());
320  hists_["efficiencyD0_" + suffix]->Fill(d0);
321  hists_["efficiencyZ0_" + suffix]->Fill(z0);
322  hists_["efficiencyCharge_" + suffix]->Fill(muon.charge());
323  }
324  }
325  }
326 
327  // Fill plots for tag and probe
328  // Muon cannot be a tag because doesn't match an hlt muon
329  if(matches[i] >= targetMuons.size()) continue;
330  for (size_t k = 0; k < targetMuons.size(); k++) {
331  if(k == i) continue;
332  if(muon.pt() < 20.0) continue;
333  Muon & theProbe = targetMuons[k];
334  if (muon.charge() != theProbe.charge() && !pairalreadyconsidered) {
335  double mass = (muon.p4() + theProbe.p4()).M();
336  if(mass > 60 && mass < 120) {
337  hists_["massVsEtaZ_denom"]->Fill(theProbe.eta());
338  hists_["massVsPtZ_denom"]->Fill(theProbe.pt());
339  hists_["massVsVertexZ_denom"]->Fill(vertices->size());
340  if(matches[k] < targetMuons.size()) {
341  hists_["massVsEtaZ_numer"]->Fill(theProbe.eta());
342  hists_["massVsPtZ_numer"]->Fill(theProbe.pt());
343  hists_["massVsVertexZ_numer"]->Fill(vertices->size());
344  }
345  pairalreadyconsidered = true;
346  }
347  if(mass > 1 && mass < 4) {
348  hists_["massVsEtaJpsi_denom"]->Fill(theProbe.eta());
349  hists_["massVsPtJpsi_denom"]->Fill(theProbe.pt());
350  hists_["massVsVertexJpsi_denom"]->Fill(vertices->size());
351  if(matches[k] < targetMuons.size()) {
352  hists_["massVsEtaJpsi_numer"]->Fill(theProbe.eta());
353  hists_["massVsPtJpsi_numer"]->Fill(theProbe.pt());
354  hists_["massVsVertexJpsi_numer"]->Fill(vertices->size());
355  }
356  pairalreadyconsidered = true;
357  }
358  }
359  } // End loop over denominator and numerator.
360  } // End loop over targetMuons.
361 
362  // Plot fake rates (efficiency for HLT objects to not get matched to RECO).
363  vector<size_t> hltMatches = matchByDeltaR(hltMuons, targetMuons,
364  plotCuts_[triggerLevel_ + "DeltaR"]);
365  for (size_t i = 0; i < hltMuons.size(); i++) {
366  TriggerObject & hltMuon = hltMuons[i];
367  bool isFake = hltMatches[i] > hltMuons.size();
368  for (size_t j = 0; j < 2; j++) {
369  string suffix = EFFICIENCY_SUFFIXES[j];
370  // If match is found, then numerator plots should not get filled
371  if (suffix == "numer" && ! isFake) continue;
372  hists_["fakerateVertex_" + suffix]->Fill(vertices->size());
373  hists_["fakerateEta_" + suffix]->Fill(hltMuon.eta());
374  hists_["fakeratePhi_" + suffix]->Fill(hltMuon.phi());
375  hists_["fakerateTurnOn_" + suffix]->Fill(hltMuon.pt());
376  } // End loop over numerator and denominator.
377  } // End loop over hltMuons.
378 
379 
380 } // End analyze() method.
reco::MuonCollection selectedMuons(const reco::MuonCollection &, const reco::BeamSpot &, bool, const StringCutObjectSelector< reco::Muon > &, double, double)
int i
Definition: DBlmapReader.cc:9
StringCutObjectSelector< reco::Muon > probeMuonSelector_
const std::string EFFICIENCY_SUFFIXES[2]
virtual TrackRef innerTrack() const
Definition: Muon.h:49
float phi() const
Definition: TriggerObject.h:60
edm::ParameterSet targetParams_
bool isTrackerMuon() const
Definition: Muon.h:212
std::map< std::string, MonitorElement * > hists_
float eta() const
Definition: TriggerObject.h:59
T eta() const
bool isStandAloneMuon() const
Definition: Muon.h:213
virtual double eta() const
momentum pseudorapidity
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
std::vector< std::string > requiredTriggers_
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:24
std::vector< size_t > matchByDeltaR(const std::vector< T1 > &, const std::vector< T2 > &, const double maxDeltaR=NOMATCH)
virtual int charge() const
electric charge
int j
Definition: DBlmapReader.cc:9
tuple allMuons
Definition: allMuons_cfi.py:3
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual TrackRef outerTrack() const
reference to Track reconstructed in the muon detector only
Definition: Muon.h:52
int k[5][pyjets_maxn]
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:127
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:83
StringCutObjectSelector< reco::Muon > targetMuonSelector_
virtual double pt() const
transverse momentum
std::map< std::string, edm::InputTag > inputTags_
trigger::TriggerObjectCollection selectedTriggerObjects(const trigger::TriggerObjectCollection &, const trigger::TriggerEvent &, const edm::ParameterSet &)
tuple mass
Definition: scaleCards.py:27
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
virtual double phi() const
momentum azimuthal angle
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
std::map< std::string, double > plotCuts_
tuple size
Write out results.
Definition: DDAxes.h:10
void HLTMuonMatchAndPlot::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)

Definition at line 102 of file HLTMuonMatchAndPlot.cc.

References book1D(), book2D(), dbe_, destination_, EFFICIENCY_SUFFIXES, hltPath_, i, DQMStore::setCurrentFolder(), and createPayload::suffix.

104 {
105 
106  TPRegexp suffixPtCut("Mu[0-9]+$");
107 
108  string baseDir = destination_;
109  if (baseDir[baseDir.size() - 1] != '/') baseDir += '/';
110  string pathSansSuffix = hltPath_;
111  if (hltPath_.rfind("_v") < hltPath_.length())
112  pathSansSuffix = hltPath_.substr(0, hltPath_.rfind("_v"));
113  dbe_->setCurrentFolder(baseDir + pathSansSuffix);
114 
115  // Form is book1D(name, binningType, title) where 'binningType' is used
116  // to fetch the bin settings from binParams_.
117  book1D("deltaR", "deltaR", ";#Deltar(reco, HLT);");
118  book1D("hltPt", "pt", ";p_{T} of HLT object");
119  book1D("hltEta", "eta", ";#eta of HLT object");
120  book1D("hltPhi", "phi", ";#phi of HLT object");
121  book1D("resolutionEta", "resolutionEta", ";#eta^{reco}-#eta^{HLT};");
122  book1D("resolutionPhi", "resolutionPhi", ";#phi^{reco}-#phi^{HLT};");
123  book1D("resolutionPt", "resolutionRel",
124  ";(p_{T}^{reco}-p_{T}^{HLT})/|p_{T}^{reco}|;");
125 
126  for (size_t i = 0; i < 2; i++) {
127 
128  string suffix = EFFICIENCY_SUFFIXES[i];
129 
130  book1D("efficiencyEta_" + suffix, "eta", ";#eta;");
131  book1D("efficiencyPhi_" + suffix, "phi", ";#phi;");
132  book1D("efficiencyTurnOn_" + suffix, "pt", ";p_{T};");
133  book1D("efficiencyD0_" + suffix, "d0", ";d0;");
134  book1D("efficiencyZ0_" + suffix, "z0", ";z0;");
135  book1D("efficiencyCharge_" + suffix, "charge", ";charge;");
136  book1D("efficiencyVertex_" + suffix, "NVertex", ";NVertex;");
137 
138  book2D("efficiencyPhiVsEta_" + suffix, "etaCoarse", "phiCoarse",
139  ";#eta;#phi");
140 
141  book1D("fakerateEta_" + suffix, "eta", ";#eta;");
142  book1D("fakerateVertex_" + suffix, "NVertex", ";NVertex;");
143  book1D("fakeratePhi_" + suffix, "phi", ";#phi;");
144  book1D("fakerateTurnOn_" + suffix, "pt", ";p_{T};");
145 
146  book1D("massVsEtaZ_" + suffix, "etaCoarse", ";#eta");
147  book1D("massVsEtaJpsi_" + suffix, "etaCoarse", ";#eta");
148  book1D("massVsPtZ_" + suffix, "ptCoarse", ";p_{T}");
149  book1D("massVsPtJpsi_" + suffix, "ptCoarse", ";p_{T}");
150  book1D("massVsVertexZ_" + suffix, "NVertex", ";NVertex");
151  book1D("massVsVertexJpsi_" + suffix, "NVertex", ";NVertex");
152 
153  }
154 
155 }
int i
Definition: DBlmapReader.cc:9
const std::string EFFICIENCY_SUFFIXES[2]
void book1D(std::string, std::string, std::string)
void book2D(std::string, std::string, std::string, std::string)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
void HLTMuonMatchAndPlot::book1D ( std::string  ,
std::string  ,
std::string   
)
private

Definition at line 555 of file HLTMuonMatchAndPlot.cc.

References binParams_, DQMStore::book1D(), dbe_, prof2calltree::edges, fillEdges(), getTH1F(), hists_, and mergeVDriftHistosByStation::name.

Referenced by beginRun().

556 {
557  /* Properly delete the array of floats that has been allocated on
558  * the heap by fillEdges. Avoid multiple copies and internal ROOT
559  * clones by simply creating the histograms directly in the DQMStore
560  * using the appropriate book1D method to handle the variable bins
561  * case. */
562 
563  size_t nBins;
564  float * edges = 0;
565  fillEdges(nBins, edges, binParams_[binningType]);
566 
567  hists_[name] = dbe_->book1D(name, title, nBins, edges);
568  if (hists_[name])
569  if (hists_[name]->getTH1F()->GetSumw2N())
570  hists_[name]->getTH1F()->Sumw2();
571 
572  if (edges)
573  delete [] edges;
574 }
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
void fillEdges(size_t &nBins, float *&edges, std::vector< double > binning)
std::map< std::string, std::vector< double > > binParams_
std::map< std::string, MonitorElement * > hists_
dictionary edges
TH1F * getTH1F(std::string name, std::string process, std::string rootfolder, DQMStore *dbe_, bool verb, bool clone)
void HLTMuonMatchAndPlot::book2D ( std::string  ,
std::string  ,
std::string  ,
std::string   
)
private

Definition at line 579 of file HLTMuonMatchAndPlot.cc.

References binParams_, DQMStore::book2D(), dbe_, fillEdges(), getTH2F(), hists_, and mergeVDriftHistosByStation::name.

Referenced by beginRun().

581 {
582 
583  /* Properly delete the arrays of floats that have been allocated on
584  * the heap by fillEdges. Avoid multiple copies and internal ROOT
585  * clones by simply creating the histograms directly in the DQMStore
586  * using the appropriate book2D method to handle the variable bins
587  * case. */
588 
589  size_t nBinsX;
590  float * edgesX = 0;
591  fillEdges(nBinsX, edgesX, binParams_[binningTypeX]);
592 
593  size_t nBinsY;
594  float * edgesY = 0;
595  fillEdges(nBinsY, edgesY, binParams_[binningTypeY]);
596 
597  hists_[name] = dbe_->book2D(name.c_str(), title.c_str(),
598  nBinsX, edgesX, nBinsY, edgesY);
599  if (hists_[name])
600  if (hists_[name]->getTH2F()->GetSumw2N())
601  hists_[name]->getTH2F()->Sumw2();
602 
603  if (edgesX)
604  delete [] edgesX;
605  if (edgesY)
606  delete [] edgesY;
607 }
void fillEdges(size_t &nBins, float *&edges, std::vector< double > binning)
std::map< std::string, std::vector< double > > binParams_
std::map< std::string, MonitorElement * > hists_
TH2F * getTH2F(std::string name, std::string process, std::string rootfolder, DQMStore *dbe_, bool verb, bool clone)
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:845
void HLTMuonMatchAndPlot::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)

Definition at line 159 of file HLTMuonMatchAndPlot.cc.

161 {
162 }
void HLTMuonMatchAndPlot::fillEdges ( size_t &  nBins,
float *&  edges,
std::vector< double >  binning 
)

Definition at line 386 of file HLTMuonMatchAndPlot.cc.

References i, and min.

Referenced by book1D(), and book2D().

387  {
388 
389  if (binning.size() < 3) {
390  LogWarning("HLTMuonVal") << "Invalid binning parameters!";
391  return;
392  }
393 
394  // Fixed-width binning.
395  if (binning.size() == 3) {
396  nBins = binning[0];
397  edges = new float[nBins + 1];
398  const double min = binning[1];
399  const double binwidth = (binning[2] - binning[1]) / nBins;
400  for (size_t i = 0; i <= nBins; i++) edges[i] = min + (binwidth * i);
401  }
402 
403  // Variable-width binning.
404  else {
405  nBins = binning.size() - 1;
406  edges = new float[nBins + 1];
407  for (size_t i = 0; i <= nBins; i++) edges[i] = binning[i];
408  }
409 
410 }
int i
Definition: DBlmapReader.cc:9
#define min(a, b)
Definition: mlp_lapack.h:161
dictionary edges
template<class T >
void HLTMuonMatchAndPlot::fillMapFromPSet ( std::map< std::string, T > &  ,
edm::ParameterSet  ,
std::string   
)

Referenced by HLTMuonMatchAndPlot().

template<class T >
void HLTMuonMatchAndPlot::fillMapFromPSet ( map< string, T > &  m,
ParameterSet  pset,
string  target 
)

Definition at line 419 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.

420  {
421 
422  // Get the ParameterSet with name 'target' from 'pset'
423  ParameterSet targetPset;
424  if (pset.existsAs<ParameterSet>(target, true)) // target is tracked
425  targetPset = pset.getParameterSet(target);
426  else if (pset.existsAs<ParameterSet>(target, false)) // target is untracked
427  targetPset = pset.getUntrackedParameterSet(target);
428 
429  // Get the parameter names from targetPset
430  vector<string> names = targetPset.getParameterNames();
431  vector<string>::const_iterator iter;
432 
433  for (iter = names.begin(); iter != names.end(); ++iter) {
434  if (targetPset.existsAs<T>(* iter, true)) // target is tracked
435  m[* iter] = targetPset.getParameter<T>(* iter);
436  else if (targetPset.existsAs<T>(* iter, false)) // target is untracked
437  m[* iter] = targetPset.getUntrackedParameter<T>(* iter);
438  }
439 
440 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:187
static const HistoName names[]
ParameterSet const & getUntrackedParameterSet(std::string const &name, ParameterSet const &defaultValue) const
std::vector< std::string > getParameterNames() const
ParameterSet const & getParameterSet(std::string const &) const
long double T
template<class T1 , class T2 >
std::vector<size_t> HLTMuonMatchAndPlot::matchByDeltaR ( const std::vector< T1 > &  ,
const std::vector< T2 > &  ,
const double  maxDeltaR = NOMATCH 
)

Referenced by analyze().

template<class T1 , class T2 >
vector<size_t> HLTMuonMatchAndPlot::matchByDeltaR ( const vector< T1 > &  collection1,
const vector< T2 > &  collection2,
const double  maxDeltaR 
)

Definition at line 447 of file HLTMuonMatchAndPlot.cc.

References deltaR(), i, j, gen::k, NOMATCH, and query::result.

449  {
450 
451  const size_t n1 = collection1.size();
452  const size_t n2 = collection2.size();
453 
454  vector<size_t> result(n1, -1);
455  vector<vector<double> > deltaRMatrix(n1, vector<double>(n2, NOMATCH));
456 
457  for (size_t i = 0; i < n1; i++)
458  for (size_t j = 0; j < n2; j++) {
459  deltaRMatrix[i][j] = deltaR(collection1[i], collection2[j]);
460  }
461 
462  // Run through the matrix n1 times to make sure we've found all matches.
463  for (size_t k = 0; k < n1; k++) {
464  size_t i_min = -1;
465  size_t j_min = -1;
466  double minDeltaR = maxDeltaR;
467  // find the smallest deltaR
468  for (size_t i = 0; i < n1; i++)
469  for (size_t j = 0; j < n2; j++)
470  if (deltaRMatrix[i][j] < minDeltaR) {
471  i_min = i;
472  j_min = j;
473  minDeltaR = deltaRMatrix[i][j];
474  }
475  // If a match has been made, save it and make those candidates unavailable.
476  if (minDeltaR < maxDeltaR) {
477  result[i_min] = j_min;
478  deltaRMatrix[i_min] = vector<double>(n2, NOMATCH);
479  for (size_t i = 0; i < n1; i++)
480  deltaRMatrix[i][j_min] = NOMATCH;
481  }
482  }
483 
484  return result;
485 
486 }
int i
Definition: DBlmapReader.cc:9
const double NOMATCH
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
int k[5][pyjets_maxn]
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
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 491 of file HLTMuonMatchAndPlot.cc.

References reco::TrackBase::dxy(), reco::TrackBase::dz(), and reco::BeamSpot::position().

Referenced by analyze().

496 {
497  // If there is no selector (recoCuts does not exists), return an empty collection.
498  if (!hasRecoCuts)
499  return MuonCollection();
500 
501  MuonCollection reducedMuons(allMuons);
502  MuonCollection::iterator iter = reducedMuons.begin();
503  while (iter != reducedMuons.end()) {
504  const Track * track = 0;
505  if (iter->isTrackerMuon()) track = & * iter->innerTrack();
506  else if (iter->isStandAloneMuon()) track = & * iter->outerTrack();
507  if (track && selector(* iter) &&
508  fabs(track->dxy(beamSpot.position())) < d0Cut &&
509  fabs(track->dz(beamSpot.position())) < z0Cut)
510  ++iter;
511  else reducedMuons.erase(iter);
512  }
513 
514  return reducedMuons;
515 
516 }
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
tuple allMuons
Definition: allMuons_cfi.py:3
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:127
const Point & position() const
position
Definition: BeamSpot.h:63
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
TriggerObjectCollection HLTMuonMatchAndPlot::selectedTriggerObjects ( const trigger::TriggerObjectCollection triggerObjects,
const trigger::TriggerEvent triggerSummary,
const edm::ParameterSet pset 
)
private

Definition at line 521 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().

525 {
526 
527  // If pset is empty, return an empty collection.
528  if (!pset.exists("hltCuts"))
529  return TriggerObjectCollection();
530 
532  (pset.getUntrackedParameter<string>("hltCuts"));
533 
534  InputTag filterTag(moduleLabels_[moduleLabels_.size() - 1], "",
536  size_t filterIndex = triggerSummary.filterIndex(filterTag);
537 
538  TriggerObjectCollection selectedObjects;
539 
540  if (filterIndex < triggerSummary.sizeFilters()) {
541  const Keys &keys = triggerSummary.filterKeys(filterIndex);
542  for (size_t j = 0; j < keys.size(); j++ ){
543  TriggerObject foundObject = triggerObjects[keys[j]];
544  if (selector(foundObject))
545  selectedObjects.push_back(foundObject);
546  }
547  }
548 
549  return selectedObjects;
550 
551 }
T getUntrackedParameter(std::string const &, T const &) const
trigger::size_type sizeFilters() const
Definition: TriggerEvent.h:137
std::vector< std::string > moduleLabels_
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:113
trigger::size_type filterIndex(const edm::InputTag &filterTag) const
find index of filter in data-member vector from filter tag
Definition: TriggerEvent.h:125
bool exists(std::string const &parameterName) const
checks if a parameter exists
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:24
int j
Definition: DBlmapReader.cc:9
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:83
std::vector< size_type > Keys

Member Data Documentation

std::map<std::string, std::vector<double> > HLTMuonMatchAndPlot::binParams_
private

Definition at line 108 of file HLTMuonMatchAndPlot.h.

Referenced by book1D(), 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(), book1D(), 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(), book1D(), and 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().

StringCutObjectSelector<reco::Muon> HLTMuonMatchAndPlot::probeMuonSelector_
private

Definition at line 128 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

edm::ParameterSet HLTMuonMatchAndPlot::probeParams_
private

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().

StringCutObjectSelector<reco::Muon> HLTMuonMatchAndPlot::targetMuonSelector_
private

Definition at line 125 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

edm::ParameterSet HLTMuonMatchAndPlot::targetParams_
private

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().