CMS 3D CMS Logo

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