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 (edm::Handle< reco::MuonCollection > &, edm::Handle< reco::BeamSpot > &, edm::Handle< reco::VertexCollection > &, edm::Handle< trigger::TriggerEvent > &, edm::Handle< edm::TriggerResults > &)
 
void beginRun (DQMStore::IBooker &, const edm::Run &, const edm::EventSetup &)
 
void endRun (const edm::Run &, const edm::EventSetup &)
 
void fillEdges (size_t &nBins, float *&edges, const std::vector< double > &binning)
 
template<class T >
void fillMapFromPSet (std::map< std::string, T > &, const edm::ParameterSet &, std::string)
 
template<class T >
void fillMapFromPSet (map< string, T > &m, const ParameterSet &pset, string target)
 
 HLTMuonMatchAndPlot (const edm::ParameterSet &, std::string, const 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 (DQMStore::IBooker &, std::string, std::string, std::string)
 
void book2D (DQMStore::IBooker &, 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 &, bool hasTriggerCuts, const StringCutObjectSelector< trigger::TriggerObject > triggerSelector)
 

Private Attributes

std::map< std::string,
std::vector< double > > 
binParams_
 
unsigned int cutMinPt_
 
std::string destination_
 
bool hasProbeRecoCuts
 
bool hasTargetRecoCuts
 
bool hasTriggerCuts_
 
std::map< std::string,
MonitorElement * > 
hists_
 
std::string hltPath_
 
std::string hltProcessName_
 
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 targetptCutJpsi_
 
double targetptCutZ_
 
double targetZ0Cut_
 
std::string triggerLevel_
 
StringCutObjectSelector
< trigger::TriggerObject
triggerSelector_
 

Detailed Description

Match reconstructed muons to HLT objects and plot efficiencies.

Note that this is not a true EDAnalyzer;

Documentation available on the CMS TWiki: https://twiki.cern.ch/twiki/bin/view/CMS/MuonHLTOfflinePerformance

Author
J. Slaunwhite, Jeff Klukas

Contanier class to handle vector of reconstructed muons matched to HLT objects used to plot efficiencies.

Note that this is not a true EDAnalyzer; rather, the intent is that one EDAnalyzer would call an instance of HLTMuonMatchAndPlotContainer

Documentation available on the CMS TWiki: https://twiki.cern.ch/twiki/bin/view/CMS/MuonHLTOfflinePerformance

Author
C. Battilana

Definition at line 60 of file HLTMuonMatchAndPlot.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 34 of file HLTMuonMatchAndPlot.cc.

References binParams_, cutMinPt_, fillMapFromPSet(), hltPath_, moduleLabels_, plotCuts_, and triggerLevel_.

36  :
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),
49  targetptCutZ_(targetParams_.getUntrackedParameter<double>("ptCut_Z",20.)),
50  targetptCutJpsi_(targetParams_.getUntrackedParameter<double>("ptCut_Jpsi",20.)),
51  probeMuonSelector_(probeParams_.getUntrackedParameter<string>("recoCuts", "")),
52  probeZ0Cut_(probeParams_.getUntrackedParameter<double>("z0Cut",0.)),
53  probeD0Cut_(probeParams_.getUntrackedParameter<double>("d0Cut",0.)),
56 {
57  // Create std::map<string, T> from ParameterSets.
58  fillMapFromPSet(binParams_, pset, "binParams");
59  fillMapFromPSet(plotCuts_, pset, "plotCuts");
60 
61  // Get the trigger level.
62  triggerLevel_ = "L3";
63  TPRegexp levelRegexp("L[1-3]");
64  size_t nModules = moduleLabels_.size();
65  TObjArray * levelArray = levelRegexp.MatchS(moduleLabels_[nModules - 1]);
66  if (levelArray->GetEntriesFast() > 0) {
67  triggerLevel_ = ((TObjString *)levelArray->At(0))->GetString();
68  }
69  delete levelArray;
70 
71  // Get the pT cut by parsing the name of the HLT path.
72  cutMinPt_ = 3;
73  TPRegexp ptRegexp("Mu([0-9]*)");
74  TObjArray * objArray = ptRegexp.MatchS(hltPath_);
75  if (objArray->GetEntriesFast() >= 2) {
76  TObjString * ptCutString = (TObjString *)objArray->At(1);
77  cutMinPt_ = atoi(ptCutString->GetString());
78  cutMinPt_ = ceil(cutMinPt_ * plotCuts_["minPtFactor"]);
79  }
80  delete objArray;
81 
82 }
T getUntrackedParameter(std::string const &, T const &) const
StringCutObjectSelector< reco::Muon > probeMuonSelector_
vector< string > vstring
Definition: ExoticaDQM.cc:8
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 fillMapFromPSet(std::map< std::string, T > &, const edm::ParameterSet &, std::string)
StringCutObjectSelector< reco::Muon > targetMuonSelector_
edm::ParameterSet probeParams_
StringCutObjectSelector< trigger::TriggerObject > triggerSelector_
std::map< std::string, double > plotCuts_

