CMS 3D CMS Logo

HLTMuonGenericRate Class Reference

Get L1/HLT efficiency/rate plots Documentation available on the CMS TWiki: httpss://twiki.cern.ch/twiki/bin/view/CMS/MuonHLTOfflinePerformance. More...

#include <HLTriggerOffline/Muon/interface/HLTMuonGenericRate.h>

List of all members.

Public Member Functions

void analyze (const edm::Event &iEvent)
void begin ()
MonitorElementbookIt (TString name, TString title, std::vector< double >)
void finish ()
 HLTMuonGenericRate (const edm::ParameterSet &pset, std::string triggerName, std::vector< std::string > moduleNames)
 Constructor.

Private Member Functions

int findGenMatch (double eta, double phi, double maxDeltaR, std::vector< MatchStruct > matches)
const reco::CandidatefindMother (const reco::Candidate *)
int findRecMatch (double eta, double phi, double maxdeltaR, std::vector< MatchStruct > matches)

Private Attributes

DQMStoredbe_
int eventNumber
MonitorElementhNumObjects
MonitorElementhNumOrphansGen
MonitorElementhNumOrphansRec
std::vector< MonitorElement * > hPassEtaGen
std::vector< MonitorElement * > hPassEtaRec
std::vector< MonitorElement * > hPassMaxPtGen
std::vector< MonitorElement * > hPassMaxPtRec
std::vector< MonitorElement * > hPassPhiGen
std::vector< MonitorElement * > hPassPhiRec
bool isIsolatedPath
bool makeNtuple
MonitorElementmeNumberOfEvents
unsigned int numHltLabels
std::string theAodL1Label
std::string theAodL2Label
std::vector< double > theEtaParameters
TFile * theFile
std::string theGenLabel
std::vector< std::string > theHltCollectionLabels
std::string theHltProcessName
std::string theL1CollectionLabel
double theL1DrCut
double theL2DrCut
double theL3DrCut
double theMaxEtaCut
std::vector< double > theMaxPtParameters
double theMinPtCut
int theMotherParticleId
std::vector< double > theNSigmas
TNtuple * theNtuple
std::string theNtupleFileName
float theNtuplePars [100]
std::string theNtuplePath
unsigned int theNumberOfObjects
std::vector< double > thePhiParameters
std::vector< double > thePtParameters
std::string theRecoLabel
std::string theTriggerName
bool useAod
bool useMuonFromGenerator
bool useMuonFromReco

Classes

struct  MatchStruct


Detailed Description

Get L1/HLT efficiency/rate plots Documentation available on the CMS TWiki: httpss://twiki.cern.ch/twiki/bin/view/CMS/MuonHLTOfflinePerformance.

Date
2009/02/11 20:21:26
Revision
1.30
Author:
M. Vander Donckt, J. Klukas (copied from J. Alcaraz)

Definition at line 45 of file HLTMuonGenericRate.h.


Constructor & Destructor Documentation

HLTMuonGenericRate::HLTMuonGenericRate ( const edm::ParameterSet pset,
std::string  triggerName,
std::vector< std::string >  moduleNames 
)

Constructor.


Member Function Documentation

void HLTMuonGenericRate::analyze ( const edm::Event iEvent  ) 

Definition at line 151 of file HLTMuonGenericRate.cc.

References funct::abs(), assign, coneSizes, reco::Particle::eta(), eta, eventNumber, edm::Handle< T >::failedToGet(), MonitorElement::Fill(), findGenMatch(), findMother(), findRecMatch(), HLTMuonGenericRate::MatchStruct::genCand, genParticles_cfi::genParticles, edm::Event::getByLabel(), MonitorElement::getTH1(), MonitorElement::getTH1F(), hNumObjects, hNumOrphansGen, hNumOrphansRec, hPassEtaGen, hPassEtaRec, hPassMaxPtGen, hPassMaxPtRec, hPassPhiGen, hPassPhiRec, isIsolatedPath, edm::Handle< T >::isValid(), j, k, bookConverter::keys, LogTrace, m, edm::match(), meNumberOfEvents, trajectoryFilterForConversions_cfi::minPt, minPtCuts, metsig::muon, n, numHltLabels, p4, reco::Particle::pdgId(), phi, reco::Particle::pt(), HLTMuonGenericRate::MatchStruct::recCand, size, StDecayID::status, reco::Particle::status(), theGenLabel, theHltCollectionLabels, theHltProcessName, theL1CollectionLabel, theNtuple, theNtuplePars, theTriggerName, trigger::TriggerL1Mu, trigger::TriggerMuon, useAod, useMuonFromReco, and dimuonsSequences_cff::veto.

