00001 00013 #include "HiggsAnalysis/Skimming/interface/HiggsToWW2LeptonsSkim.h" 00014 00015 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h" 00016 #include "DataFormats/TrackReco/interface/TrackFwd.h" 00017 #include "DataFormats/TrackReco/interface/Track.h" 00018 00019 // Muons: 00020 #include <DataFormats/TrackReco/interface/Track.h> 00021 00022 // Electrons 00023 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h" 00024 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h" 00025 00026 #include "DataFormats/Candidate/interface/Candidate.h" 00027 00028 00029 #include <iostream> 00030 00031 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00032 00033 using namespace edm; 00034 using namespace std; 00035 using namespace reco; 00036 00037 HiggsToWW2LeptonsSkim::HiggsToWW2LeptonsSkim(const edm::ParameterSet& iConfig) : 00038 nEvents_(0), nAccepted_(0) 00039 { 00040 00041 // Reconstructed objects 00042 recTrackLabel = iConfig.getParameter<edm::InputTag>("RecoTrackLabel"); 00043 theGLBMuonLabel = iConfig.getParameter<edm::InputTag>("GlobalMuonCollectionLabel"); 00044 theGsfELabel = iConfig.getParameter<edm::InputTag>("ElectronCollectionLabel"); 00045 00046 singleTrackPtMin_ = iConfig.getParameter<double>("SingleTrackPtMin"); 00047 diTrackPtMin_ = iConfig.getParameter<double>("DiTrackPtMin"); 00048 etaMin_ = iConfig.getParameter<double>("etaMin"); 00049 etaMax_ = iConfig.getParameter<double>("etaMax"); 00050 } 00051 00052 00053 HiggsToWW2LeptonsSkim::~HiggsToWW2LeptonsSkim() 00054 { 00055 } 00056 00057 void HiggsToWW2LeptonsSkim::endJob() 00058 { 00059 edm::LogVerbatim("HiggsToWW2LeptonsSkim") 00060 << "Events read " << nEvents_ 00061 << " Events accepted " << nAccepted_ 00062 << "\nEfficiency " << ((double)nAccepted_)/((double)nEvents_) 00063 << std::endl; 00064 } 00065 00066 00067 bool HiggsToWW2LeptonsSkim::filter(edm::Event& event, const edm::EventSetup& iSetup) 00068 { 00069 00070 nEvents_++; 00071 bool accepted = false; 00072 bool accepted1 = false; 00073 int nTrackOver2ndCut = 0; 00074 00075 00076 // Handle<CandidateCollection> tracks; 00077 00078 using reco::TrackCollection; 00079 00080 // Get the muon track collection from the event 00081 edm::Handle<reco::TrackCollection> muTracks; 00082 event.getByLabel(theGLBMuonLabel.label(), muTracks); 00083 00084 if ( muTracks.isValid() ) { 00085 00086 reco::TrackCollection::const_iterator muons; 00087 00088 // Loop over muon collections and count how many muons there are, 00089 // and how many are above threshold 00090 for ( muons = muTracks->begin(); muons != muTracks->end(); ++muons ) { 00091 if ( muons->eta() > etaMin_ && muons->eta() < etaMax_ ) { 00092 if ( muons->pt() > singleTrackPtMin_ ) accepted1 = true; 00093 if ( muons->pt() > diTrackPtMin_ ) nTrackOver2ndCut++; 00094 } 00095 } 00096 } 00097 00098 // Now look at electrons: 00099 00100 // Get the electron track collection from the event 00101 edm::Handle<reco::GsfElectronCollection> pTracks; 00102 00103 event.getByLabel(theGsfELabel.label(),pTracks); 00104 00105 if ( pTracks.isValid() ) { 00106 00107 const reco::GsfElectronCollection* eTracks = pTracks.product(); 00108 00109 reco::GsfElectronCollection::const_iterator electrons; 00110 00111 // Loop over electron collections and count how many muons there are, 00112 // and how many are above threshold 00113 for ( electrons = eTracks->begin(); electrons != eTracks->end(); ++electrons ) { 00114 if ( electrons->eta() > etaMin_ && electrons->eta() < etaMax_ ) { 00115 if ( electrons->pt() > singleTrackPtMin_ ) accepted1 = true; 00116 if ( electrons->pt() > diTrackPtMin_ ) nTrackOver2ndCut++; 00117 } 00118 } 00119 } 00120 00121 00122 if ( accepted1 && nTrackOver2ndCut >= 2 ) accepted = true; 00123 00124 if ( accepted ) nAccepted_++; 00125 00126 return accepted; 00127 00128 }