Member Function Documentation

void HLTMuonMatchAndPlot::analyze ( edm::Handle< reco::MuonCollection > &  allMuons,
edm::Handle< reco::BeamSpot > &  beamSpot,
edm::Handle< reco::VertexCollection > &  vertices,
edm::Handle< trigger::TriggerEvent > &  triggerSummary,
edm::Handle< edm::TriggerResults > &  triggerResults 
)

Definition at line 150 of file HLTMuonMatchAndPlot.cc.

References reco::LeafCandidate::charge(), cutMinPt_, deltaR(), reco::TrackBase::dxy(), reco::TrackBase::dz(), EFFICIENCY_SUFFIXES, eta, trigger::TriggerObject::eta(), reco::LeafCandidate::eta(), hasProbeRecoCuts, hasTargetRecoCuts, hasTriggerCuts_, hists_, HLT_FULL_cff::hltMuons, i, reco::Muon::innerTrack(), reco::Muon::isStandAloneMuon(), reco::Muon::isTrackerMuon(), j, relval_2017::k, ResonanceBuilder::mass, matchByDeltaR(), matches, metsig::muon, reco::Muon::outerTrack(), reco::LeafCandidate::p4(), phi, trigger::TriggerObject::phi(), reco::LeafCandidate::phi(), plotCuts_, probeD0Cut_, probeMuonSelector_, probeZ0Cut_, EnergyCorrector::pt, trigger::TriggerObject::pt(), reco::LeafCandidate::pt(), requiredTriggers_, selectedMuons(), selectedTriggerObjects(), findQualityFiles::size, createPayload::suffix, targetD0Cut_, targetMuonSelector_, targetptCutJpsi_, targetptCutZ_, targetZ0Cut_, triggerLevel_, and triggerSelector_.