00152 {
00153   
00154   eventNumber++;
00155   LogTrace( "HLTMuonVal" ) << "In analyze for trigger path " << 
00156     theTriggerName << ", Event:" << eventNumber;
00157 
00158   // Update event numbers
00159   meNumberOfEvents->Fill(eventNumber); 
00160 
00162   // Get all generated and reconstructed muons and create structs to hold  
00163   // matches to trigger candidates 
00164 
00165   double genMuonPt = -1;
00166   double recMuonPt = -1;
00167 
00168   std::vector<MatchStruct> genMatches;
00169   if ( useMuonFromGenerator ) {
00170     Handle<GenParticleCollection> genParticles;
00171     iEvent.getByLabel(theGenLabel, genParticles);
00172     for ( size_t i = 0; i < genParticles->size(); i++ ) {
00173       const reco::GenParticle *genParticle = &(*genParticles)[i];
00174       const Candidate *mother = findMother(genParticle);
00175       int    momId  = ( mother ) ? mother->pdgId() : 0;
00176       int    id     = genParticle->pdgId();
00177       int    status = genParticle->status();
00178       double pt     = genParticle->pt();
00179       double eta    = genParticle->eta();
00180       if ( abs(id) == 13  && status == 1 && 
00181            ( theMotherParticleId == 0 || abs(momId) == theMotherParticleId ) )
00182       {
00183         MatchStruct newMatchStruct;
00184         newMatchStruct.genCand = genParticle;
00185         genMatches.push_back(newMatchStruct);
00186         if ( pt > genMuonPt && fabs(eta) < theMaxEtaCut )
00187           genMuonPt = pt;
00188   } } }
00189 
00190   std::vector<MatchStruct> recMatches;
00191   if ( useMuonFromReco ) {
00192     Handle<reco::TrackCollection> muTracks;
00193     iEvent.getByLabel(theRecoLabel, muTracks);    
00194     reco::TrackCollection::const_iterator muon;
00195     if  ( muTracks.failedToGet() ) {
00196       LogTrace("HLTMuonVal") << "No reco tracks to compare to";
00197       useMuonFromReco = false;
00198     } else {
00199       for ( muon = muTracks->begin(); muon != muTracks->end(); ++muon ) {
00200         float pt  = muon->pt();
00201         float eta = muon->eta();
00202         MatchStruct newMatchStruct;
00203         newMatchStruct.recCand = &*muon;
00204         recMatches.push_back(newMatchStruct);
00205         if ( pt > recMuonPt && fabs(eta) < theMaxEtaCut ) 
00206           recMuonPt = pt;
00207   } } } 
00208   
00209   LogTrace("HLTMuonVal") << "genMuonPt: " << genMuonPt << ", "  
00210                          << "recMuonPt: " << recMuonPt;
00211 
00213   // Get the L1 and HLT trigger collections
00214 
00215   edm::Handle<trigger::TriggerEventWithRefs> rawTriggerEvent;
00216   edm::Handle<trigger::TriggerEvent>         aodTriggerEvent;
00217   vector<LorentzVector>                      l1Particles;
00218   vector< vector<LorentzVector> >            hltParticles(numHltLabels);
00219   vector< vector<RecoChargedCandidateRef> >  hltCands(numHltLabels);
00220   InputTag collectionTag;
00221   size_t   filterIndex;
00222 
00223 
00225 
00226   if ( !useAod ) {
00227 
00228     iEvent.getByLabel( "hltTriggerSummaryRAW", rawTriggerEvent );
00229     if ( !rawTriggerEvent.isValid() ) { 
00230       LogError("HLTMuonVal") << "No RAW trigger summary found! "
00231                              << "Change UseAod to ""True"" in configuration";
00232       return;
00233     }
00234 
00235     collectionTag = InputTag( theL1CollectionLabel, "", theHltProcessName );
00236     filterIndex   = rawTriggerEvent->filterIndex(collectionTag);
00237     vector<L1MuonParticleRef> l1Cands;
00238 
00239     if ( filterIndex < rawTriggerEvent->size() ) 
00240       rawTriggerEvent->getObjects( filterIndex, TriggerL1Mu, l1Cands );
00241     else LogTrace("HLTMuonVal") << "No L1 Collection with label " 
00242                                 << collectionTag;
00243     
00244     for ( size_t i = 0; i < l1Cands.size(); i++ ) 
00245       l1Particles.push_back( l1Cands[i]->p4() );
00246     
00247     for ( size_t i = 0; i < numHltLabels; i++ ) {
00248 
00249       collectionTag = InputTag( theHltCollectionLabels[i], 
00250                                 "", theHltProcessName );
00251       filterIndex   = rawTriggerEvent->filterIndex(collectionTag);
00252 
00253       if ( filterIndex < rawTriggerEvent->size() )
00254         rawTriggerEvent->getObjects( filterIndex, TriggerMuon, hltCands[i]);
00255       else LogTrace("HLTMuonVal") << "No HLT Collection with label " 
00256                                   << collectionTag;
00257 
00258       for ( size_t j = 0; j < hltCands[i].size(); j++ )
00259         hltParticles[i].push_back( hltCands[i][j]->p4() );
00260 
00261     } // End loop over theHltCollectionLabels
00262 
00263   } // Done getting RAW trigger summary
00264 
00265 
00267 
00268   if ( useAod ) {
00269 
00270     iEvent.getByLabel("hltTriggerSummaryAOD", aodTriggerEvent);
00271     if ( !aodTriggerEvent.isValid() ) { 
00272       LogError("HLTMuonVal") << "No AOD trigger summary found! Returning..."; 
00273       return; 
00274     }
00275 
00276     const TriggerObjectCollection objects = aodTriggerEvent->getObjects();
00277 
00278     collectionTag = InputTag( theAodL1Label, "", theHltProcessName );
00279     filterIndex   = aodTriggerEvent->filterIndex( collectionTag );
00280     if ( filterIndex < aodTriggerEvent->sizeFilters() ) {
00281       const Keys &keys = aodTriggerEvent->filterKeys( filterIndex );
00282       for ( size_t j = 0; j < keys.size(); j++ )
00283         l1Particles.push_back( objects[keys[j]].particle().p4() );
00284     } 
00285 
00286     collectionTag = InputTag( theAodL2Label, "", theHltProcessName );
00287     filterIndex   = aodTriggerEvent->filterIndex( collectionTag );
00288     if ( filterIndex < aodTriggerEvent->sizeFilters() ) {
00289       const Keys &keys = aodTriggerEvent->filterKeys( filterIndex );
00290       for ( size_t j = 0; j < keys.size(); j++ )
00291         hltParticles[0].push_back( objects[keys[j]].particle().p4() );
00292     } 
00293 
00294     collectionTag = InputTag( theHltCollectionLabels.back(), "", 
00295                               theHltProcessName );
00296     filterIndex   = aodTriggerEvent->filterIndex( collectionTag );
00297     if ( filterIndex < aodTriggerEvent->sizeFilters() ) {
00298       const Keys &keys = aodTriggerEvent->filterKeys( filterIndex );
00299       for ( size_t j = 0; j < keys.size(); j++ )
00300         hltParticles[1].push_back( objects[keys[j]].particle().p4() );
00301     } 
00302 
00303     // At this point, we should check whether the prescaled L1 and L2
00304     // triggers actually fired, and exit if not.
00305     if ( l1Particles.size() == 0 || hltParticles[0].size() == 0 ) 
00306       { LogTrace("HLTMuonVal") << "L1 or L2 didn't fire"; return; }
00307     
00308   } // Done getting AOD trigger summary
00309   
00310   hNumObjects->getTH1()->AddBinContent( 3, l1Particles.size() );
00311 
00312   for ( size_t i = 0; i < numHltLabels; i++ ) 
00313     hNumObjects->getTH1()->AddBinContent( i + 4, hltParticles[i].size() );
00314 
00316   // Initialize MatchStructs
00317 
00318   LorentzVector nullLorentzVector( 0., 0., 0., -999. );
00319 
00320   for ( size_t i = 0; i < genMatches.size(); i++ ) {
00321     genMatches[i].l1Cand = nullLorentzVector;
00322     genMatches[i].hltCands. assign( numHltLabels, nullLorentzVector );
00323     genMatches[i].hltTracks.assign( numHltLabels, false );
00324   }
00325 
00326   for ( size_t i = 0; i < recMatches.size(); i++ ) {
00327     recMatches[i].l1Cand = nullLorentzVector;
00328     recMatches[i].hltCands. assign( numHltLabels, nullLorentzVector );
00329     recMatches[i].hltTracks.assign( numHltLabels, false );
00330   }
00331 
00333   // Loop through L1 candidates, matching to gen/reco muons 
00334 
00335   unsigned int numL1Cands = 0;
00336 
00337   for ( size_t i = 0; i < l1Particles.size(); i++ ) {
00338 
00339     LorentzVector l1Cand = l1Particles[i];
00340     double eta           = l1Cand.eta();
00341     double phi           = l1Cand.phi();
00342     // L1 pt is taken from a lookup table
00343     // double ptLUT      = l1Cand->pt();  
00344 
00345     double maxDeltaR = theL1DrCut;
00346     numL1Cands++;
00347 
00348     if ( useMuonFromGenerator ){
00349       int match = findGenMatch( eta, phi, maxDeltaR, genMatches );
00350       if ( match != -1 && genMatches[match].l1Cand.E() < 0 ) 
00351         genMatches[match].l1Cand = l1Cand;
00352       else hNumOrphansGen->getTH1F()->AddBinContent( 1 );
00353     }
00354 
00355     if ( useMuonFromReco ){
00356       int match = findRecMatch( eta, phi, maxDeltaR, recMatches );
00357       if ( match != -1 && recMatches[match].l1Cand.E() < 0 ) 
00358         recMatches[match].l1Cand = l1Cand;
00359       else hNumOrphansRec->getTH1F()->AddBinContent( 1 );
00360     }
00361 
00362   } // End loop over l1Particles
00363 
00364   LogTrace("HLTMuonVal") << "Number of L1 Cands: " << numL1Cands;
00365 
00367   // Loop through HLT candidates, matching to gen/reco muons
00368 
00369   vector<unsigned int> numHltCands( numHltLabels, 0) ;
00370 
00371   for ( size_t i = 0; i < numHltLabels; i++ ) { 
00372 
00373     int triggerLevel      = ( i < ( numHltLabels / 2 ) ) ? 2 : 3;
00374     double maxDeltaR      = ( triggerLevel == 2 ) ? theL2DrCut : theL3DrCut;
00375 
00376     for ( size_t candNum = 0; candNum < hltParticles[i].size(); candNum++ ) {
00377 
00378       LorentzVector hltCand = hltParticles[i][candNum];
00379       double eta            = hltCand.eta();
00380       double phi            = hltCand.phi();
00381 
00382       numHltCands[i]++;
00383 
00384       if ( useMuonFromGenerator ){
00385         int match = findGenMatch( eta, phi, maxDeltaR, genMatches );
00386       
00387         if ( match != -1 && genMatches[match].hltCands[i].E() < 0 ) {
00388           genMatches[match].hltCands[i] = hltCand;
00389           if ( !useAod ) genMatches[match].hltTracks[i] = 
00390              &*hltCands[i][candNum];
00391         }
00392         else hNumOrphansGen->getTH1F()->AddBinContent( i + 2 );
00393       }
00394 
00395       if ( useMuonFromReco ){
00396         int match  = findRecMatch( eta, phi, maxDeltaR, recMatches );
00397         if ( match != -1 && recMatches[match].hltCands[i].E() < 0 ) 
00398           recMatches[match].hltCands[i] = hltCand;
00399         else hNumOrphansRec->getTH1F()->AddBinContent( i + 2 );
00400       }
00401 
00402       LogTrace("HLTMuonVal") << "Number of L1 Cands: " << numHltCands[i];
00403 
00404     } // End loop over HLT particles
00405 
00406   } // End loop over HLT labels
00407 
00409   // Fill ntuple
00410 
00411   if ( makeNtuple ) {
00412     Handle<reco::IsoDepositMap> caloDepMap, trackDepMap;
00413     iEvent.getByLabel("hltL2MuonIsolations",caloDepMap);
00414     iEvent.getByLabel("hltL3MuonIsolations",trackDepMap);
00415     IsoDeposit::Vetos vetos;
00416     if ( isIsolatedPath )
00417       for ( size_t i = 0; i < hltCands[2].size(); i++ ) {
00418         TrackRef tk = hltCands[2][i]->get<TrackRef>();
00419         vetos.push_back( (*trackDepMap)[tk].veto() );
00420       }
00421     for ( size_t i = 0; i < genMatches.size(); i++ ) {
00422       for ( int k = 0; k < 50; k++ ) theNtuplePars[k] = -99;
00423       theNtuplePars[0] = eventNumber;
00424       theNtuplePars[1] = (findMother(genMatches[i].genCand))->pdgId();
00425       theNtuplePars[4] = genMatches[i].genCand->pt();
00426       theNtuplePars[5] = genMatches[i].genCand->eta();
00427       theNtuplePars[6] = genMatches[i].genCand->phi();
00428       if ( genMatches[i].l1Cand.E() > 0 ) {
00429         theNtuplePars[7] = genMatches[i].l1Cand.pt();
00430         theNtuplePars[8] = genMatches[i].l1Cand.eta();
00431         theNtuplePars[9] = genMatches[i].l1Cand.phi();
00432       }
00433       for ( size_t j = 0; j < genMatches[i].hltCands.size(); j++ ) {
00434         if ( genMatches[i].hltCands[j].E() > 0 ) {
00435           if ( j == 0 ) {
00436             theNtuplePars[10] = genMatches[i].hltCands[j].pt();
00437             theNtuplePars[11] = genMatches[i].hltCands[j].eta();
00438             theNtuplePars[12] = genMatches[i].hltCands[j].phi();
00439             if ( isIsolatedPath && !useAod ) {
00440               TrackRef tk = genMatches[i].hltTracks[j]->get<TrackRef>();
00441               const IsoDeposit &dep = (*caloDepMap)[tk];
00442               for ( int m = 0; m < numCones; m++ ) {
00443                 double dr = coneSizes[m];
00444                 std::pair<double,int> depInfo = dep.depositAndCountWithin(dr);
00445                 theNtuplePars[ 16 + 4*m + 0 ] = depInfo.first;
00446                 theNtuplePars[ 16 + 4*m + 1 ] = depInfo.second;
00447           } } }
00448           if ( ( !isIsolatedPath && j == 1 ) ||
00449                (  isIsolatedPath && j == 2 ) ) {
00450             theNtuplePars[13] = genMatches[i].hltCands[j].pt();
00451             theNtuplePars[14] = genMatches[i].hltCands[j].eta();
00452             theNtuplePars[15] = genMatches[i].hltCands[j].phi();
00453             if ( isIsolatedPath ) {
00454               TrackRef tk = genMatches[i].hltTracks[j]->get<TrackRef>();
00455               const IsoDeposit &dep = (*trackDepMap)[tk];
00456               for ( int m = 0; m < numCones; m++ ) {
00457                 for ( int n = 0; n < numMinPtCuts; n++ ) {
00458                   double dr = coneSizes[m];
00459                   double minPt = minPtCuts[n];
00460                   std::pair<double,int> depInfo;
00461                   depInfo = dep.depositAndCountWithin(dr, vetos, minPt);
00462                   int currentPlace = 16 + 4*numCones + 2*numMinPtCuts*m + 2*n;
00463                   theNtuplePars[ currentPlace + 0 ] = depInfo.first;
00464                   theNtuplePars[ currentPlace + 1 ] = depInfo.second;
00465                 }
00466           } } }
00467           if ( isIsolatedPath && j == 1 ) theNtuplePars[2] = true;
00468           if ( isIsolatedPath && j == 3 ) theNtuplePars[3] = true;
00469         }
00470       }
00471       theNtuple->Fill(theNtuplePars); 
00472     } // Done filling ntuple
00473   }
00474   
00476   // Fill histograms
00477 
00478   if ( genMuonPt > 0 ) hPassMaxPtGen[0]->Fill( genMuonPt );
00479   if ( recMuonPt > 0 ) hPassMaxPtRec[0]->Fill( recMuonPt );
00480   if ( numL1Cands >= theNumberOfObjects ) {
00481     if ( genMuonPt > 0 ) hPassMaxPtGen[1]->Fill( genMuonPt );
00482     if ( recMuonPt > 0 ) hPassMaxPtRec[1]->Fill( recMuonPt );
00483   }
00484   for ( size_t i = 0; i < numHltLabels; i++ ) {
00485     if ( numHltCands[i] >= theNumberOfObjects ) {
00486       if ( genMuonPt > 0 ) hPassMaxPtGen[i+2]->Fill( genMuonPt );
00487       if ( recMuonPt > 0 ) hPassMaxPtRec[i+2]->Fill( recMuonPt );
00488     }
00489   }
00490 
00491   for ( size_t i = 0; i < genMatches.size(); i++ ) {
00492     double pt  = genMatches[i].genCand->pt();
00493     double eta = genMatches[i].genCand->eta();
00494     double phi = genMatches[i].genCand->phi();
00495     if ( pt > theMinPtCut &&  fabs(eta) < theMaxEtaCut ) {
00496       hNumObjects->getTH1()->AddBinContent(1);
00497       hPassEtaGen[0]->Fill(eta);
00498       hPassPhiGen[0]->Fill(phi);
00499       if ( genMatches[i].l1Cand.E() > 0 ) {
00500         hPassEtaGen[1]->Fill(eta);
00501         hPassPhiGen[1]->Fill(phi);
00502         bool foundAllPreviousCands = true;
00503         for ( size_t j = 0; j < genMatches[i].hltCands.size(); j++ ) {
00504           if ( foundAllPreviousCands && genMatches[i].hltCands[j].E() > 0 ) {
00505             hPassEtaGen[j+2]->Fill(eta);
00506             hPassPhiGen[j+2]->Fill(phi);
00507           } else foundAllPreviousCands = false;
00508   } } } }
00509 
00510   for ( size_t i = 0; i < recMatches.size(); i++ ) {
00511     double pt  = recMatches[i].recCand->pt();
00512     double eta = recMatches[i].recCand->eta();
00513     double phi = recMatches[i].recCand->phi();
00514     if ( pt > theMinPtCut &&  fabs(eta) < theMaxEtaCut ) {
00515       hNumObjects->getTH1()->AddBinContent(2);
00516       hPassEtaRec[0]->Fill(eta);
00517       hPassPhiRec[0]->Fill(phi);
00518       if ( recMatches[i].l1Cand.E() > 0 ) {
00519         hPassEtaRec[1]->Fill(eta);
00520         hPassPhiRec[1]->Fill(phi);
00521         bool foundAllPreviousCands = true;
00522         for ( size_t j = 0; j < recMatches[i].hltCands.size(); j++ ) {
00523           if ( foundAllPreviousCands && recMatches[i].hltCands[j].E() > 0 ) {
00524             hPassEtaRec[j+2]->Fill(eta);
00525             hPassPhiRec[j+2]->Fill(phi);
00526           } else foundAllPreviousCands = false;
00527   } } } }
00528 
00529 } // Done filling histograms

