CMS 3D CMS Logo

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