155 {
156 
157  /*
158  if(gen != 0) {
159  for(g_part = gen->begin(); g_part != gen->end(); g_part++){
160  if( abs(g_part->pdgId()) == 13) {
161  if(!( g_part->status() ==1 || (g_part->status() ==2 && abs(g_part->pdgId())==5))) continue;
162  bool GenMomExists (true);
163  bool GenGrMomExists(true);
164  if( g_part->pt() < 10.0 ) continue;
165  if( fabs(g_part->eta()) > 2.4 ) continue;
166  int gen_id= g_part->pdgId();
167  const GenParticle* gen_lept = &(*g_part);
168  // get mother of gen_lept
169  const GenParticle* gen_mom = static_cast<const GenParticle*> (gen_lept->mother());
170  if(gen_mom==NULL) GenMomExists=false;
171  int m_id=-999;
172  if(GenMomExists) m_id = gen_mom->pdgId();
173  if(m_id != gen_id || !GenMomExists);
174  else{
175  int id= m_id;
176  while(id == gen_id && GenMomExists){
177  gen_mom = static_cast<const GenParticle*> (gen_mom->mother());
178  if(gen_mom==NULL){
179  GenMomExists=false;
180  }
181  if(GenMomExists) id=gen_mom->pdgId();
182  }
183  }
184  if(GenMomExists) m_id = gen_mom->pdgId();
185  gen_leptsp.push_back(gen_lept);
186  gen_momsp.push_back(gen_mom);
187  }
188  }
189  }
190 
191 
192  std::vector<GenParticle> gen_lepts;
193  for(size_t i = 0; i < gen_leptsp.size(); i++) {
194  gen_lepts.push_back(*gen_leptsp[i]);
195  }
196 
197 
198  */
199 
200 
201 
202  // Throw out this event if it doesn't pass the required triggers.
203  for (size_t i = 0; i < requiredTriggers_.size(); i++) {
204  unsigned int triggerIndex = triggerResults->find(requiredTriggers_[i]);
205  if (triggerIndex < triggerResults->size() ||
206  !triggerResults->accept(triggerIndex))
207  return;
208  }
209 
210 
211  // Select objects based on the configuration.
214  TriggerObjectCollection allTriggerObjects = triggerSummary->getObjects();
216  selectedTriggerObjects(allTriggerObjects, * triggerSummary, hasTriggerCuts_,triggerSelector_);
217 
218  // Fill plots for HLT muons.
219  for (size_t i = 0; i < hltMuons.size(); i++) {
220  hists_["hltPt"]->Fill(hltMuons[i].pt());
221  hists_["hltEta"]->Fill(hltMuons[i].eta());
222  hists_["hltPhi"]->Fill(hltMuons[i].phi());
223  }
224 
225  // Find the best trigger object matches for the targetMuons.
226  vector<size_t> matches = matchByDeltaR(targetMuons, hltMuons,
227  plotCuts_[triggerLevel_ + "DeltaR"]);
228 
229 
230  // Fill plots for matched muons.
231  bool pairalreadyconsidered = false;
232  for (size_t i = 0; i < targetMuons.size(); i++) {
233 
234  Muon & muon = targetMuons[i];
235 
236  // Fill plots which are not efficiencies.
237  if (matches[i] < targetMuons.size()) {
238  TriggerObject & hltMuon = hltMuons[matches[i]];
239  double ptRes = (muon.pt() - hltMuon.pt()) / muon.pt();
240  double etaRes = muon.eta() - hltMuon.eta();
241  double phiRes = muon.phi() - hltMuon.phi();
242  hists_["resolutionEta"]->Fill(etaRes);
243  hists_["resolutionPhi"]->Fill(phiRes);
244  hists_["resolutionPt"]->Fill(ptRes);
245  hists_["deltaR"]->Fill(deltaR(muon, hltMuon));
246  }
247 
248  // Fill numerators and denominators for efficiency plots.
249  for (size_t j = 0; j < 2; j++) {
250 
251  string suffix = EFFICIENCY_SUFFIXES[j];
252 
253  // If no match was found, then the numerator plots don't get filled.
254  if (suffix == "numer" && matches[i] >= targetMuons.size()) continue;
255 
256  if (muon.pt() > cutMinPt_) {
257  hists_["efficiencyEta_" + suffix]->Fill(muon.eta());
258  hists_["efficiencyPhiVsEta_" + suffix]->Fill(muon.eta(), muon.phi());
259  }
260 
261  if (fabs(muon.eta()) < plotCuts_["maxEta"]) {
262  hists_["efficiencyTurnOn_" + suffix]->Fill(muon.pt());
263  }
264 
265  if (muon.pt() > cutMinPt_ && fabs(muon.eta()) < plotCuts_["maxEta"]) {
266  const Track * track = 0;
267  if (muon.isTrackerMuon()) track = & * muon.innerTrack();
268  else if (muon.isStandAloneMuon()) track = & * muon.outerTrack();
269  if (track) {
270  double d0 = track->dxy(beamSpot->position());
271  double z0 = track->dz(beamSpot->position());
272  hists_["efficiencyVertex_" + suffix]->Fill(vertices->size());
273  hists_["efficiencyPhi_" + suffix]->Fill(muon.phi());
274  hists_["efficiencyD0_" + suffix]->Fill(d0);
275  hists_["efficiencyZ0_" + suffix]->Fill(z0);
276  hists_["efficiencyCharge_" + suffix]->Fill(muon.charge());
277  }
278  }
279  }
280 
281  // Fill plots for tag and probe
282  // Muon cannot be a tag because doesn't match an hlt muon
283  if(matches[i] >= targetMuons.size()) continue;
284  for (size_t k = 0; k < targetMuons.size(); k++) {
285  if(k == i) continue;
286  Muon & theProbe = targetMuons[k];
287  if (muon.charge() != theProbe.charge() && !pairalreadyconsidered) {
288  double mass = (muon.p4() + theProbe.p4()).M();
289  if(mass > 60 && mass < 120) {
290  if(muon.pt() < targetptCutZ_) continue;
291  hists_["massVsEtaZ_denom"]->Fill(theProbe.eta());
292  hists_["massVsPtZ_denom"]->Fill(theProbe.pt());
293  hists_["massVsVertexZ_denom"]->Fill(vertices->size());
294  if(matches[k] < targetMuons.size()) {
295  hists_["massVsEtaZ_numer"]->Fill(theProbe.eta());
296  hists_["massVsPtZ_numer"]->Fill(theProbe.pt());
297  hists_["massVsVertexZ_numer"]->Fill(vertices->size());
298  }
299  pairalreadyconsidered = true;
300  }
301  if(mass > 1 && mass < 4) {
302  if(muon.pt() < targetptCutJpsi_) continue;
303  hists_["massVsEtaJpsi_denom"]->Fill(theProbe.eta());
304  hists_["massVsPtJpsi_denom"]->Fill(theProbe.pt());
305  hists_["massVsVertexJpsi_denom"]->Fill(vertices->size());
306  if(matches[k] < targetMuons.size()) {
307  hists_["massVsEtaJpsi_numer"]->Fill(theProbe.eta());
308  hists_["massVsPtJpsi_numer"]->Fill(theProbe.pt());
309  hists_["massVsVertexJpsi_numer"]->Fill(vertices->size());
310  }
311  pairalreadyconsidered = true;
312  }
313  }
314  } // End loop over denominator and numerator.
315  } // End loop over targetMuons.
316 
317  // Plot fake rates (efficiency for HLT objects to not get matched to RECO).
318  vector<size_t> hltMatches = matchByDeltaR(hltMuons, targetMuons,
319  plotCuts_[triggerLevel_ + "DeltaR"]);
320  for (size_t i = 0; i < hltMuons.size(); i++) {
321  TriggerObject & hltMuon = hltMuons[i];
322  bool isFake = hltMatches[i] > hltMuons.size();
323  for (size_t j = 0; j < 2; j++) {
324  string suffix = EFFICIENCY_SUFFIXES[j];
325  // If match is found, then numerator plots should not get filled
326  if (suffix == "numer" && ! isFake) continue;
327  hists_["fakerateVertex_" + suffix]->Fill(vertices->size());
328  hists_["fakerateEta_" + suffix]->Fill(hltMuon.eta());
329  hists_["fakeratePhi_" + suffix]->Fill(hltMuon.phi());
330  hists_["fakerateTurnOn_" + suffix]->Fill(hltMuon.pt());
331  } // End loop over numerator and denominator.
332  } // End loop over hltMuons.
333 
334 
335 } // 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:48
float phi() const
Definition: TriggerObject.h:58
bool isTrackerMuon() const
Definition: Muon.h:223
virtual double phi() const final
momentum azimuthal angle
std::map< std::string, MonitorElement * > hists_
float eta() const
Definition: TriggerObject.h:57
bool isStandAloneMuon() const
Definition: Muon.h:224
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
trigger::TriggerObjectCollection selectedTriggerObjects(const trigger::TriggerObjectCollection &, const trigger::TriggerEvent &, bool hasTriggerCuts, const StringCutObjectSelector< trigger::TriggerObject > triggerSelector)
std::vector< size_t > matchByDeltaR(const std::vector< T1 > &, const std::vector< T2 > &, const double maxDeltaR=NOMATCH)
virtual int charge() const final
electric charge
Definition: LeafCandidate.h:91
int j
Definition: DBlmapReader.cc:9
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:604
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
StringCutObjectSelector< reco::Muon > targetMuonSelector_
StringCutObjectSelector< trigger::TriggerObject > triggerSelector_
virtual double eta() const final
momentum pseudorapidity
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:586
std::map< std::string, double > plotCuts_
tuple size
Write out results.
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
virtual double pt() const final
transverse momentum
void HLTMuonMatchAndPlot::beginRun ( DQMStore::IBooker iBooker,
const edm::Run iRun,
const edm::EventSetup iSetup 
)