void HLTMuonGenericRate::begin (  ) 

Definition at line 584 of file HLTMuonGenericRate.cc.

References DQMStore::book1D(), DQMStore::bookFloat(), DQMStore::bookInt(), bookIt(), DQMStore::cd(), dbe_, MonitorElement::Fill(), MonitorElement::getTH1(), h, hNumObjects, hNumOrphansGen, hNumOrphansRec, hPassEtaGen, hPassEtaRec, hPassMaxPtGen, hPassMaxPtRec, hPassPhiGen, hPassPhiRec, level, meNumberOfEvents, MonitorElement::setBinLabel(), DQMStore::setCurrentFolder(), theHltCollectionLabels, theL1CollectionLabel, theMaxPtParameters, theTriggerName, and useMuonFromReco.

00585 {
00586 
00587   TString myLabel, newFolder;
00588   vector<TH1F*> h;
00589 
00590   if ( dbe_ ) {
00591     dbe_->cd();
00592     dbe_->setCurrentFolder("HLT/Muon");
00593 
00594     myLabel = theL1CollectionLabel;
00595     myLabel = myLabel(myLabel.Index("L1"),myLabel.Length());
00596     myLabel = myLabel(0,myLabel.Index("Filtered")+8);
00597 
00598     newFolder = "HLT/Muon/Distributions/" + theTriggerName;
00599     dbe_->setCurrentFolder( newFolder.Data() );
00600 
00601     meNumberOfEvents            = dbe_->bookInt("NumberOfEvents");
00602     MonitorElement *meMinPtCut  = dbe_->bookFloat("MinPtCut"    );
00603     MonitorElement *meMaxEtaCut = dbe_->bookFloat("MaxEtaCut"   );
00604     meMinPtCut ->Fill(theMinPtCut );
00605     meMaxEtaCut->Fill(theMaxEtaCut);
00606     
00607     vector<string> binLabels;
00608     binLabels.push_back( theL1CollectionLabel.c_str() );
00609     for ( size_t i = 0; i < theHltCollectionLabels.size(); i++ )
00610       binLabels.push_back( theHltCollectionLabels[i].c_str() );
00611 
00612     hNumObjects = dbe_->book1D( "numObjects", "Number of Objects", 7, 0, 7 );
00613     hNumObjects->setBinLabel( 1, "Gen" );
00614     hNumObjects->setBinLabel( 2, "Reco" );
00615     for ( size_t i = 0; i < binLabels.size(); i++ )
00616       hNumObjects->setBinLabel( i + 3, binLabels[i].c_str() );
00617     hNumObjects->getTH1()->LabelsDeflate("X");
00618 
00619     if ( useMuonFromGenerator ){
00620 
00621       hNumOrphansGen = dbe_->book1D( "genNumOrphans", "Number of Orphans;;Number of Objects Not Matched to a Gen #mu", 5, 0, 5 );
00622       for ( size_t i = 0; i < binLabels.size(); i++ )
00623         hNumOrphansGen->setBinLabel( i + 1, binLabels[i].c_str() );
00624       hNumOrphansGen->getTH1()->LabelsDeflate("X");
00625 
00626       hPassMaxPtGen.push_back( bookIt( "genPassMaxPt_All", "pt of Leading Gen Muon", theMaxPtParameters) );
00627       hPassMaxPtGen.push_back( bookIt( "genPassMaxPt_" + myLabel, "pt of Leading Gen Muon, if matched to " + myLabel, theMaxPtParameters) );
00628       hPassEtaGen.push_back( bookIt( "genPassEta_All", "#eta of Gen Muons", theEtaParameters) );
00629       hPassEtaGen.push_back( bookIt( "genPassEta_" + myLabel, "#eta of Gen Muons matched to " + myLabel, theEtaParameters) );
00630       hPassPhiGen.push_back( bookIt( "genPassPhi_All", "#phi of Gen Muons", thePhiParameters) );
00631       hPassPhiGen.push_back( bookIt( "genPassPhi_" + myLabel, "#phi of Gen Muons matched to " + myLabel, thePhiParameters) );
00632 
00633     }
00634 
00635     if ( useMuonFromReco ){
00636 
00637       hNumOrphansRec = dbe_->book1D( "recNumOrphans", "Number of Orphans;;Number of Objects Not Matched to a Reconstructed #mu", 5, 0, 5 );
00638       for ( size_t i = 0; i < binLabels.size(); i++ )
00639         hNumOrphansRec->setBinLabel( i + 1, binLabels[i].c_str() );
00640       hNumOrphansRec->getTH1()->LabelsDeflate("X");
00641 
00642       hPassMaxPtRec.push_back( bookIt( "recPassMaxPt_All", "pt of Leading Reco Muon" + myLabel,  theMaxPtParameters) );
00643       hPassMaxPtRec.push_back( bookIt( "recPassMaxPt_" + myLabel, "pt of Leading Reco Muon, if matched to " + myLabel,  theMaxPtParameters) );
00644       hPassEtaRec.push_back( bookIt( "recPassEta_All", "#eta of Reco Muons", theEtaParameters) );
00645       hPassEtaRec.push_back( bookIt( "recPassEta_" + myLabel, "#eta of Reco Muons matched to " + myLabel, theEtaParameters) );
00646       hPassPhiRec.push_back( bookIt( "recPassPhi_All", "#phi of Reco Muons", thePhiParameters) );
00647       hPassPhiRec.push_back( bookIt( "recPassPhi_" + myLabel, "#phi of Reco Muons matched to " + myLabel, thePhiParameters) );
00648 
00649     }
00650 
00651     for ( unsigned int i = 0; i < theHltCollectionLabels.size(); i++ ) {
00652 
00653       myLabel = theHltCollectionLabels[i];
00654       TString level = ( myLabel.Contains("L2") ) ? "L2" : "L3";
00655       myLabel = myLabel(myLabel.Index(level),myLabel.Length());
00656       myLabel = myLabel(0,myLabel.Index("Filtered")+8);
00657       
00658       if ( useMuonFromGenerator ) {
00659 
00660         hPassMaxPtGen.push_back( bookIt( "genPassMaxPt_" + myLabel, "pt of Leading Gen Muon, if matched to " + myLabel, theMaxPtParameters) );   
00661         hPassEtaGen.push_back( bookIt( "genPassEta_" + myLabel, "#eta of Gen Muons matched to " + myLabel, theEtaParameters) );
00662         hPassPhiGen.push_back( bookIt( "genPassPhi_" + myLabel, "#phi of Gen Muons matched to " + myLabel, thePhiParameters) );
00663 
00664       }
00665 
00666       if ( useMuonFromReco ) {
00667 
00668         hPassMaxPtRec.push_back( bookIt( "recPassMaxPt_" + myLabel, "pt of Leading Reco Muon, if matched to " + myLabel, theMaxPtParameters) );     
00669         hPassEtaRec.push_back( bookIt( "recPassEta_" + myLabel, "#eta of Reco Muons matched to " + myLabel, theEtaParameters) );
00670         hPassPhiRec.push_back( bookIt( "recPassPhi_" + myLabel, "#phi of Reco Muons matched to " + myLabel, thePhiParameters) );
00671       }
00672 
00673     }
00674   }
00675 
00676 }

