CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTMuonMatchAndPlot.cc
Go to the documentation of this file.
1 
7 
15 
16 #include <iostream>
17 
20 
21 using namespace std;
22 using namespace edm;
23 using namespace reco;
24 using namespace trigger;
25 using namespace l1extra;
26 
27 typedef std::vector<std::string> vstring;
28 
29 
32 
35  string hltPath,
36  const vector<string>& moduleLabels) :
37  hltProcessName_(pset.getParameter<string>("hltProcessName")),
38  destination_(pset.getUntrackedParameter<string>("destination")),
39  requiredTriggers_(pset.getUntrackedParameter<vstring>("requiredTriggers")),
40  targetParams_(pset.getParameterSet("targetParams")),
41  probeParams_(pset.getParameterSet("probeParams")),
42  hltPath_(hltPath),
43  moduleLabels_(moduleLabels),
44  hasTargetRecoCuts(targetParams_.exists("recoCuts")),
45  hasProbeRecoCuts(probeParams_.exists("recoCuts")),
46  targetMuonSelector_(targetParams_.getUntrackedParameter<string>("recoCuts", "")),
47  targetZ0Cut_(targetParams_.getUntrackedParameter<double>("z0Cut",0.)),
48  targetD0Cut_(targetParams_.getUntrackedParameter<double>("d0Cut",0.)),
49  probeMuonSelector_(probeParams_.getUntrackedParameter<string>("recoCuts", "")),
50  probeZ0Cut_(probeParams_.getUntrackedParameter<double>("z0Cut",0.)),
51  probeD0Cut_(probeParams_.getUntrackedParameter<double>("d0Cut",0.)),
52  triggerSelector_(targetParams_.getUntrackedParameter<string>("hltCuts","")),
53  hasTriggerCuts_(targetParams_.exists("hltCuts"))
54 {
55  // Create std::map<string, T> from ParameterSets.
56  fillMapFromPSet(binParams_, pset, "binParams");
57  fillMapFromPSet(plotCuts_, pset, "plotCuts");
58 
59  // Get the trigger level.
60  triggerLevel_ = "L3";
61  TPRegexp levelRegexp("L[1-3]");
62  size_t nModules = moduleLabels_.size();
63  TObjArray * levelArray = levelRegexp.MatchS(moduleLabels_[nModules - 1]);
64  if (levelArray->GetEntriesFast() > 0) {
65  triggerLevel_ = ((TObjString *)levelArray->At(0))->GetString();
66  }
67  delete levelArray;
68 
69  // Get the pT cut by parsing the name of the HLT path.
70  cutMinPt_ = 3;
71  TPRegexp ptRegexp("Mu([0-9]*)");
72  TObjArray * objArray = ptRegexp.MatchS(hltPath_);
73  if (objArray->GetEntriesFast() >= 2) {
74  TObjString * ptCutString = (TObjString *)objArray->At(1);
75  cutMinPt_ = atoi(ptCutString->GetString());
76  cutMinPt_ = ceil(cutMinPt_ * plotCuts_["minPtFactor"]);
77  }
78  delete objArray;
79 
80 }
81 
83  const edm::Run& iRun,
84  const edm::EventSetup& iSetup)
85 {
86 
87  TPRegexp suffixPtCut("Mu[0-9]+$");
88 
89  string baseDir = destination_;
90  if (baseDir[baseDir.size() - 1] != '/') baseDir += '/';
91  string pathSansSuffix = hltPath_;
92  if (hltPath_.rfind("_v") < hltPath_.length())
93  pathSansSuffix = hltPath_.substr(0, hltPath_.rfind("_v"));
94  iBooker.setCurrentFolder(baseDir + pathSansSuffix);
95 
96  // Form is book1D(name, binningType, title) where 'binningType' is used
97  // to fetch the bin settings from binParams_.
98  book1D(iBooker, "deltaR", "deltaR", ";#Deltar(reco, HLT);");
99  book1D(iBooker, "hltPt", "pt", ";p_{T} of HLT object");
100  book1D(iBooker, "hltEta", "eta", ";#eta of HLT object");
101  book1D(iBooker, "hltPhi", "phi", ";#phi of HLT object");
102  book1D(iBooker, "resolutionEta", "resolutionEta", ";#eta^{reco}-#eta^{HLT};");
103  book1D(iBooker, "resolutionPhi", "resolutionPhi", ";#phi^{reco}-#phi^{HLT};");
104  book1D(iBooker, "resolutionPt", "resolutionRel",
105  ";(p_{T}^{reco}-p_{T}^{HLT})/|p_{T}^{reco}|;");
106 
107  for (size_t i = 0; i < 2; i++) {
108 
109  string suffix = EFFICIENCY_SUFFIXES[i];
110 
111  book1D(iBooker, "efficiencyEta_" + suffix, "eta", ";#eta;");
112  book1D(iBooker, "efficiencyPhi_" + suffix, "phi", ";#phi;");
113  book1D(iBooker, "efficiencyTurnOn_" + suffix, "pt", ";p_{T};");
114  book1D(iBooker, "efficiencyD0_" + suffix, "d0", ";d0;");
115  book1D(iBooker, "efficiencyZ0_" + suffix, "z0", ";z0;");
116  book1D(iBooker, "efficiencyCharge_" + suffix, "charge", ";charge;");
117  book1D(iBooker, "efficiencyVertex_" + suffix, "NVertex", ";NVertex;");
118 
119  book2D(iBooker, "efficiencyPhiVsEta_" + suffix, "etaCoarse",
120  "phiCoarse", ";#eta;#phi");
121 
122  book1D(iBooker, "fakerateEta_" + suffix, "eta", ";#eta;");
123  book1D(iBooker, "fakerateVertex_" + suffix, "NVertex", ";NVertex;");
124  book1D(iBooker, "fakeratePhi_" + suffix, "phi", ";#phi;");
125  book1D(iBooker, "fakerateTurnOn_" + suffix, "pt", ";p_{T};");
126 
127  book1D(iBooker, "massVsEtaZ_" + suffix, "etaCoarse", ";#eta");
128  book1D(iBooker, "massVsEtaJpsi_" + suffix, "etaCoarse", ";#eta");
129  book1D(iBooker, "massVsPtZ_" + suffix, "ptCoarse", ";p_{T}");
130  book1D(iBooker, "massVsPtJpsi_" + suffix, "ptCoarse", ";p_{T}");
131  book1D(iBooker, "massVsVertexZ_" + suffix, "NVertex", ";NVertex");
132  book1D(iBooker, "massVsVertexJpsi_" + suffix, "NVertex", ";NVertex");
133 
134  }
135 
136 }
137 
138 
139 
141  const edm::EventSetup& iSetup)
142 {
143 
144 }
145 
146 
147 
151  Handle<TriggerEvent> & triggerSummary,
153 {
154 
155  /*
156  if(gen != 0) {
157  for(g_part = gen->begin(); g_part != gen->end(); g_part++){
158  if( abs(g_part->pdgId()) == 13) {
159  if(!( g_part->status() ==1 || (g_part->status() ==2 && abs(g_part->pdgId())==5))) continue;
160  bool GenMomExists (true);
161  bool GenGrMomExists(true);
162  if( g_part->pt() < 10.0 ) continue;
163  if( fabs(g_part->eta()) > 2.4 ) continue;
164  int gen_id= g_part->pdgId();
165  const GenParticle* gen_lept = &(*g_part);
166  // get mother of gen_lept
167  const GenParticle* gen_mom = static_cast<const GenParticle*> (gen_lept->mother());
168  if(gen_mom==NULL) GenMomExists=false;
169  int m_id=-999;
170  if(GenMomExists) m_id = gen_mom->pdgId();
171  if(m_id != gen_id || !GenMomExists);
172  else{
173  int id= m_id;
174  while(id == gen_id && GenMomExists){
175  gen_mom = static_cast<const GenParticle*> (gen_mom->mother());
176  if(gen_mom==NULL){
177  GenMomExists=false;
178  }
179  if(GenMomExists) id=gen_mom->pdgId();
180  }
181  }
182  if(GenMomExists) m_id = gen_mom->pdgId();
183  gen_leptsp.push_back(gen_lept);
184  gen_momsp.push_back(gen_mom);
185  }
186  }
187  }
188 
189 
190  std::vector<GenParticle> gen_lepts;
191  for(size_t i = 0; i < gen_leptsp.size(); i++) {
192  gen_lepts.push_back(*gen_leptsp[i]);
193  }
194 
195 
196  */
197 
198 
199 
200  // Throw out this event if it doesn't pass the required triggers.
201  for (size_t i = 0; i < requiredTriggers_.size(); i++) {
202  unsigned int triggerIndex = triggerResults->find(requiredTriggers_[i]);
203  if (triggerIndex < triggerResults->size() ||
204  !triggerResults->accept(triggerIndex))
205  return;
206  }
207 
208 
209  // Select objects based on the configuration.
212  TriggerObjectCollection allTriggerObjects = triggerSummary->getObjects();
214  selectedTriggerObjects(allTriggerObjects, * triggerSummary, hasTriggerCuts_,triggerSelector_);
215 
216  // Fill plots for HLT muons.
217  for (size_t i = 0; i < hltMuons.size(); i++) {
218  hists_["hltPt"]->Fill(hltMuons[i].pt());
219  hists_["hltEta"]->Fill(hltMuons[i].eta());
220  hists_["hltPhi"]->Fill(hltMuons[i].phi());
221  }
222 
223  // Find the best trigger object matches for the targetMuons.
224  vector<size_t> matches = matchByDeltaR(targetMuons, hltMuons,
225  plotCuts_[triggerLevel_ + "DeltaR"]);
226 
227 
228  // Fill plots for matched muons.
229  bool pairalreadyconsidered = false;
230  for (size_t i = 0; i < targetMuons.size(); i++) {
231 
232  Muon & muon = targetMuons[i];
233 
234  // Fill plots which are not efficiencies.
235  if (matches[i] < targetMuons.size()) {
236  TriggerObject & hltMuon = hltMuons[matches[i]];
237  double ptRes = (muon.pt() - hltMuon.pt()) / muon.pt();
238  double etaRes = muon.eta() - hltMuon.eta();
239  double phiRes = muon.phi() - hltMuon.phi();
240  hists_["resolutionEta"]->Fill(etaRes);
241  hists_["resolutionPhi"]->Fill(phiRes);
242  hists_["resolutionPt"]->Fill(ptRes);
243  hists_["deltaR"]->Fill(deltaR(muon, hltMuon));
244  }
245 
246  // Fill numerators and denominators for efficiency plots.
247  for (size_t j = 0; j < 2; j++) {
248 
249  string suffix = EFFICIENCY_SUFFIXES[j];
250 
251  // If no match was found, then the numerator plots don't get filled.
252  if (suffix == "numer" && matches[i] >= targetMuons.size()) continue;
253 
254  if (muon.pt() > cutMinPt_) {
255  hists_["efficiencyEta_" + suffix]->Fill(muon.eta());
256  hists_["efficiencyPhiVsEta_" + suffix]->Fill(muon.eta(), muon.phi());
257  }
258 
259  if (fabs(muon.eta()) < plotCuts_["maxEta"]) {
260  hists_["efficiencyTurnOn_" + suffix]->Fill(muon.pt());
261  }
262 
263  if (muon.pt() > cutMinPt_ && fabs(muon.eta()) < plotCuts_["maxEta"]) {
264  const Track * track = 0;
265  if (muon.isTrackerMuon()) track = & * muon.innerTrack();
266  else if (muon.isStandAloneMuon()) track = & * muon.outerTrack();
267  if (track) {
268  double d0 = track->dxy(beamSpot->position());
269  double z0 = track->dz(beamSpot->position());
270  hists_["efficiencyVertex_" + suffix]->Fill(vertices->size());
271  hists_["efficiencyPhi_" + suffix]->Fill(muon.phi());
272  hists_["efficiencyD0_" + suffix]->Fill(d0);
273  hists_["efficiencyZ0_" + suffix]->Fill(z0);
274  hists_["efficiencyCharge_" + suffix]->Fill(muon.charge());
275  }
276  }
277  }
278 
279  // Fill plots for tag and probe
280  // Muon cannot be a tag because doesn't match an hlt muon
281  if(matches[i] >= targetMuons.size()) continue;
282  for (size_t k = 0; k < targetMuons.size(); k++) {
283  if(k == i) continue;
284  if(muon.pt() < 20.0) continue;
285  Muon & theProbe = targetMuons[k];
286  if (muon.charge() != theProbe.charge() && !pairalreadyconsidered) {
287  double mass = (muon.p4() + theProbe.p4()).M();
288  if(mass > 60 && mass < 120) {
289  hists_["massVsEtaZ_denom"]->Fill(theProbe.eta());
290  hists_["massVsPtZ_denom"]->Fill(theProbe.pt());
291  hists_["massVsVertexZ_denom"]->Fill(vertices->size());
292  if(matches[k] < targetMuons.size()) {
293  hists_["massVsEtaZ_numer"]->Fill(theProbe.eta());
294  hists_["massVsPtZ_numer"]->Fill(theProbe.pt());
295  hists_["massVsVertexZ_numer"]->Fill(vertices->size());
296  }
297  pairalreadyconsidered = true;
298  }
299  if(mass > 1 && mass < 4) {
300  hists_["massVsEtaJpsi_denom"]->Fill(theProbe.eta());
301  hists_["massVsPtJpsi_denom"]->Fill(theProbe.pt());
302  hists_["massVsVertexJpsi_denom"]->Fill(vertices->size());
303  if(matches[k] < targetMuons.size()) {
304  hists_["massVsEtaJpsi_numer"]->Fill(theProbe.eta());
305  hists_["massVsPtJpsi_numer"]->Fill(theProbe.pt());
306  hists_["massVsVertexJpsi_numer"]->Fill(vertices->size());
307  }
308  pairalreadyconsidered = true;
309  }
310  }
311  } // End loop over denominator and numerator.
312  } // End loop over targetMuons.
313 
314  // Plot fake rates (efficiency for HLT objects to not get matched to RECO).
315  vector<size_t> hltMatches = matchByDeltaR(hltMuons, targetMuons,
316  plotCuts_[triggerLevel_ + "DeltaR"]);
317  for (size_t i = 0; i < hltMuons.size(); i++) {
318  TriggerObject & hltMuon = hltMuons[i];
319  bool isFake = hltMatches[i] > hltMuons.size();
320  for (size_t j = 0; j < 2; j++) {
321  string suffix = EFFICIENCY_SUFFIXES[j];
322  // If match is found, then numerator plots should not get filled
323  if (suffix == "numer" && ! isFake) continue;
324  hists_["fakerateVertex_" + suffix]->Fill(vertices->size());
325  hists_["fakerateEta_" + suffix]->Fill(hltMuon.eta());
326  hists_["fakeratePhi_" + suffix]->Fill(hltMuon.phi());
327  hists_["fakerateTurnOn_" + suffix]->Fill(hltMuon.pt());
328  } // End loop over numerator and denominator.
329  } // End loop over hltMuons.
330 
331 
332 } // End analyze() method.
333 
334 
335 
336 // Method to fill binning parameters from a vector of doubles.
337 void
338 HLTMuonMatchAndPlot::fillEdges(size_t & nBins, float * & edges,
339  const vector<double>& binning)
340 {
341 
342  if (binning.size() < 3) {
343  LogWarning("HLTMuonVal") << "Invalid binning parameters!";
344  return;
345  }
346 
347  // Fixed-width binning.
348  if (binning.size() == 3) {
349  nBins = binning[0];
350  edges = new float[nBins + 1];
351  const double min = binning[1];
352  const double binwidth = (binning[2] - binning[1]) / nBins;
353  for (size_t i = 0; i <= nBins; i++) edges[i] = min + (binwidth * i);
354  }
355 
356  // Variable-width binning.
357  else {
358  nBins = binning.size() - 1;
359  edges = new float[nBins + 1];
360  for (size_t i = 0; i <= nBins; i++) edges[i] = binning[i];
361  }
362 
363 }
364 
365 
366 
367 // This is an unorthodox method of getting parameters, but cleaner in my mind
368 // It looks for parameter 'target' in the pset, and assumes that all entries
369 // have the same class (T), filling the values into a std::map.
370 template <class T>
371 void
373  const ParameterSet& pset, string target)
374 {
375 
376  // Get the ParameterSet with name 'target' from 'pset'
377  ParameterSet targetPset;
378  if (pset.existsAs<ParameterSet>(target, true)) // target is tracked
379  targetPset = pset.getParameterSet(target);
380  else if (pset.existsAs<ParameterSet>(target, false)) // target is untracked
381  targetPset = pset.getUntrackedParameterSet(target);
382 
383  // Get the parameter names from targetPset
384  vector<string> names = targetPset.getParameterNames();
385  vector<string>::const_iterator iter;
386 
387  for (iter = names.begin(); iter != names.end(); ++iter) {
388  if (targetPset.existsAs<T>(* iter, true)) // target is tracked
389  m[* iter] = targetPset.getParameter<T>(* iter);
390  else if (targetPset.existsAs<T>(* iter, false)) // target is untracked
391  m[* iter] = targetPset.getUntrackedParameter<T>(* iter);
392  }
393 
394 }
395 
396 
397 
398 // A generic method to find the best deltaR matches from 2 collections.
399 template <class T1, class T2>
400 vector<size_t>
401 HLTMuonMatchAndPlot::matchByDeltaR(const vector<T1> & collection1,
402  const vector<T2> & collection2,
403  const double maxDeltaR)
404 {
405 
406  const size_t n1 = collection1.size();
407  const size_t n2 = collection2.size();
408 
409  vector<size_t> result(n1, -1);
410  vector<vector<double> > deltaRMatrix(n1, vector<double>(n2, NOMATCH));
411 
412  for (size_t i = 0; i < n1; i++)
413  for (size_t j = 0; j < n2; j++) {
414  deltaRMatrix[i][j] = deltaR(collection1[i], collection2[j]);
415  }
416 
417  // Run through the matrix n1 times to make sure we've found all matches.
418  for (size_t k = 0; k < n1; k++) {
419  size_t i_min = -1;
420  size_t j_min = -1;
421  double minDeltaR = maxDeltaR;
422  // find the smallest deltaR
423  for (size_t i = 0; i < n1; i++)
424  for (size_t j = 0; j < n2; j++)
425  if (deltaRMatrix[i][j] < minDeltaR) {
426  i_min = i;
427  j_min = j;
428  minDeltaR = deltaRMatrix[i][j];
429  }
430  // If a match has been made, save it and make those candidates unavailable.
431  if (minDeltaR < maxDeltaR) {
432  result[i_min] = j_min;
433  deltaRMatrix[i_min] = vector<double>(n2, NOMATCH);
434  for (size_t i = 0; i < n1; i++)
435  deltaRMatrix[i][j_min] = NOMATCH;
436  }
437  }
438 
439  return result;
440 
441 }
442 
443 
444 
447  const BeamSpot & beamSpot,
448  bool hasRecoCuts,
449  const StringCutObjectSelector<reco::Muon> &selector,
450  double d0Cut, double z0Cut)
451 {
452 
453  // If there is no selector (recoCuts does not exists), return an empty collection.
454  if (!hasRecoCuts)
455  return MuonCollection();
456 
457  MuonCollection reducedMuons(allMuons);
458  MuonCollection::iterator iter = reducedMuons.begin();
459  while (iter != reducedMuons.end()) {
460  const Track * track = 0;
461  if (iter->isTrackerMuon()) track = & * iter->innerTrack();
462  else if (iter->isStandAloneMuon()) track = & * iter->outerTrack();
463  if (track && selector(* iter) &&
464  fabs(track->dxy(beamSpot.position())) < d0Cut &&
465  fabs(track->dz(beamSpot.position())) < z0Cut)
466  ++iter;
467  else reducedMuons.erase(iter);
468  }
469 
470  return reducedMuons;
471 
472 }
473 
474 
475 
478  const TriggerObjectCollection & triggerObjects,
479  const TriggerEvent & triggerSummary,
480  bool hasTriggerCuts,
481  const StringCutObjectSelector<TriggerObject> triggerSelector)
482 {
483  if ( !hasTriggerCuts) return TriggerObjectCollection();
484 
485  InputTag filterTag(moduleLabels_[moduleLabels_.size() - 1], "",
487  size_t filterIndex = triggerSummary.filterIndex(filterTag);
488 
489  TriggerObjectCollection selectedObjects;
490 
491  if (filterIndex < triggerSummary.sizeFilters()) {
492  const Keys &keys = triggerSummary.filterKeys(filterIndex);
493  for (size_t j = 0; j < keys.size(); j++ ){
494  TriggerObject foundObject = triggerObjects[keys[j]];
495  if (triggerSelector(foundObject))
496  selectedObjects.push_back(foundObject);
497  }
498  }
499 
500  return selectedObjects;
501 
502 }
503 
504 
505 
507  string binningType, string title)
508 {
509 
510  /* Properly delete the array of floats that has been allocated on
511  * the heap by fillEdges. Avoid multiple copies and internal ROOT
512  * clones by simply creating the histograms directly in the DQMStore
513  * using the appropriate book1D method to handle the variable bins
514  * case. */
515 
516  size_t nBins;
517  float * edges = 0;
518  fillEdges(nBins, edges, binParams_[binningType]);
519 
520  hists_[name] = iBooker.book1D(name, title, nBins, edges);
521  if (hists_[name])
522  if (hists_[name]->getTH1F()->GetSumw2N())
523  hists_[name]->getTH1F()->Sumw2();
524 
525  if (edges)
526  delete [] edges;
527 
528 }
529 
530 
531 
532 void
534  string binningTypeX, string binningTypeY,
535  string title)
536 {
537 
538  /* Properly delete the arrays of floats that have been allocated on
539  * the heap by fillEdges. Avoid multiple copies and internal ROOT
540  * clones by simply creating the histograms directly in the DQMStore
541  * using the appropriate book2D method to handle the variable bins
542  * case. */
543 
544  size_t nBinsX;
545  float * edgesX = 0;
546  fillEdges(nBinsX, edgesX, binParams_[binningTypeX]);
547 
548  size_t nBinsY;
549  float * edgesY = 0;
550  fillEdges(nBinsY, edgesY, binParams_[binningTypeY]);
551 
552  hists_[name] = iBooker.book2D(name.c_str(), title.c_str(),
553  nBinsX, edgesX, nBinsY, edgesY);
554  if (hists_[name])
555  if (hists_[name]->getTH2F()->GetSumw2N())
556  hists_[name]->getTH2F()->Sumw2();
557 
558  if (edgesX)
559  delete [] edgesX;
560  if (edgesY)
561  delete [] edgesY;
562 
563 }
564 
565 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
reco::MuonCollection selectedMuons(const reco::MuonCollection &, const reco::BeamSpot &, bool, const StringCutObjectSelector< reco::Muon > &, double, double)
int i
Definition: DBlmapReader.cc:9
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
StringCutObjectSelector< reco::Muon > probeMuonSelector_
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:25
trigger::size_type sizeFilters() const
Definition: TriggerEvent.h:135
static const HistoName names[]
const std::string EFFICIENCY_SUFFIXES[2]
vector< string > vstring
Definition: ExoticaDQM.cc:86
virtual TrackRef innerTrack() const
Definition: Muon.h:48
float phi() const
Definition: TriggerObject.h:58
void analyze(edm::Handle< reco::MuonCollection > &, edm::Handle< reco::BeamSpot > &, edm::Handle< reco::VertexCollection > &, edm::Handle< trigger::TriggerEvent > &, edm::Handle< edm::TriggerResults > &)
bool isTrackerMuon() const
Definition: Muon.h:219
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::vector< std::string > moduleLabels_
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:111
trigger::size_type filterIndex(const edm::InputTag &filterTag) const
find index of filter in data-member vector from filter tag
Definition: TriggerEvent.h:123
std::map< std::string, std::vector< double > > binParams_
ParameterSet getUntrackedParameterSet(std::string const &name, ParameterSet const &defaultValue) const
std::map< std::string, MonitorElement * > hists_
float eta() const
Definition: TriggerObject.h:57
bool isStandAloneMuon() const
Definition: Muon.h:220
dictionary edges
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
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:22
const double NOMATCH
trigger::TriggerObjectCollection selectedTriggerObjects(const trigger::TriggerObjectCollection &, const trigger::TriggerEvent &, bool hasTriggerCuts, const StringCutObjectSelector< trigger::TriggerObject > triggerSelector)
void fillMapFromPSet(std::map< std::string, T > &, const edm::ParameterSet &, std::string)
std::vector< size_t > matchByDeltaR(const std::vector< T1 > &, const std::vector< T2 > &, const double maxDeltaR=NOMATCH)
tuple result
Definition: query.py:137
virtual int charge() const
electric charge
Definition: LeafCandidate.h:91
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
int j
Definition: DBlmapReader.cc:9
tuple allMuons
Definition: allMuons_cfi.py:3
TH1F * getTH1F(std::string name, std::string process, std::string rootfolder, DQMStore *dbe_, bool verb, bool clone)
T min(T a, T b)
Definition: MathUtil.h:58
static std::string const triggerResults
Definition: EdmProvDump.cc:40
std::vector< std::string > getParameterNames() const
void book2D(DQMStore::IBooker &, std::string, std::string, std::string, std::string)
virtual TrackRef outerTrack() const
reference to Track reconstructed in the muon detector only
Definition: Muon.h:51
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:596
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:81
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
StringCutObjectSelector< reco::Muon > targetMuonSelector_
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
ParameterSet const & getParameterSet(std::string const &) const
std::vector< size_type > Keys
TH2F * getTH2F(std::string name, std::string process, std::string rootfolder, DQMStore *dbe_, bool verb, bool clone)
StringCutObjectSelector< trigger::TriggerObject > triggerSelector_
void endRun(const edm::Run &, const edm::EventSetup &)
void book1D(DQMStore::IBooker &, std::string, std::string, std::string)
const Point & position() const
position
Definition: BeamSpot.h:62
HLTMuonMatchAndPlot(const edm::ParameterSet &, std::string, const std::vector< std::string > &)
Constructor.
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:578
long double T
virtual double phi() const
momentum azimuthal angle
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
std::map< std::string, double > plotCuts_
void fillEdges(size_t &nBins, float *&edges, const std::vector< double > &binning)
tuple size
Write out results.
Definition: Run.h:43
void beginRun(DQMStore::IBooker &, const edm::Run &, const edm::EventSetup &)