Definition at line 84 of file HLTMuonMatchAndPlot.cc.

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

87 {
88 
89  TPRegexp suffixPtCut("Mu[0-9]+$");
90 
91  string baseDir = destination_;
92  if (baseDir[baseDir.size() - 1] != '/') baseDir += '/';
93  string pathSansSuffix = hltPath_;
94  if (hltPath_.rfind("_v") < hltPath_.length())
95  pathSansSuffix = hltPath_.substr(0, hltPath_.rfind("_v"));
96  iBooker.setCurrentFolder(baseDir + pathSansSuffix);
97 
98  // Form is book1D(name, binningType, title) where 'binningType' is used
99  // to fetch the bin settings from binParams_.
100  book1D(iBooker, "deltaR", "deltaR", ";#Deltar(reco, HLT);");
101  book1D(iBooker, "hltPt", "pt", ";p_{T} of HLT object");
102  book1D(iBooker, "hltEta", "eta", ";#eta of HLT object");
103  book1D(iBooker, "hltPhi", "phi", ";#phi of HLT object");
104  book1D(iBooker, "resolutionEta", "resolutionEta", ";#eta^{reco}-#eta^{HLT};");
105  book1D(iBooker, "resolutionPhi", "resolutionPhi", ";#phi^{reco}-#phi^{HLT};");
106  book1D(iBooker, "resolutionPt", "resolutionRel",
107  ";(p_{T}^{reco}-p_{T}^{HLT})/|p_{T}^{reco}|;");
108 
109  for (size_t i = 0; i < 2; i++) {
110 
111  string suffix = EFFICIENCY_SUFFIXES[i];
112 
113  book1D(iBooker, "efficiencyEta_" + suffix, "eta", ";#eta;");
114  book1D(iBooker, "efficiencyPhi_" + suffix, "phi", ";#phi;");
115  book1D(iBooker, "efficiencyTurnOn_" + suffix, "pt", ";p_{T};");
116  book1D(iBooker, "efficiencyD0_" + suffix, "d0", ";d0;");
117  book1D(iBooker, "efficiencyZ0_" + suffix, "z0", ";z0;");
118  book1D(iBooker, "efficiencyCharge_" + suffix, "charge", ";charge;");
119  book1D(iBooker, "efficiencyVertex_" + suffix, "NVertex", ";NVertex;");
120 
121  book2D(iBooker, "efficiencyPhiVsEta_" + suffix, "etaCoarse",
122  "phiCoarse", ";#eta;#phi");
123 
124  book1D(iBooker, "fakerateEta_" + suffix, "eta", ";#eta;");
125  book1D(iBooker, "fakerateVertex_" + suffix, "NVertex", ";NVertex;");
126  book1D(iBooker, "fakeratePhi_" + suffix, "phi", ";#phi;");
127  book1D(iBooker, "fakerateTurnOn_" + suffix, "pt", ";p_{T};");
128 
129  book1D(iBooker, "massVsEtaZ_" + suffix, "etaCoarse", ";#eta");
130  book1D(iBooker, "massVsEtaJpsi_" + suffix, "etaCoarse", ";#eta");
131  book1D(iBooker, "massVsPtZ_" + suffix, "ptCoarse", ";p_{T}");
132  book1D(iBooker, "massVsPtJpsi_" + suffix, "ptCoarse", ";p_{T}");
133  book1D(iBooker, "massVsVertexZ_" + suffix, "NVertex", ";NVertex");
134  book1D(iBooker, "massVsVertexJpsi_" + suffix, "NVertex", ";NVertex");
135 
136  }
137 
138 }
int i
Definition: DBlmapReader.cc:9
const std::string EFFICIENCY_SUFFIXES[2]
void book2D(DQMStore::IBooker &, std::string, std::string, std::string, std::string)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
void book1D(DQMStore::IBooker &, std::string, std::string, std::string)
void HLTMuonMatchAndPlot::book1D ( DQMStore::IBooker ,
std::string  ,
std::string  ,
std::string   
)
private