MonitorElement* HLTMuonGenericRate::bookIt ( TString  name,
TString  title,
std::vector< double >   
)

Referenced by begin().

int HLTMuonGenericRate::findGenMatch ( double  eta,
double  phi,
double  maxDeltaR,
std::vector< MatchStruct matches 
) [private]

Referenced by analyze().

const reco::Candidate * HLTMuonGenericRate::findMother ( const reco::Candidate p  )  [private]

Definition at line 534 of file HLTMuonGenericRate.cc.

References reco::Candidate::mother(), and reco::Particle::pdgId().

Referenced by analyze().

00535 {
00536   const reco::Candidate* mother = p->mother();
00537   if ( mother ) {
00538     if ( mother->pdgId() == p->pdgId() ) return findMother(mother);
00539     else return mother;
00540   }
00541   else return 0;
00542 }

int HLTMuonGenericRate::findRecMatch ( double  eta,
double  phi,
double  maxdeltaR,
std::vector< MatchStruct matches 
) [private]

Referenced by analyze().

void HLTMuonGenericRate::finish (  ) 

Definition at line 140 of file HLTMuonGenericRate.cc.

References theFile, and theNtuple.

00141 {
00142   if ( makeNtuple ) {
00143     theFile->cd();
00144     theNtuple->Write();
00145     theFile->Close();
00146   }
00147 }


