CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/PhysicsTools/PatExamples/plugins/PatTriggerAnalyzer.cc

Go to the documentation of this file.
00001 #include "TMath.h"
00002 #include "FWCore/ServiceRegistry/interface/Service.h"
00003 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00004 #include "PhysicsTools/PatUtils/interface/TriggerHelper.h"
00005 #include "PhysicsTools/PatExamples/plugins/PatTriggerAnalyzer.h"
00006 
00007 
00008 using namespace pat;
00009 
00010 
00011 PatTriggerAnalyzer::PatTriggerAnalyzer( const edm::ParameterSet & iConfig ) :
00012   // pat::Trigger
00013   trigger_( iConfig.getParameter< edm::InputTag >( "trigger" ) ),
00014   // pat::TriggerEvent
00015   triggerEvent_( iConfig.getParameter< edm::InputTag >( "triggerEvent" ) ),
00016   // muon input collection
00017   muons_( iConfig.getParameter< edm::InputTag >( "muons" ) ),
00018   // muon match objects
00019   muonMatch_( iConfig.getParameter< std::string >( "muonMatch" ) ),
00020   // minimal id for of all trigger objects
00021   minID_( iConfig.getParameter< unsigned >( "minID" ) ),
00022   // maximal id for of all trigger objects
00023   maxID_( iConfig.getParameter< unsigned >( "maxID" ) ),
00024   histos1D_(), histos2D_()
00025 {
00026 }
00027 
00028 PatTriggerAnalyzer::~PatTriggerAnalyzer()
00029 {
00030 }
00031 
00032 void PatTriggerAnalyzer::beginJob()
00033 {
00034   edm::Service< TFileService > fileService;
00035 
00036   // pt correlation plot
00037   histos2D_[ "ptTrigCand" ] = fileService->make< TH2D >( "ptTrigCand", "Object vs. candidate p_{T} (GeV)", 60, 0., 300., 60, 0., 300. );
00038   histos2D_[ "ptTrigCand" ]->SetXTitle( "candidate p_{T} (GeV)" );
00039   histos2D_[ "ptTrigCand" ]->SetYTitle( "object p_{T} (GeV)" );
00040   // eta correlation plot
00041   histos2D_[ "etaTrigCand" ] = fileService->make< TH2D >( "etaTrigCand", "Object vs. candidate #eta", 50, -2.5, 2.5, 50, -2.5, 2.5 );
00042   histos2D_[ "etaTrigCand" ]->SetXTitle( "candidate #eta" );
00043   histos2D_[ "etaTrigCand" ]->SetYTitle( "object #eta" );
00044   // phi correlation plot
00045   histos2D_[ "phiTrigCand" ] = fileService->make< TH2D >( "phiTrigCand", "Object vs. candidate #phi", 60, -TMath::Pi(), TMath::Pi(), 60, -TMath::Pi(), TMath::Pi() );
00046   histos2D_[ "phiTrigCand" ]->SetXTitle( "candidate #phi" );
00047   histos2D_[ "phiTrigCand" ]->SetYTitle( "object #phi" );
00048   // turn-on curves
00049   histos1D_[ "turnOn" ] = fileService->make< TH1D >( "turnOn", "p_{T} (GeV) of matched candidate", 10, 0., 50.);
00050   histos1D_[ "turnOn" ]->SetXTitle( "candidate p_{T} (GeV)" );
00051   histos1D_[ "turnOn" ]->SetYTitle( "# of objects" );
00052   // mean pt for all trigger objects
00053   histos1D_[ "ptMean" ] = fileService->make< TH1D >( "ptMean", "Mean p_{T} (GeV) per trigger object type", maxID_ - minID_ + 1, minID_ - 0.5, maxID_ + 0.5);
00054   histos1D_[ "ptMean" ]->SetXTitle( "trigger object type" );
00055   histos1D_[ "ptMean" ]->SetYTitle( "mean p_{T} (GeV)" );
00056 
00057   // initialize counters for mean pt calculation
00058   for( unsigned id = minID_; id <= maxID_; ++id ){
00059     sumN_ [ id ] = 0;
00060     sumPt_[ id ] = 0.;
00061   }
00062 }
00063 
00064 void PatTriggerAnalyzer::analyze( const edm::Event & iEvent, const edm::EventSetup & iSetup )
00065 {
00066   // PAT trigger event
00067   edm::Handle< TriggerEvent > triggerEvent;
00068   iEvent.getByLabel( triggerEvent_, triggerEvent );
00069 
00070   // PAT object collection
00071   edm::Handle< MuonCollection > muons;
00072   iEvent.getByLabel( muons_, muons );
00073 
00074   // PAT trigger helper for trigger matching information
00075   const helper::TriggerMatchHelper matchHelper;
00076 
00077   /*
00078     kinematics comparison
00079   */
00080 
00081   // loop over muon references (PAT muons have been used in the matcher in task 3)
00082   for( size_t iMuon = 0; iMuon < muons->size(); ++iMuon ) {
00083     // we need all these ingedients to recieve matched trigger objects from the matchHelper
00084     const TriggerObjectRef trigRef( matchHelper.triggerMatchObject( muons, iMuon, muonMatch_, iEvent, *triggerEvent ) );
00085     // finally we can fill the histograms
00086     if ( trigRef.isAvailable() && trigRef.isNonnull() ) { // check references (necessary!)
00087       histos2D_[ "ptTrigCand" ]->Fill( muons->at( iMuon ).pt(), trigRef->pt() );
00088       histos2D_[ "etaTrigCand" ]->Fill( muons->at( iMuon ).eta(), trigRef->eta() );
00089       histos2D_[ "phiTrigCand" ]->Fill( muons->at( iMuon ).phi(), trigRef->phi() );
00090     }
00091   }
00092 
00093   /*
00094     turn-on curve
00095   */
00096 
00097   // get the trigger objects corresponding to the used matching (HLT muons)
00098   const TriggerObjectRefVector trigRefs( triggerEvent->objects( trigger::TriggerMuon ) );
00099   // loop over selected trigger objects
00100   for ( TriggerObjectRefVector::const_iterator iTrig = trigRefs.begin(); iTrig != trigRefs.end(); ++iTrig ) {
00101     // get all matched candidates for the trigger object
00102     const reco::CandidateBaseRefVector candRefs( matchHelper.triggerMatchCandidates( ( *iTrig ), muonMatch_, iEvent, *triggerEvent ) );
00103     if ( candRefs.empty() ) continue;
00104     // fill the histogram...
00105     // (only for the first match, since we resolved ambiguities in the matching configuration,
00106     // so that we have one at maximum per trigger object)
00107     reco::CandidateBaseRef muonRef( candRefs.at( 0 ) );
00108     if ( muonRef.isAvailable() && muonRef.isNonnull() ) {
00109       histos1D_[ "turnOn" ]->Fill( muonRef->pt() );
00110     }
00111   }
00112 
00113   /*
00114     mean pt
00115   */
00116 
00117   // loop over all trigger match objects from minID to maxID; have
00118   // a look to DataFormats/HLTReco/interface/TriggerTypeDefs.h to
00119   // know more about the available trrigger object id's
00120   for ( unsigned id=minID_; id<=maxID_; ++id ) {
00121     // vector of all objects for a given object id
00122     const TriggerObjectRefVector objRefs( triggerEvent->objects( id ) );
00123     // buffer the number of objects
00124     sumN_[ id ] += objRefs.size();
00125     // iterate the objects and buffer the pt of the objects
00126     for ( TriggerObjectRefVector::const_iterator iRef = objRefs.begin(); iRef != objRefs.end(); ++iRef ) {
00127       sumPt_[ id ] += ( *iRef )->pt();
00128     }
00129   }
00130 }
00131 
00132 void PatTriggerAnalyzer::endJob()
00133 {
00134   // normalize the entries for the mean pt plot
00135   for(unsigned id=minID_; id<=maxID_; ++id){
00136     if( sumN_[ id ]!=0 ) histos1D_[ "ptMean" ]->Fill( id, sumPt_[ id ]/sumN_[ id ] );
00137   }
00138 }
00139 
00140 
00141 #include "FWCore/Framework/interface/MakerMacros.h"
00142 DEFINE_FWK_MODULE( PatTriggerAnalyzer );