Definition at line 507 of file HLTMuonMatchAndPlot.cc.

References binParams_, DQMStore::IBooker::book1D(), fillEdges(), getTH1F(), hists_, and mergeVDriftHistosByStation::name.

Referenced by beginRun().

509 {
510 
511  /* Properly delete the array of floats that has been allocated on
512  * the heap by fillEdges. Avoid multiple copies and internal ROOT
513  * clones by simply creating the histograms directly in the DQMStore
514  * using the appropriate book1D method to handle the variable bins
515  * case. */
516 
517  size_t nBins;
518  float * edges = 0;
519  fillEdges(nBins, edges, binParams_[binningType]);
520 
521  hists_[name] = iBooker.book1D(name, title, nBins, edges);
522  if (hists_[name])
523  if (hists_[name]->getTH1F()->GetSumw2N())
524  hists_[name]->getTH1F()->Sumw2();
525 
526  if (edges)
527  delete [] edges;
528 
529 }
std::map< std::string, std::vector< double > > binParams_
std::map< std::string, MonitorElement * > hists_
TH1F * getTH1F(std::string name, std::string process, std::string rootfolder, DQMStore *dbe_, bool verb, bool clone)
void fillEdges(size_t &nBins, float *&edges, const std::vector< double > &binning)
void HLTMuonMatchAndPlot::book2D ( DQMStore::IBooker ,
std::string  ,
std::string  ,
std::string  ,
std::string   
)
private