Member Data Documentation

DQMStore* HLTMuonGenericRate::dbe_ [private]

Definition at line 125 of file HLTMuonGenericRate.h.

Referenced by begin().

int HLTMuonGenericRate::eventNumber [private]

Definition at line 119 of file HLTMuonGenericRate.h.

Referenced by analyze().

MonitorElement* HLTMuonGenericRate::hNumObjects [private]

Definition at line 134 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

MonitorElement* HLTMuonGenericRate::hNumOrphansGen [private]

Definition at line 135 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

MonitorElement* HLTMuonGenericRate::hNumOrphansRec [private]

Definition at line 136 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::vector<MonitorElement*> HLTMuonGenericRate::hPassEtaGen [private]

Definition at line 128 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::vector<MonitorElement*> HLTMuonGenericRate::hPassEtaRec [private]

Definition at line 131 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::vector<MonitorElement*> HLTMuonGenericRate::hPassMaxPtGen [private]

Definition at line 127 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::vector<MonitorElement*> HLTMuonGenericRate::hPassMaxPtRec [private]

Definition at line 130 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::vector<MonitorElement*> HLTMuonGenericRate::hPassPhiGen [private]

Definition at line 129 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::vector<MonitorElement*> HLTMuonGenericRate::hPassPhiRec [private]

