CMS 3D CMS Logo

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