Definition at line 534 of file HLTMuonMatchAndPlot.cc.

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

Referenced by beginRun().

537 {
538 
539  /* Properly delete the arrays of floats that have been allocated on
540  * the heap by fillEdges. Avoid multiple copies and internal ROOT
541  * clones by simply creating the histograms directly in the DQMStore
542  * using the appropriate book2D method to handle the variable bins
543  * case. */
544 
545  size_t nBinsX;
546  float * edgesX = 0;
547  fillEdges(nBinsX, edgesX, binParams_[binningTypeX]);
548 
549  size_t nBinsY;
550  float * edgesY = 0;
551  fillEdges(nBinsY, edgesY, binParams_[binningTypeY]);
552 
553  hists_[name] = iBooker.book2D(name.c_str(), title.c_str(),
554  nBinsX, edgesX, nBinsY, edgesY);
555  if (hists_[name])
556  if (hists_[name]->getTH2F()->GetSumw2N())
557  hists_[name]->getTH2F()->Sumw2();
558 
559  if (edgesX)
560  delete [] edgesX;
561  if (edgesY)
562  delete [] edgesY;
563 
564 }
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)
void fillEdges(size_t &nBins, float *&edges, const std::vector< double > &binning)
void HLTMuonMatchAndPlot::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)

Definition at line 142 of file HLTMuonMatchAndPlot.cc.

144 {
145 
146 }
void HLTMuonMatchAndPlot::fillEdges ( size_t &  nBins,
float *&  edges,
const std::vector< double > &  binning 
)

Definition at line 341 of file HLTMuonMatchAndPlot.cc.

References i, and min().

Referenced by book1D(), and book2D().

343 {
344 
345  if (binning.size() < 3) {
346  LogWarning("HLTMuonVal") << "Invalid binning parameters!";
347  return;
348  }
349 
350  // Fixed-width binning.
351  if (binning.size() == 3) {
352  nBins = binning[0];
353  edges = new float[nBins + 1];
354  const double min = binning[1];
355  const double binwidth = (binning[2] - binning[1]) / nBins;
356  for (size_t i = 0; i <= nBins; i++) edges[i] = min + (binwidth * i);
357  }
358 
359  // Variable-width binning.
360  else {
361  nBins = binning.size() - 1;
362  edges = new float[nBins + 1];
363  for (size_t i = 0; i <= nBins; i++) edges[i] = binning[i];
364  }
365 
366 }
int i
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58
template<class T >
void HLTMuonMatchAndPlot::fillMapFromPSet ( std::map< std::string, T > &  ,
const edm::ParameterSet ,
std::string   
)

Referenced by HLTMuonMatchAndPlot().

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

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