Definition at line 132 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

bool HLTMuonGenericRate::isIsolatedPath [private]

Definition at line 121 of file HLTMuonGenericRate.h.

Referenced by analyze().

bool HLTMuonGenericRate::makeNtuple [private]

Definition at line 79 of file HLTMuonGenericRate.h.

MonitorElement* HLTMuonGenericRate::meNumberOfEvents [private]

Definition at line 137 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

unsigned int HLTMuonGenericRate::numHltLabels [private]

Definition at line 120 of file HLTMuonGenericRate.h.

Referenced by analyze().

std::string HLTMuonGenericRate::theAodL1Label [private]

Definition at line 98 of file HLTMuonGenericRate.h.

std::string HLTMuonGenericRate::theAodL2Label [private]

Definition at line 99 of file HLTMuonGenericRate.h.

std::vector<double> HLTMuonGenericRate::theEtaParameters [private]

Definition at line 103 of file HLTMuonGenericRate.h.

TFile* HLTMuonGenericRate::theFile [private]

Definition at line 82 of file HLTMuonGenericRate.h.

Referenced by finish().

std::string HLTMuonGenericRate::theGenLabel [private]

Definition at line 94 of file HLTMuonGenericRate.h.

Referenced by analyze().

std::vector<std::string> HLTMuonGenericRate::theHltCollectionLabels [private]

