27 using namespace trigger;
28 using namespace l1extra;
38 const std::vector<string>& moduleLabels,
39 const std::vector<string>& stepLabels,
93 if (
hltPath_.find(
"eta2p1") != string::npos)
98 TPRegexp ptRegexp(
"Mu([0-9]+)");
99 TObjArray * regexArray = ptRegexp.MatchS(
hltPath_);
100 if (regexArray->GetEntriesFast() == 2) {
101 threshold = atoi(((TObjString *)regexArray->At(1))->GetString());
107 cutMinPt_ = ceil(threshold * 1.1) - 0.01;
110 string baseDir =
"HLT/Muon/Distributions/";
113 vector<string> sources(2);
122 for (
size_t i = 0;
i < sources.size();
i++) {
140 static int eventNumber = 0;
142 LogTrace(
"HLTMuonVal") <<
"In HLTMuonPlotter::analyze, "
143 <<
"Event: " << eventNumber;
155 {
LogError(
"HLTMuonVal") <<
"No trigger summary found";
return;}
159 vector<string> sources;
160 if (genParticles.
isValid()) sources.push_back(
"gen");
161 if ( recMuons.
isValid()) sources.push_back(
"rec");
163 for (
size_t sourceNo = 0; sourceNo < sources.size(); sourceNo++) {
165 string source = sources[sourceNo];
175 if (source ==
"gen" && genParticles.
isValid())
176 for (
size_t i = 0;
i < genParticles->size();
i++)
179 if (source ==
"rec" && recMuons.
isValid())
180 for (
size_t i = 0;
i < recMuons->size();
i++)
187 const bool isDoubleMuonPath = (
hltPath_.find(
"Double") != string::npos);
190 const size_t nStepsHlt = nSteps - 2;
191 const int nObjectsToPassPath = (isDoubleMuonPath) ? 2 : 1;
192 vector< L1MuonParticleRef > candsL1;
193 vector< vector< RecoChargedCandidateRef > > refsHlt(nStepsHlt);
194 vector< vector< const RecoChargedCandidate * > > candsHlt(nStepsHlt);
196 for (
size_t i = 0;
i < nFilters;
i++) {
197 const int hltStep =
i - 1;
199 size_t iFilter = rawTriggerEvent->filterIndex(tag);
200 if (iFilter < rawTriggerEvent->
size()) {
202 rawTriggerEvent->getObjects(iFilter,
TriggerL1Mu, candsL1);
207 else LogTrace(
"HLTMuonVal") <<
"No collection with label " <<
tag;
209 for (
size_t i = 0;
i < nStepsHlt;
i++)
210 for (
size_t j = 0;
j < refsHlt[
i].size();
j++)
211 if (refsHlt[
i][
j].isAvailable()) {
212 candsHlt[
i].push_back(& * refsHlt[
i][
j]);
215 <<
"Ref refsHlt[i][j]: product not available "
222 vector<size_t> matchesInEtaRange;
223 vector<bool> hasMatch(matches.size(),
true);
227 size_t hltStep = (
step >= 2) ?
step - 2 : 0;
228 if (nSteps == 6) hltStep=hltStep-1;
237 for (
size_t j = 0;
j < matches.size();
j++) {
240 matchesInEtaRange.push_back(
j);
242 else if (level == 1) {
243 if (matches[
j].candL1 == 0)
246 else if (level >= 2) {
247 if (matches[
j].candHlt[hltStep] == 0)
249 else if (!hasMatch[
j]) {
250 LogTrace(
"HLTMuonVal") <<
"Match found for HLT step " << hltStep
251 <<
" of " << nStepsHlt
252 <<
" without previous match!";
258 if (
std::count(hasMatch.begin(), hasMatch.end(),
true) <
262 string pre = source +
"Pass";
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();
270 if (matchesInEtaRange.size() >= 1 && j == matchesInEtaRange[0])
272 if (matchesInEtaRange.size() >= 2 && j == matchesInEtaRange[1])
275 elements_[pre +
"Eta" + post]->Fill(eta);
277 elements_[pre +
"Phi" + post]->Fill(phi);
314 vector<MatchStruct> & matches,
315 const std::vector<L1MuonParticleRef>& candsL1,
316 const std::vector< vector< const RecoChargedCandidate *> >& candsHlt)
319 set<size_t>::iterator it;
321 set<size_t> indicesL1;
322 for (
size_t i = 0;
i < candsL1.size();
i++)
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++)
330 for (
size_t i = 0;
i < matches.size();
i++) {
334 double bestDeltaR =
cutsDr_[0];
336 for (it = indicesL1.begin(); it != indicesL1.end(); it++) {
337 if (candsL1[*it].isAvailable()) {
339 candsL1[*it]->eta(), candsL1[*it]->phi());
340 if (dR < bestDeltaR) {
354 <<
"Ref candsL1[*it]: product not available "
358 if (bestMatch !=
kNull)
360 indicesL1.erase(bestMatch);
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 :
367 bestDeltaR =
cutsDr_[level - 2];
369 for (it = indicesHlt[
j].
begin(); it != indicesHlt[
j].end(); it++) {
371 candsHlt[
j][*it]->eta(), candsHlt[
j][*it]->phi());
372 if (dR < bestDeltaR) {
377 if (bestMatch !=
kNull)
379 indicesHlt[
j].erase(bestMatch);
403 string sourceUpper =
source;
404 sourceUpper[0] = toupper(sourceUpper[0]);
405 string name = source +
"Pass" + type +
"_" +
label;
408 if (type.find(
"MaxPt") != string::npos) {
409 string desc = (type ==
"MaxPt1") ?
"Leading" :
"Next-to-Leading";
410 string title =
"pT of " + desc +
" " + sourceUpper +
" Muon "+
411 "matched to " +
label;
413 float *
edges =
new float[nBins + 1];
415 h =
new TH1F(name.c_str(), title.c_str(), nBins,
edges);
419 string symbol = (type ==
"Eta") ?
"#eta" :
"#phi";
420 string title = symbol +
" of " + sourceUpper +
" Muons " +
421 "matched to " +
label;
423 int nBins = (int)params[0];
424 double min = params[1];
425 double max = params[2];
426 h =
new TH1F(name.c_str(), title.c_str(), nBins,
min,
max);
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...
bool getByToken(EDGetToken token, Handle< PROD > &result) const
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::vector< double > parametersEta_
L1MuonMatcherAlgo l1Matcher_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
std::vector< std::string > moduleLabels_
std::vector< Muon > MuonCollection
collection of Muon objects
static boost::tuple< edm::EDGetTokenT< trigger::TriggerEventWithRefs >, edm::EDGetTokenT< reco::GenParticleCollection >, edm::EDGetTokenT< reco::MuonCollection > > getTokens(const edm::ParameterSet &, edm::ConsumesCollector &&)
vector< ParameterSet > Parameters
StringCutObjectSelector< reco::Muon > * recMuonSelector_
std::vector< double > cutsDr_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
edm::EDGetTokenT< reco::MuonCollection > recMuonLabel_
MonitorElement * book1D(Args &&...args)
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
void analyze(const edm::Event &, const edm::EventSetup &)
bool insert(Storage &iStorage, ItemType *iItem, const IdTag &iIdTag)
StringCutObjectSelector< reco::GenParticle > * genMuonSelector_
double deltaR(double eta1, double eta2, double phi1, double phi2)
std::vector< std::string > stepLabels_
void findMatches(std::vector< MatchStruct > &, const std::vector< l1extra::L1MuonParticleRef > &, const std::vector< std::vector< const reco::RecoChargedCandidate * > > &)
void setCurrentFolder(const std::string &fullpath)
void bookHist(DQMStore::IBooker &, std::string, std::string, std::string, std::string)
std::vector< double > parametersPhi_
MonitorElement * bookFloat(Args &&...args)
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
tuple size
Write out results.
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity