CMS 3D CMS Logo

HLTMuonPlotter.cc
Go to the documentation of this file.
1 
14 
19 
20 using namespace std;
21 using namespace edm;
22 using namespace reco;
23 using namespace trigger;
24 
25 typedef vector<ParameterSet> Parameters;
26 
28  string hltPath,
29  const std::vector<string> &moduleLabels,
30  const std::vector<string> &stepLabels,
34  :
35 
36  l1Matcher_(pset) {
37  hltPath_ = hltPath;
38  moduleLabels_ = moduleLabels;
39  stepLabels_ = stepLabels;
40  hltProcessName_ = pset.getParameter<string>("hltProcessName");
41 
42  cutsDr_ = pset.getParameter<vector<double>>("cutsDr");
43 
44  parametersEta_ = pset.getParameter<vector<double>>("parametersEta");
45  parametersPhi_ = pset.getParameter<vector<double>>("parametersPhi");
46  parametersTurnOn_ = pset.getParameter<vector<double>>("parametersTurnOn");
47 
48  genMuonCut_ = pset.getParameter<string>("genMuonCut");
49  recMuonCut_ = pset.getParameter<string>("recMuonCut");
50 
51  genMuonSelector_ = nullptr;
52  recMuonSelector_ = nullptr;
53 
54  // set tokens
55  hltTriggerSummaryRAW_ = tokens.get<0>();
56  genParticleLabel_ = tokens.get<1>();
57  recMuonLabel_ = tokens.get<2>();
58 }
59 
61 
62 void HLTMuonPlotter::beginRun(DQMStore::IBooker &iBooker, const Run &iRun, const EventSetup &iSetup) {
63  static int runNumber = 0;
64  runNumber++;
65 
66  l1Matcher_.init(iSetup);
67 
68  cutMaxEta_ = 2.4;
69  if (hltPath_.find("eta2p1") != string::npos)
70  cutMaxEta_ = 2.1;
71 
72  // Choose a pT cut for gen/rec muons based on the pT cut in the hltPath_
73  unsigned int threshold = 0;
74  TPRegexp ptRegexp("Mu([0-9]+)");
75  TObjArray *regexArray = ptRegexp.MatchS(hltPath_);
76  if (regexArray->GetEntriesFast() == 2) {
77  threshold = atoi(((TObjString *)regexArray->At(1))->GetString());
78  }
79  delete regexArray;
80  // We select a whole number min pT cut slightly above the hltPath_'s final
81  // pt threshold, then subtract a bit to let through particle gun muons with
82  // exact integer pT:
83  cutMinPt_ = ceil(threshold * 1.1) - 0.01;
84  if (cutMinPt_ < 0.)
85  cutMinPt_ = 0.;
86 
87  string baseDir = "HLT/Muon/Distributions/";
88  iBooker.setCurrentFolder(baseDir + hltPath_);
89 
90  vector<string> sources(2);
91  sources[0] = "gen";
92  sources[1] = "rec";
93 
94  elements_["CutMinPt"] = iBooker.bookFloat("CutMinPt");
95  elements_["CutMaxEta"] = iBooker.bookFloat("CutMaxEta");
96  elements_["CutMinPt"]->Fill(cutMinPt_);
97  elements_["CutMaxEta"]->Fill(cutMaxEta_);
98 
99  for (size_t i = 0; i < sources.size(); i++) {
100  string source = sources[i];
101  for (size_t j = 0; j < stepLabels_.size(); j++) {
102  bookHist(iBooker, hltPath_, stepLabels_[j], source, "Eta");
103  bookHist(iBooker, hltPath_, stepLabels_[j], source, "Phi");
104  bookHist(iBooker, hltPath_, stepLabels_[j], source, "MaxPt1");
105  bookHist(iBooker, hltPath_, stepLabels_[j], source, "MaxPt2");
106  }
107  }
108 }
109 
110 void HLTMuonPlotter::analyze(const Event &iEvent, const EventSetup &iSetup) {
111  static int eventNumber = 0;
112  eventNumber++;
113  LogTrace("HLTMuonVal") << "In HLTMuonPlotter::analyze, "
114  << "Event: " << eventNumber;
115 
116  // cout << hltPath_ << endl;
117  // for (size_t i = 0; i < moduleLabels_.size(); i++)
118  // cout << " " << moduleLabels_[i] << endl;
119 
120  Handle<TriggerEventWithRefs> rawTriggerEvent;
121  Handle<MuonCollection> recMuons;
123 
124  iEvent.getByToken(hltTriggerSummaryRAW_, rawTriggerEvent);
125  if (rawTriggerEvent.failedToGet()) {
126  LogError("HLTMuonVal") << "No trigger summary found";
127  return;
128  }
129  iEvent.getByToken(recMuonLabel_, recMuons);
130  iEvent.getByToken(genParticleLabel_, genParticles);
131 
132  vector<string> sources;
133  if (genParticles.isValid())
134  sources.push_back("gen");
135  if (recMuons.isValid())
136  sources.push_back("rec");
137 
138  for (size_t sourceNo = 0; sourceNo < sources.size(); sourceNo++) {
139  string source = sources[sourceNo];
140 
141  // If this is the first event, initialize selectors
142  if (!genMuonSelector_)
144  if (!recMuonSelector_)
146 
147  // Make each good gen/rec muon into the base cand for a MatchStruct
148  vector<MatchStruct> matches;
149  if (source == "gen" && genParticles.isValid())
150  for (size_t i = 0; i < genParticles->size(); i++)
151  if ((*genMuonSelector_)(genParticles->at(i)))
152  matches.push_back(MatchStruct(&genParticles->at(i)));
153  if (source == "rec" && recMuons.isValid())
154  for (size_t i = 0; i < recMuons->size(); i++)
155  if ((*recMuonSelector_)(recMuons->at(i)))
156  matches.push_back(MatchStruct(&recMuons->at(i)));
157 
158  // Sort the MatchStructs by pT for later filling of turn-on curve
159  sort(matches.begin(), matches.end(), matchesByDescendingPt());
160 
161  const bool isDoubleMuonPath = (hltPath_.find("Double") != string::npos);
162  const size_t nFilters = moduleLabels_.size();
163  const size_t nSteps = stepLabels_.size();
164  const size_t nStepsHlt = nSteps - 2;
165  const int nObjectsToPassPath = (isDoubleMuonPath) ? 2 : 1;
166  l1t::MuonVectorRef candsL1;
167  vector<vector<RecoChargedCandidateRef>> refsHlt(nStepsHlt);
168  vector<vector<const RecoChargedCandidate *>> candsHlt(nStepsHlt);
169 
170  for (size_t i = 0; i < nFilters; i++) {
171  const int hltStep = i - 1;
173  size_t iFilter = rawTriggerEvent->filterIndex(tag);
174  if (iFilter < rawTriggerEvent->size()) {
175  if (i == 0)
176  rawTriggerEvent->getObjects(iFilter, TriggerL1Mu, candsL1);
177  else
178  rawTriggerEvent->getObjects(iFilter, TriggerMuon, refsHlt[hltStep]);
179  } else
180  LogTrace("HLTMuonVal") << "No collection with label " << tag;
181  }
182  for (size_t i = 0; i < nStepsHlt; i++)
183  for (size_t j = 0; j < refsHlt[i].size(); j++)
184  if (refsHlt[i][j].isAvailable()) {
185  candsHlt[i].push_back(&*refsHlt[i][j]);
186  } else {
187  LogWarning("HLTMuonPlotter") << "Ref refsHlt[i][j]: product not available " << i << " " << j;
188  }
189 
190  // Add trigger objects to the MatchStructs
191  findMatches(matches, candsL1, candsHlt);
192 
193  vector<size_t> matchesInEtaRange;
194  vector<bool> hasMatch(matches.size(), true);
195 
196  for (size_t step = 0; step < nSteps; step++) {
197  size_t hltStep = (step >= 2) ? step - 2 : 0;
198  if (nSteps == 6)
199  hltStep = hltStep - 1; // case of the tracker muon (it has no L2)
200  size_t level = 0;
201  if ((stepLabels_[step].find("L3TkIso") != string::npos) || (stepLabels_[step].find("TkTkIso") != string::npos))
202  level = 6;
203  else if ((stepLabels_[step].find("L3HcalIso") != string::npos) ||
204  (stepLabels_[step].find("TkEcalIso") != string::npos))
205  level = 5;
206  else if ((stepLabels_[step].find("L3EcalIso") != string::npos) ||
207  (stepLabels_[step].find("TkEcalIso") != string::npos))
208  level = 4;
209  else if ((stepLabels_[step].find("L3") != string::npos) || (stepLabels_[step].find("Tk") != string::npos))
210  level = 3;
211  else if (stepLabels_[step].find("L2") != string::npos)
212  level = 2;
213  else if (stepLabels_[step].find("L1") != string::npos)
214  level = 1;
215 
216  for (size_t j = 0; j < matches.size(); j++) {
217  if (level == 0) {
218  if (fabs(matches[j].candBase->eta()) < cutMaxEta_)
219  matchesInEtaRange.push_back(j);
220  } else if (level == 1) {
221  if (matches[j].candL1 == nullptr)
222  hasMatch[j] = false;
223  } else if (level >= 2) {
224  if (matches[j].candHlt[hltStep] == nullptr)
225  hasMatch[j] = false;
226  else if (!hasMatch[j]) {
227  LogTrace("HLTMuonVal") << "Match found for HLT step " << hltStep << " of " << nStepsHlt
228  << " without previous match!";
229  break;
230  }
231  }
232  }
233 
234  if (std::count(hasMatch.begin(), hasMatch.end(), true) < nObjectsToPassPath)
235  break;
236 
237  string pre = source + "Pass";
238  string post = "_" + stepLabels_[step];
239 
240  for (size_t j = 0; j < matches.size(); j++) {
241  float pt = matches[j].candBase->pt();
242  float eta = matches[j].candBase->eta();
243  float phi = matches[j].candBase->phi();
244  if (hasMatch[j]) {
245  if (!matchesInEtaRange.empty() && j == matchesInEtaRange[0])
246  elements_[pre + "MaxPt1" + post]->Fill(pt);
247  if (matchesInEtaRange.size() >= 2 && j == matchesInEtaRange[1])
248  elements_[pre + "MaxPt2" + post]->Fill(pt);
249  if (pt > cutMinPt_) {
250  elements_[pre + "Eta" + post]->Fill(eta);
251  if (fabs(eta) < cutMaxEta_)
252  elements_[pre + "Phi" + post]->Fill(phi);
253  }
254  }
255  }
256  }
257 
258  } // End loop over sources
259 }
260 
261 boost::tuple<edm::EDGetTokenT<trigger::TriggerEventWithRefs>,
266  iC.consumes<TriggerEventWithRefs>(edm::InputTag("hltTriggerSummaryRAW"));
268  iC.consumes<GenParticleCollection>(pset.getParameter<string>("genParticleLabel"));
270  iC.consumes<MuonCollection>(pset.getParameter<string>("recMuonLabel"));
271 
272  boost::tuple<edm::EDGetTokenT<trigger::TriggerEventWithRefs>,
275  myTuple(_hltTriggerSummaryRAW, _genParticleLabel, _recMuonLabel);
276 
277  return (myTuple);
278 }
279 
280 void HLTMuonPlotter::findMatches(vector<MatchStruct> &matches,
281  const l1t::MuonVectorRef &candsL1,
282  const std::vector<vector<const RecoChargedCandidate *>> &candsHlt) {
283  set<size_t>::iterator it;
284 
285  set<size_t> indicesL1;
286  for (size_t i = 0; i < candsL1.size(); i++)
287  indicesL1.insert(i);
288 
289  vector<set<size_t>> indicesHlt(candsHlt.size());
290  for (size_t i = 0; i < candsHlt.size(); i++)
291  for (size_t j = 0; j < candsHlt[i].size(); j++)
292  indicesHlt[i].insert(j);
293 
294  for (size_t i = 0; i < matches.size(); i++) {
295  const Candidate *cand = matches[i].candBase;
296 
297  double bestDeltaR = cutsDr_[0];
298  size_t bestMatch = kNull;
299  for (it = indicesL1.begin(); it != indicesL1.end(); it++) {
300  if (candsL1[*it].isAvailable()) {
301  double dR = deltaR(cand->eta(), cand->phi(), candsL1[*it]->eta(), candsL1[*it]->phi());
302  if (dR < bestDeltaR) {
303  bestMatch = *it;
304  bestDeltaR = dR;
305  }
306  // TrajectoryStateOnSurface propagated;
307  // float dR = 999., dPhi = 999.;
308  // bool isValid = l1Matcher_.match(* cand, * candsL1[*it],
309  // dR, dPhi, propagated);
310  // if (isValid && dR < bestDeltaR) {
311  // bestMatch = *it;
312  // bestDeltaR = dR;
313  // }
314  } else {
315  LogWarning("HLTMuonPlotter") << "Ref candsL1[*it]: product not available " << *it;
316  }
317  }
318 
319  if (bestMatch != kNull)
320  matches[i].candL1 = &*candsL1[bestMatch];
321  indicesL1.erase(bestMatch);
322 
323  matches[i].candHlt.assign(candsHlt.size(), nullptr);
324  for (size_t j = 0; j < candsHlt.size(); j++) {
325  size_t level = (candsHlt.size() == 4) ? (j < 2) ? 2 : 3 : (candsHlt.size() == 2) ? (j < 1) ? 2 : 3 : 2;
326  bestDeltaR = cutsDr_[level - 2];
327  bestMatch = kNull;
328  for (it = indicesHlt[j].begin(); it != indicesHlt[j].end(); it++) {
329  double dR = deltaR(cand->eta(), cand->phi(), candsHlt[j][*it]->eta(), candsHlt[j][*it]->phi());
330  if (dR < bestDeltaR) {
331  bestMatch = *it;
332  bestDeltaR = dR;
333  }
334  }
335  if (bestMatch != kNull)
336  matches[i].candHlt[j] = candsHlt[j][bestMatch];
337  indicesHlt[j].erase(bestMatch);
338  }
339 
340  // cout << " Muon: " << cand->eta() << ", ";
341  // if (matches[i].candL1) cout << matches[i].candL1->eta() << ", ";
342  // else cout << "none, ";
343  // for (size_t j = 0; j < candsHlt.size(); j++)
344  // if (matches[i].candHlt[j]) cout << matches[i].candHlt[j]->eta() <<
345  // ", "; else cout << "none, ";
346  // cout << endl;
347  }
348 }
349 
350 void HLTMuonPlotter::bookHist(DQMStore::IBooker &iBooker, string path, string label, string source, string type) {
351  string sourceUpper = source;
352  sourceUpper[0] = toupper(sourceUpper[0]);
353  string name = source + "Pass" + type + "_" + label;
354  TH1F *h;
355 
356  if (type.find("MaxPt") != string::npos) {
357  string desc = (type == "MaxPt1") ? "Leading" : "Next-to-Leading";
358  string title = "pT of " + desc + " " + sourceUpper + " Muon " + "matched to " + label;
359  const size_t nBins = parametersTurnOn_.size() - 1;
360  float *edges = new float[nBins + 1];
361  for (size_t i = 0; i < nBins + 1; i++)
362  edges[i] = parametersTurnOn_[i];
363  h = new TH1F(name.c_str(), title.c_str(), nBins, edges);
364  }
365 
366  else {
367  string symbol = (type == "Eta") ? "#eta" : "#phi";
368  string title = symbol + " of " + sourceUpper + " Muons " + "matched to " + label;
369  vector<double> params = (type == "Eta") ? parametersEta_ : parametersPhi_;
370  int nBins = (int)params[0];
371  double min = params[1];
372  double max = params[2];
373  h = new TH1F(name.c_str(), title.c_str(), nBins, min, max);
374  }
375 
376  h->Sumw2();
377  elements_[name] = iBooker.book1D(name, h);
378  delete h;
379 }
size
Write out results.
HLTMuonPlotter(const edm::ParameterSet &, std::string, const std::vector< std::string > &, const std::vector< std::string > &, const boost::tuple< edm::EDGetTokenT< trigger::TriggerEventWithRefs >, edm::EDGetTokenT< reco::GenParticleCollection >, edm::EDGetTokenT< reco::MuonCollection >> &)
type
Definition: HCALResponse.h:21
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
void init(const edm::EventSetup &iSetup)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
enum start value shifted to 81 so as to avoid clashes with PDG codes
void beginRun(DQMStore::IBooker &, const edm::Run &, const edm::EventSetup &)
std::string hltPath_
std::vector< double > parametersEta_
L1MuonMatcherAlgo l1Matcher_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::string genMuonCut_
std::vector< std::string > moduleLabels_
std::vector< MuonRef > MuonVectorRef
Definition: Muon.h:15
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
static boost::tuple< edm::EDGetTokenT< trigger::TriggerEventWithRefs >, edm::EDGetTokenT< reco::GenParticleCollection >, edm::EDGetTokenT< reco::MuonCollection > > getTokens(const edm::ParameterSet &, edm::ConsumesCollector &&)
char const * label
std::string recMuonCut_
int iEvent
Definition: GenABIO.cc:224
vector< ParameterSet > Parameters
void findMatches(std::vector< MatchStruct > &, const l1t::MuonVectorRef &candsL1, const std::vector< std::vector< const reco::RecoChargedCandidate * >> &)
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::vector< double > cutsDr_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
const unsigned int kNull
edm::EDGetTokenT< reco::MuonCollection > recMuonLabel_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
void analyze(const edm::Event &, const edm::EventSetup &)
T min(T a, T b)
Definition: MathUtil.h:58
bool insert(Storage &iStorage, ItemType *iItem, const IdTag &iIdTag)
Definition: HCMethods.h:50
bool isValid() const
Definition: HandleBase.h:74
StringCutObjectSelector< reco::GenParticle > * genMuonSelector_
#define LogTrace(id)
size_type filterIndex(const edm::InputTag &filterTag) const
index from tag
std::vector< std::string > stepLabels_
bool failedToGet() const
Definition: HandleBase.h:78
virtual double eta() const =0
momentum pseudorapidity
void bookHist(DQMStore::IBooker &, std::string, std::string, std::string, std::string)
def bestMatch(object, matchCollection)
Definition: deltar.py:138
void getObjects(size_type filter, Vids &ids, VRphoton &photons) const
extract Ref<C>s for a specific filter and of specific physics type
StringCutObjectSelector< reco::Muon > * recMuonSelector_
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
std::vector< double > parametersPhi_
step
Definition: StallMonitor.cc:94
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:105
std::map< std::string, MonitorElement * > elements_
std::string hltProcessName_
edm::EDGetTokenT< reco::GenParticleCollection > genParticleLabel_
edm::EDGetTokenT< trigger::TriggerEventWithRefs > hltTriggerSummaryRAW_
std::vector< double > parametersTurnOn_
virtual double phi() const =0
momentum azimuthal angle
static std::string const source
Definition: EdmProvDump.cc:47
Definition: Run.h:45