CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/DQMOffline/Trigger/src/HLTMuonMatchAndPlot.cc

Go to the documentation of this file.
00001 
00009 #include "DQMOffline/Trigger/interface/HLTMuonMatchAndPlot.h"
00010 
00011 #include "FWCore/Framework/interface/Run.h"
00012 #include "FWCore/Framework/interface/EventSetup.h"
00013 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00014 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "DataFormats/Common/interface/Handle.h"
00018 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
00019 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00020 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00021 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00022 #include "DataFormats/VertexReco/interface/Vertex.h"
00023 #include "DataFormats/Common/interface/TriggerResults.h"
00024 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00025 
00026 #include <iostream>
00027 
00030 
00031 using namespace std;
00032 using namespace edm;
00033 using namespace reco;
00034 using namespace trigger;
00035 using namespace l1extra;
00036 
00037 typedef std::vector<std::string> vstring;
00038 
00039 
00042 
00044 HLTMuonMatchAndPlot::HLTMuonMatchAndPlot(const ParameterSet & pset, 
00045                                          string hltPath, 
00046                                          vector<string> moduleLabels) :
00047   hltProcessName_(pset.getParameter<string>("hltProcessName")),
00048   destination_(pset.getUntrackedParameter<string>("destination")),
00049   requiredTriggers_(pset.getUntrackedParameter<vstring>("requiredTriggers")),
00050   targetParams_(pset.getParameterSet("targetParams")),
00051   probeParams_(pset.getParameterSet("probeParams")),
00052   hltPath_(hltPath),
00053   moduleLabels_(moduleLabels),
00054   hasTargetRecoCuts(targetParams_.exists("recoCuts")),
00055   hasProbeRecoCuts(probeParams_.exists("recoCuts")),
00056   targetMuonSelector_(targetParams_.getUntrackedParameter<string>("recoCuts", "")),
00057   targetZ0Cut_(targetParams_.getUntrackedParameter<double>("z0Cut",0.)),
00058   targetD0Cut_(targetParams_.getUntrackedParameter<double>("d0Cut",0.)),
00059   probeMuonSelector_(probeParams_.getUntrackedParameter<string>("recoCuts", "")),
00060   probeZ0Cut_(probeParams_.getUntrackedParameter<double>("z0Cut",0.)),
00061   probeD0Cut_(probeParams_.getUntrackedParameter<double>("d0Cut",0.))
00062 {
00063 
00064   // Create std::map<string, T> from ParameterSets. 
00065   fillMapFromPSet(inputTags_, pset, "inputTags");
00066   fillMapFromPSet(binParams_, pset, "binParams");
00067   fillMapFromPSet(plotCuts_, pset, "plotCuts");
00068 
00069   // Set HLT process name for TriggerResults and TriggerSummary.
00070   InputTag & resTag = inputTags_["triggerResults"];
00071   resTag = InputTag(resTag.label(), resTag.instance(), hltProcessName_);
00072   InputTag & sumTag = inputTags_["triggerSummary"];
00073   sumTag = InputTag(sumTag.label(), sumTag.instance(), hltProcessName_);
00074 
00075   // Prepare the DQMStore object.
00076   dbe_ = edm::Service<DQMStore>().operator->();
00077   dbe_->setVerbose(0);
00078 
00079   // Get the trigger level.
00080   triggerLevel_ = "L3";
00081   TPRegexp levelRegexp("L[1-3]");
00082   size_t nModules = moduleLabels_.size();
00083   TObjArray * levelArray = levelRegexp.MatchS(moduleLabels_[nModules - 1]);
00084   if (levelArray->GetEntriesFast() > 0) {
00085     triggerLevel_ = ((TObjString *)levelArray->At(0))->GetString();
00086   }
00087   delete levelArray;
00088 
00089   // Get the pT cut by parsing the name of the HLT path.
00090   cutMinPt_ = 3;
00091   TPRegexp ptRegexp("Mu([0-9]*)");
00092   TObjArray * objArray = ptRegexp.MatchS(hltPath_);
00093   if (objArray->GetEntriesFast() >= 2) {
00094     TObjString * ptCutString = (TObjString *)objArray->At(1);
00095     cutMinPt_ = atoi(ptCutString->GetString());
00096     cutMinPt_ = ceil(cutMinPt_ * plotCuts_["minPtFactor"]);
00097   }
00098   delete objArray;
00099 }
00100 
00101 void HLTMuonMatchAndPlot::beginRun(const edm::Run& iRun, 
00102                                    const edm::EventSetup& iSetup)
00103 {
00104 
00105   TPRegexp suffixPtCut("Mu[0-9]+$");
00106 
00107   string baseDir = destination_;
00108   if (baseDir[baseDir.size() - 1] != '/') baseDir += '/';
00109   string pathSansSuffix = hltPath_;
00110   if (hltPath_.rfind("_v") < hltPath_.length())
00111     pathSansSuffix = hltPath_.substr(0, hltPath_.rfind("_v"));
00112   dbe_->setCurrentFolder(baseDir + pathSansSuffix);
00113 
00114   // Form is book1D(name, binningType, title) where 'binningType' is used 
00115   // to fetch the bin settings from binParams_.
00116   book1D("deltaR", "deltaR", ";#Deltar(reco, HLT);");
00117   book1D("hltPt", "pt", ";p_{T} of HLT object");
00118   book1D("hltEta", "eta", ";#eta of HLT object");
00119   book1D("hltPhi", "phi", ";#phi of HLT object");
00120   book1D("resolutionEta", "resolutionEta", ";#eta^{reco}-#eta^{HLT};");
00121   book1D("resolutionPhi", "resolutionPhi", ";#phi^{reco}-#phi^{HLT};");
00122   book1D("resolutionPt", "resolutionRel", 
00123          ";(p_{T}^{reco}-p_{T}^{HLT})/|p_{T}^{reco}|;");
00124 
00125   for (size_t i = 0; i < 2; i++) {
00126 
00127     string suffix = EFFICIENCY_SUFFIXES[i];
00128 
00129     book1D("efficiencyEta_" + suffix, "eta", ";#eta;");
00130     book1D("efficiencyPhi_" + suffix, "phi", ";#phi;");
00131     book1D("efficiencyTurnOn_" + suffix, "pt", ";p_{T};");
00132     book1D("efficiencyD0_" + suffix, "d0", ";d0;");
00133     book1D("efficiencyZ0_" + suffix, "z0", ";z0;");
00134     book1D("efficiencyCharge_" + suffix, "charge", ";charge;");
00135     book1D("efficiencyVertex_" + suffix, "NVertex", ";NVertex;");
00136 
00137     book2D("efficiencyPhiVsEta_" + suffix, "etaCoarse", "phiCoarse", 
00138            ";#eta;#phi");
00139 
00140     book1D("fakerateEta_" + suffix, "eta", ";#eta;");
00141     book1D("fakerateVertex_" + suffix, "NVertex", ";NVertex;");
00142     book1D("fakeratePhi_" + suffix, "phi", ";#phi;");
00143     book1D("fakerateTurnOn_" + suffix, "pt", ";p_{T};");
00144 
00145     book1D("massVsEtaZ_" + suffix, "etaCoarse", ";#eta");
00146     book1D("massVsEtaJpsi_" + suffix, "etaCoarse", ";#eta");
00147     book1D("massVsPtZ_" + suffix, "ptCoarse", ";p_{T}");
00148     book1D("massVsPtJpsi_" + suffix, "ptCoarse", ";p_{T}");
00149     book1D("massVsVertexZ_" + suffix, "NVertex", ";NVertex");
00150     book1D("massVsVertexJpsi_" + suffix, "NVertex", ";NVertex");
00151 
00152   }
00153   
00154 }
00155 
00156 
00157 
00158 void HLTMuonMatchAndPlot::endRun(const edm::Run& iRun, 
00159                                  const edm::EventSetup& iSetup)
00160 {
00161 }
00162 
00163 
00164 
00165 void HLTMuonMatchAndPlot::analyze(const Event & iEvent,
00166                                   const edm::EventSetup& iSetup)
00167 
00168 {
00169 
00170   // Get objects from the event.
00171   Handle<MuonCollection> allMuons;
00172   iEvent.getByLabel(inputTags_["recoMuon"], allMuons);
00173 
00174   Handle<BeamSpot> beamSpot;
00175   iEvent.getByLabel(inputTags_["beamSpot"], beamSpot);
00176   
00177   Handle<TriggerEvent> triggerSummary;
00178   iEvent.getByLabel(inputTags_["triggerSummary"], triggerSummary);
00179   Handle<TriggerResults> triggerResults;
00180 
00181   edm::Handle<VertexCollection> vertices;
00182   iEvent.getByLabel("offlinePrimaryVertices", vertices);
00183 
00184   //edm::Handle<GenParticleCollection> gen;
00185   //iEvent.getByLabel("genParticles", gen);
00186   //GenParticleCollection::const_iterator g_part;
00187   //std::vector<const GenParticle *> gen_leptsp;
00188   //std::vector<const GenParticle *> gen_momsp;
00189 
00190 
00191   if(!triggerSummary.isValid()) 
00192   {
00193     edm::LogError("HLTMuonMatchAndPlot")<<"Missing triggerSummary with label " << inputTags_["triggerSummary"] <<std::endl;
00194     return;
00195   }
00196   iEvent.getByLabel(inputTags_["triggerResults"], triggerResults);
00197   if(!triggerResults.isValid()) 
00198   {
00199     edm::LogError("HLTMuonMatchAndPlot")<<"Missing triggerResults with label " << inputTags_["triggerResults"] <<std::endl;
00200     return;
00201   }
00202 
00203   /*
00204   if(gen != 0) {
00205     for(g_part = gen->begin(); g_part != gen->end(); g_part++){
00206       if( abs(g_part->pdgId()) ==  13) {
00207         if(!( g_part->status() ==1 || (g_part->status() ==2 && abs(g_part->pdgId())==5))) continue;
00208         bool GenMomExists  (true);
00209         bool GenGrMomExists(true);
00210         if( g_part->pt() < 10.0 )  continue;
00211         if( fabs(g_part->eta()) > 2.4 ) continue;
00212         int gen_id= g_part->pdgId();
00213         const GenParticle* gen_lept = &(*g_part);
00214         // get mother of gen_lept
00215         const GenParticle* gen_mom = static_cast<const GenParticle*> (gen_lept->mother());
00216         if(gen_mom==NULL) GenMomExists=false;
00217         int m_id=-999;
00218         if(GenMomExists) m_id = gen_mom->pdgId();
00219         if(m_id != gen_id || !GenMomExists);
00220         else{
00221           int id= m_id;
00222           while(id == gen_id && GenMomExists){
00223             gen_mom = static_cast<const GenParticle*> (gen_mom->mother());
00224             if(gen_mom==NULL){
00225               GenMomExists=false;
00226             }
00227             if(GenMomExists) id=gen_mom->pdgId();
00228           }
00229         }
00230         if(GenMomExists) m_id = gen_mom->pdgId();
00231         gen_leptsp.push_back(gen_lept);
00232         gen_momsp.push_back(gen_mom);
00233       }
00234     }
00235   }
00236 
00237 
00238   std::vector<GenParticle> gen_lepts;
00239   for(size_t i = 0; i < gen_leptsp.size(); i++) {
00240     gen_lepts.push_back(*gen_leptsp[i]);
00241   }
00242 
00243 
00244   */
00245 
00246 
00247 
00248   // Throw out this event if it doesn't pass the required triggers.
00249   for (size_t i = 0; i < requiredTriggers_.size(); i++) {
00250     unsigned int triggerIndex = triggerResults->find(requiredTriggers_[i]);
00251     if (triggerIndex < triggerResults->size() ||
00252         !triggerResults->accept(triggerIndex))
00253       return;
00254   }
00255 
00256   // Select objects based on the configuration.
00257   MuonCollection targetMuons = selectedMuons(* allMuons, * beamSpot, hasTargetRecoCuts, targetMuonSelector_, targetD0Cut_, targetZ0Cut_);
00258   MuonCollection probeMuons = selectedMuons(* allMuons, * beamSpot, hasProbeRecoCuts, probeMuonSelector_, probeD0Cut_, probeZ0Cut_);
00259   TriggerObjectCollection allTriggerObjects = triggerSummary->getObjects();
00260   TriggerObjectCollection hltMuons = 
00261     selectedTriggerObjects(allTriggerObjects, * triggerSummary, targetParams_);
00262 
00263   // Fill plots for HLT muons.
00264   for (size_t i = 0; i < hltMuons.size(); i++) {
00265     hists_["hltPt"]->Fill(hltMuons[i].pt());
00266     hists_["hltEta"]->Fill(hltMuons[i].eta());
00267     hists_["hltPhi"]->Fill(hltMuons[i].phi());
00268   }
00269 
00270   // Find the best trigger object matches for the targetMuons.
00271   vector<size_t> matches = matchByDeltaR(targetMuons, hltMuons, 
00272                                          plotCuts_[triggerLevel_ + "DeltaR"]);
00273 
00274 
00275   // Fill plots for matched muons.
00276   bool pairalreadyconsidered = false;
00277   for (size_t i = 0; i < targetMuons.size(); i++) {
00278 
00279     Muon & muon = targetMuons[i];
00280 
00281     // Fill plots which are not efficiencies.
00282     if (matches[i] < targetMuons.size()) {
00283       TriggerObject & hltMuon = hltMuons[matches[i]];
00284       double ptRes = (muon.pt() - hltMuon.pt()) / muon.pt();
00285       double etaRes = muon.eta() - hltMuon.eta();
00286       double phiRes = muon.phi() - hltMuon.phi();
00287       hists_["resolutionEta"]->Fill(etaRes);
00288       hists_["resolutionPhi"]->Fill(phiRes);
00289       hists_["resolutionPt"]->Fill(ptRes);
00290       hists_["deltaR"]->Fill(deltaR(muon, hltMuon));
00291     }
00292 
00293     // Fill numerators and denominators for efficiency plots.
00294     for (size_t j = 0; j < 2; j++) {
00295 
00296       string suffix = EFFICIENCY_SUFFIXES[j];
00297 
00298       // If no match was found, then the numerator plots don't get filled.
00299       if (suffix == "numer" && matches[i] >= targetMuons.size()) continue;
00300 
00301       if (muon.pt() > cutMinPt_) {
00302         hists_["efficiencyEta_" + suffix]->Fill(muon.eta());
00303         hists_["efficiencyPhiVsEta_" + suffix]->Fill(muon.eta(), muon.phi());
00304       }
00305       
00306       if (fabs(muon.eta()) < plotCuts_["maxEta"]) {
00307         hists_["efficiencyTurnOn_" + suffix]->Fill(muon.pt());
00308       }
00309       
00310       if (muon.pt() > cutMinPt_ && fabs(muon.eta()) < plotCuts_["maxEta"]) {
00311         const Track * track = 0;
00312         if (muon.isTrackerMuon()) track = & * muon.innerTrack();
00313         else if (muon.isStandAloneMuon()) track = & * muon.outerTrack();
00314         if (track) {
00315           double d0 = track->dxy(beamSpot->position());
00316           double z0 = track->dz(beamSpot->position());
00317           hists_["efficiencyVertex_" + suffix]->Fill(vertices->size());
00318           hists_["efficiencyPhi_" + suffix]->Fill(muon.phi());
00319           hists_["efficiencyD0_" + suffix]->Fill(d0);
00320           hists_["efficiencyZ0_" + suffix]->Fill(z0);
00321           hists_["efficiencyCharge_" + suffix]->Fill(muon.charge());
00322         }
00323       }
00324     }
00325       
00326     // Fill plots for tag and probe
00327     // Muon cannot be a tag because doesn't match an hlt muon     
00328     if(matches[i] >= targetMuons.size()) continue;
00329     for (size_t k = 0; k < targetMuons.size(); k++) {
00330       if(k == i) continue;
00331       if(muon.pt() < 20.0) continue;
00332       Muon & theProbe = targetMuons[k];
00333       if (muon.charge() != theProbe.charge() && !pairalreadyconsidered) {
00334         double mass = (muon.p4() + theProbe.p4()).M();
00335         if(mass > 60 && mass < 120) {
00336           hists_["massVsEtaZ_denom"]->Fill(theProbe.eta());
00337           hists_["massVsPtZ_denom"]->Fill(theProbe.pt());
00338           hists_["massVsVertexZ_denom"]->Fill(vertices->size());
00339           if(matches[k] < targetMuons.size()) {
00340             hists_["massVsEtaZ_numer"]->Fill(theProbe.eta());
00341             hists_["massVsPtZ_numer"]->Fill(theProbe.pt());
00342             hists_["massVsVertexZ_numer"]->Fill(vertices->size());
00343           }  
00344           pairalreadyconsidered = true;
00345         }
00346         if(mass > 1 && mass < 4) {
00347           hists_["massVsEtaJpsi_denom"]->Fill(theProbe.eta());
00348           hists_["massVsPtJpsi_denom"]->Fill(theProbe.pt());
00349           hists_["massVsVertexJpsi_denom"]->Fill(vertices->size());
00350           if(matches[k] < targetMuons.size()) {
00351             hists_["massVsEtaJpsi_numer"]->Fill(theProbe.eta());
00352             hists_["massVsPtJpsi_numer"]->Fill(theProbe.pt());
00353             hists_["massVsVertexJpsi_numer"]->Fill(vertices->size());
00354           }
00355           pairalreadyconsidered = true;
00356         }
00357       }
00358     } // End loop over denominator and numerator.
00359   } // End loop over targetMuons.
00360 
00361   // Plot fake rates (efficiency for HLT objects to not get matched to RECO).
00362   vector<size_t> hltMatches = matchByDeltaR(hltMuons, targetMuons,
00363                                             plotCuts_[triggerLevel_ + "DeltaR"]);
00364   for (size_t i = 0; i < hltMuons.size(); i++) {
00365     TriggerObject & hltMuon = hltMuons[i];
00366     bool isFake = hltMatches[i] > hltMuons.size();
00367     for (size_t j = 0; j < 2; j++) {
00368       string suffix = EFFICIENCY_SUFFIXES[j];
00369       // If match is found, then numerator plots should not get filled
00370       if (suffix == "numer" && ! isFake) continue;
00371       hists_["fakerateVertex_" + suffix]->Fill(vertices->size());
00372       hists_["fakerateEta_" + suffix]->Fill(hltMuon.eta());
00373       hists_["fakeratePhi_" + suffix]->Fill(hltMuon.phi());
00374       hists_["fakerateTurnOn_" + suffix]->Fill(hltMuon.pt());
00375     } // End loop over numerator and denominator.
00376   } // End loop over hltMuons.
00377   
00378 
00379 } // End analyze() method.
00380 
00381 
00382 
00383 // Method to fill binning parameters from a vector of doubles.
00384 void 
00385 HLTMuonMatchAndPlot::fillEdges(size_t & nBins, float * & edges, 
00386                                vector<double> binning) {
00387 
00388   if (binning.size() < 3) {
00389     LogWarning("HLTMuonVal") << "Invalid binning parameters!"; 
00390     return;
00391   }
00392 
00393   // Fixed-width binning.
00394   if (binning.size() == 3) {
00395     nBins = binning[0];
00396     edges = new float[nBins + 1];
00397     const double min = binning[1];
00398     const double binwidth = (binning[2] - binning[1]) / nBins;
00399     for (size_t i = 0; i <= nBins; i++) edges[i] = min + (binwidth * i);
00400   } 
00401 
00402   // Variable-width binning.
00403   else {
00404     nBins = binning.size() - 1;
00405     edges = new float[nBins + 1];
00406     for (size_t i = 0; i <= nBins; i++) edges[i] = binning[i];
00407   }
00408 
00409 }
00410 
00411 
00412 
00413 // This is an unorthodox method of getting parameters, but cleaner in my mind
00414 // It looks for parameter 'target' in the pset, and assumes that all entries
00415 // have the same class (T), filling the values into a std::map.
00416 template <class T>
00417 void 
00418 HLTMuonMatchAndPlot::fillMapFromPSet(map<string, T> & m, 
00419                                      ParameterSet pset, string target) {
00420 
00421   // Get the ParameterSet with name 'target' from 'pset'
00422   ParameterSet targetPset;
00423   if (pset.existsAs<ParameterSet>(target, true))        // target is tracked
00424     targetPset = pset.getParameterSet(target);
00425   else if (pset.existsAs<ParameterSet>(target, false))  // target is untracked
00426     targetPset = pset.getUntrackedParameterSet(target);
00427   
00428   // Get the parameter names from targetPset
00429   vector<string> names = targetPset.getParameterNames();
00430   vector<string>::const_iterator iter;
00431   
00432   for (iter = names.begin(); iter != names.end(); ++iter) {
00433     if (targetPset.existsAs<T>(* iter, true))           // target is tracked
00434       m[* iter] = targetPset.getParameter<T>(* iter);
00435     else if (targetPset.existsAs<T>(* iter, false))     // target is untracked
00436       m[* iter] = targetPset.getUntrackedParameter<T>(* iter);
00437   }
00438 
00439 }
00440 
00441 
00442 
00443 // A generic method to find the best deltaR matches from 2 collections.
00444 template <class T1, class T2> 
00445 vector<size_t> 
00446 HLTMuonMatchAndPlot::matchByDeltaR(const vector<T1> & collection1, 
00447                                    const vector<T2> & collection2,
00448                                    const double maxDeltaR) {
00449 
00450   const size_t n1 = collection1.size();
00451   const size_t n2 = collection2.size();
00452 
00453   vector<size_t> result(n1, -1);
00454   vector<vector<double> > deltaRMatrix(n1, vector<double>(n2, NOMATCH));
00455 
00456   for (size_t i = 0; i < n1; i++)
00457     for (size_t j = 0; j < n2; j++) {
00458       deltaRMatrix[i][j] = deltaR(collection1[i], collection2[j]);
00459     }
00460 
00461   // Run through the matrix n1 times to make sure we've found all matches.
00462   for (size_t k = 0; k < n1; k++) {
00463     size_t i_min = -1;
00464     size_t j_min = -1;
00465     double minDeltaR = maxDeltaR;
00466     // find the smallest deltaR
00467     for (size_t i = 0; i < n1; i++)
00468       for (size_t j = 0; j < n2; j++)
00469         if (deltaRMatrix[i][j] < minDeltaR) {
00470           i_min = i;
00471           j_min = j;
00472           minDeltaR = deltaRMatrix[i][j];
00473         }
00474     // If a match has been made, save it and make those candidates unavailable.
00475     if (minDeltaR < maxDeltaR) {
00476       result[i_min] = j_min;
00477       deltaRMatrix[i_min] = vector<double>(n2, NOMATCH);
00478       for (size_t i = 0; i < n1; i++)
00479         deltaRMatrix[i][j_min] = NOMATCH;
00480     }
00481   }
00482 
00483   return result;
00484 
00485 }
00486 
00487 
00488 
00489 MuonCollection
00490 HLTMuonMatchAndPlot::selectedMuons(const MuonCollection & allMuons, 
00491                                    const BeamSpot & beamSpot,
00492                                    bool hasRecoCuts,
00493                                    const StringCutObjectSelector<reco::Muon> &selector,
00494                                    double d0Cut, double z0Cut)
00495 {
00496   // If there is no selector (recoCuts does not exists), return an empty collection. 
00497   if (!hasRecoCuts)
00498     return MuonCollection();
00499 
00500   MuonCollection reducedMuons(allMuons);
00501   MuonCollection::iterator iter = reducedMuons.begin();
00502   while (iter != reducedMuons.end()) {
00503     const Track * track = 0;
00504     if (iter->isTrackerMuon()) track = & * iter->innerTrack();
00505     else if (iter->isStandAloneMuon()) track = & * iter->outerTrack();
00506     if (track && selector(* iter) &&
00507         fabs(track->dxy(beamSpot.position())) < d0Cut &&
00508         fabs(track->dz(beamSpot.position())) < z0Cut)
00509       ++iter;
00510     else reducedMuons.erase(iter);
00511   }
00512 
00513   return reducedMuons;
00514 
00515 }
00516 
00517 
00518 
00519 TriggerObjectCollection
00520 HLTMuonMatchAndPlot::selectedTriggerObjects(
00521   const TriggerObjectCollection & triggerObjects,
00522   const TriggerEvent & triggerSummary,
00523   const ParameterSet & pset) 
00524 {
00525 
00526   // If pset is empty, return an empty collection.
00527   if (!pset.exists("hltCuts"))
00528     return TriggerObjectCollection();
00529 
00530   StringCutObjectSelector<TriggerObject> selector
00531     (pset.getUntrackedParameter<string>("hltCuts"));
00532 
00533   InputTag filterTag(moduleLabels_[moduleLabels_.size() - 1], "", 
00534                      hltProcessName_);
00535   size_t filterIndex = triggerSummary.filterIndex(filterTag);
00536 
00537   TriggerObjectCollection selectedObjects;
00538 
00539   if (filterIndex < triggerSummary.sizeFilters()) {
00540     const Keys &keys = triggerSummary.filterKeys(filterIndex);
00541     for (size_t j = 0; j < keys.size(); j++ ){
00542       TriggerObject foundObject = triggerObjects[keys[j]];
00543       if (selector(foundObject))
00544         selectedObjects.push_back(foundObject);
00545     }
00546   }
00547 
00548   return selectedObjects;
00549 
00550 }
00551 
00552 
00553 
00554 void HLTMuonMatchAndPlot::book1D(string name, string binningType, string title)
00555 {
00556   /* Properly delete the array of floats that has been allocated on
00557    * the heap by fillEdges.  Avoid multiple copies and internal ROOT
00558    * clones by simply creating the histograms directly in the DQMStore
00559    * using the appropriate book1D method to handle the variable bins
00560    * case. */ 
00561 
00562   size_t nBins; 
00563   float * edges = 0; 
00564   fillEdges(nBins, edges, binParams_[binningType]);
00565 
00566   hists_[name] = dbe_->book1D(name, title, nBins, edges);
00567   if (hists_[name])
00568     if (hists_[name]->getTH1F()->GetSumw2N())
00569       hists_[name]->getTH1F()->Sumw2();
00570 
00571   if (edges)
00572     delete [] edges;
00573 }
00574 
00575 
00576 
00577 void
00578 HLTMuonMatchAndPlot::book2D(string name, string binningTypeX, 
00579                             string binningTypeY, string title) 
00580 {
00581   
00582   /* Properly delete the arrays of floats that have been allocated on
00583    * the heap by fillEdges.  Avoid multiple copies and internal ROOT
00584    * clones by simply creating the histograms directly in the DQMStore
00585    * using the appropriate book2D method to handle the variable bins
00586    * case. */ 
00587 
00588   size_t  nBinsX;
00589   float * edgesX = 0;
00590   fillEdges(nBinsX, edgesX, binParams_[binningTypeX]);
00591 
00592   size_t  nBinsY;
00593   float * edgesY = 0;
00594   fillEdges(nBinsY, edgesY, binParams_[binningTypeY]);
00595 
00596   hists_[name] = dbe_->book2D(name.c_str(), title.c_str(),
00597                               nBinsX, edgesX, nBinsY, edgesY);
00598   if (hists_[name])
00599     if (hists_[name]->getTH2F()->GetSumw2N())
00600       hists_[name]->getTH2F()->Sumw2();
00601 
00602   if (edgesX)
00603     delete [] edgesX;
00604   if (edgesY)
00605     delete [] edgesY;
00606 }
00607