CMS 3D CMS Logo

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