00001
00008 #include "HLTriggerOffline/Muon/interface/HLTMuonValidator.h"
00009 #include "FWCore/ServiceRegistry/interface/Service.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00013 #include "DataFormats/Math/interface/deltaPhi.h"
00014 #include "DataFormats/Math/interface/deltaR.h"
00015 #include "DataFormats/MuonReco/interface/Muon.h"
00016 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00017
00018 #include "DataFormats/MuonSeed/interface/L2MuonTrajectorySeed.h"
00019 #include "DataFormats/MuonSeed/interface/L2MuonTrajectorySeedCollection.h"
00020 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
00021 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
00022
00023
00024
00025 using namespace std;
00026 using namespace edm;
00027 using namespace reco;
00028 using namespace trigger;
00029 using namespace l1extra;
00030
00031
00032
00033 typedef vector<ParameterSet> Parameters;
00034
00035
00036
00037 HLTMuonValidator::HLTMuonValidator(const ParameterSet & pset) :
00038 l1Matcher_(pset)
00039 {
00040
00041 hltProcessName_ = pset.getParameter< string >("hltProcessName");
00042 hltPathsToCheck_ = pset.getParameter< vector<string> >("hltPathsToCheck");
00043
00044 genParticleLabel_ = pset.getParameter<string>("genParticleLabel");
00045 recMuonLabel_ = pset.getParameter<string>( "recMuonLabel");
00046 l1CandLabel_ = pset.getParameter<string>( "l1CandLabel");
00047 l2CandLabel_ = pset.getParameter<string>( "l2CandLabel");
00048 l3CandLabel_ = pset.getParameter<string>( "l3CandLabel");
00049
00050 cutsDr_ = pset.getParameter< vector<double> >("cutsDr" );
00051
00052 parametersEta_ = pset.getParameter< vector<double> >("parametersEta");
00053 parametersPhi_ = pset.getParameter< vector<double> >("parametersPhi");
00054 parametersTurnOn_ = pset.getParameter< vector<double> >("parametersTurnOn");
00055
00056 genMuonCut_ = pset.getParameter<string>("genMuonCut");
00057 recMuonCut_ = pset.getParameter<string>("recMuonCut");
00058
00059 genMuonSelector_ = 0;
00060 recMuonSelector_ = 0;
00061
00062 dbe_ = Service<DQMStore>().operator->();
00063 dbe_->setVerbose(0);
00064
00065 }
00066
00067
00068
00069 void
00070 HLTMuonValidator::beginJob()
00071 {
00072 }
00073
00074
00075
00076 void
00077 HLTMuonValidator::beginRun(const Run & iRun, const EventSetup & iSetup)
00078 {
00079
00080 static int runNumber = 0;
00081 runNumber++;
00082
00083 l1Matcher_.init(iSetup);
00084
00085 bool changedConfig;
00086 if (!hltConfig_.init(iRun, iSetup, hltProcessName_, changedConfig)) {
00087 LogError("HLTMuonVal") << "Initialization of HLTConfigProvider failed!!";
00088 return;
00089 }
00090
00091 hltPaths_.clear();
00092 filterLabels_.clear();
00093 vector<string> validTriggerNames = hltConfig_.triggerNames();
00094 for (size_t i = 0; i < hltPathsToCheck_.size(); i++) {
00095 TPRegexp pattern(hltPathsToCheck_[i]);
00096 for (size_t j = 0; j < validTriggerNames.size(); j++)
00097 if (TString(validTriggerNames[j]).Contains(pattern))
00098 hltPaths_.insert(validTriggerNames[j]);
00099
00100 if (TString("NoFilters").Contains(pattern))
00101 hltPaths_.insert("NoFilters");
00102 }
00103
00104 initializeHists();
00105
00106 }
00107
00108
00109
00110 void
00111 HLTMuonValidator::initializeHists()
00112 {
00113
00114 vector<string> sources(2);
00115 sources[0] = "gen";
00116 sources[1] = "rec";
00117
00118 set<string>::iterator iPath;
00119 TPRegexp suffixPtCut("[0-9]+$");
00120
00121 for (iPath = hltPaths_.begin(); iPath != hltPaths_.end(); iPath++) {
00122
00123 string path = * iPath;
00124
00125 if (path == "NoFilters") {
00126 filterLabels_[path].push_back("hltL1sL1SingleMuOpenL1SingleMu0");
00127
00128 filterLabels_[path].push_back("");
00129 filterLabels_[path].push_back("");
00130 }
00131
00132 vector<string> moduleLabels;
00133 if (path != "NoFilters") moduleLabels = hltConfig_.moduleLabels(path);
00134
00135 for (size_t i = 0; i < moduleLabels.size(); i++)
00136 if (moduleLabels[i].find("Filtered") != string::npos)
00137 filterLabels_[path].push_back(moduleLabels[i]);
00138
00139
00140
00141 double cutMaxEta = (TString(path).Contains(kLooseL1Requirement)) ? 2.4:2.1;
00142
00143
00144 unsigned int index = TString(path).Index(suffixPtCut);
00145 unsigned int threshold = 3;
00146 if (index < path.length()) threshold = atoi(path.substr(index).c_str());
00147
00148
00149
00150 double cutMinPt = ceil(threshold * 1.1) - 0.01;
00151 if (cutMinPt < 0. || path == "NoFilters") cutMinPt = 0.;
00152 cutsMinPt_[path] = cutMinPt;
00153
00154 string baseDir = "HLT/Muon/Distributions/";
00155 dbe_->setCurrentFolder(baseDir + path);
00156
00157 if (dbe_->get(baseDir + path + "/CutMinPt") == 0) {
00158
00159 elements_[path + "_" + "CutMinPt" ] = dbe_->bookFloat("CutMinPt" );
00160 elements_[path + "_" + "CutMaxEta"] = dbe_->bookFloat("CutMaxEta");
00161 elements_[path + "_" + "CutMinPt" ]->Fill(cutMinPt);
00162 elements_[path + "_" + "CutMaxEta"]->Fill(cutMaxEta);
00163
00164
00165 const int nFilters = filterLabels_[path].size();
00166 stepLabels_[path].push_back("All");
00167 stepLabels_[path].push_back("L1");
00168 if (nFilters == 2) {
00169 stepLabels_[path].push_back("L2");
00170 }
00171 if (nFilters == 3) {
00172 stepLabels_[path].push_back("L2");
00173 stepLabels_[path].push_back("L3");
00174 }
00175 if (nFilters == 5) {
00176 stepLabels_[path].push_back("L2");
00177 stepLabels_[path].push_back("L2Iso");
00178 stepLabels_[path].push_back("L3");
00179 stepLabels_[path].push_back("L3Iso");
00180 }
00181
00182
00183
00184
00185
00186
00187
00188 for (size_t i = 0; i < sources.size(); i++) {
00189 string source = sources[i];
00190 for (size_t j = 0; j < stepLabels_[path].size(); j++) {
00191 bookHist(path, stepLabels_[path][j], source, "Eta");
00192 bookHist(path, stepLabels_[path][j], source, "Phi");
00193 bookHist(path, stepLabels_[path][j], source, "MaxPt1");
00194 bookHist(path, stepLabels_[path][j], source, "MaxPt2");
00195 }
00196 }
00197
00198 }
00199
00200 }
00201
00202 }
00203
00204
00205
00206 void
00207 HLTMuonValidator::analyze(const Event & iEvent, const EventSetup & iSetup)
00208 {
00209
00210 static int eventNumber = 0;
00211 eventNumber++;
00212 LogTrace("HLTMuonVal") << "In HLTMuonValidator::analyze, "
00213 << "Event: " << eventNumber;
00214
00215 Handle< TriggerEventWithRefs> rawTriggerEvent;
00216 Handle< MuonCollection> recMuons;
00217 Handle< GenParticleCollection> genParticles;
00218
00219 iEvent.getByLabel("hltTriggerSummaryRAW", rawTriggerEvent);
00220 if (rawTriggerEvent.failedToGet())
00221 {LogError("HLTMuonVal") << "No trigger summary found"; return;}
00222 iEvent.getByLabel( recMuonLabel_, recMuons );
00223 iEvent.getByLabel(genParticleLabel_, genParticles );
00224
00225 vector<string> sources;
00226 if (genParticles.isValid()) sources.push_back("gen");
00227 if ( recMuons.isValid()) sources.push_back("rec");
00228
00229 for (size_t sourceNo = 0; sourceNo < sources.size(); sourceNo++) {
00230
00231 string source = sources[sourceNo];
00232
00233
00234 if (!genMuonSelector_) genMuonSelector_ =
00235 new StringCutObjectSelector<reco::GenParticle>(genMuonCut_);
00236 if (!recMuonSelector_) recMuonSelector_ =
00237 new StringCutObjectSelector<reco::Muon >(recMuonCut_);
00238
00239
00240 vector<MatchStruct> matches;
00241 if (source == "gen" && genParticles.isValid())
00242 for (size_t i = 0; i < genParticles->size(); i++)
00243 if ((*genMuonSelector_)(genParticles->at(i)))
00244 matches.push_back(MatchStruct(& genParticles->at(i)));
00245 if (source == "rec" && recMuons.isValid())
00246 for (size_t i = 0; i < recMuons->size(); i++)
00247 if ((*recMuonSelector_)(recMuons->at(i)))
00248 matches.push_back(MatchStruct(& recMuons->at(i)));
00249
00250
00251 sort(matches.begin(), matches.end(), matchesByDescendingPt());
00252
00253 set<string>::iterator iPath;
00254 for (iPath = hltPaths_.begin(); iPath != hltPaths_.end(); iPath++)
00255 analyzePath(iEvent, * iPath, source, matches, rawTriggerEvent);
00256
00257 }
00258
00259 }
00260
00261
00262
00263 void
00264 HLTMuonValidator::analyzePath(const Event & iEvent,
00265 const string & path,
00266 const string & source,
00267 vector<MatchStruct> matches,
00268 Handle<TriggerEventWithRefs> rawTriggerEvent)
00269 {
00270
00271 const bool skipFilters = (path == "NoFilters");
00272
00273 const float maxEta = elements_[path + "_" + "CutMaxEta"]->getFloatValue();
00274 const bool isDoubleMuonPath = (path.find("Double") != string::npos);
00275 const size_t nFilters = filterLabels_[path].size();
00276 const size_t nSteps = stepLabels_[path].size();
00277 const size_t nStepsHlt = nSteps - 2;
00278 const int nObjectsToPassPath = (isDoubleMuonPath) ? 2 : 1;
00279 vector< L1MuonParticleRef > candsL1;
00280 vector< vector< RecoChargedCandidateRef > > refsHlt(nStepsHlt);
00281 vector< vector< const RecoChargedCandidate * > > candsHlt(nStepsHlt);
00282
00283 for (size_t i = 0; i < nFilters; i++) {
00284 const int hltStep = i - 1;
00285 InputTag tag = InputTag(filterLabels_[path][i], "", hltProcessName_);
00286 size_t iFilter = rawTriggerEvent->filterIndex(tag);
00287 if (iFilter < rawTriggerEvent->size()) {
00288 if (i == 0)
00289 rawTriggerEvent->getObjects(iFilter, TriggerL1Mu, candsL1);
00290 else
00291 rawTriggerEvent->getObjects(iFilter, TriggerMuon,
00292 refsHlt[hltStep]);
00293 }
00294 else if (!skipFilters)
00295 LogTrace("HLTMuonVal") << "No collection with label " << tag;
00296 }
00297 if (skipFilters) {
00298 Handle<RecoChargedCandidateCollection> handleCandsL2;
00299 Handle<RecoChargedCandidateCollection> handleCandsL3;
00300 iEvent.getByLabel(l2CandLabel_, handleCandsL2);
00301 iEvent.getByLabel(l3CandLabel_, handleCandsL3);
00302 if (handleCandsL2.isValid())
00303 for (size_t i = 0; i < handleCandsL2->size(); i++)
00304 candsHlt[0].push_back(& handleCandsL2->at(i));
00305 if (handleCandsL3.isValid())
00306 for (size_t i = 0; i < handleCandsL3->size(); i++)
00307 candsHlt[1].push_back(& handleCandsL3->at(i));
00308 }
00309 else for (size_t i = 0; i < nStepsHlt; i++)
00310 for (size_t j = 0; j < refsHlt[i].size(); j++)
00311 candsHlt[i].push_back(& * refsHlt[i][j]);
00312
00313
00314 findMatches(matches, candsL1, candsHlt);
00315
00316 vector<size_t> matchesInEtaRange;
00317 vector<bool> hasMatch(matches.size(), true);
00318
00319 for (size_t step = 0; step < nSteps; step++) {
00320
00321 const size_t hltStep = (step >= 2) ? step - 2 : 0;
00322 const size_t level = (step == 1) ? 1 :
00323 (step == 2) ? 2 :
00324 (step == 3) ? ((nStepsHlt == 4) ? 2 : 3) :
00325 (step >= 4) ? 3 :
00326 0;
00327
00328 for (size_t j = 0; j < matches.size(); j++) {
00329 if (level == 0) {
00330 if (fabs(matches[j].candBase->eta()) < maxEta)
00331 matchesInEtaRange.push_back(j);
00332 }
00333 else if (level == 1) {
00334 if (matches[j].candL1 == 0)
00335 hasMatch[j] = false;
00336 }
00337 else if (level >= 2) {
00338 if (matches[j].candHlt[hltStep] == 0)
00339 hasMatch[j] = false;
00340 else if (!hasMatch[j]) {
00341 LogTrace("HLTMuonVal") << "Match found for HLT step " << hltStep
00342 << " of " << nStepsHlt
00343 << " without previous match!";
00344 break;
00345 }
00346 }
00347 }
00348
00349 if (std::count(hasMatch.begin(), hasMatch.end(), true) <
00350 nObjectsToPassPath)
00351 break;
00352
00353 string pre = path + "_" + source + "Pass";
00354 string post = "_" + stepLabels_[path][step];
00355
00356 for (size_t j = 0; j < matches.size(); j++) {
00357 float pt = matches[j].candBase->pt();
00358 float eta = matches[j].candBase->eta();
00359 float phi = matches[j].candBase->phi();
00360 if (hasMatch[j]) {
00361 if (matchesInEtaRange.size() >= 1 && j == matchesInEtaRange[0])
00362 elements_[pre + "MaxPt1" + post]->Fill(pt);
00363 if (matchesInEtaRange.size() >= 2 && j == matchesInEtaRange[1])
00364 elements_[pre + "MaxPt2" + post]->Fill(pt);
00365 if(fabs(eta) < maxEta && pt > cutsMinPt_[path]) {
00366 elements_[pre + "Eta" + post]->Fill(eta);
00367 elements_[pre + "Phi" + post]->Fill(phi);
00368 }
00369 }
00370 }
00371
00372 }
00373
00374 }
00375
00376
00377
00378 void
00379 HLTMuonValidator::findMatches(
00380 vector<MatchStruct> & matches,
00381 vector<L1MuonParticleRef> candsL1,
00382 vector< vector< const RecoChargedCandidate *> > candsHlt)
00383 {
00384
00385 set<size_t>::iterator it;
00386
00387 set<size_t> indicesL1;
00388 for (size_t i = 0; i < candsL1.size(); i++)
00389 indicesL1.insert(i);
00390
00391 vector< set<size_t> > indicesHlt(candsHlt.size());
00392 for (size_t i = 0; i < candsHlt.size(); i++)
00393 for (size_t j = 0; j < candsHlt[i].size(); j++)
00394 indicesHlt[i].insert(j);
00395
00396 for (size_t i = 0; i < matches.size(); i++) {
00397
00398 const Candidate * cand = matches[i].candBase;
00399
00400 double bestDeltaR = cutsDr_[0];
00401 size_t bestMatch = kNull;
00402 for (it = indicesL1.begin(); it != indicesL1.end(); it++) {
00403 double dR = deltaR(cand->eta(), cand->phi(),
00404 candsL1[*it]->eta(), candsL1[*it]->phi());
00405 if (dR < bestDeltaR) {
00406 bestMatch = *it;
00407 bestDeltaR = dR;
00408 }
00409
00410
00411
00412
00413
00414
00415
00416
00417 }
00418 if (bestMatch != kNull)
00419 matches[i].candL1 = & * candsL1[bestMatch];
00420 indicesL1.erase(bestMatch);
00421
00422 matches[i].candHlt.assign(candsHlt.size(), 0);
00423 for (size_t j = 0; j < candsHlt.size(); j++) {
00424 size_t level = (candsHlt.size() == 4) ? (j < 2) ? 2 : 3 :
00425 (candsHlt.size() == 2) ? (j < 1) ? 2 : 3 :
00426 2;
00427 bestDeltaR = cutsDr_[level - 2];
00428 bestMatch = kNull;
00429 for (it = indicesHlt[j].begin(); it != indicesHlt[j].end(); it++) {
00430 double dR = deltaR(cand->eta(), cand->phi(),
00431 candsHlt[j][*it]->eta(), candsHlt[j][*it]->phi());
00432 if (dR < bestDeltaR) {
00433 bestMatch = *it;
00434 bestDeltaR = dR;
00435 }
00436 }
00437 if (bestMatch != kNull)
00438 matches[i].candHlt[j] = candsHlt[j][bestMatch];
00439 indicesHlt[j].erase(bestMatch);
00440 }
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450 }
00451
00452 }
00453
00454
00455
00456
00457 void
00458 HLTMuonValidator::bookHist(string path, string label,
00459 string source, string type)
00460 {
00461
00462 string sourceUpper = source;
00463 sourceUpper[0] = toupper(sourceUpper[0]);
00464 string name = source + "Pass" + type + "_" + label;
00465 string rootName = path + "_" + name;
00466 TH1F * h;
00467
00468 if (type.find("MaxPt") != string::npos) {
00469 string desc = (type == "MaxPt1") ? "Leading" : "Next-to-Leading";
00470 string title = "pT of " + desc + " " + sourceUpper + " Muon "+
00471 "matched to " + label;
00472 const size_t nBins = parametersTurnOn_.size() - 1;
00473 float * edges = new float[nBins + 1];
00474 for (size_t i = 0; i < nBins + 1; i++) edges[i] = parametersTurnOn_[i];
00475 h = new TH1F(rootName.c_str(), title.c_str(), nBins, edges);
00476 }
00477
00478 else {
00479 string symbol = (type == "Eta") ? "#eta" : "#phi";
00480 string title = symbol + " of " + sourceUpper + " Muons " +
00481 "matched to " + label;
00482 vector<double> params = (type == "Eta") ? parametersEta_ : parametersPhi_;
00483 int nBins = (int)params[0];
00484 double min = params[1];
00485 double max = params[2];
00486 h = new TH1F(rootName.c_str(), title.c_str(), nBins, min, max);
00487 }
00488
00489 h->Sumw2();
00490 elements_[rootName] = dbe_->book1D(name, h);
00491 delete h;
00492
00493 }
00494
00495
00496
00497 DEFINE_FWK_MODULE(HLTMuonValidator);