377 {
378 
379  // Get the ParameterSet with name 'target' from 'pset'
380  ParameterSet targetPset;
381  if (pset.existsAs<ParameterSet>(target, true)) // target is tracked
382  targetPset = pset.getParameterSet(target);
383  else if (pset.existsAs<ParameterSet>(target, false)) // target is untracked
384  targetPset = pset.getUntrackedParameterSet(target);
385 
386  // Get the parameter names from targetPset
387  vector<string> names = targetPset.getParameterNames();
388  vector<string>::const_iterator iter;
389 
390  for (iter = names.begin(); iter != names.end(); ++iter) {
391  if (targetPset.existsAs<T>(* iter, true)) // target is tracked
392  m[* iter] = targetPset.getParameter<T>(* iter);
393  else if (targetPset.existsAs<T>(* iter, false)) // target is untracked
394  m[* iter] = targetPset.getUntrackedParameter<T>(* iter);
395  }
396 
397 }
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:186
static const HistoName names[]
ParameterSet 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 404 of file HLTMuonMatchAndPlot.cc.

References deltaR(), i, j, relval_2017::k, HLT_25ns10e33_v2_cff::minDeltaR, NOMATCH, and mps_fire::result.

407 {
408 
409  const size_t n1 = collection1.size();
410  const size_t n2 = collection2.size();
411 
412  vector<size_t> result(n1, -1);
413  vector<vector<double> > deltaRMatrix(n1, vector<double>(n2, NOMATCH));
414 
415  for (size_t i = 0; i < n1; i++)
416  for (size_t j = 0; j < n2; j++) {
417  deltaRMatrix[i][j] = deltaR(collection1[i], collection2[j]);
418  }
419 
420  // Run through the matrix n1 times to make sure we've found all matches.
421  for (size_t k = 0; k < n1; k++) {
422  size_t i_min = -1;
423  size_t j_min = -1;
424  double minDeltaR = maxDeltaR;
425  // find the smallest deltaR
426  for (size_t i = 0; i < n1; i++)
427  for (size_t j = 0; j < n2; j++)
428  if (deltaRMatrix[i][j] < minDeltaR) {
429  i_min = i;
430  j_min = j;
431  minDeltaR = deltaRMatrix[i][j];
432  }
433  // If a match has been made, save it and make those candidates unavailable.
434  if (minDeltaR < maxDeltaR) {
435  result[i_min] = j_min;
436  deltaRMatrix[i_min] = vector<double>(n2, NOMATCH);
437  for (size_t i = 0; i < n1; i++)
438  deltaRMatrix[i][j_min] = NOMATCH;
439  }
440  }
441 
442  return result;
443 
444 }
int i
Definition: DBlmapReader.cc:9
tuple result
Definition: mps_fire.py:84
const double NOMATCH
int j
Definition: DBlmapReader.cc:9
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 449 of file HLTMuonMatchAndPlot.cc.

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

Referenced by analyze().

454 {
455 
456  // If there is no selector (recoCuts does not exists), return an empty collection.
457  if (!hasRecoCuts)
458  return MuonCollection();
459 
460  MuonCollection reducedMuons;
461  for (auto const& mu : allMuons){
462  const Track * track = 0;
463  if (mu.isTrackerMuon()) track = & * mu.innerTrack();
464  else if (mu.isStandAloneMuon()) track = & * mu.outerTrack();
465  if (track && selector(mu) &&
466  fabs(track->dxy(beamSpot.position())) < d0Cut &&
467  fabs(track->dz(beamSpot.position())) < z0Cut)
468  reducedMuons.push_back(mu);
469  }
470 
471  return reducedMuons;
472 
473 }
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
tuple allMuons
Definition: allMuons_cfi.py:3
const int mu
Definition: Constants.h:22
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:604
const Point & position() const
position
Definition: BeamSpot.h:62
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:586
TriggerObjectCollection HLTMuonMatchAndPlot::selectedTriggerObjects ( const trigger::TriggerObjectCollection triggerObjects,
const trigger::TriggerEvent triggerSummary,
bool  hasTriggerCuts,
const StringCutObjectSelector< trigger::TriggerObject triggerSelector 
)
private

Definition at line 478 of file HLTMuonMatchAndPlot.cc.

References trigger::TriggerEvent::filterIndex(), trigger::TriggerEvent::filterKeys(), hltProcessName_, j, relativeConstraints::keys, moduleLabels_, and trigger::TriggerEvent::sizeFilters().

Referenced by analyze().