Definition at line 89 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

std::string HLTMuonGenericRate::theHltProcessName [private]

Definition at line 86 of file HLTMuonGenericRate.h.

Referenced by analyze().

std::string HLTMuonGenericRate::theL1CollectionLabel [private]

Definition at line 88 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

double HLTMuonGenericRate::theL1DrCut [private]

Definition at line 108 of file HLTMuonGenericRate.h.

double HLTMuonGenericRate::theL2DrCut [private]

Definition at line 109 of file HLTMuonGenericRate.h.

double HLTMuonGenericRate::theL3DrCut [private]

Definition at line 110 of file HLTMuonGenericRate.h.

double HLTMuonGenericRate::theMaxEtaCut [private]

Definition at line 107 of file HLTMuonGenericRate.h.

std::vector<double> HLTMuonGenericRate::theMaxPtParameters [private]

Definition at line 101 of file HLTMuonGenericRate.h.

Referenced by begin().

double HLTMuonGenericRate::theMinPtCut [private]

Definition at line 106 of file HLTMuonGenericRate.h.

int HLTMuonGenericRate::theMotherParticleId [private]

Definition at line 111 of file HLTMuonGenericRate.h.

std::vector<double> HLTMuonGenericRate::theNSigmas [private]

Definition at line 112 of file HLTMuonGenericRate.h.

