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