483 {
484  if ( !hasTriggerCuts) return TriggerObjectCollection();
485 
486  InputTag filterTag(moduleLabels_[moduleLabels_.size() - 1], "",
488  size_t filterIndex = triggerSummary.filterIndex(filterTag);
489 
490  TriggerObjectCollection selectedObjects;
491 
492  if (filterIndex < triggerSummary.sizeFilters()) {
493  const Keys &keys = triggerSummary.filterKeys(filterIndex);
494  for (size_t j = 0; j < keys.size(); j++ ){
495  TriggerObject foundObject = triggerObjects[keys[j]];
496  if (triggerSelector(foundObject))
497  selectedObjects.push_back(foundObject);
498  }
499  }
500 
501  return selectedObjects;
502 
503 }
trigger::size_type sizeFilters() const
Definition: TriggerEvent.h:135
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
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:22
int j
Definition: DBlmapReader.cc:9
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:81
std::vector< size_type > Keys

Member Data Documentation

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

Definition at line 106 of file HLTMuonMatchAndPlot.h.

Referenced by book1D(), book2D(), and HLTMuonMatchAndPlot().

unsigned int HLTMuonMatchAndPlot::cutMinPt_
private

Definition at line 113 of file HLTMuonMatchAndPlot.h.

Referenced by analyze(), and HLTMuonMatchAndPlot().

std::string HLTMuonMatchAndPlot::destination_
private

Definition at line 104 of file HLTMuonMatchAndPlot.h.

Referenced by beginRun().

bool HLTMuonMatchAndPlot::hasProbeRecoCuts
private

Definition at line 120 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

bool HLTMuonMatchAndPlot::hasTargetRecoCuts
private

Definition at line 119 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

bool HLTMuonMatchAndPlot::hasTriggerCuts_
private

Definition at line 132 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

std::map<std::string, MonitorElement *> HLTMuonMatchAndPlot::hists_
private

Definition at line 116 of file HLTMuonMatchAndPlot.h.

Referenced by analyze(), book1D(), and book2D().

std::string HLTMuonMatchAndPlot::hltPath_
private

Definition at line 114 of file HLTMuonMatchAndPlot.h.

Referenced by beginRun(), and HLTMuonMatchAndPlot().

std::string HLTMuonMatchAndPlot::hltProcessName_
private

Definition at line 103 of file HLTMuonMatchAndPlot.h.

Referenced by selectedTriggerObjects().

std::vector<std::string> HLTMuonMatchAndPlot::moduleLabels_
private

Definition at line 115 of file HLTMuonMatchAndPlot.h.

Referenced by HLTMuonMatchAndPlot(), and selectedTriggerObjects().

std::map<std::string, double> HLTMuonMatchAndPlot::plotCuts_
private

Definition at line 107 of file HLTMuonMatchAndPlot.h.

Referenced by analyze(), and HLTMuonMatchAndPlot().

double HLTMuonMatchAndPlot::probeD0Cut_
private

Definition at line 129 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

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

Definition at line 127 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

edm::ParameterSet HLTMuonMatchAndPlot::probeParams_
private

Definition at line 109 of file HLTMuonMatchAndPlot.h.

double HLTMuonMatchAndPlot::probeZ0Cut_
private

Definition at line 128 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

std::vector<std::string> HLTMuonMatchAndPlot::requiredTriggers_
private

Definition at line 105 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

double HLTMuonMatchAndPlot::targetD0Cut_
private

Definition at line 124 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

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

Definition at line 122 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

edm::ParameterSet HLTMuonMatchAndPlot::targetParams_
private

Definition at line 108 of file HLTMuonMatchAndPlot.h.

double HLTMuonMatchAndPlot::targetptCutJpsi_
private

Definition at line 126 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

double HLTMuonMatchAndPlot::targetptCutZ_
private

Definition at line 125 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

double HLTMuonMatchAndPlot::targetZ0Cut_
private

Definition at line 123 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().

std::string HLTMuonMatchAndPlot::triggerLevel_
private

Definition at line 112 of file HLTMuonMatchAndPlot.h.

Referenced by analyze(), and HLTMuonMatchAndPlot().

StringCutObjectSelector<trigger::TriggerObject> HLTMuonMatchAndPlot::triggerSelector_
private

Definition at line 131 of file HLTMuonMatchAndPlot.h.

Referenced by analyze().