TNtuple* HLTMuonGenericRate::theNtuple [private]

Definition at line 81 of file HLTMuonGenericRate.h.

Referenced by analyze(), and finish().

std::string HLTMuonGenericRate::theNtupleFileName [private]

Definition at line 114 of file HLTMuonGenericRate.h.

float HLTMuonGenericRate::theNtuplePars[100] [private]

Definition at line 80 of file HLTMuonGenericRate.h.

Referenced by analyze().

std::string HLTMuonGenericRate::theNtuplePath [private]

Definition at line 115 of file HLTMuonGenericRate.h.

unsigned int HLTMuonGenericRate::theNumberOfObjects [private]

Definition at line 90 of file HLTMuonGenericRate.h.

std::vector<double> HLTMuonGenericRate::thePhiParameters [private]

Definition at line 104 of file HLTMuonGenericRate.h.

std::vector<double> HLTMuonGenericRate::thePtParameters [private]

Definition at line 102 of file HLTMuonGenericRate.h.

std::string HLTMuonGenericRate::theRecoLabel [private]

Definition at line 95 of file HLTMuonGenericRate.h.

std::string HLTMuonGenericRate::theTriggerName [private]

Definition at line 87 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().

bool HLTMuonGenericRate::useAod [private]

Definition at line 97 of file HLTMuonGenericRate.h.

Referenced by analyze().

bool HLTMuonGenericRate::useMuonFromGenerator [private]

Definition at line 92 of file HLTMuonGenericRate.h.

bool HLTMuonGenericRate::useMuonFromReco [private]

Definition at line 93 of file HLTMuonGenericRate.h.

Referenced by analyze(), and begin().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:24:36 2009 for CMSSW by  doxygen 1.5.4